Go to the documentation of this file.
35 if (interpolatorPtr_.empty())
38 tableSamplesPtr_.reset(
new scalarField(table_.size()));
42 tableSamples[i] = table_[i].first();
51 return interpolatorPtr_();
66 dict.lookupOrDefault<
word>(
"outOfBounds",
"clamp")
71 dict.lookupOrDefault<
word>(
"interpolationScheme",
"linear")
83 boundsHandling_(tbl.boundsHandling_),
84 interpolationScheme_(tbl.interpolationScheme_),
86 dimensions_(tbl.dimensions_),
87 tableSamplesPtr_(tbl.tableSamplesPtr_),
88 interpolatorPtr_(tbl.interpolatorPtr_)
107 word enumName(
"warn");
144 if (
bound ==
"error")
148 else if (
bound ==
"warn")
152 else if (
bound ==
"clamp")
156 else if (
bound ==
"repeat")
163 <<
"bad outOfBounds specifier " <<
bound <<
" using 'warn'"
179 boundsHandling_ =
bound;
191 <<
"Table for entry " << this->name_ <<
" is invalid (empty)"
196 scalar prevValue = table_[0].first();
198 for (
label i = 1; i <
n; ++i)
200 const scalar currValue = table_[i].first();
203 if (currValue <= prevValue)
206 <<
"out-of-order value: " << currValue <<
" at index " << i
209 prevValue = currValue;
221 if (
x < table_[0].first())
223 switch (boundsHandling_)
228 <<
"value (" <<
x <<
") underflow"
235 <<
"value (" <<
x <<
") underflow" <<
nl
242 xDash = table_[0].first();
249 scalar span = table_.last().first() - table_[0].first();
250 xDash = fmod(
x - table_[0].first(), span) + table_[0].first();
271 if (
x > table_.last().first())
273 switch (boundsHandling_)
278 <<
"value (" <<
x <<
") overflow"
285 <<
"value (" <<
x <<
") overflow" <<
nl
292 xDash = table_.last().first();
299 scalar span = table_.last().first() - table_[0].first();
300 xDash = fmod(
x - table_[0].first(), span) + table_[0].first();
319 scalar value = table_[i].first();
323 tableSamplesPtr_.clear();
324 interpolatorPtr_.clear();
333 if (checkMinBounds(
x, xDash))
335 return table_[0].second();
338 if (checkMaxBounds(xDash, xDash))
340 return table_.last().second();
344 interpolator().valueWeights(xDash, currentIndices_, currentWeights_);
346 Type t = currentWeights_[0]*table_[currentIndices_[0]].second();
347 for (
label i = 1; i < currentIndices_.size(); i++)
349 t += currentWeights_[i]*table_[currentIndices_[i]].second();
360 interpolator().integrationWeights(x1, x2, currentIndices_, currentWeights_);
362 Type
sum = currentWeights_[0]*table_[currentIndices_[0]].second();
363 for (
label i = 1; i < currentIndices_.size(); i++)
365 sum += currentWeights_[i]*table_[currentIndices_[i]].second();
383 const scalar x1,
const scalar x2
390 this->integrate(x2, x1)
403 fld[i] = table_[i].first();
418 fld[i] = table_[i].second();
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
virtual tmp< scalarField > x() const
Return the reference values.
A class for handling words, derived from string.
virtual tmp< Field< Type > > y() const
Return the dependent values.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
virtual dimensioned< Type > dimValue(const scalar x) const
Return dimensioned constant value.
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Ostream & endl(Ostream &os)
Add newline and flush stream.
boundsHandling
Enumeration for handling out-of-bound values.
bool checkMinBounds(const scalar x, scalar &xDash) const
Check minimum table bounds.
virtual Type integrate(const scalar x1, const scalar x2) const
Integrate between two (scalar) values.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
TableBase(const word &name, const dictionary &dict)
Construct from dictionary - note table is not populated.
virtual ~TableBase()
Destructor.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Abstract base class for interpolating in 1D.
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){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool checkMaxBounds(const scalar x, scalar &xDash) const
Check maximum table bounds.
errorManipArg< error, int > exit(error &err, const int errNo=1)
boundsHandling outOfBounds(const boundsHandling &bound)
Set the out-of-bounds handling from enum, return previous setting.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Traits class for primitives.
word boundsHandlingToWord(const boundsHandling &bound) const
Return the out-of-bounds handling as a word.
volScalarField scalarField(fieldObject, mesh)
virtual Type value(const scalar x) const
Return Table value.
void check() const
Check the table for size and consistency.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
boundsHandling wordToBoundsHandling(const word &bound) const
Return the out-of-bounds handling as an enumeration.
virtual void convertTimeBase(const Time &t)
Convert time.
virtual dimensioned< Type > dimIntegrate(const scalar x1, const scalar x2) const
Integrate between two values and return dimensioned type.
const interpolationWeights & interpolator() const
Return (demand driven) interpolator.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
#define WarningInFunction
Report a warning using Foam::Warning.
word name(const complex &)
Return a string representation of a complex.
Base class for table with bounds handling, interpolation and integration.