scalarFieldField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 FieldField<T> for scalar.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "scalarFieldField.H"
30 
31 #define TEMPLATE template<template<class> class Field>
32 #include "FieldFieldFunctionsM.C"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 template<template<class> class Field>
42 void stabilise
43 (
45  const FieldField<Field, scalar>& f1,
46  const scalar s
47 )
48 {
49  forAll(f, i)
50  {
51  stabilise(f[i], f1[i], s);
52  }
53 }
54 
55 template<template<class> class Field>
57 (
58  const FieldField<Field, scalar>& f1,
59  const scalar s
60 )
61 {
63  (
65  );
66  stabilise(tf(), f1, s);
67  return tf;
68 }
69 
70 template<template<class> class Field>
72 (
73  const tmp<FieldField<Field, scalar>>& tf1,
74  const scalar s
75 )
76 {
78  stabilise(tf(), tf(), s);
79  return tf;
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, add)
86 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, subtract)
87 
88 BINARY_OPERATOR(scalar, scalar, scalar, *, multiply)
89 BINARY_OPERATOR(scalar, scalar, scalar, /, divide)
90 
91 BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, divide)
92 
93 BINARY_FUNCTION(scalar, scalar, scalar, pow)
94 BINARY_TYPE_FUNCTION(scalar, scalar, scalar, pow)
95 
96 BINARY_FUNCTION(scalar, scalar, scalar, atan2)
97 BINARY_TYPE_FUNCTION(scalar, scalar, scalar, atan2)
98 
99 
100 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
101 
102 UNARY_FUNCTION(scalar, scalar, pow3)
103 UNARY_FUNCTION(scalar, scalar, pow4)
104 UNARY_FUNCTION(scalar, scalar, pow5)
105 UNARY_FUNCTION(scalar, scalar, pow6)
106 UNARY_FUNCTION(scalar, scalar, pow025)
107 UNARY_FUNCTION(scalar, scalar, sqrt)
108 UNARY_FUNCTION(scalar, scalar, cbrt)
109 UNARY_FUNCTION(scalar, scalar, sign)
110 UNARY_FUNCTION(scalar, scalar, pos)
111 UNARY_FUNCTION(scalar, scalar, pos0)
112 UNARY_FUNCTION(scalar, scalar, neg)
113 UNARY_FUNCTION(scalar, scalar, neg0)
114 UNARY_FUNCTION(scalar, scalar, posPart)
115 UNARY_FUNCTION(scalar, scalar, negPart)
116 UNARY_FUNCTION(scalar, scalar, exp)
117 UNARY_FUNCTION(scalar, scalar, log)
118 UNARY_FUNCTION(scalar, scalar, log10)
119 UNARY_FUNCTION(scalar, scalar, sin)
120 UNARY_FUNCTION(scalar, scalar, cos)
121 UNARY_FUNCTION(scalar, scalar, tan)
122 UNARY_FUNCTION(scalar, scalar, asin)
123 UNARY_FUNCTION(scalar, scalar, acos)
124 UNARY_FUNCTION(scalar, scalar, atan)
125 UNARY_FUNCTION(scalar, scalar, sinh)
126 UNARY_FUNCTION(scalar, scalar, cosh)
127 UNARY_FUNCTION(scalar, scalar, tanh)
128 UNARY_FUNCTION(scalar, scalar, asinh)
129 UNARY_FUNCTION(scalar, scalar, acosh)
130 UNARY_FUNCTION(scalar, scalar, atanh)
131 UNARY_FUNCTION(scalar, scalar, erf)
132 UNARY_FUNCTION(scalar, scalar, erfc)
133 UNARY_FUNCTION(scalar, scalar, lgamma)
134 UNARY_FUNCTION(scalar, scalar, j0)
135 UNARY_FUNCTION(scalar, scalar, j1)
136 UNARY_FUNCTION(scalar, scalar, y0)
137 UNARY_FUNCTION(scalar, scalar, y1)
138 
139 
140 #define BesselFunc(func) \
141  \
142 template<template<class> class Field> \
143 void func \
144 ( \
145  FieldField<Field, scalar>& res, \
146  const int n, \
147  const FieldField<Field, scalar>& sf \
148 ) \
149 { \
150  forAll(res, i) \
151  { \
152  func(res[i], n, sf[i]); \
153  } \
154 } \
155  \
156 template<template<class> class Field> \
157 tmp<FieldField<Field, scalar>> func \
158 ( \
159  const int n, \
160  const FieldField<Field, scalar>& sf \
161 ) \
162 { \
163  tmp<FieldField<Field, scalar>> tRes \
164  ( \
165  FieldField<Field, scalar>::NewCalculatedType(sf) \
166  ); \
167  func(tRes(), n, sf); \
168  return tRes; \
169 } \
170  \
171 template<template<class> class Field> \
172 tmp<FieldField<Field, scalar>> func \
173 ( \
174  const int n, \
175  const tmp<FieldField<Field, scalar>>& tsf \
176 ) \
177 { \
178  tmp<FieldField<Field, scalar>> tRes(New(tsf)); \
179  func(tRes(), n, tsf()); \
180  tsf.clear(); \
181  return tRes; \
182 }
183 
186 
187 #undef BesselFunc
188 
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 } // End namespace Foam
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 #include "undefFieldFunctionsM.H"
197 
198 // ************************************************************************* //
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionedScalar log(const dimensionedScalar &ds)
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
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)
const tensorField & tf
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
Generic field type.
Definition: FieldField.H:51
dimensionedScalar j1(const dimensionedScalar &ds)
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
Pre-declare SubField and related Field type.
Definition: Field.H:56
dimensionedScalar asinh(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
Specialisation of FieldField<T> for scalar.
dimensionedScalar y1(const dimensionedScalar &ds)
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 BesselFunc(func)
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensionedScalar atan(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
A class for managing temporary objects.
Definition: PtrList.H:53
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
dimensionedScalar log10(const dimensionedScalar &ds)
Namespace for OpenFOAM.
dimensionedScalar negPart(const dimensionedScalar &ds)