ops.H
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-2020 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 InClass
25  Foam::Pstream
26 
27 Description
28  Combination-Reduction operation for a parallel run.
29 
30  The information from all nodes is collected on the master node,
31  combined using the given combination function and the result is
32  broadcast to all nodes
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef ops_H
37 #define ops_H
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 #define EqOp(opName, op) \
47  \
48 template<class T1, class T2> \
49 class opName##Op2 \
50 { \
51 public: \
52  \
53  void operator()(T1& x, const T2& y) const \
54  { \
55  op; \
56  } \
57 }; \
58  \
59 template<class T> \
60 class opName##Op \
61 { \
62 public: \
63  \
64  void operator()(T& x, const T& y) const \
65  { \
66  op; \
67  } \
68 };
69 
70 EqOp(eq, x = y)
71 EqOp(plusEq, x += y)
72 EqOp(minusEq, x -= y)
73 EqOp(multiplyEq, x *= y)
74 EqOp(divideEq, x /= y)
75 EqOp(eqSqr, x = sqr(y))
76 EqOp(eqMag, x = mag(y))
77 EqOp(plusEqMagSqr, x += magSqr(y))
78 EqOp(maxEq, x = max(x, y))
79 EqOp(minEq, x = min(x, y))
80 EqOp(minMagSqrEq, x = (magSqr(x)<=magSqr(y) ? x : y))
81 EqOp(maxMagSqrEq, x = (magSqr(x)>=magSqr(y) ? x : y))
82 EqOp(andEq, x = (x && y))
83 EqOp(orEq, x = (x || y))
84 EqOp(notEq, x = (x != y))
85 
86 EqOp(eqMinus, x = -y)
87 
88 EqOp(nopEq, (void)x)
89 
90 #undef EqOp
91 
92 
93 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
94 
95 #if __GNUC__
96 #define WARNRETURN __attribute__((warn_unused_result))
97 #else
98 #define WARNRETURN
99 #endif
100 
101 #define Op(opName, op) \
102  \
103  template<class T, class T1, class T2> \
104  class opName##Op3 \
105  { \
106  public: \
107  \
108  T operator()(const T1& x, const T2& y) const WARNRETURN \
109  { \
110  return op; \
111  } \
112  }; \
113  \
114  template<class T1, class T2> \
115  class opName##Op2 \
116  { \
117  public: \
118  \
119  T1 operator()(const T1& x, const T2& y) const WARNRETURN \
120  { \
121  return op; \
122  } \
123  }; \
124  \
125  template<class T> \
126  class opName##Op \
127  { \
128  public: \
129  \
130  T operator()(const T& x, const T& y) const WARNRETURN \
131  { \
132  return op; \
133  } \
134  };
135 
136 
137 #define weightedOp(opName, op) \
138  \
139  template<class Type, class CombineOp> \
140  class opName##WeightedOp \
141  { \
142  const CombineOp& cop_; \
143  \
144  public: \
145  \
146  opName##WeightedOp(const CombineOp& cop) \
147  : \
148  cop_(cop) \
149  {} \
150  \
151  void operator() \
152  ( \
153  Type& x, \
154  const label index, \
155  const Type& y, \
156  const scalar weight \
157  ) const \
158  { \
159  cop_(x, op); \
160  } \
161  }; \
162 
163 
164 Op(sum, x + y)
165 
166 Op(plus, x + y)
167 Op(minus, x - y)
168 Op(multiply, x * y)
169 Op(divide, x / y)
170 Op(cmptMultiply, cmptMultiply(x, y))
171 Op(cmptPow, cmptPow(x, y))
172 Op(cmptDivide, cmptDivide(x, y))
173 Op(stabilise, stabilise(x, y))
174 Op(max, max(x, y))
175 Op(min, min(x, y))
176 Op(minMagSqr, (magSqr(x)<=magSqr(y) ? x : y))
177 Op(maxMagSqr, (magSqr(x)>=magSqr(y) ? x : y))
178 Op(minMod, minMod(x, y))
179 Op(and, x && y)
180 Op(or, x || y)
181 Op(not, x != y)
182 Op(eqEq, x == y)
183 Op(less, x < y)
184 Op(lessEq, x <= y)
185 Op(greater, x > y)
186 Op(greaterEq, x >= y)
187 
188 weightedOp(multiply, (weight*y))
189 
190 #undef Op
191 #undef weightedOp
192 #undef WARNRETURN
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 } // End namespace Foam
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 #endif
201 
202 // ************************************************************************* //
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Type minMagSqr(const UList< Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
#define Op(opName, op)
Definition: ops.H:101
static bool less(const vector &x, const vector &y)
To compare normals.
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
scalar y
#define EqOp(opName, op)
Definition: ops.H:46
Scalar cmptPow(const Scalar s1, const Scalar s2)
Definition: Scalar.H:293
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Scalar minMod(const Scalar s1, const Scalar s2)
Definition: Scalar.H:227
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Type maxMagSqr(const UList< Type > &f)
#define weightedOp(opName, op)
Definition: ops.H:137
dimensioned< scalar > mag(const dimensioned< Type > &)
Namespace for OpenFOAM.