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