1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "GeometricFieldReuseFunctions.H"
29 template<class Type, template<class> class PatchField, class GeoMesh>
30 #include "GeometricFieldFunctionsM.C"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 // * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * //
39 template<class Type, template<class> class PatchField, class GeoMesh>
44 typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
48 const GeometricField<Type, PatchField, GeoMesh>& gf,
52 component(gcf.internalField(), gf.internalField(), d);
53 component(gcf.boundaryField(), gf.boundaryField(), d);
57 template<class Type, template<class> class PatchField, class GeoMesh>
60 GeometricField<Type, PatchField, GeoMesh>& gf,
61 const GeometricField<Type, PatchField, GeoMesh>& gf1
64 T(gf.internalField(), gf1.internalField());
65 T(gf.boundaryField(), gf1.boundaryField());
69 template<class Type, template<class> class PatchField, class GeoMesh, int r>
72 GeometricField<typename powProduct<Type, r>::type, PatchField, GeoMesh>& gf,
73 const GeometricField<Type, PatchField, GeoMesh>& gf1
76 pow(gf.internalField(), gf1.internalField(), r);
77 pow(gf.boundaryField(), gf1.boundaryField(), r);
80 template<class Type, template<class> class PatchField, class GeoMesh, int r>
81 tmp<GeometricField<typename powProduct<Type, r>::type, PatchField, GeoMesh> >
84 const GeometricField<Type, PatchField, GeoMesh>& gf,
85 typename powProduct<Type, r>::type
88 typedef typename powProduct<Type, r>::type powProductType;
90 tmp<GeometricField<powProductType, PatchField, GeoMesh> > tPow
92 new GeometricField<powProductType, PatchField, GeoMesh>
96 "pow(" + gf.name() + ',' + name(r) + ')',
103 pow(gf.dimensions(), r)
107 pow<Type, r, PatchField, GeoMesh>(tPow(), gf);
113 template<class Type, template<class> class PatchField, class GeoMesh, int r>
114 tmp<GeometricField<typename powProduct<Type, r>::type, PatchField, GeoMesh> >
117 const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf,
118 typename powProduct<Type, r>::type
121 typedef typename powProduct<Type, r>::type powProductType;
123 const GeometricField<Type, PatchField, GeoMesh>& gf = tgf();
125 tmp<GeometricField<powProductType, PatchField, GeoMesh> > tPow
127 new GeometricField<powProductType, PatchField, GeoMesh>
131 "pow(" + gf.name() + ',' + name(r) + ')',
138 pow(gf.dimensions(), r)
142 pow<Type, r, PatchField, GeoMesh>(tPow(), gf);
150 template<class Type, template<class> class PatchField, class GeoMesh>
154 <typename outerProduct<Type, Type>::type, PatchField, GeoMesh>& gf,
155 const GeometricField<Type, PatchField, GeoMesh>& gf1
158 sqr(gf.internalField(), gf1.internalField());
159 sqr(gf.boundaryField(), gf1.boundaryField());
162 template<class Type, template<class> class PatchField, class GeoMesh>
167 typename outerProduct<Type, Type>::type,
172 sqr(const GeometricField<Type, PatchField, GeoMesh>& gf)
174 typedef typename outerProduct<Type, Type>::type outerProductType;
176 tmp<GeometricField<outerProductType, PatchField, GeoMesh> > tSqr
178 new GeometricField<outerProductType, PatchField, GeoMesh>
182 "sqr(" + gf.name() + ')',
198 template<class Type, template<class> class PatchField, class GeoMesh>
203 typename outerProduct<Type, Type>::type,
208 sqr(const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf)
210 typedef typename outerProduct<Type, Type>::type outerProductType;
212 const GeometricField<Type, PatchField, GeoMesh>& gf = tgf();
214 tmp<GeometricField<outerProductType, PatchField, GeoMesh> > tSqr
216 new GeometricField<outerProductType, PatchField, GeoMesh>
220 "sqr(" + gf.name() + ')',
239 template<class Type, template<class> class PatchField, class GeoMesh>
242 GeometricField<scalar, PatchField, GeoMesh>& gsf,
243 const GeometricField<Type, PatchField, GeoMesh>& gf
246 magSqr(gsf.internalField(), gf.internalField());
247 magSqr(gsf.boundaryField(), gf.boundaryField());
250 template<class Type, template<class> class PatchField, class GeoMesh>
251 tmp<GeometricField<scalar, PatchField, GeoMesh> > magSqr
253 const GeometricField<Type, PatchField, GeoMesh>& gf
256 tmp<GeometricField<scalar, PatchField, GeoMesh> > tMagSqr
258 new GeometricField<scalar, PatchField, GeoMesh>
262 "magSqr(" + gf.name() + ')',
273 magSqr(tMagSqr(), gf);
278 template<class Type, template<class> class PatchField, class GeoMesh>
279 tmp<GeometricField<scalar, PatchField, GeoMesh> > magSqr
281 const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf
284 const GeometricField<Type, PatchField, GeoMesh>& gf = tgf();
286 tmp<GeometricField<scalar, PatchField, GeoMesh> > tMagSqr
288 new GeometricField<scalar, PatchField, GeoMesh>
292 "magSqr(" + gf.name() + ')',
303 magSqr(tMagSqr(), gf);
311 template<class Type, template<class> class PatchField, class GeoMesh>
314 GeometricField<scalar, PatchField, GeoMesh>& gsf,
315 const GeometricField<Type, PatchField, GeoMesh>& gf
318 mag(gsf.internalField(), gf.internalField());
319 mag(gsf.boundaryField(), gf.boundaryField());
322 template<class Type, template<class> class PatchField, class GeoMesh>
323 tmp<GeometricField<scalar, PatchField, GeoMesh> > mag
325 const GeometricField<Type, PatchField, GeoMesh>& gf
328 tmp<GeometricField<scalar, PatchField, GeoMesh> > tMag
330 new GeometricField<scalar, PatchField, GeoMesh>
334 "mag(" + gf.name() + ')',
350 template<class Type, template<class> class PatchField, class GeoMesh>
351 tmp<GeometricField<scalar, PatchField, GeoMesh> > mag
353 const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf
356 const GeometricField<Type, PatchField, GeoMesh>& gf = tgf();
358 tmp<GeometricField<scalar, PatchField, GeoMesh> > tMag
360 new GeometricField<scalar, PatchField, GeoMesh>
364 "mag(" + gf.name() + ')',
383 template<class Type, template<class> class PatchField, class GeoMesh>
388 typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
392 const GeometricField<Type, PatchField, GeoMesh>& gf
395 cmptAv(gcf.internalField(), gf.internalField());
396 cmptAv(gcf.boundaryField(), gf.boundaryField());
399 template<class Type, template<class> class PatchField, class GeoMesh>
404 typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
409 cmptAv(const GeometricField<Type, PatchField, GeoMesh>& gf)
411 typedef typename GeometricField<Type, PatchField, GeoMesh>::cmptType
414 tmp<GeometricField<cmptType, PatchField, GeoMesh> > CmptAv
416 new GeometricField<scalar, PatchField, GeoMesh>
420 "cmptAv(" + gf.name() + ')',
431 cmptAv(CmptAv(), gf);
436 template<class Type, template<class> class PatchField, class GeoMesh>
441 typename GeometricField<Type, PatchField, GeoMesh>::cmptType,
446 cmptAv(const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf)
448 typedef typename GeometricField<Type, PatchField, GeoMesh>::cmptType
451 const GeometricField<Type, PatchField, GeoMesh>& gf = tgf();
453 tmp<GeometricField<cmptType, PatchField, GeoMesh> > CmptAv
455 new GeometricField<scalar, PatchField, GeoMesh>
459 "cmptAv(" + gf.name() + ')',
470 cmptAv(CmptAv(), gf);
478 #define UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(returnType, func, gFunc) \
480 template<class Type, template<class> class PatchField, class GeoMesh> \
481 dimensioned<returnType> func \
483 const GeometricField<Type, PatchField, GeoMesh>& gf \
486 return dimensioned<Type> \
488 #func "(" + gf.name() + ')', \
490 Foam::func(gFunc(gf.internalField()), gFunc(gf.boundaryField())) \
494 template<class Type, template<class> class PatchField, class GeoMesh> \
495 dimensioned<returnType> func \
497 const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf1 \
500 dimensioned<returnType> res = func(tgf1()); \
505 UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(Type, max, gMax)
506 UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(Type, min, gMin)
508 #undef UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY
511 #define UNARY_REDUCTION_FUNCTION(returnType, func, gFunc) \
513 template<class Type, template<class> class PatchField, class GeoMesh> \
514 dimensioned<returnType> func \
516 const GeometricField<Type, PatchField, GeoMesh>& gf \
519 return dimensioned<Type> \
521 #func "(" + gf.name() + ')', \
523 gFunc(gf.internalField()) \
527 template<class Type, template<class> class PatchField, class GeoMesh> \
528 dimensioned<returnType> func \
530 const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf1 \
533 dimensioned<returnType> res = func(tgf1()); \
538 UNARY_REDUCTION_FUNCTION(Type, sum, gSum)
539 UNARY_REDUCTION_FUNCTION(scalar, sumMag, gSumMag)
540 UNARY_REDUCTION_FUNCTION(Type, average, gAverage)
542 #undef UNARY_REDUCTION_FUNCTION
545 BINARY_FUNCTION(Type, Type, Type, max)
546 BINARY_FUNCTION(Type, Type, Type, min)
547 BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
548 BINARY_FUNCTION(Type, Type, Type, cmptDivide)
550 BINARY_TYPE_FUNCTION(Type, Type, Type, max)
551 BINARY_TYPE_FUNCTION(Type, Type, Type, min)
552 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
553 BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
556 // * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * //
558 UNARY_OPERATOR(Type, Type, -, negate, transform)
560 #ifndef __INTEL_COMPILER
561 BINARY_OPERATOR(Type, Type, scalar, *, '*', multiply)
562 BINARY_OPERATOR(Type, scalar, Type, *, '*', multiply)
564 BINARY_OPERATOR(Type, Type, scalar, /, '|', divide)
566 BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, '*', multiply)
567 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, '*', multiply)
569 BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, '|', divide)
572 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
574 #define PRODUCT_OPERATOR(product, op, opFunc) \
577 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
581 <typename product<Type1, Type2>::type, PatchField, GeoMesh>& gf, \
582 const GeometricField<Type1, PatchField, GeoMesh>& gf1, \
583 const GeometricField<Type2, PatchField, GeoMesh>& gf2 \
586 Foam::opFunc(gf.internalField(), gf1.internalField(), gf2.internalField());\
587 Foam::opFunc(gf.boundaryField(), gf1.boundaryField(), gf2.boundaryField());\
591 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
594 GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
598 const GeometricField<Type1, PatchField, GeoMesh>& gf1, \
599 const GeometricField<Type2, PatchField, GeoMesh>& gf2 \
602 typedef typename product<Type1, Type2>::type productType; \
603 tmp<GeometricField<productType, PatchField, GeoMesh> > tRes \
605 new GeometricField<productType, PatchField, GeoMesh> \
609 '(' + gf1.name() + #op + gf2.name() + ')', \
616 gf1.dimensions() op gf2.dimensions() \
620 Foam::opFunc(tRes(), gf1, gf2); \
626 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
629 GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
633 const GeometricField<Type1, PatchField, GeoMesh>& gf1, \
634 const tmp<GeometricField<Type2, PatchField, GeoMesh> >& tgf2 \
637 typedef typename product<Type1, Type2>::type productType; \
639 const GeometricField<Type2, PatchField, GeoMesh>& gf2 = tgf2(); \
641 tmp<GeometricField<productType, PatchField, GeoMesh> > tRes = \
642 reuseTmpGeometricField<productType, Type2, PatchField, GeoMesh>::New \
645 '(' + gf1.name() + #op + gf2.name() + ')', \
646 gf1.dimensions() op gf2.dimensions() \
649 Foam::opFunc(tRes(), gf1, gf2); \
651 reuseTmpGeometricField<productType, Type2, PatchField, GeoMesh> \
658 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
661 GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
665 const tmp<GeometricField<Type1, PatchField, GeoMesh> >& tgf1, \
666 const GeometricField<Type2, PatchField, GeoMesh>& gf2 \
669 typedef typename product<Type1, Type2>::type productType; \
671 const GeometricField<Type1, PatchField, GeoMesh>& gf1 = tgf1(); \
673 tmp<GeometricField<productType, PatchField, GeoMesh> > tRes = \
674 reuseTmpGeometricField<productType, Type1, PatchField, GeoMesh>::New \
677 '(' + gf1.name() + #op + gf2.name() + ')', \
678 gf1.dimensions() op gf2.dimensions() \
681 Foam::opFunc(tRes(), gf1, gf2); \
683 reuseTmpGeometricField<productType, Type1, PatchField, GeoMesh> \
690 <class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
693 GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
697 const tmp<GeometricField<Type1, PatchField, GeoMesh> >& tgf1, \
698 const tmp<GeometricField<Type2, PatchField, GeoMesh> >& tgf2 \
701 typedef typename product<Type1, Type2>::type productType; \
703 const GeometricField<Type1, PatchField, GeoMesh>& gf1 = tgf1(); \
704 const GeometricField<Type2, PatchField, GeoMesh>& gf2 = tgf2(); \
706 tmp<GeometricField<productType, PatchField, GeoMesh> > tRes = \
707 reuseTmpTmpGeometricField \
708 <productType, Type1, Type1, Type2, PatchField, GeoMesh>::New \
712 '(' + gf1.name() + #op + gf2.name() + ')', \
713 gf1.dimensions() op gf2.dimensions() \
716 Foam::opFunc(tRes(), gf1, gf2); \
718 reuseTmpTmpGeometricField \
719 <productType, Type1, Type1, Type2, PatchField, GeoMesh> \
720 ::clear(tgf1, tgf2); \
726 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
730 <typename product<Type, Form>::type, PatchField, GeoMesh>& gf, \
731 const GeometricField<Type, PatchField, GeoMesh>& gf1, \
732 const dimensioned<Form>& dvs \
735 Foam::opFunc(gf.internalField(), gf1.internalField(), dvs.value()); \
736 Foam::opFunc(gf.boundaryField(), gf1.boundaryField(), dvs.value()); \
740 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
741 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh> > \
744 const GeometricField<Type, PatchField, GeoMesh>& gf1, \
745 const dimensioned<Form>& dvs \
748 typedef typename product<Type, Form>::type productType; \
750 tmp<GeometricField<productType, PatchField, GeoMesh> > tRes \
752 new GeometricField<productType, PatchField, GeoMesh> \
756 '(' + gf1.name() + #op + dvs.name() + ')', \
763 gf1.dimensions() op dvs.dimensions() \
767 Foam::opFunc(tRes(), gf1, dvs); \
777 class Type, template<class> class PatchField, \
780 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh> > \
783 const GeometricField<Type, PatchField, GeoMesh>& gf1, \
784 const VectorSpace<Form,Cmpt,nCmpt>& vs \
787 return gf1 op dimensioned<Form>(static_cast<const Form&>(vs)); \
792 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
793 tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh> > \
796 const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf1, \
797 const dimensioned<Form>& dvs \
800 typedef typename product<Type, Form>::type productType; \
802 const GeometricField<Type, PatchField, GeoMesh>& gf1 = tgf1(); \
804 tmp<GeometricField<productType, PatchField, GeoMesh> > tRes = \
805 reuseTmpGeometricField<productType, Type, PatchField, GeoMesh>::New \
808 '(' + gf1.name() + #op + dvs.name() + ')', \
809 gf1.dimensions() op dvs.dimensions() \
812 Foam::opFunc(tRes(), gf1, dvs); \
814 reuseTmpGeometricField<productType, Type, PatchField, GeoMesh> \
825 class Type, template<class> class PatchField, \
828 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh> > \
831 const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf1, \
832 const VectorSpace<Form,Cmpt,nCmpt>& vs \
835 return tgf1 op dimensioned<Form>(static_cast<const Form&>(vs)); \
840 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
844 <typename product<Form, Type>::type, PatchField, GeoMesh>& gf, \
845 const dimensioned<Form>& dvs, \
846 const GeometricField<Type, PatchField, GeoMesh>& gf1 \
849 Foam::opFunc(gf.internalField(), dvs.value(), gf1.internalField()); \
850 Foam::opFunc(gf.boundaryField(), dvs.value(), gf1.boundaryField()); \
854 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
855 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh> > \
858 const dimensioned<Form>& dvs, \
859 const GeometricField<Type, PatchField, GeoMesh>& gf1 \
862 typedef typename product<Form, Type>::type productType; \
863 tmp<GeometricField<productType, PatchField, GeoMesh> > tRes \
865 new GeometricField<productType, PatchField, GeoMesh> \
869 '(' + dvs.name() + #op + gf1.name() + ')', \
876 dvs.dimensions() op gf1.dimensions() \
880 Foam::opFunc(tRes(), dvs, gf1); \
890 class Type, template<class> class PatchField, \
893 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh> > \
896 const VectorSpace<Form,Cmpt,nCmpt>& vs, \
897 const GeometricField<Type, PatchField, GeoMesh>& gf1 \
900 return dimensioned<Form>(static_cast<const Form&>(vs)) op gf1; \
904 <class Form, class Type, template<class> class PatchField, class GeoMesh> \
905 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh> > \
908 const dimensioned<Form>& dvs, \
909 const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf1 \
912 typedef typename product<Form, Type>::type productType; \
914 const GeometricField<Type, PatchField, GeoMesh>& gf1 = tgf1(); \
916 tmp<GeometricField<productType, PatchField, GeoMesh> > tRes = \
917 reuseTmpGeometricField<productType, Type, PatchField, GeoMesh>::New \
920 '(' + dvs.name() + #op + gf1.name() + ')', \
921 dvs.dimensions() op gf1.dimensions() \
924 Foam::opFunc(tRes(), dvs, gf1); \
926 reuseTmpGeometricField<productType, Type, PatchField, GeoMesh> \
937 class Type, template<class> class PatchField, \
940 tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh> > \
943 const VectorSpace<Form,Cmpt,nCmpt>& vs, \
944 const tmp<GeometricField<Type, PatchField, GeoMesh> >& tgf1 \
947 return dimensioned<Form>(static_cast<const Form&>(vs)) op tgf1; \
950 PRODUCT_OPERATOR(typeOfSum, +, add)
951 PRODUCT_OPERATOR(typeOfSum, -, subtract)
953 PRODUCT_OPERATOR(outerProduct, *, outer)
954 PRODUCT_OPERATOR(crossProduct, ^, cross)
955 PRODUCT_OPERATOR(innerProduct, &, dot)
956 PRODUCT_OPERATOR(scalarProduct, &&, dotdot)
958 #undef PRODUCT_OPERATOR
961 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
963 } // End namespace Foam
965 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
967 #include "undefFieldFunctionsM.H"
969 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //