sampledCuttingPlane.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "sampledCuttingPlane.H"
30 #include "dictionary.H"
31 #include "fvMesh.H"
32 #include "volFields.H"
33 #include "volPointInterpolation.H"
34 #include "cartesianCS.H"
36 #include "PtrList.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(sampledCuttingPlane, 0);
44  (
45  sampledSurface,
46  sampledCuttingPlane,
47  word,
48  cuttingPlane
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 Foam::plane Foam::sampledCuttingPlane::definePlane
56 (
57  const polyMesh& mesh,
58  const dictionary& dict
59 )
60 {
61  plane pln(dict);
62 
63  // Optional (cartesian) coordinate transform.
64  // - with registry to allow lookup from globally defined systems
65 
66  auto csysPtr = coordinateSystem::NewIfPresent(mesh, dict);
67 
68  if (!csysPtr)
69  {
70  csysPtr = coordinateSystem::NewIfPresent(dict, "transform");
71  }
72 
73  // Make plane relative to the Cartesian coordinate system
74  if (csysPtr)
75  {
76  coordSystem::cartesian cs(csysPtr());
77 
78  const point orig = cs.globalPosition(pln.origin());
79  const vector norm = cs.globalVector(pln.normal());
80 
81  DebugInfo
82  << "plane "
83  << " origin:" << pln.origin()
84  << " normal:" << pln.normal()
85  << " =>"
86  << " origin:" << orig << " normal:" << norm
87  << endl;
88 
89  // Reassign the plane
90  pln = plane(orig, norm);
91  }
92 
93  return pln;
94 }
95 
96 
97 void Foam::sampledCuttingPlane::checkBoundsIntersection
98 (
99  const plane& pln,
100  const boundBox& meshBb
101 ) const
102 {
103  // Verify specified bounding box
104  const boundBox& clipBb = isoParams_.getClipBounds();
105 
106  if (clipBb.good())
107  {
108  // Bounding box does not overlap with (global) mesh!
109  if (!clipBb.overlaps(meshBb))
110  {
112  << nl
113  << name() << " : "
114  << "Bounds " << clipBb
115  << " do not overlap the mesh bounding box " << meshBb
116  << nl << endl;
117  }
118 
119  // Plane does not intersect the bounding box
120  if (!clipBb.intersects(pln))
121  {
123  << nl
124  << name() << " : "
125  << "Plane "<< pln << " does not intersect the bounds "
126  << clipBb
127  << nl << endl;
128  }
129  }
130 
131  // Plane does not intersect the (global) mesh!
132  if (!meshBb.intersects(pln))
133  {
135  << nl
136  << name() << " : "
137  << "Plane "<< pln << " does not intersect the mesh bounds "
138  << meshBb
139  << nl << endl;
140  }
141 }
142 
143 
144 void Foam::sampledCuttingPlane::setDistanceFields(const plane& pln)
145 {
146  volScalarField& cellDistance = cellDistancePtr_();
147 
148  // Get mesh from volField,
149  // so automatically submesh or baseMesh
150 
151  const fvMesh& mesh = cellDistance.mesh();
152 
153  // Distance to cell centres
154  // ~~~~~~~~~~~~~~~~~~~~~~~~
155 
156  // Internal field
157  {
158  const auto& cc = mesh.cellCentres();
159  scalarField& fld = cellDistance.primitiveFieldRef();
160 
161  forAll(cc, i)
162  {
163  fld[i] = pln.signedDistance(cc[i]);
164  }
165  }
166 
167  // Patch fields
168  {
169  volScalarField::Boundary& cellDistanceBf =
170  cellDistance.boundaryFieldRef();
171 
172  forAll(cellDistanceBf, patchi)
173  {
174  if
175  (
176  isA<emptyFvPatchScalarField>
177  (
178  cellDistanceBf[patchi]
179  )
180  )
181  {
182  cellDistanceBf.set
183  (
184  patchi,
185  new calculatedFvPatchScalarField
186  (
187  mesh.boundary()[patchi],
188  cellDistance
189  )
190  );
191 
192  const polyPatch& pp = mesh.boundary()[patchi].patch();
193  pointField::subField cc = pp.patchSlice(mesh.faceCentres());
194 
195  fvPatchScalarField& fld = cellDistanceBf[patchi];
196  fld.setSize(pp.size());
197  forAll(fld, i)
198  {
199  fld[i] = pln.signedDistance(cc[i]);
200  }
201  }
202  else
203  {
204  // Other side cell centres?
205  const pointField& cc = mesh.C().boundaryField()[patchi];
206  fvPatchScalarField& fld = cellDistanceBf[patchi];
207 
208  forAll(fld, i)
209  {
210  fld[i] = pln.signedDistance(cc[i]);
211  }
212  }
213  }
214  }
215 
216 
217  // On processor patches the mesh.C() will already be the cell centre
218  // on the opposite side so no need to swap cellDistance.
219 
220  // Distance to points
221  pointDistance_.resize(mesh.nPoints());
222  {
223  const pointField& pts = mesh.points();
224 
225  forAll(pointDistance_, i)
226  {
227  pointDistance_[i] = pln.signedDistance(pts[i]);
228  }
229  }
230 }
231 
232 
233 void Foam::sampledCuttingPlane::combineSurfaces
234 (
235  PtrList<isoSurfaceBase>& isoSurfPtrs
236 )
237 {
238  isoSurfacePtr_.reset(nullptr);
239 
240  // Already checked previously for ALGO_POINT, but do it again
241  // - ALGO_POINT still needs fields (for interpolate)
242  // The others can do straight transfer
243  if
244  (
245  isoParams_.algorithm() == isoSurfaceParams::ALGO_POINT
246  && isoSurfPtrs.size() == 1
247  )
248  {
249  // Shift ownership from list to autoPtr
250  isoSurfacePtr_.reset(isoSurfPtrs.release(0));
251  }
252  else if (isoSurfPtrs.size() == 1)
253  {
254  autoPtr<isoSurfaceBase> ptr(isoSurfPtrs.release(0));
255  auto& surf = *ptr;
256 
257  surface_.transfer(static_cast<meshedSurface&>(surf));
258  meshCells_.transfer(surf.meshCells());
259  }
260  else
261  {
262  // Combine faces with point offsets
263  //
264  // Note: use points().size() from surface, not nPoints()
265  // since there may be uncompacted dangling nodes
266 
267  label nFaces = 0, nPoints = 0;
268 
269  for (const auto& surf : isoSurfPtrs)
270  {
271  nFaces += surf.size();
272  nPoints += surf.points().size();
273  }
274 
275  faceList newFaces(nFaces);
276  pointField newPoints(nPoints);
277  meshCells_.resize(nFaces);
278 
279  surfZoneList newZones(isoSurfPtrs.size());
280 
281  nFaces = 0;
282  nPoints = 0;
283  forAll(isoSurfPtrs, surfi)
284  {
285  autoPtr<isoSurfaceBase> ptr(isoSurfPtrs.release(surfi));
286  auto& surf = *ptr;
287 
288  SubList<face> subFaces(newFaces, surf.size(), nFaces);
289  SubList<point> subPoints(newPoints, surf.points().size(), nPoints);
290  SubList<label> subCells(meshCells_, surf.size(), nFaces);
291 
292  newZones[surfi] = surfZone
293  (
295  subFaces.size(), // size
296  nFaces, // start
297  surfi // index
298  );
299 
300  subFaces = surf.surfFaces();
301  subPoints = surf.points();
302  subCells = surf.meshCells();
303 
304  if (nPoints)
305  {
306  for (face& f : subFaces)
307  {
308  for (label& pointi : f)
309  {
310  pointi += nPoints;
311  }
312  }
313  }
314 
315  nFaces += subFaces.size();
316  nPoints += subPoints.size();
317  }
318 
319  meshedSurface combined
320  (
321  std::move(newPoints),
322  std::move(newFaces),
323  newZones
324  );
325 
326  surface_.transfer(combined);
327  }
328 
329  // Addressing into the full mesh
330  if (subMeshPtr_ && meshCells_.size())
331  {
332  meshCells_ =
333  UIndirectList<label>(subMeshPtr_->cellMap(), meshCells_);
334  }
335 }
336 
337 
338 void Foam::sampledCuttingPlane::createGeometry()
339 {
340  if (debug)
341  {
342  Pout<< "sampledCuttingPlane::createGeometry :updating geometry."
343  << endl;
344  }
345 
346  // Clear any previously stored topologies
347  surface_.clear();
348  meshCells_.clear();
349  isoSurfacePtr_.reset(nullptr);
350 
351  // Clear derived data
353 
354  // Clear any stored fields
355  pointDistance_.clear();
356  cellDistancePtr_.clear();
357 
358  const bool hasCellZones =
359  (-1 != mesh().cellZones().findIndex(zoneNames_));
360 
361  const fvMesh& fvm = static_cast<const fvMesh&>(this->mesh());
362 
363  // Geometry
364  if
365  (
366  simpleSubMesh_
367  && isoParams_.algorithm() != isoSurfaceParams::ALGO_POINT
368  )
369  {
370  subMeshPtr_.reset(nullptr);
371 
372  // Handle cell zones as inverse (blocked) selection
373  if (!ignoreCellsPtr_)
374  {
375  ignoreCellsPtr_.reset(new bitSet);
376 
377  if (hasCellZones)
378  {
379  bitSet select(mesh().cellZones().selection(zoneNames_));
380 
381  if (select.any() && !select.all())
382  {
383  // From selection to blocking
384  select.flip();
385 
386  *ignoreCellsPtr_ = std::move(select);
387  }
388  }
389  }
390  }
391  else
392  {
393  // A standard subMesh treatment
394 
395  if (ignoreCellsPtr_)
396  {
397  ignoreCellsPtr_->clearStorage();
398  }
399  else
400  {
401  ignoreCellsPtr_.reset(new bitSet);
402  }
403 
404  // Get sub-mesh if any
405  if (!subMeshPtr_ && hasCellZones)
406  {
407  const label exposedPatchi =
408  mesh().boundaryMesh().findPatchID(exposedPatchName_);
409 
410  bitSet cellsToSelect(mesh().cellZones().selection(zoneNames_));
411 
412  DebugInfo
413  << "Allocating subset of size "
414  << cellsToSelect.count() << " with exposed faces into patch "
415  << exposedPatchi << endl;
416 
417 
418  // If we will use a fvMeshSubset so can apply bounds as well to make
419  // the initial selection smaller.
420 
421  const boundBox& clipBb = isoParams_.getClipBounds();
422  if (clipBb.good() && cellsToSelect.any())
423  {
424  const auto& cellCentres = fvm.C();
425 
426  for (const label celli : cellsToSelect)
427  {
428  const point& cc = cellCentres[celli];
429 
430  if (!clipBb.contains(cc))
431  {
432  cellsToSelect.unset(celli);
433  }
434  }
435 
436  DebugInfo
437  << "Bounded subset of size "
438  << cellsToSelect.count() << endl;
439  }
440 
441  subMeshPtr_.reset
442  (
443  new fvMeshSubset(fvm, cellsToSelect, exposedPatchi)
444  );
445  }
446  }
447 
448 
449  // Select either the submesh or the underlying mesh
450  const fvMesh& mesh =
451  (
452  subMeshPtr_
453  ? subMeshPtr_->subMesh()
454  : fvm
455  );
456 
457  checkBoundsIntersection(plane_, mesh.bounds());
458 
459 
460  // Distance to cell centres
461  // ~~~~~~~~~~~~~~~~~~~~~~~~
462 
463  cellDistancePtr_.reset
464  (
465  new volScalarField
466  (
467  IOobject
468  (
469  "cellDistance",
470  mesh.time().timeName(),
471  mesh.time(),
475  ),
476  mesh,
478  )
479  );
480  const volScalarField& cellDistance = cellDistancePtr_();
481 
482  setDistanceFields(plane_);
483 
484  if (debug)
485  {
486  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
487  cellDistance.write();
488  pointScalarField pointDist
489  (
490  IOobject
491  (
492  "pointDistance",
493  mesh.time().timeName(),
494  mesh.time(),
498  ),
501  );
502  pointDist.primitiveFieldRef() = pointDistance_;
503 
504  Pout<< "Writing point distance:" << pointDist.objectPath() << endl;
505  pointDist.write();
506  }
507 
508 
509  // Create surfaces for each offset
510 
511  PtrList<isoSurfaceBase> isoSurfPtrs(offsets_.size());
512 
513  forAll(offsets_, surfi)
514  {
515  isoSurfPtrs.set
516  (
517  surfi,
519  (
520  isoParams_,
521  cellDistance,
522  pointDistance_,
523  offsets_[surfi],
524  *ignoreCellsPtr_
525  )
526  );
527  }
528 
529  // And flatten
530  combineSurfaces(isoSurfPtrs);
531 
532 
533  // Discard fields if not required by an iso-surface
534  if (!isoSurfacePtr_)
535  {
536  cellDistancePtr_.reset(nullptr);
537  pointDistance_.clear();
538  }
539 
540  if (debug)
541  {
542  print(Pout, debug);
543  Pout<< endl;
544  }
545 }
546 
547 
548 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
549 
551 (
552  const word& name,
553  const polyMesh& mesh,
554  const dictionary& dict
555 )
556 :
557  sampledSurface(name, mesh, dict),
558  plane_(definePlane(mesh, dict)),
559  offsets_(),
560  isoParams_
561  (
562  dict,
563  isoSurfaceParams::ALGO_TOPO,
564  isoSurfaceParams::filterType::DIAGCELL
565  ),
566  average_(dict.getOrDefault("average", false)),
567  simpleSubMesh_(dict.getOrDefault("simpleSubMesh", false)),
568  zoneNames_(),
569  exposedPatchName_(),
570  needsUpdate_(true),
571 
572  surface_(),
573  meshCells_(),
574  isoSurfacePtr_(nullptr),
575 
576  subMeshPtr_(nullptr),
577  ignoreCellsPtr_(nullptr),
578  cellDistancePtr_(nullptr),
579  pointDistance_()
580 {
581  dict.readIfPresent("offsets", offsets_);
582 
583  if (offsets_.empty())
584  {
585  offsets_.resize(1);
586  offsets_.first() = Zero;
587  }
588 
589  if (offsets_.size() > 1)
590  {
591  const label nOrig = offsets_.size();
592 
593  inplaceUniqueSort(offsets_);
594 
595  if (nOrig != offsets_.size())
596  {
598  << "Removed non-unique offsets" << nl;
599  }
600  }
601 
602  if (isoParams_.algorithm() == isoSurfaceParams::ALGO_POINT)
603  {
604  // Not possible for ALGO_POINT
605  simpleSubMesh_ = false;
606 
607  // Not possible for ALGO_POINT
608  if (offsets_.size() > 1)
609  {
611  << "Multiple offsets with iso-surface (point) not supported"
612  << " since needs original interpolators." << nl
613  << exit(FatalIOError);
614  }
615  }
616 
617 
618  // Zones
619 
620  if (!dict.readIfPresent("zones", zoneNames_) && dict.found("zone"))
621  {
622  zoneNames_.resize(1);
623  dict.readEntry("zone", zoneNames_.first());
624  }
625 
626  if (-1 != mesh.cellZones().findIndex(zoneNames_))
627  {
628  dict.readIfPresent("exposedPatchName", exposedPatchName_);
629 
630  DebugInfo
631  << "Restricting to cellZone(s) " << flatOutput(zoneNames_)
632  << " with exposed internal faces into patch "
633  << mesh.boundaryMesh().findPatchID(exposedPatchName_) << endl;
634  }
635 }
636 
637 
638 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
641 {
642  return needsUpdate_;
643 }
644 
645 
647 {
648  if (debug)
649  {
650  Pout<< "sampledCuttingPlane::expire :"
651  << " needsUpdate:" << needsUpdate_ << endl;
652  }
653 
654  surface_.clear();
655  meshCells_.clear();
656  isoSurfacePtr_.reset(nullptr);
657 
658  // Clear derived data
660 
661  // Already marked as expired
662  if (needsUpdate_)
663  {
664  return false;
665  }
666 
667  needsUpdate_ = true;
668  return true;
669 }
670 
671 
673 {
674  if (debug)
675  {
676  Pout<< "sampledCuttingPlane::update :"
677  << " needsUpdate:" << needsUpdate_ << endl;
678  }
679 
680  if (!needsUpdate_)
681  {
682  return false;
683  }
684 
685  createGeometry();
686 
687  needsUpdate_ = false;
688  return true;
689 }
690 
691 
694 (
695  const interpolation<scalar>& sampler
696 ) const
697 {
698  return sampleOnFaces(sampler);
699 }
700 
701 
704 (
705  const interpolation<vector>& sampler
706 ) const
707 {
708  return sampleOnFaces(sampler);
709 }
710 
711 
714 (
715  const interpolation<sphericalTensor>& sampler
716 ) const
717 {
718  return sampleOnFaces(sampler);
719 }
720 
721 
724 (
725  const interpolation<symmTensor>& sampler
726 ) const
727 {
728  return sampleOnFaces(sampler);
729 }
730 
731 
734 (
735  const interpolation<tensor>& sampler
736 ) const
737 {
738  return sampleOnFaces(sampler);
739 }
740 
741 
744 (
745  const interpolation<scalar>& interpolator
746 ) const
747 {
748  return sampleOnPoints(interpolator);
749 }
750 
751 
754 (
755  const interpolation<vector>& interpolator
756 ) const
757 {
758  return sampleOnPoints(interpolator);
759 }
760 
761 
764 (
765  const interpolation<sphericalTensor>& interpolator
766 ) const
767 {
768  return sampleOnPoints(interpolator);
769 }
770 
771 
774 (
775  const interpolation<symmTensor>& interpolator
776 ) const
777 {
778  return sampleOnPoints(interpolator);
779 }
780 
781 
784 (
785  const interpolation<tensor>& interpolator
786 ) const
787 {
788  return sampleOnPoints(interpolator);
789 }
790 
791 
792 void Foam::sampledCuttingPlane::print(Ostream& os, int level) const
793 {
794  os << "sampledCuttingPlane: " << name() << " :"
795  << " plane:" << plane_
796  << " offsets:" << flatOutput(offsets_);
797 
798  if (level)
799  {
800  os << " faces:" << faces().size()
801  << " points:" << points().size();
802  }
803 }
804 
805 
806 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
bool interpolate() const noexcept
Same as isPointData()
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
virtual bool expire()
Mark the surface as needing an update.
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
label nPoints() const noexcept
Number of mesh points.
label findIndex(const wordRe &key) const
Zone index for the first match, return -1 if not found.
Definition: ZoneMesh.C:601
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
Ostream & print(Ostream &os, UIntType value, char off='0', char on='1')
Print 0/1 bits in the (unsigned) integral type.
Definition: BitOps.H:323
List< bool > select(const label n, const labelUList &locations)
Construct a selection list of bools (all false) with the given pre-size, subsequently add specified l...
Definition: BitOps.C:134
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:90
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
Macros for easy insertion into run-time selection tables.
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
fvPatchField< scalar > fvPatchScalarField
virtual bool needsUpdate() const
Does the surface need an update?
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
MeshedSurface< face > meshedSurface
static word defaultName(const label n=-1)
Default zone name: "zone" or "zoneN".
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Vector< scalar > vector
Definition: vector.H:57
void inplaceUniqueSort(ListType &input)
Inplace sorting and removal of duplicates.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
const vectorField & cellCentres() const
#define DebugInfo
Report an information message using Foam::Info.
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
algorithmType algorithm() const noexcept
Get current algorithm.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const vectorField & faceCentres() const
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
sampledCuttingPlane(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
List< surfZone > surfZoneList
List of surfZone.
Definition: surfZoneList.H:32
Nothing to be read.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
SubField< vector > subField
Declare type of subField.
Definition: Field.H:128
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:616
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
const volVectorField & C() const
Return cell centres as volVectorField.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
virtual bool update()
Update the surface as required.
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
virtual void print(Ostream &os, int level=0) const
Print information.
static autoPtr< isoSurfaceBase > New(const isoSurfaceParams &params, const volScalarField &cellValues, const scalarField &pointValues, const scalar iso, const bitSet &ignoreCells=bitSet())
Create for specified algorithm type.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
const pointField & pts
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127