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-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 
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 resultType;
94 
96 
97  pow<Type, r>(tres.ref(), f);
98  return tres;
99 }
100 
101 template<template<class> class Field, class Type, direction r>
103 pow
104 (
105  const tmp<FieldField<Field, Type>>& tf, typename powProduct<Type, r>::type
106 )
107 {
108  typedef typename powProduct<Type, r>::type resultType;
109 
111 
112  pow<Type, r>(tres.ref(), tf());
113  tf.clear();
114  return tres;
115 }
116 
117 
118 template<template<class> class Field, class Type>
119 void sqr
120 (
121  FieldField<Field, typename outerProduct<Type, Type>::type>& f,
122  const FieldField<Field, Type>& vf
123 )
124 {
125  const label loopLen = (f).size();
126 
127  for (label i = 0; i < loopLen; ++i)
128  {
129  sqr(f[i], vf[i]);
130  }
131 }
132 
133 template<template<class> class Field, class Type>
135 sqr(const FieldField<Field, Type>& f)
136 {
137  typedef typename outerProduct<Type, Type>::type resultType;
138 
141  sqr(tres.ref(), f);
142  return tres;
143 }
144 
145 template<template<class> class Field, class Type>
147 sqr(const tmp<FieldField<Field, Type>>& tf)
148 {
149  typedef typename outerProduct<Type, Type>::type resultType;
150 
152 
153  sqr(tres.ref(), tf());
154  tf.clear();
155  return tres;
156 }
157 
158 
159 template<template<class> class Field, class Type>
160 void magSqr
161 (
162  FieldField<Field, typename typeOfMag<Type>::type>& sf,
163  const FieldField<Field, Type>& f
164 )
165 {
166  const label loopLen = (sf).size();
167 
168  for (label i = 0; i < loopLen; ++i)
169  {
170  magSqr(sf[i], f[i]);
171  }
172 }
173 
174 template<template<class> class Field, class Type>
176 magSqr(const FieldField<Field, Type>& f)
177 {
178  typedef typename typeOfMag<Type>::type resultType;
179 
182  magSqr(tres.ref(), f);
183  return tres;
184 }
185 
186 template<template<class> class Field, class Type>
188 magSqr(const tmp<FieldField<Field, Type>>& tf)
189 {
190  typedef typename typeOfMag<Type>::type resultType;
191 
193 
194  magSqr(tres.ref(), tf());
195  tf.clear();
196  return tres;
197 }
198 
199 
200 template<template<class> class Field, class Type>
201 void mag
202 (
203  FieldField<Field, typename typeOfMag<Type>::type>& sf,
204  const FieldField<Field, Type>& f
205 )
206 {
207  const label loopLen = (sf).size();
208 
209  for (label i = 0; i < loopLen; ++i)
210  {
211  mag(sf[i], f[i]);
212  }
213 }
214 
215 template<template<class> class Field, class Type>
217 mag(const FieldField<Field, Type>& f)
218 {
219  typedef typename typeOfMag<Type>::type resultType;
220 
223  mag(tres.ref(), f);
224  return tres;
225 }
226 
227 template<template<class> class Field, class Type>
229 mag(const tmp<FieldField<Field, Type>>& tf)
230 {
231  typedef typename typeOfMag<Type>::type resultType;
232 
234 
235  mag(tres.ref(), tf());
236  tf.clear();
237  return tres;
238 }
239 
240 
241 template<template<class> class Field, class Type>
242 void cmptMax
243 (
244  FieldField<Field, typename FieldField<Field, Type>::cmptType>& cf,
245  const FieldField<Field, Type>& f
246 )
247 {
248  const label loopLen = (cf).size();
249 
250  for (label i = 0; i < loopLen; ++i)
251  {
252  cmptMax(cf[i], f[i]);
253  }
254 }
255 
256 template<template<class> class Field, class Type>
257 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptMax
258 (
259  const FieldField<Field, Type>& f
260 )
261 {
262  typedef typename FieldField<Field, Type>::cmptType resultType;
263 
266  cmptMax(tres.ref(), f);
267  return tres;
268 }
269 
270 template<template<class> class Field, class Type>
271 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptMax
272 (
273  const tmp<FieldField<Field, Type>>& tf
274 )
275 {
276  typedef typename FieldField<Field, Type>::cmptType resultType;
277 
279 
280  cmptMax(tres.ref(), tf());
281  tf.clear();
282  return tres;
283 }
284 
285 
286 template<template<class> class Field, class Type>
287 void cmptMin
288 (
289  FieldField<Field, typename FieldField<Field, Type>::cmptType>& cf,
290  const FieldField<Field, Type>& f
291 )
292 {
293  const label loopLen = (cf).size();
294 
295  for (label i = 0; i < loopLen; ++i)
296  {
297  cmptMin(cf[i], f[i]);
298  }
299 }
300 
301 template<template<class> class Field, class Type>
302 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptMin
303 (
304  const FieldField<Field, Type>& f
305 )
306 {
307  typedef typename FieldField<Field, Type>::cmptType resultType;
308 
311  cmptMin(tres.ref(), f);
312  return tres;
313 }
314 
315 template<template<class> class Field, class Type>
316 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptMin
317 (
318  const tmp<FieldField<Field, Type>>& tf
319 )
320 {
321  typedef typename FieldField<Field, Type>::cmptType resultType;
322 
324 
325  cmptMin(tres.ref(), tf());
326  tf.clear();
327  return tres;
328 }
329 
330 
331 template<template<class> class Field, class Type>
332 void cmptAv
333 (
334  FieldField<Field, typename FieldField<Field, Type>::cmptType>& cf,
335  const FieldField<Field, Type>& f
336 )
337 {
338  const label loopLen = (cf).size();
339 
340  for (label i = 0; i < loopLen; ++i)
341  {
342  cmptAv(cf[i], f[i]);
343  }
344 }
345 
346 template<template<class> class Field, class Type>
347 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptAv
348 (
349  const FieldField<Field, Type>& f
350 )
351 {
352  typedef typename FieldField<Field, Type>::cmptType resultType;
353 
356  cmptAv(tres.ref(), f);
357  return tres;
358 }
359 
360 template<template<class> class Field, class Type>
361 tmp<FieldField<Field, typename FieldField<Field, Type>::cmptType>> cmptAv
362 (
363  const tmp<FieldField<Field, Type>>& tf
364 )
365 {
366  typedef typename FieldField<Field, Type>::cmptType resultType;
367 
369 
370  cmptAv(tres.ref(), tf());
371  tf.clear();
372  return tres;
373 }
374 
375 
376 template<template<class> class Field, class Type>
377 void cmptMag
378 (
379  FieldField<Field, Type>& cf,
380  const FieldField<Field, Type>& f
381 )
382 {
383  const label loopLen = (cf).size();
384 
385  for (label i = 0; i < loopLen; ++i)
386  {
387  cmptMag(cf[i], f[i]);
388  }
389 }
390 
391 template<template<class> class Field, class Type>
392 tmp<FieldField<Field, Type>> cmptMag
393 (
394  const FieldField<Field, Type>& f
395 )
396 {
399  cmptMag(tres.ref(), f);
400  return tres;
401 }
402 
403 template<template<class> class Field, class Type>
404 tmp<FieldField<Field, Type>> cmptMag
405 (
406  const tmp<FieldField<Field, Type>>& tf
407 )
408 {
409  tmp<FieldField<Field, Type>> tres(New(tf));
410  cmptMag(tres.ref(), tf());
411  tf.clear();
412  return tres;
413 }
414 
415 
416 #define TMP_UNARY_FUNCTION(ReturnType, Func) \
417  \
418 template<template<class> class Field, class Type> \
419 ReturnType Func(const tmp<FieldField<Field, Type>>& tf1) \
420 { \
421  ReturnType res = Func(tf1()); \
422  tf1.clear(); \
423  return res; \
424 }
425 
426 template<template<class> class Field, class Type>
427 Type max(const FieldField<Field, Type>& f)
428 {
429  Type result = pTraits<Type>::min;
430 
431  const label loopLen = (f).size();
432 
433  for (label i = 0; i < loopLen; ++i)
434  {
435  if (f[i].size())
436  {
437  result = max(max(f[i]), result);
438  }
439 
440  }
441 
442  return result;
443 }
444 
445 TMP_UNARY_FUNCTION(Type, max)
446 
447 
448 template<template<class> class Field, class Type>
449 Type min(const FieldField<Field, Type>& f)
450 {
451  Type result = pTraits<Type>::max;
452 
453  const label loopLen = (f).size();
454 
455  for (label i = 0; i < loopLen; ++i)
456  {
457  if (f[i].size())
458  {
459  result = min(min(f[i]), result);
460  }
461  }
462 
463  return result;
464 }
465 
466 TMP_UNARY_FUNCTION(Type, min)
467 
468 
469 template<template<class> class Field, class Type>
470 Type sum(const FieldField<Field, Type>& f)
471 {
472  Type result = Zero;
473 
474  const label loopLen = (f).size();
475 
476  for (label i = 0; i < loopLen; ++i)
477  {
478  result += sum(f[i]);
479  }
481  return result;
482 }
483 
484 TMP_UNARY_FUNCTION(Type, sum)
485 
486 template<template<class> class Field, class Type>
487 typename typeOfMag<Type>::type sumMag(const FieldField<Field, Type>& f)
488 {
489  typedef typename typeOfMag<Type>::type resultType;
490 
491  resultType result = Zero;
492 
493  const label loopLen = (f).size();
494 
495  for (label i = 0; i < loopLen; ++i)
496  {
497  result += sumMag(f[i]);
498  }
500  return result;
501 }
502 
504 
505 template<template<class> class Field, class Type>
506 Type average(const FieldField<Field, Type>& f)
507 {
508  label n = 0;
509 
510  const label loopLen = (f).size();
511 
512  for (label i = 0; i < loopLen; ++i)
513  {
514  n += f[i].size();
515  }
516 
517  if (n)
518  {
519  Type avrg = sum(f)/n;
520 
521  return avrg;
522  }
525  << "empty fieldField, returning zero" << endl;
526 
527  return Zero;
528 }
529 
531 
532 
533 template<template<class> class Field, class Type>
534 MinMax<Type> minMax(const FieldField<Field, Type>& f)
535 {
536  MinMax<Type> result;
537 
538  const label loopLen = (f).size();
539 
540  for (label i = 0; i < loopLen; ++i)
541  {
542  result += minMax(f[i]);
543  }
545  return result;
546 }
547 
548 TMP_UNARY_FUNCTION(MinMax<Type>, minMax)
549 
550 template<template<class> class Field, class Type>
551 scalarMinMax minMaxMag(const FieldField<Field, Type>& f)
552 {
553  scalarMinMax result;
554 
555  const label loopLen = (f).size();
556 
557  for (label i = 0; i < loopLen; ++i)
558  {
559  result += minMaxMag(f[i]);
560  }
561 
562  return result;
563 }
564 
566 
567 
568 // With reduction on ReturnType
569 #define G_UNARY_FUNCTION(ReturnType, gFunc, Func, rFunc) \
570  \
571 template<template<class> class Field, class Type> \
572 ReturnType gFunc(const FieldField<Field, Type>& f) \
573 { \
574  ReturnType res = Func(f); \
575  reduce(res, rFunc##Op<ReturnType>()); \
576  return res; \
577 } \
578 TMP_UNARY_FUNCTION(ReturnType, gFunc)
580 G_UNARY_FUNCTION(Type, gMax, max, max)
581 G_UNARY_FUNCTION(Type, gMin, min, min)
582 G_UNARY_FUNCTION(Type, gSum, sum, sum)
586 G_UNARY_FUNCTION(typename typeOfMag<Type>::type, gSumMag, sumMag, sum)
587 
588 #undef G_UNARY_FUNCTION
589 
590 
591 template<template<class> class Field, class Type>
593 {
594  label n = 0;
595 
596  const label loopLen = (f).size();
597 
598  for (label i = 0; i < loopLen; ++i)
599  {
600  n += f[i].size();
601  }
602 
603  reduce(n, sumOp<label>());
604 
605  if (n)
606  {
607  Type avrg = gSum(f)/n;
608 
609  return avrg;
610  }
613  << "Empty FieldField, returning zero" << endl;
614 
615  return Zero;
616 }
617 
619 
620 #undef TMP_UNARY_FUNCTION
621 
622 
623 BINARY_FUNCTION(Type, Type, Type, max)
624 BINARY_FUNCTION(Type, Type, Type, min)
625 BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
626 BINARY_FUNCTION(Type, Type, Type, cmptDivide)
627 
628 BINARY_TYPE_FUNCTION(Type, Type, Type, max)
629 BINARY_TYPE_FUNCTION(Type, Type, Type, min)
630 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
631 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
632 
633 // Note: works with zero_one through implicit conversion to MinMax
634 BINARY_TYPE_FUNCTION_FS(Type, Type, MinMax<Type>, clamp)
635 
636 
637 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
638 
639 TERNARY_FUNCTION(Type, Type, Type, scalar, lerp)
640 TERNARY_TYPE_FUNCTION_FFS(Type, Type, Type, scalar, lerp)
641 
642 
643 /* * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * */
644 
645 UNARY_OPERATOR(Type, Type, -, negate)
646 
647 BINARY_OPERATOR(Type, Type, scalar, *, multiply)
648 BINARY_OPERATOR(Type, scalar, Type, *, multiply)
649 BINARY_OPERATOR(Type, Type, scalar, /, divide)
650 
651 BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, multiply)
652 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, multiply)
653 
654 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, divide)
655 
656 
657 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
658 
659 #define PRODUCT_OPERATOR(product, Op, OpFunc) \
660  \
661 template \
662 < \
663  template<class> class Field1, \
664  template<class> class Field2, \
665  class Type1, \
666  class Type2 \
667 > \
668 void OpFunc \
669 ( \
670  FieldField<Field1, typename product<Type1, Type2>::type>& f, \
671  const FieldField<Field1, Type1>& f1, \
672  const FieldField<Field2, Type2>& f2 \
673 ) \
674 { \
675  const label loopLen = (f).size(); \
676  \
677  for (label i = 0; i < loopLen; ++i) \
678  { \
679  OpFunc(f[i], f1[i], f2[i]); \
680  } \
681 } \
682  \
683 template \
684 < \
685  template<class> class Field1, \
686  template<class> class Field2, \
687  class Type1, \
688  class Type2 \
689 > \
690 tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
691 operator Op \
692 ( \
693  const FieldField<Field1, Type1>& f1, \
694  const FieldField<Field2, Type2>& f2 \
695 ) \
696 { \
697  typedef typename product<Type1, Type2>::type resultType; \
698  auto tres = FieldField<Field1, resultType>::NewCalculatedType(f1); \
699  OpFunc(tres.ref(), f1, f2); \
700  return tres; \
701 } \
702  \
703 template<template<class> class Field, class Type1, class Type2> \
704 tmp<FieldField<Field, typename product<Type1, Type2>::type>> \
705 operator Op \
706 ( \
707  const FieldField<Field, Type1>& f1, \
708  const tmp<FieldField<Field, Type2>>& tf2 \
709 ) \
710 { \
711  typedef typename product<Type1, Type2>::type resultType; \
712  auto tres = reuseTmpFieldField<Field, resultType, Type2>::New(tf2); \
713  OpFunc(tres.ref(), f1, tf2()); \
714  tf2.clear(); \
715  return tres; \
716 } \
717  \
718 template \
719 < \
720  template<class> class Field1, \
721  template<class> class Field2, \
722  class Type1, \
723  class Type2 \
724 > \
725 tmp<FieldField<Field, typename product<Type1, Type2>::type>> \
726 operator Op \
727 ( \
728  const FieldField<Field1, Type1>& f1, \
729  const tmp<FieldField<Field2, Type2>>& tf2 \
730 ) \
731 { \
732  typedef typename product<Type1, Type2>::type resultType; \
733  auto tres = FieldField<Field1, resultType>::NewCalculatedType(f1); \
734  OpFunc(tres.ref(), f1, tf2()); \
735  tf2.clear(); \
736  return tres; \
737 } \
738  \
739 template \
740 < \
741  template<class> class Field1, \
742  template<class> class Field2, \
743  class Type1, \
744  class Type2 \
745 > \
746 tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
747 operator Op \
748 ( \
749  const tmp<FieldField<Field1, Type1>>& tf1, \
750  const FieldField<Field2, Type2>& f2 \
751 ) \
752 { \
753  typedef typename product<Type1, Type2>::type resultType; \
754  auto tres = reuseTmpFieldField<Field1, resultType, Type1>::New(tf1); \
755  OpFunc(tres.ref(), tf1(), f2); \
756  tf1.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<Field1, typename product<Type1, Type2>::type>> \
768 operator Op \
769 ( \
770  const tmp<FieldField<Field1, Type1>>& tf1, \
771  const tmp<FieldField<Field2, Type2>>& tf2 \
772 ) \
773 { \
774  typedef typename product<Type1, Type2>::type resultType; \
775  auto tres \
776  ( \
777  reuseTmpTmpFieldField<Field1, resultType, Type1, Type1, Type2>::New \
778  (tf1, tf2) \
779  ); \
780  OpFunc(tres.ref(), tf1(), tf2()); \
781  tf1.clear(); \
782  tf2.clear(); \
783  return tres; \
784 } \
785  \
786 template \
787 < \
788  template<class> class Field, \
789  class Type, \
790  class Form, \
791  class Cmpt, \
792  direction nCmpt \
793 > \
794 void OpFunc \
795 ( \
796  FieldField<Field, typename product<Type, Form>::type>& result, \
797  const FieldField<Field, Type>& f1, \
798  const VectorSpace<Form,Cmpt,nCmpt>& vs \
799 ) \
800 { \
801  const label loopLen = (result).size(); \
802  \
803  for (label i = 0; i < loopLen; ++i) \
804  { \
805  OpFunc(result[i], f1[i], vs); \
806  } \
807 } \
808  \
809 template \
810 < \
811  template<class> class Field, \
812  class Type, \
813  class Form, \
814  class Cmpt, \
815  direction nCmpt \
816 > \
817 tmp<FieldField<Field, typename product<Type, Form>::type>> \
818 operator Op \
819 ( \
820  const FieldField<Field, Type>& f1, \
821  const VectorSpace<Form,Cmpt,nCmpt>& vs \
822 ) \
823 { \
824  typedef typename product<Type, Form>::type resultType; \
825  auto tres = FieldField<Field, resultType>::NewCalculatedType(f1); \
826  OpFunc(tres.ref(), f1, static_cast<const Form&>(vs)); \
827  return tres; \
828 } \
829  \
830 template \
831 < \
832  template<class> class Field, \
833  class Type, \
834  class Form, \
835  class Cmpt, \
836  direction nCmpt \
837 > \
838 tmp<FieldField<Field, typename product<Type, Form>::type>> \
839 operator Op \
840 ( \
841  const tmp<FieldField<Field, Type>>& tf1, \
842  const VectorSpace<Form,Cmpt,nCmpt>& vs \
843 ) \
844 { \
845  typedef typename product<Type, Form>::type resultType; \
846  auto tres = reuseTmpFieldField<Field, resultType, Type>::New(tf1); \
847  OpFunc(tres.ref(), tf1(), static_cast<const Form&>(vs)); \
848  tf1.clear(); \
849  return tres; \
850 } \
851  \
852 template \
853 < \
854  template<class> class Field, \
855  class Type, \
856  class Form, \
857  class Cmpt, \
858  direction nCmpt \
859 > \
860 void OpFunc \
861 ( \
862  FieldField<Field, typename product<Form, Type>::type>& result, \
863  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
864  const FieldField<Field, Type>& f1 \
865 ) \
866 { \
867  const label loopLen = (result).size(); \
868  \
869  for (label i = 0; i < loopLen; ++i) \
870  { \
871  OpFunc(result[i], vs, f1[i]); \
872  } \
873 } \
874  \
875 template \
876 < \
877  template<class> class Field, \
878  class Type, \
879  class Form, \
880  class Cmpt, \
881  direction nCmpt \
882 > \
883 tmp<FieldField<Field, typename product<Form, Type>::type>> \
884 operator Op \
885 ( \
886  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
887  const FieldField<Field, Type>& f1 \
888 ) \
889 { \
890  typedef typename product<Form, Type>::type resultType; \
891  auto tres = FieldField<Field, resultType>::NewCalculatedType(f1); \
892  OpFunc(tres.ref(), static_cast<const Form&>(vs), f1); \
893  return tres; \
894 } \
895  \
896 template \
897 < \
898  template<class> class Field, \
899  class Type, \
900  class Form, \
901  class Cmpt, \
902  direction nCmpt \
903 > \
904 tmp<FieldField<Field, typename product<Form, Type>::type>> \
905 operator Op \
906 ( \
907  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
908  const tmp<FieldField<Field, Type>>& tf1 \
909 ) \
910 { \
911  typedef typename product<Form, Type>::type resultType; \
912  auto tres = reuseTmpFieldField<Field, resultType, Type>::New(tf1); \
913  OpFunc(tres.ref(), static_cast<const Form&>(vs), tf1()); \
914  tf1.clear(); \
915  return tres; \
916 }
920 
925 
926 #undef PRODUCT_OPERATOR
927 
928 
929 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
930 
931 } // End namespace Foam
932 
933 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
934 
935 #include "undefFieldFunctionsM.H"
936 
937 // ************************************************************************* //
dimensioned< typename typeOfMag< Type >::type > sumMag(const DimensionedField< Type, GeoMesh > &f1)
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
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)
uint8_t direction
Definition: direction.H:46
Inter-processor communication reduction functions.
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:97
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)
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:531
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
#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)
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 > &)
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:54
#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:799
Type gSum(const FieldField< Field, Type > &f)
Generic templated field type.
Definition: Field.H:62
Type max(const tmp< FieldField< Field, Type >> &tf1)
#define TMP_UNARY_FUNCTION(ReturnType, Func)
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
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)
#define G_UNARY_FUNCTION(ReturnType, gFunc, Func, rFunc)
scalarMinMax gMinMaxMag(const FieldField< Field, Type > &f)
labelList f(nPoints)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
symmTypeOfRank< typename pTraits< arg1 >::cmptType, arg2 *direction(pTraits< arg1 >::rank) >::type type
Definition: products.H:176
#define PRODUCT_OPERATOR(product, Op, OpFunc)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
pTraits< Type >::cmptType cmptType
Component type.
Definition: FieldField.H:83
dimensioned< scalarMinMax > minMaxMag(const DimensionedField< Type, GeoMesh > &f1)
#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.
BINARY_TYPE_FUNCTION_FS(Type, Type, MinMax< Type >, clip)
label n
static tmp< FieldField< Field, TypeR > > New(const FieldField< Field, Type1 > &f1)
Pass-through to NewCalculatedType.
#define UNARY_OPERATOR(ReturnType, Type1, Op, OpFunc, Dfunc)
#define TERNARY_TYPE_FUNCTION_FFS(ReturnType, Type1, Type2, Type3, Func)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
#define TERNARY_FUNCTION(ReturnType, Type1, Type2, Type3, Func)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
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
pTraits< typename pTraits< arg1 >::cmptType >::magType type
Definition: products.H:96
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127