modifyMesh.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-2021 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 Application
25  modifyMesh
26 
27 Description
28  Manipulates mesh elements.
29 
30  Actions are:
31  (boundary)points:
32  - move
33 
34  (boundary)edges:
35  - split and move introduced point
36 
37  (boundary)faces:
38  - split(triangulate) and move introduced point
39 
40  edges:
41  - collapse
42 
43  cells:
44  - split into polygonal base pyramids around newly introduced mid
45  point
46 
47  Is a bit of a loose collection of mesh change drivers.
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #include "argList.H"
52 #include "Time.H"
53 #include "polyMesh.H"
54 #include "polyTopoChange.H"
55 #include "mapPolyMesh.H"
56 #include "boundaryCutter.H"
57 #include "cellSplitter.H"
58 #include "edgeCollapser.H"
59 #include "meshTools.H"
60 #include "Pair.H"
61 #include "globalIndex.H"
62 #include "systemDict.H"
63 
64 using namespace Foam;
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 // Locate point on patch. Returns (mesh) point label.
69 label findPoint(const primitivePatch& pp, const point& nearPoint)
70 {
71  const pointField& points = pp.points();
72  const labelList& meshPoints = pp.meshPoints();
73 
74  // Find nearest and next nearest
75  scalar minDistSqr = great;
76  label minI = -1;
77 
78  scalar almostMinDistSqr = great;
79  label almostMinI = -1;
80 
81  forAll(meshPoints, i)
82  {
83  label pointi = meshPoints[i];
84 
85  scalar distSqr = magSqr(nearPoint - points[pointi]);
86 
87  if (distSqr < minDistSqr)
88  {
89  almostMinDistSqr = minDistSqr;
90  almostMinI = minI;
91 
92  minDistSqr = distSqr;
93  minI = pointi;
94  }
95  else if (distSqr < almostMinDistSqr)
96  {
97  almostMinDistSqr = distSqr;
98  almostMinI = pointi;
99  }
100  }
101 
102 
103  // Decide if nearPoint unique enough.
104  Info<< "Found to point " << nearPoint << nl
105  << " nearest point : " << minI
106  << " distance " << Foam::sqrt(minDistSqr)
107  << " at " << points[minI] << nl
108  << " next nearest point : " << almostMinI
109  << " distance " << Foam::sqrt(almostMinDistSqr)
110  << " at " << points[almostMinI] << endl;
111 
112  if (almostMinDistSqr < 4*minDistSqr)
113  {
114  Info<< "Next nearest too close to nearest. Aborting" << endl;
115 
116  return -1;
117  }
118  else
119  {
120  return minI;
121  }
122 }
123 
124 
125 // Locate edge on patch. Return mesh edge label.
127 (
128  const primitiveMesh& mesh,
129  const primitivePatch& pp,
130  const point& nearPoint
131 )
132 {
133  const pointField& localPoints = pp.localPoints();
134  const pointField& points = pp.points();
135  const labelList& meshPoints = pp.meshPoints();
136  const edgeList& edges = pp.edges();
137 
138  // Find nearest and next nearest
139  scalar minDist = great;
140  label minI = -1;
141 
142  scalar almostMinDist = great;
143  label almostMinI = -1;
144 
145  forAll(edges, edgeI)
146  {
147  const edge& e = edges[edgeI];
148 
149  pointHit pHit(e.line(localPoints).nearestDist(nearPoint));
150 
151  if (pHit.hit())
152  {
153  if (pHit.distance() < minDist)
154  {
155  almostMinDist = minDist;
156  almostMinI = minI;
157 
158  minDist = pHit.distance();
159  minI = meshTools::findEdge
160  (
161  mesh,
162  meshPoints[e[0]],
163  meshPoints[e[1]]
164  );
165  }
166  else if (pHit.distance() < almostMinDist)
167  {
168  almostMinDist = pHit.distance();
169  almostMinI = meshTools::findEdge
170  (
171  mesh,
172  meshPoints[e[0]],
173  meshPoints[e[1]]
174  );
175  }
176  }
177  }
178 
179  if (minI == -1)
180  {
181  Info<< "Did not find edge close to point " << nearPoint << endl;
182 
183  return -1;
184  }
185 
186 
187  // Decide if nearPoint unique enough.
188  Info<< "Found to point " << nearPoint << nl
189  << " nearest edge : " << minI
190  << " distance " << minDist << " endpoints "
191  << mesh.edges()[minI].line(points) << nl
192  << " next nearest edge : " << almostMinI
193  << " distance " << almostMinDist << " endpoints "
194  << mesh.edges()[almostMinI].line(points) << nl
195  << endl;
196 
197  if (almostMinDist < 2*minDist)
198  {
199  Info<< "Next nearest too close to nearest. Aborting" << endl;
200 
201  return -1;
202  }
203  else
204  {
205  return minI;
206  }
207 }
208 
209 
210 // Find face on patch. Return mesh face label.
211 label findFace
212 (
213  const primitiveMesh& mesh,
214  const primitivePatch& pp,
215  const point& nearPoint
216 )
217 {
218  const pointField& points = pp.points();
219 
220  // Find nearest and next nearest
221  scalar minDist = great;
222  label minI = -1;
223 
224  scalar almostMinDist = great;
225  label almostMinI = -1;
226 
227  forAll(pp, patchFacei)
228  {
229  pointHit pHit(pp[patchFacei].nearestPoint(nearPoint, points));
230 
231  if (pHit.hit())
232  {
233  if (pHit.distance() < minDist)
234  {
235  almostMinDist = minDist;
236  almostMinI = minI;
237 
238  minDist = pHit.distance();
239  minI = patchFacei + mesh.nInternalFaces();
240  }
241  else if (pHit.distance() < almostMinDist)
242  {
243  almostMinDist = pHit.distance();
244  almostMinI = patchFacei + mesh.nInternalFaces();
245  }
246  }
247  }
248 
249  if (minI == -1)
250  {
251  Info<< "Did not find face close to point " << nearPoint << endl;
252 
253  return -1;
254  }
255 
256 
257  // Decide if nearPoint unique enough.
258  Info<< "Found to point " << nearPoint << nl
259  << " nearest face : " << minI
260  << " distance " << minDist
261  << " to face centre " << mesh.faceCentres()[minI] << nl
262  << " next nearest face : " << almostMinI
263  << " distance " << almostMinDist
264  << " to face centre " << mesh.faceCentres()[almostMinI] << nl
265  << endl;
266 
267  if (almostMinDist < 2*minDist)
268  {
269  Info<< "Next nearest too close to nearest. Aborting" << endl;
270 
271  return -1;
272  }
273  else
274  {
275  return minI;
276  }
277 }
278 
279 
280 // Find cell with cell centre close to given point.
281 label findCell(const primitiveMesh& mesh, const point& nearPoint)
282 {
283  label celli = mesh.findCell(nearPoint);
284 
285  if (celli != -1)
286  {
287  scalar distToCcSqr = magSqr(nearPoint - mesh.cellCentres()[celli]);
288 
289  const labelList& cPoints = mesh.cellPoints()[celli];
290 
291  label minI = -1;
292  scalar minDistSqr = great;
293 
294  forAll(cPoints, i)
295  {
296  label pointi = cPoints[i];
297 
298  scalar distSqr = magSqr(nearPoint - mesh.points()[pointi]);
299 
300  if (distSqr < minDistSqr)
301  {
302  minDistSqr = distSqr;
303  minI = pointi;
304  }
305  }
306 
307  // Decide if nearPoint unique enough.
308  Info<< "Found to point " << nearPoint << nl
309  << " nearest cell : " << celli
310  << " distance " << Foam::sqrt(distToCcSqr)
311  << " to cell centre " << mesh.cellCentres()[celli] << nl
312  << " nearest mesh point : " << minI
313  << " distance " << Foam::sqrt(minDistSqr)
314  << " to " << mesh.points()[minI] << nl
315  << endl;
316 
317  if (minDistSqr < 4*distToCcSqr)
318  {
319  Info<< "Mesh point too close to nearest cell centre. Aborting"
320  << endl;
321 
322  celli = -1;
323  }
324  }
325 
326  return celli;
327 }
328 
329 
330 int main(int argc, char *argv[])
331 {
332  #include "addOverwriteOption.H"
333 
334  #include "setRootCase.H"
335  #include "createTime.H"
337  #include "createPolyMesh.H"
338  const word oldInstance = mesh.pointsInstance();
339 
340  const bool overwrite = args.optionFound("overwrite");
341 
342  const dictionary dict(systemDict("modifyMeshDict", args, mesh));
343 
344  // Read all from the dictionary.
345  List<Pair<point>> pointsToMove(dict.lookup("pointsToMove"));
346  List<Pair<point>> edgesToSplit(dict.lookup("edgesToSplit"));
347  List<Pair<point>> facesToTriangulate
348  (
349  dict.lookup("facesToTriangulate")
350  );
351 
352  bool cutBoundary =
353  (
354  pointsToMove.size()
355  || edgesToSplit.size()
356  || facesToTriangulate.size()
357  );
358 
359  List<Pair<point>> edgesToCollapse(dict.lookup("edgesToCollapse"));
360 
361  bool collapseEdge = edgesToCollapse.size();
362 
363  List<Pair<point>> cellsToPyramidise(dict.lookup("cellsToSplit"));
364 
365  bool cellsToSplit = cellsToPyramidise.size();
366 
367  // List<Tuple2<pointField,point>>
368  // cellsToCreate(dict.lookup("cellsToCreate"));
369 
370  Info<< "Read from " << dict.name() << nl
371  << " Boundary cutting module:" << nl
372  << " points to move :" << pointsToMove.size() << nl
373  << " edges to split :" << edgesToSplit.size() << nl
374  << " faces to triangulate:" << facesToTriangulate.size() << nl
375  << " Cell splitting module:" << nl
376  << " cells to split :" << cellsToPyramidise.size() << nl
377  << " Edge collapsing module:" << nl
378  << " edges to collapse :" << edgesToCollapse.size() << nl
379  //<< " cells to create :" << cellsToCreate.size() << nl
380  << endl;
381 
382  if
383  (
384  (cutBoundary && collapseEdge)
385  || (cutBoundary && cellsToSplit)
386  || (collapseEdge && cellsToSplit)
387  )
388  {
390  << "Used more than one mesh modifying module "
391  << "(boundary cutting, cell splitting, edge collapsing)" << nl
392  << "Please do them in separate passes." << exit(FatalError);
393  }
394 
395 
396 
397  // Get calculating engine for all of outside
398  const SubList<face> outsideFaces
399  (
400  mesh.faces(),
401  mesh.nFaces() - mesh.nInternalFaces(),
402  mesh.nInternalFaces()
403  );
404 
405  primitivePatch allBoundary(outsideFaces, mesh.points());
406 
407 
408  // Look up mesh labels and convert to input for boundaryCutter.
409 
410  bool validInputs = true;
411 
412 
413  Info<< nl << "Looking up points to move ..." << nl << endl;
414  Map<point> pointToPos(pointsToMove.size());
415  forAll(pointsToMove, i)
416  {
417  const Pair<point>& pts = pointsToMove[i];
418 
419  label pointi = findPoint(allBoundary, pts.first());
420 
421  if (pointi == -1 || !pointToPos.insert(pointi, pts.second()))
422  {
423  Info<< "Could not insert mesh point " << pointi
424  << " for input point " << pts.first() << nl
425  << "Perhaps the point is already marked for moving?" << endl;
426  validInputs = false;
427  }
428  }
429 
430 
431  Info<< nl << "Looking up edges to split ..." << nl << endl;
432  Map<List<point>> edgeToCuts(edgesToSplit.size());
433  forAll(edgesToSplit, i)
434  {
435  const Pair<point>& pts = edgesToSplit[i];
436 
437  label edgeI = findEdge(mesh, allBoundary, pts.first());
438 
439  if
440  (
441  edgeI == -1
442  || !edgeToCuts.insert(edgeI, List<point>(1, pts.second()))
443  )
444  {
445  Info<< "Could not insert mesh edge " << edgeI
446  << " for input point " << pts.first() << nl
447  << "Perhaps the edge is already marked for cutting?" << endl;
448 
449  validInputs = false;
450  }
451  }
452 
453 
454  Info<< nl << "Looking up faces to triangulate ..." << nl << endl;
455  Map<point> faceToDecompose(facesToTriangulate.size());
456  forAll(facesToTriangulate, i)
457  {
458  const Pair<point>& pts = facesToTriangulate[i];
459 
460  label facei = findFace(mesh, allBoundary, pts.first());
461 
462  if (facei == -1 || !faceToDecompose.insert(facei, pts.second()))
463  {
464  Info<< "Could not insert mesh face " << facei
465  << " for input point " << pts.first() << nl
466  << "Perhaps the face is already marked for splitting?" << endl;
467 
468  validInputs = false;
469  }
470  }
471 
472 
473 
474  Info<< nl << "Looking up cells to convert to pyramids around"
475  << " cell centre ..." << nl << endl;
476  Map<point> cellToPyrCentre(cellsToPyramidise.size());
477  forAll(cellsToPyramidise, i)
478  {
479  const Pair<point>& pts = cellsToPyramidise[i];
480 
481  label celli = findCell(mesh, pts.first());
482 
483  if (celli == -1 || !cellToPyrCentre.insert(celli, pts.second()))
484  {
485  Info<< "Could not insert mesh cell " << celli
486  << " for input point " << pts.first() << nl
487  << "Perhaps the cell is already marked for splitting?" << endl;
488 
489  validInputs = false;
490  }
491  }
492 
493 
494  Info<< nl << "Looking up edges to collapse ..." << nl << endl;
495  Map<point> edgeToPos(edgesToCollapse.size());
496  forAll(edgesToCollapse, i)
497  {
498  const Pair<point>& pts = edgesToCollapse[i];
499 
500  label edgeI = findEdge(mesh, allBoundary, pts.first());
501 
502  if (edgeI == -1 || !edgeToPos.insert(edgeI, pts.second()))
503  {
504  Info<< "Could not insert mesh edge " << edgeI
505  << " for input point " << pts.first() << nl
506  << "Perhaps the edge is already marked for collapsing?" << endl;
507 
508  validInputs = false;
509  }
510  }
511 
512 
513 
514  if (!validInputs)
515  {
516  Info<< nl << "There was a problem in one of the inputs in the"
517  << " dictionary. Not modifying mesh." << endl;
518  }
519  else if (cellToPyrCentre.size())
520  {
521  Info<< nl << "All input cells located. Modifying mesh." << endl;
522 
523  // Mesh change engine
524  cellSplitter cutter(mesh);
525 
526  // Topo change container
527  polyTopoChange meshMod(mesh);
528 
529  // Insert commands into meshMod
530  cutter.setRefinement(cellToPyrCentre, meshMod);
531 
532  // Do changes
533  autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
534 
535  if (morphMap().hasMotionPoints())
536  {
537  mesh.movePoints(morphMap().preMotionPoints());
538  }
539 
540  cutter.updateMesh(morphMap());
541 
542  if (!overwrite)
543  {
544  runTime++;
545  }
546  else
547  {
548  mesh.setInstance(oldInstance);
549  }
550 
551  // Write resulting mesh
552  Info<< "Writing modified mesh to time " << runTime.timeName() << endl;
553  mesh.write();
554  }
555  else if (edgeToPos.size())
556  {
557  Info<< nl << "All input edges located. Modifying mesh." << endl;
558 
559  // Mesh change engine
560  edgeCollapser cutter(mesh);
561 
562  const edgeList& edges = mesh.edges();
563  const pointField& points = mesh.points();
564 
565  pointField newPoints(points);
566 
568  Map<point> collapsePointToLocation(mesh.nPoints());
569 
570  // Get new positions and construct collapse network
571  forAllConstIter(Map<point>, edgeToPos, iter)
572  {
573  label edgeI = iter.key();
574  const edge& e = edges[edgeI];
575 
576  collapseEdge[edgeI] = true;
577  collapsePointToLocation.set(e[1], points[e[0]]);
578 
579  newPoints[e[0]] = iter();
580  }
581 
582  // Move master point to destination.
583  mesh.movePoints(newPoints);
584 
585  List<pointEdgeCollapse> allPointInfo;
586  const globalIndex globalPoints(mesh.nPoints());
587  labelList pointPriority(mesh.nPoints(), 0);
588 
589  cutter.consistentCollapse
590  (
591  globalPoints,
592  pointPriority,
593  collapsePointToLocation,
594  collapseEdge,
595  allPointInfo
596  );
597 
598  // Topo change container
599  polyTopoChange meshMod(mesh);
600 
601  // Insert
602  cutter.setRefinement(allPointInfo, meshMod);
603 
604  // Do changes
605  autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
606 
607  if (morphMap().hasMotionPoints())
608  {
609  mesh.movePoints(morphMap().preMotionPoints());
610  }
611 
612  // Not implemented yet:
613  // cutter.updateMesh(morphMap());
614 
615 
616  if (!overwrite)
617  {
618  runTime++;
619  }
620  else
621  {
622  mesh.setInstance(oldInstance);
623  }
624 
625  // Write resulting mesh
626  Info<< "Writing modified mesh to time " << runTime.timeName() << endl;
627  mesh.write();
628  }
629  else
630  {
631  Info<< nl << "All input points located. Modifying mesh." << endl;
632 
633  // Mesh change engine
634  boundaryCutter cutter(mesh);
635 
636  // Topo change container
637  polyTopoChange meshMod(mesh);
638 
639  // Insert commands into meshMod
640  cutter.setRefinement
641  (
642  pointToPos,
643  edgeToCuts,
644  Map<labelPair>(0), // Faces to split diagonally
645  faceToDecompose, // Faces to triangulate
646  meshMod
647  );
648 
649  // Do changes
650  autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
651 
652  if (morphMap().hasMotionPoints())
653  {
654  mesh.movePoints(morphMap().preMotionPoints());
655  }
656 
657  cutter.updateMesh(morphMap());
658 
659  if (!overwrite)
660  {
661  runTime++;
662  }
663  else
664  {
665  mesh.setInstance(oldInstance);
666  }
667 
668  // Write resulting mesh
669  Info<< "Writing modified mesh to time " << runTime.timeName() << endl;
670  mesh.write();
671  }
672 
673 
674  Info<< "\nEnd\n" << endl;
675  return 0;
676 }
677 
678 
679 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
label findCell(const point &location) const
Find cell enclosing this location (-1 if not in mesh)
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void off()
Switch the function objects off.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
label nInternalFaces() const
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
label nFaces() const
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
Does polyTopoChanges to remove edges. Can remove faces due to edge collapse but can not remove cells ...
Definition: edgeCollapser.H:66
virtual const pointField & points() const =0
Return mesh points.
label collapseEdge(triSurface &surf, const scalar minLen)
Keep collapsing all edges < minLen.
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:337
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const Field< PointType > & localPoints() const
Return pointField of points in patch.
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
scalar minDist(const List< pointIndexHit > &hitList)
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
Does pyramidal decomposition of selected cells. So all faces will become base of pyramid with as top ...
Definition: cellSplitter.H:57
A list of faces which address into the list of points.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
A List obtained as a section of another List.
Definition: SubList.H:53
linePointRef line(const pointField &) const
Return edge line.
Definition: edgeI.H:187
tmp< scalarField > movePoints(const pointField &p, const pointField &oldP)
Move points, returns volumes swept by faces in motion.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
A class for handling words, derived from string.
Definition: word.H:59
const Field< PointType > & points() const
Return reference to global points.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const vectorField & cellCentres() const
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:260
const Type & second() const
Return second.
Definition: Pair.H:110
label nEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const vectorField & faceCentres() const
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:100
A bit-packed bool list.
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:410
Direct mesh changes based on v1.3 polyTopoChange syntax.
virtual const faceList & faces() const =0
Return faces.
messageStream Info
label nPoints() const
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Foam::argList args(argc, argv)
const Type & first() const
Return first.
Definition: Pair.H:98
Does modifications to boundary faces.
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition: lineI.H:95
const labelListList & cellPoints() const
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49