scalarField.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Typedef
25  Foam::scalarField
26 
27 Description
28  Specialisation of Field<T> for scalar.
29 
30 SourceFiles
31  scalarField.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef scalarField_H
36 #define scalarField_H
37 
38 #include "Field.H"
39 #include "scalar.H"
40 
41 #define TEMPLATE
42 #include "FieldFunctionsM.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 typedef Field<scalar> scalarField;
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 template<>
54 tmp<scalarField> scalarField::component(const direction) const;
55 
56 void component
57 (
58  scalarField& sf,
59  const UList<scalar>& f,
60  const direction
61 );
62 
63 
64 template<>
65 void scalarField::replace(const direction, const UList<scalar>& sf);
66 
67 template<>
68 void scalarField::replace(const direction, const scalar& s);
69 
70 
71 void stabilise(scalarField& Res, const UList<scalar>& sf, const scalar s);
72 
73 tmp<scalarField> stabilise(const UList<scalar>&, const scalar s);
74 
75 tmp<scalarField> stabilise(const tmp<scalarField>&, const scalar s);
76 
77 
78 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79 
80 template<>
81 scalar sumProd(const UList<scalar>& f1, const UList<scalar>& f2);
82 
83 
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, add)
87 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, subtract)
88 
89 BINARY_OPERATOR(scalar, scalar, scalar, *, multiply)
90 BINARY_OPERATOR(scalar, scalar, scalar, /, divide)
91 
92 BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, divide)
93 
94 BINARY_FUNCTION(scalar, scalar, scalar, pow)
95 BINARY_TYPE_FUNCTION(scalar, scalar, scalar, pow)
96 BINARY_FUNCTION(scalar, scalar, label, integerPow)
98 BINARY_FUNCTION(scalar, scalar, label, integerRoot)
100 
101 BINARY_FUNCTION(scalar, scalar, scalar, atan2)
102 BINARY_TYPE_FUNCTION(scalar, scalar, scalar, atan2)
103 
104 
105 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
106 
107 UNARY_FUNCTION(scalar, scalar, pow3)
108 UNARY_FUNCTION(scalar, scalar, pow4)
109 UNARY_FUNCTION(scalar, scalar, pow5)
110 UNARY_FUNCTION(scalar, scalar, pow6)
111 UNARY_FUNCTION(scalar, scalar, pow025)
112 UNARY_FUNCTION(scalar, scalar, sqrt)
113 UNARY_FUNCTION(scalar, scalar, cbrt)
114 UNARY_FUNCTION(scalar, scalar, sign)
115 UNARY_FUNCTION(scalar, scalar, pos)
116 UNARY_FUNCTION(scalar, scalar, pos0)
117 UNARY_FUNCTION(scalar, scalar, neg)
118 UNARY_FUNCTION(scalar, scalar, neg0)
119 UNARY_FUNCTION(scalar, scalar, posPart)
120 UNARY_FUNCTION(scalar, scalar, negPart)
121 UNARY_FUNCTION(scalar, scalar, exp)
122 UNARY_FUNCTION(scalar, scalar, log)
123 UNARY_FUNCTION(scalar, scalar, log10)
124 UNARY_FUNCTION(scalar, scalar, sin)
125 UNARY_FUNCTION(scalar, scalar, cos)
126 UNARY_FUNCTION(scalar, scalar, tan)
127 UNARY_FUNCTION(scalar, scalar, asin)
128 UNARY_FUNCTION(scalar, scalar, acos)
129 UNARY_FUNCTION(scalar, scalar, atan)
130 UNARY_FUNCTION(scalar, scalar, sinh)
131 UNARY_FUNCTION(scalar, scalar, cosh)
132 UNARY_FUNCTION(scalar, scalar, tanh)
133 UNARY_FUNCTION(scalar, scalar, asinh)
134 UNARY_FUNCTION(scalar, scalar, acosh)
135 UNARY_FUNCTION(scalar, scalar, atanh)
136 UNARY_FUNCTION(scalar, scalar, erf)
137 UNARY_FUNCTION(scalar, scalar, erfc)
138 UNARY_FUNCTION(scalar, scalar, lgamma)
139 UNARY_FUNCTION(scalar, scalar, j0)
140 UNARY_FUNCTION(scalar, scalar, j1)
141 UNARY_FUNCTION(scalar, scalar, y0)
142 UNARY_FUNCTION(scalar, scalar, y1)
143 
144 UNARY_FUNCTION(scalar, scalar, degToRad)
145 UNARY_FUNCTION(scalar, scalar, radToDeg)
146 UNARY_FUNCTION(scalar, scalar, atmToPa)
147 UNARY_FUNCTION(scalar, scalar, paToAtm)
148 
149 #define BesselFunc(func) \
150 void func(scalarField& Res, const int n, const UList<scalar>& sf); \
151 tmp<scalarField> func(const int n, const UList<scalar>&); \
152 tmp<scalarField> func(const int n, const tmp<scalarField>&);
153 
154 BesselFunc(jn)
155 BesselFunc(yn)
156 
157 #undef BesselFunc
158 
159 
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 
162 } // End namespace Foam
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 #include "undefFieldFunctionsM.H"
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #endif
171 
172 // ************************************************************************* //
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BINARY_TYPE_FUNCTION_FS(ReturnType, Type1, Type2, Func)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BINARY_TYPE_FUNCTION(ReturnType, Type1, Type2, Func)
High performance macro functions for Field<Type> algebra. These expand using either array element acc...
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:467
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:455
volScalarField sf(fieldObject, mesh)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
scalar integerRoot(const scalar x, const label e)
Compute the power of the number x to the reciprocal integer 1/e.
Definition: scalarI.H:55
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar tan(const dimensionedScalar &ds)
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
dimensionedScalar cosh(const dimensionedScalar &ds)
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
scalar integerPow(const scalar x, const label e)
Compute the power of the number x to the integer e.
Definition: scalarI.H:30
dimensionedScalar log10(const dimensionedScalar &ds)
UNARY_FUNCTION(Type, Type, cmptMag, cmptMag)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
scalar sumProd(const UList< Type > &f1, const UList< Type > &f2)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
scalar atmToPa(const scalar atm)
Conversion from atm to Pa.
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar atan(const dimensionedScalar &ds)
scalar paToAtm(const scalar pa)
Conversion from atm to Pa.
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar posPart(const dimensionedScalar &ds)
uint8_t direction
Definition: direction.H:45
dimensionedScalar acos(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
dimensionedScalar pow025(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)
labelList f(nPoints)
#define BesselFunc(func)
Definition: scalarField.H:149