Go to the documentation of this file.
35 #undef defineExpressionMethod
36 #define defineExpressionMethod(Type, Var) \
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>();
70 if (ok && fieldPtr_ !=
nullptr)
72 delete static_cast<Field<Type>*
>(fieldPtr_);
82 inline bool Foam::expressions::exprResult::readChecked
85 const dictionary&
dict,
90 const bool ok = isType<Type>();
97 const Type val(
dict.get<Type>(key));
100 fieldPtr_ =
new Field<Type>(size_, val);
107 fieldPtr_ =
new Field<Type>(key,
dict, size_);
110 isUniform_ = uniform;
118 bool Foam::expressions::exprResult::getUniformChecked
133 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
141 if (limits.mag() > SMALL)
144 <<
"Different min/max values: " << limits
145 <<
" Using the average " << avg <<
nl;
149 result.setResult(avg, size);
156 bool Foam::expressions::exprResult::plusEqChecked
161 const bool ok = isType<Type>();
165 *
static_cast<Field<Type>*
>(fieldPtr_)
166 += *
static_cast<const Field<Type>*
>(
b.fieldPtr_);
174 bool Foam::expressions::exprResult::multiplyEqChecked
179 const bool ok = isType<Type>();
183 *
static_cast<Field<Type>*
>(fieldPtr_) *=
b;
246 return (!valType_.empty() && fieldPtr_ !=
nullptr);
258 const bool isPointVal
261 return isPointVal_ == isPointVal;
286 return objectPtr_.valid();
305 target().setResultImpl(val, isPointVal);
318 target().setResultImpl(val, isPointVal);
323 void Foam::expressions::exprResult::setResultImpl
334 isPointVal_ = isPointVal;
345 void Foam::expressions::exprResult::setResultImpl
356 isPointVal_ = isPointVal;
369 target().setObjectResultImpl(ap.
ptr());
376 target().setObjectResultImpl(ap.ptr());
381 void Foam::expressions::exprResult::setObjectResultImpl(
T* ptr)
389 valType_ = ptr->typeName;
390 objectPtr_.reset(ptr);
395 void Foam::expressions::exprResult::setObjectResultImpl(
autoPtr<T>& o)
397 setObjectResultImpl(o.
ptr());
402 void Foam::expressions::exprResult::setObjectResultImpl(
autoPtr<T>&& o)
404 setObjectResultImpl(o.ptr());
415 target().setResultImpl(fldPtr, isPointVal);
420 void Foam::expressions::exprResult::setResultImpl
430 isPointVal_ = isPointVal;
433 size_ = fldPtr->size();
446 target().setResultImpl(val, size);
451 void Foam::expressions::exprResult::setResultImpl
475 target().setSingleValueImpl(val);
480 bool Foam::expressions::exprResult::writeSingleValueChecked(
Ostream& os)
const
487 if (this->size() <= 0)
491 os << single_.get<Type>();
496 os << pTraits<Type>::zero;
501 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
511 bool Foam::expressions::exprResult::writeValueFieldChecked(Ostream& os)
const
518 if (this->size() <= 0)
522 const Type& val = single_.get<Type>();
523 os.writeEntry(
"value", val);
528 const Field<Type>
fld;
529 fld.writeEntry(
"value", os);
534 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
538 os.writeEntry(
"value",
fld.first());
542 fld.writeEntry(
"value", os);
551 bool Foam::expressions::exprResult::writeEntryChecked
562 if (this->size() <= 0)
564 if (isUniform_ && is_contiguous<Type>::value)
566 const Type& val = single_.get<Type>();
570 os.writeKeyword(keyword);
578 const Field<Type>
fld;
579 fld.writeEntry(keyword, os);
584 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
586 if (isUniform_ && is_contiguous<Type>::value)
590 os.writeKeyword(keyword);
597 fld.writeEntry(keyword, os);
606 bool Foam::expressions::exprResult::setAverageValueChecked(
const bool parRun)
613 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(fieldPtr_);
617 isUniform_ = (limits.mag() <= SMALL);
619 const Type avg = limits.centre();
628 bool Foam::expressions::exprResult::duplicateFieldChecked(
const void* ptr)
637 deleteChecked<Type>();
640 const Field<Type>&
fld = *
static_cast<const Field<Type>*
>(ptr);
643 fieldPtr_ =
new Field<Type>(
fld);
650 void Foam::expressions::exprResult::setSingleValueImpl(
const Type& val)
662 valType_ = pTraits<Type>::typeName;
663 fieldPtr_ =
new Field<Type>(size_, val);
677 <<
" is different from the stored result type "
678 << valType_ <<
nl <<
nl
682 if (fieldPtr_ ==
nullptr)
685 <<
"Cannot create tmp from nullptr." <<
nl
686 <<
"This error message should never appear!!" <<
nl
719 <<
" is different from the stored result type "
720 << valType_ <<
nl <<
nl
724 if (fieldPtr_ ==
nullptr)
727 <<
"Cannot return reference from nullptr." <<
nl
728 <<
"This error message should never appear!!" <<
nl
732 return *
static_cast<const Field<Type>*
>(fieldPtr_);
740 return const_cast<Field<Type>&
>(this->cref<Type>());
748 return const_cast<Field<Type>&
>(this->cref<Type>());
762 <<
" is different from the stored result type "
763 << valType_ <<
nl <<
nl
767 Type* ptr =
dynamic_cast<Type*
>(objectPtr_.get());
784 objectPtr_.release();
792 template<
template<
class>
class BinaryOp,
class Type>
795 const BinaryOp<Type>& bop,
803 <<
" is different from the stored result type "
804 << valType_ <<
nl <<
nl
808 Type result = initial;
812 for (
const Type& val :
fld)
814 result = bop(result, val);
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
bool isObject() const
True if the object pointer is being used.
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.
Type gAverage(const FieldField< Field, Type > &f)
const word & valueType() const
Basic type for the field or single value.
tmp< Field< Type > > getResult(bool cacheCopy=false)
void setResult(Field< Type > *, bool isPointVal=false)
Set result field, taking ownership of the pointer.
bool isBool() const
True if valueType is a bool.
A polymorphic field/result from evaluating an expression.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Generic templated field type.
#define DebugInFunction
Report an information message using Foam::Info.
#define defineExpressionMethod(Type, Var)
Field< Type > & ref()
Return non-const reference to the field.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void setSingleValue(const Type &val)
Set single-value uniform result.
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
T * ptr() noexcept
Same as release().
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))
tmp< Type > getObjectResult(bool cacheCopy=false)
bool hasValue() const
Has a value?
bool isUniform() const
True if single, uniform value.
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.
Field< Type > & getRef() const
Return non-const reference to the field, casting away constness.
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.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
exprResult()
Construct null.
Traits class for primitives.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
bool isPointValue(const bool isPointVal=true) const
True if representing point values, or test if same as isPointVal.
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,...
#define WarningInFunction
Report a warning using Foam::Warning.
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
void setObjectResult(autoPtr< Type > &o)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)