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-2021 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 UNARY_FUNCTION(Type, Type, cmptMag, cmptMag);
426 
427 
428 #define UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(returnType, func, gFunc) \
429  \
430 template<class Type, template<class> class PatchField, class GeoMesh> \
431 dimensioned<returnType> func \
432 ( \
433  const GeometricField<Type, PatchField, GeoMesh>& gf \
434 ) \
435 { \
436  return dimensioned<Type> \
437  ( \
438  #func "(" + gf.name() + ')', \
439  gf.dimensions(), \
440  Foam::func(gFunc(gf.primitiveField()), gFunc(gf.boundaryField())) \
441  ); \
442 } \
443  \
444 template<class Type, template<class> class PatchField, class GeoMesh> \
445 dimensioned<returnType> func \
446 ( \
447  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1 \
448 ) \
449 { \
450  dimensioned<returnType> res = func(tgf1()); \
451  tgf1.clear(); \
452  return res; \
453 }
454 
457 
458 #undef UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY
459 
460 
461 #define UNARY_REDUCTION_FUNCTION(returnType, func, gFunc) \
462  \
463 template<class Type, template<class> class PatchField, class GeoMesh> \
464 dimensioned<returnType> func \
465 ( \
466  const GeometricField<Type, PatchField, GeoMesh>& gf \
467 ) \
468 { \
469  return dimensioned<Type> \
470  ( \
471  #func "(" + gf.name() + ')', \
472  gf.dimensions(), \
473  gFunc(gf.primitiveField()) \
474  ); \
475 } \
476  \
477 template<class Type, template<class> class PatchField, class GeoMesh> \
478 dimensioned<returnType> func \
479 ( \
480  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1 \
481 ) \
482 { \
483  dimensioned<returnType> res = func(tgf1()); \
484  tgf1.clear(); \
485  return res; \
486 }
487 
491 
492 #undef UNARY_REDUCTION_FUNCTION
493 
494 
495 BINARY_FUNCTION(Type, Type, Type, max)
496 BINARY_FUNCTION(Type, Type, Type, min)
497 BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
498 BINARY_FUNCTION(Type, Type, Type, cmptDivide)
499 
500 BINARY_TYPE_FUNCTION(Type, Type, Type, max)
501 BINARY_TYPE_FUNCTION(Type, Type, Type, min)
502 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
503 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
504 
505 
506 // * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * //
507 
508 UNARY_OPERATOR(Type, Type, -, negate, transform)
509 
510 BINARY_OPERATOR(Type, Type, scalar, *, '*', multiply)
511 BINARY_OPERATOR(Type, scalar, Type, *, '*', multiply)
512 BINARY_OPERATOR(Type, Type, scalar, /, '|', divide)
513 
514 BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, '*', multiply)
515 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, '*', multiply)
516 
517 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, '|', divide)
518 
519 
520 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
521 
522 #define PRODUCT_OPERATOR(product, op, opFunc) \
523  \
524 template \
525 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
526 void opFunc \
527 ( \
528  GeometricField \
529  <typename product<Type1, Type2>::type, PatchField, GeoMesh>& gf, \
530  const GeometricField<Type1, PatchField, GeoMesh>& gf1, \
531  const GeometricField<Type2, PatchField, GeoMesh>& gf2 \
532 ) \
533 { \
534  Foam::opFunc \
535  ( \
536  gf.primitiveFieldRef(), \
537  gf1.primitiveField(), \
538  gf2.primitiveField() \
539  ); \
540  Foam::opFunc \
541  ( \
542  gf.boundaryFieldRef(), \
543  gf1.boundaryField(), \
544  gf2.boundaryField() \
545  ); \
546 } \
547  \
548 template \
549 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
550 tmp \
551 < \
552  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
553 > \
554 operator op \
555 ( \
556  const GeometricField<Type1, PatchField, GeoMesh>& gf1, \
557  const GeometricField<Type2, PatchField, GeoMesh>& gf2 \
558 ) \
559 { \
560  typedef typename product<Type1, Type2>::type productType; \
561  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes \
562  ( \
563  GeometricField<productType, PatchField, GeoMesh>::New \
564  ( \
565  '(' + gf1.name() + #op + gf2.name() + ')', \
566  gf1.mesh(), \
567  gf1.dimensions() op gf2.dimensions() \
568  ) \
569  ); \
570  \
571  Foam::opFunc(tRes.ref(), gf1, gf2); \
572  \
573  return tRes; \
574 } \
575  \
576 template \
577 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
578 tmp \
579 < \
580  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
581 > \
582 operator op \
583 ( \
584  const GeometricField<Type1, PatchField, GeoMesh>& gf1, \
585  const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tgf2 \
586 ) \
587 { \
588  typedef typename product<Type1, Type2>::type productType; \
589  \
590  const GeometricField<Type2, PatchField, GeoMesh>& gf2 = tgf2(); \
591  \
592  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
593  reuseTmpGeometricField<productType, Type2, PatchField, GeoMesh>::New \
594  ( \
595  tgf2, \
596  '(' + gf1.name() + #op + gf2.name() + ')', \
597  gf1.dimensions() op gf2.dimensions() \
598  ); \
599  \
600  Foam::opFunc(tRes.ref(), gf1, gf2); \
601  \
602  tgf2.clear(); \
603  \
604  return tRes; \
605 } \
606  \
607 template \
608 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
609 tmp \
610 < \
611  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
612 > \
613 operator op \
614 ( \
615  const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tgf1, \
616  const GeometricField<Type2, PatchField, GeoMesh>& gf2 \
617 ) \
618 { \
619  typedef typename product<Type1, Type2>::type productType; \
620  \
621  const GeometricField<Type1, PatchField, GeoMesh>& gf1 = tgf1(); \
622  \
623  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
624  reuseTmpGeometricField<productType, Type1, PatchField, GeoMesh>::New \
625  ( \
626  tgf1, \
627  '(' + gf1.name() + #op + gf2.name() + ')', \
628  gf1.dimensions() op gf2.dimensions() \
629  ); \
630  \
631  Foam::opFunc(tRes.ref(), gf1, gf2); \
632  \
633  tgf1.clear(); \
634  \
635  return tRes; \
636 } \
637  \
638 template \
639 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
640 tmp \
641 < \
642  GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
643 > \
644 operator op \
645 ( \
646  const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tgf1, \
647  const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tgf2 \
648 ) \
649 { \
650  typedef typename product<Type1, Type2>::type productType; \
651  \
652  const GeometricField<Type1, PatchField, GeoMesh>& gf1 = tgf1(); \
653  const GeometricField<Type2, PatchField, GeoMesh>& gf2 = tgf2(); \
654  \
655  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
656  reuseTmpTmpGeometricField \
657  <productType, Type1, Type2, PatchField, GeoMesh>::New \
658  ( \
659  tgf1, \
660  tgf2, \
661  '(' + gf1.name() + #op + gf2.name() + ')', \
662  gf1.dimensions() op gf2.dimensions() \
663  ); \
664  \
665  Foam::opFunc(tRes.ref(), gf1, gf2); \
666  \
667  tgf1.clear(); \
668  tgf2.clear(); \
669  \
670  return tRes; \
671 } \
672  \
673 template \
674 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
675 void opFunc \
676 ( \
677  GeometricField \
678  <typename product<Type, Form>::type, PatchField, GeoMesh>& gf, \
679  const GeometricField<Type, PatchField, GeoMesh>& gf1, \
680  const dimensioned<Form>& dvs \
681 ) \
682 { \
683  Foam::opFunc(gf.primitiveFieldRef(), gf1.primitiveField(), dvs.value()); \
684  Foam::opFunc(gf.boundaryFieldRef(), gf1.boundaryField(), dvs.value()); \
685 } \
686  \
687 template \
688 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
689 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
690 operator op \
691 ( \
692  const GeometricField<Type, PatchField, GeoMesh>& gf1, \
693  const dimensioned<Form>& dvs \
694 ) \
695 { \
696  typedef typename product<Type, Form>::type productType; \
697  \
698  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes \
699  ( \
700  GeometricField<productType, PatchField, GeoMesh>::New \
701  ( \
702  '(' + gf1.name() + #op + dvs.name() + ')', \
703  gf1.mesh(), \
704  gf1.dimensions() op dvs.dimensions() \
705  ) \
706  ); \
707  \
708  Foam::opFunc(tRes.ref(), gf1, dvs); \
709  \
710  return tRes; \
711 } \
712  \
713 template \
714 < \
715  class Form, \
716  class Cmpt, \
717  direction nCmpt, \
718  class Type, template<class> class PatchField, \
719  class GeoMesh \
720 > \
721 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
722 operator op \
723 ( \
724  const GeometricField<Type, PatchField, GeoMesh>& gf1, \
725  const VectorSpace<Form,Cmpt,nCmpt>& vs \
726 ) \
727 { \
728  return gf1 op dimensioned<Form>(static_cast<const Form&>(vs)); \
729 } \
730  \
731  \
732 template \
733 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
734 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
735 operator op \
736 ( \
737  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1, \
738  const dimensioned<Form>& dvs \
739 ) \
740 { \
741  typedef typename product<Type, Form>::type productType; \
742  \
743  const GeometricField<Type, PatchField, GeoMesh>& gf1 = tgf1(); \
744  \
745  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
746  reuseTmpGeometricField<productType, Type, PatchField, GeoMesh>::New \
747  ( \
748  tgf1, \
749  '(' + gf1.name() + #op + dvs.name() + ')', \
750  gf1.dimensions() op dvs.dimensions() \
751  ); \
752  \
753  Foam::opFunc(tRes.ref(), gf1, dvs); \
754  \
755  tgf1.clear(); \
756  \
757  return tRes; \
758 } \
759  \
760 template \
761 < \
762  class Form, \
763  class Cmpt, \
764  direction nCmpt, \
765  class Type, template<class> class PatchField, \
766  class GeoMesh \
767 > \
768 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
769 operator op \
770 ( \
771  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1, \
772  const VectorSpace<Form,Cmpt,nCmpt>& vs \
773 ) \
774 { \
775  return tgf1 op dimensioned<Form>(static_cast<const Form&>(vs)); \
776 } \
777  \
778  \
779 template \
780 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
781 void opFunc \
782 ( \
783  GeometricField \
784  <typename product<Form, Type>::type, PatchField, GeoMesh>& gf, \
785  const dimensioned<Form>& dvs, \
786  const GeometricField<Type, PatchField, GeoMesh>& gf1 \
787 ) \
788 { \
789  Foam::opFunc(gf.primitiveFieldRef(), dvs.value(), gf1.primitiveField()); \
790  Foam::opFunc(gf.boundaryFieldRef(), dvs.value(), gf1.boundaryField()); \
791 } \
792  \
793 template \
794 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
795 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
796 operator op \
797 ( \
798  const dimensioned<Form>& dvs, \
799  const GeometricField<Type, PatchField, GeoMesh>& gf1 \
800 ) \
801 { \
802  typedef typename product<Form, Type>::type productType; \
803  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes \
804  ( \
805  GeometricField<productType, PatchField, GeoMesh>::New \
806  ( \
807  '(' + dvs.name() + #op + gf1.name() + ')', \
808  gf1.mesh(), \
809  dvs.dimensions() op gf1.dimensions() \
810  ) \
811  ); \
812  \
813  Foam::opFunc(tRes.ref(), dvs, gf1); \
814  \
815  return tRes; \
816 } \
817  \
818 template \
819 < \
820  class Form, \
821  class Cmpt, \
822  direction nCmpt, \
823  class Type, template<class> class PatchField, \
824  class GeoMesh \
825 > \
826 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
827 operator op \
828 ( \
829  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
830  const GeometricField<Type, PatchField, GeoMesh>& gf1 \
831 ) \
832 { \
833  return dimensioned<Form>(static_cast<const Form&>(vs)) op gf1; \
834 } \
835  \
836 template \
837 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
838 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
839 operator op \
840 ( \
841  const dimensioned<Form>& dvs, \
842  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1 \
843 ) \
844 { \
845  typedef typename product<Form, Type>::type productType; \
846  \
847  const GeometricField<Type, PatchField, GeoMesh>& gf1 = tgf1(); \
848  \
849  tmp<GeometricField<productType, PatchField, GeoMesh>> tRes = \
850  reuseTmpGeometricField<productType, Type, PatchField, GeoMesh>::New \
851  ( \
852  tgf1, \
853  '(' + dvs.name() + #op + gf1.name() + ')', \
854  dvs.dimensions() op gf1.dimensions() \
855  ); \
856  \
857  Foam::opFunc(tRes.ref(), dvs, gf1); \
858  \
859  tgf1.clear(); \
860  \
861  return tRes; \
862 } \
863  \
864 template \
865 < \
866  class Form, \
867  class Cmpt, \
868  direction nCmpt, \
869  class Type, template<class> class PatchField, \
870  class GeoMesh \
871 > \
872 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
873 operator op \
874 ( \
875  const VectorSpace<Form,Cmpt,nCmpt>& vs, \
876  const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf1 \
877 ) \
878 { \
879  return dimensioned<Form>(static_cast<const Form&>(vs)) op tgf1; \
880 }
881 
884 
889 
890 #undef PRODUCT_OPERATOR
891 
892 
893 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
894 
895 } // End namespace Foam
896 
897 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
898 
899 #include "undefFieldFunctionsM.H"
900 
901 // ************************************************************************* //
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)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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:315
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
#define UNARY_REDUCTION_FUNCTION(returnType, func, gFunc)
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)
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:275
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
UNARY_FUNCTION(Type, Type, cmptMag, cmptMag)
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:483