All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
cellShapeControlMesh.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) 2012-2022 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 \*---------------------------------------------------------------------------*/
25 
26 #include "cellShapeControlMesh.H"
28 #include "pointFields.H"
29 #include "pointIOField.H"
30 #include "scalarIOField.H"
31 #include "triadIOField.H"
32 #include "tetPointRef.H"
33 #include "plane.H"
34 #include "transform.H"
35 #include "meshTools.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(cellShapeControlMesh, 0);
42 
43  word cellShapeControlMesh::meshSubDir = "cellShapeControlMesh";
44 }
45 
46 
47 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
48 
49 //Foam::tensor Foam::cellShapeControlMesh::requiredAlignment
50 //(
51 // const Foam::point& pt,
52 // const searchableSurfaces& allGeometry,
53 // const conformationSurfaces& geometryToConformTo
54 //) const
55 //{
56 // pointIndexHit surfHit;
57 // label hitSurface;
58 //
59 // geometryToConformTo.findSurfaceNearest
60 // (
61 // pt,
62 // sqr(great),
63 // surfHit,
64 // hitSurface
65 // );
66 //
67 // if (!surfHit.hit())
68 // {
69 // FatalErrorInFunction
70 // << "findSurfaceNearest did not find a hit across the surfaces."
71 // << exit(FatalError) << endl;
72 // }
73 //
74 // // Primary alignment
75 //
76 // vectorField norm(1);
77 //
78 // allGeometry[hitSurface].getNormal
79 // (
80 // List<pointIndexHit>(1, surfHit),
81 // norm
82 // );
83 //
84 // const vector np = norm[0];
85 //
86 // // Generate equally spaced 'spokes' in a circle normal to the
87 // // direction from the vertex to the closest point on the surface
88 // // and look for a secondary intersection.
89 //
90 // const vector d = surfHit.hitPoint() - pt;
91 //
92 // const tensor Rp = rotationTensor(vector(0,0,1), np);
93 //
94 // const label s = 36;//foamyHexMeshControls().alignmentSearchSpokes();
95 //
96 // scalar closestSpokeHitDistance = great;
97 //
98 // pointIndexHit closestSpokeHit;
99 //
100 // label closestSpokeSurface = -1;
101 //
102 // const scalar spanMag = geometryToConformTo.globalBounds().mag();
103 //
104 // for (label i = 0; i < s; i++)
105 // {
106 // vector spoke
107 // (
108 // Foam::cos(i*constant::mathematical::twoPi/s),
109 // Foam::sin(i*constant::mathematical::twoPi/s),
110 // 0
111 // );
112 //
113 // spoke *= spanMag;
114 //
115 // spoke = Rp & spoke;
116 //
117 // pointIndexHit spokeHit;
118 //
119 // label spokeSurface = -1;
120 //
121 // // internal spoke
122 //
123 // geometryToConformTo.findSurfaceNearestIntersection
124 // (
125 // pt,
126 // pt + spoke,
127 // spokeHit,
128 // spokeSurface
129 // );
130 //
131 // if (spokeHit.hit())
132 // {
133 // scalar spokeHitDistance = mag
134 // (
135 // spokeHit.hitPoint() - pt
136 // );
137 //
138 // if (spokeHitDistance < closestSpokeHitDistance)
139 // {
140 // closestSpokeHit = spokeHit;
141 // closestSpokeSurface = spokeSurface;
142 // closestSpokeHitDistance = spokeHitDistance;
143 // }
144 // }
145 //
146 // // external spoke
147 //
148 // Foam::point mirrorPt = pt + 2*d;
149 //
150 // geometryToConformTo.findSurfaceNearestIntersection
151 // (
152 // mirrorPt,
153 // mirrorPt + spoke,
154 // spokeHit,
155 // spokeSurface
156 // );
157 //
158 // if (spokeHit.hit())
159 // {
160 // scalar spokeHitDistance = mag
161 // (
162 // spokeHit.hitPoint() - mirrorPt
163 // );
164 //
165 // if (spokeHitDistance < closestSpokeHitDistance)
166 // {
167 // closestSpokeHit = spokeHit;
168 // closestSpokeSurface = spokeSurface;
169 // closestSpokeHitDistance = spokeHitDistance;
170 // }
171 // }
172 // }
173 //
174 // if (closestSpokeSurface == -1)
175 // {
181 //
182 // return I;
183 // }
184 //
185 // // Auxiliary alignment generated by spoke intersection normal.
186 //
187 // allGeometry[closestSpokeSurface].getNormal
188 // (
189 // List<pointIndexHit>(1, closestSpokeHit),
190 // norm
191 // );
192 //
193 // const vector& na = norm[0];
194 //
195 // // Secondary alignment
196 // vector ns = np ^ na;
197 //
198 // if (mag(ns) < small)
199 // {
200 // FatalErrorInFunction
201 // << "Parallel normals detected in spoke search." << nl
202 // << "point: " << pt << nl
203 // << "closest surface point: " << surfHit.hitPoint() << nl
204 // << "closest spoke hit: " << closestSpokeHit.hitPoint() << nl
205 // << "np: " << surfHit.hitPoint() + np << nl
206 // << "ns: " << closestSpokeHit.hitPoint() + na << nl
207 // << exit(FatalError);
208 // }
209 //
210 // ns /= mag(ns);
211 //
212 // tensor Rs = rotationTensor((Rp & vector(0,1,0)), ns);
213 //
214 // return (Rs & Rp);
215 //}
216 
217 
219 {
220  label nRemoved = 0;
221  for
222  (
223  CellSizeDelaunay::Finite_vertices_iterator vit =
224  finite_vertices_begin();
225  vit != finite_vertices_end();
226  ++vit
227  )
228  {
229  std::list<Vertex_handle> verts;
230  adjacent_vertices(vit, std::back_inserter(verts));
231 
232  bool removePt = true;
233  for
234  (
235  std::list<Vertex_handle>::iterator aVit = verts.begin();
236  aVit != verts.end();
237  ++aVit
238  )
239  {
240  Vertex_handle avh = *aVit;
241 
242  scalar diff =
243  mag(avh->targetCellSize() - vit->targetCellSize())
244  /max(vit->targetCellSize(), 1e-6);
245 
246  if (diff > 0.05)
247  {
248  removePt = false;
249  }
250  }
251 
252  if (removePt)
253  {
254  remove(vit);
255  nRemoved++;
256  }
257  }
258 
259  return nRemoved;
260 }
261 
262 
264 {
265  tmp<pointField> tcellCentres(new pointField(number_of_finite_cells()));
266  pointField& cellCentres = tcellCentres.ref();
267 
268  label count = 0;
269  for
270  (
271  CellSizeDelaunay::Finite_cells_iterator c = finite_cells_begin();
272  c != finite_cells_end();
273  ++c
274  )
275  {
276  if (c->hasFarPoint())
277  {
278  continue;
279  }
280 
281  scalarList bary;
283 
284  const Foam::point centre = topoint
285  (
286  CGAL::centroid<baseK>
287  (
288  c->vertex(0)->point(),
289  c->vertex(1)->point(),
290  c->vertex(2)->point(),
291  c->vertex(3)->point()
292  )
293  );
294 
295  cellCentres[count++] = centre;
296  }
297 
298  cellCentres.resize(count);
299 
300  return tcellCentres;
301 }
302 
303 
305 {
306  OFstream str
307  (
308  "refinementTriangulation_"
310  + ".obj"
311  );
312 
313  label count = 0;
314 
315  Info<< "Write refinementTriangulation" << endl;
316 
317  for
318  (
319  CellSizeDelaunay::Finite_edges_iterator e = finite_edges_begin();
320  e != finite_edges_end();
321  ++e
322  )
323  {
324  Cell_handle c = e->first;
325  Vertex_handle vA = c->vertex(e->second);
326  Vertex_handle vB = c->vertex(e->third);
327 
328  // Don't write far edges
329  if (vA->farPoint() || vB->farPoint())
330  {
331  continue;
332  }
333 
334  // Don't write unowned edges
335  if (vA->referred() && vB->referred())
336  {
337  continue;
338  }
339 
340  pointFromPoint p1 = topoint(vA->point());
341  pointFromPoint p2 = topoint(vB->point());
342 
343  meshTools::writeOBJ(str, p1, p2, count);
344  }
345 
346  if (is_valid())
347  {
348  Info<< " Triangulation is valid" << endl;
349  }
350  else
351  {
353  << "Triangulation is not valid"
354  << abort(FatalError);
355  }
356 }
357 
358 
359 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
360 
362 :
363  DistributedDelaunayMesh<CellSizeDelaunay>
364  (
365  runTime,
366  meshSubDir
367  ),
368  runTime_(runTime)
369 {
370  if (this->vertexCount())
371  {
372  fvMesh mesh
373  (
374  IOobject
375  (
376  meshSubDir,
377  runTime.timeName(),
378  runTime,
381  ),
382  false
383  );
384 
385  if (mesh.nPoints() == this->vertexCount())
386  {
387  pointScalarField sizes
388  (
389  IOobject
390  (
391  "sizes",
392  runTime.timeName(),
393  meshSubDir,
394  runTime,
397  false
398  ),
399  pointMesh::New(mesh)
400  );
401 
402  triadIOField alignments
403  (
404  IOobject
405  (
406  "alignments",
407  mesh.time().timeName(),
408  meshSubDir,
409  mesh.time(),
412  false
413  )
414  );
415 
416  if
417  (
418  sizes.size() == this->vertexCount()
419  && alignments.size() == this->vertexCount()
420  )
421  {
422  for
423  (
424  Finite_vertices_iterator vit = finite_vertices_begin();
425  vit != finite_vertices_end();
426  ++vit
427  )
428  {
429  vit->targetCellSize() = sizes[vit->index()];
430  vit->alignment() = alignments[vit->index()];
431  }
432  }
433  else
434  {
436  << "Cell size point field is not the same size as the "
437  << "mesh."
438  << abort(FatalError);
439  }
440  }
441  }
442 }
443 
444 
445 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
446 
448 {}
449 
450 
451 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
452 
454 (
455  const Foam::point& pt,
456  barycentric& bary,
457  Cell_handle& ch
458 ) const
459 {
460  // Use the previous cell handle as a hint on where to start searching
461  // Giving a hint causes strange errors...
462  ch = locate(toPoint(pt));
463 
464  if (dimension() > 2 && !is_infinite(ch))
465  {
466  oldCellHandle_ = ch;
467 
468  tetPointRef tet
469  (
470  topoint(ch->vertex(0)->point()),
471  topoint(ch->vertex(1)->point()),
472  topoint(ch->vertex(2)->point()),
473  topoint(ch->vertex(3)->point())
474  );
475 
476  bary = tet.pointToBarycentric(pt);
477  }
478 }
479 
480 
482 {
483  DynamicList<Foam::point> pts(number_of_vertices());
484 
485  for
486  (
487  Finite_vertices_iterator vit = finite_vertices_begin();
488  vit != finite_vertices_end();
489  ++vit
490  )
491  {
492  if (vit->real())
493  {
494  pts.append(topoint(vit->point()));
495  }
496  }
497 
498  boundBox bb(pts);
499 
500  return bb;
501 }
502 
503 
505 (
506  const backgroundMeshDecomposition& decomposition
507 )
508 {
509  DynamicList<Foam::point> points(number_of_vertices());
510  DynamicList<scalar> sizes(number_of_vertices());
511  DynamicList<tensor> alignments(number_of_vertices());
512 
513  DynamicList<Vb> farPts(8);
514 
515  for
516  (
517  Finite_vertices_iterator vit = finite_vertices_begin();
518  vit != finite_vertices_end();
519  ++vit
520  )
521  {
522  if (vit->real())
523  {
524  points.append(topoint(vit->point()));
525  sizes.append(vit->targetCellSize());
526  alignments.append(vit->alignment());
527  }
528  else if (vit->farPoint())
529  {
530  farPts.append
531  (
532  Vb
533  (
534  vit->point(),
535  -1,
536  Vb::vtFar,
538  )
539  );
540 
541  farPts.last().targetCellSize() = vit->targetCellSize();
542  farPts.last().alignment() = vit->alignment();
543  }
544  }
545 
546  autoPtr<distributionMap> mapDist =
548  (
549  decomposition,
550  points
551  );
552 
553  mapDist().distribute(sizes);
554  mapDist().distribute(alignments);
555 
556  // Reset the entire tessellation
558 
559 
560  // Internal points have to be inserted first
561  DynamicList<Vb> verticesToInsert(points.size());
562 
563 
564  forAll(farPts, ptI)
565  {
566  verticesToInsert.append(farPts[ptI]);
567  }
568 
569 
570  forAll(points, pI)
571  {
572  verticesToInsert.append
573  (
574  Vb
575  (
576  toPoint(points[pI]),
577  -1,
580  )
581  );
582 
583  verticesToInsert.last().targetCellSize() = sizes[pI];
584  verticesToInsert.last().alignment() = alignments[pI];
585  }
586 
587  Info<< nl << " Inserting distributed background tessellation..." << endl;
588 
589  this->rangeInsertWithInfo
590  (
591  verticesToInsert.begin(),
592  verticesToInsert.end(),
593  true
594  );
595 
596  sync(decomposition.procBounds());
597 
598  Info<< " Total number of vertices after redistribution "
599  << returnReduce(label(number_of_vertices()), sumOp<label>()) << endl;
600 }
601 
602 
604 {
605  tensorField alignmentsTmp(number_of_vertices(), Zero);
606 
607  label count = 0;
608  for
609  (
610  Finite_vertices_iterator vit = finite_vertices_begin();
611  vit != finite_vertices_end();
612  ++vit
613  )
614  {
615  alignmentsTmp[count++] = vit->alignment();
616  }
617 
618  return alignmentsTmp;
619 }
620 
621 
623 {
624  Info<< "Writing " << meshSubDir << endl;
625 
626  // Reindex the cells
627  label cellCount = 0;
628  for
629  (
630  Finite_cells_iterator cit = finite_cells_begin();
631  cit != finite_cells_end();
632  ++cit
633  )
634  {
635  if (!cit->hasFarPoint() && !is_infinite(cit))
636  {
637  cit->cellIndex() = cellCount++;
638  }
639  }
640 
642  labelList cellMap;
643 
644  autoPtr<polyMesh> meshPtr = DelaunayMesh<CellSizeDelaunay>::createMesh
645  (
646  meshSubDir,
647  vertexMap,
648  cellMap
649  );
650  const polyMesh& mesh = meshPtr();
651 
652  pointScalarField sizes
653  (
654  IOobject
655  (
656  "sizes",
657  mesh.time().timeName(),
658  meshSubDir,
659  mesh.time(),
662  ),
663  pointMesh::New(mesh),
664  scalar(0)
665  );
666 
667  triadIOField alignments
668  (
669  IOobject
670  (
671  "alignments",
672  mesh.time().timeName(),
673  meshSubDir,
674  mesh.time(),
677  ),
678  sizes.size()
679  );
680 
681  // Write alignments
682 // OFstream str(runTime_.path()/"alignments.obj");
683 
684  for
685  (
686  Finite_vertices_iterator vit = finite_vertices_begin();
687  vit != finite_vertices_end();
688  ++vit
689  )
690  {
691  if (!vit->farPoint())
692  {
693  // Populate sizes
694  sizes[vertexMap[labelPair(vit->index(), vit->procIndex())]] =
695  vit->targetCellSize();
696 
697  alignments[vertexMap[labelPair(vit->index(), vit->procIndex())]] =
698  vit->alignment();
699 
700 // // Write alignments
701 // const tensor& alignment = vit->alignment();
702 // pointFromPoint pt = topoint(vit->point());
703 //
704 // if
705 // (
706 // alignment.x() == triad::unset[0]
707 // || alignment.y() == triad::unset[0]
708 // || alignment.z() == triad::unset[0]
709 // )
710 // {
711 // Info<< "Bad alignment = " << vit->info();
712 //
713 // vit->alignment() = tensor::I;
714 //
715 // Info<< "New alignment = " << vit->info();
716 //
717 // continue;
718 // }
719 //
720 // meshTools::writeOBJ(str, pt, alignment.x() + pt);
721 // meshTools::writeOBJ(str, pt, alignment.y() + pt);
722 // meshTools::writeOBJ(str, pt, alignment.z() + pt);
723  }
724  }
725 
726  mesh.write();
727  sizes.write();
728  alignments.write();
729 }
730 
731 
733 (
734  const autoPtr<backgroundMeshDecomposition>& decomposition
735 ) const
736 {
737  // Loop over all the tets and estimate the cell count in each one
738 
739  scalar cellCount = 0;
740 
741  for
742  (
743  Finite_cells_iterator cit = finite_cells_begin();
744  cit != finite_cells_end();
745  ++cit
746  )
747  {
748  if (!cit->hasFarPoint() && !is_infinite(cit))
749  {
750  // TODO: Check if tet centre is on the processor..
751  CGAL::Tetrahedron_3<baseK> tet
752  (
753  cit->vertex(0)->point(),
754  cit->vertex(1)->point(),
755  cit->vertex(2)->point(),
756  cit->vertex(3)->point()
757  );
758 
759  const Foam::point centre(topoint(CGAL::centroid(tet)));
760 
761  if
762  (
764  && !decomposition().positionOnThisProcessor(centre)
765  )
766  {
767  continue;
768  }
769 
770  scalar volume = CGAL::to_double(tet.volume());
771 
772  scalar averagedPointCellSize = 0;
773  // scalar averagedPointCellSize = 1;
774 
775  // Get an average volume by averaging the cell size of the vertices
776  for (label vI = 0; vI < 4; ++vI)
777  {
778  averagedPointCellSize += cit->vertex(vI)->targetCellSize();
779  // averagedPointCellSize *= cit->vertex(vI)->targetCellSize();
780  }
781 
782  averagedPointCellSize /= 4;
783  // averagedPointCellSize = ::sqrt(averagedPointCellSize);
784 
785 // if (averagedPointCellSize < small)
786 // {
787 // Pout<< "Volume = " << volume << endl;
788 //
789 // for (label vI = 0; vI < 4; ++vI)
790 // {
791 // Pout<< "Point " << vI
792 // << ", point = " << topoint(cit->vertex(vI)->point())
793 // << ", size = " << cit->vertex(vI)->targetCellSize()
794 // << endl;
795 // }
796 // }
797 
798  cellCount += volume/pow(averagedPointCellSize, 3);
799  }
800  }
801 
802  return cellCount;
803 }
804 
805 
806 // ************************************************************************* //
tetrahedron< point, const point & > tetPointRef
Definition: tetPointRef.H:44
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
~cellShapeControlMesh()
Destructor.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
pointFromPoint topoint(const Point &P)
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
tmp< pointField > cellCentres() const
Get the centres of all the tets.
static fvMesh * meshPtr
Definition: globalFoam.H:52
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
HashTable< label, labelPair, FixedList< label, 2 >::Hash<> > labelTolabelPairHashTable
Definition: DelaunayMesh.H:90
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool distribute(const boundBox &bb)
Holds information (coordinate and normal) regarding nearest wall point.
void barycentricCoords(const Foam::point &pt, barycentric &bary, Cell_handle &ch) const
Calculate and return the barycentric coordinates for.
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:45
fvMesh & mesh
const dimensionedScalar c
Speed of light in a vacuum.
tensorField dumpAlignments() const
An indexed form of CGAL::Triangulation_vertex_base_3<K> used to keep track of the Delaunay vertices i...
Definition: indexedVertex.H:51
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const pointField & points
3D tensor transformation operations.
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
void reset()
Clear the entire triangulation.
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
PointFrompoint toPoint(const Foam::point &p)
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
void distribute(const backgroundMeshDecomposition &decomposition)
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
errorManip< error > abort(error &err)
Definition: errorManip.H:131
cellShapeControlMesh(const Time &runTime)
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
CGAL::Delaunay_triangulation_3< K, Tds, FastLocator > CellSizeDelaunay
static word meshSubDir
Return the mesh sub-directory name (usually "cellShapeControlMesh")
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
CellSizeDelaunay::Cell_handle Cell_handle
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
boundBox bounds() const
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
label estimateCellCount(const autoPtr< backgroundMeshDecomposition > &decomposition) const
A class for managing temporary objects.
Definition: PtrList.H:53
IOField< triad > triadIOField
triadField with IO.
Definition: triadIOField.H:42
autoPtr< polyMesh > createMesh(const fileName &name, labelTolabelPairHashTable &vertexMap, labelList &cellMap, const bool writeDelaunayData=true) const
Create an fvMesh from the triangulation.
Namespace for OpenFOAM.