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