scalarField.H
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
27Typedef
28 Foam::scalarField
29
30Description
31 Specialisation of Field<T> for scalar.
32
33SourceFiles
34 scalarField.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef scalarField_H
39#define scalarField_H
40
41#include "Field.H"
42#include "scalar.H"
43
44#define TEMPLATE
45#include "FieldFunctionsM.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52typedef Field<scalar> scalarField;
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56template<>
57tmp<scalarField> scalarField::component(const direction) const;
58
59void component
60(
61 scalarField& sf,
62 const UList<scalar>& f,
63 const direction
64);
65
66template<>
67void scalarField::replace(const direction, const UList<scalar>& sf);
68
69template<>
70void scalarField::replace(const direction, const scalar& s);
71
72
73void stabilise(scalarField& Res, const UList<scalar>& sf, const scalar s);
74tmp<scalarField> stabilise(const UList<scalar>&, const scalar s);
75tmp<scalarField> stabilise(const tmp<scalarField>&, const scalar s);
76
77
78// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79
80//- Sum product for float
81template<>
82float sumProd(const UList<float>& f1, const UList<float>& f2);
83
84//- Sum product for double
85template<>
86double sumProd(const UList<double>& f1, const UList<double>& f2);
87
88
89// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
90
91BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, add)
92BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, subtract)
93
94BINARY_OPERATOR(scalar, scalar, scalar, *, multiply)
95BINARY_OPERATOR(scalar, scalar, scalar, /, divide)
96
97BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, divide)
98
99BINARY_FUNCTION(scalar, scalar, scalar, pow)
100BINARY_TYPE_FUNCTION(scalar, scalar, scalar, pow)
101
102BINARY_FUNCTION(scalar, scalar, scalar, atan2)
103BINARY_TYPE_FUNCTION(scalar, scalar, scalar, atan2)
104
105BINARY_FUNCTION(scalar, scalar, scalar, hypot)
106BINARY_TYPE_FUNCTION(scalar, scalar, scalar, hypot)
107
108
109// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110
111UNARY_FUNCTION(scalar, scalar, pow3)
112UNARY_FUNCTION(scalar, scalar, pow4)
113UNARY_FUNCTION(scalar, scalar, pow5)
114UNARY_FUNCTION(scalar, scalar, pow6)
115UNARY_FUNCTION(scalar, scalar, pow025)
116UNARY_FUNCTION(scalar, scalar, sqrt)
117UNARY_FUNCTION(scalar, scalar, cbrt)
118UNARY_FUNCTION(scalar, scalar, sign)
119UNARY_FUNCTION(scalar, scalar, pos)
120UNARY_FUNCTION(scalar, scalar, pos0)
121UNARY_FUNCTION(scalar, scalar, neg)
122UNARY_FUNCTION(scalar, scalar, neg0)
123UNARY_FUNCTION(scalar, scalar, posPart)
124UNARY_FUNCTION(scalar, scalar, negPart)
125UNARY_FUNCTION(scalar, scalar, exp)
126UNARY_FUNCTION(scalar, scalar, log)
127UNARY_FUNCTION(scalar, scalar, log10)
128UNARY_FUNCTION(scalar, scalar, sin)
129UNARY_FUNCTION(scalar, scalar, cos)
130UNARY_FUNCTION(scalar, scalar, tan)
131UNARY_FUNCTION(scalar, scalar, asin)
132UNARY_FUNCTION(scalar, scalar, acos)
133UNARY_FUNCTION(scalar, scalar, atan)
134UNARY_FUNCTION(scalar, scalar, sinh)
135UNARY_FUNCTION(scalar, scalar, cosh)
136UNARY_FUNCTION(scalar, scalar, tanh)
137UNARY_FUNCTION(scalar, scalar, asinh)
138UNARY_FUNCTION(scalar, scalar, acosh)
139UNARY_FUNCTION(scalar, scalar, atanh)
140UNARY_FUNCTION(scalar, scalar, erf)
141UNARY_FUNCTION(scalar, scalar, erfc)
142UNARY_FUNCTION(scalar, scalar, lgamma)
143UNARY_FUNCTION(scalar, scalar, j0)
144UNARY_FUNCTION(scalar, scalar, j1)
145UNARY_FUNCTION(scalar, scalar, y0)
146UNARY_FUNCTION(scalar, scalar, y1)
147
148UNARY_FUNCTION(scalar, scalar, degToRad)
149UNARY_FUNCTION(scalar, scalar, radToDeg)
150UNARY_FUNCTION(scalar, scalar, atmToPa)
151UNARY_FUNCTION(scalar, scalar, barToPa)
152UNARY_FUNCTION(scalar, scalar, paToAtm)
153UNARY_FUNCTION(scalar, scalar, paToBar)
154
155#define BesselFunc(func) \
156void func(scalarField& Res, const int n, const UList<scalar>& sf); \
157tmp<scalarField> func(const int n, const UList<scalar>&); \
158tmp<scalarField> func(const int n, const tmp<scalarField>&);
159
162
163#undef BesselFunc
164
165
166// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167
168} // End namespace Foam
169
170// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171
172#include "undefFieldFunctionsM.H"
173
174// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175
176#endif
177
178// ************************************************************************* //
#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)
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:557
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)
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)
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)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
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)
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)
dimensionedScalar sqrt(const dimensionedScalar &ds)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar hypot(const dimensionedScalar &x, const dimensionedScalar &y)
complex sumProd(const UList< complex > &f1, const UList< complex > &f2)
Sum product.
Definition: complexField.C:191
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
uint8_t direction
Definition: direction.H:56
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 divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &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 atmToPa(const scalar atm) noexcept
Conversion from atm to Pa.
dimensionedScalar asinh(const dimensionedScalar &ds)
labelList f(nPoints)