pointMapper.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) 2011-2016 OpenFOAM Foundation
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 "pointMapper.H"
29#include "demandDrivenData.H"
30#include "pointMesh.H"
31#include "mapPolyMesh.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35void Foam::pointMapper::calcAddressing() const
36{
37 if
38 (
39 directAddrPtr_
40 || interpolationAddrPtr_
41 || weightsPtr_
42 || insertedPointLabelsPtr_
43 )
44 {
46 << "Addressing already calculated."
47 << abort(FatalError);
48 }
49
50 if (direct())
51 {
52 // Direct addressing, no weights
53
54 directAddrPtr_ = new labelList(mpm_.pointMap());
55 labelList& directAddr = *directAddrPtr_;
56
57 // Not necessary to resize the list as there are no retired points
58 // directAddr.setSize(pMesh_.size());
59
60 insertedPointLabelsPtr_ = new labelList(pMesh_.size());
61 labelList& insertedPoints = *insertedPointLabelsPtr_;
62
63 label nInsertedPoints = 0;
64
65 forAll(directAddr, pointi)
66 {
67 if (directAddr[pointi] < 0)
68 {
69 // Found inserted point
70 directAddr[pointi] = 0;
71 insertedPoints[nInsertedPoints] = pointi;
72 nInsertedPoints++;
73 }
74 }
75
76 insertedPoints.setSize(nInsertedPoints);
77 }
78 else
79 {
80 // Interpolative addressing
81
82 interpolationAddrPtr_ = new labelListList(pMesh_.size());
83 labelListList& addr = *interpolationAddrPtr_;
84
85 weightsPtr_ = new scalarListList(pMesh_.size());
86 scalarListList& w = *weightsPtr_;
87
88 // Points created from other points (i.e. points merged into it).
89 const List<objectMap>& cfc = mpm_.pointsFromPointsMap();
90
91 forAll(cfc, cfcI)
92 {
93 // Get addressing
94 const labelList& mo = cfc[cfcI].masterObjects();
95
96 label pointi = cfc[cfcI].index();
97
98 if (addr[pointi].size())
99 {
101 << "Master point " << pointi
102 << " mapped from points " << mo
103 << " already destination of mapping." << abort(FatalError);
104 }
105
106 // Map from masters, uniform weights
107 addr[pointi] = mo;
108 w[pointi] = scalarList(mo.size(), 1.0/mo.size());
109 }
110
111
112 // Do mapped points. Note that can already be set from pointsFromPoints
113 // so check if addressing size still zero.
114
115 const labelList& cm = mpm_.pointMap();
116
117 forAll(cm, pointi)
118 {
119 if (cm[pointi] > -1 && addr[pointi].empty())
120 {
121 // Mapped from a single point
122 addr[pointi] = labelList(1, cm[pointi]);
123 w[pointi] = scalarList(1, scalar(1));
124 }
125 }
126
127 // Grab inserted points (for them the size of addressing is still zero)
128
129 insertedPointLabelsPtr_ = new labelList(pMesh_.size());
130 labelList& insertedPoints = *insertedPointLabelsPtr_;
131
132 label nInsertedPoints = 0;
133
134 forAll(addr, pointi)
135 {
136 if (addr[pointi].empty())
137 {
138 // Mapped from a dummy point. Take point 0 with weight 1.
139 addr[pointi] = labelList(1, Zero);
140 w[pointi] = scalarList(1, scalar(1));
141
142 insertedPoints[nInsertedPoints] = pointi;
143 nInsertedPoints++;
144 }
145 }
146
147 insertedPoints.setSize(nInsertedPoints);
148 }
149}
150
151
152void Foam::pointMapper::clearOut()
153{
154 deleteDemandDrivenData(directAddrPtr_);
155 deleteDemandDrivenData(interpolationAddrPtr_);
156 deleteDemandDrivenData(weightsPtr_);
157 deleteDemandDrivenData(insertedPointLabelsPtr_);
158}
159
160
161// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
162
164:
165 pMesh_(pMesh),
166 mpm_(mpm),
167 insertedPoints_(true),
168 direct_(false),
169 directAddrPtr_(nullptr),
170 interpolationAddrPtr_(nullptr),
171 weightsPtr_(nullptr),
172 insertedPointLabelsPtr_(nullptr)
173{
174 // Check for possibility of direct mapping
175 if (mpm_.pointsFromPointsMap().empty())
176 {
177 direct_ = true;
178 }
179 else
180 {
181 direct_ = false;
182 }
183
184 // Check for inserted points
185 if (direct_ && (mpm_.pointMap().empty() || min(mpm_.pointMap()) > -1))
186 {
187 insertedPoints_ = false;
188 }
189 else
190 {
191 //Check if there are inserted points with no owner
192
193 // Make a copy of the point map, add the entries for points from points
194 // and check for left-overs
195 labelList cm(pMesh_.size(), -1);
196
197 const List<objectMap>& cfc = mpm_.pointsFromPointsMap();
198
199 forAll(cfc, cfcI)
200 {
201 cm[cfc[cfcI].index()] = 0;
202 }
203
204 if (min(cm) < 0)
205 {
206 insertedPoints_ = true;
207 }
208 }
209}
210
211
212// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
213
215{
216 clearOut();
217}
218
219
220// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221
222Foam::label Foam::pointMapper::size() const
223{
224 return mpm_.pointMap().size();
225}
226
227
229{
230 return mpm_.nOldPoints();
231}
232
233
235{
236 if (!direct())
237 {
239 << "Requested direct addressing for an interpolative mapper."
240 << abort(FatalError);
241 }
242
243 if (!insertedObjects())
244 {
245 // No inserted points. Re-use pointMap
246 return mpm_.pointMap();
247 }
248 else
249 {
250 if (!directAddrPtr_)
251 {
252 calcAddressing();
253 }
254
255 return *directAddrPtr_;
256 }
257}
258
259
261{
262 if (direct())
263 {
265 << "Requested interpolative addressing for a direct mapper."
266 << abort(FatalError);
267 }
268
269 if (!interpolationAddrPtr_)
270 {
271 calcAddressing();
272 }
273
274 return *interpolationAddrPtr_;
275}
276
277
279{
280 if (direct())
281 {
283 << "Requested interpolative weights for a direct mapper."
284 << abort(FatalError);
285 }
286
287 if (!weightsPtr_)
288 {
289 calcAddressing();
290 }
291
292 return *weightsPtr_;
293}
294
295
297{
298 if (!insertedPointLabelsPtr_)
299 {
300 if (!insertedObjects())
301 {
302 // There are no inserted points
303 insertedPointLabelsPtr_ = new labelList(0);
304 }
305 else
306 {
307 calcAddressing();
308 }
309 }
310
311 return *insertedPointLabelsPtr_;
312}
313
314
315// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
316
317
318// * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
319
320
321// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
322
323
324// ************************************************************************* //
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
const List< objectMap > & pointsFromPointsMap() const
Points originating from points.
Definition: mapPolyMesh.H:402
const labelList & pointMap() const
Old point map.
Definition: mapPolyMesh.H:396
This object provides mapping and fill-in information for point data between the two meshes after the ...
Definition: pointMapper.H:60
virtual const labelListList & addressing() const
Return interpolated addressing.
Definition: pointMapper.C:260
virtual const scalarListList & weights() const
Return interpolation weights.
Definition: pointMapper.C:278
virtual const labelUList & directAddressing() const
Return direct addressing.
Definition: pointMapper.C:234
virtual label size() const
Return size.
Definition: pointMapper.C:222
const labelList & insertedObjectLabels() const
Return list of inserted points.
Definition: pointMapper.C:296
virtual ~pointMapper()
Destructor.
Definition: pointMapper.C:214
virtual label sizeBeforeMapping() const
Return size before mapping.
Definition: pointMapper.C:228
virtual bool direct() const
Is the mapping direct.
Definition: pointMapper.H:128
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:55
static label size(const Mesh &mesh)
Return size. Number of points.
Definition: pointMesh.H:102
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
List< label > labelList
A List of labels.
Definition: List.H:66
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
List< scalarList > scalarListList
A List of scalarList.
Definition: scalarList.H:66
void deleteDemandDrivenData(DataPtr &dataPtr)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333