GeometricFieldFunctions.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) 2018-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 \
32  template<class Type, template<class> class PatchField, class GeoMesh>
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
41 
42 template<class Type, template<class> class PatchField, class GeoMesh>
43 void component
44 (
46  <
48  PatchField,
49  GeoMesh
50  >& result,
52  const direction d
53 )
54 {
55  component(result.primitiveFieldRef(), f1.primitiveField(), d);
57  result.oriented() = f1.oriented();
58 }
59 
60 
61 template<class Type, template<class> class PatchField, class GeoMesh>
62 void T
63 (
64  GeometricField<Type, PatchField, GeoMesh>& result,
65  const GeometricField<Type, PatchField, GeoMesh>& f1
66 )
67 {
68  T(result.primitiveFieldRef(), f1.primitiveField());
69  T(result.boundaryFieldRef(), f1.boundaryField());
70  result.oriented() = f1.oriented();
71 }
72 
73 
74 template
75 <
76  class Type,
77  template<class> class PatchField,
78  class GeoMesh,
79  direction r
80 >
81 void pow
82 (
84  <typename powProduct<Type, r>::type, PatchField, GeoMesh>& result,
86 )
87 {
88  pow(result.primitiveFieldRef(), f1.primitiveField(), r);
89  pow(result.boundaryFieldRef(), f1.boundaryField(), r);
90  result.oriented() = pow(f1.oriented(), r);
91 }
92 
93 
94 template
95 <
96  class Type,
97  template<class> class PatchField,
98  class GeoMesh,
99  direction r
100 >
102 pow
103 (
106 )
107 {
108  typedef typename powProduct<Type, r>::type resultType;
109 
110  auto tres =
112  (
113  f1,
114  "pow(" + f1.name() + ',' + Foam::name(r) + ')',
115  pow(f1.dimensions(), r)
116  );
117 
118  pow<Type, r, PatchField, GeoMesh>(tres.ref(), f1);
119 
120  return tres;
121 }
122 
123 
124 template
125 <
126  class Type,
127  template<class> class PatchField,
128  class GeoMesh,
129  direction r
130 >
132 pow
133 (
136 )
137 {
138  typedef typename powProduct<Type, r>::type resultType;
139 
140  const auto& f1 = tf1();
141 
142  auto tres =
144  (
145  tf1,
146  "pow(" + f1.name() + ',' + Foam::name(r) + ')',
147  pow(f1.dimensions(), r)
148  );
149 
150  pow<Type, r, PatchField, GeoMesh>(tres.ref(), f1);
151 
152  tf1.clear();
153  return tres;
154 }
155 
156 
157 template<class Type, template<class> class PatchField, class GeoMesh>
158 void sqr
159 (
160  GeometricField
161  <
162  typename outerProduct<Type, Type>::type, PatchField, GeoMesh
163  >& result,
164  const GeometricField<Type, PatchField, GeoMesh>& f1
165 )
166 {
167  sqr(result.primitiveFieldRef(), f1.primitiveField());
168  sqr(result.boundaryFieldRef(), f1.boundaryField());
169  result.oriented() = sqr(f1.oriented());
170 }
171 
172 
173 template<class Type, template<class> class PatchField, class GeoMesh>
174 tmp
175 <
177  <
179  PatchField,
180  GeoMesh
181  >
182 >
184 {
185  typedef typename outerProduct<Type, Type>::type resultType;
186 
187  auto tres =
189  (
190  f1,
191  "sqr(" + f1.name() + ')',
192  sqr(f1.dimensions())
193  );
194 
195  sqr(tres.ref(), f1);
196 
197  return tres;
198 }
199 
200 
201 template<class Type, template<class> class PatchField, class GeoMesh>
202 tmp
203 <
205  <
207  PatchField,
208  GeoMesh
209  >
210 >
212 {
213  typedef typename outerProduct<Type, Type>::type resultType;
214 
215  const auto& f1 = tf1();
216 
217  auto tres =
219  (
220  tf1,
221  "sqr(" + f1.name() + ')',
222  sqr(f1.dimensions())
223  );
224 
225  sqr(tres.ref(), f1);
226 
227  tf1.clear();
228  return tres;
229 }
230 
231 
232 template<class Type, template<class> class PatchField, class GeoMesh>
233 void magSqr
234 (
235  GeometricField<typename typeOfMag<Type>::type, PatchField, GeoMesh>& result,
236  const GeometricField<Type, PatchField, GeoMesh>& f1
237 )
238 {
239  magSqr(result.primitiveFieldRef(), f1.primitiveField());
240  magSqr(result.boundaryFieldRef(), f1.boundaryField());
241  result.oriented() = magSqr(f1.oriented());
242 }
243 
244 
245 template<class Type, template<class> class PatchField, class GeoMesh>
247 magSqr
248 (
249  const GeometricField<Type, PatchField, GeoMesh>& f1
250 )
251 {
252  typedef typename typeOfMag<Type>::type resultType;
253 
254  auto tres =
256  (
257  f1,
258  "magSqr(" + f1.name() + ')',
259  sqr(f1.dimensions())
260  );
261 
262  magSqr(tres.ref(), f1);
264  return tres;
265 }
266 
267 template<class Type, template<class> class PatchField, class GeoMesh>
269 magSqr
270 (
271  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1
272 )
273 {
274  auto tres = magSqr(tf1.cref());
275  tf1.clear();
277  return tres;
278 }
279 
280 
281 template<class Type, template<class> class PatchField, class GeoMesh>
282 void mag
283 (
284  GeometricField<typename typeOfMag<Type>::type, PatchField, GeoMesh>& result,
285  const GeometricField<Type, PatchField, GeoMesh>& f1
286 )
287 {
288  mag(result.primitiveFieldRef(), f1.primitiveField());
289  mag(result.boundaryFieldRef(), f1.boundaryField());
290  result.oriented() = mag(f1.oriented());
291 }
292 
293 
294 template<class Type, template<class> class PatchField, class GeoMesh>
296 mag
297 (
298  const GeometricField<Type, PatchField, GeoMesh>& f1
299 )
300 {
301  typedef typename typeOfMag<Type>::type resultType;
302 
303  auto tres =
305  (
306  f1,
307  "mag(" + f1.name() + ')',
308  f1.dimensions()
309  );
310 
311  mag(tres.ref(), f1);
313  return tres;
314 }
315 
316 template<class Type, template<class> class PatchField, class GeoMesh>
318 mag
319 (
320  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1
321 )
322 {
323  auto tres = mag(tf1.cref());
324  tf1.clear();
326  return tres;
327 }
328 
329 
330 template<class Type, template<class> class PatchField, class GeoMesh>
331 void cmptAv
332 (
333  GeometricField
334  <
336  PatchField,
337  GeoMesh
338  >& result,
339  const GeometricField<Type, PatchField, GeoMesh>& f1
340 )
341 {
342  cmptAv(result.primitiveFieldRef(), f1.primitiveField());
343  cmptAv(result.boundaryFieldRef(), f1.boundaryField());
344  result.oriented() = cmptAv(f1.oriented());
345 }
346 
347 template<class Type, template<class> class PatchField, class GeoMesh>
348 tmp
349 <
351  <
353  PatchField,
354  GeoMesh
355  >
356 >
358 {
359  typedef typename
361 
362  auto tres =
364  (
365  f1,
366  "cmptAv(" + f1.name() + ')',
367  f1.dimensions()
368  );
369 
370  cmptAv(tres.ref(), f1);
371 
372  return tres;
373 }
374 
375 
376 template<class Type, template<class> class PatchField, class GeoMesh>
377 tmp
378 <
380  <
382  PatchField,
383  GeoMesh
384  >
385 >
387 {
388  auto tres = cmptAv(tf1.cref());
389  tf1.clear();
390 
391  return tres;
392 }
393 
394 
395 #define UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(ReturnType, Func, BinaryOp) \
396  \
397 template<class Type, template<class> class PatchField, class GeoMesh> \
398 dimensioned<ReturnType> Func \
399 ( \
400  const GeometricField<Type, PatchField, GeoMesh>& f1 \
401 ) \
402 { \
403  return dimensioned<ReturnType> \
404  ( \
405  #Func "(" + f1.name() + ')', \
406  f1.dimensions(), \
407  returnReduce \
408  ( \
409  Foam::Func \
410  ( \
411  Foam::Func(f1.primitiveField()), \
412  Foam::Func(f1.boundaryField()) \
413  ), \
414  BinaryOp<ReturnType>() \
415  ) \
416  ); \
417 } \
418  \
419 template<class Type, template<class> class PatchField, class GeoMesh> \
420 dimensioned<ReturnType> Func \
421 ( \
422  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1 \
423 ) \
424 { \
425  dimensioned<ReturnType> res = Func(tf1()); \
426  tf1.clear(); \
427  return res; \
428 }
429 
432 UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(MinMax<Type>, minMax, minMaxOp)
434 
435 #undef UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY
436 
437 
438 #define UNARY_REDUCTION_FUNCTION(ReturnType, Func, gFunc) \
439  \
440 template<class Type, template<class> class PatchField, class GeoMesh> \
441 dimensioned<ReturnType> Func \
442 ( \
443  const GeometricField<Type, PatchField, GeoMesh>& f1 \
444 ) \
445 { \
446  return dimensioned<ReturnType> \
447  ( \
448  #Func "(" + f1.name() + ')', \
449  f1.dimensions(), \
450  gFunc(f1.primitiveField()) \
451  ); \
452 } \
453  \
454 template<class Type, template<class> class PatchField, class GeoMesh> \
455 dimensioned<ReturnType> Func \
456 ( \
457  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1 \
458 ) \
459 { \
460  dimensioned<ReturnType> res = Func(tf1()); \
461  tf1.clear(); \
462  return res; \
463 }
464 
468 
469 #undef UNARY_REDUCTION_FUNCTION
470 
471 
472 BINARY_FUNCTION(Type, Type, Type, max)
473 BINARY_FUNCTION(Type, Type, Type, min)
474 BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
475 BINARY_FUNCTION(Type, Type, Type, cmptDivide)
476 
477 BINARY_TYPE_FUNCTION(Type, Type, Type, max)
478 BINARY_TYPE_FUNCTION(Type, Type, Type, min)
479 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
480 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
481 
483 // ------------------------------------------------------------------------- //
484 
485 // Clamp Methods
486 
487 template<class Type, template<class> class PatchField, class GeoMesh>
488 void clamp
489 (
490  GeometricField<Type, PatchField, GeoMesh>& result,
491  const GeometricField<Type, PatchField, GeoMesh>& f1,
492  const Foam::zero_one
493 )
494 {
496 
497  clamp(result.primitiveFieldRef(), f1.primitiveField(), range);
498  clamp(result.boundaryFieldRef(), f1.boundaryField(), range);
499  result.oriented() = f1.oriented();
500 }
501 
502 template<class Type, template<class> class PatchField, class GeoMesh>
503 tmp<GeometricField<Type, PatchField, GeoMesh>>
504 clamp
505 (
506  const GeometricField<Type, PatchField, GeoMesh>& f1,
507  const Foam::zero_one
508 )
509 {
510  auto tres =
512  (
513  f1,
514  "clamp01(" + f1.name() + ')',
515  f1.dimensions()
516  );
517 
518  clamp(tres.ref(), f1, Foam::zero_one{});
519 
520  return tres;
521 }
522 
523 
524 template<class Type, template<class> class PatchField, class GeoMesh>
525 tmp<GeometricField<Type, PatchField, GeoMesh>>
526 clamp
527 (
528  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1,
529  const Foam::zero_one
530 )
531 {
532  const auto& f1 = tf1();
533 
534  auto tres =
536  (
537  tf1,
538  "clamp01(" + f1.name() + ')',
539  f1.dimensions()
540  );
541 
542  clamp(tres.ref(), f1, Foam::zero_one{});
543 
544  tf1.clear();
545  return tres;
546 }
547 
548 BINARY_TYPE_FUNCTION_FS(Type, Type, MinMax<Type>, clamp)
549 
550 
551 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
552 
553 TERNARY_FUNCTION(Type, Type, Type, scalar, lerp)
554 TERNARY_TYPE_FUNCTION_FFS(Type, Type, Type, scalar, lerp)
555 
556 
557 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
558 
559 UNARY_OPERATOR(Type, Type, -, negate, transform)
560 
561 BINARY_OPERATOR(Type, Type, scalar, *, '*', multiply)
562 BINARY_OPERATOR(Type, scalar, Type, *, '*', multiply)
563 BINARY_OPERATOR(Type, Type, scalar, /, '|', divide)
564 
565 BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, '*', multiply)
566 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, '*', multiply)
567 
568 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, '|', divide)
569 
570 
571 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
572 
573 #define PRODUCT_OPERATOR(product, Op, OpFunc) \
574  \
575 template \
576 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
577 void OpFunc \
578 ( \
579  GeometricField \
580  <typename product<Type1, Type2>::type, PatchField, GeoMesh>& result, \
581  const GeometricField<Type1, PatchField, GeoMesh>& f1, \
582  const GeometricField<Type2, PatchField, GeoMesh>& f2 \
583 ) \
584 { \
585  Foam::OpFunc \
586  ( \
587  result.primitiveFieldRef(), \
588  f1.primitiveField(), \
589  f2.primitiveField() \
590  ); \
591  Foam::OpFunc \
592  ( \
593  result.boundaryFieldRef(), \
594  f1.boundaryField(), \
595  f2.boundaryField() \
596  ); \
597  \
598  result.oriented() = (f1.oriented() Op f2.oriented()); \
599 } \
600  \
601  \
602 template \
603 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
604 tmp \
605 < \
606  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
607 > \
608 operator Op \
609 ( \
610  const GeometricField<Type1, PatchField, GeoMesh>& f1, \
611  const GeometricField<Type2, PatchField, GeoMesh>& f2 \
612 ) \
613 { \
614  typedef typename product<Type1, Type2>::type resultType; \
615  \
616  auto tres = \
617  reuseTmpGeometricField<resultType, Type1, PatchField, GeoMesh>::New \
618  ( \
619  f1, \
620  '(' + f1.name() + #Op + f2.name() + ')', \
621  (f1.dimensions() Op f2.dimensions()) \
622  ); \
623  \
624  Foam::OpFunc(tres.ref(), f1, f2); \
625  \
626  return tres; \
627 } \
628  \
629  \
630 template \
631 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
632 tmp \
633 < \
634  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
635 > \
636 operator Op \
637 ( \
638  const GeometricField<Type1, PatchField, GeoMesh>& f1, \
639  const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tf2 \
640 ) \
641 { \
642  typedef typename product<Type1, Type2>::type resultType; \
643  \
644  const auto& f2 = tf2(); \
645  \
646  auto tres = \
647  reuseTmpGeometricField<resultType, Type2, PatchField, GeoMesh>::New \
648  ( \
649  tf2, \
650  '(' + f1.name() + #Op + f2.name() + ')', \
651  (f1.dimensions() Op f2.dimensions()) \
652  ); \
653  \
654  Foam::OpFunc(tres.ref(), f1, f2); \
655  \
656  tf2.clear(); \
657  return tres; \
658 } \
659  \
660 template \
661 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
662 tmp \
663 < \
664  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
665 > \
666 operator Op \
667 ( \
668  const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tf1, \
669  const GeometricField<Type2, PatchField, GeoMesh>& f2 \
670 ) \
671 { \
672  typedef typename product<Type1, Type2>::type resultType; \
673  \
674  const auto& f1 = tf1(); \
675  \
676  auto tres = \
677  reuseTmpGeometricField<resultType, Type1, PatchField, GeoMesh>::New \
678  ( \
679  tf1, \
680  '(' + f1.name() + #Op + f2.name() + ')', \
681  (f1.dimensions() Op f2.dimensions()) \
682  ); \
683  \
684  Foam::OpFunc(tres.ref(), f1, f2); \
685  \
686  tf1.clear(); \
687  return tres; \
688 } \
689  \
690 template \
691 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
692 tmp \
693 < \
694  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
695 > \
696 operator Op \
697 ( \
698  const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tf1, \
699  const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tf2 \
700 ) \
701 { \
702  typedef typename product<Type1, Type2>::type resultType; \
703  \
704  const auto& f1 = tf1(); \
705  const auto& f2 = tf2(); \
706  \
707  auto tres = \
708  reuseTmpTmpGeometricField \
709  <resultType, Type1, Type1, Type2, PatchField, GeoMesh>::New \
710  ( \
711  tf1, \
712  tf2, \
713  '(' + f1.name() + #Op + f2.name() + ')', \
714  (f1.dimensions() Op f2.dimensions()) \
715  ); \
716  \
717  Foam::OpFunc(tres.ref(), f1, f2); \
718  \
719  tf1.clear(); \
720  tf2.clear(); \
721  return tres; \
722 } \
723  \
724 template \
725 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
726 void OpFunc \
727 ( \
728  GeometricField \
729  <typename product<Type, Form>::type, PatchField, GeoMesh>& result, \
730  const GeometricField<Type, PatchField, GeoMesh>& f1, \
731  const dimensioned<Form>& dvs \
732 ) \
733 { \
734  Foam::OpFunc(result.primitiveFieldRef(), f1.primitiveField(), dvs.value());\
735  Foam::OpFunc(result.boundaryFieldRef(), f1.boundaryField(), dvs.value()); \
736  result.oriented() = f1.oriented(); \
737 } \
738  \
739 template \
740 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
741 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
742 operator Op \
743 ( \
744  const GeometricField<Type, PatchField, GeoMesh>& f1, \
745  const dimensioned<Form>& dvs \
746 ) \
747 { \
748  typedef typename product<Type, Form>::type resultType; \
749  \
750  auto tres = \
751  reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
752  ( \
753  f1, \
754  '(' + f1.name() + #Op + dvs.name() + ')', \
755  (f1.dimensions() Op dvs.dimensions()) \
756  ); \
757  \
758  Foam::OpFunc(tres.ref(), f1, dvs); \
759  \
760  return tres; \
761 } \
762  \
763 template \
764 < \
765  class Form, \
766  class Cmpt, \
767  direction nCmpt, \
768  class Type, template<class> class PatchField, \
769  class GeoMesh \
770 > \
771 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
772 operator Op \
773 ( \
774  const GeometricField<Type, PatchField, GeoMesh>& f1, \
775  const VectorSpace<Form,Cmpt,nCmpt>& vs \
776 ) \
777 { \
778  return f1 Op dimensioned<Form>(static_cast<const Form&>(vs)); \
779 } \
780  \
781  \
782 template \
783 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
784 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
785 operator Op \
786 ( \
787  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1, \
788  const dimensioned<Form>& dvs \
789 ) \
790 { \
791  typedef typename product<Type, Form>::type resultType; \
792  \
793  const auto& f1 = tf1(); \
794  \
795  auto tres = \
796  reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
797  ( \
798  tf1, \
799  '(' + f1.name() + #Op + dvs.name() + ')', \
800  (f1.dimensions() Op dvs.dimensions()) \
801  ); \
802  \
803  Foam::OpFunc(tres.ref(), f1, dvs); \
804  \
805  tf1.clear(); \
806  return tres; \
807 } \
808  \
809 template \
810 < \
811  class Form, \
812  class Cmpt, \
813  direction nCmpt, \
814  class Type, template<class> class PatchField, \
815  class GeoMesh \
816 > \
817 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
818 operator Op \
819 ( \
820  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1, \
821  const VectorSpace<Form,Cmpt,nCmpt>& vs \
822 ) \
823 { \
824  return tf1 Op dimensioned<Form>(static_cast<const Form&>(vs)); \
825 } \
826  \
827  \
828 template \
829 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
830 void OpFunc \
831 ( \
832  GeometricField \
833  <typename product<Form, Type>::type, PatchField, GeoMesh>& result, \
834  const dimensioned<Form>& dvs, \
835  const GeometricField<Type, PatchField, GeoMesh>& f2 \
836 ) \
837 { \
838  Foam::OpFunc(result.primitiveFieldRef(), dvs.value(), f2.primitiveField());\
839  Foam::OpFunc(result.boundaryFieldRef(), dvs.value(), f2.boundaryField()); \
840  result.oriented() = f2.oriented(); \
841 } \
842  \
843 template \
844 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
845 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
846 operator Op \
847 ( \
848  const dimensioned<Form>& dvs, \
849  const GeometricField<Type, PatchField, GeoMesh>& f2 \
850 ) \
851 { \
852  typedef typename product<Form, Type>::type resultType; \
853  \
854  auto tres = \
855  reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
856  ( \
857  f2, \
858  '(' + dvs.name() + #Op + f2.name() + ')', \
859  (dvs.dimensions() Op f2.dimensions()) \
860  ); \
861  \
862  Foam::OpFunc(tres.ref(), dvs, f2); \
863  \
864  return tres; \
865 } \
866  \
867 template \
868 < \
869  class Form, \
870  class Cmpt, \
871  direction nCmpt, \
872  class Type, template<class> class PatchField, \
873  class GeoMesh \
874 > \
875 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
876 operator Op \
877 ( \
878  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
879  const GeometricField<Type, PatchField, GeoMesh>& f2 \
880 ) \
881 { \
882  return dimensioned<Form>(static_cast<const Form&>(vs)) Op f2; \
883 } \
884  \
885 template \
886 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
887 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
888 operator Op \
889 ( \
890  const dimensioned<Form>& dvs, \
891  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf2 \
892 ) \
893 { \
894  typedef typename product<Form, Type>::type resultType; \
895  \
896  const auto& f2 = tf2(); \
897  \
898  auto tres = \
899  reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
900  ( \
901  tf2, \
902  '(' + dvs.name() + #Op + f2.name() + ')', \
903  (dvs.dimensions() Op f2.dimensions()) \
904  ); \
905  \
906  Foam::OpFunc(tres.ref(), dvs, f2); \
907  \
908  tf2.clear(); \
909  return tres; \
910 } \
911  \
912 template \
913 < \
914  class Form, \
915  class Cmpt, \
916  direction nCmpt, \
917  class Type, template<class> class PatchField, \
918  class GeoMesh \
919 > \
920 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
921 operator Op \
922 ( \
923  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
924  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf2 \
925 ) \
926 { \
927  return dimensioned<Form>(static_cast<const Form&>(vs)) Op tf2; \
928 }
933 
938 
939 #undef PRODUCT_OPERATOR
940 
941 
942 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
943 
944 } // End namespace Foam
945 
946 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
947 
948 #include "undefFieldFunctionsM.H"
949 
950 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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)
uint8_t direction
Definition: direction.H:48
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
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)
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:162
Field< Type >::cmptType cmptType
Component type of the field elements.
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)
orientedType oriented() const noexcept
Return oriented type.
#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
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
#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)
scalar range
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
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 INVALID.
Definition: exprTraits.C:52
static tmp< GeometricField< TypeR, PatchField, GeoMesh > > New(const GeometricField< Type1, PatchField, GeoMesh > &f1, const word &name, const dimensionSet &dimensions)
Pass-through to New GeometricField.
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 > &)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
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)
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition: pTraits.H:84
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
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
dimensioned< scalarMinMax > minMaxMag(const DimensionedField< Type, GeoMesh > &f1)
Type gAverage(const FieldField< Field, Type > &f)
#define PRODUCT_OPERATOR(product, Op, OpFunc)
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
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)
#define UNARY_REDUCTION_FUNCTION(ReturnType, Func, gFunc)
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)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
#define UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(ReturnType, Func, BinaryOp)
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