Loading...
Searching...
No Matches
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.
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<>
47{
48 return *this;
49}
50
51void component(scalarField& sf, const UList<scalar>& f, const direction)
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{
71 // std::transform
73 (
74 scalar, res, =, ::Foam::stabilise, scalar, s, scalar, sf
75 )
76}
77
78tmp<scalarField> stabilise(const UList<scalar>& sf, const scalar s)
80 auto tresult = tmp<scalarField>::New(sf.size());
81 stabilise(tresult.ref(), sf, s);
82 return tresult;
83}
84
85tmp<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;
92
93
94// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
95
96template<>
97float 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)
105 return result;
106}
107
108
109template<>
110double 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
124BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, add)
125BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, subtract)
126
127BINARY_OPERATOR(scalar, scalar, scalar, *, multiply)
128BINARY_OPERATOR(scalar, scalar, scalar, /, divide)
129
130BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, divide)
131
132BINARY_FUNCTION(scalar, scalar, scalar, pow)
133BINARY_TYPE_FUNCTION(scalar, scalar, scalar, pow)
134
135BINARY_FUNCTION(scalar, scalar, scalar, atan2)
136BINARY_TYPE_FUNCTION(scalar, scalar, scalar, atan2)
137
138BINARY_FUNCTION(scalar, scalar, scalar, hypot)
139BINARY_TYPE_FUNCTION(scalar, scalar, scalar, hypot)
140
141// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142
143UNARY_FUNCTION(scalar, scalar, pow3)
144UNARY_FUNCTION(scalar, scalar, pow4)
145UNARY_FUNCTION(scalar, scalar, pow5)
146UNARY_FUNCTION(scalar, scalar, pow6)
147UNARY_FUNCTION(scalar, scalar, pow025)
148UNARY_FUNCTION(scalar, scalar, sqrt)
149UNARY_FUNCTION(scalar, scalar, cbrt)
150UNARY_FUNCTION(scalar, scalar, sign)
151UNARY_FUNCTION(scalar, scalar, pos)
152UNARY_FUNCTION(scalar, scalar, pos0)
153UNARY_FUNCTION(scalar, scalar, neg)
154UNARY_FUNCTION(scalar, scalar, neg0)
155UNARY_FUNCTION(scalar, scalar, posPart)
156UNARY_FUNCTION(scalar, scalar, negPart)
157UNARY_FUNCTION(scalar, scalar, exp)
158UNARY_FUNCTION(scalar, scalar, log)
159UNARY_FUNCTION(scalar, scalar, log10)
160UNARY_FUNCTION(scalar, scalar, sin)
161UNARY_FUNCTION(scalar, scalar, cos)
162UNARY_FUNCTION(scalar, scalar, tan)
163UNARY_FUNCTION(scalar, scalar, asin)
164UNARY_FUNCTION(scalar, scalar, acos)
165UNARY_FUNCTION(scalar, scalar, atan)
166UNARY_FUNCTION(scalar, scalar, sinh)
167UNARY_FUNCTION(scalar, scalar, cosh)
168UNARY_FUNCTION(scalar, scalar, tanh)
169UNARY_FUNCTION(scalar, scalar, asinh)
170UNARY_FUNCTION(scalar, scalar, acosh)
171UNARY_FUNCTION(scalar, scalar, atanh)
172UNARY_FUNCTION(scalar, scalar, erf)
173UNARY_FUNCTION(scalar, scalar, erfc)
174UNARY_FUNCTION(scalar, scalar, lgamma)
175UNARY_FUNCTION(scalar, scalar, j0)
176UNARY_FUNCTION(scalar, scalar, j1)
177UNARY_FUNCTION(scalar, scalar, y0)
178UNARY_FUNCTION(scalar, scalar, y1)
179
180UNARY_FUNCTION(scalar, scalar, degToRad)
181UNARY_FUNCTION(scalar, scalar, radToDeg)
182UNARY_FUNCTION(scalar, scalar, atmToPa)
183UNARY_FUNCTION(scalar, scalar, barToPa)
184UNARY_FUNCTION(scalar, scalar, paToAtm)
185UNARY_FUNCTION(scalar, scalar, paToBar)
186
187
188#define BesselFunc(func) \
189void 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 \
194tmp<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 \
201tmp<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
211
212#undef BesselFunc
213
214
215// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216
217} // End namespace Foam
218
219// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220
221#include "undefFieldFunctionsM.H"
222
223// ************************************************************************* //
#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 BesselFunc(func)
#define TFOR_ALL_S_OP_F_OP_F(typeS, s, OP1, typeF1, f1, OP2, typeF2, f2)
Definition FieldM.H:540
#define TFOR_ALL_F_OP_FUNC_S_F(typeF1, f1, OP, FUNC, typeS, s, typeF2, f2)
Definition FieldM.H:277
void replace(const direction, const UList< cmptType > &)
tmp< Field< cmptType > > component(const direction) const
Definition scalarField.C:40
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:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
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.
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar asin(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 exp(const dimensionedScalar &ds)
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
constexpr scalar barToPa(const scalar bar) noexcept
Conversion from bar to Pa.
dimensionedScalar lgamma(const dimensionedScalar &ds)
void subtract(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
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)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar log10(const dimensionedScalar &ds)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar negPart(const dimensionedScalar &ds)
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar acosh(const dimensionedScalar &ds)
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar sqrt(const dimensionedScalar &ds)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
dimensionedScalar hypot(const dimensionedScalar &x, const dimensionedScalar &y)
complex sumProd(const UList< complex > &f1, const UList< complex > &f2)
Sum product.
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
uint8_t direction
Definition direction.H:49
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar atan(const dimensionedScalar &ds)
constexpr scalar paToBar(const scalar pa) noexcept
Conversion from Pa to bar.
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar posPart(const dimensionedScalar &ds)
constexpr scalar paToAtm(const scalar pa) noexcept
Conversion from Pa to atm.
dimensionedScalar acos(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
constexpr scalar atmToPa(const scalar atm) noexcept
Conversion from atm to Pa.
dimensionedScalar asinh(const dimensionedScalar &ds)
labelList f(nPoints)
dict add("bounds", meshBb)
Unit conversion functions.