ops.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 | 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 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(eqMag, x = mag(y))
76 EqOp(plusEqMagSqr, x += magSqr(y))
77 EqOp(maxEq, x = max(x, y))
78 EqOp(minEq, x = min(x, y))
79 EqOp(minMagSqrEq, x = (magSqr(x)<=magSqr(y) ? x : y))
80 EqOp(maxMagSqrEq, x = (magSqr(x)>=magSqr(y) ? x : y))
81 EqOp(andEq, x = (x && y))
82 EqOp(orEq, x = (x || y))
83 
84 EqOp(eqMinus, x = -y)
85 
86 EqOp(nopEq, (void)x)
87 
88 #undef EqOp
89 
90 
91 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
92 
93 #if __GNUC__
94 #define WARNRETURN __attribute__((warn_unused_result))
95 #else
96 #define WARNRETURN
97 #endif
98 
99 #define Op(opName, op) \
100  \
101  template<class T, class T1, class T2> \
102  class opName##Op3 \
103  { \
104  public: \
105  \
106  T operator()(const T1& x, const T2& y) const WARNRETURN \
107  { \
108  return op; \
109  } \
110  }; \
111  \
112  template<class T1, class T2> \
113  class opName##Op2 \
114  { \
115  public: \
116  \
117  T1 operator()(const T1& x, const T2& y) const WARNRETURN \
118  { \
119  return op; \
120  } \
121  }; \
122  \
123  template<class T> \
124  class opName##Op \
125  { \
126  public: \
127  \
128  T operator()(const T& x, const T& y) const WARNRETURN \
129  { \
130  return op; \
131  } \
132  };
133 
134 
135 #define weightedOp(opName, op) \
136  \
137  template<class Type, class CombineOp> \
138  class opName##WeightedOp \
139  { \
140  const CombineOp& cop_; \
141  \
142  public: \
143  \
144  opName##WeightedOp(const CombineOp& cop) \
145  : \
146  cop_(cop) \
147  {} \
148  \
149  void operator() \
150  ( \
151  Type& x, \
152  const label index, \
153  const Type& y, \
154  const scalar weight \
155  ) const \
156  { \
157  cop_(x, op); \
158  } \
159  }; \
160 
161 
162 Op(sum, x + y)
163 
164 Op(plus, x + y)
165 Op(minus, x - y)
172 Op(max, max(x, y))
173 Op(min, min(x, y))
177 Op(and, x && y)
178 Op(or, x || y)
179 Op(eqEq, x == y)
180 Op(less, x < y)
181 Op(lessEq, x <= y)
182 Op(greater, x > y)
183 Op(greaterEq, x >= y)
184 
186 
187 #undef Op
188 #undef weightedOp
189 #undef WARNRETURN
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 } // End namespace Foam
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 #endif
198 
199 // ************************************************************************* //
EqOp
#define EqOp(opName, op)
Definition: ops.H:46
Foam::cmptPow
Scalar cmptPow(const Scalar s1, const Scalar s2)
Definition: Scalar.H:249
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::less
static bool less(const vector &x, const vector &y)
To compare normals.
Definition: meshRefinementRefine.C:49
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::divide
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
Foam::maxMagSqr
Type maxMagSqr(const UList< Type > &f)
Definition: FieldFunctions.C:362
weightedOp
#define weightedOp(opName, op)
Definition: ops.H:135
Foam::stabilise
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Definition: DimensionedScalarField.C:40
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::minMagSqr
Type minMagSqr(const UList< Type > &f)
Definition: FieldFunctions.C:389
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::cmptDivide
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::minMod
Scalar minMod(const Scalar s1, const Scalar s2)
Definition: Scalar.H:167
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::multiply
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
Op
#define Op(opName, op)
Definition: ops.H:99
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
y
scalar y
Definition: LISASMDCalcMethod1.H:14