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