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-------------------------------------------------------------------------------
11License
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
27Description
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
40namespace Foam
41{
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45template<>
46tmp<scalarField> scalarField::component(const direction) const
47{
48 return *this;
49}
50
52{
53 sf = f;
54}
55
56template<>
57void scalarField::replace(const direction, const UList<scalar>& sf)
58{
59 *this = sf;
60}
61
62template<>
63void scalarField::replace(const direction, const scalar& s)
64{
65 *this = s;
66}
67
68
69void stabilise(scalarField& res, const UList<scalar>& sf, const scalar s)
70{
72 (
73 scalar, res, =, ::Foam::stabilise, scalar, s, scalar, sf
74 )
75}
76
77tmp<scalarField> stabilise(const UList<scalar>& sf, const scalar s)
78{
79 auto tresult = tmp<scalarField>::New(sf.size());
80 stabilise(tresult.ref(), sf, s);
81 return tresult;
82}
83
85{
86 tmp<scalarField> tresult = New(tsf);
87 stabilise(tresult.ref(), tsf(), s);
88 tsf.clear();
89 return tresult;
90}
91
92
93// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
94
95template<>
96float sumProd(const UList<float>& f1, const UList<float>& f2)
97{
98 float result = 0.0;
99 if (f1.size() && (f1.size() == f2.size()))
100 {
101 TFOR_ALL_S_OP_F_OP_F(float, result, +=, float, f1, *, float, f2)
102 }
103 return result;
104}
105
106
107template<>
108double sumProd(const UList<double>& f1, const UList<double>& f2)
109{
110 double result = 0.0;
111 if (f1.size() && (f1.size() == f2.size()))
112 {
113 TFOR_ALL_S_OP_F_OP_F(double, result, +=, double, f1, *, double, f2)
114 }
115 return result;
116}
117
118
119// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120
121BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, add)
122BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, subtract)
123
124BINARY_OPERATOR(scalar, scalar, scalar, *, multiply)
125BINARY_OPERATOR(scalar, scalar, scalar, /, divide)
126
127BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, divide)
128
129BINARY_FUNCTION(scalar, scalar, scalar, pow)
130BINARY_TYPE_FUNCTION(scalar, scalar, scalar, pow)
131
132BINARY_FUNCTION(scalar, scalar, scalar, atan2)
133BINARY_TYPE_FUNCTION(scalar, scalar, scalar, atan2)
134
135BINARY_FUNCTION(scalar, scalar, scalar, hypot)
136BINARY_TYPE_FUNCTION(scalar, scalar, scalar, hypot)
137
138// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139
140UNARY_FUNCTION(scalar, scalar, pow3)
141UNARY_FUNCTION(scalar, scalar, pow4)
142UNARY_FUNCTION(scalar, scalar, pow5)
143UNARY_FUNCTION(scalar, scalar, pow6)
144UNARY_FUNCTION(scalar, scalar, pow025)
145UNARY_FUNCTION(scalar, scalar, sqrt)
146UNARY_FUNCTION(scalar, scalar, cbrt)
147UNARY_FUNCTION(scalar, scalar, sign)
148UNARY_FUNCTION(scalar, scalar, pos)
149UNARY_FUNCTION(scalar, scalar, pos0)
150UNARY_FUNCTION(scalar, scalar, neg)
151UNARY_FUNCTION(scalar, scalar, neg0)
152UNARY_FUNCTION(scalar, scalar, posPart)
153UNARY_FUNCTION(scalar, scalar, negPart)
154UNARY_FUNCTION(scalar, scalar, exp)
155UNARY_FUNCTION(scalar, scalar, log)
156UNARY_FUNCTION(scalar, scalar, log10)
157UNARY_FUNCTION(scalar, scalar, sin)
158UNARY_FUNCTION(scalar, scalar, cos)
159UNARY_FUNCTION(scalar, scalar, tan)
160UNARY_FUNCTION(scalar, scalar, asin)
161UNARY_FUNCTION(scalar, scalar, acos)
162UNARY_FUNCTION(scalar, scalar, atan)
163UNARY_FUNCTION(scalar, scalar, sinh)
164UNARY_FUNCTION(scalar, scalar, cosh)
165UNARY_FUNCTION(scalar, scalar, tanh)
166UNARY_FUNCTION(scalar, scalar, asinh)
167UNARY_FUNCTION(scalar, scalar, acosh)
168UNARY_FUNCTION(scalar, scalar, atanh)
169UNARY_FUNCTION(scalar, scalar, erf)
170UNARY_FUNCTION(scalar, scalar, erfc)
171UNARY_FUNCTION(scalar, scalar, lgamma)
172UNARY_FUNCTION(scalar, scalar, j0)
173UNARY_FUNCTION(scalar, scalar, j1)
174UNARY_FUNCTION(scalar, scalar, y0)
175UNARY_FUNCTION(scalar, scalar, y1)
176
177UNARY_FUNCTION(scalar, scalar, degToRad)
178UNARY_FUNCTION(scalar, scalar, radToDeg)
179UNARY_FUNCTION(scalar, scalar, atmToPa)
180UNARY_FUNCTION(scalar, scalar, barToPa)
181UNARY_FUNCTION(scalar, scalar, paToAtm)
182UNARY_FUNCTION(scalar, scalar, paToBar)
183
184
185#define BesselFunc(func) \
186void func(scalarField& res, const int n, const UList<scalar>& sf) \
187{ \
188 TFOR_ALL_F_OP_FUNC_S_F(scalar, res, =, ::Foam::func, int, n, scalar, sf) \
189} \
190 \
191tmp<scalarField> func(const int n, const UList<scalar>& sf) \
192{ \
193 auto tresult = tmp<scalarField>::New(sf.size()); \
194 func(tresult.ref(), n, sf); \
195 return tresult; \
196} \
197 \
198tmp<scalarField> func(const int n, const tmp<scalarField>& tsf) \
199{ \
200 tmp<scalarField> tresult = New(tsf); \
201 func(tresult.ref(), n, tsf()); \
202 tsf.clear(); \
203 return tresult; \
204}
205
208
209#undef BesselFunc
210
211
212// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213
214} // End namespace Foam
215
216// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217
218#include "undefFieldFunctionsM.H"
219
220// ************************************************************************* //
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
#define BINARY_TYPE_FUNCTION(ReturnType, Type1, Type2, Func)
#define TFOR_ALL_S_OP_F_OP_F(typeS, s, OP1, typeF1, f1, OP2, typeF2, f2)
Definition: FieldM.H:389
#define TFOR_ALL_F_OP_FUNC_S_F(typeF1, f1, OP, FUNC, typeS, s, typeF2, f2)
Definition: FieldM.H:219
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
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))
Namespace for OpenFOAM.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
complex sumProd(const UList< complex > &f1, const UList< complex > &f2)
Sum product.
Definition: complexField.C:191
uint8_t direction
Definition: direction.H:56
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
labelList f(nPoints)
#define BesselFunc(func)
Definition: scalarField.C:185
dict add("bounds", meshBb)
Unit conversion functions.