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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2022 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 
29 #include "PstreamReduceOps.H"
31 
32 #define TEMPLATE template<template<class> class Field, class Type>
33 #include "FieldFieldFunctionsM.C"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 /* * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * */
41 
42 template<template<class> class Field, class Type>
43 void component
44 (
47  const direction d
48 )
49 {
50  const label loopLen = (sf).size();
51 
52  for (label i = 0; i < loopLen; ++i)
53  {
54  component(sf[i], f[i], d);
55  }
56 }
57 
58 
59 template<template<class> class Field, class Type>
60 void T(FieldField<Field, Type>& f1, const FieldField<Field, Type>& f2)
61 {
62  const label loopLen = (f1).size();
63 
64  for (label i = 0; i < loopLen; ++i)
65  {
66  T(f1[i], f2[i]);
67  }
68 }
69 
70 
71 template<template<class> class Field, class Type, direction r>
72 void pow
73 (
74  FieldField<Field, typename powProduct<Type, r>::type>& f,
75  const FieldField<Field, Type>& vf
76 )
77 {
78  const label loopLen = (f).size();
79 
80  for (label i = 0; i < loopLen; ++i)
81  {
82  pow(f[i], vf[i]);
83  }
84 }
85 
86 template<template<class> class Field, class Type, direction r>
88 pow
89 (
90  const FieldField<Field, Type>& f, typename powProduct<Type, r>::type
91 )
92 {
93  typedef typename powProduct<Type, r>::type powProductType;
94 
95  auto tres
96  (
98  );
99 
100  pow<Type, r>(tres.ref(), f);
101  return tres;
102 }
103 
104 template<template<class> class Field, class Type, direction r>
106 pow
107 (
108  const tmp<FieldField<Field, Type>>& tf, typename powProduct<Type, r>::type
109 )
110 {
111  typedef typename powProduct<Type, r>::type powProductType;
112 
113  auto tres
114  (
116  );
117 
118  pow<Type, r>(tres.ref(), tf());
119  tf.clear();
120  return tres;
121 }
122 
123 
124 template<template<class> class Field, class Type>
125 void sqr
126 (
127  FieldField<Field, typename outerProduct<Type, Type>::type>& f,
128  const FieldField<Field, Type>& vf
129 )
130 {
131  const label loopLen = (f).size();
132 
133  for (label i = 0; i < loopLen; ++i)
134  {
135  sqr(f[i], vf[i]);
136  }
137 }
138 
139 template<template<class> class Field, class Type>
141 sqr(const FieldField<Field, Type>& f)
142 {
143  typedef typename outerProduct<Type, Type>::type outerProductType;
144  tmp<FieldField<Field, outerProductType>> tres
145  (
147  );
148  sqr(tres.ref(), f);
149  return tres;
150 }
151 
152 template<template<class> class Field, class Type>
154 sqr(const tmp<FieldField<Field, Type>>& tf)
155 {
156  typedef typename outerProduct<Type, Type>::type outerProductType;
157 
158  auto tres
159  (
161  );
162 
163  sqr(tres.ref(), tf());
164  tf.clear();
165  return tres;
166 }
167 
168 
169 template<template<class> class Field, class Type>
170 void magSqr
171 (
172  FieldField<Field, typename typeOfMag<Type>::type>& sf,
173  const FieldField<Field, Type>& f
174 )
175 {
176  const label loopLen = (sf).size();
177 
178  for (label i = 0; i < loopLen; ++i)
179  {
180  magSqr(sf[i], f[i]);
181  }
182 }
183 
184 template<template<class> class Field, class Type>
186 magSqr(const FieldField<Field, Type>& f)
187 {
188  typedef typename typeOfMag<Type>::type magType;
189 
190  auto tres
191  (
193  );
195  magSqr(tres.ref(), f);
196  return tres;
197 }
198 
199 template<template<class> class Field, class Type>
201 magSqr(const tmp<FieldField<Field, Type>>& tf)
202 {
203  typedef typename typeOfMag<Type>::type magType;
204 
205  auto tres
206  (
208  );
209 
210  magSqr(tres.ref(), tf());
211  tf.clear();
212  return tres;
213 }
214 
215 
216 template<template<class> class Field, class Type>
217 void mag
218 (
219  FieldField<Field, typename typeOfMag<Type>::type>& sf,
220  const FieldField<Field, Type>& f
221 )
222 {
223  const label loopLen = (sf).size();
224 
225  for (label i = 0; i < loopLen; ++i)
226  {
227  mag(sf[i], f[i]);
228  }
229 }
230 
231 template<template<class> class Field, class Type>
233 mag(const FieldField<Field, Type>& f)
234 {
235  typedef typename typeOfMag<Type>::type magType;
236 
237  auto tres
238  (
240  );
242  mag(tres.ref(), f);
243  return tres;
244 }
245 
246 template<template<class> class Field, class Type>
248 mag(const tmp<FieldField<Field, Type>>& tf)
249 {
250  typedef typename typeOfMag<Type>::type magType;
251 
252  auto tres
253  (
255  );
256 
257  mag(tres.ref(), tf());
258  tf.clear();
259  return tres;
260 }
261 
262 
263 template<template<class> class Field, class Type>
264 void cmptMax
265 (
266  FieldField<Field, typename FieldField<Field, Type>::cmptType>& cf,
267  const FieldField<Field, Type>& f
268 )
269 {
270  const label loopLen = (cf).size();
271 
272  for (label i = 0; i < loopLen; ++i)
273  {
274  cmptMax(cf[i], f[i]);
275  }
276 }
277 
278 template<template<class> class Field, class Type>
279 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptMax
280 (
281  const FieldField<Field, Type>& f
282 )
283 {
284  typedef typename FieldField<Field, Type>::cmptType cmptType;
285 
286  auto tres
287  (
289  );
291  cmptMax(tres.ref(), f);
292  return tres;
293 }
294 
295 template<template<class> class Field, class Type>
296 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptMax
297 (
298  const tmp<FieldField<Field, Type>>& tf
299 )
300 {
301  typedef typename FieldField<Field, Type>::cmptType cmptType;
302 
303  auto tres
304  (
306  );
307 
308  cmptMax(tres.ref(), tf());
309  tf.clear();
310  return tres;
311 }
312 
313 
314 template<template<class> class Field, class Type>
315 void cmptMin
316 (
317  FieldField<Field, typename FieldField<Field, Type>::cmptType>& cf,
318  const FieldField<Field, Type>& f
319 )
320 {
321  const label loopLen = (cf).size();
322 
323  for (label i = 0; i < loopLen; ++i)
324  {
325  cmptMin(cf[i], f[i]);
326  }
327 }
328 
329 template<template<class> class Field, class Type>
330 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptMin
331 (
332  const FieldField<Field, Type>& f
333 )
334 {
335  typedef typename FieldField<Field, Type>::cmptType cmptType;
336 
337  auto tres
338  (
340  );
342  cmptMin(tres.ref(), f);
343  return tres;
344 }
345 
346 template<template<class> class Field, class Type>
347 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptMin
348 (
349  const tmp<FieldField<Field, Type>>& tf
350 )
351 {
352  typedef typename FieldField<Field, Type>::cmptType cmptType;
353 
354  auto tres
355  (
357  );
358 
359  cmptMin(tres.ref(), tf());
360  tf.clear();
361  return tres;
362 }
363 
364 
365 template<template<class> class Field, class Type>
366 void cmptAv
367 (
368  FieldField<Field, typename FieldField<Field, Type>::cmptType>& cf,
369  const FieldField<Field, Type>& f
370 )
371 {
372  const label loopLen = (cf).size();
373 
374  for (label i = 0; i < loopLen; ++i)
375  {
376  cmptAv(cf[i], f[i]);
377  }
378 }
379 
380 template<template<class> class Field, class Type>
381 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptAv
382 (
383  const FieldField<Field, Type>& f
384 )
385 {
386  typedef typename FieldField<Field, Type>::cmptType cmptType;
387 
388  auto tres
389  (
391  );
393  cmptAv(tres.ref(), f);
394  return tres;
395 }
396 
397 template<template<class> class Field, class Type>
398 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptAv
399 (
400  const tmp<FieldField<Field, Type>>& tf
401 )
402 {
403  typedef typename FieldField<Field, Type>::cmptType cmptType;
404 
405  auto tres
406  (
408  );
409 
410  cmptAv(tres.ref(), tf());
411  tf.clear();
412  return tres;
413 }
414 
415 
416 template<template<class> class Field, class Type>
417 void cmptMag
418 (
419  FieldField<Field, Type>& cf,
420  const FieldField<Field, Type>& f
421 )
422 {
423  const label loopLen = (cf).size();
424 
425  for (label i = 0; i < loopLen; ++i)
426  {
427  cmptMag(cf[i], f[i]);
428  }
429 }
430 
431 template<template<class> class Field, class Type>
432 tmp<FieldField<Field, Type>> cmptMag
433 (
434  const FieldField<Field, Type>& f
435 )
436 {
437  auto tres
438  (
440  );
442  cmptMag(tres.ref(), f);
443  return tres;
444 }
445 
446 template<template<class> class Field, class Type>
447 tmp<FieldField<Field, Type>> cmptMag
448 (
449  const tmp<FieldField<Field, Type>>& tf
450 )
451 {
452  tmp<FieldField<Field, Type>> tres(New(tf));
453  cmptMag(tres.ref(), tf());
454  tf.clear();
455  return tres;
456 }
457 
458 
459 #define TMP_UNARY_FUNCTION(returnType, func) \
460  \
461 template<template<class> class Field, class Type> \
462 returnType func(const tmp<FieldField<Field, Type>>& tf1) \
463 { \
464  returnType res = func(tf1()); \
465  tf1.clear(); \
466  return res; \
467 }
468 
469 template<template<class> class Field, class Type>
470 Type max(const FieldField<Field, Type>& f)
471 {
472  Type result = pTraits<Type>::min;
473 
474  const label loopLen = (f).size();
475 
476  for (label i = 0; i < loopLen; ++i)
477  {
478  if (f[i].size())
479  {
480  result = max(max(f[i]), result);
481  }
482 
483  }
484 
485  return result;
486 }
487 
488 TMP_UNARY_FUNCTION(Type, max)
489 
490 
491 template<template<class> class Field, class Type>
492 Type min(const FieldField<Field, Type>& f)
493 {
494  Type result = pTraits<Type>::max;
495 
496  const label loopLen = (f).size();
497 
498  for (label i = 0; i < loopLen; ++i)
499  {
500  if (f[i].size())
501  {
502  result = min(min(f[i]), result);
503  }
504  }
505 
506  return result;
507 }
508 
509 TMP_UNARY_FUNCTION(Type, min)
510 
511 
512 template<template<class> class Field, class Type>
513 Type sum(const FieldField<Field, Type>& f)
514 {
515  Type result = Zero;
516 
517  const label loopLen = (f).size();
518 
519  for (label i = 0; i < loopLen; ++i)
520  {
521  result += sum(f[i]);
522  }
524  return result;
525 }
526 
527 TMP_UNARY_FUNCTION(Type, sum)
528 
529 template<template<class> class Field, class Type>
530 typename typeOfMag<Type>::type sumMag(const FieldField<Field, Type>& f)
531 {
532  typedef typename typeOfMag<Type>::type magType;
533 
534  magType result = Zero;
535 
536  const label loopLen = (f).size();
537 
538  for (label i = 0; i < loopLen; ++i)
539  {
540  result += sumMag(f[i]);
541  }
543  return result;
544 }
545 
547 
548 template<template<class> class Field, class Type>
549 Type average(const FieldField<Field, Type>& f)
550 {
551  label n = 0;
552 
553  const label loopLen = (f).size();
554 
555  for (label i = 0; i < loopLen; ++i)
556  {
557  n += f[i].size();
558  }
559 
560  if (n)
561  {
562  Type avrg = sum(f)/n;
563 
564  return avrg;
565  }
568  << "empty fieldField, returning zero" << endl;
569 
570  return Zero;
571 }
572 
574 
575 
576 template<template<class> class Field, class Type>
577 MinMax<Type> minMax(const FieldField<Field, Type>& f)
578 {
579  MinMax<Type> result;
580 
581  const label loopLen = (f).size();
582 
583  for (label i = 0; i < loopLen; ++i)
584  {
585  result += minMax(f[i]);
586  }
588  return result;
589 }
590 
591 TMP_UNARY_FUNCTION(MinMax<Type>, minMax)
592 
593 template<template<class> class Field, class Type>
594 scalarMinMax minMaxMag(const FieldField<Field, Type>& f)
595 {
596  scalarMinMax result;
597 
598  const label loopLen = (f).size();
599 
600  for (label i = 0; i < loopLen; ++i)
601  {
602  result += minMaxMag(f[i]);
603  }
604 
605  return result;
606 }
607 
609 
610 
611 // With reduction on ReturnType
612 #define G_UNARY_FUNCTION(ReturnType, gFunc, func, rFunc) \
613  \
614 template<template<class> class Field, class Type> \
615 ReturnType gFunc(const FieldField<Field, Type>& f) \
616 { \
617  ReturnType res = func(f); \
618  reduce(res, rFunc##Op<ReturnType>()); \
619  return res; \
620 } \
621 TMP_UNARY_FUNCTION(ReturnType, gFunc)
623 G_UNARY_FUNCTION(Type, gMax, max, max)
624 G_UNARY_FUNCTION(Type, gMin, min, min)
625 G_UNARY_FUNCTION(Type, gSum, sum, sum)
629 G_UNARY_FUNCTION(typename typeOfMag<Type>::type, gSumMag, sumMag, sum)
630 
631 #undef G_UNARY_FUNCTION
632 
633 
634 template<template<class> class Field, class Type>
636 {
637  label n = 0;
638 
639  const label loopLen = (f).size();
640 
641  for (label i = 0; i < loopLen; ++i)
642  {
643  n += f[i].size();
644  }
645 
646  reduce(n, sumOp<label>());
647 
648  if (n)
649  {
650  Type avrg = gSum(f)/n;
651 
652  return avrg;
653  }
656  << "empty fieldField, returning zero" << endl;
657 
658  return Zero;
659 }
660 
662 
663 #undef TMP_UNARY_FUNCTION
664 
665 
666 BINARY_FUNCTION(Type, Type, Type, max)
667 BINARY_FUNCTION(Type, Type, Type, min)
668 BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
669 BINARY_FUNCTION(Type, Type, Type, cmptDivide)
670 
671 BINARY_TYPE_FUNCTION(Type, Type, Type, max)
672 BINARY_TYPE_FUNCTION(Type, Type, Type, min)
673 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
674 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
675 
676 BINARY_TYPE_FUNCTION_FS(Type, Type, MinMax<Type>, clip)
677 
678 
679 /* * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * */
680 
681 UNARY_OPERATOR(Type, Type, -, negate)
682 
683 BINARY_OPERATOR(Type, Type, scalar, *, multiply)
684 BINARY_OPERATOR(Type, scalar, Type, *, multiply)
685 BINARY_OPERATOR(Type, Type, scalar, /, divide)
686 
687 BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, multiply)
688 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, multiply)
689 
690 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, divide)
691 
692 
693 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
694 
695 #define PRODUCT_OPERATOR(product, op, opFunc) \
696  \
697 template \
698 < \
699  template<class> class Field1, \
700  template<class> class Field2, \
701  class Type1, \
702  class Type2 \
703 > \
704 void opFunc \
705 ( \
706  FieldField<Field1, typename product<Type1, Type2>::type>& f, \
707  const FieldField<Field1, Type1>& f1, \
708  const FieldField<Field2, Type2>& f2 \
709 ) \
710 { \
711  const label loopLen = (f).size(); \
712  \
713  for (label i = 0; i < loopLen; ++i) \
714  { \
715  opFunc(f[i], f1[i], f2[i]); \
716  } \
717 } \
718  \
719 template \
720 < \
721  template<class> class Field1, \
722  template<class> class Field2, \
723  class Type1, \
724  class Type2 \
725 > \
726 tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
727 operator op \
728 ( \
729  const FieldField<Field1, Type1>& f1, \
730  const FieldField<Field2, Type2>& f2 \
731 ) \
732 { \
733  typedef typename product<Type1, Type2>::type productType; \
734  tmp<FieldField<Field1, productType>> tres \
735  ( \
736  FieldField<Field1, productType>::NewCalculatedType(f1) \
737  ); \
738  opFunc(tres.ref(), f1, f2); \
739  return tres; \
740 } \
741  \
742 template<template<class> class Field, class Type1, class Type2> \
743 tmp<FieldField<Field, typename product<Type1, Type2>::type>> \
744 operator op \
745 ( \
746  const FieldField<Field, Type1>& f1, \
747  const tmp<FieldField<Field, Type2>>& tf2 \
748 ) \
749 { \
750  typedef typename product<Type1, Type2>::type productType; \
751  tmp<FieldField<Field, productType>> tres \
752  ( \
753  reuseTmpFieldField<Field, productType, Type2>::New(tf2) \
754  ); \
755  opFunc(tres.ref(), f1, tf2()); \
756  tf2.clear(); \
757  return tres; \
758 } \
759  \
760 template \
761 < \
762  template<class> class Field1, \
763  template<class> class Field2, \
764  class Type1, \
765  class Type2 \
766 > \
767 tmp<FieldField<Field, typename product<Type1, Type2>::type>> \
768 operator op \
769 ( \
770  const FieldField<Field1, Type1>& f1, \
771  const tmp<FieldField<Field2, Type2>>& tf2 \
772 ) \
773 { \
774  typedef typename product<Type1, Type2>::type productType; \
775  tmp<FieldField<Field1, productType>> tres \
776  ( \
777  FieldField<Field1, productType>::NewCalculatedType(f1) \
778  ); \
779  opFunc(tres.ref(), f1, tf2()); \
780  tf2.clear(); \
781  return tres; \
782 } \
783  \
784 template \
785 < \
786  template<class> class Field1, \
787  template<class> class Field2, \
788  class Type1, \
789  class Type2 \
790 > \
791 tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
792 operator op \
793 ( \
794  const tmp<FieldField<Field1, Type1>>& tf1, \
795  const FieldField<Field2, Type2>& f2 \
796 ) \
797 { \
798  typedef typename product<Type1, Type2>::type productType; \
799  tmp<FieldField<Field1, productType>> tres \
800  ( \
801  reuseTmpFieldField<Field1, productType, Type1>::New(tf1) \
802  ); \
803  opFunc(tres.ref(), tf1(), f2); \
804  tf1.clear(); \
805  return tres; \
806 } \
807  \
808 template \
809 < \
810  template<class> class Field1, \
811  template<class> class Field2, \
812  class Type1, \
813  class Type2 \
814 > \
815 tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
816 operator op \
817 ( \
818  const tmp<FieldField<Field1, Type1>>& tf1, \
819  const tmp<FieldField<Field2, Type2>>& tf2 \
820 ) \
821 { \
822  typedef typename product<Type1, Type2>::type productType; \
823  tmp<FieldField<Field1, productType>> tres \
824  ( \
825  reuseTmpTmpFieldField<Field1, productType, Type1, Type1, Type2>::New \
826  (tf1, tf2) \
827  ); \
828  opFunc(tres.ref(), tf1(), tf2()); \
829  tf1.clear(); \
830  tf2.clear(); \
831  return tres; \
832 } \
833  \
834 template \
835 < \
836  template<class> class Field, \
837  class Type, \
838  class Form, \
839  class Cmpt, \
840  direction nCmpt \
841 > \
842 void opFunc \
843 ( \
844  FieldField<Field, typename product<Type, Form>::type>& f, \
845  const FieldField<Field, Type>& f1, \
846  const VectorSpace<Form,Cmpt,nCmpt>& vs \
847 ) \
848 { \
849  const label loopLen = (f).size(); \
850  \
851  for (label i = 0; i < loopLen; ++i) \
852  { \
853  opFunc(f[i], f1[i], vs); \
854  } \
855 } \
856  \
857 template \
858 < \
859  template<class> class Field, \
860  class Type, \
861  class Form, \
862  class Cmpt, \
863  direction nCmpt \
864 > \
865 tmp<FieldField<Field, typename product<Type, Form>::type>> \
866 operator op \
867 ( \
868  const FieldField<Field, Type>& f1, \
869  const VectorSpace<Form,Cmpt,nCmpt>& vs \
870 ) \
871 { \
872  typedef typename product<Type, Form>::type productType; \
873  tmp<FieldField<Field, productType>> tres \
874  ( \
875  FieldField<Field, productType>::NewCalculatedType(f1) \
876  ); \
877  opFunc(tres.ref(), f1, static_cast<const Form&>(vs)); \
878  return tres; \
879 } \
880  \
881 template \
882 < \
883  template<class> class Field, \
884  class Type, \
885  class Form, \
886  class Cmpt, \
887  direction nCmpt \
888 > \
889 tmp<FieldField<Field, typename product<Type, Form>::type>> \
890 operator op \
891 ( \
892  const tmp<FieldField<Field, Type>>& tf1, \
893  const VectorSpace<Form,Cmpt,nCmpt>& vs \
894 ) \
895 { \
896  typedef typename product<Type, Form>::type productType; \
897  tmp<FieldField<Field, productType>> tres \
898  ( \
899  reuseTmpFieldField<Field, productType, Type>::New(tf1) \
900  ); \
901  opFunc(tres.ref(), tf1(), static_cast<const Form&>(vs)); \
902  tf1.clear(); \
903  return tres; \
904 } \
905  \
906 template \
907 < \
908  template<class> class Field, \
909  class Type, \
910  class Form, \
911  class Cmpt, \
912  direction nCmpt \
913 > \
914 void opFunc \
915 ( \
916  FieldField<Field, typename product<Form, Type>::type>& f, \
917  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
918  const FieldField<Field, Type>& f1 \
919 ) \
920 { \
921  const label loopLen = (f).size(); \
922  \
923  for (label i = 0; i < loopLen; ++i) \
924  { \
925  opFunc(f[i], vs, f1[i]); \
926  } \
927 } \
928  \
929 template \
930 < \
931  template<class> class Field, \
932  class Type, \
933  class Form, \
934  class Cmpt, \
935  direction nCmpt \
936 > \
937 tmp<FieldField<Field, typename product<Form, Type>::type>> \
938 operator op \
939 ( \
940  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
941  const FieldField<Field, Type>& f1 \
942 ) \
943 { \
944  typedef typename product<Form, Type>::type productType; \
945  tmp<FieldField<Field, productType>> tres \
946  ( \
947  FieldField<Field, productType>::NewCalculatedType(f1) \
948  ); \
949  opFunc(tres.ref(), static_cast<const Form&>(vs), f1); \
950  return tres; \
951 } \
952  \
953 template \
954 < \
955  template<class> class Field, \
956  class Type, \
957  class Form, \
958  class Cmpt, \
959  direction nCmpt \
960 > \
961 tmp<FieldField<Field, typename product<Form, Type>::type>> \
962 operator op \
963 ( \
964  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
965  const tmp<FieldField<Field, Type>>& tf1 \
966 ) \
967 { \
968  typedef typename product<Form, Type>::type productType; \
969  tmp<FieldField<Field, productType>> tres \
970  ( \
971  reuseTmpFieldField<Field, productType, Type>::New(tf1) \
972  ); \
973  opFunc(tres.ref(), static_cast<const Form&>(vs), tf1()); \
974  tf1.clear(); \
975  return tres; \
976 }
980 
985 
986 #undef PRODUCT_OPERATOR
987 
988 
989 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
990 
991 } // End namespace Foam
992 
993 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
994 
995 #include "undefFieldFunctionsM.H"
996 
997 // ************************************************************************* //
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
type
Types of root.
Definition: Roots.H:52
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
uint8_t direction
Definition: direction.H:46
Inter-processor communication reduction functions.
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:111
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type gMin(const FieldField< Field, Type > &f)
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)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
A traits class, which is primarily used for primitives.
Definition: pTraits.H:50
#define BINARY_TYPE_OPERATOR_FS(ReturnType, Type1, Type2, Op, OpName, OpFunc)
The magnitude type for given argument.
Definition: products.H:92
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)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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)
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:54
#define G_UNARY_FUNCTION(ReturnType, gFunc, func, rFunc)
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
Type gSum(const FieldField< Field, Type > &f)
Generic templated field type.
Definition: Field.H:61
dimensionSet clip(const dimensionSet &ds1, const dimensionSet &ds2)
Definition: dimensionSet.C:271
#define TMP_UNARY_FUNCTION(returnType, func)
Type max(const tmp< FieldField< Field, Type >> &tf1)
#define BINARY_TYPE_FUNCTION_FS(ReturnType, Type1, Type2, Func)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
static tmp< FieldField< Field, Type > > NewCalculatedType(const FieldField< Field, Type2 > &ff)
Return a pointer to a new calculatedFvPatchFieldField created on.
Definition: FieldField.C:191
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
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)
labelList f(nPoints)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
dimensioned< typename typeOfMag< Type >::type > sumMag(const DimensionedField< Type, GeoMesh > &df)
symmTypeOfRank< typename pTraits< arg1 >::cmptType, arg2 *direction(pTraits< arg1 >::rank) >::type type
Definition: products.H:176
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:83
#define WarningInFunction
Report a warning using Foam::Warning.
Type gAverage(const FieldField< Field, Type > &f)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
label n
#define UNARY_OPERATOR(ReturnType, Type1, Op, OpFunc, Dfunc)
static tmp< FieldField< Field, TypeR > > New(const tmp< FieldField< Field, Type1 >> &tf1)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensioned< scalarMinMax > minMaxMag(const DimensionedField< Type, GeoMesh > &df)
Namespace for OpenFOAM.
pTraits< typename pTraits< arg1 >::cmptType >::magType type
Definition: products.H:96
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157