FieldFieldFunctions.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 | Copyright (C) 2011-2016 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 \*---------------------------------------------------------------------------*/
25 
26 #include "PstreamReduceOps.H"
28 
29 #define TEMPLATE template<template<class> class Field, class Type>
30 #include "FieldFieldFunctionsM.C"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 /* * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * */
38 
39 template<template<class> class Field, class Type>
40 void component
41 (
43  const FieldField<Field, Type>& f,
44  const direction d
45 )
46 {
47  forAll(sf, i)
48  {
49  component(sf[i], f[i], d);
50  }
51 }
52 
53 
54 template<template<class> class Field, class Type>
56 {
57  forAll(f1, i)
58  {
59  T(f1[i], f2[i]);
60  }
61 }
62 
63 
64 template<template<class> class Field, class Type, direction r>
65 void pow
66 (
68  const FieldField<Field, Type>& vf
69 )
70 {
71  forAll(f, i)
72  {
73  pow(f[i], vf[i]);
74  }
75 }
76 
77 template<template<class> class Field, class Type, direction r>
79 pow
80 (
82 )
83 {
84  typedef typename powProduct<Type, r>::type powProductType;
86  (
88  );
89  pow<Type, r>(tRes.ref(), f);
90  return tRes;
91 }
92 
93 template<template<class> class Field, class Type, direction r>
95 pow
96 (
98 )
99 {
100  typedef typename powProduct<Type, r>::type powProductType;
102  (
104  );
105  pow<Type, r>(tRes.ref(), tf());
106  tf.clear();
107  return tRes;
108 }
109 
110 
111 template<template<class> class Field, class Type>
112 void sqr
113 (
115  const FieldField<Field, Type>& vf
116 )
117 {
118  forAll(f, i)
119  {
120  sqr(f[i], vf[i]);
121  }
122 }
123 
124 template<template<class> class Field, class Type>
127 {
128  typedef typename outerProduct<Type, Type>::type outerProductType;
130  (
132  );
133  sqr(tRes.ref(), f);
134  return tRes;
135 }
136 
137 template<template<class> class Field, class Type>
140 {
141  typedef typename outerProduct<Type, Type>::type outerProductType;
143  (
145  );
146  sqr(tRes.ref(), tf());
147  tf.clear();
148  return tRes;
149 }
150 
151 
152 template<template<class> class Field, class Type>
154 {
155  forAll(sf, i)
156  {
157  magSqr(sf[i], f[i]);
158  }
159 }
160 
161 template<template<class> class Field, class Type>
163 {
165  (
167  );
168 
169  magSqr(tRes.ref(), f);
170  return tRes;
171 }
172 
173 template<template<class> class Field, class Type>
175 {
177  (
179  );
180 
181  magSqr(tRes.ref(), tf());
182  tf.clear();
183  return tRes;
184 }
185 
186 
187 template<template<class> class Field, class Type>
189 {
190  forAll(sf, i)
191  {
192  mag(sf[i], f[i]);
193  }
194 }
195 
196 template<template<class> class Field, class Type>
198 {
200  (
202  );
203 
204  mag(tRes.ref(), f);
205  return tRes;
206 }
207 
208 template<template<class> class Field, class Type>
210 {
212  (
214  );
215 
216  mag(tRes.ref(), tf());
217  tf.clear();
218  return tRes;
219 }
220 
221 
222 template<template<class> class Field, class Type>
223 void cmptMax
224 (
226  const FieldField<Field, Type>& f
227 )
228 {
229  forAll(cf, i)
230  {
231  cmptMax(cf[i], f[i]);
232  }
233 }
234 
235 template<template<class> class Field, class Type>
237 (
238  const FieldField<Field, Type>& f
239 )
240 {
241  typedef typename FieldField<Field, Type>::cmptType cmptType;
243  (
245  );
246  cmptMax(tRes.ref(), f);
247  return tRes;
248 }
249 
250 template<template<class> class Field, class Type>
252 (
254 )
255 {
256  typedef typename FieldField<Field, Type>::cmptType cmptType;
258  (
260  );
261  cmptMax(tRes.ref(), tf());
262  tf.clear();
263  return tRes;
264 }
265 
266 
267 template<template<class> class Field, class Type>
268 void cmptMin
269 (
271  const FieldField<Field, Type>& f
272 )
273 {
274  forAll(cf, i)
275  {
276  cmptMin(cf[i], f[i]);
277  }
278 }
279 
280 template<template<class> class Field, class Type>
282 (
283  const FieldField<Field, Type>& f
284 )
285 {
286  typedef typename FieldField<Field, Type>::cmptType cmptType;
288  (
290  );
291  cmptMin(tRes.ref(), f);
292  return tRes;
293 }
294 
295 template<template<class> class Field, class Type>
297 (
299 )
300 {
301  typedef typename FieldField<Field, Type>::cmptType cmptType;
303  (
305  );
306  cmptMin(tRes.ref(), tf());
307  tf.clear();
308  return tRes;
309 }
310 
311 
312 template<template<class> class Field, class Type>
313 void cmptAv
314 (
316  const FieldField<Field, Type>& f
317 )
318 {
319  forAll(cf, i)
320  {
321  cmptAv(cf[i], f[i]);
322  }
323 }
324 
325 template<template<class> class Field, class Type>
327 (
328  const FieldField<Field, Type>& f
329 )
330 {
331  typedef typename FieldField<Field, Type>::cmptType cmptType;
333  (
335  );
336  cmptAv(tRes.ref(), f);
337  return tRes;
338 }
339 
340 template<template<class> class Field, class Type>
342 (
344 )
345 {
346  typedef typename FieldField<Field, Type>::cmptType cmptType;
348  (
350  );
351  cmptAv(tRes.ref(), tf());
352  tf.clear();
353  return tRes;
354 }
355 
356 
357 template<template<class> class Field, class Type>
358 void cmptMag
359 (
361  const FieldField<Field, Type>& f
362 )
363 {
364  forAll(cf, i)
365  {
366  cmptMag(cf[i], f[i]);
367  }
368 }
369 
370 template<template<class> class Field, class Type>
372 (
373  const FieldField<Field, Type>& f
374 )
375 {
377  (
379  );
380  cmptMag(tRes.ref(), f);
381  return tRes;
382 }
383 
384 template<template<class> class Field, class Type>
386 (
388 )
389 {
391  cmptMag(tRes.ref(), tf());
392  tf.clear();
393  return tRes;
394 }
395 
396 
397 #define TMP_UNARY_FUNCTION(returnType, func) \
398  \
399 template<template<class> class Field, class Type> \
400 returnType func(const tmp<FieldField<Field, Type>>& tf1) \
401 { \
402  returnType res = func(tf1()); \
403  tf1.clear(); \
404  return res; \
405 }
406 
407 template<template<class> class Field, class Type>
409 {
410  label i = 0;
411  while (i < f.size() && !f[i].size()) i++;
412 
413  if (i < f.size())
414  {
415  Type Max(max(f[i]));
416 
417  for (label j=i+1; j<f.size(); j++)
418  {
419  if (f[j].size())
420  {
421  Max = max(max(f[j]), Max);
422  }
423  }
424 
425  return Max;
426  }
427  else
428  {
429  return pTraits<Type>::min;
430  }
431 }
432 
434 
435 template<template<class> class Field, class Type>
436 Type min(const FieldField<Field, Type>& f)
437 {
438  label i = 0;
439  while (i < f.size() && !f[i].size()) i++;
440 
441  if (i < f.size())
442  {
443  label i = 0;
444  while (!f[i].size()) i++;
445 
446  Type Min(min(f[i]));
447 
448  for (label j=i+1; j<f.size(); j++)
449  {
450  if (f[j].size())
451  {
452  Min = min(min(f[j]), Min);
453  }
454  }
455 
456  return Min;
457  }
458  else
459  {
460  return pTraits<Type>::max;
461  }
462 }
463 
465 
466 template<template<class> class Field, class Type>
467 Type sum(const FieldField<Field, Type>& f)
468 {
469  if (f.size())
470  {
471  Type Sum = Zero;
472 
473  forAll(f, i)
474  {
475  Sum += sum(f[i]);
476  }
477 
478  return Sum;
479  }
480  else
481  {
482  return Zero;
483  }
484 }
485 
487 
488 template<template<class> class Field, class Type>
489 scalar sumMag(const FieldField<Field, Type>& f)
490 {
491  if (f.size())
492  {
493  scalar SumMag = 0.0;
494 
495  forAll(f, i)
496  {
497  SumMag += sumMag(f[i]);
498  }
499 
500  return SumMag;
501  }
502  else
503  {
504  return 0.0;
505  }
506 }
507 
508 TMP_UNARY_FUNCTION(scalar, sumMag)
509 
510 template<template<class> class Field, class Type>
511 Type average(const FieldField<Field, Type>& f)
512 {
513  if (f.size())
514  {
515  label n = 0;
516 
517  forAll(f, i)
518  {
519  n += f[i].size();
520  }
521 
522  if (n == 0)
523  {
525  << "empty fieldField, returning zero" << endl;
526 
527  return Zero;
528  }
529 
530  Type avrg = sum(f)/n;
531 
532  return avrg;
533  }
534  else
535  {
537  << "empty fieldField, returning zero" << endl;
538 
539  return Zero;
540  }
541 }
542 
544 
545 
546 #define G_UNARY_FUNCTION(returnType, gFunc, func, rFunc) \
547  \
548 template<template<class> class Field, class Type> \
549 returnType gFunc(const FieldField<Field, Type>& f) \
550 { \
551  returnType res = func(f); \
552  reduce(res, rFunc##Op<Type>()); \
553  return res; \
554 } \
555 TMP_UNARY_FUNCTION(returnType, gFunc)
556 
558 G_UNARY_FUNCTION(Type, gMin, min, min)
559 G_UNARY_FUNCTION(Type, gSum, sum, sum)
560 G_UNARY_FUNCTION(scalar, gSumMag, sumMag, sum)
561 
562 #undef G_UNARY_FUNCTION
563 
564 
565 template<template<class> class Field, class Type>
567 {
568  label n = 0;
569 
570  forAll(f, i)
571  {
572  n += f[i].size();
573  }
574 
575  reduce(n, sumOp<label>());
576 
577  if (n > 0)
578  {
579  Type avrg = gSum(f)/n;
580 
581  return avrg;
582  }
583  else
584  {
586  << "empty fieldField, returning zero" << endl;
587 
588  return Zero;
589  }
590 }
591 
593 
594 #undef TMP_UNARY_FUNCTION
595 
596 
597 BINARY_FUNCTION(Type, Type, Type, max)
598 BINARY_FUNCTION(Type, Type, Type, min)
599 BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
600 BINARY_FUNCTION(Type, Type, Type, cmptDivide)
601 
602 BINARY_TYPE_FUNCTION(Type, Type, Type, max)
603 BINARY_TYPE_FUNCTION(Type, Type, Type, min)
604 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
605 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
606 
607 
608 /* * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * */
609 
610 UNARY_OPERATOR(Type, Type, -, negate)
611 
612 BINARY_OPERATOR(Type, Type, scalar, *, multiply)
613 BINARY_OPERATOR(Type, scalar, Type, *, multiply)
614 BINARY_OPERATOR(Type, Type, scalar, /, divide)
615 
616 BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, multiply)
617 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, multiply)
618 
619 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, divide)
620 
621 
622 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
623 
624 #define PRODUCT_OPERATOR(product, op, opFunc) \
625  \
626 template \
627 < \
628  template<class> class Field1, \
629  template<class> class Field2, \
630  class Type1, \
631  class Type2 \
632 > \
633 void opFunc \
634 ( \
635  FieldField<Field1, typename product<Type1, Type2>::type>& f, \
636  const FieldField<Field1, Type1>& f1, \
637  const FieldField<Field2, Type2>& f2 \
638 ) \
639 { \
640  forAll(f, i) \
641  { \
642  opFunc(f[i], f1[i], f2[i]); \
643  } \
644 } \
645  \
646 template \
647 < \
648  template<class> class Field1, \
649  template<class> class Field2, \
650  class Type1, \
651  class Type2 \
652 > \
653 tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
654 operator op \
655 ( \
656  const FieldField<Field1, Type1>& f1, \
657  const FieldField<Field2, Type2>& f2 \
658 ) \
659 { \
660  typedef typename product<Type1, Type2>::type productType; \
661  tmp<FieldField<Field1, productType>> tRes \
662  ( \
663  FieldField<Field1, productType>::NewCalculatedType(f1) \
664  ); \
665  opFunc(tRes.ref(), f1, f2); \
666  return tRes; \
667 } \
668  \
669 template<template<class> class Field, class Type1, class Type2> \
670 tmp<FieldField<Field, typename product<Type1, Type2>::type>> \
671 operator op \
672 ( \
673  const FieldField<Field, Type1>& f1, \
674  const tmp<FieldField<Field, Type2>>& tf2 \
675 ) \
676 { \
677  typedef typename product<Type1, Type2>::type productType; \
678  tmp<FieldField<Field, productType>> tRes \
679  ( \
680  reuseTmpFieldField<Field, productType, Type2>::New(tf2) \
681  ); \
682  opFunc(tRes.ref(), f1, tf2()); \
683  tf2.clear(); \
684  return tRes; \
685 } \
686  \
687 template \
688 < \
689  template<class> class Field1, \
690  template<class> class Field2, \
691  class Type1, \
692  class Type2 \
693 > \
694 tmp<FieldField<Field, typename product<Type1, Type2>::type>> \
695 operator op \
696 ( \
697  const FieldField<Field1, Type1>& f1, \
698  const tmp<FieldField<Field2, Type2>>& tf2 \
699 ) \
700 { \
701  typedef typename product<Type1, Type2>::type productType; \
702  tmp<FieldField<Field1, productType>> tRes \
703  ( \
704  FieldField<Field1, productType>::NewCalculatedType(f1) \
705  ); \
706  opFunc(tRes.ref(), f1, tf2()); \
707  tf2.clear(); \
708  return tRes; \
709 } \
710  \
711 template \
712 < \
713  template<class> class Field1, \
714  template<class> class Field2, \
715  class Type1, \
716  class Type2 \
717 > \
718 tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
719 operator op \
720 ( \
721  const tmp<FieldField<Field1, Type1>>& tf1, \
722  const FieldField<Field2, Type2>& f2 \
723 ) \
724 { \
725  typedef typename product<Type1, Type2>::type productType; \
726  tmp<FieldField<Field1, productType>> tRes \
727  ( \
728  reuseTmpFieldField<Field1, productType, Type1>::New(tf1) \
729  ); \
730  opFunc(tRes.ref(), tf1(), f2); \
731  tf1.clear(); \
732  return tRes; \
733 } \
734  \
735 template \
736 < \
737  template<class> class Field1, \
738  template<class> class Field2, \
739  class Type1, \
740  class Type2 \
741 > \
742 tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
743 operator op \
744 ( \
745  const tmp<FieldField<Field1, Type1>>& tf1, \
746  const tmp<FieldField<Field2, Type2>>& tf2 \
747 ) \
748 { \
749  typedef typename product<Type1, Type2>::type productType; \
750  tmp<FieldField<Field1, productType>> tRes \
751  ( \
752  reuseTmpTmpFieldField<Field1, productType, Type1, Type1, Type2>::New \
753  (tf1, tf2) \
754  ); \
755  opFunc(tRes.ref(), tf1(), tf2()); \
756  tf1.clear(); \
757  tf2.clear(); \
758  return tRes; \
759 } \
760  \
761 template \
762 < \
763  template<class> class Field, \
764  class Type, \
765  class Form, \
766  class Cmpt, \
767  direction nCmpt \
768 > \
769 void opFunc \
770 ( \
771  FieldField<Field, typename product<Type, Form>::type>& f, \
772  const FieldField<Field, Type>& f1, \
773  const VectorSpace<Form,Cmpt,nCmpt>& vs \
774 ) \
775 { \
776  forAll(f, i) \
777  { \
778  opFunc(f[i], f1[i], vs); \
779  } \
780 } \
781  \
782 template \
783 < \
784  template<class> class Field, \
785  class Type, \
786  class Form, \
787  class Cmpt, \
788  direction nCmpt \
789 > \
790 tmp<FieldField<Field, typename product<Type, Form>::type>> \
791 operator op \
792 ( \
793  const FieldField<Field, Type>& f1, \
794  const VectorSpace<Form,Cmpt,nCmpt>& vs \
795 ) \
796 { \
797  typedef typename product<Type, Form>::type productType; \
798  tmp<FieldField<Field, productType>> tRes \
799  ( \
800  FieldField<Field, productType>::NewCalculatedType(f1) \
801  ); \
802  opFunc(tRes.ref(), f1, static_cast<const Form&>(vs)); \
803  return tRes; \
804 } \
805  \
806 template \
807 < \
808  template<class> class Field, \
809  class Type, \
810  class Form, \
811  class Cmpt, \
812  direction nCmpt \
813 > \
814 tmp<FieldField<Field, typename product<Type, Form>::type>> \
815 operator op \
816 ( \
817  const tmp<FieldField<Field, Type>>& tf1, \
818  const VectorSpace<Form,Cmpt,nCmpt>& vs \
819 ) \
820 { \
821  typedef typename product<Type, Form>::type productType; \
822  tmp<FieldField<Field, productType>> tRes \
823  ( \
824  reuseTmpFieldField<Field, productType, Type>::New(tf1) \
825  ); \
826  opFunc(tRes.ref(), tf1(), static_cast<const Form&>(vs)); \
827  tf1.clear(); \
828  return tRes; \
829 } \
830  \
831 template \
832 < \
833  template<class> class Field, \
834  class Type, \
835  class Form, \
836  class Cmpt, \
837  direction nCmpt \
838 > \
839 void opFunc \
840 ( \
841  FieldField<Field, typename product<Form, Type>::type>& f, \
842  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
843  const FieldField<Field, Type>& f1 \
844 ) \
845 { \
846  forAll(f, i) \
847  { \
848  opFunc(f[i], vs, f1[i]); \
849  } \
850 } \
851  \
852 template \
853 < \
854  template<class> class Field, \
855  class Type, \
856  class Form, \
857  class Cmpt, \
858  direction nCmpt \
859 > \
860 tmp<FieldField<Field, typename product<Form, Type>::type>> \
861 operator op \
862 ( \
863  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
864  const FieldField<Field, Type>& f1 \
865 ) \
866 { \
867  typedef typename product<Form, Type>::type productType; \
868  tmp<FieldField<Field, productType>> tRes \
869  ( \
870  FieldField<Field, productType>::NewCalculatedType(f1) \
871  ); \
872  opFunc(tRes.ref(), static_cast<const Form&>(vs), f1); \
873  return tRes; \
874 } \
875  \
876 template \
877 < \
878  template<class> class Field, \
879  class Type, \
880  class Form, \
881  class Cmpt, \
882  direction nCmpt \
883 > \
884 tmp<FieldField<Field, typename product<Form, Type>::type>> \
885 operator op \
886 ( \
887  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
888  const tmp<FieldField<Field, Type>>& tf1 \
889 ) \
890 { \
891  typedef typename product<Form, Type>::type productType; \
892  tmp<FieldField<Field, productType>> tRes \
893  ( \
894  reuseTmpFieldField<Field, productType, Type>::New(tf1) \
895  ); \
896  opFunc(tRes.ref(), static_cast<const Form&>(vs), tf1()); \
897  tf1.clear(); \
898  return tRes; \
899 }
900 
903 
908 
909 #undef PRODUCT_OPERATOR
910 
911 
912 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
913 
914 } // End namespace Foam
915 
916 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
917 
918 #include "undefFieldFunctionsM.H"
919 
920 // ************************************************************************* //
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
#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)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
uint8_t direction
Definition: direction.H:46
Inter-processor communication reduction functions.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Traits class for primitives.
Definition: pTraits.H:50
#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:90
#define BINARY_TYPE_FUNCTION(ReturnType, Type1, Type2, Func)
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)
const tensorField & tf
scalar f1
Definition: createFields.H:28
void outer(FieldField< Field1, typename outerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Generic field type.
Definition: FieldField.H:51
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
Type gSum(const FieldField< Field, Type > &f)
Pre-declare SubField and related Field type.
Definition: Field.H:57
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:356
#define TMP_UNARY_FUNCTION(returnType, func)
static const zero Zero
Definition: zero.H:91
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
dimensioned< scalar > magSqr(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)
labelList f(nPoints)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField sf(fieldObject, mesh)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
symmTypeOfRank< typename pTraits< arg1 >::cmptType, arg2 *direction(pTraits< arg1 >::rank) >::type type
Definition: products.H:136
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
#define PRODUCT_OPERATOR(product, op, opFunc)
void negate(FieldField< Field, Type > &res, const FieldField< Field, Type > &f)
pTraits< Type >::cmptType cmptType
Component type.
Definition: FieldField.H:82
#define WarningInFunction
Report a warning using Foam::Warning.
Type gAverage(const FieldField< Field, Type > &f)
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
#define G_UNARY_FUNCTION(returnType, gFunc, func, rFunc)
#define UNARY_OPERATOR(ReturnType, Type1, Op, OpFunc, Dfunc)
A class for managing temporary objects.
Definition: PtrList.H:54
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
dimensioned< scalar > sumMag(const DimensionedField< Type, GeoMesh > &df)