38 #define checkField(df1, df2, op) \ 39 if (&(df1).mesh() != &(df2).mesh()) \ 41 FatalErrorInFunction \ 42 << "different mesh for fields " \ 43 << (df1).name() << " and " << (df2).name() \ 44 << " during operatrion " << op \ 45 << abort(FatalError); \ 51 template<
class Type,
class GeoMesh>
65 if (field.
size() && field.
size() != GeoMesh::size(mesh))
68 <<
"size of field = " << field.
size()
69 <<
" is not the same as the size of mesh = " 70 << GeoMesh::size(mesh)
76 template<
class Type,
class GeoMesh>
82 const bool checkIOFlags
97 template<
class Type,
class GeoMesh>
103 const bool checkIOFlags
118 template<
class Type,
class GeoMesh>
127 dimensions_(df.dimensions_)
131 template<
class Type,
class GeoMesh>
141 dimensions_(df.dimensions_)
145 template<
class Type,
class GeoMesh>
154 dimensions_(move(df.dimensions_))
158 template<
class Type,
class GeoMesh>
171 dimensions_(tdf().dimensions_)
177 template<
class Type,
class GeoMesh>
187 dimensions_(df.dimensions_)
191 template<
class Type,
class GeoMesh>
202 dimensions_(df.dimensions_)
206 template<
class Type,
class GeoMesh>
216 dimensions_(df.dimensions_)
220 template<
class Type,
class GeoMesh>
231 dimensions_(df.dimensions_)
235 template<
class Type,
class GeoMesh>
249 dimensions_(tdf().dimensions_)
255 template<
class Type,
class GeoMesh>
266 template<
class Type,
class GeoMesh>
275 const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
284 mesh.thisDb().time().timeName(),
299 template<
class Type,
class GeoMesh>
308 const bool cacheTmp = mesh.thisDb().cacheTemporaryObject(name);
317 mesh.thisDb().time().timeName(),
332 template<
class Type,
class GeoMesh>
363 template<
class Type,
class GeoMesh>
371 const bool cacheTmp = tdf().db().cacheTemporaryObject(newName);
396 template<
class Type,
class GeoMesh>
399 db().cacheTemporaryObject(*
this);
405 template<
class Type,
class GeoMesh>
409 <
typename DimensionedField<Type, GeoMesh>::cmptType,
GeoMesh>
432 template<
class Type,
class GeoMesh>
437 <
typename DimensionedField<Type, GeoMesh>::cmptType,
GeoMesh>& df
444 template<
class Type,
class GeoMesh>
451 <
typename DimensionedField<Type, GeoMesh>::cmptType,
GeoMesh>
460 template<
class Type,
class GeoMesh>
480 template<
class Type,
class GeoMesh>
485 this->
name() +
".average()",
494 template<
class Type,
class GeoMesh>
504 this->
name() +
".weightedAverage(weights)",
512 template<
class Type,
class GeoMesh>
519 tweightField.clear();
526 template<
class Type,
class GeoMesh>
527 void DimensionedField<Type, GeoMesh>::operator=
536 <<
"attempted assignment to self" 547 template<
class Type,
class GeoMesh>
548 void DimensionedField<Type, GeoMesh>::operator=
557 <<
"attempted assignment to self" 568 template<
class Type,
class GeoMesh>
569 void DimensionedField<Type, GeoMesh>::operator=
580 <<
"attempted assignment to self" 592 template<
class Type,
class GeoMesh>
593 void DimensionedField<Type, GeoMesh>::operator=
603 template<
class Type,
class GeoMesh>
610 #define COMPUTED_ASSIGNMENT(TYPE, op) \ 612 template<class Type, class GeoMesh> \ 613 void DimensionedField<Type, GeoMesh>::operator op \ 615 const DimensionedField<TYPE, GeoMesh>& df \ 618 checkField(*this, df, #op); \ 620 dimensions_ op df.dimensions(); \ 621 Field<Type>::operator op(df); \ 624 template<class Type, class GeoMesh> \ 625 void DimensionedField<Type, GeoMesh>::operator op \ 627 const tmp<DimensionedField<TYPE, GeoMesh>>& tdf \ 630 operator op(tdf()); \ 634 template<class Type, class GeoMesh> \ 635 void DimensionedField<Type, GeoMesh>::operator op \ 637 const dimensioned<TYPE>& dt \ 640 dimensions_ op dt.dimensions(); \ 641 Field<Type>::operator op(dt.value()); \ 649 #undef COMPUTED_ASSIGNMENT bool cacheTemporaryObject(const word &name) const
Return true if given name is in the cacheTemporaryObjects set.
const word & name() const
Return name.
#define checkField(df1, df2, op)
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
T & ref() const
Return non-const reference or generate a fatal error.
tmp< DimensionedField< Type, GeoMesh > > clone() const
Clone.
void size(const label)
Override size to be inconsistent with allocated storage.
tmp< DimensionedField< Type, GeoMesh > > T() const
Return the field transpose (only defined for second rank tensors)
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &) const
Calculate and return weighted average.
Generic dimensioned Type class.
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
const dimensionSet & dimensions() const
Return dimensions.
Dimension set for the base types.
void operator=(const DimensionedField< Type, GeoMesh > &)
Type gSum(const FieldField< Field, Type > &f)
Pre-declare SubField and related Field type.
void clear()
Clear the list, i.e. set size to zero.
A class for handling words, derived from string.
const fileName & local() const
const Type & value() const
Return const reference to value.
virtual ~DimensionedField()
Destructor.
Foam::pointMesh ::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
errorManip< error > abort(error &err)
#define COMPUTED_ASSIGNMENT(TYPE, op)
void operator=(const Field< Type > &)
DimensionedField(const IOobject &, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Construct from components.
void replace(const direction, const DimensionedField< cmptType, GeoMesh > &)
Replace a component field of the field.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
dimensioned< Type > average() const
Calculate and return arithmetic average.
const fileName & instance() const
Type gAverage(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return const reference to dimensions.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
tmp< DimensionedField< cmptType, GeoMesh > > component(const direction) const
Return a component field of the field.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
A class for managing temporary objects.
const objectRegistry & db() const
Return the local objectRegistry.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)