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