Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37namespace Foam
38{
39
40// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
41
42template<class Type, template<class> class PatchField, class GeoMesh>
43void component
44(
46 <
48 PatchField,
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
65template<class Type, template<class> class PatchField, class GeoMesh>
66void T
67(
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
82template
83<
84 class Type,
85 template<class> class PatchField,
86 class GeoMesh,
87 direction r
88>
89void 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);
99 result.correctLocalBoundaryConditions();
101 {
102 result.boundaryField().check();
103 }
104}
105
106
107template
109 class Type,
110 template<class> class PatchField,
111 class GeoMesh,
112 direction r
113>
115pow
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
137template
139 class Type,
140 template<class> class PatchField,
141 class GeoMesh,
142 direction r
143>
145pow
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);
165 tf1.clear();
166 return tres;
167}
168
169
170template<class Type, template<class> class PatchField, class GeoMesh>
171void 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());
183 result.correctLocalBoundaryConditions();
185 {
186 result.boundaryField().check();
187 }
188}
189
190
191template<class Type, template<class> class PatchField, class GeoMesh>
192tmp
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
219template<class Type, template<class> class PatchField, class GeoMesh>
220tmp
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);
245 tf1.clear();
246 return tres;
247}
248
249
250template<class Type, template<class> class PatchField, class GeoMesh>
251void 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());
260 result.correctLocalBoundaryConditions();
262 {
263 result.boundaryField().check();
264 }
265}
266
267
268template<class Type, template<class> class PatchField, class GeoMesh>
270magSqr
271(
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);
286
287 return tres;
288}
289
290template<class Type, template<class> class PatchField, class GeoMesh>
291tmp<GeometricField<typename typeOfMag<Type>::type, PatchField, GeoMesh>>
292magSqr
293(
294 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1
295)
296{
297 auto tres = magSqr(tf1.cref());
298 tf1.clear();
299
300 return tres;
301}
302
303
304template<class Type, template<class> class PatchField, class GeoMesh>
305void 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());
314 result.correctLocalBoundaryConditions();
316 {
317 result.boundaryField().check();
318 }
319}
320
321
322template<class Type, template<class> class PatchField, class GeoMesh>
324mag
325(
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);
340
341 return tres;
342}
343
344template<class Type, template<class> class PatchField, class GeoMesh>
345tmp<GeometricField<typename typeOfMag<Type>::type, PatchField, GeoMesh>>
346mag
347(
348 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1
349)
350{
351 auto tres = mag(tf1.cref());
352 tf1.clear();
353
354 return tres;
355}
356
357
358template<class Type, template<class> class PatchField, class GeoMesh>
359void 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
379template<class Type, template<class> class PatchField, class GeoMesh>
380tmp
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
408template<class Type, template<class> class PatchField, class GeoMesh>
409tmp
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 \
429template<class Type, template<class> class PatchField, class GeoMesh> \
430dimensioned<ReturnType> Func \
431( \
432 const GeometricField<Type, PatchField, GeoMesh>& f1, \
433 const label comm \
434) \
435{ \
436 return dimensioned<ReturnType> \
437 ( \
438 #Func "(" + f1.name() + ')', \
439 f1.dimensions(), \
440 returnReduce \
441 ( \
442 Foam::Func \
443 ( \
444 Foam::Func(f1.primitiveField()), \
445 Foam::Func(f1.boundaryField()) \
446 ), \
447 BinaryOp<ReturnType>(), \
448 UPstream::msgType(), \
449 comm \
450 ) \
451 ); \
452} \
453 \
454template<class Type, template<class> class PatchField, class GeoMesh> \
455dimensioned<ReturnType> Func \
456( \
457 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1, \
458 const label comm \
459) \
460{ \
461 dimensioned<ReturnType> res = Func(tf1(), comm); \
462 tf1.clear(); \
463 return res; \
464}
465
466
471
472#undef UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY
473
474
475// Forward to DimensionedField directly (same name)
476#define UNARY_REDUCTION_FUNCTION(ReturnType, Func) \
477 \
478template<class Type, template<class> class PatchField, class GeoMesh> \
479dimensioned<ReturnType> Func \
480( \
481 const GeometricField<Type, PatchField, GeoMesh>& f1, \
482 const label comm \
483) \
484{ \
485 return Func(f1.internalField(), comm); \
486} \
487 \
488template<class Type, template<class> class PatchField, class GeoMesh> \
489dimensioned<ReturnType> Func \
490( \
491 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1, \
492 const label comm \
493) \
494{ \
495 auto result = Func(tf1(), comm); \
496 tf1.clear(); \
497 return result; \
498}
499
504#undef UNARY_REDUCTION_FUNCTION
507BINARY_FUNCTION(Type, Type, Type, max)
508BINARY_FUNCTION(Type, Type, Type, min)
509BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
510BINARY_FUNCTION(Type, Type, Type, cmptDivide)
511
512BINARY_TYPE_FUNCTION(Type, Type, Type, max)
513BINARY_TYPE_FUNCTION(Type, Type, Type, min)
514BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
515BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
517
518// ------------------------------------------------------------------------- //
519
520// Clamp Methods
521
522template<class Type, template<class> class PatchField, class GeoMesh>
523void clamp
524(
525 GeometricField<Type, PatchField, GeoMesh>& result,
526 const GeometricField<Type, PatchField, GeoMesh>& f1,
528)
529{
531
532 clamp(result.primitiveFieldRef(), f1.primitiveField(), range);
533 clamp(result.boundaryFieldRef(), f1.boundaryField(), range);
534 result.oriented() = f1.oriented();
535 result.correctLocalBoundaryConditions();
538 result.boundaryField().check();
539 }
540}
541
542template<class Type, template<class> class PatchField, class GeoMesh>
544clamp
545(
547 Foam::zero_one
548)
549{
550 auto tres =
552 (
553 f1,
554 "clamp01(" + f1.name() + ')',
555 f1.dimensions()
556 );
557
558 clamp(tres.ref(), f1, Foam::zero_one{});
560 return tres;
561}
562
563
564template<class Type, template<class> class PatchField, class GeoMesh>
565tmp<GeometricField<Type, PatchField, GeoMesh>>
566clamp
567(
568 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1,
570)
571{
572 const auto& f1 = tf1();
573
574 auto tres =
576 (
577 tf1,
578 "clamp01(" + f1.name() + ')',
579 f1.dimensions()
580 );
582 clamp(tres.ref(), f1, Foam::zero_one{});
583
584 tf1.clear();
585 return tres;
589
590
591// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
593TERNARY_FUNCTION(Type, Type, Type, scalar, lerp)
594TERNARY_TYPE_FUNCTION_FFS(Type, Type, Type, scalar, lerp)
597// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
600
601BINARY_OPERATOR(Type, Type, scalar, *, '*', multiply)
602BINARY_OPERATOR(Type, scalar, Type, *, '*', multiply)
603BINARY_OPERATOR(Type, Type, scalar, /, '|', divide)
604
605BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, '*', multiply)
606BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, '*', multiply)
607
608BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, '|', divide)
609
610
611// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
612
613#define PRODUCT_OPERATOR(product, Op, OpFunc) \
614 \
615template \
616<class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
617void OpFunc \
618( \
619 GeometricField \
620 <typename product<Type1, Type2>::type, PatchField, GeoMesh>& result, \
621 const GeometricField<Type1, PatchField, GeoMesh>& f1, \
622 const GeometricField<Type2, PatchField, GeoMesh>& f2 \
623) \
624{ \
625 Foam::OpFunc \
626 ( \
627 result.primitiveFieldRef(), \
628 f1.primitiveField(), \
629 f2.primitiveField() \
630 ); \
631 Foam::OpFunc \
632 ( \
633 result.boundaryFieldRef(), \
634 f1.boundaryField(), \
635 f2.boundaryField() \
636 ); \
637 \
638 result.oriented() = (f1.oriented() Op f2.oriented()); \
639 result.correctLocalBoundaryConditions(); \
640 if (GeometricBoundaryField<Type1, PatchField, GeoMesh>::debug) \
641 { \
642 result.boundaryField().check(); \
643 } \
644} \
645 \
646 \
647template \
648<class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
649tmp \
650< \
651 GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
652> \
653operator Op \
654( \
655 const GeometricField<Type1, PatchField, GeoMesh>& f1, \
656 const GeometricField<Type2, PatchField, GeoMesh>& f2 \
657) \
658{ \
659 typedef typename product<Type1, Type2>::type resultType; \
660 \
661 auto tres = \
662 reuseTmpGeometricField<resultType, Type1, PatchField, GeoMesh>::New \
663 ( \
664 f1, \
665 '(' + f1.name() + #Op + f2.name() + ')', \
666 (f1.dimensions() Op f2.dimensions()) \
667 ); \
668 \
669 Foam::OpFunc(tres.ref(), f1, f2); \
670 \
671 return tres; \
672} \
673 \
674 \
675template \
676<class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
677tmp \
678< \
679 GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
680> \
681operator Op \
682( \
683 const GeometricField<Type1, PatchField, GeoMesh>& f1, \
684 const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tf2 \
685) \
686{ \
687 typedef typename product<Type1, Type2>::type resultType; \
688 \
689 const auto& f2 = tf2(); \
690 \
691 auto tres = \
692 reuseTmpGeometricField<resultType, Type2, PatchField, GeoMesh>::New \
693 ( \
694 tf2, \
695 '(' + f1.name() + #Op + f2.name() + ')', \
696 (f1.dimensions() Op f2.dimensions()) \
697 ); \
698 \
699 Foam::OpFunc(tres.ref(), f1, f2); \
700 \
701 tf2.clear(); \
702 return tres; \
703} \
704 \
705template \
706<class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
707tmp \
708< \
709 GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
710> \
711operator Op \
712( \
713 const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tf1, \
714 const GeometricField<Type2, PatchField, GeoMesh>& f2 \
715) \
716{ \
717 typedef typename product<Type1, Type2>::type resultType; \
718 \
719 const auto& f1 = tf1(); \
720 \
721 auto tres = \
722 reuseTmpGeometricField<resultType, Type1, PatchField, GeoMesh>::New \
723 ( \
724 tf1, \
725 '(' + f1.name() + #Op + f2.name() + ')', \
726 (f1.dimensions() Op f2.dimensions()) \
727 ); \
728 \
729 Foam::OpFunc(tres.ref(), f1, f2); \
730 \
731 tf1.clear(); \
732 return tres; \
733} \
734 \
735template \
736<class Type1, class Type2, template<class> class PatchField, class GeoMesh> \
737tmp \
738< \
739 GeometricField<typename product<Type1, Type2>::type, PatchField, GeoMesh> \
740> \
741operator Op \
742( \
743 const tmp<GeometricField<Type1, PatchField, GeoMesh>>& tf1, \
744 const tmp<GeometricField<Type2, PatchField, GeoMesh>>& tf2 \
745) \
746{ \
747 typedef typename product<Type1, Type2>::type resultType; \
748 \
749 const auto& f1 = tf1(); \
750 const auto& f2 = tf2(); \
751 \
752 auto tres = \
753 reuseTmpTmpGeometricField \
754 <resultType, Type1, Type1, Type2, PatchField, GeoMesh>::New \
755 ( \
756 tf1, \
757 tf2, \
758 '(' + f1.name() + #Op + f2.name() + ')', \
759 (f1.dimensions() Op f2.dimensions()) \
760 ); \
761 \
762 Foam::OpFunc(tres.ref(), f1, f2); \
763 \
764 tf1.clear(); \
765 tf2.clear(); \
766 return tres; \
767} \
768 \
769template \
770<class Form, class Type, template<class> class PatchField, class GeoMesh> \
771void OpFunc \
772( \
773 GeometricField \
774 <typename product<Type, Form>::type, PatchField, GeoMesh>& result, \
775 const GeometricField<Type, PatchField, GeoMesh>& f1, \
776 const dimensioned<Form>& dvs \
777) \
778{ \
779 Foam::OpFunc(result.primitiveFieldRef(), f1.primitiveField(), dvs.value());\
780 Foam::OpFunc(result.boundaryFieldRef(), f1.boundaryField(), dvs.value()); \
781 result.oriented() = f1.oriented(); \
782 if (GeometricBoundaryField<Type, PatchField, GeoMesh>::debug) \
783 { \
784 result.boundaryField().check(); \
785 } \
786} \
787 \
788template \
789<class Form, class Type, template<class> class PatchField, class GeoMesh> \
790tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
791operator Op \
792( \
793 const GeometricField<Type, PatchField, GeoMesh>& f1, \
794 const dimensioned<Form>& dvs \
795) \
796{ \
797 typedef typename product<Type, Form>::type resultType; \
798 \
799 auto tres = \
800 reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
801 ( \
802 f1, \
803 '(' + f1.name() + #Op + dvs.name() + ')', \
804 (f1.dimensions() Op dvs.dimensions()) \
805 ); \
806 \
807 Foam::OpFunc(tres.ref(), f1, dvs); \
808 \
809 return tres; \
810} \
811 \
812template \
813< \
814 class Form, \
815 class Cmpt, \
816 direction nCmpt, \
817 class Type, template<class> class PatchField, \
818 class GeoMesh \
819> \
820tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
821operator Op \
822( \
823 const GeometricField<Type, PatchField, GeoMesh>& f1, \
824 const VectorSpace<Form,Cmpt,nCmpt>& vs \
825) \
826{ \
827 return f1 Op dimensioned<Form>(static_cast<const Form&>(vs)); \
828} \
829 \
830 \
831template \
832<class Form, class Type, template<class> class PatchField, class GeoMesh> \
833tmp<GeometricField<typename product<Type, Form>::type, PatchField, GeoMesh>> \
834operator Op \
835( \
836 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1, \
837 const dimensioned<Form>& dvs \
838) \
839{ \
840 typedef typename product<Type, Form>::type resultType; \
841 \
842 const auto& f1 = tf1(); \
843 \
844 auto tres = \
845 reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
846 ( \
847 tf1, \
848 '(' + f1.name() + #Op + dvs.name() + ')', \
849 (f1.dimensions() Op dvs.dimensions()) \
850 ); \
851 \
852 Foam::OpFunc(tres.ref(), f1, dvs); \
853 \
854 tf1.clear(); \
855 return tres; \
856} \
857 \
858template \
859< \
860 class Form, \
861 class Cmpt, \
862 direction nCmpt, \
863 class Type, template<class> class PatchField, \
864 class GeoMesh \
865> \
866tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
867operator Op \
868( \
869 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf1, \
870 const VectorSpace<Form,Cmpt,nCmpt>& vs \
871) \
872{ \
873 return tf1 Op dimensioned<Form>(static_cast<const Form&>(vs)); \
874} \
875 \
876 \
877template \
878<class Form, class Type, template<class> class PatchField, class GeoMesh> \
879void OpFunc \
880( \
881 GeometricField \
882 <typename product<Form, Type>::type, PatchField, GeoMesh>& result, \
883 const dimensioned<Form>& dvs, \
884 const GeometricField<Type, PatchField, GeoMesh>& f2 \
885) \
886{ \
887 Foam::OpFunc(result.primitiveFieldRef(), dvs.value(), f2.primitiveField());\
888 Foam::OpFunc(result.boundaryFieldRef(), dvs.value(), f2.boundaryField()); \
889 result.oriented() = f2.oriented(); \
890 if (GeometricBoundaryField<Type, PatchField, GeoMesh>::debug) \
891 { \
892 result.boundaryField().check(); \
893 } \
894} \
895 \
896template \
897<class Form, class Type, template<class> class PatchField, class GeoMesh> \
898tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
899operator Op \
900( \
901 const dimensioned<Form>& dvs, \
902 const GeometricField<Type, PatchField, GeoMesh>& f2 \
903) \
904{ \
905 typedef typename product<Form, Type>::type resultType; \
906 \
907 auto tres = \
908 reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
909 ( \
910 f2, \
911 '(' + dvs.name() + #Op + f2.name() + ')', \
912 (dvs.dimensions() Op f2.dimensions()) \
913 ); \
914 \
915 Foam::OpFunc(tres.ref(), dvs, f2); \
916 \
917 return tres; \
918} \
919 \
920template \
921< \
922 class Form, \
923 class Cmpt, \
924 direction nCmpt, \
925 class Type, template<class> class PatchField, \
926 class GeoMesh \
927> \
928tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
929operator Op \
930( \
931 const VectorSpace<Form,Cmpt,nCmpt>& vs, \
932 const GeometricField<Type, PatchField, GeoMesh>& f2 \
933) \
934{ \
935 return dimensioned<Form>(static_cast<const Form&>(vs)) Op f2; \
936} \
937 \
938template \
939<class Form, class Type, template<class> class PatchField, class GeoMesh> \
940tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
941operator Op \
942( \
943 const dimensioned<Form>& dvs, \
944 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf2 \
945) \
946{ \
947 typedef typename product<Form, Type>::type resultType; \
948 \
949 const auto& f2 = tf2(); \
950 \
951 auto tres = \
952 reuseTmpGeometricField<resultType, Type, PatchField, GeoMesh>::New \
953 ( \
954 tf2, \
955 '(' + dvs.name() + #Op + f2.name() + ')', \
956 (dvs.dimensions() Op f2.dimensions()) \
957 ); \
958 \
959 Foam::OpFunc(tres.ref(), dvs, f2); \
960 \
961 tf2.clear(); \
962 return tres; \
963} \
964 \
965template \
966< \
967 class Form, \
968 class Cmpt, \
969 direction nCmpt, \
970 class Type, template<class> class PatchField, \
971 class GeoMesh \
972> \
973tmp<GeometricField<typename product<Form, Type>::type, PatchField, GeoMesh>> \
974operator Op \
975( \
976 const VectorSpace<Form,Cmpt,nCmpt>& vs, \
977 const tmp<GeometricField<Type, PatchField, GeoMesh>>& tf2 \
978) \
979{ \
980 return dimensioned<Form>(static_cast<const Form&>(vs)) Op tf2; \
986
991
992#undef PRODUCT_OPERATOR
993
994
995// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
996
997} // End namespace Foam
998
999// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1000
1001#include "undefFieldFunctionsM.H"
1002
1003// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#define BINARY_FUNCTION(ReturnType, Type1, Type2, Func)
#define BINARY_TYPE_FUNCTION_FS(ReturnType, Type1, Type2, Func)
#define UNARY_OPERATOR(ReturnType, Type1, Op, OpFunc, Dfunc)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define TERNARY_FUNCTION(ReturnType, Type1, Type2, Type3, Func)
#define TERNARY_TYPE_FUNCTION_FFS(ReturnType, Type1, Type2, Type3, Func)
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BINARY_TYPE_OPERATOR_FS(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BINARY_TYPE_FUNCTION(ReturnType, Type1, Type2, Func)
#define UNARY_REDUCTION_FUNCTION(ReturnType, Func, gFunc)
#define UNARY_REDUCTION_FUNCTION_WITH_BOUNDARY(ReturnType, Func, BinaryOp)
scalar range
const dimensionSet & dimensions() const noexcept
Return dimensions.
orientedType oriented() const noexcept
Return oriented type.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition GeoMesh.H:46
bool check() const
Helper: check if field has been evaluated. See instantiations.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
Field< Type >::cmptType cmptType
Component type of the field elements.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition MinMax.H:106
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
symmTypeOfRank< typenamepTraits< arg1 >::cmptType, arg2 *direction(pTraits< arg1 >::rank)>::type type
Definition products.H:176
A class for managing temporary objects.
Definition tmp.H:75
The magnitude type for given argument.
Definition products.H:93
pTraits< typenamepTraits< arg1 >::cmptType >::magType type
Definition products.H:96
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
const volScalarField & T
#define PRODUCT_OPERATOR(product, op, opFunc)
Namespace for OpenFOAM.
dimensioned< scalarMinMax > minMaxMag(const DimensionedField< Type, GeoMesh > &f1, const label comm)
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)
void subtract(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition hashSets.C:54
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh > > cmptAv(const DimensionedField< Type, GeoMesh > &f1)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
uint8_t direction
Definition direction.H:49
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void dotdot(FieldField< Field1, typename scalarProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
void cross(FieldField< Field1, typename crossProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
void negate(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1)
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
void outer(FieldField< Field1, typename outerProduct< 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)
dimensioned< typename typeOfMag< Type >::type > sumMag(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static tmp< GeometricField< TypeR, PatchField, GeoMesh > > New(const GeometricField< Type1, PatchField, GeoMesh > &f1, const word &name, const dimensionSet &dimensions)
Pass-through to New GeometricField.