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-2022 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(eqNeg, x = -y)
72 EqOp(eqSqr, x = sqr(y))
73 EqOp(eqMag, x = mag(y))
74 
75 EqOp(plusEq, x += y)
76 EqOp(minusEq, x -= y)
77 EqOp(multiplyEq, x *= y)
78 EqOp(divideEq, x /= y)
79 EqOp(plusEqMagSqr, x += magSqr(y))
80 EqOp(maxEq, x = max(x, y))
81 EqOp(minEq, x = min(x, y))
82 EqOp(minMagSqrEq, x = (magSqr(x) <= magSqr(y) ? x : y))
83 EqOp(maxMagSqrEq, x = (magSqr(x) >= magSqr(y) ? x : y))
84 EqOp(andEq, x = (x && y))
85 EqOp(orEq, x = (x || y))
86 EqOp(notEq, x = (x != y))
87 
88 EqOp(nopEq, (void)x)
89 
90 #undef EqOp
91 
92 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93 
94 #if __GNUC__
95 #define WARNRETURN __attribute__((warn_unused_result))
96 #else
97 #define WARNRETURN
98 #endif
99 
100 #define Op(opName, op) \
101  \
102  template<class T, class T1, class T2> \
103  class opName##Op3 \
104  { \
105  public: \
106  \
107  T operator()(const T1& x, const T2& y) const WARNRETURN \
108  { \
109  return op; \
110  } \
111  }; \
112  \
113  template<class T1, class T2> \
114  class opName##Op2 \
115  { \
116  public: \
117  \
118  T1 operator()(const T1& x, const T2& y) const WARNRETURN \
119  { \
120  return op; \
121  } \
122  }; \
123  \
124  template<class T> \
125  class opName##Op \
126  { \
127  public: \
128  \
129  T operator()(const T& x, const T& y) const WARNRETURN \
130  { \
131  return op; \
132  } \
133  };
134 
135 Op(sum, x + y)
136 Op(plus, x + y)
137 Op(minus, x - y)
144 Op(max, max(x, y))
145 Op(min, min(x, y))
146 Op(minMagSqr, (magSqr(x) <= magSqr(y) ? x : y))
147 Op(maxMagSqr, (magSqr(x) >= magSqr(y) ? x : y))
149 Op(and, x && y)
150 Op(or, x || y)
151 Op(not, x != y)
152 Op(equal, x == y)
153 Op(less, x < y)
154 Op(lessEq, x <= y)
155 Op(greater, x > y)
156 Op(greaterEq, x >= y)
157 
158 #undef Op
159 #undef WARNRETURN
160 
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 
163 } // End namespace Foam
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 #endif
168 
169 // ************************************************************************* //
scalar y
Namespace for OpenFOAM.
bool equal(const T &s1, const T &s2)
Definition: doubleFloat.H:62
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Type maxMagSqr(const UList< Type > &f)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
Scalar minMod(const Scalar s1, const Scalar s2)
Definition: Scalar.H:227
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static bool less(const vector &x, const vector &y)
To compare normals.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Scalar cmptPow(const Scalar s1, const Scalar s2)
Definition: Scalar.H:293
Type minMagSqr(const UList< Type > &f)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
#define Op(opName, op)
Definition: ops.H:100
#define EqOp(opName, op)
Definition: ops.H:46