fvSurfaceMapper.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  FV surface mapper.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fvSurfaceMapper.H"
30 #include "fvMesh.H"
31 #include "mapPolyMesh.H"
32 #include "faceMapper.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 void Foam::fvSurfaceMapper::calcAddressing() const
37 {
38  if
39  (
40  directAddrPtr_
41  || interpolationAddrPtr_
42  || weightsPtr_
43  || insertedObjectLabelsPtr_
44  )
45  {
47  << "Addressing already calculated"
48  << abort(FatalError);
49  }
50 
51  // Mapping
52 
53  const label oldNInternal = faceMap_.nOldInternalFaces();
54 
55  // Assemble the maps
56  if (direct())
57  {
58  // Direct mapping - slice to size
59  directAddrPtr_ =
60  new labelList
61  (
63  );
64  labelList& addr = *directAddrPtr_;
65 
66  // Adjust for creation of an internal face from a boundary face
67  forAll(addr, facei)
68  {
69  if (addr[facei] > oldNInternal)
70  {
71  addr[facei] = 0;
72  }
73  }
74  }
75  else
76  {
77  // Interpolative mapping - slice to size
78  interpolationAddrPtr_ =
79  new labelListList
80  (
82  );
83  labelListList& addr = *interpolationAddrPtr_;
84 
85  weightsPtr_ =
86  new scalarListList
87  (
88  scalarListList::subList(faceMap_.weights(), size())
89  );
90  scalarListList& w = *weightsPtr_;
91 
92  // Adjust for creation of an internal face from a boundary face
93  forAll(addr, facei)
94  {
95  if (max(addr[facei]) >= oldNInternal)
96  {
97  addr[facei] = labelList(1, label(0));
98  w[facei] = scalarList(1, 1.0);
99  }
100  }
101  }
102 
103  // Inserted objects
104 
105  // If there are, assemble the labels
106  if (insertedObjects())
107  {
108  const labelList& insFaces = faceMap_.insertedObjectLabels();
109 
110  insertedObjectLabelsPtr_ = new labelList(insFaces.size());
111  labelList& ins = *insertedObjectLabelsPtr_;
112 
113  label nIns = 0;
114 
115  forAll(insFaces, facei)
116  {
117  // If the face is internal, keep it here
118  if (insFaces[facei] < size())
119  {
120  ins[nIns] = insFaces[facei];
121  nIns++;
122  }
123  }
124 
125  ins.setSize(nIns);
126  }
127  else
128  {
129  // No inserted objects
130  insertedObjectLabelsPtr_ = new labelList(0);
131  }
132 }
133 
134 
135 void Foam::fvSurfaceMapper::clearOut()
136 {
137  deleteDemandDrivenData(directAddrPtr_);
138  deleteDemandDrivenData(interpolationAddrPtr_);
139  deleteDemandDrivenData(weightsPtr_);
140 
141  deleteDemandDrivenData(insertedObjectLabelsPtr_);
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
146 
147 Foam::fvSurfaceMapper::fvSurfaceMapper
148 (
149  const fvMesh& mesh,
150  const faceMapper& fMapper
151 )
152 :
153  mesh_(mesh),
154  faceMap_(fMapper),
155  directAddrPtr_(NULL),
156  interpolationAddrPtr_(NULL),
157  weightsPtr_(NULL),
158  insertedObjectLabelsPtr_(NULL)
159 {}
160 
161 
162 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
163 
165 {
166  clearOut();
167 }
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
173 {
174  if (!direct())
175  {
177  << "Requested direct addressing for an interpolative mapper."
178  << abort(FatalError);
179  }
180 
181  if (!directAddrPtr_)
182  {
183  calcAddressing();
184  }
185 
186  return *directAddrPtr_;
187 }
188 
189 
191 {
192  if (direct())
193  {
195  << "Requested interpolative addressing for a direct mapper."
196  << abort(FatalError);
197  }
198 
199  if (!interpolationAddrPtr_)
200  {
201  calcAddressing();
202  }
203 
204  return *interpolationAddrPtr_;
205 }
206 
207 
209 {
210  if (direct())
211  {
213  << "Requested interpolative weights for a direct mapper."
214  << abort(FatalError);
215  }
216 
217  if (!weightsPtr_)
218  {
219  calcAddressing();
220  }
221 
222  return *weightsPtr_;
223 }
224 
225 
227 {
228  if (!insertedObjectLabelsPtr_)
229  {
230  calcAddressing();
231  }
232 
233  return *insertedObjectLabelsPtr_;
234 }
235 
236 
237 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
238 
239 
240 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
241 
242 
243 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
244 
245 
246 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
virtual label nOldInternalFaces() const
Return number of old internalFaces.
Definition: faceMapper.C:388
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
This object provides mapping and fill-in information for face data between the two meshes after the t...
Definition: faceMapper.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual bool insertedObjects() const
Are there any inserted faces.
virtual bool direct() const
Is the mapping direct.
virtual label size() const
Return size.
SubList< label > subList
Declare type of subList.
Definition: List.H:155
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
virtual ~fvSurfaceMapper()
Destructor.
virtual const scalarListList & weights() const
Return interpolaion weights.
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faceMapper.C:327
List< scalarList > scalarListList
Definition: scalarList.H:51
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
virtual const labelUList & directAddressing() const
Return direct addressing.
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faceMapper.C:301
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual const labelListList & addressing() const
Return interpolated addressing.
virtual const scalarListList & weights() const
Return interpolaion weights.
Definition: faceMapper.C:345
void deleteDemandDrivenData(DataPtr &dataPtr)
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
Definition: faceMapper.C:363