DimensionedFieldFunctions.C
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-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 \*---------------------------------------------------------------------------*/
28 
30 
31 #define TEMPLATE template<class Type, class GeoMesh>
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
40 
41 template<class Type, class GeoMesh, direction r>
43 pow
44 (
47 )
48 {
49  typedef typename powProduct<Type, r>::type resultType;
50 
51  auto tres =
53  (
54  f1,
55  "pow(" + f1.name() + ',' + Foam::name(r) + ')',
56  pow(f1.dimensions(), r)
57  );
58 
59  pow<Type, r, GeoMesh>(tres.ref().field(), f1.field());
60 
61  return tres;
62 }
63 
64 
65 template<class Type, class GeoMesh, direction r>
67 pow
68 (
69  const tmp<DimensionedField<Type, GeoMesh>>& tf1,
71 )
72 {
73  typedef typename powProduct<Type, r>::type resultType;
74 
75  const auto& f1 = tf1();
76 
77  auto tres =
79  (
80  tf1,
81  "pow(" + f1.name() + ',' + Foam::name(r) + ')',
82  pow(f1.dimensions(), r)
83  );
84 
85  pow<Type, r, GeoMesh>(tres.ref().field(), f1.field());
86 
87  tf1.clear();
88  return tres;
89 }
90 
91 
92 template<class Type, class GeoMesh>
94 sqr(const DimensionedField<Type, GeoMesh>& f1)
95 {
96  typedef typename outerProduct<Type, Type>::type resultType;
97 
98  auto tres =
100  (
101  f1,
102  "sqr(" + f1.name() + ')',
103  sqr(f1.dimensions())
104  );
105 
106  sqr(tres.ref().field(), f1.field());
107 
108  return tres;
109 }
110 
111 template<class Type, class GeoMesh>
113 sqr(const tmp<DimensionedField<Type, GeoMesh>>& tf1)
114 {
115  typedef typename outerProduct<Type, Type>::type resultType;
116 
117  const auto& f1 = tf1();
118 
119  auto tres =
121  (
122  tf1,
123  "sqr(" + f1.name() + ')',
124  sqr(f1.dimensions())
125  );
126 
127  sqr(tres.ref().field(), f1.field());
128 
129  tf1.clear();
130  return tres;
131 }
132 
133 
134 template<class Type, class GeoMesh>
136 magSqr(const DimensionedField<Type, GeoMesh>& f1)
137 {
138  typedef typename typeOfMag<Type>::type resultType;
139 
140  auto tres =
142  (
143  f1,
144  "magSqr(" + f1.name() + ')',
145  sqr(f1.dimensions())
146  );
147 
148  magSqr(tres.ref().field(), f1.field());
149 
150  return tres;
151 }
152 
153 template<class Type, class GeoMesh>
155 magSqr(const tmp<DimensionedField<Type, GeoMesh>>& tf1)
156 {
157  typedef typename typeOfMag<Type>::type resultType;
158 
159  const auto& f1 = tf1();
160 
161  auto tres =
163  (
164  tf1,
165  "magSqr(" + f1.name() + ')',
166  sqr(f1.dimensions())
167  );
168 
169  magSqr(tres.ref().field(), f1.field());
170 
171  tf1.clear();
172  return tres;
173 }
174 
175 
176 template<class Type, class GeoMesh>
178 mag(const DimensionedField<Type, GeoMesh>& f1)
179 {
180  typedef typename typeOfMag<Type>::type resultType;
181 
182  auto tres =
184  (
185  f1,
186  "mag(" + f1.name() + ')',
187  f1.dimensions()
188  );
189 
190  mag(tres.ref().field(), f1.field());
191 
192  return tres;
193 }
194 
195 template<class Type, class GeoMesh>
197 mag(const tmp<DimensionedField<Type, GeoMesh>>& tf1)
198 {
199  typedef typename typeOfMag<Type>::type resultType;
200 
201  const auto& f1 = tf1();
202 
203  auto tres =
205  (
206  tf1,
207  "mag(" + f1.name() + ')',
208  f1.dimensions()
209  );
210 
211  mag(tres.ref().field(), f1.field());
212 
213  tf1.clear();
214  return tres;
215 }
216 
217 
218 template<class Type, class GeoMesh>
220 <
222  <
224  >
225 >
227 {
228  typedef typename DimensionedField<Type, GeoMesh>::cmptType resultType;
229 
230  auto tres =
232  (
233  f1,
234  "cmptAv(" + f1.name() + ')',
235  f1.dimensions()
236  );
237 
238  cmptAv(tres.ref().field(), f1.field());
239 
240  return tres;
241 }
242 
243 
244 template<class Type, class GeoMesh>
246 <
248  <
250  >
251 >
253 {
254  typedef typename DimensionedField<Type, GeoMesh>::cmptType resultType;
255 
256  const auto& f1 = tf1();
257 
258  auto tres =
260  (
261  tf1,
262  "cmptAv(" + f1.name() + ')',
263  f1.dimensions()
264  );
265 
266  cmptAv(tres.ref().field(), f1.field());
267 
268  tf1.clear();
269  return tres;
270 }
271 
272 
273 #define UNARY_REDUCTION_FUNCTION(ReturnType, Func, gFunc) \
274  \
275 template<class Type, class GeoMesh> \
276 dimensioned<ReturnType> Func \
277 ( \
278  const DimensionedField<Type, GeoMesh>& f1 \
279 ) \
280 { \
281  return dimensioned<ReturnType> \
282  ( \
283  #Func "(" + f1.name() + ')', \
284  f1.dimensions(), \
285  gFunc(f1.field()) \
286  ); \
287 } \
288  \
289 template<class Type, class GeoMesh> \
290 dimensioned<ReturnType> Func \
291 ( \
292  const tmp<DimensionedField<Type, GeoMesh>>& tf1 \
293 ) \
294 { \
295  dimensioned<ReturnType> res = Func(tf1()); \
296  tf1.clear(); \
297  return res; \
298 }
299 
306 
308 
309 #undef UNARY_REDUCTION_FUNCTION
310 
311 
312 BINARY_FUNCTION(Type, Type, Type, max)
313 BINARY_FUNCTION(Type, Type, Type, min)
314 BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
315 BINARY_FUNCTION(Type, Type, Type, cmptDivide)
316 
317 BINARY_TYPE_FUNCTION(Type, Type, Type, max)
318 BINARY_TYPE_FUNCTION(Type, Type, Type, min)
319 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
320 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
321 
323 // ------------------------------------------------------------------------- //
324 
325 // Clamp Methods
326 
327 template<class Type, class GeoMesh>
328 void clamp
329 (
330  DimensionedField<Type, GeoMesh>& result,
331  const DimensionedField<Type, GeoMesh>& f1,
332  const Foam::zero_one
333 )
334 {
336 
337  clamp(result.field(), f1.field(), range);
338  result.oriented() = f1.oriented();
339 }
340 
341 template<class Type, class GeoMesh>
342 tmp<DimensionedField<Type, GeoMesh>>
343 clamp
344 (
345  const DimensionedField<Type, GeoMesh>& f1,
346  const Foam::zero_one
347 )
348 {
349  auto tres =
351  (
352  f1,
353  "clamp01(" + f1.name() + ')',
354  f1.dimensions()
355  );
356 
357  clamp(tres.ref(), f1, Foam::zero_one{});
358 
359  return tres;
360 }
361 
362 
363 template<class Type, class GeoMesh>
364 tmp<DimensionedField<Type, GeoMesh>>
365 clamp
366 (
367  const tmp<DimensionedField<Type, GeoMesh>>& tf1,
368  const Foam::zero_one
369 )
370 {
371  const auto& f1 = tf1();
372 
373  auto tres =
375  (
376  tf1,
377  "clamp01(" + f1.name() + ')',
378  f1.dimensions()
379  );
380 
381  clamp(tres.field(), f1, Foam::zero_one{});
382 
383  tf1.clear();
384  return tres;
385 }
386 
387 BINARY_TYPE_FUNCTION_FS(Type, Type, MinMax<Type>, clamp)
388 
389 
390 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
391 
392 TERNARY_FUNCTION(Type, Type, Type, scalar, lerp)
393 TERNARY_TYPE_FUNCTION_FFS(Type, Type, Type, scalar, lerp)
394 
395 
396 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
397 
398 UNARY_OPERATOR(Type, Type, -, negate, transform)
399 
400 BINARY_OPERATOR(Type, Type, scalar, *, '*', multiply)
401 BINARY_OPERATOR(Type, scalar, Type, *, '*', multiply)
402 BINARY_OPERATOR(Type, Type, scalar, /, '|', divide)
403 
404 BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, '*', multiply)
405 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, '*', multiply)
406 
407 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, '|', divide)
408 
409 
410 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
411 
412 #define PRODUCT_OPERATOR(product, Op, OpFunc) \
413  \
414 template<class Type1, class Type2, class GeoMesh> \
415 tmp<DimensionedField<typename product<Type1, Type2>::type, GeoMesh>> \
416 operator Op \
417 ( \
418  const DimensionedField<Type1, GeoMesh>& f1, \
419  const DimensionedField<Type2, GeoMesh>& f2 \
420 ) \
421 { \
422  typedef typename product<Type1, Type2>::type resultType; \
423  \
424  auto tres = \
425  reuseTmpDimensionedField<resultType, Type1, GeoMesh>::New \
426  ( \
427  f1, \
428  '(' + f1.name() + #Op + f2.name() + ')', \
429  (f1.dimensions() Op f2.dimensions()) \
430  ); \
431  \
432  Foam::OpFunc(tres.ref().field(), f1.field(), f2.field()); \
433  \
434  return tres; \
435 } \
436  \
437  \
438 template<class Type1, class Type2, class GeoMesh> \
439 tmp<DimensionedField<typename product<Type1, Type2>::type, GeoMesh>> \
440 operator Op \
441 ( \
442  const DimensionedField<Type1, GeoMesh>& f1, \
443  const tmp<DimensionedField<Type2, GeoMesh>>& tf2 \
444 ) \
445 { \
446  typedef typename product<Type1, Type2>::type resultType; \
447  \
448  const auto& f2 = tf2(); \
449  \
450  auto tres = \
451  reuseTmpDimensionedField<resultType, Type2, GeoMesh>::New \
452  ( \
453  tf2, \
454  '(' + f1.name() + #Op + f2.name() + ')', \
455  (f1.dimensions() Op f2.dimensions()) \
456  ); \
457  \
458  Foam::OpFunc(tres.ref().field(), f1.field(), f2.field()); \
459  \
460  tf2.clear(); \
461  return tres; \
462 } \
463  \
464  \
465 template<class Type1, class Type2, class GeoMesh> \
466 tmp<DimensionedField<typename product<Type1, Type2>::type, GeoMesh>> \
467 operator Op \
468 ( \
469  const tmp<DimensionedField<Type1, GeoMesh>>& tf1, \
470  const DimensionedField<Type2, GeoMesh>& f2 \
471 ) \
472 { \
473  typedef typename product<Type1, Type2>::type resultType; \
474  \
475  const auto& f1 = tf1(); \
476  \
477  auto tres = \
478  reuseTmpDimensionedField<resultType, Type1, GeoMesh>::New \
479  ( \
480  tf1, \
481  '(' + f1.name() + #Op + f2.name() + ')', \
482  (f1.dimensions() Op f2.dimensions()) \
483  ); \
484  \
485  Foam::OpFunc(tres.ref().field(), f1.field(), f2.field()); \
486  \
487  tf1.clear(); \
488  return tres; \
489 } \
490  \
491  \
492 template<class Type1, class Type2, class GeoMesh> \
493 tmp<DimensionedField<typename product<Type1, Type2>::type, GeoMesh>> \
494 operator Op \
495 ( \
496  const tmp<DimensionedField<Type1, GeoMesh>>& tf1, \
497  const tmp<DimensionedField<Type2, GeoMesh>>& tf2 \
498 ) \
499 { \
500  typedef typename product<Type1, Type2>::type resultType; \
501  \
502  const auto& f1 = tf1(); \
503  const auto& f2 = tf2(); \
504  \
505  auto tres = \
506  reuseTmpTmpDimensionedField \
507  <resultType, Type1, Type1, Type2, GeoMesh>::New \
508  ( \
509  tf1, \
510  tf2, \
511  '(' + f1.name() + #Op + f2.name() + ')', \
512  (f1.dimensions() Op f2.dimensions()) \
513  ); \
514  \
515  Foam::OpFunc(tres.ref().field(), f1.field(), f2.field()); \
516  \
517  tf1.clear(); \
518  tf2.clear(); \
519  return tres; \
520 } \
521  \
522  \
523 template<class Form, class Type, class GeoMesh> \
524 tmp<DimensionedField<typename product<Type, Form>::type, GeoMesh>> \
525 operator Op \
526 ( \
527  const DimensionedField<Type, GeoMesh>& f1, \
528  const dimensioned<Form>& dvs \
529 ) \
530 { \
531  typedef typename product<Type, Form>::type resultType; \
532  \
533  auto tres = \
534  reuseTmpDimensionedField<resultType, Type, GeoMesh>::New \
535  ( \
536  f1, \
537  '(' + f1.name() + #Op + dvs.name() + ')', \
538  (f1.dimensions() Op dvs.dimensions()) \
539  ); \
540  \
541  Foam::OpFunc(tres.ref().field(), f1.field(), dvs.value()); \
542  \
543  return tres; \
544 } \
545  \
546  \
547 template<class Form, class Cmpt, direction nCmpt, class Type, class GeoMesh> \
548 tmp<DimensionedField<typename product<Form, Type>::type, GeoMesh>> \
549 operator Op \
550 ( \
551  const DimensionedField<Type, GeoMesh>& f1, \
552  const VectorSpace<Form,Cmpt,nCmpt>& vs \
553 ) \
554 { \
555  return f1 Op dimensioned<Form>(static_cast<const Form&>(vs)); \
556 } \
557  \
558  \
559 template<class Form, class Type, class GeoMesh> \
560 tmp<DimensionedField<typename product<Type, Form>::type, GeoMesh>> \
561 operator Op \
562 ( \
563  const tmp<DimensionedField<Type, GeoMesh>>& tf1, \
564  const dimensioned<Form>& dvs \
565 ) \
566 { \
567  typedef typename product<Type, Form>::type resultType; \
568  \
569  const auto& f1 = tf1(); \
570  \
571  auto tres = \
572  reuseTmpDimensionedField<resultType, Type, GeoMesh>::New \
573  ( \
574  tf1, \
575  '(' + f1.name() + #Op + dvs.name() + ')', \
576  (f1.dimensions() Op dvs.dimensions()) \
577  ); \
578  \
579  Foam::OpFunc(tres.ref().field(), f1.field(), dvs.value()); \
580  \
581  tf1.clear(); \
582  return tres; \
583 } \
584  \
585  \
586 template<class Form, class Cmpt, direction nCmpt, class Type, class GeoMesh> \
587 tmp<DimensionedField<typename product<Form, Type>::type, GeoMesh>> \
588 operator Op \
589 ( \
590  const tmp<DimensionedField<Type, GeoMesh>>& tf1, \
591  const VectorSpace<Form,Cmpt,nCmpt>& vs \
592 ) \
593 { \
594  return tf1 Op dimensioned<Form>(static_cast<const Form&>(vs)); \
595 } \
596  \
597  \
598 template<class Form, class Type, class GeoMesh> \
599 tmp<DimensionedField<typename product<Form, Type>::type, GeoMesh>> \
600 operator Op \
601 ( \
602  const dimensioned<Form>& dvs, \
603  const DimensionedField<Type, GeoMesh>& f2 \
604 ) \
605 { \
606  typedef typename product<Form, Type>::type resultType; \
607  \
608  auto tres = \
609  reuseTmpDimensionedField<resultType, Type, GeoMesh>::New \
610  ( \
611  f2, \
612  '(' + dvs.name() + #Op + f2.name() + ')', \
613  (dvs.dimensions() Op f2.dimensions()) \
614  ); \
615  \
616  Foam::OpFunc(tres.ref().field(), dvs.value(), f2.field()); \
617  \
618  return tres; \
619 } \
620  \
621  \
622 template<class Form, class Cmpt, direction nCmpt, class Type, class GeoMesh> \
623 tmp<DimensionedField<typename product<Form, Type>::type, GeoMesh>> \
624 operator Op \
625 ( \
626  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
627  const DimensionedField<Type, GeoMesh>& f2 \
628 ) \
629 { \
630  return dimensioned<Form>(static_cast<const Form&>(vs)) Op f2; \
631 } \
632  \
633  \
634 template<class Form, class Type, class GeoMesh> \
635 tmp<DimensionedField<typename product<Form, Type>::type, GeoMesh>> \
636 operator Op \
637 ( \
638  const dimensioned<Form>& dvs, \
639  const tmp<DimensionedField<Type, GeoMesh>>& tf2 \
640 ) \
641 { \
642  typedef typename product<Form, Type>::type resultType; \
643  \
644  const auto& f2 = tf2(); \
645  \
646  auto tres = \
647  reuseTmpDimensionedField<resultType, Type, GeoMesh>::New \
648  ( \
649  tf2, \
650  '(' + dvs.name() + #Op + f2.name() + ')', \
651  (dvs.dimensions() Op f2.dimensions()) \
652  ); \
653  \
654  Foam::OpFunc(tres.ref().field(), dvs.value(), f2.field()); \
655  \
656  tf2.clear(); \
657  return tres; \
658 } \
659  \
660  \
661 template<class Form, class Cmpt, direction nCmpt, class Type, class GeoMesh> \
662 tmp<DimensionedField<typename product<Form, Type>::type, GeoMesh>> \
663 operator Op \
664 ( \
665  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
666  const tmp<DimensionedField<Type, GeoMesh>>& tf2 \
667 ) \
668 { \
669  return dimensioned<Form>(static_cast<const Form&>(vs)) Op tf2; \
670 }
675 
680 
681 #undef PRODUCT_OPERATOR
682 
683 
684 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
685 
686 } // End namespace Foam
687 
688 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
689 
690 #include "undefFieldFunctionsM.H"
691 
692 // ************************************************************************* //
dimensioned< typename typeOfMag< Type >::type > sumMag(const DimensionedField< Type, GeoMesh > &f1)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
type
Types of root.
Definition: Roots.H:52
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &f1)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type gMin(const FieldField< Field, Type > &f)
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
A min/max value pair with additional methods. In addition to conveniently storing values...
Definition: HashSet.H:72
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define UNARY_REDUCTION_FUNCTION(ReturnType, Func, gFunc)
#define BINARY_TYPE_OPERATOR_FS(ReturnType, Type1, Type2, Op, OpName, OpFunc)
void cross(FieldField< Field1, typename crossProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:118
#define BINARY_TYPE_FUNCTION(ReturnType, Type1, Type2, Func)
Field< Type >::cmptType cmptType
Component type of the field elements.
void dotdot(FieldField< Field1, typename scalarProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
void outer(FieldField< Field1, typename outerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
scalar range
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
#define PRODUCT_OPERATOR(product, Op, OpFunc)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:54
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
void negate(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
static tmp< DimensionedField< TypeR, GeoMesh > > New(const DimensionedField< Type1, GeoMesh > &f1, const word &name, const dimensionSet &dimensions)
Pass-through to New DimensionedField.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Type gMax(const FieldField< Field, Type > &f)
scalarMinMax gMinMaxMag(const FieldField< Field, Type > &f)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
const Field< Type > & field() const noexcept
Return const-reference to the field values.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition: pTraits.H:65
symmTypeOfRank< typename pTraits< arg1 >::cmptType, arg2 *direction(pTraits< arg1 >::rank) >::type type
Definition: products.H:176
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void clear()
Clear exponents - resets to be dimensionless.
Definition: dimensionSet.C:136
dimensioned< scalarMinMax > minMaxMag(const DimensionedField< Type, GeoMesh > &f1)
Type gAverage(const FieldField< Field, Type > &f)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
BINARY_TYPE_FUNCTION_FS(Type, Type, MinMax< Type >, clip)
#define UNARY_OPERATOR(ReturnType, Type1, Op, OpFunc, Dfunc)
#define TERNARY_TYPE_FUNCTION_FFS(ReturnType, Type1, Type2, Type3, Func)
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
#define TERNARY_FUNCTION(ReturnType, Type1, Type2, Type3, Func)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
const dimensionSet & dimensions() const noexcept
Return dimensions.
pTraits< typename pTraits< arg1 >::cmptType >::magType type
Definition: products.H:96