sampledTriSurfaceMesh.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-2015 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 "Tuple2.H"
29 #include "globalIndex.H"
30 #include "treeDataCell.H"
31 #include "treeDataFace.H"
32 #include "meshTools.H"
33 
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(sampledTriSurfaceMesh, 0);
42  (
43  sampledSurface,
44  sampledTriSurfaceMesh,
45  word
46  );
47 
48  template<>
50  {
51  "cells",
52  "insideCells",
53  "boundaryFaces"
54  };
55 
58 
59 
60  //- Private class for finding nearest
61  // Comprising:
62  // - global index
63  // - sqr(distance)
65 
67  {
68 
69  public:
70 
71  void operator()(nearInfo& x, const nearInfo& y) const
72  {
73  if (y.first() < x.first())
74  {
75  x = y;
76  }
77  }
78  };
79 }
80 
81 
82 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
83 
86 {
87  // Variant of meshSearch::boundaryTree() that only does non-coupled
88  // boundary faces.
89 
90  if (!boundaryTreePtr_.valid())
91  {
92  // all non-coupled boundary faces (not just walls)
94 
95  labelList bndFaces(mesh().nFaces()-mesh().nInternalFaces());
96  label bndI = 0;
97  forAll(patches, patchI)
98  {
99  const polyPatch& pp = patches[patchI];
100  if (!pp.coupled())
101  {
102  forAll(pp, i)
103  {
104  bndFaces[bndI++] = pp.start()+i;
105  }
106  }
107  }
108  bndFaces.setSize(bndI);
109 
110 
111  treeBoundBox overallBb(mesh().points());
112  Random rndGen(123456);
113  overallBb = overallBb.extend(rndGen, 1e-4);
114  overallBb.min() -= point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
115  overallBb.max() += point(ROOTVSMALL, ROOTVSMALL, ROOTVSMALL);
116 
117  boundaryTreePtr_.reset
118  (
120  (
121  treeDataFace // all information needed to search faces
122  (
123  false, // do not cache bb
124  mesh(),
125  bndFaces // boundary faces only
126  ),
127  overallBb, // overall search domain
128  8, // maxLevel
129  10, // leafsize
130  3.0 // duplicity
131  )
132  );
133  }
134 
135  return boundaryTreePtr_();
136 }
137 
138 
140 {
141  // Find the cells the triangles of the surface are in.
142  // Does approximation by looking at the face centres only
143  const pointField& fc = surface_.faceCentres();
144 
145  List<nearInfo> nearest(fc.size());
146 
147  // Global numbering for cells/faces - only used to uniquely identify local
148  // elements
149  globalIndex globalCells
150  (
151  (sampleSource_ == cells || sampleSource_ == insideCells)
152  ? mesh().nCells()
153  : mesh().nFaces()
154  );
155 
156  forAll(nearest, i)
157  {
158  nearest[i].first() = GREAT;
159  nearest[i].second() = labelMax;
160  }
161 
162  if (sampleSource_ == cells)
163  {
164  // Search for nearest cell
165 
166  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
167 
168  forAll(fc, triI)
169  {
171  (
172  fc[triI],
173  sqr(GREAT)
174  );
175  if (nearInfo.hit())
176  {
177  nearest[triI].first() = magSqr(nearInfo.hitPoint()-fc[triI]);
178  nearest[triI].second() = globalCells.toGlobal(nearInfo.index());
179  }
180  }
181  }
182  else if (sampleSource_ == insideCells)
183  {
184  // Search for cell containing point
185 
186  const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
187 
188  forAll(fc, triI)
189  {
190  if (cellTree.bb().contains(fc[triI]))
191  {
192  label index = cellTree.findInside(fc[triI]);
193  if (index != -1)
194  {
195  nearest[triI].first() = 0.0;
196  nearest[triI].second() = globalCells.toGlobal(index);
197  }
198  }
199  }
200  }
201  else
202  {
203  // Search for nearest boundaryFace
204 
206  //const indexedOctree<treeDataFace>& bTree = meshSearcher.boundaryTree()
207  //- Search on all non-coupled boundary faces
208  const indexedOctree<treeDataFace>& bTree = nonCoupledboundaryTree();
209 
210  forAll(fc, triI)
211  {
213  (
214  fc[triI],
215  sqr(GREAT)
216  );
217  if (nearInfo.hit())
218  {
219  nearest[triI].first() = magSqr(nearInfo.hitPoint()-fc[triI]);
220  nearest[triI].second() = globalCells.toGlobal
221  (
222  bTree.shapes().faceLabels()[nearInfo.index()]
223  );
224  }
225  }
226  }
227 
228 
229  // See which processor has the nearest. Mark and subset
230  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
231 
234 
235  labelList cellOrFaceLabels(fc.size(), -1);
236 
237  label nFound = 0;
238  forAll(nearest, triI)
239  {
240  if (nearest[triI].second() == labelMax)
241  {
242  // Not found on any processor. How to map?
243  }
244  else if (globalCells.isLocal(nearest[triI].second()))
245  {
246  cellOrFaceLabels[triI] = globalCells.toLocal
247  (
248  nearest[triI].second()
249  );
250  nFound++;
251  }
252  }
253 
254 
255  if (debug)
256  {
257  Pout<< "Local out of faces:" << cellOrFaceLabels.size()
258  << " keeping:" << nFound << endl;
259  }
260 
261  // Now subset the surface. Do not use triSurface::subsetMesh since requires
262  // original surface to be in compact numbering.
263 
264  const triSurface& s = surface_;
265 
266  // Compact to original triangle
267  labelList faceMap(s.size());
268  // Compact to original points
269  labelList pointMap(s.points().size());
270  // From original point to compact points
271  labelList reversePointMap(s.points().size(), -1);
272 
273  {
274  label newPointI = 0;
275  label newFaceI = 0;
276 
277  forAll(s, faceI)
278  {
279  if (cellOrFaceLabels[faceI] != -1)
280  {
281  faceMap[newFaceI++] = faceI;
282 
283  const triSurface::FaceType& f = s[faceI];
284  forAll(f, fp)
285  {
286  if (reversePointMap[f[fp]] == -1)
287  {
288  pointMap[newPointI] = f[fp];
289  reversePointMap[f[fp]] = newPointI++;
290  }
291  }
292  }
293  }
294  faceMap.setSize(newFaceI);
295  pointMap.setSize(newPointI);
296  }
297 
298 
299  // Subset cellOrFaceLabels
300  cellOrFaceLabels = UIndirectList<label>(cellOrFaceLabels, faceMap)();
301 
302  // Store any face per point (without using pointFaces())
303  labelList pointToFace(pointMap.size());
304 
305  // Create faces and points for subsetted surface
306  faceList& faces = this->storedFaces();
307  faces.setSize(faceMap.size());
308  forAll(faceMap, i)
309  {
310  const triFace& f = s[faceMap[i]];
311  triFace newF
312  (
313  reversePointMap[f[0]],
314  reversePointMap[f[1]],
315  reversePointMap[f[2]]
316  );
317  faces[i] = newF.triFaceFace();
318 
319  forAll(newF, fp)
320  {
321  pointToFace[newF[fp]] = i;
322  }
323  }
324 
325  this->storedPoints() = pointField(s.points(), pointMap);
326 
327  if (debug)
328  {
329  print(Pout);
330  Pout<< endl;
331  }
332 
333 
334 
335  // Collect the samplePoints and sampleElements
336  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
337 
339  {
340  samplePoints_.setSize(pointMap.size());
341  sampleElements_.setSize(pointMap.size(), -1);
342 
343  if (sampleSource_ == cells)
344  {
345  // samplePoints_ : per surface point a location inside the cell
346  // sampleElements_ : per surface point the cell
347 
348  forAll(points(), pointI)
349  {
350  const point& pt = points()[pointI];
351  label cellI = cellOrFaceLabels[pointToFace[pointI]];
352  sampleElements_[pointI] = cellI;
353 
354  // Check if point inside cell
355  if
356  (
357  mesh().pointInCell
358  (
359  pt,
360  sampleElements_[pointI],
361  meshSearcher.decompMode()
362  )
363  )
364  {
365  samplePoints_[pointI] = pt;
366  }
367  else
368  {
369  // Find nearest point on faces of cell
370  const cell& cFaces = mesh().cells()[cellI];
371 
372  scalar minDistSqr = VGREAT;
373 
374  forAll(cFaces, i)
375  {
376  const face& f = mesh().faces()[cFaces[i]];
377  pointHit info = f.nearestPoint(pt, mesh().points());
378  if (info.distance() < minDistSqr)
379  {
380  minDistSqr = info.distance();
381  samplePoints_[pointI] = info.rawPoint();
382  }
383  }
384  }
385  }
386  }
387  else if (sampleSource_ == insideCells)
388  {
389  // samplePoints_ : per surface point a location inside the cell
390  // sampleElements_ : per surface point the cell
391 
392  forAll(points(), pointI)
393  {
394  const point& pt = points()[pointI];
395  label cellI = cellOrFaceLabels[pointToFace[pointI]];
396  sampleElements_[pointI] = cellI;
397  samplePoints_[pointI] = pt;
398  }
399  }
400  else
401  {
402  // samplePoints_ : per surface point a location on the boundary
403  // sampleElements_ : per surface point the boundary face containing
404  // the location
405 
406  forAll(points(), pointI)
407  {
408  const point& pt = points()[pointI];
409  label faceI = cellOrFaceLabels[pointToFace[pointI]];
410  sampleElements_[pointI] = faceI;
411  samplePoints_[pointI] = mesh().faces()[faceI].nearestPoint
412  (
413  pt,
414  mesh().points()
415  ).rawPoint();
416  }
417  }
418  }
419  else
420  {
421  // if sampleSource_ == cells:
422  // samplePoints_ : n/a
423  // sampleElements_ : per surface triangle the cell
424  // if sampleSource_ == insideCells:
425  // samplePoints_ : n/a
426  // sampleElements_ : -1 or per surface triangle the cell
427  // else:
428  // samplePoints_ : n/a
429  // sampleElements_ : per surface triangle the boundary face
430  samplePoints_.clear();
431  sampleElements_.transfer(cellOrFaceLabels);
432  }
433 
434 
435  if (debug)
436  {
437  this->clearOut();
438  OFstream str(mesh().time().path()/"surfaceToMesh.obj");
439  Info<< "Dumping correspondence from local surface (points or faces)"
440  << " to mesh (cells or faces) to " << str.name() << endl;
441  label vertI = 0;
442 
444  {
445  if (sampleSource_ == cells || sampleSource_ == insideCells)
446  {
447  forAll(samplePoints_, pointI)
448  {
449  meshTools::writeOBJ(str, points()[pointI]);
450  vertI++;
451 
452  meshTools::writeOBJ(str, samplePoints_[pointI]);
453  vertI++;
454 
455  label cellI = sampleElements_[pointI];
456  meshTools::writeOBJ(str, mesh().cellCentres()[cellI]);
457  vertI++;
458  str << "l " << vertI-2 << ' ' << vertI-1 << ' ' << vertI
459  << nl;
460  }
461  }
462  else
463  {
464  forAll(samplePoints_, pointI)
465  {
466  meshTools::writeOBJ(str, points()[pointI]);
467  vertI++;
468 
469  meshTools::writeOBJ(str, samplePoints_[pointI]);
470  vertI++;
471 
472  label faceI = sampleElements_[pointI];
473  meshTools::writeOBJ(str, mesh().faceCentres()[faceI]);
474  vertI++;
475  str << "l " << vertI-2 << ' ' << vertI-1 << ' ' << vertI
476  << nl;
477  }
478  }
479  }
480  else
481  {
482  if (sampleSource_ == cells || sampleSource_ == insideCells)
483  {
484  forAll(sampleElements_, triI)
485  {
486  meshTools::writeOBJ(str, faceCentres()[triI]);
487  vertI++;
488 
489  label cellI = sampleElements_[triI];
490  meshTools::writeOBJ(str, mesh().cellCentres()[cellI]);
491  vertI++;
492  str << "l " << vertI-1 << ' ' << vertI << nl;
493  }
494  }
495  else
496  {
497  forAll(sampleElements_, triI)
498  {
499  meshTools::writeOBJ(str, faceCentres()[triI]);
500  vertI++;
501 
502  label faceI = sampleElements_[triI];
503  meshTools::writeOBJ(str, mesh().faceCentres()[faceI]);
504  vertI++;
505  str << "l " << vertI-1 << ' ' << vertI << nl;
506  }
507  }
508  }
509  }
510 
511 
512  needsUpdate_ = false;
513  return true;
514 }
515 
516 
517 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
518 
520 (
521  const word& name,
522  const polyMesh& mesh,
523  const word& surfaceName,
524  const samplingSource sampleSource
525 )
526 :
528  surface_
529  (
530  IOobject
531  (
532  surfaceName,
533  mesh.time().constant(), // instance
534  "triSurface", // local
535  mesh, // registry
538  false
539  )
540  ),
541  sampleSource_(sampleSource),
542  needsUpdate_(true),
543  sampleElements_(0),
544  samplePoints_(0)
545 {}
546 
547 
549 (
550  const word& name,
551  const polyMesh& mesh,
552  const dictionary& dict
553 )
554 :
556  surface_
557  (
558  IOobject
559  (
560  dict.lookup("surface"),
561  mesh.time().constant(), // instance
562  "triSurface", // local
563  mesh, // registry
566  false
567  )
568  ),
569  sampleSource_(samplingSourceNames_[dict.lookup("source")]),
570  needsUpdate_(true),
571  sampleElements_(0),
572  samplePoints_(0)
573 {}
574 
575 
577 (
578  const word& name,
579  const polyMesh& mesh,
580  const triSurface& surface,
581  const word& sampleSourceName
582 )
583 :
585  surface_
586  (
587  IOobject
588  (
589  name,
590  mesh.time().constant(), // instance
591  "triSurface", // local
592  mesh, // registry
595  false
596  ),
597  surface
598  ),
599  sampleSource_(samplingSourceNames_[sampleSourceName]),
600  needsUpdate_(true),
601  sampleElements_(0),
602  samplePoints_(0)
603 {}
604 
605 
606 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
607 
609 {}
610 
611 
612 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
613 
615 {
616  return needsUpdate_;
617 }
618 
619 
621 {
622  // already marked as expired
623  if (needsUpdate_)
624  {
625  return false;
626  }
627 
630 
631  boundaryTreePtr_.clear();
632  sampleElements_.clear();
633  samplePoints_.clear();
634 
635  needsUpdate_ = true;
636  return true;
637 }
638 
639 
641 {
642  if (!needsUpdate_)
643  {
644  return false;
645  }
646 
647  Info <<endl << "current STL model name is : "<<surface_.objectRegistry::name() << endl<<endl;
648  // Calculate surface and mesh overlap bounding box
649  treeBoundBox bb
650  (
651  surface_.triSurface::points(),
652  surface_.triSurface::meshPoints()
653  );
654  bb.min() = max(bb.min(), mesh().bounds().min());
655  bb.max() = min(bb.max(), mesh().bounds().max());
656 
657  // Extend a bit
658  const vector span(bb.span());
659 
660  bb.min() -= 0.5*span;
661  bb.max() += 0.5*span;
662 
663  bb.inflate(1e-6);
664 
665  // Mesh search engine, no triangulation of faces.
666  meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
667 
668  return update(meshSearcher);
669 }
670 
671 
673 {
674  if (!needsUpdate_)
675  {
676  return false;
677  }
678 
679  // Mesh search engine on subset, no triangulation of faces.
680  meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
681 
682  return update(meshSearcher);
683 }
684 
685 
687 (
688  const volScalarField& vField
689 ) const
690 {
691  return sampleField(vField);
692 }
693 
694 
696 (
697  const volVectorField& vField
698 ) const
699 {
700  return sampleField(vField);
701 }
702 
704 (
705  const volSphericalTensorField& vField
706 ) const
707 {
708  return sampleField(vField);
709 }
710 
711 
713 (
714  const volSymmTensorField& vField
715 ) const
716 {
717  return sampleField(vField);
718 }
719 
720 
722 (
723  const volTensorField& vField
724 ) const
725 {
726  return sampleField(vField);
727 }
728 
729 
731 (
732  const interpolation<scalar>& interpolator
733 ) const
734 {
735  return interpolateField(interpolator);
736 }
737 
738 
740 (
741  const interpolation<vector>& interpolator
742 ) const
743 {
744  return interpolateField(interpolator);
745 }
746 
748 (
749  const interpolation<sphericalTensor>& interpolator
750 ) const
751 {
752  return interpolateField(interpolator);
753 }
754 
755 
757 (
758  const interpolation<symmTensor>& interpolator
759 ) const
760 {
761  return interpolateField(interpolator);
762 }
763 
764 
766 (
767  const interpolation<tensor>& interpolator
768 ) const
769 {
770  return interpolateField(interpolator);
771 }
772 
773 
775 {
776  os << "sampledTriSurfaceMesh: " << name() << " :"
777  << " surface:" << surface_.objectRegistry::name()
778  << " faces:" << faces().size()
779  << " points:" << points().size();
780 }
781 
782 
783 // ************************************************************************* //
Foam::sampledTriSurfaceMesh::points
virtual const pointField & points() const
Points of surface.
Definition: sampledTriSurfaceMesh.H:207
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::Random
Simple random number generator.
Definition: Random.H:49
Foam::sampledTriSurfaceMesh::print
virtual void print(Ostream &) const
Write.
Definition: sampledTriSurfaceMesh.C:774
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
meshTools.H
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:57
print
void print(const char *msg, Ostream &os, const PtrList< GeoField > &flds)
Definition: foamToVTK.C:164
Foam::nearestEqOp
Definition: sampledTriSurfaceMesh.C:66
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
Foam::sampledTriSurfaceMesh::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledTriSurfaceMesh.C:614
Foam::sampledTriSurfaceMesh::sampledTriSurfaceMesh
sampledTriSurfaceMesh(const word &name, const polyMesh &mesh, const word &surfaceName, const samplingSource sampleSource)
Construct from components.
Definition: sampledTriSurfaceMesh.C:520
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
clear
UEqn clear()
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::sampledTriSurfaceMesh::boundaryTreePtr_
autoPtr< indexedOctree< treeDataFace > > boundaryTreePtr_
Search tree for all non-coupled boundary faces.
Definition: sampledTriSurfaceMesh.H:120
Foam::sampledTriSurfaceMesh::update
virtual bool update()
Update the surface as required.
Definition: sampledTriSurfaceMesh.C:640
Tuple2.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::PointHit
This class describes the interaction of a face and a point. It carries the info of a successful hit a...
Definition: PointHit.H:51
Foam::sampledTriSurfaceMesh::samplingSource
samplingSource
Types of communications.
Definition: sampledTriSurfaceMesh.H:93
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::boundBox::inflate
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Definition: boundBox.C:210
globalIndex.H
Foam::labelMax
static const label labelMax
Definition: label.H:62
Foam::PointHit::rawPoint
const Point & rawPoint() const
Return point with no checking.
Definition: PointHit.H:158
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:322
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
Foam::indexedOctree::shapes
const Type & shapes() const
Reference to shape.
Definition: indexedOctree.H:449
Foam::triFace::triFaceFace
face triFaceFace() const
Return triangle as a face.
Definition: triFaceI.H:140
Foam::sampledTriSurfaceMesh::samplingSourceNames_
static const NamedEnum< samplingSource, 3 > samplingSourceNames_
Definition: sampledTriSurfaceMesh.H:108
Foam::globalIndex::isLocal
bool isLocal(const label i) const
Is on local processor.
Definition: globalIndexI.H:95
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::nearInfo
Tuple2< scalar, label > nearInfo
Private class for finding nearest.
Definition: sampledTriSurfaceMesh.C:64
sampledTriSurfaceMesh.H
treeDataFace.H
Foam::boundBox::span
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
Foam::meshSearch::cellTree
const indexedOctree< treeDataCell > & cellTree() const
Get (demand driven) reference to octree holding all cells.
Definition: meshSearch.C:605
Foam::meshSearch::decompMode
polyMesh::cellDecomposition decompMode() const
Definition: meshSearch.H:199
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::PointHit::distance
scalar distance() const
Return distance to hit.
Definition: PointHit.H:139
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
Foam::treeBoundBox::extend
treeBoundBox extend(Random &, const scalar s) const
Return slightly wider bounding box.
Definition: treeBoundBoxI.H:317
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::indexedOctree< Foam::treeDataFace >
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:77
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::indexedOctree::findInside
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
Definition: indexedOctree.C:2828
Foam::interpolation< scalar >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::sampledTriSurfaceMesh::sample
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledTriSurfaceMesh.C:687
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::sampledTriSurfaceMesh::nonCoupledboundaryTree
const indexedOctree< treeDataFace > & nonCoupledboundaryTree() const
Get tree of all non-coupled boundary faces.
Definition: sampledTriSurfaceMesh.C:85
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
treeDataCell.H
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
s
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::indexedOctree::bb
const treeBoundBox & bb() const
Top bounding box.
Definition: indexedOctree.H:468
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:282
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::globalIndex::toLocal
label toLocal(const label i) const
From global to local on current processor.
Definition: globalIndexI.H:117
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::treeBoundBox::contains
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:395
meshSearch.H
Foam::triFace
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:68
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::treeDataFace
Encapsulation of data needed to search for faces.
Definition: treeDataFace.H:58
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePaths.H:130
Foam::pointToFace
A topoSetSource to select faces based on use of points.
Definition: pointToFace.H:49
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:379
Foam::surface
Definition: surface.H:55
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::sampledSurface::clearGeom
virtual void clearGeom() const
Definition: sampledSurface.C:42
Foam::nearestEqOp::operator()
void operator()(nearInfo &x, const nearInfo &y) const
Definition: sampledTriSurfaceMesh.C:71
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
Foam::sampledTriSurfaceMesh::~sampledTriSurfaceMesh
virtual ~sampledTriSurfaceMesh()
Destructor.
Definition: sampledTriSurfaceMesh.C:608
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::indexedOctree::findNearest
void findNearest(const label nodeI, const linePointRef &ln, treeBoundBox &tightest, label &nearestShapeI, point &linePoint, point &nearestPoint, const FindNearestOp &fnOp) const
Find nearest point to line.
Definition: indexedOctree.C:590
Foam::polyMesh::FACE_PLANES
@ FACE_PLANES
Definition: polyMesh.H:100
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::sampledSurface::mesh
const polyMesh & mesh() const
Access to the underlying mesh.
Definition: sampledSurface.H:244
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::OFstream::name
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:118
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::Tuple2
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:47
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:417
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::sampledTriSurfaceMesh::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledTriSurfaceMesh.C:620
Foam::sampledSurface::interpolate
bool interpolate() const
Interpolation requested for surface.
Definition: sampledSurface.H:256
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
rndGen
cachedRandom rndGen(label(0), -1)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
y
scalar y
Definition: LISASMDCalcMethod1.H:14