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 | Copyright (C) 2011-2015 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 Description
25  Specialisation of Field<T> for scalar.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "scalarField.H"
30 #include "unitConversion.H"
31 
32 #define TEMPLATE
33 #include "FieldFunctionsM.C"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<>
43 tmp<scalarField> scalarField::component(const direction) const
44 {
45  return *this;
46 }
47 
49 {
50  sf = f;
51 }
52 
53 template<>
54 void scalarField::replace(const direction, const UList<scalar>& sf)
55 {
56  *this = sf;
57 }
58 
59 template<>
60 void scalarField::replace(const direction, const scalar& s)
61 {
62  *this = s;
63 }
64 
65 
66 void stabilise(scalarField& res, const UList<scalar>& sf, const scalar s)
67 {
69  (
70  scalar, res, =, ::Foam::stabilise, scalar, s, scalar, sf
71  )
72 }
73 
75 {
76  tmp<scalarField> tRes(new scalarField(sf.size()));
77  stabilise(tRes(), sf, s);
78  return tRes;
79 }
80 
81 tmp<scalarField> stabilise(const tmp<scalarField>& tsf, const scalar s)
82 {
84  stabilise(tRes(), tsf(), s);
86  return tRes;
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 
92 template<>
93 scalar sumProd(const UList<scalar>& f1, const UList<scalar>& f2)
94 {
95  if (f1.size() && (f1.size() == f2.size()))
96  {
97  scalar SumProd = 0.0;
98  TFOR_ALL_S_OP_F_OP_F(scalar, SumProd, +=, scalar, f1, *, scalar, f2)
99  return SumProd;
100  }
101  else
102  {
103  return 0.0;
104  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, add)
111 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, subtract)
112 
113 BINARY_OPERATOR(scalar, scalar, scalar, *, multiply)
114 BINARY_OPERATOR(scalar, scalar, scalar, /, divide)
115 
116 BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, divide)
117 
118 BINARY_FUNCTION(scalar, scalar, scalar, pow)
119 BINARY_TYPE_FUNCTION(scalar, scalar, scalar, pow)
120 
121 BINARY_FUNCTION(scalar, scalar, scalar, atan2)
122 BINARY_TYPE_FUNCTION(scalar, scalar, scalar, atan2)
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 UNARY_FUNCTION(scalar, scalar, pow3)
127 UNARY_FUNCTION(scalar, scalar, pow4)
128 UNARY_FUNCTION(scalar, scalar, pow5)
129 UNARY_FUNCTION(scalar, scalar, pow6)
130 UNARY_FUNCTION(scalar, scalar, pow025)
131 UNARY_FUNCTION(scalar, scalar, sqrt)
132 UNARY_FUNCTION(scalar, scalar, cbrt)
133 UNARY_FUNCTION(scalar, scalar, sign)
134 UNARY_FUNCTION(scalar, scalar, pos)
135 UNARY_FUNCTION(scalar, scalar, neg)
136 UNARY_FUNCTION(scalar, scalar, posPart)
137 UNARY_FUNCTION(scalar, scalar, negPart)
138 UNARY_FUNCTION(scalar, scalar, exp)
139 UNARY_FUNCTION(scalar, scalar, log)
140 UNARY_FUNCTION(scalar, scalar, log10)
141 UNARY_FUNCTION(scalar, scalar, sin)
142 UNARY_FUNCTION(scalar, scalar, cos)
143 UNARY_FUNCTION(scalar, scalar, tan)
144 UNARY_FUNCTION(scalar, scalar, asin)
145 UNARY_FUNCTION(scalar, scalar, acos)
146 UNARY_FUNCTION(scalar, scalar, atan)
147 UNARY_FUNCTION(scalar, scalar, sinh)
148 UNARY_FUNCTION(scalar, scalar, cosh)
149 UNARY_FUNCTION(scalar, scalar, tanh)
150 UNARY_FUNCTION(scalar, scalar, asinh)
151 UNARY_FUNCTION(scalar, scalar, acosh)
152 UNARY_FUNCTION(scalar, scalar, atanh)
153 UNARY_FUNCTION(scalar, scalar, erf)
154 UNARY_FUNCTION(scalar, scalar, erfc)
155 UNARY_FUNCTION(scalar, scalar, lgamma)
156 UNARY_FUNCTION(scalar, scalar, j0)
157 UNARY_FUNCTION(scalar, scalar, j1)
158 UNARY_FUNCTION(scalar, scalar, y0)
159 UNARY_FUNCTION(scalar, scalar, y1)
160 
161 UNARY_FUNCTION(scalar, scalar, degToRad)
162 UNARY_FUNCTION(scalar, scalar, radToDeg)
163 UNARY_FUNCTION(scalar, scalar, atmToPa)
164 UNARY_FUNCTION(scalar, scalar, paToAtm)
165 
166 
167 #define BesselFunc(func) \
168 void func(scalarField& res, const int n, const UList<scalar>& sf) \
169 { \
170  TFOR_ALL_F_OP_FUNC_S_F(scalar, res, =, ::Foam::func, int, n, scalar, sf) \
171 } \
172  \
173 tmp<scalarField> func(const int n, const UList<scalar>& sf) \
174 { \
175  tmp<scalarField> tRes(new scalarField(sf.size())); \
176  func(tRes(), n, sf); \
177  return tRes; \
178 } \
179  \
180 tmp<scalarField> func(const int n, const tmp<scalarField>& tsf) \
181 { \
182  tmp<scalarField> tRes = reuseTmp<scalar, scalar>::New(tsf); \
183  func(tRes(), n, tsf()); \
184  reuseTmp<scalar, scalar>::clear(tsf); \
185  return tRes; \
186 }
187 
190 
191 #undef BesselFunc
192 
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 } // End namespace Foam
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 #include "undefFieldFunctionsM.H"
201 
202 // ************************************************************************* //
BINARY_TYPE_OPERATOR_SF
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
Definition: DimensionedFieldFunctionsM.C:531
BINARY_OPERATOR
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
Definition: DimensionedFieldFunctionsM.C:417
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::subtract
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Definition: FieldFieldFunctions.C:871
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:41
Foam::tan
dimensionedScalar tan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:257
UNARY_FUNCTION
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
Definition: DimensionedFieldFunctionsM.C:30
Foam::reuseTmp::clear
static void clear(const tmp< Field< Type1 > > &tf1)
Definition: FieldReuseFunctions.H:46
Foam::cosh
dimensionedScalar cosh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
Foam::y1
dimensionedScalar y1(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:273
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
scalarField.H
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
Foam::jn
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
Definition: dimensionedScalar.C:296
Foam::posPart
dimensionedScalar posPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:212
Foam::atan2
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
Definition: dimensionedScalar.C:303
BINARY_TYPE_OPERATOR
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
Definition: DimensionedFieldFunctionsM.C:684
undefFieldFunctionsM.H
unitConversion.H
Unit conversion functions.
BesselFunc
#define BesselFunc(func)
Definition: scalarField.C:167
Foam::paToAtm
scalar paToAtm(const scalar pa)
Conversion from atm to Pa.
Definition: unitConversion.H:63
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
BINARY_FUNCTION
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
Definition: DimensionedFieldFunctionsM.C:142
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:179
Foam::erf
dimensionedScalar erf(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:267
Foam::divide
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:131
Foam::atanh
dimensionedScalar atanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:266
Foam::lgamma
dimensionedScalar lgamma(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:269
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:98
Foam::pow6
dimensionedScalar pow6(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:120
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::stabilise
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Definition: DimensionedScalarField.C:40
Foam::reuseTmp::New
static tmp< Field< TypeR > > New(const tmp< Field< Type1 > > &tf1)
Definition: FieldReuseFunctions.H:41
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:254
Foam::y0
dimensionedScalar y0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:272
Foam::erfc
dimensionedScalar erfc(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
f1
scalar f1
Definition: createFields.H:28
Foam::Field::replace
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Foam::asinh
dimensionedScalar asinh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::add
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Definition: FieldFieldFunctions.C:870
Foam::pow5
dimensionedScalar pow5(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:109
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
TFOR_ALL_S_OP_F_OP_F
#define TFOR_ALL_S_OP_F_OP_F(typeS, s, OP1, typeF1, f1, OP2, typeF2, f2)
Definition: FieldM.H:376
sf
volScalarField sf(fieldObject, mesh)
Foam::Field::component
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:644
s
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
BINARY_TYPE_FUNCTION
#define BINARY_TYPE_FUNCTION(ReturnType, Type1, Type2, Func)
Definition: DimensionedFieldFunctionsM.C:410
Foam::yn
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
Definition: dimensionedScalar.C:297
Foam::acosh
dimensionedScalar acosh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::negPart
dimensionedScalar negPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:223
TFOR_ALL_F_OP_FUNC_S_F
#define TFOR_ALL_F_OP_FUNC_S_F(typeF1, f1, OP, FUNC, typeS, s, typeF2, f2)
Definition: FieldM.H:210
f
labelList f(nPoints)
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::j0
dimensionedScalar j0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:270
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:259
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
FieldFunctionsM.C
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::atan
dimensionedScalar atan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:260
Foam::radToDeg
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
Definition: unitConversion.H:51
Foam::multiply
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
Foam::UList::size
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:153
Foam::j1
dimensionedScalar j1(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:271
Foam::sumProd
scalar sumProd(const UList< Type > &f1, const UList< Type > &f2)
Definition: FieldFunctions.C:416
Foam::atmToPa
scalar atmToPa(const scalar atm)
Conversion from atm to Pa.
Definition: unitConversion.H:57
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:201
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:258
Foam::degToRad
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Definition: unitConversion.H:45
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190
Foam::sinh
dimensionedScalar sinh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261