54bool Foam::sampledThresholdCellFaces::updateGeometry()
const
56 const fvMesh& fvm =
static_cast<const fvMesh&
>(
mesh());
59 if (fvm.time().timeIndex() == prevTimeIndex_)
64 prevTimeIndex_ = fvm.time().timeIndex();
68 const auto* cellFldPtr = fvm.findObject<
volScalarField>(fieldName_);
79 <<
"Reading " << fieldName_
80 <<
" from time " << fvm.time().timeName()
86 autoPtr<volScalarField> fieldReadPtr;
97 fvm.time().timeName(),
108 (fieldReadPtr ? *fieldReadPtr : *cellFldPtr);
112 thresholdCellFaces surf
115 cellFld.primitiveField(),
121 mySurface.transfer(
static_cast<Mesh&
>(surf));
122 meshCells_.transfer(surf.meshCells());
129 Pout<<
"sampledThresholdCellFaces::updateGeometry() : constructed"
131 <<
" field : " << fieldName_ <<
nl
132 <<
" lowerLimit : " << lowerThreshold_ <<
nl
133 <<
" upperLimit : " << upperThreshold_ <<
nl
136 <<
" cut cells : " << meshCells_.size() <<
endl;
153 fieldName_(
dict.get<
word>(
"field")),
154 lowerThreshold_(
dict.getOrDefault<scalar>(
"lowerLimit", -VGREAT)),
155 upperThreshold_(
dict.getOrDefault<scalar>(
"upperLimit", VGREAT)),
156 triangulate_(
dict.getOrDefault(
"triangulate", false)),
160 if (!
dict.found(
"lowerLimit") && !
dict.found(
"upperLimit"))
163 <<
"require at least one of 'lowerLimit' or 'upperLimit'" << endl
164 << abort(FatalError);
175 return fvm.time().timeIndex() != prevTimeIndex_;
182 if (prevTimeIndex_ == -1)
195 return updateGeometry();
204 return sampleOnFaces(sampler);
213 return sampleOnFaces(sampler);
222 return sampleOnFaces(sampler);
231 return sampleOnFaces(sampler);
240 return sampleOnFaces(sampler);
249 return sampleOnPoints(interpolator);
258 return sampleOnPoints(interpolator);
267 return sampleOnPoints(interpolator);
276 return sampleOnPoints(interpolator);
285 return sampleOnPoints(interpolator);
291 os <<
"sampledThresholdCellFaces: " <<
name() <<
" :"
292 <<
" field:" << fieldName_
293 <<
" lowerLimit:" << lowerThreshold_
294 <<
" upperLimit:" << upperThreshold_;
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
bool get(const label i) const
void size(const label n)
Older name for setAddressableSize.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh data needed to do the Finite Volume discretisation.
Abstract base class for volume field interpolation.
Mesh consisting of general polyhedral cells.
An abstract class for surfaces with sampling.
sampledSurface(const word &name, std::nullptr_t)
Construct null.
const word & name() const noexcept
Name of surface.
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
bool interpolate() const noexcept
Same as isPointData().
A sampledSurface defined by the cell faces corresponding to a threshold value.
virtual const pointField & points() const
Points of surface.
sampledThresholdCellFaces(const word &name, const polyMesh &, const dictionary &)
Construct from dictionary.
virtual void print(Ostream &os, int level=0) const
Print information.
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
virtual const faceList & faces() const
Faces of 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.
Selects the mesh cell faces specified by a threshold value. Non-triangulated by default.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
OBJstream os(runTime.globalPath()/outputName)
#define InfoInFunction
Report an information message using Foam::Info.
Namespace of functions to calculate implicit derivatives returning a matrix.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a).