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-2018 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 "IOdictionary.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 
331 
332 int main(int argc, char *argv[])
333 {
334  #include "addOverwriteOption.H"
335 
336  #include "setRootCase.H"
337  #include "createTime.H"
339  #include "createPolyMesh.H"
340  const word oldInstance = mesh.pointsInstance();
341 
342  const bool overwrite = args.optionFound("overwrite");
343 
344  Info<< "Reading modifyMeshDict\n" << endl;
345 
347  (
348  IOobject
349  (
350  "modifyMeshDict",
351  runTime.system(),
352  mesh,
355  )
356  );
357 
358  // Read all from the dictionary.
359  List<Pair<point>> pointsToMove(dict.lookup("pointsToMove"));
360  List<Pair<point>> edgesToSplit(dict.lookup("edgesToSplit"));
361  List<Pair<point>> facesToTriangulate
362  (
363  dict.lookup("facesToTriangulate")
364  );
365 
366  bool cutBoundary =
367  (
368  pointsToMove.size()
369  || edgesToSplit.size()
370  || facesToTriangulate.size()
371  );
372 
373  List<Pair<point>> edgesToCollapse(dict.lookup("edgesToCollapse"));
374 
375  bool collapseEdge = edgesToCollapse.size();
376 
377  List<Pair<point>> cellsToPyramidise(dict.lookup("cellsToSplit"));
378 
379  bool cellsToSplit = cellsToPyramidise.size();
380 
381  // List<Tuple2<pointField,point>>
382  // cellsToCreate(dict.lookup("cellsToCreate"));
383 
384  Info<< "Read from " << dict.name() << nl
385  << " Boundary cutting module:" << nl
386  << " points to move :" << pointsToMove.size() << nl
387  << " edges to split :" << edgesToSplit.size() << nl
388  << " faces to triangulate:" << facesToTriangulate.size() << nl
389  << " Cell splitting module:" << nl
390  << " cells to split :" << cellsToPyramidise.size() << nl
391  << " Edge collapsing module:" << nl
392  << " edges to collapse :" << edgesToCollapse.size() << nl
393  //<< " cells to create :" << cellsToCreate.size() << nl
394  << endl;
395 
396  if
397  (
398  (cutBoundary && collapseEdge)
399  || (cutBoundary && cellsToSplit)
400  || (collapseEdge && cellsToSplit)
401  )
402  {
404  << "Used more than one mesh modifying module "
405  << "(boundary cutting, cell splitting, edge collapsing)" << nl
406  << "Please do them in separate passes." << exit(FatalError);
407  }
408 
409 
410 
411  // Get calculating engine for all of outside
412  const SubList<face> outsideFaces
413  (
414  mesh.faces(),
415  mesh.nFaces() - mesh.nInternalFaces(),
416  mesh.nInternalFaces()
417  );
418 
419  primitivePatch allBoundary(outsideFaces, mesh.points());
420 
421 
422  // Look up mesh labels and convert to input for boundaryCutter.
423 
424  bool validInputs = true;
425 
426 
427  Info<< nl << "Looking up points to move ..." << nl << endl;
428  Map<point> pointToPos(pointsToMove.size());
429  forAll(pointsToMove, i)
430  {
431  const Pair<point>& pts = pointsToMove[i];
432 
433  label pointi = findPoint(allBoundary, pts.first());
434 
435  if (pointi == -1 || !pointToPos.insert(pointi, pts.second()))
436  {
437  Info<< "Could not insert mesh point " << pointi
438  << " for input point " << pts.first() << nl
439  << "Perhaps the point is already marked for moving?" << endl;
440  validInputs = false;
441  }
442  }
443 
444 
445  Info<< nl << "Looking up edges to split ..." << nl << endl;
446  Map<List<point>> edgeToCuts(edgesToSplit.size());
447  forAll(edgesToSplit, i)
448  {
449  const Pair<point>& pts = edgesToSplit[i];
450 
451  label edgeI = findEdge(mesh, allBoundary, pts.first());
452 
453  if
454  (
455  edgeI == -1
456  || !edgeToCuts.insert(edgeI, List<point>(1, pts.second()))
457  )
458  {
459  Info<< "Could not insert mesh edge " << edgeI
460  << " for input point " << pts.first() << nl
461  << "Perhaps the edge is already marked for cutting?" << endl;
462 
463  validInputs = false;
464  }
465  }
466 
467 
468  Info<< nl << "Looking up faces to triangulate ..." << nl << endl;
469  Map<point> faceToDecompose(facesToTriangulate.size());
470  forAll(facesToTriangulate, i)
471  {
472  const Pair<point>& pts = facesToTriangulate[i];
473 
474  label facei = findFace(mesh, allBoundary, pts.first());
475 
476  if (facei == -1 || !faceToDecompose.insert(facei, pts.second()))
477  {
478  Info<< "Could not insert mesh face " << facei
479  << " for input point " << pts.first() << nl
480  << "Perhaps the face is already marked for splitting?" << endl;
481 
482  validInputs = false;
483  }
484  }
485 
486 
487 
488  Info<< nl << "Looking up cells to convert to pyramids around"
489  << " cell centre ..." << nl << endl;
490  Map<point> cellToPyrCentre(cellsToPyramidise.size());
491  forAll(cellsToPyramidise, i)
492  {
493  const Pair<point>& pts = cellsToPyramidise[i];
494 
495  label celli = findCell(mesh, pts.first());
496 
497  if (celli == -1 || !cellToPyrCentre.insert(celli, pts.second()))
498  {
499  Info<< "Could not insert mesh cell " << celli
500  << " for input point " << pts.first() << nl
501  << "Perhaps the cell is already marked for splitting?" << endl;
502 
503  validInputs = false;
504  }
505  }
506 
507 
508  Info<< nl << "Looking up edges to collapse ..." << nl << endl;
509  Map<point> edgeToPos(edgesToCollapse.size());
510  forAll(edgesToCollapse, i)
511  {
512  const Pair<point>& pts = edgesToCollapse[i];
513 
514  label edgeI = findEdge(mesh, allBoundary, pts.first());
515 
516  if (edgeI == -1 || !edgeToPos.insert(edgeI, pts.second()))
517  {
518  Info<< "Could not insert mesh edge " << edgeI
519  << " for input point " << pts.first() << nl
520  << "Perhaps the edge is already marked for collaping?" << endl;
521 
522  validInputs = false;
523  }
524  }
525 
526 
527 
528  if (!validInputs)
529  {
530  Info<< nl << "There was a problem in one of the inputs in the"
531  << " dictionary. Not modifying mesh." << endl;
532  }
533  else if (cellToPyrCentre.size())
534  {
535  Info<< nl << "All input cells located. Modifying mesh." << endl;
536 
537  // Mesh change engine
538  cellSplitter cutter(mesh);
539 
540  // Topo change container
541  polyTopoChange meshMod(mesh);
542 
543  // Insert commands into meshMod
544  cutter.setRefinement(cellToPyrCentre, meshMod);
545 
546  // Do changes
547  autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
548 
549  if (morphMap().hasMotionPoints())
550  {
551  mesh.movePoints(morphMap().preMotionPoints());
552  }
553 
554  cutter.updateMesh(morphMap());
555 
556  if (!overwrite)
557  {
558  runTime++;
559  }
560  else
561  {
562  mesh.setInstance(oldInstance);
563  }
564 
565  // Write resulting mesh
566  Info<< "Writing modified mesh to time " << runTime.timeName() << endl;
567  mesh.write();
568  }
569  else if (edgeToPos.size())
570  {
571  Info<< nl << "All input edges located. Modifying mesh." << endl;
572 
573  // Mesh change engine
574  edgeCollapser cutter(mesh);
575 
576  const edgeList& edges = mesh.edges();
577  const pointField& points = mesh.points();
578 
579  pointField newPoints(points);
580 
582  Map<point> collapsePointToLocation(mesh.nPoints());
583 
584  // Get new positions and construct collapse network
585  forAllConstIter(Map<point>, edgeToPos, iter)
586  {
587  label edgeI = iter.key();
588  const edge& e = edges[edgeI];
589 
590  collapseEdge[edgeI] = true;
591  collapsePointToLocation.set(e[1], points[e[0]]);
592 
593  newPoints[e[0]] = iter();
594  }
595 
596  // Move master point to destination.
597  mesh.movePoints(newPoints);
598 
599  List<pointEdgeCollapse> allPointInfo;
600  const globalIndex globalPoints(mesh.nPoints());
601  labelList pointPriority(mesh.nPoints(), 0);
602 
603  cutter.consistentCollapse
604  (
605  globalPoints,
606  pointPriority,
607  collapsePointToLocation,
608  collapseEdge,
609  allPointInfo
610  );
611 
612  // Topo change container
613  polyTopoChange meshMod(mesh);
614 
615  // Insert
616  cutter.setRefinement(allPointInfo, meshMod);
617 
618  // Do changes
619  autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
620 
621  if (morphMap().hasMotionPoints())
622  {
623  mesh.movePoints(morphMap().preMotionPoints());
624  }
625 
626  // Not implemented yet:
627  // cutter.updateMesh(morphMap());
628 
629 
630  if (!overwrite)
631  {
632  runTime++;
633  }
634  else
635  {
636  mesh.setInstance(oldInstance);
637  }
638 
639  // Write resulting mesh
640  Info<< "Writing modified mesh to time " << runTime.timeName() << endl;
641  mesh.write();
642  }
643  else
644  {
645  Info<< nl << "All input points located. Modifying mesh." << endl;
646 
647  // Mesh change engine
648  boundaryCutter cutter(mesh);
649 
650  // Topo change container
651  polyTopoChange meshMod(mesh);
652 
653  // Insert commands into meshMod
654  cutter.setRefinement
655  (
656  pointToPos,
657  edgeToCuts,
658  Map<labelPair>(0), // Faces to split diagonally
659  faceToDecompose, // Faces to triangulate
660  meshMod
661  );
662 
663  // Do changes
664  autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
665 
666  if (morphMap().hasMotionPoints())
667  {
668  mesh.movePoints(morphMap().preMotionPoints());
669  }
670 
671  cutter.updateMesh(morphMap());
672 
673  if (!overwrite)
674  {
675  runTime++;
676  }
677  else
678  {
679  mesh.setInstance(oldInstance);
680  }
681 
682  // Write resulting mesh
683  Info<< "Writing modified mesh to time " << runTime.timeName() << endl;
684  mesh.write();
685  }
686 
687 
688  Info<< "\nEnd\n" << endl;
689  return 0;
690 }
691 
692 
693 // ************************************************************************* //
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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:163
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
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:626
const Field< PointType > & localPoints() const
Return pointField of points in patch.
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
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
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
dynamicFvMesh & mesh
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:111
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
const word & system() const
Return system name.
Definition: TimePaths.H:114
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:265
const Type & second() const
Return second.
Definition: Pair.H:99
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)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const Type & first() const
Return first.
Definition: Pair.H:87
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:583
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49