splitCells.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-2024 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  splitCells
26 
27 Description
28  Utility to split cells with flat faces.
29 
30  Uses a geometric cut with a plane dividing the edge angle into two so
31  might produce funny cells. For hexes it will use by default a cut from
32  edge onto opposite edge (i.e. purely topological).
33 
34  Options:
35  - split cells from cellSet only
36  - use geometric cut for hexes as well
37 
38  The angle is the angle between two faces sharing an edge as seen from
39  inside each cell. So a cube will have all angles 90. If you want
40  to split cells with cell centre outside use e.g. angle 200
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #include "argList.H"
45 #include "Time.H"
46 #include "polyTopoChange.H"
47 #include "polyTopoChangeMap.H"
48 #include "polyMesh.H"
49 #include "cellCuts.H"
50 #include "cellSet.H"
51 #include "cellModeller.H"
52 #include "meshCutter.H"
53 #include "geomCellLooper.H"
54 #include "plane.H"
55 #include "edgeVertex.H"
56 #include "meshTools.H"
57 #include "ListOps.H"
58 
59 using namespace Foam;
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 
64 labelList pack(const boolList& lst)
65 {
66  labelList packedLst(lst.size());
67  label packedI = 0;
68 
69  forAll(lst, i)
70  {
71  if (lst[i])
72  {
73  packedLst[packedI++] = i;
74  }
75  }
76  packedLst.setSize(packedI);
77 
78  return packedLst;
79 }
80 
81 
82 scalarField pack(const boolList& lst, const scalarField& elems)
83 {
84  scalarField packedElems(lst.size());
85  label packedI = 0;
86 
87  forAll(lst, i)
88  {
89  if (lst[i])
90  {
91  packedElems[packedI++] = elems[i];
92  }
93  }
94  packedElems.setSize(packedI);
95 
96  return packedElems;
97 }
98 
99 
100 // Given sin and cos of max angle between normals calculate whether f0 and f1
101 // on celli make larger angle. Uses sinAngle only for quadrant detection.
102 bool largerAngle
103 (
104  const primitiveMesh& mesh,
105  const scalar cosAngle,
106  const scalar sinAngle,
107 
108  const label celli,
109  const label f0, // face label
110  const label f1,
111 
112  const vector& n0, // normal at f0
113  const vector& n1
114 )
115 {
116  const labelList& own = mesh.faceOwner();
117 
118  bool sameFaceOrder = !((own[f0] == celli) ^ (own[f1] == celli));
119 
120  // Get cos between faceArea vectors. Correct so flat angle (180 degrees)
121  // gives -1.
122  scalar normalCosAngle = n0 & n1;
123 
124  if (sameFaceOrder)
125  {
126  normalCosAngle = -normalCosAngle;
127  }
128 
129 
130  // Get cos between faceCentre and normal vector to determine in
131  // which quadrant angle is. (Is correct for unwarped faces only!)
132  // Correct for non-outwards pointing normal.
133  vector c1c0(mesh.faceCentres()[f1] - mesh.faceCentres()[f0]);
134  c1c0 /= mag(c1c0) + vSmall;
135 
136  scalar fcCosAngle = n0 & c1c0;
137 
138  if (own[f0] != celli)
139  {
140  fcCosAngle = -fcCosAngle;
141  }
142 
143  if (sinAngle < 0.0)
144  {
145  // Looking for concave angles (quadrant 3 or 4)
146  if (fcCosAngle <= 0)
147  {
148  // Angle is convex so smaller.
149  return false;
150  }
151  else
152  {
153  if (normalCosAngle < cosAngle)
154  {
155  return false;
156  }
157  else
158  {
159  return true;
160  }
161  }
162  }
163  else
164  {
165  // Looking for convex angles (quadrant 1 or 2)
166  if (fcCosAngle > 0)
167  {
168  // Concave angle
169  return true;
170  }
171  else
172  {
173  // Convex. Check cos of normal vectors.
174  if (normalCosAngle > cosAngle)
175  {
176  return false;
177  }
178  else
179  {
180  return true;
181  }
182  }
183  }
184 }
185 
186 
187 // Split hex (and hex only) along edgeI creating two prisms
188 bool splitHex
189 (
190  const polyMesh& mesh,
191  const label celli,
192  const label edgeI,
193 
194  DynamicList<label>& cutCells,
195  DynamicList<labelList>& cellLoops,
196  DynamicList<scalarField>& cellEdgeWeights
197 )
198 {
199  // cut handling functions
200  edgeVertex ev(mesh);
201 
202  const edgeList& edges = mesh.edges();
203  const faceList& faces = mesh.faces();
204 
205  const edge& e = edges[edgeI];
206 
207  // Get faces on the side, i.e. faces not using edge but still using one of
208  // the edge endpoints.
209 
210  label leftI = -1;
211  label rightI = -1;
212  label leftFp = -1;
213  label rightFp = -1;
214 
215  const cell& cFaces = mesh.cells()[celli];
216 
217  forAll(cFaces, i)
218  {
219  label facei = cFaces[i];
220 
221  const face& f = faces[facei];
222 
223  label fp0 = findIndex(f, e[0]);
224  label fp1 = findIndex(f, e[1]);
225 
226  if (fp0 == -1)
227  {
228  if (fp1 != -1)
229  {
230  // Face uses e[1] but not e[0]
231  rightI = facei;
232  rightFp = fp1;
233 
234  if (leftI != -1)
235  {
236  // Have both faces so exit
237  break;
238  }
239  }
240  }
241  else
242  {
243  if (fp1 != -1)
244  {
245  // Face uses both e[1] and e[0]
246  }
247  else
248  {
249  leftI = facei;
250  leftFp = fp0;
251 
252  if (rightI != -1)
253  {
254  break;
255  }
256  }
257  }
258  }
259 
260  if (leftI == -1 || rightI == -1)
261  {
263  << " rightI:" << rightI << abort(FatalError);
264  }
265 
266  // Walk two vertices further on faces.
267 
268  const face& leftF = faces[leftI];
269 
270  label leftV = leftF[(leftFp + 2) % leftF.size()];
271 
272  const face& rightF = faces[rightI];
273 
274  label rightV = rightF[(rightFp + 2) % rightF.size()];
275 
276  labelList loop(4);
277  loop[0] = ev.vertToEVert(e[0]);
278  loop[1] = ev.vertToEVert(leftV);
279  loop[2] = ev.vertToEVert(rightV);
280  loop[3] = ev.vertToEVert(e[1]);
281 
282  scalarField loopWeights(4);
283  loopWeights[0] = -great;
284  loopWeights[1] = -great;
285  loopWeights[2] = -great;
286  loopWeights[3] = -great;
287 
288  cutCells.append(celli);
289  cellLoops.append(loop);
290  cellEdgeWeights.append(loopWeights);
291 
292  return true;
293 }
294 
295 
296 // Split celli along edgeI with a plane along halfNorm direction.
297 bool splitCell
298 (
299  const polyMesh& mesh,
300  const geomCellLooper& cellCutter,
301 
302  const label celli,
303  const label edgeI,
304  const vector& halfNorm,
305 
306  const boolList& vertIsCut,
307  const boolList& edgeIsCut,
308  const scalarField& edgeWeight,
309 
310  DynamicList<label>& cutCells,
311  DynamicList<labelList>& cellLoops,
312  DynamicList<scalarField>& cellEdgeWeights
313 )
314 {
315  const edge& e = mesh.edges()[edgeI];
316 
317  vector eVec = e.vec(mesh.points());
318  eVec /= mag(eVec);
319 
320  vector planeN = eVec ^ halfNorm;
321 
322  // Slightly tilt plane to make it not cut edges exactly
323  // halfway on fully regular meshes (since we want cuts
324  // to be snapped to vertices)
325  planeN += 0.01*halfNorm;
326 
327  planeN /= mag(planeN);
328 
329  // Define plane through edge
330  plane cutPlane(mesh.points()[e.start()], planeN);
331 
332  labelList loop;
333  scalarField loopWeights;
334 
335  if
336  (
337  cellCutter.cut
338  (
339  cutPlane,
340  celli,
341  vertIsCut,
342  edgeIsCut,
343  edgeWeight,
344  loop,
345  loopWeights
346  )
347  )
348  {
349  // Did manage to cut cell. Copy into overall list.
350  cutCells.append(celli);
351  cellLoops.append(loop);
352  cellEdgeWeights.append(loopWeights);
353 
354  return true;
355  }
356  else
357  {
358  return false;
359  }
360 }
361 
362 
363 // Collects cuts for all cells in cellSet
364 void collectCuts
365 (
366  const polyMesh& mesh,
367  const geomCellLooper& cellCutter,
368  const bool geometry,
369  const scalar minCos,
370  const scalar minSin,
371  const cellSet& cellsToCut,
372 
373  DynamicList<label>& cutCells,
374  DynamicList<labelList>& cellLoops,
375  DynamicList<scalarField>& cellEdgeWeights
376 )
377 {
378  // Get data from mesh
379  const cellShapeList& cellShapes = mesh.cellShapes();
380  const labelList& own = mesh.faceOwner();
381  const labelListList& cellEdges = mesh.cellEdges();
382  const vectorField& faceAreas = mesh.faceAreas();
383 
384  // Hex shape
385  const cellModel& hex = *(cellModeller::lookup("hex"));
386 
387  // cut handling functions
388  edgeVertex ev(mesh);
389 
390 
391  // Cut information per mesh entity
392  boolList vertIsCut(mesh.nPoints(), false);
393  boolList edgeIsCut(mesh.nEdges(), false);
394  scalarField edgeWeight(mesh.nEdges(), -great);
395 
396  forAllConstIter(cellSet, cellsToCut, iter)
397  {
398  const label celli = iter.key();
399  const labelList& cEdges = cellEdges[celli];
400 
401  forAll(cEdges, i)
402  {
403  label edgeI = cEdges[i];
404 
405  label f0, f1;
406  meshTools::getEdgeFaces(mesh, celli, edgeI, f0, f1);
407 
408  vector n0 = faceAreas[f0];
409  n0 /= mag(n0);
410 
411  vector n1 = faceAreas[f1];
412  n1 /= mag(n1);
413 
414  if
415  (
416  largerAngle
417  (
418  mesh,
419  minCos,
420  minSin,
421 
422  celli,
423  f0,
424  f1,
425  n0,
426  n1
427  )
428  )
429  {
430  bool splitOk = false;
431 
432  if (!geometry && cellShapes[celli].model() == hex)
433  {
434  splitOk =
435  splitHex
436  (
437  mesh,
438  celli,
439  edgeI,
440 
441  cutCells,
442  cellLoops,
443  cellEdgeWeights
444  );
445  }
446  else
447  {
448  vector halfNorm;
449 
450  if ((own[f0] == celli) ^ (own[f1] == celli))
451  {
452  // Opposite owner orientation
453  halfNorm = 0.5*(n0 - n1);
454  }
455  else
456  {
457  // Faces have same owner or same neighbour so
458  // normals point in same direction
459  halfNorm = 0.5*(n0 + n1);
460  }
461 
462  splitOk =
463  splitCell
464  (
465  mesh,
466  cellCutter,
467  celli,
468  edgeI,
469  halfNorm,
470 
471  vertIsCut,
472  edgeIsCut,
473  edgeWeight,
474 
475  cutCells,
476  cellLoops,
477  cellEdgeWeights
478  );
479  }
480 
481 
482  if (splitOk)
483  {
484  // Update cell/edge/vertex wise info.
485  label index = cellLoops.size() - 1;
486  const labelList& loop = cellLoops[index];
487  const scalarField& loopWeights = cellEdgeWeights[index];
488 
489  forAll(loop, i)
490  {
491  label cut = loop[i];
492 
493  if (ev.isEdge(cut))
494  {
495  edgeIsCut[ev.getEdge(cut)] = true;
496  edgeWeight[ev.getEdge(cut)] = loopWeights[i];
497  }
498  else
499  {
500  vertIsCut[ev.getVertex(cut)] = true;
501  }
502  }
503 
504  // Stop checking edges for this cell.
505  break;
506  }
507  }
508  }
509  }
510 
511  cutCells.shrink();
512  cellLoops.shrink();
513  cellEdgeWeights.shrink();
514 }
515 
516 
517 
518 int main(int argc, char *argv[])
519 {
521  (
522  "split cells with flat faces"
523  );
524  #include "addOverwriteOption.H"
526  argList::validArgs.append("edgeAngle [0..360]");
527 
529  (
530  "set",
531  "name",
532  "split cells from specified cellSet only"
533  );
535  (
536  "geometry",
537  "use geometric cut for hexes as well"
538  );
540  (
541  "tol",
542  "scalar", "edge snap tolerance (default 0.2)"
543  );
544 
545  #include "setRootCase.H"
547  #include "createPolyMesh.H"
548  const word oldInstance = mesh.pointsInstance();
549 
550  const scalar featureAngle = degToRad(args.argRead<scalar>(1));
551  const scalar minCos = Foam::cos(featureAngle);
552  const scalar minSin = Foam::sin(featureAngle);
553 
554  const bool readSet = args.optionFound("set");
555  const bool geometry = args.optionFound("geometry");
556  const bool overwrite = args.optionFound("overwrite");
557 
558  const scalar edgeTol = args.optionLookupOrDefault("tol", 0.2);
559 
560  Info<< "Trying to split cells with internal angles > feature angle\n" << nl
561  << "featureAngle : " << radToDeg(featureAngle) << nl
562  << "edge snapping tol : " << edgeTol << nl;
563  if (readSet)
564  {
565  Info<< "candidate cells : cellSet " << args["set"] << nl;
566  }
567  else
568  {
569  Info<< "candidate cells : all cells" << nl;
570  }
571  if (geometry)
572  {
573  Info<< "hex cuts : geometric; using edge tolerance" << nl;
574  }
575  else
576  {
577  Info<< "hex cuts : topological; cut to opposite edge" << nl;
578  }
579  Info<< endl;
580 
581 
582  // Cell circumference cutter
583  geomCellLooper cellCutter(mesh);
584  // Snap all edge cuts close to endpoints to vertices.
586 
587  // Candidate cells to cut
588  cellSet cellsToCut(mesh, "cellsToCut", mesh.nCells()/100);
589 
590  if (readSet)
591  {
592  // Read cells to cut from cellSet
593  cellSet cells(mesh, args["set"]);
594 
595  cellsToCut = cells;
596  }
597 
598  while (true)
599  {
600  if (!readSet)
601  {
602  // Try all cells for cutting
603  for (label celli = 0; celli < mesh.nCells(); celli++)
604  {
605  cellsToCut.insert(celli);
606  }
607  }
608 
609 
610  // Cut information per cut cell
611  DynamicList<label> cutCells(mesh.nCells()/10 + 10);
612  DynamicList<labelList> cellLoops(mesh.nCells()/10 + 10);
613  DynamicList<scalarField> cellEdgeWeights(mesh.nCells()/10 + 10);
614 
615  collectCuts
616  (
617  mesh,
618  cellCutter,
619  geometry,
620  minCos,
621  minSin,
622  cellsToCut,
623 
624  cutCells,
625  cellLoops,
626  cellEdgeWeights
627  );
628 
629  cellSet cutSet(mesh, "cutSet", cutCells.size());
630  forAll(cutCells, i)
631  {
632  cutSet.insert(cutCells[i]);
633  }
634 
635  // Gets cuts across cells from cuts through edges.
636  Info<< "Writing " << cutSet.size() << " cells to cut to cellSet "
637  << cutSet.instance()/cutSet.local()/cutSet.name()
638  << endl << endl;
639  cutSet.write();
640 
641  // Analyze cuts for clashes.
642  cellCuts cuts
643  (
644  mesh,
645  cutCells, // cells candidate for cutting
646  cellLoops,
647  cellEdgeWeights
648  );
649 
650  Info<< "Actually cut cells:" << cuts.nLoops() << nl << endl;
651 
652  if (cuts.nLoops() == 0)
653  {
654  break;
655  }
656 
657  // Remove cut cells from cellsToCut (Note:only relevant if -readSet)
658  forAll(cuts.cellLoops(), celli)
659  {
660  if (cuts.cellLoops()[celli].size())
661  {
662  // Info<< "Removing cut cell " << celli << " from wishlist"
663  // << endl;
664  cellsToCut.erase(celli);
665  }
666  }
667 
668  // At least some cells are cut.
669  polyTopoChange meshMod(mesh);
670 
671  // Cutting engine
672  meshCutter cutter(mesh);
673 
674  // Insert mesh refinement into polyTopoChange.
675  cutter.setRefinement(cuts, meshMod);
676 
677  // Do all changes
678  Info<< "Morphing ..." << endl;
679 
680  if (!overwrite)
681  {
682  runTime++;
683  }
684 
685  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
686 
687  // Update stored labels on meshCutter
688  cutter.topoChange(map());
689 
690  // Update cellSet
691  cellsToCut.topoChange(map());
692 
693  Info<< "Remaining:" << cellsToCut.size() << endl;
694 
695  // Write resulting mesh
696  if (overwrite)
697  {
698  mesh.setInstance(oldInstance);
699  }
700 
701  Info<< "Writing refined morphMesh to time " << runTime.name()
702  << endl;
703 
704  mesh.write();
705  }
706 
707  Info<< "End\n" << endl;
708 
709  return 0;
710 }
711 
712 
713 // ************************************************************************* //
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:445
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
virtual Ostream & write(const char)=0
Write character.
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:183
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:243
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Description of cuts across cells.
Definition: cellCuts.H:111
Maps a geometry to a set of cell primitives, which enables geometric cell data to be calculated witho...
Definition: cellModel.H:65
static const cellModel * lookup(const word &)
Look up a model by name and return a pointer to the model or nullptr.
Definition: cellModeller.C:100
A collection of cell labels.
Definition: cellSet.H:51
virtual void topoChange(const polyTopoChangeMap &map)
Update any stored data for new labels.
Definition: cellSet.C:225
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
Combines edge or vertex in single label. Used to specify cuts across cell circumference.
Definition: edgeVertex.H:53
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Implementation of cellLooper. Does pure geometric cut through cell.
static void setSnapTol(const scalar tol)
virtual bool cut(const vector &refDir, const label celli, const boolList &vertIsCut, const boolList &edgeIsCut, const scalarField &edgeWeight, labelList &loop, scalarField &loopWeights) const
Create cut along circumference of celli. Gets current mesh cuts.
Cuts (splits) cells.
Definition: meshCutter.H:131
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:61
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:965
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
Direct mesh changes based on v1.3 polyTopoChange syntax.
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
label nEdges() const
label nCells() const
label nPoints() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const labelListList & cellEdges() const
const vectorField & faceCentres() const
const cellShapeList & cellShapes() const
Return cell shapes.
virtual const labelList & faceOwner() const =0
Face face-owner addressing.
const vectorField & faceAreas() const
const cellList & cells() const
virtual bool write(const bool write=true) const
Write using setting from DB.
Description of cell after splitting. Contains cellLabel and pointers to cells it split in....
Definition: splitCell.H:52
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const cellShapeList & cellShapes
const cellShapeList & cells
void getEdgeFaces(const primitiveMesh &, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition: meshTools.C:458
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
scalar radToDeg(const scalar rad)
Convert radians to degrees.
scalar degToRad(const scalar deg)
Convert degrees to radians.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar sin(const dimensionedScalar &ds)
IOstream & hex(IOstream &io)
Definition: IOstream.H:561
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
static const char nl
Definition: Ostream.H:266
dimensionedScalar cos(const dimensionedScalar &ds)
labelList f(nPoints)
Foam::argList args(argc, argv)