28 #define TEMPLATE template<class GeoMesh, template<class> class PrimitiveField>
33 template<class> class PrimitiveField1, \
34 template<class> class PrimitiveField2 \
45 template<
class GeoMesh,
template<
class>
class PrimitiveField>
56 "stabilise(" + dsf.
name() +
',' + ds.
name() +
')',
68 template<
class GeoMesh,
template<
class>
class PrimitiveField>
80 "stabilise(" + dsf.
name() +
',' + ds.
name() +
')',
104 template<
class GeoMesh,
template<
class>
class PrimitiveField>
114 <<
"Base field is not dimensionless: " << dsf1.
dimensions()
121 <<
"Exponent field is not dimensionless: " << dsf2.
dimensions()
129 "pow(" + dsf1.
name() +
',' + dsf2.
name() +
')',
137 tPow.
ref().primitiveFieldRef(),
149 template<
class>
class PrimitiveField1,
150 template<
class>
class PrimitiveField2
163 <<
"Base field is not dimensionless: " << dsf1.
dimensions()
170 <<
"Exponent field is not dimensionless: " << dsf2.
dimensions()
177 "pow(" + dsf1.
name() +
',' + dsf2.
name() +
')',
183 tPow.
ref().primitiveFieldRef(),
197 template<
class>
class PrimitiveField1,
198 template<
class>
class PrimitiveField2
211 <<
"Base field is not dimensionless: " << dsf1.
dimensions()
218 <<
"Exponent field is not dimensionless: " << dsf2.
dimensions()
225 "pow(" + dsf1.
name() +
',' + dsf2.
name() +
')',
231 tPow.
ref().primitiveFieldRef(),
245 template<
class>
class PrimitiveField1,
246 template<
class>
class PrimitiveField2
260 <<
"Base field is not dimensionless: " << dsf1.
dimensions()
267 <<
"Exponent field is not dimensionless: " << dsf2.
dimensions()
284 "pow(" + dsf1.
name() +
',' + dsf2.
name() +
')',
290 tPow.
ref().primitiveFieldRef(),
302 template<
class GeoMesh,
template<
class>
class PrimitiveField>
312 <<
"Exponent is not dimensionless: " << ds.
dimensions()
320 "pow(" + dsf.
name() +
',' + ds.
name() +
')',
332 template<
class GeoMesh,
template<
class>
class PrimitiveField>
342 <<
"Exponent is not dimensionless: " << ds.
dimensions()
351 "pow(" + dsf.
name() +
',' + ds.
name() +
')',
363 template<
class GeoMesh,
template<
class>
class PrimitiveField>
374 template<
class GeoMesh,
template<
class>
class PrimitiveField>
385 template<
class GeoMesh,
template<
class>
class PrimitiveField>
395 <<
"Base scalar is not dimensionless: " << ds.
dimensions()
402 <<
"Exponent field is not dimensionless: " << dsf.
dimensions()
410 "pow(" + ds.
name() +
',' + dsf.
name() +
')',
422 template<
class GeoMesh,
template<
class>
class PrimitiveField>
434 <<
"Base scalar is not dimensionless: " << ds.
dimensions()
441 <<
"Exponent field is not dimensionless: " << dsf.
dimensions()
448 "pow(" + ds.
name() +
',' + dsf.
name() +
')',
459 template<
class GeoMesh,
template<
class>
class PrimitiveField>
469 template<
class GeoMesh,
template<
class>
class PrimitiveField>
485 template<
class>
class PrimitiveField1,
486 template<
class>
class PrimitiveField2
498 "atan2(" + dsf1.
name() +
',' + dsf2.
name() +
')',
506 tAtan2.
ref().primitiveFieldRef(),
518 template<
class>
class PrimitiveField1,
519 template<
class>
class PrimitiveField2
532 "atan2(" + dsf1.
name() +
',' + dsf2.
name() +
')',
538 tAtan2.
ref().primitiveFieldRef(),
552 template<
class>
class PrimitiveField1,
553 template<
class>
class PrimitiveField2
566 "atan2(" + dsf1.
name() +
',' + dsf2.
name() +
')',
572 tAtan2.
ref().primitiveFieldRef(),
585 template<
class>
class PrimitiveField1,
586 template<
class>
class PrimitiveField2
610 "atan2(" + dsf1.
name() +
',' + dsf2.
name() +
')',
616 tAtan2.
ref().primitiveFieldRef(),
628 template<
class GeoMesh,
template<
class>
class PrimitiveField>
639 "atan2(" + dsf.
name() +
',' + ds.
name() +
')',
650 template<
class GeoMesh,
template<
class>
class PrimitiveField>
662 "atan2(" + dsf.
name() +
',' + ds.
name() +
')',
673 template<
class GeoMesh,
template<
class>
class PrimitiveField>
683 template<
class GeoMesh,
template<
class>
class PrimitiveField>
694 template<
class GeoMesh,
template<
class>
class PrimitiveField>
705 "atan2(" + ds.
name() +
',' + dsf.
name() +
')',
717 template<
class GeoMesh,
template<
class>
class PrimitiveField>
729 "atan2(" + ds.
name() +
',' + dsf.
name() +
')',
740 template<
class GeoMesh,
template<
class>
class PrimitiveField>
750 template<
class GeoMesh,
template<
class>
class PrimitiveField>
804 #define BesselFunc(func) \
806 template<class GeoMesh, template<class> class PrimitiveField> \
807 tmp<DimensionedField<scalar, GeoMesh, Field>> func \
810 const DimensionedField<scalar, GeoMesh, PrimitiveField>& dsf \
813 if (!dsf.dimensions().dimensionless()) \
815 FatalErrorInFunction \
816 << "dsf not dimensionless" \
817 << abort(FatalError); \
820 tmp<DimensionedField<scalar, GeoMesh, Field>> tFunc \
822 DimensionedField<scalar, GeoMesh, Field>::New \
824 #func "(" + name(n) + ',' + dsf.name() + ')', \
830 func(tFunc.ref().primitiveFieldRef(), n, dsf.primitiveField()); \
835 template<class GeoMesh, template<class> class PrimitiveField> \
836 tmp<DimensionedField<scalar, GeoMesh, Field>> func \
839 const tmp<DimensionedField<scalar, GeoMesh, PrimitiveField>>& tdsf \
842 const DimensionedField<scalar, GeoMesh, PrimitiveField>& dsf = tdsf(); \
844 if (!dsf.dimensions().dimensionless()) \
846 FatalErrorInFunction \
847 << " : dsf not dimensionless" \
848 << abort(FatalError); \
851 tmp<DimensionedField<scalar, GeoMesh, Field>> tFunc \
856 #func "(" + name(n) + ',' + dsf.name() + ')', \
861 func(tFunc.ref().primitiveFieldRef(), n, dsf.primitiveField()); \
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
Scalar specific part of the implementation of DimensionedField.
#define BINARY_OPERATOR(Template, Type, Type1, Type2, op, opFunc)
#define UNARY_FUNCTION(Template, Type, Type1, func)
#define BINARY_TYPE_OPERATOR_SF(TYPE, op, opFunc)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
const PrimitiveField< Type > & primitiveField() const
Return a const-reference to the primitive field.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
const word & name() const
Return name.
bool dimensionless() const
Return true if it is dimensionless.
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.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void subtract(LagrangianPatchField< typename typeOfSum< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar tan(const dimensionedScalar &ds)
void divide(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1, const LagrangianPatchField< scalar > &f2)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
void pow025(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
const dimensionSet dimless
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensionedScalar log10(const dimensionedScalar &ds)
tmp< DimensionedField< scalar, GeoMesh, Field > > stabilise(const DimensionedField< scalar, GeoMesh, PrimitiveField > &dsf, const dimensioned< scalar > &ds)
void pow4(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void pow6(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
void multiply(LagrangianPatchField< Type > &f, const LagrangianPatchField< scalar > &f1, const LagrangianPatchField< Type > &f2)
void pow5(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar neg(const dimensionedScalar &ds)
void cbrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar atanh(const dimensionedScalar &ds)
dimensionedScalar neg0(const dimensionedScalar &ds)
void add(LagrangianPatchField< typename typeOfSum< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)
dimensionedScalar atan(const dimensionedScalar &ds)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionSet trans(const dimensionSet &)
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar j0(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)