sampledIsoSurface.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 "sampledIsoSurface.H"
27 #include "dictionary.H"
28 #include "volFields.H"
29 #include "volPointInterpolation.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(sampledIsoSurface, 0);
38  (
39  sampledSurface,
40  sampledIsoSurface,
41  word,
42  isoSurface
43  );
44 }
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
49 {
50  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
51 
52  // Get volField
53  // ~~~~~~~~~~~~
54 
56  {
57  if (debug)
58  {
59  Info<< "sampledIsoSurface::getIsoFields() : lookup volField "
60  << isoField_ << endl;
61  }
62  storedVolFieldPtr_.clear();
64  }
65  else
66  {
67  // Bit of a hack. Read field and store.
68 
69  if (debug)
70  {
71  Info<< "sampledIsoSurface::getIsoFields() : checking "
72  << isoField_ << " for same time " << fvm.time().timeName()
73  << endl;
74  }
75 
76  if
77  (
78  storedVolFieldPtr_.empty()
79  || (fvm.time().timeName() != storedVolFieldPtr_().instance())
80  )
81  {
82  if (debug)
83  {
84  Info<< "sampledIsoSurface::getIsoFields() : reading volField "
85  << isoField_ << " from time " << fvm.time().timeName()
86  << endl;
87  }
88 
89  IOobject vfHeader
90  (
91  isoField_,
92  fvm.time().timeName(),
93  fvm,
96  false
97  );
98 
99  if (vfHeader.headerOk())
100  {
101  storedVolFieldPtr_.reset
102  (
103  new volScalarField
104  (
105  vfHeader,
106  fvm
107  )
108  );
109  volFieldPtr_ = storedVolFieldPtr_.operator->();
110  }
111  else
112  {
114  << "Cannot find isosurface field " << isoField_
115  << " in database or directory " << vfHeader.path()
116  << exit(FatalError);
117  }
118  }
119  }
120 
121 
122  // Get pointField
123  // ~~~~~~~~~~~~~~
124 
125  // In case of multiple iso values we don't want to calculate multiple e.g.
126  // "volPointInterpolate(p)" so register it and re-use it. This is the
127  // same as the 'cache' functionality from volPointInterpolate but
128  // unfortunately that one does not guarantee that the field pointer
129  // remain: e.g. some other
130  // functionObject might delete the cached version.
131  // (volPointInterpolation::interpolate with cache=false deletes any
132  // registered one or if mesh.changing())
133 
134  if (!subMeshPtr_.valid())
135  {
136  const word pointFldName =
137  "volPointInterpolate_"
138  + type()
139  + "("
140  + isoField_
141  + ')';
142 
143  if (fvm.foundObject<pointScalarField>(pointFldName))
144  {
145  if (debug)
146  {
147  Info<< "sampledIsoSurface::getIsoFields() :"
148  << " lookup pointField " << pointFldName << endl;
149  }
151  (
152  pointFldName
153  );
154 
155  if (!pfld.upToDate(*volFieldPtr_))
156  {
157  if (debug)
158  {
159  Info<< "sampledIsoSurface::getIsoFields() :"
160  << " updating pointField " << pointFldName << endl;
161  }
162  // Update the interpolated value
164  (
165  *volFieldPtr_,
166  const_cast<pointScalarField&>(pfld)
167  );
168  }
169 
170  pointFieldPtr_ = &pfld;
171  }
172  else
173  {
174  // Not in registry. Interpolate.
175 
176  if (debug)
177  {
178  Info<< "sampledIsoSurface::getIsoFields() :"
179  << " creating and storing pointField "
180  << pointFldName << " for time "
181  << fvm.time().timeName() << endl;
182  }
183 
185  (
187  (
188  *volFieldPtr_,
189  pointFldName,
190  false
191  )
192  );
193  pointFieldPtr_ = tpfld.ptr();
194  const_cast<pointScalarField*>(pointFieldPtr_)->store();
195  }
196 
197 
198  // If averaging redo the volField. Can only be done now since needs the
199  // point field.
200  if (average_)
201  {
202  storedVolFieldPtr_.reset
203  (
204  pointAverage(*pointFieldPtr_).ptr()
205  );
206  volFieldPtr_ = storedVolFieldPtr_.operator->();
207  }
208 
209 
210  if (debug)
211  {
212  Info<< "sampledIsoSurface::getIsoFields() : volField "
213  << volFieldPtr_->name() << " min:" << min(*volFieldPtr_).value()
214  << " max:" << max(*volFieldPtr_).value() << endl;
215  Info<< "sampledIsoSurface::getIsoFields() : pointField "
216  << pointFieldPtr_->name()
217  << " min:" << gMin(pointFieldPtr_->internalField())
218  << " max:" << gMax(pointFieldPtr_->internalField()) << endl;
219  }
220  }
221  else
222  {
223  // Get subMesh variants
224  const fvMesh& subFvm = subMeshPtr_().subMesh();
225 
226  // Either lookup on the submesh or subset the whole-mesh volField
227 
228  if (subFvm.foundObject<volScalarField>(isoField_))
229  {
230  if (debug)
231  {
232  Info<< "sampledIsoSurface::getIsoFields() :"
233  << " submesh lookup volField "
234  << isoField_ << endl;
235  }
236  storedVolSubFieldPtr_.clear();
238  }
239  else
240  {
241  if (debug)
242  {
243  Info<< "sampledIsoSurface::getIsoFields() : "
244  << "subsetting volField " << isoField_ << endl;
245  }
247  (
248  subMeshPtr_().interpolate
249  (
250  *volFieldPtr_
251  ).ptr()
252  );
253  storedVolSubFieldPtr_->checkOut();
254  volSubFieldPtr_ = storedVolSubFieldPtr_.operator->();
255  }
256 
257 
258  // Pointfield on submesh
259 
260  const word pointFldName =
261  "volPointInterpolate_"
262  + type()
263  + "("
264  + volSubFieldPtr_->name()
265  + ')';
266 
267  if (subFvm.foundObject<pointScalarField>(pointFldName))
268  {
269  if (debug)
270  {
271  Info<< "sampledIsoSurface::getIsoFields() :"
272  << " submesh lookup pointField " << pointFldName << endl;
273  }
274  const pointScalarField& pfld = subFvm.lookupObject<pointScalarField>
275  (
276  pointFldName
277  );
278 
279  if (!pfld.upToDate(*volSubFieldPtr_))
280  {
281  if (debug)
282  {
283  Info<< "sampledIsoSurface::getIsoFields() :"
284  << " updating submesh pointField "
285  << pointFldName << endl;
286  }
287  // Update the interpolated value
289  (
291  const_cast<pointScalarField&>(pfld)
292  );
293  }
294 
295  pointFieldPtr_ = &pfld;
296  }
297  else
298  {
299  if (debug)
300  {
301  Info<< "sampledIsoSurface::getIsoFields() :"
302  << " interpolating submesh volField "
303  << volSubFieldPtr_->name()
304  << " to get submesh pointField " << pointFldName << endl;
305  }
307  (
309  (
310  subFvm
312  );
313  pointSubFieldPtr_ = tpfld.ptr();
314  const_cast<pointScalarField*>(pointSubFieldPtr_)->store();
315  }
316 
317 
318  // If averaging redo the volField. Can only be done now since needs the
319  // point field.
320  if (average_)
321  {
323  (
324  pointAverage(*pointSubFieldPtr_).ptr()
325  );
326  volSubFieldPtr_ = storedVolSubFieldPtr_.operator->();
327  }
328 
329 
330  if (debug)
331  {
332  Info<< "sampledIsoSurface::getIsoFields() : volSubField "
333  << volSubFieldPtr_->name()
334  << " min:" << min(*volSubFieldPtr_).value()
335  << " max:" << max(*volSubFieldPtr_).value() << endl;
336  Info<< "sampledIsoSurface::getIsoFields() : pointSubField "
337  << pointSubFieldPtr_->name()
338  << " min:" << gMin(pointSubFieldPtr_->internalField())
339  << " max:" << gMax(pointSubFieldPtr_->internalField()) << endl;
340  }
341  }
342 }
343 
344 
346 {
347  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
348 
349  // No update needed
350  if (fvm.time().timeIndex() == prevTimeIndex_)
351  {
352  return false;
353  }
354 
355  // Get any subMesh
356  if (zoneID_.index() != -1 && !subMeshPtr_.valid())
357  {
359 
360  // Patch to put exposed internal faces into
361  const label exposedPatchI = patches.findPatchID(exposedPatchName_);
362 
363  if (debug)
364  {
365  Info<< "Allocating subset of size "
366  << mesh().cellZones()[zoneID_.index()].size()
367  << " with exposed faces into patch "
368  << patches[exposedPatchI].name() << endl;
369  }
370 
371  subMeshPtr_.reset
372  (
373  new fvMeshSubset(fvm)
374  );
375  subMeshPtr_().setLargeCellSubset
376  (
377  labelHashSet(mesh().cellZones()[zoneID_.index()]),
378  exposedPatchI
379  );
380  }
381 
382 
383  prevTimeIndex_ = fvm.time().timeIndex();
384  getIsoFields();
385 
386  // Clear any stored topo
387  surfPtr_.clear();
388  facesPtr_.clear();
389 
390  // Clear derived data
391  clearGeom();
392 
393  if (subMeshPtr_.valid())
394  {
395  const volScalarField& vfld = *volSubFieldPtr_;
396 
397  surfPtr_.reset
398  (
399  new isoSurface
400  (
401  vfld,
402  *pointSubFieldPtr_,
403  isoVal_,
404  regularise_,
405  bounds_,
406  mergeTol_
407  )
408  );
409  }
410  else
411  {
412  const volScalarField& vfld = *volFieldPtr_;
413 
414  surfPtr_.reset
415  (
416  new isoSurface
417  (
418  vfld,
419  *pointFieldPtr_,
420  isoVal_,
421  regularise_,
422  bounds_,
423  mergeTol_
424  )
425  );
426  }
427 
428 
429  if (debug)
430  {
431  Pout<< "sampledIsoSurface::updateGeometry() : constructed iso:"
432  << nl
433  << " regularise : " << regularise_ << nl
434  << " average : " << average_ << nl
435  << " isoField : " << isoField_ << nl
436  << " isoValue : " << isoVal_ << nl;
437  if (subMeshPtr_.valid())
438  {
439  Pout<< " zone size : " << subMeshPtr_().subMesh().nCells()
440  << nl;
441  }
442  Pout<< " points : " << points().size() << nl
443  << " tris : " << surface().size() << nl
444  << " cut cells : " << surface().meshCells().size()
445  << endl;
446  }
447 
448  return true;
449 }
450 
451 
452 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
453 
455 (
456  const word& name,
457  const polyMesh& mesh,
458  const dictionary& dict
459 )
460 :
462  isoField_(dict.lookup("isoField")),
463  isoVal_(readScalar(dict.lookup("isoValue"))),
464  bounds_(dict.lookupOrDefault("bounds", boundBox::greatBox)),
465  mergeTol_(dict.lookupOrDefault("mergeTol", 1e-6)),
466  regularise_(dict.lookupOrDefault("regularise", true)),
467  average_(dict.lookupOrDefault("average", false)),
468  zoneID_(dict.lookupOrDefault("zone", word::null), mesh.cellZones()),
469  exposedPatchName_(word::null),
470  surfPtr_(NULL),
471  facesPtr_(NULL),
472  prevTimeIndex_(-1),
473  storedVolFieldPtr_(NULL),
474  volFieldPtr_(NULL),
475  pointFieldPtr_(NULL)
476 {
478  {
480  (
481  dict
482  ) << "Non-interpolated iso surface not supported since triangles"
483  << " span across cells." << exit(FatalIOError);
484  }
485 
486  if (zoneID_.index() != -1)
487  {
488  dict.lookup("exposedPatchName") >> exposedPatchName_;
489 
490  if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
491  {
493  (
494  dict
495  ) << "Cannot find patch " << exposedPatchName_
496  << " in which to put exposed faces." << endl
497  << "Valid patches are " << mesh.boundaryMesh().names()
498  << exit(FatalIOError);
499  }
500 
501  if (debug && zoneID_.index() != -1)
502  {
503  Info<< "Restricting to cellZone " << zoneID_.name()
504  << " with exposed internal faces into patch "
505  << exposedPatchName_ << endl;
506  }
507  }
508 }
509 
510 
511 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
512 
514 {}
515 
516 
517 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
518 
520 {
521  const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
522 
523  return fvm.time().timeIndex() != prevTimeIndex_;
524 }
525 
526 
528 {
529  surfPtr_.clear();
530  facesPtr_.clear();
531  subMeshPtr_.clear();
532 
533  // Clear derived data
534  clearGeom();
535 
536  // already marked as expired
537  if (prevTimeIndex_ == -1)
538  {
539  return false;
540  }
541 
542  // force update
543  prevTimeIndex_ = -1;
544  return true;
545 }
546 
547 
549 {
550  return updateGeometry();
551 }
552 
553 
555 (
556  const volScalarField& vField
557 ) const
558 {
559  return sampleField(vField);
560 }
561 
562 
564 (
565  const volVectorField& vField
566 ) const
567 {
568  return sampleField(vField);
569 }
570 
571 
573 (
574  const volSphericalTensorField& vField
575 ) const
576 {
577  return sampleField(vField);
578 }
579 
580 
582 (
583  const volSymmTensorField& vField
584 ) const
585 {
586  return sampleField(vField);
587 }
588 
589 
591 (
592  const volTensorField& vField
593 ) const
594 {
595  return sampleField(vField);
596 }
597 
598 
600 (
601  const interpolation<scalar>& interpolator
602 ) const
603 {
604  return interpolateField(interpolator);
605 }
606 
607 
609 (
610  const interpolation<vector>& interpolator
611 ) const
612 {
613  return interpolateField(interpolator);
614 }
615 
617 (
618  const interpolation<sphericalTensor>& interpolator
619 ) const
620 {
621  return interpolateField(interpolator);
622 }
623 
624 
626 (
627  const interpolation<symmTensor>& interpolator
628 ) const
629 {
630  return interpolateField(interpolator);
631 }
632 
633 
635 (
636  const interpolation<tensor>& interpolator
637 ) const
638 {
639  return interpolateField(interpolator);
640 }
641 
642 
644 {
645  os << "sampledIsoSurface: " << name() << " :"
646  << " field :" << isoField_
647  << " value :" << isoVal_;
648 }
649 
650 
651 // ************************************************************************* //
volFields.H
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
Foam::sampledIsoSurface::updateGeometry
bool updateGeometry() const
Create iso surface (if time has changed)
Definition: sampledIsoSurface.C:345
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(fvPatch, cyclicAMIFvPatch, polyPatch, cyclicPeriodicAMI)
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::sampledIsoSurface::sample
virtual tmp< scalarField > sample(const volScalarField &) const
Sample field on surface.
Definition: sampledIsoSurface.C:555
Foam::sampledIsoSurface::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledIsoSurface.C:519
Foam::sampledSurface::pointAverage
tmp< GeometricField< Type, fvPatchField, volMesh > > pointAverage(const GeometricField< Type, pointPatchField, pointMesh > &pfld) const
Interpolate from points to cell centre.
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
Foam::sampledIsoSurface::pointSubFieldPtr_
const pointScalarField * pointSubFieldPtr_
Cached pointfield.
Definition: sampledIsoSurface.H:111
sampledIsoSurface.H
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::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
Foam::sampledIsoSurface::storedVolFieldPtr_
autoPtr< volScalarField > storedVolFieldPtr_
Cached volfield.
Definition: sampledIsoSurface.H:95
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::FatalIOError
IOerror FatalIOError
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::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::sampledIsoSurface::update
virtual bool update()
Update the surface as required.
Definition: sampledIsoSurface.C:548
Foam::sampledIsoSurface::print
virtual void print(Ostream &) const
Write.
Definition: sampledIsoSurface.C:643
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation >::New
static const volPointInterpolation & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:527
Foam::IOobject::headerOk
bool headerOk()
Read and check header info.
Definition: IOobject.C:439
Foam::isoSurface
A surface formed by the iso value. After "Regularised Marching Tetrahedra: improved iso-surface extra...
Definition: isoSurface.H:85
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::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::sampledIsoSurface::volSubFieldPtr_
const volScalarField * volSubFieldPtr_
Definition: sampledIsoSurface.H:108
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:77
Foam::interpolation< scalar >
Foam::sampledIsoSurface::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledIsoSurface.C:527
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::sampledIsoSurface::isoField_
const word isoField_
Field to get isoSurface of.
Definition: sampledIsoSurface.H:60
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::tmp::ptr
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:146
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::volPointInterpolation::interpolate
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Definition: volPointInterpolate.C:502
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::sampledIsoSurface::getIsoFields
void getIsoFields() const
Get fields needed to recreate iso surface.
Definition: sampledIsoSurface.C:48
volPointInterpolation.H
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sampledIsoSurface::average_
const Switch average_
Whether to recalculate cell values as average of point values.
Definition: sampledIsoSurface.H:75
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::boundBox::greatBox
static const boundBox greatBox
A very large boundBox: min/max == -/+ VGREAT.
Definition: boundBox.H:76
Foam::sampledIsoSurface::pointFieldPtr_
const pointScalarField * pointFieldPtr_
Cached pointfield.
Definition: sampledIsoSurface.H:99
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
points
const pointField & points
Definition: gmvOutputHeader.H:1
dictionary.H
Foam::surface
Definition: surface.H:55
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::sampledIsoSurface::storedVolSubFieldPtr_
autoPtr< volScalarField > storedVolSubFieldPtr_
Cached volfield.
Definition: sampledIsoSurface.H:107
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
patches
patches[0]
Definition: createSingleCellMesh.H:36
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::sampledIsoSurface::subMeshPtr_
autoPtr< fvMeshSubset > subMeshPtr_
Cached submesh.
Definition: sampledIsoSurface.H:104
Foam::sampledSurface::mesh
const polyMesh & mesh() const
Access to the underlying mesh.
Definition: sampledSurface.H:244
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::TimeState::timeIndex
label timeIndex() const
Return current time index.
Definition: TimeState.C:73
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::IOobject::path
fileName path() const
Return complete path.
Definition: IOobject.C:293
Foam::sampledSurface::interpolate
bool interpolate() const
Interpolation requested for surface.
Definition: sampledSurface.H:256
Foam::sampledIsoSurface::~sampledIsoSurface
virtual ~sampledIsoSurface()
Destructor.
Definition: sampledIsoSurface.C:513
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::sampledIsoSurface::sampledIsoSurface
sampledIsoSurface(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: sampledIsoSurface.C:455
Foam::sampledIsoSurface::volFieldPtr_
const volScalarField * volFieldPtr_
Definition: sampledIsoSurface.H:96