37 #define checkFieldAssignment(df1, df2) \
41 static_cast<const regIOobject*>(&df1) \
42 == static_cast<const regIOobject*>(&df2) \
45 FatalErrorInFunction \
46 << "attempted assignment to self for field " \
47 << (df1).name() << abort(FatalError); \
51 #define checkFieldOperation(df1, df2, op) \
53 if (&(df1).mesh() != &(df2).mesh()) \
55 FatalErrorInFunction \
56 << "different mesh for fields " \
57 << (df1).name() << " and " << (df2).name() \
58 << " during operation " << op \
59 << abort(FatalError); \
65 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
71 const PrimitiveField<Type>& field
75 PrimitiveField<Type>(field),
80 if (field.size() && field.size() != GeoMesh::size(
mesh))
83 <<
"size of field = " << field.size()
84 <<
" is not the same as the size of mesh = "
85 << GeoMesh::size(
mesh)
91 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
97 const bool checkIOFlags
116 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
122 const bool checkIOFlags
126 PrimitiveField<Type>(
GeoMesh::size(
mesh), dt.value()),
138 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
145 PrimitiveField<Type>(df),
148 dimensions_(df.dimensions_)
152 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
159 PrimitiveField<Type>(move(df)),
162 dimensions_(move(df.dimensions_))
166 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
167 template<
template<
class>
class PrimitiveField2>
174 PrimitiveField<Type>(df),
181 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
182 template<
template<
class>
class PrimitiveField2>
190 PrimitiveField<Type>(df, reuse),
197 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
211 dimensions_(tdf().dimensions_)
217 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
218 template<
template<
class>
class PrimitiveField2>
223 const bool checkIOFlags
227 PrimitiveField<Type>(df),
230 dimensions_(df.dimensions_)
232 if (!checkIOFlags || !readIfPresent())
234 copyOldTimes(io, df);
239 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
240 template<
template<
class>
class PrimitiveField2>
246 const bool checkIOFlags
250 PrimitiveField<Type>(df, reuse),
253 dimensions_(df.dimensions_)
262 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
267 const bool checkIOFlags
278 dimensions_(tdf().dimensions_)
289 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
290 template<
template<
class>
class PrimitiveField2>
298 PrimitiveField<Type>(df),
301 dimensions_(df.dimensions_)
307 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
308 template<
template<
class>
class PrimitiveField2>
317 PrimitiveField<Type>(df, reuse),
320 dimensions_(df.dimensions_)
324 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
339 dimensions_(tdf().dimensions_)
345 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
356 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
363 const PrimitiveField<Type>& field
390 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
423 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
456 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
457 template<
template<
class>
class PrimitiveField2>
488 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
496 const bool cacheTmp = tdf().db().cacheTemporaryObject(newName);
521 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
524 db().cacheTemporaryObject(*
this);
530 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
531 PrimitiveField<Type>&
540 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
571 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
572 template<
template<
class>
class PrimitiveField2>
584 PrimitiveField<Type>::replace(d, df);
588 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
589 template<
template<
class>
class PrimitiveField2>
609 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
629 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
635 this->
name() +
".average()",
644 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
645 template<
template<
class>
class PrimitiveField2>
656 this->
name() +
".weightedAverage(weights)",
658 gSum(weightField*primitiveField())/
gSum(weightField)
664 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
665 template<
template<
class>
class PrimitiveField2>
673 tweightField.clear();
678 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
679 template<
template<
class>
class PrimitiveField2>
692 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
706 PrimitiveField<Type>::transfer(tdf.ref());
717 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
718 template<
template<
class>
class PrimitiveField2>
732 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
742 PrimitiveField<Type>::operator=(df);
746 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
756 PrimitiveField<Type>::operator=(move(df));
760 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
786 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
787 template<
template<
class>
class PrimitiveField2>
796 PrimitiveField<Type>::operator=(df);
800 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
801 template<
template<
class>
class PrimitiveField2>
812 PrimitiveField<Type>::operator=(df);
817 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
824 PrimitiveField<Type>::operator=(dt.
value());
828 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
831 PrimitiveField<Type>::operator=(
Zero);
835 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
836 template<
template<
class>
class PrimitiveField2>
846 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
847 template<
template<
class>
class PrimitiveField2>
853 this->operator=(tdf);
857 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
867 template<
class Type,
class GeoMesh,
template<
class>
class PrimitiveField>
870 this->operator=(
Zero);
874 #define COMPUTED_ASSIGNMENT(TYPE, op) \
876 template<class Type, class GeoMesh, template<class> class PrimitiveField> \
877 template<template<class> class PrimitiveField2> \
878 void DimensionedField<Type, GeoMesh, PrimitiveField>::operator op \
880 const DimensionedField<TYPE, GeoMesh, PrimitiveField2>& df \
883 checkFieldOperation(*this, df, #op); \
885 dimensions_ op df.dimensions(); \
886 PrimitiveField<Type>::operator op(df); \
889 template<class Type, class GeoMesh, template<class> class PrimitiveField> \
890 template<template<class> class PrimitiveField2> \
891 void DimensionedField<Type, GeoMesh, PrimitiveField>::operator op \
893 const tmp<DimensionedField<TYPE, GeoMesh, PrimitiveField2>>& tdf \
896 operator op(tdf()); \
900 template<class Type, class GeoMesh, template<class> class PrimitiveField> \
901 void DimensionedField<Type, GeoMesh, PrimitiveField>::operator op \
903 const dimensioned<TYPE>& dt \
906 dimensions_ op dt.dimensions(); \
907 PrimitiveField<Type>::operator op(dt.value()); \
915 #undef COMPUTED_ASSIGNMENT
920 #undef checkFieldAssignment
921 #undef checkFieldOperation
#define checkFieldOperation(df1, df2, op)
#define COMPUTED_ASSIGNMENT(TYPE, op)
#define checkFieldAssignment(df1, df2)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
PrimitiveField< Type >::cmptType cmptType
Component type of the elements of the field.
friend class DimensionedField
Declare friendship with other dimensioned fields.
GeoMesh::Mesh Mesh
Type of mesh on which this DimensionedField is instantiated.
const dimensionSet & dimensions() const
Return dimensions.
PrimitiveField< Type > & primitiveFieldRef()
Return a reference to the internal field.
const Mesh & mesh() const
Return mesh.
const PrimitiveField< Type > & primitiveField() const
Return a const-reference to the primitive field.
Pre-declare SubField and related Field type.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const fileName & local() const
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
const objectRegistry & db() const
Return the local objectRegistry.
Class to add into field types to provide old-time storage and retrieval.
void copyOldTimes(const IOobject &io, const OtherOldTime< OtherPrimitiveField > &)
Copy the old-times from the given field.
Dimension set for the base types.
Generic dimensioned Type class.
const dimensionSet & dimensions() const
Return const reference to dimensions.
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
const Time & time() const
Return time.
bool cacheTemporaryObject(const word &name) const
Return true if given name is in the cacheTemporaryObjects set.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
A class for handling words, derived from string.
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
errorManip< error > abort(error &err)
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
Type gAverage(const FieldField< Field, Type > &f)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
points setSize(newPointi)
conserve primitiveFieldRef()+