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