Go to the documentation of this file.
45 sampledThresholdCellFaces,
54 bool 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
134 <<
" point : " <<
points().size() <<
nl
135 <<
" faces : " <<
faces().size() <<
nl
136 <<
" cut cells : " << meshCells_.size() <<
endl;
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
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_;
static autoPtr< T > New(Args &&... args)
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
virtual void print(Ostream &os, int level=0) const
Ostream & endl(Ostream &os)
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
sampledThresholdCellFaces(const word &name, const polyMesh &, const dictionary &)
virtual bool needsUpdate() const
Mesh consisting of general polyhedral cells.
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
label timeIndex() const noexcept
bool interpolate() const noexcept
void transfer(List< T > &list)
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
An abstract class for surfaces with sampling.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
errorManip< error > abort(error &err)
#define FatalErrorInFunction
const polyMesh & mesh() const noexcept
virtual void clearGeom() const
virtual const faceList & faces() const
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
word name(const expressions::valueTypeCode typeCode)
const Time & time() const
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
defineTypeNameAndDebug(combustionModel, 0)
virtual const pointField & points() const