meshSurfaceMapperMapVertices.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6  \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh is free software; you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation; either version 3 of the License, or (at your
14  option) any later version.
15 
16  cfMesh 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 cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "demandDrivenData.H"
30 #include "meshSurfaceMapper.H"
31 #include "meshSurfacePartitioner.H"
32 #include "meshOctree.H"
33 #include "triSurf.H"
34 #include "helperFunctionsPar.H"
35 #include "helperFunctions.H"
36 
37 #include <map>
38 
39 # ifdef USE_OMP
40 #include <omp.h>
41 # endif
42 
43 //#define DEBUGMapping
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 // Private member functions
52 
54 {
55  if( !Pstream::parRun() )
56  return;
57 
58  std::map<label, labelLongList> exchangeData;
59  const DynList<label>& neiProcs = surfaceEngine_.bpNeiProcs();
60  forAll(neiProcs, i)
61  exchangeData.insert(std::make_pair(neiProcs[i], labelLongList()));
62 
63  const VRWGraph& bpAtProcs = surfaceEngine_.bpAtProcs();
64  const labelList& globalPointLabel =
66  const Map<label>& globalToLocal =
68 
69  boolList selectedNode(bpAtProcs.size(), false);
70 
71  forAll(selNodes, i)
72  {
73  const label bpI = selNodes[i];
74 
75  selectedNode[bpI] = true;
76 
77  forAllRow(bpAtProcs, bpI, procI)
78  {
79  const label neiProc = bpAtProcs(bpI, procI);
80  if( neiProc == Pstream::myProcNo() )
81  continue;
82 
83  exchangeData[neiProc].append(globalPointLabel[bpI]);
84  }
85  }
86 
87  //- exchange data
88  labelLongList receivedData;
89  help::exchangeMap(exchangeData, receivedData);
90 
91  forAll(receivedData, i)
92  {
93  if( !selectedNode[globalToLocal[receivedData[i]]] )
94  {
95  selectedNode[globalToLocal[receivedData[i]]] = true;
96  const_cast<labelLongList&>(selNodes).append
97  (
98  globalToLocal[receivedData[i]]
99  );
100  }
101  }
102 }
103 
105 {
106  if( !Pstream::parRun() )
107  return;
108 
109  std::map<label, LongList<parMapperHelper> > exchangeData;
110  const DynList<label>& neiProcs = surfaceEngine_.bpNeiProcs();
111  forAll(neiProcs, i)
112  exchangeData.insert
113  (
114  std::make_pair(neiProcs[i], LongList<parMapperHelper>())
115  );
116 
117  const VRWGraph& bpAtProcs = surfaceEngine_.bpAtProcs();
118  const labelList& globalPointLabel =
120  const Map<label>& globalToLocal =
122 
123  Map<label> bpToList(parN.size());
124 
125  forAll(parN, i)
126  {
127  const label bpI = parN[i].globalLabel();
128  bpToList.insert(bpI, i);
129 
130  forAllRow(bpAtProcs, bpI, procI)
131  {
132  const label neiProc = bpAtProcs(bpI, procI);
133  if( neiProc == Pstream::myProcNo() )
134  continue;
135 
136  exchangeData[neiProc].append
137  (
139  (
140  parN[i].coordinates(),
141  parN[i].movingDistance(),
142  globalPointLabel[bpI],
143  parN[i].pointPatch()
144  )
145  );
146  }
147  }
148 
149  //- exchange data
150  LongList<parMapperHelper> receivedData;
151  help::exchangeMap(exchangeData, receivedData);
152 
153  //- select the point with the smallest moving distance
155  forAll(receivedData, i)
156  {
157  const parMapperHelper& ph = receivedData[i];
158 
159  const label bpI = globalToLocal[ph.globalLabel()];
160 
161  parMapperHelper& phOrig = parN[bpToList[bpI]];
162  if( phOrig.movingDistance() < ph.movingDistance() )
163  {
164  surfModifier.moveBoundaryVertexNoUpdate(bpI, ph.coordinates());
165  phOrig = ph;
166  }
167  }
168 }
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 void meshSurfaceMapper::mapNodeToPatch(const label bpI, const label patchI)
173 {
174  label patch, nt;
175  point mapPoint;
176  scalar dSq;
177 
179  const labelList& bPoints = surfaceEngine_.boundaryPoints();
180  const point p = points[bPoints[bpI]];
181 
182  if( patchI < 0 )
183  {
184  meshOctree_.findNearestSurfacePoint(mapPoint, dSq, nt, patch, p);
185  }
186  else
187  {
189  (
190  mapPoint,
191  dSq,
192  nt,
193  patchI,
194  p
195  );
196  }
197 
199  surfModifier.moveBoundaryVertex(bpI, mapPoint);
200 }
201 
203 {
204  Info << "Mapping vertices onto surface" << endl;
205 
207  forAll(nodesToMap, i)
208  nodesToMap[i] = i;
209 
210  mapVerticesOntoSurface(nodesToMap);
211 
212  Info << "Finished mapping vertices onto surface" << endl;
213 }
214 
216 {
217  const labelList& boundaryPoints = surfaceEngine_.boundaryPoints();
219 
220  const VRWGraph* bpAtProcsPtr(NULL);
221  if( Pstream::parRun() )
222  bpAtProcsPtr = &surfaceEngine_.bpAtProcs();
223 
224  meshSurfaceEngineModifier surfaceModifier(surfaceEngine_);
225  LongList<parMapperHelper> parallelBndNodes;
226 
227  # ifdef USE_OMP
228  const label size = nodesToMap.size();
229  # pragma omp parallel for if( size > 1000 ) shared(parallelBndNodes) \
230  schedule(dynamic, Foam::max(1, size / (3 * omp_get_max_threads())))
231  # endif
232  forAll(nodesToMap, i)
233  {
234  const label bpI = nodesToMap[i];
235 
236  # ifdef DEBUGMapping
237  Info << nl << "Mapping vertex " << bpI << " with coordinates "
238  << points[boundaryPoints[bpI]] << endl;
239  # endif
240 
241  label patch, nt;
242  point mapPoint;
243  scalar dSq;
244 
246  (
247  mapPoint,
248  dSq,
249  nt,
250  patch,
251  points[boundaryPoints[bpI]]
252  );
253 
254  surfaceModifier.moveBoundaryVertexNoUpdate(bpI, mapPoint);
255 
256  if( bpAtProcsPtr && bpAtProcsPtr->sizeOfRow(bpI) )
257  {
258  # ifdef USE_OMP
259  # pragma omp critical
260  # endif
261  parallelBndNodes.append
262  (
264  (
265  mapPoint,
266  dSq,
267  bpI,
268  patch
269  )
270  );
271  }
272 
273  # ifdef DEBUGMapping
274  Info << "Mapped point " << points[boundaryPoints[bpI]] << endl;
275  # endif
276  }
277 
278  //- make sure that the points are at the nearest location on the surface
279  mapToSmallestDistance(parallelBndNodes);
280 
281  //- re-calculate face normals, point normals, etc.
282  surfaceModifier.updateGeometry(nodesToMap);
283 }
284 
286 {
287  Info << "Mapping vertices with respect to surface patches" << endl;
288 
290  forAll(nodesToMap, i)
291  nodesToMap[i] = i;
292 
293  mapVerticesOntoSurfacePatches(nodesToMap);
294 
295  Info << "Finished mapping vertices with respect to surface patches" << endl;
296 }
297 
299 (
300  const labelLongList& nodesToMap
301 )
302 {
303  const meshSurfacePartitioner& mPart = meshPartitioner();
304  const labelHashSet& cornerPoints = mPart.corners();
305  const labelHashSet& edgePoints = mPart.edgePoints();
306  const VRWGraph& pointPatches = mPart.pointPatches();
307 
308  boolList treatedPoint(surfaceEngine_.boundaryPoints().size(), false);
309 
310  //- find corner and edge points
311  labelLongList selectedCorners, selectedEdges;
312  forAll(nodesToMap, i)
313  {
314  if( cornerPoints.found(nodesToMap[i]) )
315  {
316  treatedPoint[nodesToMap[i]] = true;
317  selectedCorners.append(nodesToMap[i]);
318  }
319  else if( edgePoints.found(nodesToMap[i]) )
320  {
321  treatedPoint[nodesToMap[i]] = true;
322  selectedEdges.append(nodesToMap[i]);
323  }
324  }
325 
326  //- map the remaining selected points
327  const labelList& bPoints = surfaceEngine_.boundaryPoints();
328  const pointFieldPMG& points = surfaceEngine_.points();
329 
330  const VRWGraph* bpAtProcsPtr(NULL);
331  if( Pstream::parRun() )
332  bpAtProcsPtr = &surfaceEngine_.bpAtProcs();
333 
334  meshSurfaceEngineModifier surfaceModifier(surfaceEngine_);
335  LongList<parMapperHelper> parallelBndNodes;
336 
337  # ifdef USE_OMP
338  const label size = nodesToMap.size();
339  # pragma omp parallel for if( size > 1000 ) shared(parallelBndNodes) \
340  schedule(dynamic, Foam::max(1, size / (3 * omp_get_max_threads())))
341  # endif
342  forAll(nodesToMap, nI)
343  {
344  const label bpI = nodesToMap[nI];
345 
346  if( treatedPoint[bpI] )
347  continue;
348 
349  const point& p = points[bPoints[bpI]];
350  point mapPoint;
351  scalar dSq;
352  label nt;
353  meshOctree_.findNearestSurfacePointInRegion
354  (
355  mapPoint,
356  dSq,
357  nt,
358  pointPatches(bpI, 0),
359  p
360  );
361 
362  surfaceModifier.moveBoundaryVertexNoUpdate(bpI, mapPoint);
363 
364  if( bpAtProcsPtr && bpAtProcsPtr->sizeOfRow(bpI) )
365  {
366  # ifdef USE_OMP
367  # pragma omp critical
368  # endif
369  {
370  parallelBndNodes.append
371  (
373  (
374  mapPoint,
375  dSq,
376  bpI,
377  -1
378  )
379  );
380  }
381  }
382 
383  # ifdef DEBUGMapping
384  Info << "Mapped point " << points[boundaryPoints[bpI]] << endl;
385  # endif
386  }
387 
388  //- map vertices at inter-processor boundaries to the nearest location
389  //- on the surface
390  mapToSmallestDistance(parallelBndNodes);
391 
392  //- update face normals, point normals, etc.
393  surfaceModifier.updateGeometry(nodesToMap);
394 
395  //- map edge nodes
396  mapEdgeNodes(selectedEdges);
397 
398  //- map corner vertices
399  mapCorners(selectedCorners);
400 }
401 
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 
404 } // End namespace Foam
405 
406 // ************************************************************************* //
Foam::meshSurfaceEngine::bpAtProcs
const VRWGraph & bpAtProcs() const
processors which contain the vertex
Definition: meshSurfaceEngineI.H:451
Foam::LongList::append
void append(const T &e)
Append an element at the end of the list.
Definition: LongListI.H:265
Foam::meshSurfaceMapper::surfaceEngine_
const meshSurfaceEngine & surfaceEngine_
mesh surface
Definition: meshSurfaceMapper.H:63
triSurf.H
Foam::parMapperHelper::globalLabel
label globalLabel() const
return global label
Definition: parMapperHelper.H:102
Foam::parMapperHelper
Definition: parMapperHelper.H:51
p
p
Definition: pEqn.H:62
meshSurfaceMapper.H
helperFunctionsPar.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::meshSurfaceEngine::globalBoundaryPointLabel
const labelList & globalBoundaryPointLabel() const
global boundary point label
Definition: meshSurfaceEngineI.H:410
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::meshOctree::findNearestSurfacePoint
void findNearestSurfacePoint(point &nearest, scalar &distSq, label &nearestTriangle, label &region, const point &p) const
find nearest surface point for vertex and its region
Definition: meshOctreeFindNearestSurfacePoint.C:44
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::meshSurfaceEngine::bpNeiProcs
const DynList< label > & bpNeiProcs() const
communication matrix for sending point data
Definition: meshSurfaceEngineI.H:470
Foam::Map< label >
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::LongList::size
label size() const
Size of the active part of the list.
Definition: LongListI.H:203
Foam::pointPatch
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
Foam::HashSet< label, Hash< label > >
Foam::meshSurfacePartitioner::pointPatches
const VRWGraph & pointPatches() const
Definition: meshSurfacePartitioner.H:123
meshOctree.H
Foam::parMapperHelper::movingDistance
const scalar & movingDistance() const
return moving distance
Definition: parMapperHelper.H:96
Foam::LongList< label >
Foam::meshSurfaceEngineModifier
Definition: meshSurfaceEngineModifier.H:48
Foam::meshSurfaceMapper::mapToSmallestDistance
void mapToSmallestDistance(LongList< parMapperHelper > &)
Definition: meshSurfaceMapperMapVertices.C:104
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::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::meshSurfaceEngineModifier::moveBoundaryVertexNoUpdate
void moveBoundaryVertexNoUpdate(const label bpI, const point &newP)
relocate the selected boundary vertex
Definition: meshSurfaceEngineModifier.C:68
Foam::meshSurfaceMapper::mapNodeToPatch
void mapNodeToPatch(const label bpI, const label patchI=-1)
Definition: meshSurfaceMapperMapVertices.C:172
Foam::VRWGraph::size
label size() const
Returns the number of rows.
Definition: VRWGraphI.H:122
forAllRow
#define forAllRow(graph, rowI, index)
Definition: VRWGraph.H:277
Foam::meshSurfaceEngineModifier::moveBoundaryVertex
void moveBoundaryVertex(const label bpI, const point &newP)
relocate the selected boundary vertex and update geometry data
Definition: meshSurfaceEngineModifier.C:77
Foam::meshSurfaceEngine::points
const pointFieldPMG & points() const
Definition: meshSurfaceEngineI.H:49
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Foam::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::meshSurfaceMapper::mapVerticesOntoSurfacePatches
void mapVerticesOntoSurfacePatches()
Definition: meshSurfaceMapperMapVertices.C:285
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::meshSurfaceEngineModifier::updateGeometry
void updateGeometry(const labelLongList &)
Definition: meshSurfaceEngineModifier.C:234
Foam::DynList< label >
Foam::meshSurfacePartitioner
Definition: meshSurfacePartitioner.H:52
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::parMapperHelper::coordinates
const point & coordinates() const
return point coordinates
Definition: parMapperHelper.H:90
Foam::meshSurfacePartitioner::edgePoints
const labelHashSet & edgePoints() const
return labels of edge points (from the list of boundary points)
Definition: meshSurfacePartitioner.H:135
helperFunctions.H
Foam::Vector< scalar >
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::meshSurfaceMapper::selectNodesAtParallelBnd
void selectNodesAtParallelBnd(const labelLongList &)
check if nodes at parallel boundaries are selected at all processors
Definition: meshSurfaceMapperMapVertices.C:53
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::help::exchangeMap
void exchangeMap(const std::map< label, ListType > &m, LongList< T > &data, const Pstream::commsTypes commsType)
Definition: helperFunctionsPar.C:129
meshSurfaceEngineModifier.H
Foam::meshSurfaceEngine::boundaryPoints
const labelList & boundaryPoints() const
Definition: meshSurfaceEngineI.H:84
Foam::meshSurfaceMapper::mapVerticesOntoSurface
void mapVerticesOntoSurface()
Definition: meshSurfaceMapperMapVertices.C:202
Foam::meshSurfaceEngine::globalToLocalBndPointAddressing
const Map< label > & globalToLocalBndPointAddressing() const
global point label to local label. Only for processors points
Definition: meshSurfaceEngineI.H:431
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::meshOctree::findNearestSurfacePointInRegion
void findNearestSurfacePointInRegion(point &nearest, scalar &distSq, label &nearestTriangle, const label region, const point &p) const
find nearest surface point for vertex in a given region
Definition: meshOctreeFindNearestSurfacePoint.C:131
Foam::meshSurfacePartitioner::corners
const labelHashSet & corners() const
return labels of corner points (from the list of boundary points)
Definition: meshSurfacePartitioner.H:129
Foam::pointFieldPMG
Definition: pointFieldPMG.H:50
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::labelLongList
LongList< label > labelLongList
Definition: labelLongList.H:46
Foam::meshSurfaceMapper::meshOctree_
const meshOctree & meshOctree_
reference to the octree
Definition: meshSurfaceMapper.H:66
meshSurfacePartitioner.H