GeometricFieldFunctions.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 
28 #define TEMPLATE \
29  template<class Type, template<class> class PatchField, class GeoMesh>
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * //
38 
39 template<class Type, template<class> class PatchField, class GeoMesh>
40 void component
41 (
43  <
44  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
45  PatchField,
46  GeoMesh
47  >& gcf,
49  const direction d
50 )
51 {
52  component(gcf.primitiveFieldRef(), gf.primitiveField(), d);
53  component(gcf.boundaryFieldRef(), gf.boundaryField(), d);
54 }
55 
56 
57 template<class Type, template<class> class PatchField, class GeoMesh>
58 void T
59 (
62 )
63 {
64  T(gf.primitiveFieldRef(), gf1.primitiveField());
65  T(gf.boundaryFieldRef(), gf1.boundaryField());
66 }
67 
68 
69 template
70 <
71  class Type,
72  template<class> class PatchField,
73  class GeoMesh,
74  direction r
75 >
76 void pow
77 (
78  GeometricField<typename powProduct<Type, r>::type, PatchField, GeoMesh>& gf,
80 )
81 {
82  pow(gf.primitiveFieldRef(), gf1.primitiveField(), r);
83  pow(gf.boundaryFieldRef(), gf1.boundaryField(), r);
84 }
85 
86 template
87 <
88  class Type,
89  template<class> class PatchField,
90  class GeoMesh,
91  direction r
92 >
93 tmp<GeometricField<typename powProduct<Type, r>::type, PatchField, GeoMesh>>
94 pow
95 (
98 )
99 {
100  typedef typename powProduct<Type, r>::type powProductType;
101 
103  (
105  (
106  "pow(" + gf.name() + ',' + name(r) + ')',
107  gf.mesh(),
108  pow(gf.dimensions(), r)
109  )
110  );
111 
112  pow<Type, r, PatchField, GeoMesh>(tPow.ref(), gf);
113 
114  return tPow;
115 }
116 
117 
118 template
119 <
120  class Type,
121  template<class> class PatchField,
122  class GeoMesh,
123  direction r
124 >
125 tmp<GeometricField<typename powProduct<Type, r>::type, PatchField, GeoMesh>>
126 pow
127 (
130 )
131 {
132  typedef typename powProduct<Type, r>::type powProductType;
133 
135 
137  (
139  (
140  "pow(" + gf.name() + ',' + name(r) + ')',
141  gf.mesh(),
142  pow(gf.dimensions(), r)
143  )
144  );
145 
146  pow<Type, r, PatchField, GeoMesh>(tPow.ref(), gf);
147 
148  tgf.clear();
149 
150  return tPow;
151 }
152 
153 
154 template<class Type, template<class> class PatchField, class GeoMesh>
155 void sqr
156 (
158  <typename outerProduct<Type, Type>::type, PatchField, GeoMesh>& gf,
160 )
161 {
162  sqr(gf.primitiveFieldRef(), gf1.primitiveField());
163  sqr(gf.boundaryFieldRef(), gf1.boundaryField());
164 }
165 
166 template<class Type, template<class> class PatchField, class GeoMesh>
167 tmp
168 <
170  <
171  typename outerProduct<Type, Type>::type,
172  PatchField,
173  GeoMesh
174  >
175 >
177 {
178  typedef typename outerProduct<Type, Type>::type outerProductType;
179 
181  (
183  (
184  "sqr(" + gf.name() + ')',
185  gf.mesh(),
186  sqr(gf.dimensions())
187  )
188  );
189 
190  sqr(tSqr.ref(), gf);
191 
192  return tSqr;
193 }
194 
195 template<class Type, template<class> class PatchField, class GeoMesh>
196 tmp
197 <
199  <
200  typename outerProduct<Type, Type>::type,
201  PatchField,
202  GeoMesh
203  >
204 >
206 {
207  typedef typename outerProduct<Type, Type>::type outerProductType;
208 
210 
212  (
214  (
215  "sqr(" + gf.name() + ')',
216  gf.mesh(),
217  sqr(gf.dimensions())
218  )
219  );
220 
221  sqr(tSqr.ref(), gf);
222 
223  tgf.clear();
224 
225  return tSqr;
226 }
227 
228 
229 template<class Type, template<class> class PatchField, class GeoMesh>
230 void magSqr
231 (
234 )
235 {
238 }
239 
240 template<class Type, template<class> class PatchField, class GeoMesh>
242 (
244 )
245 {
247  (
249  (
250  "magSqr(" + gf.name() + ')',
251  gf.mesh(),
252  sqr(gf.dimensions())
253  )
254  );
255 
256  magSqr(tMagSqr.ref(), gf);
257 
258  return tMagSqr;
259 }
260 
261 template<class Type, template<class> class PatchField, class GeoMesh>
263 (
265 )
266 {
268 
270  (
272  (
273  "magSqr(" + gf.name() + ')',
274  gf.mesh(),
275  sqr(gf.dimensions())
276  )
277  );
278 
279  magSqr(tMagSqr.ref(), gf);
280 
281  tgf.clear();
282 
283  return tMagSqr;
284 }
285 
286 
287 template<class Type, template<class> class PatchField, class GeoMesh>
288 void mag
289 (
292 )
293 {
294  mag(gsf.primitiveFieldRef(), gf.primitiveField());
295  mag(gsf.boundaryFieldRef(), gf.boundaryField());
296 }
297 
298 template<class Type, template<class> class PatchField, class GeoMesh>
300 (
302 )
303 {
305  (
307  (
308  "mag(" + gf.name() + ')',
309  gf.mesh(),
310  gf.dimensions()
311  )
312  );
313 
314  mag(tMag.ref(), gf);
315 
316  return tMag;
317 }
318 
319 template<class Type, template<class> class PatchField, class GeoMesh>
321 (
323 )
324 {
326 
328  (
330  (
331  "mag(" + gf.name() + ')',
332  gf.mesh(),
333  gf.dimensions()
334  )
335  );
336 
337  mag(tMag.ref(), gf);
338 
339  tgf.clear();
340 
341  return tMag;
342 }
343 
344 
345 template<class Type, template<class> class PatchField, class GeoMesh>
346 void cmptAv
347 (
349  <
350  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
351  PatchField,
352  GeoMesh
353  >& gcf,
355 )
356 {
357  cmptAv(gcf.primitiveFieldRef(), gf.primitiveField());
358  cmptAv(gcf.boundaryFieldRef(), gf.boundaryField());
359 }
360 
361 template<class Type, template<class> class PatchField, class GeoMesh>
362 tmp
363 <
365  <
366  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
367  PatchField,
368  GeoMesh
369  >
370 >
372 {
373  typedef typename GeometricField<Type, PatchField, GeoMesh>::cmptType
374  cmptType;
375 
377  (
379  (
380  "cmptAv(" + gf.name() + ')',
381  gf.mesh(),
382  gf.dimensions()
383  )
384  );
385 
386  cmptAv(CmptAv.ref(), gf);
387 
388  return CmptAv;
389 }
390 
391 template<class Type, template<class> class PatchField, class GeoMesh>
392 tmp
393 <
395  <
396  typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
397  PatchField,
398  GeoMesh
399  >
400 >
402 {
403  typedef typename GeometricField<Type, PatchField, GeoMesh>::cmptType
404  cmptType;
405 
407 
409  (
411  (
412  "cmptAv(" + gf.name() + ')',
413  gf.mesh(),
414  gf.dimensions()
415  )
416  );
417 
418  cmptAv(CmptAv.ref(), gf);
419 
420  tgf.clear();
421 
422  return CmptAv;
423 }
424 
425 
426 #define UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(returnType, func, gFunc) \
427  \
428 template<class Type, template<class> class PatchField, class GeoMesh> \
429 dimensioned<returnType> func \
430 ( \
431  const GeometricField<Type, PatchField, GeoMesh>& gf \
432 ) \
433 { \
434  return dimensioned<Type> \
435  ( \
436  #func "(" + gf.name() + ')', \
437  gf.dimensions(), \
438  Foam::func(gFunc(gf.primitiveField()), gFunc(gf.boundaryField())) \
439  ); \
440 } \
441  \
442 template<class Type, template<class> class PatchField, class GeoMesh> \
443 dimensioned<returnType> func \
444 ( \
445  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1 \
446 ) \
447 { \
448  dimensioned<returnType> res = func(tgf1()); \
449  tgf1.clear(); \
450  return res; \
451 }
452 
455 
456 #undef UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY
457 
458 
459 #define UNARY_REDUCTION_FUNCTION(returnType, func, gFunc) \
460  \
461 template<class Type, template<class> class PatchField, class GeoMesh> \
462 dimensioned<returnType> func \
463 ( \
464  const GeometricField<Type, PatchField, GeoMesh>& gf \
465 ) \
466 { \
467  return dimensioned<Type> \
468  ( \
469  #func "(" + gf.name() + ')', \
470  gf.dimensions(), \
471  gFunc(gf.primitiveField()) \
472  ); \
473 } \
474  \
475 template<class Type, template<class> class PatchField, class GeoMesh> \
476 dimensioned<returnType> func \
477 ( \
478  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1 \
479 ) \
480 { \
481  dimensioned<returnType> res = func(tgf1()); \
482  tgf1.clear(); \
483  return res; \
484 }
485 
489 
490 #undef UNARY_REDUCTION_FUNCTION
491 
492 
493 BINARY_FUNCTION(Type, Type, Type, max)
494 BINARY_FUNCTION(Type, Type, Type, min)
495 BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
496 BINARY_FUNCTION(Type, Type, Type, cmptDivide)
497 
498 BINARY_TYPE_FUNCTION(Type, Type, Type, max)
499 BINARY_TYPE_FUNCTION(Type, Type, Type, min)
500 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
501 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
502 
503 
504 // * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * //
505 
506 UNARY_OPERATOR(Type, Type, -, negate, transform)
507 
508 BINARY_OPERATOR(Type, Type, scalar, *, '*', multiply)
509 BINARY_OPERATOR(Type, scalar, Type, *, '*', multiply)
510 BINARY_OPERATOR(Type, Type, scalar, /, '|', divide)
511 
512 BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, '*', multiply)
513 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, '*', multiply)
514 
515 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, '|', divide)
516 
517 
518 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
519 
520 #define PRODUCT_OPERATOR(product, op, opFunc) \
521  \
522 template \
523 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
524 void opFunc \
525 ( \
526  GeometricField \
527  <typename product<Type1, Type2>::type, PatchField, GeoMesh>& gf, \
528  const GeometricField<Type1, PatchField, GeoMesh>& gf1, \
529  const GeometricField<Type2, PatchField, GeoMesh>& gf2 \
530 ) \
531 { \
532  Foam::opFunc \
533  ( \
534  gf.primitiveFieldRef(), \
535  gf1.primitiveField(), \
536  gf2.primitiveField() \
537  ); \
538  Foam::opFunc \
539  ( \
540  gf.boundaryFieldRef(), \
541  gf1.boundaryField(), \
542  gf2.boundaryField() \
543  ); \
544 } \
545  \
546 template \
547 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
548 tmp \
549 < \
550  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
551 > \
552 operator op \
553 ( \
554  const GeometricField<Type1, PatchField, GeoMesh>& gf1, \
555  const GeometricField<Type2, PatchField, GeoMesh>& gf2 \
556 ) \
557 { \
558  typedef typename product<Type1, Type2>::type productType; \
559  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes \
560  ( \
561  GeometricField<productType, PatchField, GeoMesh>::New \
562  ( \
563  '(' + gf1.name() + #op + gf2.name() + ')', \
564  gf1.mesh(), \
565  gf1.dimensions() op gf2.dimensions() \
566  ) \
567  ); \
568  \
569  Foam::opFunc(tRes.ref(), gf1, gf2); \
570  \
571  return tRes; \
572 } \
573  \
574 template \
575 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
576 tmp \
577 < \
578  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
579 > \
580 operator op \
581 ( \
582  const GeometricField<Type1, PatchField, GeoMesh>& gf1, \
583  const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tgf2 \
584 ) \
585 { \
586  typedef typename product<Type1, Type2>::type productType; \
587  \
588  const GeometricField<Type2, PatchField, GeoMesh>& gf2 = tgf2(); \
589  \
590  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
591  reuseTmpGeometricField<productType, Type2, PatchField, GeoMesh>::New \
592  ( \
593  tgf2, \
594  '(' + gf1.name() + #op + gf2.name() + ')', \
595  gf1.dimensions() op gf2.dimensions() \
596  ); \
597  \
598  Foam::opFunc(tRes.ref(), gf1, gf2); \
599  \
600  tgf2.clear(); \
601  \
602  return tRes; \
603 } \
604  \
605 template \
606 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
607 tmp \
608 < \
609  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
610 > \
611 operator op \
612 ( \
613  const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tgf1, \
614  const GeometricField<Type2, PatchField, GeoMesh>& gf2 \
615 ) \
616 { \
617  typedef typename product<Type1, Type2>::type productType; \
618  \
619  const GeometricField<Type1, PatchField, GeoMesh>& gf1 = tgf1(); \
620  \
621  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
622  reuseTmpGeometricField<productType, Type1, PatchField, GeoMesh>::New \
623  ( \
624  tgf1, \
625  '(' + gf1.name() + #op + gf2.name() + ')', \
626  gf1.dimensions() op gf2.dimensions() \
627  ); \
628  \
629  Foam::opFunc(tRes.ref(), gf1, gf2); \
630  \
631  tgf1.clear(); \
632  \
633  return tRes; \
634 } \
635  \
636 template \
637 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
638 tmp \
639 < \
640  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
641 > \
642 operator op \
643 ( \
644  const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tgf1, \
645  const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tgf2 \
646 ) \
647 { \
648  typedef typename product<Type1, Type2>::type productType; \
649  \
650  const GeometricField<Type1, PatchField, GeoMesh>& gf1 = tgf1(); \
651  const GeometricField<Type2, PatchField, GeoMesh>& gf2 = tgf2(); \
652  \
653  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
654  reuseTmpTmpGeometricField \
655  <productType, Type1, Type2, PatchField, GeoMesh>::New \
656  ( \
657  tgf1, \
658  tgf2, \
659  '(' + gf1.name() + #op + gf2.name() + ')', \
660  gf1.dimensions() op gf2.dimensions() \
661  ); \
662  \
663  Foam::opFunc(tRes.ref(), gf1, gf2); \
664  \
665  tgf1.clear(); \
666  tgf2.clear(); \
667  \
668  return tRes; \
669 } \
670  \
671 template \
672 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
673 void opFunc \
674 ( \
675  GeometricField \
676  <typename product<Type, Form>::type, PatchField, GeoMesh>& gf, \
677  const GeometricField<Type, PatchField, GeoMesh>& gf1, \
678  const dimensioned<Form>& dvs \
679 ) \
680 { \
681  Foam::opFunc(gf.primitiveFieldRef(), gf1.primitiveField(), dvs.value()); \
682  Foam::opFunc(gf.boundaryFieldRef(), gf1.boundaryField(), dvs.value()); \
683 } \
684  \
685 template \
686 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
687 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
688 operator op \
689 ( \
690  const GeometricField<Type, PatchField, GeoMesh>& gf1, \
691  const dimensioned<Form>& dvs \
692 ) \
693 { \
694  typedef typename product<Type, Form>::type productType; \
695  \
696  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes \
697  ( \
698  GeometricField<productType, PatchField, GeoMesh>::New \
699  ( \
700  '(' + gf1.name() + #op + dvs.name() + ')', \
701  gf1.mesh(), \
702  gf1.dimensions() op dvs.dimensions() \
703  ) \
704  ); \
705  \
706  Foam::opFunc(tRes.ref(), gf1, dvs); \
707  \
708  return tRes; \
709 } \
710  \
711 template \
712 < \
713  class Form, \
714  class Cmpt, \
715  direction nCmpt, \
716  class Type, template<class> class PatchField, \
717  class GeoMesh \
718 > \
719 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
720 operator op \
721 ( \
722  const GeometricField<Type, PatchField, GeoMesh>& gf1, \
723  const VectorSpace<Form,Cmpt,nCmpt>& vs \
724 ) \
725 { \
726  return gf1 op dimensioned<Form>(static_cast<const Form&>(vs)); \
727 } \
728  \
729  \
730 template \
731 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
732 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
733 operator op \
734 ( \
735  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1, \
736  const dimensioned<Form>& dvs \
737 ) \
738 { \
739  typedef typename product<Type, Form>::type productType; \
740  \
741  const GeometricField<Type, PatchField, GeoMesh>& gf1 = tgf1(); \
742  \
743  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
744  reuseTmpGeometricField<productType, Type, PatchField, GeoMesh>::New \
745  ( \
746  tgf1, \
747  '(' + gf1.name() + #op + dvs.name() + ')', \
748  gf1.dimensions() op dvs.dimensions() \
749  ); \
750  \
751  Foam::opFunc(tRes.ref(), gf1, dvs); \
752  \
753  tgf1.clear(); \
754  \
755  return tRes; \
756 } \
757  \
758 template \
759 < \
760  class Form, \
761  class Cmpt, \
762  direction nCmpt, \
763  class Type, template<class> class PatchField, \
764  class GeoMesh \
765 > \
766 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
767 operator op \
768 ( \
769  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1, \
770  const VectorSpace<Form,Cmpt,nCmpt>& vs \
771 ) \
772 { \
773  return tgf1 op dimensioned<Form>(static_cast<const Form&>(vs)); \
774 } \
775  \
776  \
777 template \
778 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
779 void opFunc \
780 ( \
781  GeometricField \
782  <typename product<Form, Type>::type, PatchField, GeoMesh>& gf, \
783  const dimensioned<Form>& dvs, \
784  const GeometricField<Type, PatchField, GeoMesh>& gf1 \
785 ) \
786 { \
787  Foam::opFunc(gf.primitiveFieldRef(), dvs.value(), gf1.primitiveField()); \
788  Foam::opFunc(gf.boundaryFieldRef(), dvs.value(), gf1.boundaryField()); \
789 } \
790  \
791 template \
792 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
793 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
794 operator op \
795 ( \
796  const dimensioned<Form>& dvs, \
797  const GeometricField<Type, PatchField, GeoMesh>& gf1 \
798 ) \
799 { \
800  typedef typename product<Form, Type>::type productType; \
801  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes \
802  ( \
803  GeometricField<productType, PatchField, GeoMesh>::New \
804  ( \
805  '(' + dvs.name() + #op + gf1.name() + ')', \
806  gf1.mesh(), \
807  dvs.dimensions() op gf1.dimensions() \
808  ) \
809  ); \
810  \
811  Foam::opFunc(tRes.ref(), dvs, gf1); \
812  \
813  return tRes; \
814 } \
815  \
816 template \
817 < \
818  class Form, \
819  class Cmpt, \
820  direction nCmpt, \
821  class Type, template<class> class PatchField, \
822  class GeoMesh \
823 > \
824 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
825 operator op \
826 ( \
827  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
828  const GeometricField<Type, PatchField, GeoMesh>& gf1 \
829 ) \
830 { \
831  return dimensioned<Form>(static_cast<const Form&>(vs)) op gf1; \
832 } \
833  \
834 template \
835 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
836 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
837 operator op \
838 ( \
839  const dimensioned<Form>& dvs, \
840  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1 \
841 ) \
842 { \
843  typedef typename product<Form, Type>::type productType; \
844  \
845  const GeometricField<Type, PatchField, GeoMesh>& gf1 = tgf1(); \
846  \
847  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
848  reuseTmpGeometricField<productType, Type, PatchField, GeoMesh>::New \
849  ( \
850  tgf1, \
851  '(' + dvs.name() + #op + gf1.name() + ')', \
852  dvs.dimensions() op gf1.dimensions() \
853  ); \
854  \
855  Foam::opFunc(tRes.ref(), dvs, gf1); \
856  \
857  tgf1.clear(); \
858  \
859  return tRes; \
860 } \
861  \
862 template \
863 < \
864  class Form, \
865  class Cmpt, \
866  direction nCmpt, \
867  class Type, template<class> class PatchField, \
868  class GeoMesh \
869 > \
870 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
871 operator op \
872 ( \
873  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
874  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1 \
875 ) \
876 { \
877  return dimensioned<Form>(static_cast<const Form&>(vs)) op tgf1; \
878 }
879 
882 
887 
888 #undef PRODUCT_OPERATOR
889 
890 
891 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
892 
893 } // End namespace Foam
894 
895 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
896 
897 #include "undefFieldFunctionsM.H"
898 
899 // ************************************************************************* //
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
#define UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(returnType, func, gFunc)
scalar gSumMag(const FieldField< Field, Type > &f)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
const word & name() const
Return name.
Definition: IOobject.H:303
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
#define UNARY_REDUCTION_FUNCTION(returnType, func, gFunc)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMin(const FieldField< Field, Type > &f)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
uint8_t direction
Definition: direction.H:45
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#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)
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Generic GeometricField class.
#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)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
#define PRODUCT_OPERATOR(product, op, opFunc)
const dimensionSet & dimensions() const
Return dimensions.
Type gSum(const FieldField< Field, Type > &f)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
const Mesh & mesh() const
Return mesh.
Type gMax(const FieldField< Field, Type > &f)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
symmTypeOfRank< typename pTraits< arg1 >::cmptType, arg2 *direction(pTraits< arg1 >::rank) >::type type
Definition: products.H:136
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void negate(FieldField< Field, Type > &res, const FieldField< Field, Type > &f)
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Type gAverage(const FieldField< Field, Type > &f)
dimensioned< scalar > mag(const dimensioned< Type > &)
#define UNARY_OPERATOR(ReturnType, Type1, Op, OpFunc, Dfunc)
A class for managing temporary objects.
Definition: PtrList.H:53
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477