sampledTriSurfaceMesh.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-2023 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 "sampledTriSurfaceMesh.H"
27 #include "meshSearch.H"
28 #include "treeDataCell.H"
29 #include "treeDataFace.H"
30 #include "OBJstream.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace sampledSurfaces
38 {
41 
42  //- Private class for finding nearest
43  // Comprising:
44  // - global index
45  // - sqr(distance)
47 
49  {
50 
51  public:
52 
53  void operator()(nearInfo& x, const nearInfo& y) const
54  {
55  if (y.first() < x.first())
56  {
57  x = y;
58  }
59  }
60  };
61 }
62 }
63 
64 template<>
65 const char* Foam::NamedEnum
66 <
68  3
69 >::names[] =
70 {
71  "cells",
72  "insideCells",
73  "boundaryFaces"
74 };
75 
76 const Foam::NamedEnum
77 <
79  3
80 > Foam::sampledSurfaces::triSurfaceMesh::samplingSourceNames_;
81 
82 
83 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
84 
86 Foam::sampledSurfaces::triSurfaceMesh::nonCoupledboundaryTree() const
87 {
88  // Variant of meshSearch::boundaryTree() that only does non-coupled
89  // boundary faces.
90 
91  if (!boundaryTreePtr_.valid())
92  {
93  // all non-coupled boundary faces (not just walls)
95 
96  labelList bndFaces(mesh().nFaces()-mesh().nInternalFaces());
97  label bndI = 0;
99  {
100  const polyPatch& pp = patches[patchi];
101  if (!pp.coupled())
102  {
103  forAll(pp, i)
104  {
105  bndFaces[bndI++] = pp.start()+i;
106  }
107  }
108  }
109  bndFaces.setSize(bndI);
110 
111 
112  treeBoundBox overallBb(mesh().points());
113  overallBb = overallBb.extend(1e-4);
114 
115  boundaryTreePtr_.reset
116  (
117  new indexedOctree<treeDataFace>
118  (
119  treeDataFace // all information needed to search faces
120  (
121  false, // do not cache bb
122  mesh(),
123  bndFaces // boundary faces only
124  ),
125  overallBb, // overall search domain
126  8, // maxLevel
127  10, // leafsize
128  3.0 // duplicity
129  )
130  );
131  }
132 
133  return boundaryTreePtr_();
134 }
135 
136 
138 (
139  const meshSearch& meshSearcher
140 )
141 {
142  // Find the cells the triangles of the surface are in.
143  // Does approximation by looking at the face centres only
144  const pointField& fc = surface_.faceCentres();
145 
146  List<nearInfo> nearest(fc.size());
147 
148  // Global numbering for cells/faces - only used to uniquely identify local
149  // elements
150  globalIndex globalCells
151  (
152  (sampleSource_ == cells || sampleSource_ == insideCells)
153  ? mesh().nCells()
154  : mesh().nFaces()
155  );
156 
157  forAll(nearest, i)
158  {
159  nearest[i].first() = great;
160  nearest[i].second() = labelMax;
161  }
162 
163  if (sampleSource_ == cells)
164  {
165  // Search for nearest cell
166 
167  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
168 
169  forAll(fc, triI)
170  {
171  pointIndexHit nearInfo = cellTree.findNearest
172  (
173  fc[triI],
174  sqr(great)
175  );
176  if (nearInfo.hit())
177  {
178  nearest[triI].first() = magSqr(nearInfo.hitPoint()-fc[triI]);
179  nearest[triI].second() = globalCells.toGlobal(nearInfo.index());
180  }
181  }
182  }
183  else if (sampleSource_ == insideCells)
184  {
185  // Search for cell containing point
186 
187  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
188 
189  forAll(fc, triI)
190  {
191  if (cellTree.bb().contains(fc[triI]))
192  {
193  label index = cellTree.findInside(fc[triI]);
194  if (index != -1)
195  {
196  nearest[triI].first() = 0.0;
197  nearest[triI].second() = globalCells.toGlobal(index);
198  }
199  }
200  }
201  }
202  else
203  {
204  // Search on all non-coupled boundary faces
205 
206  const indexedOctree<treeDataFace>& bTree = nonCoupledboundaryTree();
207 
208  forAll(fc, triI)
209  {
210  pointIndexHit nearInfo = bTree.findNearest
211  (
212  fc[triI],
213  sqr(great)
214  );
215  if (nearInfo.hit())
216  {
217  nearest[triI].first() = magSqr(nearInfo.hitPoint()-fc[triI]);
218  nearest[triI].second() = globalCells.toGlobal
219  (
220  bTree.shapes().faceLabels()[nearInfo.index()]
221  );
222  }
223  }
224  }
225 
226 
227  // See which processor has the nearest. Mark and subset
228  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
229 
230  Pstream::listCombineGather(nearest, nearestEqOp());
232 
233  labelList cellOrFaceLabels(fc.size(), -1);
234 
235  label nFound = 0;
236  forAll(nearest, triI)
237  {
238  if (nearest[triI].second() == labelMax)
239  {
240  // Not found on any processor. How to map?
241  }
242  else if (globalCells.isLocal(nearest[triI].second()))
243  {
244  cellOrFaceLabels[triI] = globalCells.toLocal
245  (
246  nearest[triI].second()
247  );
248  nFound++;
249  }
250  }
251 
252 
253  if (debug)
254  {
255  Pout<< "Local out of faces:" << cellOrFaceLabels.size()
256  << " keeping:" << nFound << endl;
257  }
258 
259  // Now subset the surface. Do not use triSurface::subsetMesh since requires
260  // original surface to be in compact numbering.
261 
262  const triSurface& s = surface_;
263 
264  // Compact to original triangle
265  labelList faceMap(s.size());
266  // Compact to original points
267  labelList pointMap(s.points().size());
268  // From original point to compact points
269  labelList reversePointMap(s.points().size(), -1);
270 
271  {
272  label newPointi = 0;
273  label newFacei = 0;
274 
275  forAll(s, facei)
276  {
277  if (cellOrFaceLabels[facei] != -1)
278  {
279  faceMap[newFacei++] = facei;
280 
281  const triSurface::FaceType& f = s[facei];
282  forAll(f, fp)
283  {
284  if (reversePointMap[f[fp]] == -1)
285  {
286  pointMap[newPointi] = f[fp];
287  reversePointMap[f[fp]] = newPointi++;
288  }
289  }
290  }
291  }
292  faceMap.setSize(newFacei);
293  pointMap.setSize(newPointi);
294  }
295 
296 
297  // Subset cellOrFaceLabels
298  cellOrFaceLabels = UIndirectList<label>(cellOrFaceLabels, faceMap)();
299 
300  // Store any face per point (without using pointFaces())
301  labelList pointToFace(pointMap.size());
302 
303  // Create faces and points for subsetted surface
304  faceList& faces = this->storedFaces();
305  faces.setSize(faceMap.size());
306  forAll(faceMap, i)
307  {
308  const triFace& f = s[faceMap[i]];
309  triFace newF
310  (
311  reversePointMap[f[0]],
312  reversePointMap[f[1]],
313  reversePointMap[f[2]]
314  );
315  faces[i] = newF.triFaceFace();
316 
317  forAll(newF, fp)
318  {
319  pointToFace[newF[fp]] = i;
320  }
321  }
322 
323  this->storedPoints() = pointField(s.points(), pointMap);
324 
325  if (debug)
326  {
327  print(Pout);
328  Pout<< endl;
329  }
330 
331 
332 
333  // Collect the samplePoints and sampleElements
334  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
335 
337  {
338  samplePoints_.setSize(pointMap.size());
339  sampleElements_.setSize(pointMap.size(), -1);
340 
341  if (sampleSource_ == cells)
342  {
343  // samplePoints_ : per surface point a location inside the cell
344  // sampleElements_ : per surface point the cell
345 
346  forAll(points(), pointi)
347  {
348  const point& pt = points()[pointi];
349  label celli = cellOrFaceLabels[pointToFace[pointi]];
350  sampleElements_[pointi] = celli;
351 
352  // Check if point inside cell
353  if
354  (
355  mesh().pointInCell
356  (
357  pt,
358  sampleElements_[pointi],
359  meshSearcher.decompMode()
360  )
361  )
362  {
363  samplePoints_[pointi] = pt;
364  }
365  else
366  {
367  // Find nearest point on faces of cell
368  const cell& cFaces = mesh().cells()[celli];
369 
370  scalar minDistSqr = vGreat;
371 
372  forAll(cFaces, i)
373  {
374  const face& f = mesh().faces()[cFaces[i]];
375  pointHit info = f.nearestPoint(pt, mesh().points());
376  if (info.distance() < minDistSqr)
377  {
378  minDistSqr = info.distance();
379  samplePoints_[pointi] = info.rawPoint();
380  }
381  }
382  }
383  }
384  }
385  else if (sampleSource_ == insideCells)
386  {
387  // samplePoints_ : per surface point a location inside the cell
388  // sampleElements_ : per surface point the cell
389 
390  forAll(points(), pointi)
391  {
392  const point& pt = points()[pointi];
393  label celli = cellOrFaceLabels[pointToFace[pointi]];
394  sampleElements_[pointi] = celli;
395  samplePoints_[pointi] = pt;
396  }
397  }
398  else
399  {
400  // samplePoints_ : per surface point a location on the boundary
401  // sampleElements_ : per surface point the boundary face containing
402  // the location
403 
404  forAll(points(), pointi)
405  {
406  const point& pt = points()[pointi];
407  label facei = cellOrFaceLabels[pointToFace[pointi]];
408  sampleElements_[pointi] = facei;
409  samplePoints_[pointi] = mesh().faces()[facei].nearestPoint
410  (
411  pt,
412  mesh().points()
413  ).rawPoint();
414  }
415  }
416  }
417  else
418  {
419  // if sampleSource_ == cells:
420  // samplePoints_ : n/a
421  // sampleElements_ : per surface triangle the cell
422  // if sampleSource_ == insideCells:
423  // samplePoints_ : n/a
424  // sampleElements_ : -1 or per surface triangle the cell
425  // else:
426  // samplePoints_ : n/a
427  // sampleElements_ : per surface triangle the boundary face
428  samplePoints_.clear();
429  sampleElements_.transfer(cellOrFaceLabels);
430  }
431 
432 
433  if (debug)
434  {
435  this->clearOut();
436 
437  OBJstream str(mesh().time().path()/"surfaceToMesh.obj");
438 
439  Info<< "Dumping correspondence from local surface (points or faces)"
440  << " to mesh (cells or faces) to " << str.name() << endl;
441 
443  {
444  if (sampleSource_ == cells || sampleSource_ == insideCells)
445  {
446  forAll(samplePoints_, pointi)
447  {
448  const label celli = sampleElements_[pointi];
449 
450  str.write
451  (
453  (
454  points()[pointi],
455  samplePoints_[pointi],
456  mesh().cellCentres()[celli]
457  )
458  );
459  }
460  }
461  else
462  {
463  forAll(samplePoints_, pointi)
464  {
465  const label facei = sampleElements_[pointi];
466 
467  str.write
468  (
470  (
471  points()[pointi],
472  samplePoints_[pointi],
473  mesh().faceCentres()[facei]
474  )
475  );
476  }
477  }
478  }
479  else
480  {
481  if (sampleSource_ == cells || sampleSource_ == insideCells)
482  {
483  forAll(sampleElements_, triI)
484  {
485  const label celli = sampleElements_[triI];
486 
487  str.write
488  (
490  (
491  faceCentres()[triI],
492  mesh().cellCentres()[celli]
493  )
494  );
495  }
496  }
497  else
498  {
499  forAll(sampleElements_, triI)
500  {
501  const label facei = sampleElements_[triI];
502 
503  str.write
504  (
506  (
507  faceCentres()[triI],
508  mesh().faceCentres()[facei]
509  )
510  );
511  }
512  }
513  }
514  }
515 
516  needsUpdate_ = false;
517 
518  return true;
519 }
520 
521 
522 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
523 
525 (
526  const word& name,
527  const polyMesh& mesh,
528  const word& surfaceName,
529  const samplingSource sampleSource
530 )
531 :
532  sampledSurface(name, mesh),
533  surface_
534  (
535  IOobject
536  (
537  surfaceName,
538  mesh.time().constant(),
539  searchableSurface::geometryDir(mesh.time()),
540  mesh,
541  IOobject::MUST_READ,
542  IOobject::NO_WRITE,
543  false
544  )
545  ),
546  sampleSource_(sampleSource),
547  needsUpdate_(true),
548  sampleElements_(0),
549  samplePoints_(0)
550 {}
551 
552 
554 (
555  const word& name,
556  const polyMesh& mesh,
557  const dictionary& dict
558 )
559 :
560  sampledSurface(name, mesh, dict),
561  surface_
562  (
563  IOobject
564  (
565  dict.lookup("surface"),
566  mesh.time().constant(),
567  searchableSurface::geometryDir(mesh.time()),
568  mesh,
569  IOobject::MUST_READ,
570  IOobject::NO_WRITE,
571  false
572  )
573  ),
574  sampleSource_(samplingSourceNames_[dict.lookup("source")]),
575  needsUpdate_(true),
576  sampleElements_(0),
577  samplePoints_(0)
578 {}
579 
580 
582 (
583  const word& name,
584  const polyMesh& mesh,
585  const triSurface& surface,
586  const word& sampleSourceName
587 )
588 :
589  sampledSurface(name, mesh),
590  surface_
591  (
592  IOobject
593  (
594  name,
595  mesh.time().constant(),
596  searchableSurface::geometryDir(mesh.time()),
597  mesh,
598  IOobject::NO_READ,
599  IOobject::NO_WRITE,
600  false
601  ),
602  surface
603  ),
604  sampleSource_(samplingSourceNames_[sampleSourceName]),
605  needsUpdate_(true),
606  sampleElements_(0),
607  samplePoints_(0)
608 {}
609 
610 
611 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
612 
614 {}
615 
616 
617 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
618 
620 {
621  return needsUpdate_;
622 }
623 
624 
626 {
627  // already marked as expired
628  if (needsUpdate_)
629  {
630  return false;
631  }
632 
635 
636  boundaryTreePtr_.clear();
637  sampleElements_.clear();
638  samplePoints_.clear();
639 
640  needsUpdate_ = true;
641  return true;
642 }
643 
644 
646 {
647  if (!needsUpdate_)
648  {
649  return false;
650  }
651 
652  // Mesh search engine, no triangulation of faces.
653  meshSearch meshSearcher(mesh(), polyMesh::FACE_PLANES);
654 
655  return update(meshSearcher);
656 }
657 
658 
660 {
661  if (!needsUpdate_)
662  {
663  return false;
664  }
665 
666  // Mesh search engine on subset, no triangulation of faces.
667  meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
668 
669  return update(meshSearcher);
670 }
671 
672 
675 (
676  const volScalarField& vField
677 ) const
678 {
679  return sampleField(vField);
680 }
681 
682 
685 (
686  const volVectorField& vField
687 ) const
688 {
689  return sampleField(vField);
690 }
691 
694 (
695  const volSphericalTensorField& vField
696 ) const
697 {
698  return sampleField(vField);
699 }
700 
701 
704 (
705  const volSymmTensorField& vField
706 ) const
707 {
708  return sampleField(vField);
709 }
710 
711 
714 (
715  const volTensorField& vField
716 ) const
717 {
718  return sampleField(vField);
719 }
720 
721 
724 (
725  const interpolation<scalar>& interpolator
726 ) const
727 {
728  return interpolateField(interpolator);
729 }
730 
731 
734 (
735  const interpolation<vector>& interpolator
736 ) const
737 {
738  return interpolateField(interpolator);
739 }
740 
743 (
744  const interpolation<sphericalTensor>& interpolator
745 ) const
746 {
747  return interpolateField(interpolator);
748 }
749 
750 
753 (
754  const interpolation<symmTensor>& interpolator
755 ) const
756 {
757  return interpolateField(interpolator);
758 }
759 
760 
763 (
764  const interpolation<tensor>& interpolator
765 ) const
766 {
767  return interpolateField(interpolator);
768 }
769 
770 
772 {
773  os << "triSurfaceMesh: " << name() << " :"
774  << " surface:" << surface_.objectRegistry::name()
775  << " faces:" << faces().size()
776  << " points:" << points().size();
777 }
778 
779 
780 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
std::remove_reference< ::Foam::List< labelledTri > >::type::value_type FaceType
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for interpolation.
Definition: interpolation.H:55
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:58
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:290
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
An abstract class for surfaces with sampling.
virtual void clearGeom() const
bool interpolate() const
Interpolation requested for surface.
const polyMesh & mesh() const
Access to the underlying mesh.
void operator()(nearInfo &x, const nearInfo &y) const
A sampledSurface from a triSurfaceMesh. It samples on the points/triangles of the triSurface....
triSurfaceMesh(const word &name, const polyMesh &mesh, const word &surfaceName, const samplingSource sampleSource)
Construct from components.
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
virtual bool expire()
Mark the surface as needing an update.
virtual bool needsUpdate() const
Does the surface need an update?
virtual bool update()
Update the surface as required.
virtual void print(Ostream &) const
Write.
virtual const pointField & points() const
Points of surface.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
A class for managing temporary objects.
Definition: tmp.H:55
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
Triangulated surface description with patch information.
Definition: triSurface.H:69
A class for handling words, derived from string.
Definition: word.H:62
label patchi
const pointField & points
const cellShapeList & cells
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const fvPatchList & patches
tUEqn clear()
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
defineTypeNameAndDebug(cutPlane, 0)
addToRunTimeSelectionTable(sampledSurface, cutPlane, word)
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
vector point
Point is a vector.
Definition: point.H:41
PointHit< point > pointHit
Definition: pointHit.H:41
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
static const label labelMax
Definition: label.H:62
List< face > faceList
Definition: faceListFwd.H:41
dimensioned< scalar > magSqr(const dimensioned< Type > &)
label newPointi
Definition: readKivaGrid.H:495
face triFace(3)
labelList f(nPoints)
dictionary dict