scalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Description
28  Specialisation of Field<T> for scalar.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "scalarField.H"
33 #include "unitConversion.H"
34 
35 #define TEMPLATE
36 #include "FieldFunctionsM.C"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 template<>
47 {
48  return *this;
49 }
50 
51 void component(scalarField& sf, const UList<scalar>& f, const direction)
52 {
53  sf = f;
54 }
55 
56 template<>
57 void scalarField::replace(const direction, const UList<scalar>& sf)
58 {
59  *this = sf;
60 }
61 
62 template<>
63 void scalarField::replace(const direction, const scalar& s)
64 {
65  *this = s;
66 }
67 
68 
69 void stabilise(scalarField& res, const UList<scalar>& sf, const scalar s)
70 {
71  // std::transform
73  (
74  scalar, res, =, ::Foam::stabilise, scalar, s, scalar, sf
75  )
76 }
77 
78 tmp<scalarField> stabilise(const UList<scalar>& sf, const scalar s)
79 {
80  auto tresult = tmp<scalarField>::New(sf.size());
81  stabilise(tresult.ref(), sf, s);
82  return tresult;
83 }
84 
85 tmp<scalarField> stabilise(const tmp<scalarField>& tsf, const scalar s)
86 {
87  tmp<scalarField> tresult = New(tsf);
88  stabilise(tresult.ref(), tsf(), s);
89  tsf.clear();
90  return tresult;
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
95 
96 template<>
97 float sumProd(const UList<float>& f1, const UList<float>& f2)
98 {
99  float result = 0.0;
100  if (f1.size() && (f1.size() == f2.size()))
101  {
102  // std::inner_product
103  TFOR_ALL_S_OP_F_OP_F(float, result, +=, float, f1, *, float, f2)
104  }
105  return result;
106 }
107 
108 
109 template<>
110 double sumProd(const UList<double>& f1, const UList<double>& f2)
111 {
112  double result = 0.0;
113  if (f1.size() && (f1.size() == f2.size()))
114  {
115  // std::inner_product
116  TFOR_ALL_S_OP_F_OP_F(double, result, +=, double, f1, *, double, f2)
117  }
118  return result;
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123 
124 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, add)
125 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, subtract)
126 
127 BINARY_OPERATOR(scalar, scalar, scalar, *, multiply)
128 BINARY_OPERATOR(scalar, scalar, scalar, /, divide)
129 
130 BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, divide)
131 
132 BINARY_FUNCTION(scalar, scalar, scalar, pow)
133 BINARY_TYPE_FUNCTION(scalar, scalar, scalar, pow)
134 
135 BINARY_FUNCTION(scalar, scalar, scalar, atan2)
136 BINARY_TYPE_FUNCTION(scalar, scalar, scalar, atan2)
137 
138 BINARY_FUNCTION(scalar, scalar, scalar, hypot)
139 BINARY_TYPE_FUNCTION(scalar, scalar, scalar, hypot)
140 
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142 
143 UNARY_FUNCTION(scalar, scalar, pow3)
144 UNARY_FUNCTION(scalar, scalar, pow4)
145 UNARY_FUNCTION(scalar, scalar, pow5)
146 UNARY_FUNCTION(scalar, scalar, pow6)
147 UNARY_FUNCTION(scalar, scalar, pow025)
148 UNARY_FUNCTION(scalar, scalar, sqrt)
149 UNARY_FUNCTION(scalar, scalar, cbrt)
150 UNARY_FUNCTION(scalar, scalar, sign)
151 UNARY_FUNCTION(scalar, scalar, pos)
152 UNARY_FUNCTION(scalar, scalar, pos0)
153 UNARY_FUNCTION(scalar, scalar, neg)
154 UNARY_FUNCTION(scalar, scalar, neg0)
155 UNARY_FUNCTION(scalar, scalar, posPart)
156 UNARY_FUNCTION(scalar, scalar, negPart)
157 UNARY_FUNCTION(scalar, scalar, exp)
158 UNARY_FUNCTION(scalar, scalar, log)
159 UNARY_FUNCTION(scalar, scalar, log10)
160 UNARY_FUNCTION(scalar, scalar, sin)
161 UNARY_FUNCTION(scalar, scalar, cos)
162 UNARY_FUNCTION(scalar, scalar, tan)
163 UNARY_FUNCTION(scalar, scalar, asin)
164 UNARY_FUNCTION(scalar, scalar, acos)
165 UNARY_FUNCTION(scalar, scalar, atan)
166 UNARY_FUNCTION(scalar, scalar, sinh)
167 UNARY_FUNCTION(scalar, scalar, cosh)
168 UNARY_FUNCTION(scalar, scalar, tanh)
169 UNARY_FUNCTION(scalar, scalar, asinh)
170 UNARY_FUNCTION(scalar, scalar, acosh)
171 UNARY_FUNCTION(scalar, scalar, atanh)
172 UNARY_FUNCTION(scalar, scalar, erf)
173 UNARY_FUNCTION(scalar, scalar, erfc)
174 UNARY_FUNCTION(scalar, scalar, lgamma)
175 UNARY_FUNCTION(scalar, scalar, j0)
176 UNARY_FUNCTION(scalar, scalar, j1)
177 UNARY_FUNCTION(scalar, scalar, y0)
178 UNARY_FUNCTION(scalar, scalar, y1)
179 
180 UNARY_FUNCTION(scalar, scalar, degToRad)
181 UNARY_FUNCTION(scalar, scalar, radToDeg)
182 UNARY_FUNCTION(scalar, scalar, atmToPa)
183 UNARY_FUNCTION(scalar, scalar, barToPa)
184 UNARY_FUNCTION(scalar, scalar, paToAtm)
185 UNARY_FUNCTION(scalar, scalar, paToBar)
186 
187 
188 #define BesselFunc(func) \
189 void func(scalarField& res, const int n, const UList<scalar>& sf) \
190 { \
191  TFOR_ALL_F_OP_FUNC_S_F(scalar, res, =, ::Foam::func, int, n, scalar, sf) \
192 } \
193  \
194 tmp<scalarField> func(const int n, const UList<scalar>& sf) \
195 { \
196  auto tresult = tmp<scalarField>::New(sf.size()); \
197  func(tresult.ref(), n, sf); \
198  return tresult; \
199 } \
200  \
201 tmp<scalarField> func(const int n, const tmp<scalarField>& tsf) \
202 { \
203  tmp<scalarField> tresult = New(tsf); \
204  func(tresult.ref(), n, tsf()); \
205  tsf.clear(); \
206  return tresult; \
207 }
208 
209 BesselFunc(jn)
210 BesselFunc(yn)
211 
212 #undef BesselFunc
213 
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 } // End namespace Foam
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 #include "undefFieldFunctionsM.H"
222 
223 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
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)
uint8_t direction
Definition: direction.H:46
constexpr scalar barToPa(const scalar bar) noexcept
Conversion from bar to Pa.
dimensionedScalar log(const dimensionedScalar &ds)
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
Unit conversion functions.
constexpr scalar paToAtm(const scalar pa) noexcept
Conversion from Pa to atm.
#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)
#define BINARY_TYPE_FUNCTION(ReturnType, Type1, Type2, Func)
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)
#define BesselFunc(func)
Definition: scalarField.C:182
dimensionedScalar j1(const dimensionedScalar &ds)
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
tmp< scalarField > jn(const int n, const tmp< scalarField > &tsf)
Definition: scalarField.C:203
tmp< scalarField > yn(const int n, const tmp< scalarField > &tsf)
Definition: scalarField.C:204
void component(scalarField &sf, const UList< scalar > &f, const direction)
Definition: scalarField.C:45
constexpr scalar paToBar(const scalar pa) noexcept
Conversion from Pa to bar.
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
double sumProd(const UList< double > &f1, const UList< double > &f2)
Sum product for double.
Definition: scalarField.C:104
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar asinh(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
dimensionedScalar y1(const dimensionedScalar &ds)
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:608
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)
#define TFOR_ALL_S_OP_F_OP_F(typeS, s, OP1, typeF1, f1, OP2, typeF2, f2)
Definition: FieldM.H:540
labelList f(nPoints)
#define TFOR_ALL_F_OP_FUNC_S_F(typeF1, f1, OP, FUNC, typeS, s, typeF2, f2)
Definition: FieldM.H:277
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
constexpr scalar atmToPa(const scalar atm) noexcept
Conversion from atm to Pa.
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
dimensionedScalar atan(const dimensionedScalar &ds)
tmp< scalarField > stabilise(const tmp< scalarField > &tsf, const scalar s)
Definition: scalarField.C:79
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
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)
dimensionedScalar hypot(const dimensionedScalar &x, const dimensionedScalar &y)
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
dimensionedScalar log10(const dimensionedScalar &ds)
Namespace for OpenFOAM.
dimensionedScalar negPart(const dimensionedScalar &ds)