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);
56  component(result.boundaryFieldRef(), f1.boundaryField(), d);
57  result.oriented() = f1.oriented();
59  {
60  result.boundaryField().check();
61  }
62 }
63 
64 
65 template<class Type, template<class> class PatchField, class GeoMesh>
66 void T
67 (
68  GeometricField<Type, PatchField, GeoMesh>& result,
69  const GeometricField<Type, PatchField, GeoMesh>& f1
70 )
71 {
72  T(result.primitiveFieldRef(), f1.primitiveField());
73  T(result.boundaryFieldRef(), f1.boundaryField());
74  result.oriented() = f1.oriented();
76  {
77  result.boundaryField().check();
78  }
79 }
80 
81 
82 template
83 <
84  class Type,
85  template<class> class PatchField,
86  class GeoMesh,
87  direction r
88 >
89 void pow
90 (
92  <typename powProduct<Type, r>::type, PatchField, GeoMesh>& result,
94 )
95 {
96  pow(result.primitiveFieldRef(), f1.primitiveField(), r);
97  pow(result.boundaryFieldRef(), f1.boundaryField(), r);
98  result.oriented() = pow(f1.oriented(), r);
101  {
102  result.boundaryField().check();
103  }
104 }
105 
106 
107 template
108 <
109  class Type,
110  template<class> class PatchField,
111  class GeoMesh,
112  direction r
113 >
115 pow
116 (
119 )
120 {
121  typedef typename powProduct<Type, r>::type resultType;
122 
123  auto tres =
125  (
126  f1,
127  "pow(" + f1.name() + ',' + Foam::name(r) + ')',
128  pow(f1.dimensions(), r)
129  );
130 
131  pow<Type, r, PatchField, GeoMesh>(tres.ref(), f1);
132 
133  return tres;
134 }
135 
136 
137 template
138 <
139  class Type,
140  template<class> class PatchField,
141  class GeoMesh,
142  direction r
143 >
145 pow
146 (
149 )
150 {
151  typedef typename powProduct<Type, r>::type resultType;
152 
153  const auto& f1 = tf1();
154 
155  auto tres =
157  (
158  tf1,
159  "pow(" + f1.name() + ',' + Foam::name(r) + ')',
160  pow(f1.dimensions(), r)
161  );
162 
163  pow<Type, r, PatchField, GeoMesh>(tres.ref(), f1);
164 
165  tf1.clear();
166  return tres;
167 }
168 
169 
170 template<class Type, template<class> class PatchField, class GeoMesh>
171 void sqr
172 (
173  GeometricField
174  <
175  typename outerProduct<Type, Type>::type, PatchField, GeoMesh
176  >& result,
177  const GeometricField<Type, PatchField, GeoMesh>& f1
178 )
179 {
180  sqr(result.primitiveFieldRef(), f1.primitiveField());
181  sqr(result.boundaryFieldRef(), f1.boundaryField());
182  result.oriented() = sqr(f1.oriented());
185  {
186  result.boundaryField().check();
187  }
188 }
189 
190 
191 template<class Type, template<class> class PatchField, class GeoMesh>
192 tmp
193 <
195  <
197  PatchField,
198  GeoMesh
199  >
200 >
202 {
203  typedef typename outerProduct<Type, Type>::type resultType;
204 
205  auto tres =
207  (
208  f1,
209  "sqr(" + f1.name() + ')',
210  sqr(f1.dimensions())
211  );
212 
213  sqr(tres.ref(), f1);
214 
215  return tres;
216 }
217 
218 
219 template<class Type, template<class> class PatchField, class GeoMesh>
220 tmp
221 <
223  <
225  PatchField,
226  GeoMesh
227  >
228 >
230 {
231  typedef typename outerProduct<Type, Type>::type resultType;
232 
233  const auto& f1 = tf1();
234 
235  auto tres =
237  (
238  tf1,
239  "sqr(" + f1.name() + ')',
240  sqr(f1.dimensions())
241  );
242 
243  sqr(tres.ref(), f1);
244 
245  tf1.clear();
246  return tres;
247 }
248 
249 
250 template<class Type, template<class> class PatchField, class GeoMesh>
251 void magSqr
252 (
253  GeometricField<typename typeOfMag<Type>::type, PatchField, GeoMesh>& result,
254  const GeometricField<Type, PatchField, GeoMesh>& f1
255 )
256 {
257  magSqr(result.primitiveFieldRef(), f1.primitiveField());
258  magSqr(result.boundaryFieldRef(), f1.boundaryField());
259  result.oriented() = magSqr(f1.oriented());
262  {
263  result.boundaryField().check();
264  }
265 }
266 
267 
268 template<class Type, template<class> class PatchField, class GeoMesh>
270 magSqr
271 (
272  const GeometricField<Type, PatchField, GeoMesh>& f1
273 )
274 {
275  typedef typename typeOfMag<Type>::type resultType;
276 
277  auto tres =
279  (
280  f1,
281  "magSqr(" + f1.name() + ')',
282  sqr(f1.dimensions())
283  );
284 
285  magSqr(tres.ref(), f1);
287  return tres;
288 }
289 
290 template<class Type, template<class> class PatchField, class GeoMesh>
292 magSqr
293 (
294  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1
295 )
296 {
297  auto tres = magSqr(tf1.cref());
298  tf1.clear();
300  return tres;
301 }
302 
303 
304 template<class Type, template<class> class PatchField, class GeoMesh>
305 void mag
306 (
307  GeometricField<typename typeOfMag<Type>::type, PatchField, GeoMesh>& result,
308  const GeometricField<Type, PatchField, GeoMesh>& f1
309 )
310 {
311  mag(result.primitiveFieldRef(), f1.primitiveField());
312  mag(result.boundaryFieldRef(), f1.boundaryField());
313  result.oriented() = mag(f1.oriented());
316  {
317  result.boundaryField().check();
318  }
319 }
320 
321 
322 template<class Type, template<class> class PatchField, class GeoMesh>
324 mag
325 (
326  const GeometricField<Type, PatchField, GeoMesh>& f1
327 )
328 {
329  typedef typename typeOfMag<Type>::type resultType;
330 
331  auto tres =
333  (
334  f1,
335  "mag(" + f1.name() + ')',
336  f1.dimensions()
337  );
338 
339  mag(tres.ref(), f1);
341  return tres;
342 }
343 
344 template<class Type, template<class> class PatchField, class GeoMesh>
346 mag
347 (
348  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1
349 )
350 {
351  auto tres = mag(tf1.cref());
352  tf1.clear();
354  return tres;
355 }
356 
357 
358 template<class Type, template<class> class PatchField, class GeoMesh>
359 void cmptAv
360 (
361  GeometricField
362  <
364  PatchField,
365  GeoMesh
366  >& result,
367  const GeometricField<Type, PatchField, GeoMesh>& f1
368 )
369 {
370  cmptAv(result.primitiveFieldRef(), f1.primitiveField());
371  cmptAv(result.boundaryFieldRef(), f1.boundaryField());
372  result.oriented() = cmptAv(f1.oriented());
374  {
375  result.boundaryField().check();
376  }
377 }
378 
379 template<class Type, template<class> class PatchField, class GeoMesh>
380 tmp
381 <
383  <
385  PatchField,
386  GeoMesh
387  >
388 >
390 {
391  typedef typename
393 
394  auto tres =
396  (
397  f1,
398  "cmptAv(" + f1.name() + ')',
399  f1.dimensions()
400  );
401 
402  cmptAv(tres.ref(), f1);
403 
404  return tres;
405 }
406 
407 
408 template<class Type, template<class> class PatchField, class GeoMesh>
409 tmp
410 <
412  <
414  PatchField,
415  GeoMesh
416  >
417 >
419 {
420  auto tres = cmptAv(tf1.cref());
421  tf1.clear();
422 
423  return tres;
424 }
425 
426 
427 #define UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(ReturnType, Func, BinaryOp) \
428  \
429 template<class Type, template<class> class PatchField, class GeoMesh> \
430 dimensioned<ReturnType> Func \
431 ( \
432  const GeometricField<Type, PatchField, GeoMesh>& f1 \
433 ) \
434 { \
435  return dimensioned<ReturnType> \
436  ( \
437  #Func "(" + f1.name() + ')', \
438  f1.dimensions(), \
439  returnReduce \
440  ( \
441  Foam::Func \
442  ( \
443  Foam::Func(f1.primitiveField()), \
444  Foam::Func(f1.boundaryField()) \
445  ), \
446  BinaryOp<ReturnType>() \
447  ) \
448  ); \
449 } \
450  \
451 template<class Type, template<class> class PatchField, class GeoMesh> \
452 dimensioned<ReturnType> Func \
453 ( \
454  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1 \
455 ) \
456 { \
457  dimensioned<ReturnType> res = Func(tf1()); \
458  tf1.clear(); \
459  return res; \
460 }
461 
464 UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(MinMax<Type>, minMax, minMaxOp)
466 
467 #undef UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY
468 
469 
470 #define UNARY_REDUCTION_FUNCTION(ReturnType, Func, gFunc) \
471  \
472 template<class Type, template<class> class PatchField, class GeoMesh> \
473 dimensioned<ReturnType> Func \
474 ( \
475  const GeometricField<Type, PatchField, GeoMesh>& f1 \
476 ) \
477 { \
478  return dimensioned<ReturnType> \
479  ( \
480  #Func "(" + f1.name() + ')', \
481  f1.dimensions(), \
482  gFunc(f1.primitiveField()) \
483  ); \
484 } \
485  \
486 template<class Type, template<class> class PatchField, class GeoMesh> \
487 dimensioned<ReturnType> Func \
488 ( \
489  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1 \
490 ) \
491 { \
492  dimensioned<ReturnType> res = Func(tf1()); \
493  tf1.clear(); \
494  return res; \
495 }
496 
500 
501 #undef UNARY_REDUCTION_FUNCTION
502 
503 
504 BINARY_FUNCTION(Type, Type, Type, max)
505 BINARY_FUNCTION(Type, Type, Type, min)
506 BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
507 BINARY_FUNCTION(Type, Type, Type, cmptDivide)
508 
509 BINARY_TYPE_FUNCTION(Type, Type, Type, max)
510 BINARY_TYPE_FUNCTION(Type, Type, Type, min)
511 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
512 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
513 
515 // ------------------------------------------------------------------------- //
516 
517 // Clamp Methods
518 
519 template<class Type, template<class> class PatchField, class GeoMesh>
520 void clamp
521 (
522  GeometricField<Type, PatchField, GeoMesh>& result,
523  const GeometricField<Type, PatchField, GeoMesh>& f1,
524  const Foam::zero_one
525 )
526 {
528 
529  clamp(result.primitiveFieldRef(), f1.primitiveField(), range);
530  clamp(result.boundaryFieldRef(), f1.boundaryField(), range);
531  result.oriented() = f1.oriented();
532  result.correctLocalBoundaryConditions();
534  {
535  result.boundaryField().check();
536  }
537 }
538 
539 template<class Type, template<class> class PatchField, class GeoMesh>
540 tmp<GeometricField<Type, PatchField, GeoMesh>>
541 clamp
542 (
543  const GeometricField<Type, PatchField, GeoMesh>& f1,
544  const Foam::zero_one
545 )
546 {
547  auto tres =
549  (
550  f1,
551  "clamp01(" + f1.name() + ')',
552  f1.dimensions()
553  );
554 
555  clamp(tres.ref(), f1, Foam::zero_one{});
556 
557  return tres;
558 }
559 
560 
561 template<class Type, template<class> class PatchField, class GeoMesh>
562 tmp<GeometricField<Type, PatchField, GeoMesh>>
563 clamp
564 (
565  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1,
566  const Foam::zero_one
567 )
568 {
569  const auto& f1 = tf1();
570 
571  auto tres =
573  (
574  tf1,
575  "clamp01(" + f1.name() + ')',
576  f1.dimensions()
577  );
578 
579  clamp(tres.ref(), f1, Foam::zero_one{});
580 
581  tf1.clear();
582  return tres;
583 }
584 
585 BINARY_TYPE_FUNCTION_FS(Type, Type, MinMax<Type>, clamp)
586 
587 
588 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
589 
590 TERNARY_FUNCTION(Type, Type, Type, scalar, lerp)
591 TERNARY_TYPE_FUNCTION_FFS(Type, Type, Type, scalar, lerp)
592 
593 
594 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
595 
596 UNARY_OPERATOR(Type, Type, -, negate, transform)
597 
598 BINARY_OPERATOR(Type, Type, scalar, *, '*', multiply)
599 BINARY_OPERATOR(Type, scalar, Type, *, '*', multiply)
600 BINARY_OPERATOR(Type, Type, scalar, /, '|', divide)
601 
602 BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, '*', multiply)
603 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, '*', multiply)
604 
605 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, '|', divide)
606 
607 
608 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
609 
610 #define PRODUCT_OPERATOR(product, Op, OpFunc) \
611  \
612 template \
613 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
614 void OpFunc \
615 ( \
616  GeometricField \
617  <typename product<Type1, Type2>::type, PatchField, GeoMesh>& result, \
618  const GeometricField<Type1, PatchField, GeoMesh>& f1, \
619  const GeometricField<Type2, PatchField, GeoMesh>& f2 \
620 ) \
621 { \
622  Foam::OpFunc \
623  ( \
624  result.primitiveFieldRef(), \
625  f1.primitiveField(), \
626  f2.primitiveField() \
627  ); \
628  Foam::OpFunc \
629  ( \
630  result.boundaryFieldRef(), \
631  f1.boundaryField(), \
632  f2.boundaryField() \
633  ); \
634  \
635  result.oriented() = (f1.oriented() Op f2.oriented()); \
636  result.correctLocalBoundaryConditions(); \
637  if (GeometricBoundaryField<Type1, PatchField, GeoMesh>::debug) \
638  { \
639  result.boundaryField().check(); \
640  } \
641 } \
642  \
643  \
644 template \
645 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
646 tmp \
647 < \
648  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
649 > \
650 operator Op \
651 ( \
652  const GeometricField<Type1, PatchField, GeoMesh>& f1, \
653  const GeometricField<Type2, PatchField, GeoMesh>& f2 \
654 ) \
655 { \
656  typedef typename product<Type1, Type2>::type resultType; \
657  \
658  auto tres = \
659  reuseTmpGeometricField<resultType, Type1, PatchField, GeoMesh>::New \
660  ( \
661  f1, \
662  '(' + f1.name() + #Op + f2.name() + ')', \
663  (f1.dimensions() Op f2.dimensions()) \
664  ); \
665  \
666  Foam::OpFunc(tres.ref(), f1, f2); \
667  \
668  return tres; \
669 } \
670  \
671  \
672 template \
673 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
674 tmp \
675 < \
676  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
677 > \
678 operator Op \
679 ( \
680  const GeometricField<Type1, PatchField, GeoMesh>& f1, \
681  const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tf2 \
682 ) \
683 { \
684  typedef typename product<Type1, Type2>::type resultType; \
685  \
686  const auto& f2 = tf2(); \
687  \
688  auto tres = \
689  reuseTmpGeometricField<resultType, Type2, PatchField, GeoMesh>::New \
690  ( \
691  tf2, \
692  '(' + f1.name() + #Op + f2.name() + ')', \
693  (f1.dimensions() Op f2.dimensions()) \
694  ); \
695  \
696  Foam::OpFunc(tres.ref(), f1, f2); \
697  \
698  tf2.clear(); \
699  return tres; \
700 } \
701  \
702 template \
703 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
704 tmp \
705 < \
706  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
707 > \
708 operator Op \
709 ( \
710  const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tf1, \
711  const GeometricField<Type2, PatchField, GeoMesh>& f2 \
712 ) \
713 { \
714  typedef typename product<Type1, Type2>::type resultType; \
715  \
716  const auto& f1 = tf1(); \
717  \
718  auto tres = \
719  reuseTmpGeometricField<resultType, Type1, PatchField, GeoMesh>::New \
720  ( \
721  tf1, \
722  '(' + f1.name() + #Op + f2.name() + ')', \
723  (f1.dimensions() Op f2.dimensions()) \
724  ); \
725  \
726  Foam::OpFunc(tres.ref(), f1, f2); \
727  \
728  tf1.clear(); \
729  return tres; \
730 } \
731  \
732 template \
733 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
734 tmp \
735 < \
736  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
737 > \
738 operator Op \
739 ( \
740  const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tf1, \
741  const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tf2 \
742 ) \
743 { \
744  typedef typename product<Type1, Type2>::type resultType; \
745  \
746  const auto& f1 = tf1(); \
747  const auto& f2 = tf2(); \
748  \
749  auto tres = \
750  reuseTmpTmpGeometricField \
751  <resultType, Type1, Type1, Type2, PatchField, GeoMesh>::New \
752  ( \
753  tf1, \
754  tf2, \
755  '(' + f1.name() + #Op + f2.name() + ')', \
756  (f1.dimensions() Op f2.dimensions()) \
757  ); \
758  \
759  Foam::OpFunc(tres.ref(), f1, f2); \
760  \
761  tf1.clear(); \
762  tf2.clear(); \
763  return tres; \
764 } \
765  \
766 template \
767 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
768 void OpFunc \
769 ( \
770  GeometricField \
771  <typename product<Type, Form>::type, PatchField, GeoMesh>& result, \
772  const GeometricField<Type, PatchField, GeoMesh>& f1, \
773  const dimensioned<Form>& dvs \
774 ) \
775 { \
776  Foam::OpFunc(result.primitiveFieldRef(), f1.primitiveField(), dvs.value());\
777  Foam::OpFunc(result.boundaryFieldRef(), f1.boundaryField(), dvs.value()); \
778  result.oriented() = f1.oriented(); \
779  if (GeometricBoundaryField<Type, PatchField, GeoMesh>::debug) \
780  { \
781  result.boundaryField().check(); \
782  } \
783 } \
784  \
785 template \
786 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
787 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
788 operator Op \
789 ( \
790  const GeometricField<Type, PatchField, GeoMesh>& f1, \
791  const dimensioned<Form>& dvs \
792 ) \
793 { \
794  typedef typename product<Type, Form>::type resultType; \
795  \
796  auto tres = \
797  reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
798  ( \
799  f1, \
800  '(' + f1.name() + #Op + dvs.name() + ')', \
801  (f1.dimensions() Op dvs.dimensions()) \
802  ); \
803  \
804  Foam::OpFunc(tres.ref(), f1, dvs); \
805  \
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 GeometricField<Type, PatchField, GeoMesh>& f1, \
821  const VectorSpace<Form,Cmpt,nCmpt>& vs \
822 ) \
823 { \
824  return f1 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 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
831 operator Op \
832 ( \
833  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1, \
834  const dimensioned<Form>& dvs \
835 ) \
836 { \
837  typedef typename product<Type, Form>::type resultType; \
838  \
839  const auto& f1 = tf1(); \
840  \
841  auto tres = \
842  reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
843  ( \
844  tf1, \
845  '(' + f1.name() + #Op + dvs.name() + ')', \
846  (f1.dimensions() Op dvs.dimensions()) \
847  ); \
848  \
849  Foam::OpFunc(tres.ref(), f1, dvs); \
850  \
851  tf1.clear(); \
852  return tres; \
853 } \
854  \
855 template \
856 < \
857  class Form, \
858  class Cmpt, \
859  direction nCmpt, \
860  class Type, template<class> class PatchField, \
861  class GeoMesh \
862 > \
863 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
864 operator Op \
865 ( \
866  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1, \
867  const VectorSpace<Form,Cmpt,nCmpt>& vs \
868 ) \
869 { \
870  return tf1 Op dimensioned<Form>(static_cast<const Form&>(vs)); \
871 } \
872  \
873  \
874 template \
875 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
876 void OpFunc \
877 ( \
878  GeometricField \
879  <typename product<Form, Type>::type, PatchField, GeoMesh>& result, \
880  const dimensioned<Form>& dvs, \
881  const GeometricField<Type, PatchField, GeoMesh>& f2 \
882 ) \
883 { \
884  Foam::OpFunc(result.primitiveFieldRef(), dvs.value(), f2.primitiveField());\
885  Foam::OpFunc(result.boundaryFieldRef(), dvs.value(), f2.boundaryField()); \
886  result.oriented() = f2.oriented(); \
887  if (GeometricBoundaryField<Type, PatchField, GeoMesh>::debug) \
888  { \
889  result.boundaryField().check(); \
890  } \
891 } \
892  \
893 template \
894 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
895 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
896 operator Op \
897 ( \
898  const dimensioned<Form>& dvs, \
899  const GeometricField<Type, PatchField, GeoMesh>& f2 \
900 ) \
901 { \
902  typedef typename product<Form, Type>::type resultType; \
903  \
904  auto tres = \
905  reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
906  ( \
907  f2, \
908  '(' + dvs.name() + #Op + f2.name() + ')', \
909  (dvs.dimensions() Op f2.dimensions()) \
910  ); \
911  \
912  Foam::OpFunc(tres.ref(), dvs, f2); \
913  \
914  return tres; \
915 } \
916  \
917 template \
918 < \
919  class Form, \
920  class Cmpt, \
921  direction nCmpt, \
922  class Type, template<class> class PatchField, \
923  class GeoMesh \
924 > \
925 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
926 operator Op \
927 ( \
928  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
929  const GeometricField<Type, PatchField, GeoMesh>& f2 \
930 ) \
931 { \
932  return dimensioned<Form>(static_cast<const Form&>(vs)) Op f2; \
933 } \
934  \
935 template \
936 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
937 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
938 operator Op \
939 ( \
940  const dimensioned<Form>& dvs, \
941  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf2 \
942 ) \
943 { \
944  typedef typename product<Form, Type>::type resultType; \
945  \
946  const auto& f2 = tf2(); \
947  \
948  auto tres = \
949  reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
950  ( \
951  tf2, \
952  '(' + dvs.name() + #Op + f2.name() + ')', \
953  (dvs.dimensions() Op f2.dimensions()) \
954  ); \
955  \
956  Foam::OpFunc(tres.ref(), dvs, f2); \
957  \
958  tf2.clear(); \
959  return tres; \
960 } \
961  \
962 template \
963 < \
964  class Form, \
965  class Cmpt, \
966  direction nCmpt, \
967  class Type, template<class> class PatchField, \
968  class GeoMesh \
969 > \
970 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
971 operator Op \
972 ( \
973  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
974  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf2 \
975 ) \
976 { \
977  return dimensioned<Form>(static_cast<const Form&>(vs)) Op tf2; \
978 }
983 
988 
989 #undef PRODUCT_OPERATOR
990 
991 
992 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
993 
994 } // End namespace Foam
995 
996 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
997 
998 #include "undefFieldFunctionsM.H"
999 
1000 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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:46
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:195
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 expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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)
Generic GeometricBoundaryField class.
Definition: areaFieldsFwd.H:46
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition: pTraits.H:65
symmTypeOfRank< typename pTraits< arg1 >::cmptType, arg2 *direction(pTraits< arg1 >::rank) >::type type
Definition: products.H:176
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void clear()
Clear exponents - resets to be dimensionless.
Definition: dimensionSet.C:136
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)
static int debug
Enable debug.
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)
void correctLocalBoundaryConditions()
Correct boundary conditions after a purely local operation. Is.
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