Fixed URL for libccmio-2.6.1 (bug report #5 by Thomas Oliveira)
[foam-extend-3.2.git] / src / mesh / cfMesh / utilities / surfaceTools / meshSurfaceMapper / meshSurfaceMapperMapVertices.C
blobc54e430a7965615d5e5a2c98c3990c1cd8c65a4c
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.
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.
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.
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/>.
24 Description
26 \*---------------------------------------------------------------------------*/
28 #include "demandDrivenData.H"
29 #include "meshSurfaceEngineModifier.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"
37 #include <map>
39 # ifdef USE_OMP
40 #include <omp.h>
41 # endif
43 //#define DEBUGMapping
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 namespace Foam
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 // Private member functions
53 void meshSurfaceMapper::selectNodesAtParallelBnd(const labelLongList& selNodes)
55     if( !Pstream::parRun() )
56         return;
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()));
63     const VRWGraph& bpAtProcs = surfaceEngine_.bpAtProcs();
64     const labelList& globalPointLabel =
65         surfaceEngine_.globalBoundaryPointLabel();
66     const Map<label>& globalToLocal =
67         surfaceEngine_.globalToLocalBndPointAddressing();
69     boolList selectedNode(bpAtProcs.size(), false);
71     forAll(selNodes, i)
72     {
73         const label bpI = selNodes[i];
75         selectedNode[bpI] = true;
77         forAllRow(bpAtProcs, bpI, procI)
78         {
79             const label neiProc = bpAtProcs(bpI, procI);
80             if( neiProc == Pstream::myProcNo() )
81                 continue;
83             exchangeData[neiProc].append(globalPointLabel[bpI]);
84         }
85     }
87     //- exchange data
88     labelLongList receivedData;
89     help::exchangeMap(exchangeData, receivedData);
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     }
104 void meshSurfaceMapper::mapToSmallestDistance(LongList<parMapperHelper>& parN)
106     if( !Pstream::parRun() )
107         return;
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         );
117     const VRWGraph& bpAtProcs = surfaceEngine_.bpAtProcs();
118     const labelList& globalPointLabel =
119         surfaceEngine_.globalBoundaryPointLabel();
120     const Map<label>& globalToLocal =
121         surfaceEngine_.globalToLocalBndPointAddressing();
123     Map<label> bpToList(parN.size());
125     forAll(parN, i)
126     {
127         const label bpI = parN[i].globalLabel();
128         bpToList.insert(bpI, i);
130         forAllRow(bpAtProcs, bpI, procI)
131         {
132             const label neiProc = bpAtProcs(bpI, procI);
133             if( neiProc == Pstream::myProcNo() )
134                 continue;
136             exchangeData[neiProc].append
137             (
138                 parMapperHelper
139                 (
140                     parN[i].coordinates(),
141                     parN[i].movingDistance(),
142                     globalPointLabel[bpI],
143                     parN[i].pointPatch()
144                 )
145             );
146         }
147     }
149     //- exchange data
150     LongList<parMapperHelper> receivedData;
151     help::exchangeMap(exchangeData, receivedData);
153     //- select the point with the smallest moving distance
154     meshSurfaceEngineModifier surfModifier(surfaceEngine_);
155     forAll(receivedData, i)
156     {
157         const parMapperHelper& ph = receivedData[i];
159         const label bpI = globalToLocal[ph.globalLabel()];
161         parMapperHelper& phOrig = parN[bpToList[bpI]];
162         if( phOrig.movingDistance() < ph.movingDistance() )
163         {
164             surfModifier.moveBoundaryVertexNoUpdate(bpI, ph.coordinates());
165             phOrig = ph;
166         }
167     }
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 void meshSurfaceMapper::mapNodeToPatch(const label bpI, const label patchI)
174     label patch, nt;
175     point mapPoint;
176     scalar dSq;
178     const pointFieldPMG& points = surfaceEngine_.points();
179     const labelList& bPoints = surfaceEngine_.boundaryPoints();
180     const point p = points[bPoints[bpI]];
182     if( patchI < 0 )
183     {
184         meshOctree_.findNearestSurfacePoint(mapPoint, dSq, nt, patch, p);
185     }
186     else
187     {
188         meshOctree_.findNearestSurfacePointInRegion
189         (
190             mapPoint,
191             dSq,
192             nt,
193             patchI,
194             p
195         );
196     }
198     meshSurfaceEngineModifier surfModifier(surfaceEngine_);
199     surfModifier.moveBoundaryVertex(bpI, mapPoint);
202 void meshSurfaceMapper::mapVerticesOntoSurface()
204     Info << "Mapping vertices onto surface" << endl;
206     labelLongList nodesToMap(surfaceEngine_.boundaryPoints().size());
207     forAll(nodesToMap, i)
208         nodesToMap[i] = i;
210     mapVerticesOntoSurface(nodesToMap);
212     Info << "Finished mapping vertices onto surface" << endl;
215 void meshSurfaceMapper::mapVerticesOntoSurface(const labelLongList& nodesToMap)
217     const labelList& boundaryPoints = surfaceEngine_.boundaryPoints();
218     const pointFieldPMG& points = surfaceEngine_.points();
220     const VRWGraph* bpAtProcsPtr(NULL);
221     if( Pstream::parRun() )
222         bpAtProcsPtr = &surfaceEngine_.bpAtProcs();
224     meshSurfaceEngineModifier surfaceModifier(surfaceEngine_);
225     LongList<parMapperHelper> parallelBndNodes;
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];
236         # ifdef DEBUGMapping
237         Info << nl << "Mapping vertex " << bpI << " with coordinates "
238             << points[boundaryPoints[bpI]] << endl;
239         # endif
241         label patch, nt;
242         point mapPoint;
243         scalar dSq;
245         meshOctree_.findNearestSurfacePoint
246         (
247             mapPoint,
248             dSq,
249             nt,
250             patch,
251             points[boundaryPoints[bpI]]
252         );
254         surfaceModifier.moveBoundaryVertexNoUpdate(bpI, mapPoint);
256         if( bpAtProcsPtr && bpAtProcsPtr->sizeOfRow(bpI) )
257         {
258             # ifdef USE_OMP
259             # pragma omp critical
260             # endif
261             parallelBndNodes.append
262             (
263                 parMapperHelper
264                 (
265                     mapPoint,
266                     dSq,
267                     bpI,
268                     patch
269                 )
270             );
271         }
273         # ifdef DEBUGMapping
274         Info << "Mapped point " << points[boundaryPoints[bpI]] << endl;
275         # endif
276     }
278     //- make sure that the points are at the nearest location on the surface
279     mapToSmallestDistance(parallelBndNodes);
281     //- re-calculate face normals, point normals, etc.
282     surfaceModifier.updateGeometry(nodesToMap);
285 void meshSurfaceMapper::mapVerticesOntoSurfacePatches()
287     Info << "Mapping vertices with respect to surface patches" << endl;
289     labelLongList nodesToMap(surfaceEngine_.boundaryPoints().size());
290     forAll(nodesToMap, i)
291         nodesToMap[i] = i;
293     mapVerticesOntoSurfacePatches(nodesToMap);
295     Info << "Finished mapping vertices with respect to surface patches" << endl;
298 void meshSurfaceMapper::mapVerticesOntoSurfacePatches
300     const labelLongList& nodesToMap
303     const meshSurfacePartitioner& mPart = meshPartitioner();
304     const labelHashSet& cornerPoints = mPart.corners();
305     const labelHashSet& edgePoints = mPart.edgePoints();
306     const VRWGraph& pointPatches = mPart.pointPatches();
308     boolList treatedPoint(surfaceEngine_.boundaryPoints().size(), false);
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     }
326     //- map the remaining selected points
327     const labelList& bPoints = surfaceEngine_.boundaryPoints();
328     const pointFieldPMG& points = surfaceEngine_.points();
330     const VRWGraph* bpAtProcsPtr(NULL);
331     if( Pstream::parRun() )
332         bpAtProcsPtr = &surfaceEngine_.bpAtProcs();
334     meshSurfaceEngineModifier surfaceModifier(surfaceEngine_);
335     LongList<parMapperHelper> parallelBndNodes;
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];
346         if( treatedPoint[bpI] )
347             continue;
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         );
362         surfaceModifier.moveBoundaryVertexNoUpdate(bpI, mapPoint);
364         if( bpAtProcsPtr && bpAtProcsPtr->sizeOfRow(bpI) )
365         {
366             # ifdef USE_OMP
367             # pragma omp critical
368             # endif
369             {
370                 parallelBndNodes.append
371                 (
372                     parMapperHelper
373                     (
374                         mapPoint,
375                         dSq,
376                         bpI,
377                         -1
378                     )
379                 );
380             }
381         }
383         # ifdef DEBUGMapping
384         Info << "Mapped point " << points[boundaryPoints[bpI]] << endl;
385         # endif
386     }
388     //- map vertices at inter-processor boundaries to the nearest location
389     //- on the surface
390     mapToSmallestDistance(parallelBndNodes);
392     //- update face normals, point normals, etc.
393     surfaceModifier.updateGeometry(nodesToMap);
395     //- map edge nodes
396     mapEdgeNodes(selectedEdges);
398     //- map corner vertices
399     mapCorners(selectedCorners);
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
404 } // End namespace Foam
406 // ************************************************************************* //