faAreaMapper.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 Wikki Ltd
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Description
27  FV surface mapper.
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "faAreaMapper.H"
32 #include "mapPolyMesh.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 void Foam::faAreaMapper::calcAddressing() const
37 {
38  if
39  (
40  newFaceLabelsPtr_
41  || newFaceLabelsMapPtr_
42  || directAddrPtr_
43  || interpolationAddrPtr_
44  || weightsPtr_
45  || insertedObjectLabelsPtr_
46  )
47  {
49  << "Addressing already calculated"
50  << abort(FatalError);
51  }
52 
53  // Mapping
54 
55  const label oldNInternal = mpm_.nOldInternalFaces();
56 
57  hasUnmapped_ = false;
58 
59  // Calculate new face labels
60 
61  // Copy old face labels
62  const labelList& oldFaces = mesh_.faceLabels();
63 
64  // Prepare a list of new face labels and (preliminary) addressing
65  // Note: dimensioned to number of boundary faces of polyMesh
66  newFaceLabelsPtr_ = new labelList(mesh_().nBoundaryFaces(), -1);
67  labelList& newFaceLabels = *newFaceLabelsPtr_;
68 
69  newFaceLabelsMapPtr_ = new labelList(mesh_().nBoundaryFaces(), -1);
70  labelList& newFaceLabelsMap = *newFaceLabelsMapPtr_;
71  label nNewFaces = 0;
72 
73  Info<< "Old face list size: " << oldFaces.size()
74  << " estimated new size " << newFaceLabels.size() << endl;
75 
76  // Get reverse face map
77  const labelList& reverseFaceMap = mpm_.reverseFaceMap();
78 
79  // Pick up live old faces
80  forAll(oldFaces, faceI)
81  {
82  if (reverseFaceMap[oldFaces[faceI]] > -1)
83  {
84  // Face is live, add it and record addressing
85  newFaceLabels[nNewFaces] = reverseFaceMap[oldFaces[faceI]];
86  newFaceLabelsMap[nNewFaces] = faceI;
87 
88  nNewFaces++;
89  }
90  }
91 
92  // Assemble the maps
93  if (direct())
94  {
95  Info<< "Direct"<< endl;
96  // Direct mapping: no further faces to add. Resize list
97  newFaceLabels.setSize(nNewFaces);
98 
99  directAddrPtr_ = new labelList(newFaceLabels.size());
100  labelList& addr = *directAddrPtr_;
101 
102  // Adjust for creation of a boundary face from an internal face
103  forAll(addr, faceI)
104  {
105  if (newFaceLabelsMap[faceI] < oldNInternal)
106  {
107  addr[faceI] = 0;
108  }
109  else
110  {
111  addr[faceI] = newFaceLabelsMap[faceI];
112  }
113  }
114  }
115  else
116  {
117  // There are further faces to add. Prepare interpolation addressing
118  // and weights to full size
119  interpolationAddrPtr_ = new labelListList(newFaceLabels.size());
120  labelListList& addr = *interpolationAddrPtr_;
121 
122  weightsPtr_ = new scalarListList(newFaceLabels.size());
123  scalarListList& w = *weightsPtr_;
124 
125  // Insert single addressing and weights
126  for (label addrI = 0; addrI < nNewFaces; ++addrI)
127  {
128  addr[addrI] = labelList(1, newFaceLabelsMap[addrI]);
129  w[addrI] = scalarList(1, scalar(1));
130  }
131 
132  // Pick up faces from points, edges and faces where the origin
133  // Only map from faces which were previously in the faMesh, using
134  // fast lookup
135 
136  // Set of faces previously in the mesh
137  labelHashSet oldFaceLookup(oldFaces);
138 
139  // Go through faces-from lists and add the ones where all
140  // old face labels belonged to the faMesh
141 
142  const List<objectMap>& ffp = mpm_.facesFromPointsMap();
143 
144  forAll(ffp, ffpI)
145  {
146  // Get addressing
147  const labelList& mo = ffp[ffpI].masterObjects();
148 
149  // Check if master objects are in faMesh
150  labelList validMo(mo.size());
151  label nValidMo = 0;
152 
153  forAll(mo, moI)
154  {
155  if (oldFaceLookup.found(mo[moI]))
156  {
157  validMo[nValidMo] = oldFaceLookup[mo[moI]];
158  nValidMo++;
159  }
160  }
161 
162  if (nValidMo > 0)
163  {
164  // Some objects found: add face and interpolation to list
165  newFaceLabels[nNewFaces] = ffp[ffpI].index();
166 
167  // No old face available
168  newFaceLabelsMap[nNewFaces] = -1;
169 
170  // Map from masters, uniform weights
171  addr[nNewFaces] = validMo;
172  w[nNewFaces] = scalarList(validMo.size(), 1.0/validMo.size());
173 
174  nNewFaces++;
175  }
176  }
177 
178  const List<objectMap>& ffe = mpm_.facesFromEdgesMap();
179 
180  forAll(ffe, ffeI)
181  {
182  // Get addressing
183  const labelList& mo = ffe[ffeI].masterObjects();
184 
185  // Check if master objects are in faMesh
186  labelList validMo(mo.size());
187  label nValidMo = 0;
188 
189  forAll(mo, moI)
190  {
191  if (oldFaceLookup.found(mo[moI]))
192  {
193  validMo[nValidMo] = oldFaceLookup[mo[moI]];
194  nValidMo++;
195  }
196  }
197 
198  if (nValidMo > 0)
199  {
200  // Some objects found: add face and interpolation to list
201  newFaceLabels[nNewFaces] = ffe[ffeI].index();
202 
203  // No old face available
204  newFaceLabelsMap[nNewFaces] = -1;
205 
206  // Map from masters, uniform weights
207  addr[nNewFaces] = validMo;
208  w[nNewFaces] = scalarList(validMo.size(), 1.0/validMo.size());
209 
210  nNewFaces++;
211  }
212  }
213 
214  const List<objectMap>& fff = mpm_.facesFromFacesMap();
215 
216  forAll(fff, fffI)
217  {
218  // Get addressing
219  const labelList& mo = fff[fffI].masterObjects();
220 
221  // Check if master objects are in faMesh
222  labelList validMo(mo.size());
223  label nValidMo = 0;
224 
225  forAll(mo, moI)
226  {
227  if (oldFaceLookup.found(mo[moI]))
228  {
229  validMo[nValidMo] = oldFaceLookup[mo[moI]];
230  nValidMo++;
231  }
232  }
233 
234  if (nValidMo > 0)
235  {
236  // Some objects found: add face and interpolation to list
237  newFaceLabels[nNewFaces] = fff[fffI].index();
238 
239  // No old face available
240  newFaceLabelsMap[nNewFaces] = -1;
241 
242  // Map from masters, uniform weights
243  addr[nNewFaces] = validMo;
244  w[nNewFaces] = scalarList(validMo.size(), 1.0/validMo.size());
245 
246  nNewFaces++;
247  }
248  }
249 
250  // All faces collected. Reset sizes of lists
251  newFaceLabels.setSize(nNewFaces);
252  newFaceLabelsMap.setSize(nNewFaces);
253  addr.setSize(nNewFaces);
254  w.setSize(nNewFaces);
255  Info<< "addr: " << addr << nl
256  << "w: " << w << endl;
257  }
258 
259  // Inserted objects cannot appear in the new faMesh as they have no master
260  // HJ, 10/Aug/2011
261  insertedObjectLabelsPtr_ = new labelList(0);
262 }
263 
264 
265 void Foam::faAreaMapper::clearOut()
266 {
267  deleteDemandDrivenData(newFaceLabelsPtr_);
268  deleteDemandDrivenData(newFaceLabelsMapPtr_);
269 
270  deleteDemandDrivenData(directAddrPtr_);
271  deleteDemandDrivenData(interpolationAddrPtr_);
272  deleteDemandDrivenData(weightsPtr_);
273 
274  deleteDemandDrivenData(insertedObjectLabelsPtr_);
275  hasUnmapped_ = false;
276 }
277 
278 
279 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
280 
281 Foam::faAreaMapper::faAreaMapper
282 (
283  const faMesh& mesh,
284  const mapPolyMesh& mpm
285 )
286 :
287  mesh_(mesh),
288  mpm_(mpm),
289  direct_(false),
290  hasUnmapped_(false),
291  sizeBeforeMapping_(mesh.nFaces()),
292  newFaceLabelsPtr_(nullptr),
293  newFaceLabelsMapPtr_(nullptr),
294  directAddrPtr_(nullptr),
295  interpolationAddrPtr_(nullptr),
296  weightsPtr_(nullptr),
297  insertedObjectLabelsPtr_(nullptr)
298 {
299  // Check for possibility of direct mapping
300  if
301  (
302  mpm_.facesFromPointsMap().empty()
303  && mpm_.facesFromEdgesMap().empty()
304  && mpm_.facesFromFacesMap().empty()
305  )
306  {
307  direct_ = true;
308  }
309  else
310  {
311  direct_ = false;
312  }
313 
314  // Inserted objects not supported: no master
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
319 
321 {
322  clearOut();
323 }
324 
325 
326 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
327 
329 {
330  if (!newFaceLabelsPtr_)
331  {
332  calcAddressing();
333  }
334 
335  return *newFaceLabelsPtr_;
336 }
337 
338 
340 {
341  if (!newFaceLabelsMapPtr_)
342  {
343  calcAddressing();
344  }
345 
346  return *newFaceLabelsMapPtr_;
347 }
348 
349 
351 {
352  if (!direct())
353  {
355  << "Requested direct addressing for an interpolative mapper."
356  << abort(FatalError);
357  }
358 
359  if (!directAddrPtr_)
360  {
361  calcAddressing();
362  }
363 
364  return *directAddrPtr_;
365 }
366 
367 
369 {
370  if (direct())
371  {
373  << "Requested interpolative addressing for a direct mapper."
374  << abort(FatalError);
375  }
376 
377  if (!interpolationAddrPtr_)
378  {
379  calcAddressing();
380  }
381 
382  return *interpolationAddrPtr_;
383 }
384 
385 
387 {
388  if (direct())
389  {
391  << "Requested interpolative weights for a direct mapper."
392  << abort(FatalError);
393  }
394 
395  if (!weightsPtr_)
396  {
397  calcAddressing();
398  }
399 
400  return *weightsPtr_;
401 }
402 
403 
405 {
406  if (!insertedObjectLabelsPtr_)
407  {
408  calcAddressing();
409  }
410 
411  return *insertedObjectLabelsPtr_;
412 }
413 
414 
415 // ************************************************************************* //
Foam::mapPolyMesh::facesFromPointsMap
const List< objectMap > & facesFromPointsMap() const
Faces inflated from points.
Definition: mapPolyMesh.H:415
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::scalarListList
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
faAreaMapper.H
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
mapPolyMesh.H
Foam::faAreaMapper::addressing
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: faAreaMapper.C:368
Foam::faMesh::faceLabels
const labelList & faceLabels() const
Return faMesh face labels.
Definition: faMesh.H:424
Foam::faAreaMapper::newFaceLabels
const labelList & newFaceLabels() const
Return new face labels.
Definition: faAreaMapper.C:328
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::faAreaMapper::weights
virtual const scalarListList & weights() const
Return interpolation weights.
Definition: faAreaMapper.C:386
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::faAreaMapper::direct
virtual bool direct() const
Is the mapping direct.
Definition: faAreaMapper.H:153
Foam::label
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:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::faAreaMapper::insertedObjectLabels
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
Definition: faAreaMapper.C:404
Foam::faAreaMapper::newFaceLabelsMap
const labelList & newFaceLabelsMap() const
Return new face labels map.
Definition: faAreaMapper.C:339
Foam::mapPolyMesh::facesFromEdgesMap
const List< objectMap > & facesFromEdgesMap() const
Faces inflated from edges.
Definition: mapPolyMesh.H:421
Foam::mapPolyMesh::facesFromFacesMap
const List< objectMap > & facesFromFacesMap() const
Faces originating from faces.
Definition: mapPolyMesh.H:427
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:500
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::faAreaMapper::~faAreaMapper
virtual ~faAreaMapper()
Destructor.
Definition: faAreaMapper.C:320
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::List< label >
Foam::UList< label >
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::faMesh
Finite area mesh. Used for 2-D non-Euclidian finite area method.
Definition: faMesh.H:77
Foam::mapPolyMesh::nOldInternalFaces
label nOldInternalFaces() const
Number of old internal faces.
Definition: mapPolyMesh.H:374
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:415
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::faAreaMapper::directAddressing
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: faAreaMapper.C:350