meshSurfaceMapperPremapVertices.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 "meshOctree.H"
32 #include "triSurf.H"
33 #include "refLabelledPoint.H"
34 #include "refLabelledPointScalar.H"
35 #include "helperFunctions.H"
36 #include "meshSurfaceOptimizer.H"
37 
38 # ifdef USE_OMP
39 #include <omp.h>
40 # endif
41 
42 //#define DEBUGMapping
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
52 {
53  Info << "Smoothing mesh surface before mapping." << endl;
54 
55  const labelList& boundaryPoints = surfaceEngine_.boundaryPoints();
57  const vectorField& faceCentres = surfaceEngine_.faceCentres();
58  const VRWGraph& pointFaces = surfaceEngine_.pointFaces();
59  const VRWGraph& pointInFace = surfaceEngine_.pointInFaces();
60  const labelList& bp = surfaceEngine_.bp();
62 
63  const triSurf& surf = meshOctree_.surface();
64 
65  List<labelledPointScalar> preMapPositions(boundaryPoints.size());
66  List<DynList<scalar, 6> > faceCentreDistances(bFaces.size());
67 
68  # ifdef USE_OMP
69  # pragma omp parallel for schedule(dynamic, 20)
70  # endif
71  forAll(bFaces, bfI)
72  {
73  const point& c = faceCentres[bfI];
74  const face& bf = bFaces[bfI];
75 
76  faceCentreDistances[bfI].setSize(bf.size());
77 
78  forAll(bf, pI)
79  {
80  faceCentreDistances[bfI][pI] = magSqr(points[bf[pI]] - c);
81  }
82  }
83 
84  for(label iterI=0;iterI<nIterations;++iterI)
85  {
86  //- find patches in the vicinity of a boundary face
87  List<DynList<label> > boundaryPointPatches(boundaryPoints.size());
88  # ifdef USE_OMP
89  # pragma omp parallel for schedule(dynamic, 50)
90  # endif
91  forAll(bFaces, bfI)
92  {
93  const face& bf = bFaces[bfI];
94  scalar boxSize(0.0);
95  forAll(bf, pI)
96  {
97  boxSize =
98  Foam::max
99  (
100  boxSize,
101  mag(faceCentres[bfI] - points[bf[pI]])
102  );
103  }
104 
105  const boundBox bb
106  (
107  faceCentres[bfI] - vector(boxSize, boxSize, boxSize),
108  faceCentres[bfI] + vector(boxSize, boxSize, boxSize)
109  );
110 
111  DynList<label> containedLeaves;
112  meshOctree_.findLeavesContainedInBox(bb, containedLeaves);
113 
115  forAll(containedLeaves, clI)
116  {
117  DynList<label> ct;
118  meshOctree_.containedTriangles(containedLeaves[clI], ct);
119 
120  forAll(ct, i)
121  patches.appendIfNotIn(surf[ct[i]].region());
122  }
123 
124  scalar metric(VGREAT);
125  label bestPatch(-1);
126  forAll(patches, ptchI)
127  {
128  const scalar m = faceMetricInPatch(bfI, patches[ptchI]);
129 
130  if( m < metric )
131  {
132  metric = m;
133  bestPatch = patches[ptchI];
134  }
135  }
136 
137  forAll(bf, pI)
138  boundaryPointPatches[bp[bf[pI]]].appendIfNotIn(bestPatch);
139  }
140 
141  //- use the shrinking laplace first
142  # ifdef USE_OMP
143  # pragma omp parallel for schedule(dynamic, 40)
144  # endif
145  forAll(pointFaces, bpI)
146  {
147  labelledPointScalar lp(bpI, vector::zero, 0.0);
148 
149  const point& p = points[boundaryPoints[bpI]];
150 
151  forAllRow(pointFaces, bpI, pfI)
152  {
153  const label bfI = pointFaces(bpI, pfI);
154  const point& fc = faceCentres[pointFaces(bpI, pfI)];
155  const label pos = pointInFace(bpI, pfI);
156  const scalar w
157  (
158  max(magSqr(p - fc) / faceCentreDistances[bfI][pos], SMALL)
159  );
160  lp.coordinates() += w * faceCentres[bfI];
161  lp.scalarValue() += w;
162  }
163 
164  preMapPositions[bpI] = lp;
165  }
166 
167  //- pointer needed in case of parallel calculation
168  const VRWGraph* bpAtProcsPtr(NULL);
169 
170  if( Pstream::parRun() )
171  {
172  const VRWGraph& bpAtProcs = surfaceEngine_.bpAtProcs();
173  bpAtProcsPtr = &bpAtProcs;
174  const labelList& globalPointLabel =
176  const Map<label>& globalToLocal =
178 
179  //- collect data to be sent to other processors
180  std::map<label, LongList<labelledPointScalar> > exchangeData;
182  exchangeData.insert
183  (
184  std::make_pair
185  (
188  )
189  );
190 
191  forAllConstIter(Map<label>, globalToLocal, it)
192  {
193  const label bpI = it();
194 
195  forAllRow(bpAtProcs, bpI, procI)
196  {
197  const label neiProc = bpAtProcs(bpI, procI);
198 
199  if( neiProc == Pstream::myProcNo() )
200  continue;
201 
202  exchangeData[neiProc].append
203  (
205  (
206  globalPointLabel[bpI],
207  preMapPositions[bpI].coordinates(),
208  preMapPositions[bpI].scalarValue()
209  )
210  );
211  }
212  }
213 
214  //- exchange data with other processors
215  LongList<labelledPointScalar> receivedData;
216  help::exchangeMap(exchangeData, receivedData);
217 
218  //- combine collected data with the available data
219  forAll(receivedData, i)
220  {
221  const labelledPointScalar& lps = receivedData[i];
222 
223  const label bpI = globalToLocal[lps.pointLabel()];
224 
225  labelledPointScalar& lp = preMapPositions[bpI];
226  lp.coordinates() += lps.coordinates();
227  lp.scalarValue() += lps.scalarValue();
228  }
229  }
230 
231  //- create the surface modifier and move the surface points
232  meshSurfaceEngineModifier surfaceModifier(surfaceEngine_);
233  LongList<parMapperHelper> parallelBndNodes;
234 
235  # ifdef USE_OMP
236  # pragma omp parallel for schedule(dynamic, 50)
237  # endif
238  forAll(boundaryPoints, bpI)
239  {
240  labelledPointScalar& lps = preMapPositions[bpI];
241 
242  lps.coordinates() /= lps.scalarValue();
243 
244  const point& p = points[boundaryPoints[bpI]];
245 
246  //label patch, nearestTri;
247  point pMap = p;
248  scalar dSq;
249 
250  if( boundaryPointPatches[bpI].size() == 1 )
251  {
252  label nt;
254  (
255  pMap,
256  dSq,
257  nt,
258  boundaryPointPatches[bpI][0],
259  lps.coordinates()
260  );
261  }
262  else
263  {
265  (
266  pMap,
267  dSq,
268  lps.coordinates(),
269  boundaryPointPatches[bpI]
270  );
271  }
272 
273  const point newP = p + 0.5 * (pMap - p);
274 
275  surfaceModifier.moveBoundaryVertexNoUpdate(bpI, newP);
276 
277  if( bpAtProcsPtr && bpAtProcsPtr->sizeOfRow(bpI) )
278  {
279  # ifdef USE_OMP
280  # pragma omp critical
281  # endif
282  parallelBndNodes.append
283  (
285  (
286  newP,
287  dSq,
288  bpI,
289  -1
290  )
291  );
292  }
293  }
294 
295  //- make sure that the vertices at inter-processor boundaries
296  //- are mapped onto the same location
297  mapToSmallestDistance(parallelBndNodes);
298 
299  //- update the surface geometry of the
300  surfaceModifier.updateGeometry();
301 
303 
304  surfaceModifier.updateGeometry();
305  }
306 
307  Info << "Finished smoothing mesh surface before mapping." << endl;
308 }
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
312 } // End namespace Foam
313 
314 // ************************************************************************* //
Foam::meshSurfaceEngine::faceCentres
const vectorField & faceCentres() const
Definition: meshSurfaceEngineI.H:277
Foam::meshSurfaceEngine::bpAtProcs
const VRWGraph & bpAtProcs() const
processors which contain the vertex
Definition: meshSurfaceEngineI.H:451
Foam::meshSurfaceEngine::pointInFaces
const VRWGraph & pointInFaces() const
Definition: meshSurfaceEngineI.H:181
Foam::meshSurfaceEngine::bp
const labelList & bp() const
Definition: meshSurfaceEngineI.H:64
Foam::LongList::append
void append(const T &e)
Append an element at the end of the list.
Definition: LongListI.H:265
Foam::labelledPointScalar::scalarValue
const scalar & scalarValue() const
return scalar value
Definition: labelledPointScalar.H:109
Foam::meshSurfaceMapper::surfaceEngine_
const meshSurfaceEngine & surfaceEngine_
mesh surface
Definition: meshSurfaceMapper.H:63
triSurf.H
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::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::parMapperHelper
Definition: parMapperHelper.H:51
p
p
Definition: pEqn.H:62
meshSurfaceMapper.H
Foam::meshSurfaceMapper::faceMetricInPatch
scalar faceMetricInPatch(const label bfI, const label patchI) const
calculate face metric
Definition: meshSurfaceMapperCornersAndEdges.C:142
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::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::meshOctree::containedTriangles
void containedTriangles(const label, DynList< label > &) const
Definition: meshOctreeI.H:81
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::labelledPointScalar::coordinates
const point & coordinates() const
return point coordinates
Definition: labelledPointScalar.H:98
Foam::meshSurfaceEngine::bpNeiProcs
const DynList< label > & bpNeiProcs() const
communication matrix for sending point data
Definition: meshSurfaceEngineI.H:470
Foam::meshSurfaceEngine::pointFaces
const VRWGraph & pointFaces() const
Definition: meshSurfaceEngineI.H:162
Foam::Map< label >
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::meshOctree::findLeavesContainedInBox
void findLeavesContainedInBox(const boundBox &, DynList< const meshOctreeCube *, 256 > &) const
find leaves contained in the box
Definition: meshOctreeInsideCalculations.C:131
meshOctree.H
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::LongList
Definition: LongList.H:55
Foam::meshSurfaceEngineModifier
Definition: meshSurfaceEngineModifier.H:48
Foam::meshSurfaceOptimizer::untangleSurface
bool untangleSurface(const labelLongList &activeBoundaryPoints, const label nAdditionalLayers=0)
runs a surface smoother on the selected boundary points
Definition: meshSurfaceOptimizerOptimizeSurface.C:327
Foam::meshSurfaceEngine::boundaryFaces
const faceList::subList & boundaryFaces() const
Definition: meshSurfaceEngineI.H:103
Foam::meshSurfaceMapper::mapToSmallestDistance
void mapToSmallestDistance(LongList< parMapperHelper > &)
Definition: meshSurfaceMapperMapVertices.C:104
refLabelledPoint.H
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::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
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::meshOctree::findNearestPointToPatches
bool findNearestPointToPatches(point &nearest, scalar &distSq, const point &p, const DynList< label > &patches, const scalar tol=1e-4) const
find the nearest vertex to the given patches
Definition: meshOctreeFindNearestSurfacePoint.C:535
forAllRow
#define forAllRow(graph, rowI, index)
Definition: VRWGraph.H:277
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::labelledPointScalar
Definition: labelledPointScalar.H:50
Foam::meshSurfaceMapper::preMapVertices
void preMapVertices(const label nIterations=3)
Definition: meshSurfaceMapperPremapVertices.C:51
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::meshSurfaceEngineModifier::updateGeometry
void updateGeometry(const labelLongList &)
Definition: meshSurfaceEngineModifier.C:234
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::DynList< label >
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
helperFunctions.H
Foam::Vector< scalar >
Foam::meshSurfaceOptimizer
Definition: meshSurfaceOptimizer.H:63
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::labelledPointScalar::pointLabel
label pointLabel() const
return point label
Definition: labelledPointScalar.H:87
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::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::meshSurfaceEngine::boundaryPoints
const labelList & boundaryPoints() const
Definition: meshSurfaceEngineI.H:84
patches
patches[0]
Definition: createSingleCellMesh.H:36
refLabelledPointScalar.H
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
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::meshOctree::surface
const triSurf & surface() const
return a reference to the surface
Definition: meshOctreeI.H:130
Foam::UList::size
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
Foam::pointFieldPMG
Definition: pointFieldPMG.H:50
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::triSurf
Definition: triSurf.H:59
meshSurfaceOptimizer.H
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::meshSurfaceMapper::meshOctree_
const meshOctree & meshOctree_
reference to the octree
Definition: meshSurfaceMapper.H:66
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190