Go to the documentation of this file.
35 #undef defineExpressionMethod
36 #define defineExpressionMethod(Type, Member) \
38 inline const Type& exprResult::singleValue::get<Type>() const \
44 inline const Type& exprResult::singleValue::set(const Type& val) \
58 #undef defineExpressionMethod
67 inline bool Foam::expressions::exprResult::deleteChecked()
69 const bool ok = isType<Type>();
71 if (ok && fieldPtr_ !=
nullptr)
73 delete static_cast<Field<Type>*
>(fieldPtr_);
83 inline bool Foam::expressions::exprResult::readChecked
86 const dictionary&
dict,
91 const bool ok = isType<Type>();
99 const Type val(
dict.get<Type>(
key));
102 fieldPtr_ =
new Field<Type>(size_, val);
109 fieldPtr_ =
new Field<Type>(
key,
dict, size_);
112 isUniform_ = uniform;
120 bool Foam::expressions::exprResult::getUniformChecked
135 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
143 if (limits.mag() > SMALL)
146 <<
"Different min/max values: " << limits
147 <<
" Using the average " << avg <<
nl;
151 result.setResult(avg, size);
158 bool Foam::expressions::exprResult::plusEqChecked
163 const bool ok = isType<Type>();
167 *
static_cast<Field<Type>*
>(fieldPtr_)
168 += *
static_cast<const Field<Type>*
>(
b.fieldPtr_);
176 bool Foam::expressions::exprResult::multiplyEqChecked
181 const bool ok = isType<Type>();
185 *
static_cast<Field<Type>*
>(fieldPtr_) *=
b;
230 return (!valType_.empty() && fieldPtr_ !=
nullptr);
243 const bool wantPointData
246 return isPointData_ == wantPointData;
266 if (!isUniform_ || !isType<Type>())
271 return single_.get<Type>();
296 target().setResultImpl(val, wantPointData);
309 target().setResultImpl(val, wantPointData);
314 void Foam::expressions::exprResult::setResultImpl
325 isPointData_ = wantPointData;
336 void Foam::expressions::exprResult::setResultImpl
347 isPointData_ = wantPointData;
364 target().setResultImpl(fldPtr, wantPointData);
369 void Foam::expressions::exprResult::setResultImpl
378 isPointData_ = wantPointData;
380 if (fldPtr !=
nullptr)
382 size_ = fldPtr->size();
396 target().setResultImpl(val, size);
401 void Foam::expressions::exprResult::setResultImpl
411 isPointData_ =
false;
425 target().setSingleValueImpl(val);
430 bool Foam::expressions::exprResult::writeSingleValueChecked(
Ostream&
os)
const
437 if (this->size() <= 0)
441 os << single_.get<Type>();
446 os << pTraits<Type>::zero;
451 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
461 bool Foam::expressions::exprResult::writeFieldChecked
472 if (this->size() <= 0)
476 const Type& val = single_.get<Type>();
483 os.writeEntry(keyword, val);
491 os << pTraits<Type>::zero;
495 Field<Type>().writeEntry(keyword,
os);
501 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
511 os.writeEntry(keyword,
fld.first());
515 fld.writeEntry(keyword,
os);
525 bool Foam::expressions::exprResult::writeEntryChecked
536 if (this->size() <= 0)
538 if (isUniform_ && is_contiguous<Type>::value)
540 const Type& val = single_.get<Type>();
544 os.writeKeyword(keyword);
552 const Field<Type>
fld;
553 fld.writeEntry(keyword,
os);
558 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
560 if (isUniform_ && is_contiguous<Type>::value)
564 os.writeKeyword(keyword);
571 fld.writeEntry(keyword,
os);
580 bool Foam::expressions::exprResult::setAverageValueChecked(
const bool parRun)
587 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
591 isUniform_ = (limits.mag() <= SMALL);
593 const Type avg = limits.centre();
602 bool Foam::expressions::exprResult::duplicateFieldChecked(
const void* ptr)
611 deleteChecked<Type>();
614 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(ptr);
617 fieldPtr_ =
new Field<Type>(
fld);
624 void Foam::expressions::exprResult::setSingleValueImpl(
const Type& val)
631 isPointData_ =
false;
636 valType_ = pTraits<Type>::typeName;
637 fieldPtr_ =
new Field<Type>(size_, val);
651 <<
" is different from the stored result type "
652 << valType_ <<
nl <<
nl
656 if (fieldPtr_ ==
nullptr)
659 <<
"Cannot create tmp from nullptr." <<
nl
660 <<
"This error message should never appear!!" <<
nl
693 <<
" is different from the stored result type "
694 << valType_ <<
nl <<
nl
698 if (fieldPtr_ ==
nullptr)
701 <<
"Cannot return reference from nullptr." <<
nl
702 <<
"This error message should never appear!!" <<
nl
706 return *
static_cast<const Field<Type>*
>(fieldPtr_);
714 return const_cast<Field<Type>&
>(this->cref<Type>());
722 return const_cast<Field<Type>&
>(this->cref<Type>());
726 template<
template<
class>
class BinaryOp,
class Type>
729 const BinaryOp<Type>& bop,
737 <<
" is different from the stored result type "
738 << valType_ <<
nl <<
nl
742 Type result = initial;
746 for (
const Type& val :
fld)
748 result = bop(result, val);
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
A class for handling words, derived from Foam::string.
label size() const
The field or object size.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Type getReduced(const BinaryOp< Type > &bop, const Type &initial=pTraits< Type >::zero)
Get a reduced result.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
bool isPointData(const bool wantPointData=true) const
Type gAverage(const FieldField< Field, Type > &f)
tmp< Field< Type > > getResult(bool cacheCopy=false)
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
const Type & value() const
Return const reference to value.
Field< Type > & constCast() const
Return non-const reference to the field, casting away constness.
A polymorphic field/result from evaluating an expression.
const word & valueType() const noexcept
Basic type for the field or single value.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Generic templated field type.
#define defineExpressionMethod(Type, Member)
#define DebugInFunction
Report an information message using Foam::Info.
Field< Type > & ref()
Return non-const reference to the field.
void setSingleValue(const Type &val)
Set single-value uniform result.
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
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))
bool hasValue() const
Has a value?
bool isUniform() const
True if single, uniform value.
OBJstream os(runTime.globalPath()/outputName)
Generic dimensioned Type class.
const Field< Type > & cref() const
Return const reference to the field.
bool isType() const
True if valueType corresponds to the given Type.
Vector< scalar > vector
A scalar version of the templated Vector.
errorManipArg< error, int > exit(error &err, const int errNo=1)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
exprResult()
Default construct.
A traits class, which is primarily used for primitives.
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
void setResult(Field< Type > *, bool wantPointData=false)
Set result field, taking ownership of the pointer.
#define WarningInFunction
Report a warning using Foam::Warning.
bool is_bool() const
True if valueType is a bool.
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)