31 #define TEMPLATE template<template<class> class PatchField, class GeoMesh> 41 template<
template<
class>
class PatchField,
class GeoMesh>
54 template<
template<
class>
class PatchField,
class GeoMesh>
55 tmp<GeometricField<scalar, PatchField, GeoMesh>>
stabilise 57 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
58 const dimensioned<scalar>& ds
63 "stabilise(" + gsf.name() +
',' + ds.name() +
')',
65 ds.dimensions() + gsf.dimensions()
74 template<
template<
class>
class PatchField,
class GeoMesh>
75 tmp<GeometricField<scalar, PatchField, GeoMesh>>
stabilise 77 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
78 const dimensioned<scalar>& ds
81 const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf();
83 tmp<GeometricField<scalar, PatchField, GeoMesh>> tRes
88 "stabilise(" + gsf.name() +
',' + ds.name() +
')',
89 ds.dimensions() + gsf.dimensions()
114 template<template<class> class PatchField, class
GeoMesh>
122 pow(Pow.primitiveFieldRef(), gsf1.primitiveField(), gsf2.primitiveField());
123 pow(Pow.boundaryFieldRef(), gsf1.boundaryField(), gsf2.boundaryField());
127 template<
template<
class>
class PatchField,
class GeoMesh>
128 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 130 const GeometricField<scalar, PatchField, GeoMesh>& f1,
131 const GeometricField<scalar, PatchField, GeoMesh>& f2
137 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
141 <<
"pow() expects dimensionless parameters, but found" <<
nl;
143 if (!f1.dimensions().dimensionless())
146 <<
" Base field dimensions: " << f1.dimensions() <<
nl;
148 if (!f2.dimensions().dimensionless())
151 <<
" Exponent field dimensions: " << f2.dimensions() <<
nl;
158 "pow(" + f1.name() +
',' + f2.name() +
')',
163 pow(tresult.ref(), f1, f2);
169 template<
template<
class>
class PatchField,
class GeoMesh>
170 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 172 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf1,
173 const GeometricField<scalar, PatchField, GeoMesh>& f2
176 const auto& f1 = tf1();
181 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
185 <<
"pow() expects dimensionless parameters, but found" <<
nl;
187 if (!f1.dimensions().dimensionless())
190 <<
" Base field dimensions: " << f1.dimensions() <<
nl;
192 if (!f2.dimensions().dimensionless())
195 <<
" Exponent field dimensions: " << f2.dimensions() <<
nl;
200 tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
205 "pow(" + f1.name() +
',' + f2.name() +
')',
210 pow(tresult.ref(), f1, f2);
218 template<
template<
class>
class PatchField,
class GeoMesh>
219 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 221 const GeometricField<scalar, PatchField, GeoMesh>& f1,
222 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf2
225 const auto& f2 = tf2();
230 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
234 <<
"pow() expects dimensionless parameters, but found" <<
nl;
236 if (!f1.dimensions().dimensionless())
239 <<
" Base field dimensions: " << f1.dimensions() <<
nl;
241 if (!f2.dimensions().dimensionless())
244 <<
" Exponent field dimensions: " << f2.dimensions() <<
nl;
249 tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
254 "pow(" + f1.name() +
',' + f2.name() +
')',
259 pow(tresult.ref(), f1, f2);
266 template<
template<
class>
class PatchField,
class GeoMesh>
267 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 269 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf1,
270 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf2
273 const auto& f1 = tf1();
274 const auto& f2 = tf2();
279 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
283 <<
"pow() expects dimensionless parameters, but found" <<
nl;
285 if (!f1.dimensions().dimensionless())
288 <<
" Base field dimensions: " << f1.dimensions() <<
nl;
290 if (!f2.dimensions().dimensionless())
293 <<
" Exponent field dimensions: " << f2.dimensions() <<
nl;
298 tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
300 reuseTmpTmpGeometricField
301 <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::
New 305 "pow(" + f1.name() +
',' + f2.name() +
')',
310 pow(tresult.ref(), f1, f2);
319 template<
template<
class>
class PatchField,
class GeoMesh>
322 GeometricField<scalar, PatchField, GeoMesh>& tPow,
323 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
324 const dimensioned<scalar>& ds
327 pow(tPow.primitiveFieldRef(), gsf.primitiveField(), ds.value());
328 pow(tPow.boundaryFieldRef(), gsf.boundaryField(), ds.value());
332 template<
template<
class>
class PatchField,
class GeoMesh>
333 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 335 const GeometricField<scalar, PatchField, GeoMesh>& f1,
342 && (!ds.dimensions().dimensionless())
346 <<
"pow() expects dimensionless parameters, but found" <<
nl 347 <<
" Exponent dimensions: " << ds.dimensions() <<
nl 353 "pow(" + f1.name() +
',' + ds.name() +
')',
355 pow(f1.dimensions(), ds)
358 pow(tresult.ref(), f1, ds);
363 template<
template<
class>
class PatchField,
class GeoMesh>
364 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 366 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf1,
370 const auto& f1 = tf1();
375 && (!ds.dimensions().dimensionless())
379 <<
"pow() expects dimensionless parameters, but found" <<
nl 380 <<
" Exponent dimensions: " << ds.dimensions() <<
nl 384 tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
389 "pow(" + f1.name() +
',' + ds.name() +
')',
390 pow(f1.dimensions(), ds)
394 pow(tresult.ref(), f1, ds);
401 template<
template<
class>
class PatchField,
class GeoMesh>
402 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 404 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
411 template<
template<
class>
class PatchField,
class GeoMesh>
412 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 414 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
422 template<
template<
class>
class PatchField,
class GeoMesh>
425 GeometricField<scalar, PatchField, GeoMesh>& tPow,
426 const dimensioned<scalar>& ds,
427 const GeometricField<scalar, PatchField, GeoMesh>& gsf
430 pow(tPow.primitiveFieldRef(), ds.value(), gsf.primitiveField());
431 pow(tPow.boundaryFieldRef(), ds.value(), gsf.boundaryField());
435 template<
template<
class>
class PatchField,
class GeoMesh>
436 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 439 const GeometricField<scalar, PatchField, GeoMesh>& f2
445 && (!ds.dimensions().dimensionless() || !f2.dimensions().dimensionless())
449 <<
"pow() expects dimensionless parameters, but found" <<
nl;
451 if (!ds.dimensions().dimensionless())
454 <<
" Base scalar dimensions: " << ds.dimensions() <<
nl;
456 if (!f2.dimensions().dimensionless())
459 <<
" Exponent field dimensions: " << f2.dimensions() <<
nl;
466 "pow(" + ds.name() +
',' + f2.name() +
')',
471 pow(tresult.ref(), ds, f2);
477 template<
template<
class>
class PatchField,
class GeoMesh>
478 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 481 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf2
484 const auto& f2 = tf2();
489 && (!ds.dimensions().dimensionless() || !f2.dimensions().dimensionless())
493 <<
"pow() expects dimensionless parameters, but found" <<
nl;
495 if (!ds.dimensions().dimensionless())
498 <<
" Base scalar dimensions: " << ds.dimensions() <<
nl;
500 if (!f2.dimensions().dimensionless())
503 <<
" Exponent field dimensions: " << f2.dimensions() <<
nl;
508 tmp<GeometricField<scalar, PatchField, GeoMesh>> tresult
513 "pow(" + ds.name() +
',' + f2.name() +
')',
518 pow(tresult.ref(), ds, f2);
525 template<
template<
class>
class PatchField,
class GeoMesh>
526 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 529 const GeometricField<scalar, PatchField, GeoMesh>& gsf
535 template<
template<
class>
class PatchField,
class GeoMesh>
536 tmp<GeometricField<scalar, PatchField, GeoMesh>>
pow 539 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf
548 template<
template<
class>
class PatchField,
class GeoMesh>
551 GeometricField<scalar, PatchField, GeoMesh>& Atan2,
552 const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
553 const GeometricField<scalar, PatchField, GeoMesh>& gsf2
558 Atan2.primitiveFieldRef(),
559 gsf1.primitiveField(),
560 gsf2.primitiveField()
564 Atan2.boundaryFieldRef(),
565 gsf1.boundaryField(),
571 template<
template<
class>
class PatchField,
class GeoMesh>
572 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 574 const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
575 const GeometricField<scalar, PatchField, GeoMesh>& gsf2
580 "atan2(" + gsf1.name() +
',' + gsf2.name() +
')',
582 atan2(gsf1.dimensions(), gsf2.dimensions())
591 template<
template<
class>
class PatchField,
class GeoMesh>
592 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 594 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf1,
595 const GeometricField<scalar, PatchField, GeoMesh>& gsf2
598 const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
600 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
605 "atan2(" + gsf1.name() +
',' + gsf2.name() +
')',
606 atan2(gsf1.dimensions(), gsf2.dimensions())
618 template<
template<
class>
class PatchField,
class GeoMesh>
619 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 621 const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
622 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf2
625 const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
627 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
632 "atan2(" + gsf1.name() +
',' + gsf2.name() +
')',
633 atan2( gsf1.dimensions(), gsf2.dimensions())
644 template<
template<
class>
class PatchField,
class GeoMesh>
645 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 647 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf1,
648 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf2
651 const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
652 const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
654 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
656 reuseTmpTmpGeometricField
657 <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::
New 661 "atan2(" + gsf1.name() +
',' + gsf2.name() +
')',
662 atan2(gsf1.dimensions(), gsf2.dimensions())
675 template<
template<
class>
class PatchField,
class GeoMesh>
678 GeometricField<scalar, PatchField, GeoMesh>& tAtan2,
679 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
680 const dimensioned<scalar>& ds
683 atan2(tAtan2.primitiveFieldRef(), gsf.primitiveField(), ds.value());
684 atan2(tAtan2.boundaryFieldRef(), gsf.boundaryField(), ds.value());
688 template<
template<
class>
class PatchField,
class GeoMesh>
689 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 691 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
697 "atan2(" + gsf.name() +
',' + ds.name() +
')',
699 atan2(gsf.dimensions(), ds)
707 template<
template<
class>
class PatchField,
class GeoMesh>
708 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 710 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
714 const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf();
716 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
721 "atan2(" + gsf.name() +
',' + ds.name() +
')',
722 atan2(gsf.dimensions(), ds)
726 atan2(tAtan2.ref(), gsf, ds);
733 template<
template<
class>
class PatchField,
class GeoMesh>
734 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 736 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
743 template<
template<
class>
class PatchField,
class GeoMesh>
744 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 746 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
754 template<
template<
class>
class PatchField,
class GeoMesh>
757 GeometricField<scalar, PatchField, GeoMesh>& tAtan2,
758 const dimensioned<scalar>& ds,
759 const GeometricField<scalar, PatchField, GeoMesh>& gsf
762 atan2(tAtan2.primitiveFieldRef(), ds.value(), gsf.primitiveField());
763 atan2(tAtan2.boundaryFieldRef(), ds.value(), gsf.boundaryField());
767 template<
template<
class>
class PatchField,
class GeoMesh>
768 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 771 const GeometricField<scalar, PatchField, GeoMesh>& gsf
776 "atan2(" + ds.name() +
',' + gsf.name() +
')',
778 atan2(ds, gsf.dimensions())
781 atan2(tAtan2.ref(), ds, gsf);
787 template<
template<
class>
class PatchField,
class GeoMesh>
788 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 791 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf
794 const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf();
796 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
801 "atan2(" + ds.name() +
',' + gsf.name() +
')',
802 atan2(ds, gsf.dimensions())
806 atan2(tAtan2.ref(), ds, gsf);
813 template<
template<
class>
class PatchField,
class GeoMesh>
814 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 817 const GeometricField<scalar, PatchField, GeoMesh>& gsf
823 template<
template<
class>
class PatchField,
class GeoMesh>
824 tmp<GeometricField<scalar, PatchField, GeoMesh>>
atan2 827 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf
877 #define BesselFunc(func) \ 879 template<template<class> class PatchField, class GeoMesh> \ 882 GeometricField<scalar, PatchField, GeoMesh>& gsf, \ 884 const GeometricField<scalar, PatchField, GeoMesh>& gsf1 \ 887 func(gsf.primitiveFieldRef(), n, gsf1.primitiveField()); \ 888 func(gsf.boundaryFieldRef(), n, gsf1.boundaryField()); \ 891 template<template<class> class PatchField, class GeoMesh> \ 892 tmp<GeometricField<scalar, PatchField, GeoMesh>> func \ 895 const GeometricField<scalar, PatchField, GeoMesh>& gsf \ 898 if (dimensionSet::checking() && !gsf.dimensions().dimensionless()) \ 900 FatalErrorInFunction \ 901 << "Field is not dimensionless: " << gsf.dimensions() << nl \ 902 << abort(FatalError); \ 905 auto tFunc = GeometricField<scalar, PatchField, GeoMesh>::New \ 907 #func "(" + gsf.name() + ')', \ 912 func(tFunc.ref(), n, gsf); \ 917 template<template<class> class PatchField, class GeoMesh> \ 918 tmp<GeometricField<scalar, PatchField, GeoMesh>> func \ 921 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf \ 924 const auto& gsf = tgsf(); \ 926 if (dimensionSet::checking() && !gsf.dimensions().dimensionless()) \ 928 FatalErrorInFunction \ 929 << "Field is not dimensionless: " << gsf.dimensions() << nl \ 930 << abort(FatalError); \ 933 tmp<GeometricField<scalar, PatchField, GeoMesh>> tFunc \ 938 #func "(" + gsf.name() + ')', \ 943 func(tFunc.ref(), n, gsf); \ dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
const Type & value() const noexcept
Return const reference to value.
dimensionedScalar acos(const dimensionedScalar &ds)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
constexpr char nl
The newline '\n' character (0x0a)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar y0(const dimensionedScalar &ds)
Generic GeometricField class.
const dimensionSet dimless
Dimensionless.
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)
Scalar specific part of the implementation of GeometricField.
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionSet trans(const dimensionSet &ds)
Check the argument is dimensionless (for transcendental functions)
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
static bool checking() noexcept
True if dimension checking is enabled (the usual default)
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar atan(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
dimensionedScalar pow6(const dimensionedScalar &ds)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
dimensionedScalar cosh(const dimensionedScalar &ds)
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=PatchField< Type >::calculatedType())
Return tmp field from name, mesh, dimensions and patch type.
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
dimensionedScalar log10(const dimensionedScalar &ds)
dimensionedScalar negPart(const dimensionedScalar &ds)