Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33void Foam::faAreaMapper::calcAddressing() const
34{
35 if
36 (
37 newFaceLabelsPtr_
38 || newFaceLabelsMapPtr_
39 || directAddrPtr_
40 || interpAddrPtr_
41 || weightsPtr_
42 || insertedObjectsPtr_
43 )
44 {
46 << "Addressing already calculated"
47 << abort(FatalError);
48 }
49
50 // Mapping
51
52 const label oldNInternal = mpm_.nOldInternalFaces();
53
54 hasUnmapped_ = false;
55
56 // Calculate new face labels
57
58 // Copy old face labels
59 const labelList& oldFaces = mesh_.faceLabels();
60
61 // Prepare a list of new face labels and (preliminary) addressing
62 // Note: dimensioned to number of boundary faces of polyMesh
63 newFaceLabelsPtr_ = std::make_unique<labelList>
64 (
65 mesh_.mesh().nBoundaryFaces(),
66 -1
67 );
68 auto& newFaceLabels = *newFaceLabelsPtr_;
69
70 newFaceLabelsMapPtr_ = std::make_unique<labelList>
71 (
72 mesh_.mesh().nBoundaryFaces(),
73 -1
74 );
75 auto& newFaceLabelsMap = *newFaceLabelsMapPtr_;
76 label nNewFaces = 0;
77
78 Info<< "Old face list size: " << oldFaces.size()
79 << " estimated new size " << newFaceLabels.size() << endl;
80
81 // Get reverse face map
82 const labelList& reverseFaceMap = mpm_.reverseFaceMap();
83
84 // Pick up live old faces
85 forAll(oldFaces, faceI)
86 {
87 if (reverseFaceMap[oldFaces[faceI]] > -1)
88 {
89 // Face is live, add it and record addressing
90 newFaceLabels[nNewFaces] = reverseFaceMap[oldFaces[faceI]];
91 newFaceLabelsMap[nNewFaces] = faceI;
92
93 ++nNewFaces;
94 }
95 }
96
97 // Assemble the maps
98 if (direct())
99 {
100 Info<< "Direct"<< endl;
101 // Direct mapping: no further faces to add. Resize list
102 newFaceLabels.resize(nNewFaces);
103
104 directAddrPtr_ = std::make_unique<labelList>
105 (
107 );
108 auto& addr = *directAddrPtr_;
109
110 // Adjust for creation of a boundary face from an internal face
111 forAll(addr, facei)
112 {
113 if (addr[facei] < oldNInternal)
114 {
115 addr[facei] = 0;
116 }
117 }
118 }
119 else
120 {
121 // There are further faces to add. Prepare interpolation addressing
122 // and weights to full size
123 interpAddrPtr_ = std::make_unique<labelListList>
124 (
125 newFaceLabels.size()
126 );
127 auto& addr = *interpAddrPtr_;
128
129 weightsPtr_ = std::make_unique<scalarListList>(addr.size());
130 auto& wght = *weightsPtr_;
131
132 // Insert single addressing and weights
133 for (label addrI = 0; addrI < nNewFaces; ++addrI)
134 {
135 addr[addrI].resize(1, newFaceLabelsMap[addrI]);
136 wght[addrI].resize(1, 1.0);
137 }
138
139 // Pick up faces from points, edges and faces where the origin
140 // Only map from faces which were previously in the faMesh, using
141 // fast lookup
142
143 // Set of faces previously in the mesh
144 const labelHashSet oldFaceLookup(oldFaces);
145
146 // Check if master objects are in faMesh
147 DynamicList<label> validMo(128);
148
149 const auto addCheckedObjects = [&](const List<objectMap>& maps)
150 {
151 for (const objectMap& map : maps)
152 {
153 // Get target index, addressing
154 const label facei = map.index();
155 const labelList& mo = map.masterObjects();
156 if (mo.empty()) continue; // safety
157
158 validMo.clear();
159 validMo.reserve(mo.size());
160
161 for (const label obji : mo)
162 {
163 if (oldFaceLookup.contains(obji))
164 {
165 validMo.push_back(obji);
166 }
167 }
168
169 if (validMo.size())
170 {
171 // Some objects found: add face and interpolation to list
172 newFaceLabels[nNewFaces] = facei;
173
174 // No old face available
175 newFaceLabelsMap[nNewFaces] = -1;
176
177 // Map from masters, uniform weights
178 addr[nNewFaces] = validMo;
179 wght[nNewFaces] =
180 scalarList(validMo.size(), 1.0/validMo.size());
181
182 ++nNewFaces;
183 }
184 }
185 };
186
187
188 // Go through faces-from lists and add the ones where all
189 // old face labels belonged to the faMesh
190
191 {
192 addCheckedObjects(mpm_.facesFromPointsMap());
193 addCheckedObjects(mpm_.facesFromEdgesMap());
194 addCheckedObjects(mpm_.facesFromFacesMap());
195 }
196
197 // All faces collected. Reset sizes of lists
198 newFaceLabels.resize(nNewFaces);
199 newFaceLabelsMap.resize(nNewFaces);
200 addr.resize(nNewFaces);
201 wght.resize(nNewFaces);
202
203 Info<< "addr: " << addr << nl
204 << "wght: " << wght << endl;
205 }
206
207 // Inserted objects cannot appear in the new faMesh as they have no master
208 // HJ, 10/Aug/2011
209 insertedObjectsPtr_ = std::make_unique<labelList>();
210}
211
212
213// void Foam::faAreaMapper::clearOut()
214// {
215// newFaceLabelsPtr_.reset(nullptr);
216// newFaceLabelsMapPtr_.reset(nullptr);
217//
218// directAddrPtr_.reset(nullptr);
219// interpAddrPtr_.reset(nullptr);
220// weightsPtr_.reset(nullptr);
221//
222// insertedObjectsPtr_.reset(nullptr);
223// hasUnmapped_ = false;
224// }
225
226
227// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
228
230(
231 const faMesh& mesh,
232 const mapPolyMesh& mpm
233)
234:
235 mesh_(mesh),
236 mpm_(mpm),
237 sizeBeforeMapping_(mesh.nFaces()),
238 direct_
239 (
240 // Mapping without interpolation?
241 mpm.facesFromPointsMap().empty()
242 && mpm.facesFromEdgesMap().empty()
243 && mpm.facesFromFacesMap().empty()
244 ),
245 hasUnmapped_(false)
247 // Inserted objects not supported: no master
248}
249
250
251// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
254{}
255
256
257// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
258
260{
261 if (!newFaceLabelsPtr_)
262 {
263 calcAddressing();
264 }
265
266 return *newFaceLabelsPtr_;
267}
268
269
271{
272 if (!newFaceLabelsMapPtr_)
273 {
274 calcAddressing();
275 }
276
277 return *newFaceLabelsMapPtr_;
278}
279
280
282{
283 if (!direct())
284 {
286 << "Requested direct addressing for an interpolative mapper."
287 << abort(FatalError);
288 }
289
290 if (!directAddrPtr_)
291 {
292 calcAddressing();
293 }
294
295 return *directAddrPtr_;
296}
297
298
300{
301 if (direct())
302 {
304 << "Requested interpolative addressing for a direct mapper."
305 << abort(FatalError);
306 }
307
308 if (!interpAddrPtr_)
309 {
310 calcAddressing();
311 }
312
313 return *interpAddrPtr_;
314}
315
316
318{
319 if (direct())
320 {
322 << "Requested interpolative weights for a direct mapper."
323 << abort(FatalError);
324 }
325
326 if (!weightsPtr_)
327 {
328 calcAddressing();
329 }
330
331 return *weightsPtr_;
332}
333
334
336{
337 if (!insertedObjectsPtr_)
338 {
339 calcAddressing();
340 }
341
342 return *insertedObjectsPtr_;
343}
344
345
346// ************************************************************************* //
SubList< label > subList
Definition List.H:129
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
virtual const labelListList & addressing() const
Return interpolated addressing.
virtual const scalarListList & weights() const
Return interpolation weights.
faAreaMapper(const faAreaMapper &)=delete
No copy construct.
virtual const labelUList & directAddressing() const
Return direct addressing.
virtual const labelList & insertedObjectLabels() const
Return list of inserted faces.
const labelList & newFaceLabelsMap() const
Return new face labels map.
virtual ~faAreaMapper()
Destructor.
const labelList & newFaceLabels() const
Return new face labels.
virtual bool direct() const
Is the mapping direct.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition faMesh.H:140
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299