FieldFieldFunctions.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) 2019 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
29#include "PstreamReduceOps.H"
31
32#define TEMPLATE template<template<class> class Field, class Type>
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39
40/* * * * * * * * * * * * * * * * Global functions * * * * * * * * * * * * * */
41
42template<template<class> class Field, class Type>
44(
47 const direction d
48)
49{
50 forAll(sf, i)
51 {
52 component(sf[i], f[i], d);
53 }
54}
55
56
57template<template<class> class Field, class Type>
59{
60 forAll(f1, i)
61 {
62 T(f1[i], f2[i]);
63 }
64}
65
66
67template<template<class> class Field, class Type, direction r>
68void pow
69(
72)
73{
74 forAll(f, i)
75 {
76 pow(f[i], vf[i]);
77 }
78}
79
80template<template<class> class Field, class Type, direction r>
81tmp<FieldField<Field, typename powProduct<Type, r>::type>>
83(
85)
86{
87 typedef typename powProduct<Type, r>::type powProductType;
88
89 auto tres
90 (
92 );
93
94 pow<Type, r>(tres.ref(), f);
95 return tres;
96}
97
98template<template<class> class Field, class Type, direction r>
99tmp<FieldField<Field, typename powProduct<Type, r>::type>>
101(
103)
104{
105 typedef typename powProduct<Type, r>::type powProductType;
106
107 auto tres
108 (
110 );
111
112 pow<Type, r>(tres.ref(), tf());
113 tf.clear();
114 return tres;
115}
116
117
118template<template<class> class Field, class Type>
119void sqr
120(
123)
124{
125 forAll(f, i)
126 {
127 sqr(f[i], vf[i]);
128 }
129}
130
131template<template<class> class Field, class Type>
132tmp<FieldField<Field, typename outerProduct<Type, Type>::type>>
134{
135 typedef typename outerProduct<Type, Type>::type outerProductType;
137 (
139 );
140 sqr(tres.ref(), f);
141 return tres;
142}
143
144template<template<class> class Field, class Type>
145tmp<FieldField<Field, typename outerProduct<Type, Type>::type>>
147{
148 typedef typename outerProduct<Type, Type>::type outerProductType;
149
150 auto tres
151 (
153 );
154
155 sqr(tres.ref(), tf());
156 tf.clear();
157 return tres;
158}
159
160
161template<template<class> class Field, class Type>
163(
166)
167{
168 forAll(sf, i)
169 {
170 magSqr(sf[i], f[i]);
171 }
172}
173
174template<template<class> class Field, class Type>
175tmp<FieldField<Field, typename typeOfMag<Type>::type>>
177{
178 typedef typename typeOfMag<Type>::type magType;
179
180 auto tres
181 (
183 );
184
185 magSqr(tres.ref(), f);
186 return tres;
187}
188
189template<template<class> class Field, class Type>
190tmp<FieldField<Field, typename typeOfMag<Type>::type>>
192{
193 typedef typename typeOfMag<Type>::type magType;
194
195 auto tres
196 (
198 );
199
200 magSqr(tres.ref(), tf());
201 tf.clear();
202 return tres;
203}
204
205
206template<template<class> class Field, class Type>
207void mag
208(
211)
212{
213 forAll(sf, i)
214 {
215 mag(sf[i], f[i]);
216 }
217}
218
219template<template<class> class Field, class Type>
220tmp<FieldField<Field, typename typeOfMag<Type>::type>>
222{
223 typedef typename typeOfMag<Type>::type magType;
224
225 auto tres
226 (
228 );
229
230 mag(tres.ref(), f);
231 return tres;
232}
233
234template<template<class> class Field, class Type>
235tmp<FieldField<Field, typename typeOfMag<Type>::type>>
237{
238 typedef typename typeOfMag<Type>::type magType;
239
240 auto tres
241 (
243 );
244
245 mag(tres.ref(), tf());
246 tf.clear();
247 return tres;
248}
249
250
251template<template<class> class Field, class Type>
253(
256)
257{
258 forAll(cf, i)
259 {
260 cmptMax(cf[i], f[i]);
261 }
262}
263
264template<template<class> class Field, class Type>
266(
268)
269{
270 typedef typename FieldField<Field, Type>::cmptType cmptType;
271
272 auto tres
273 (
275 );
276
277 cmptMax(tres.ref(), f);
278 return tres;
279}
280
281template<template<class> class Field, class Type>
283(
285)
286{
287 typedef typename FieldField<Field, Type>::cmptType cmptType;
288
289 auto tres
290 (
292 );
293
294 cmptMax(tres.ref(), tf());
295 tf.clear();
296 return tres;
297}
298
299
300template<template<class> class Field, class Type>
302(
305)
306{
307 forAll(cf, i)
308 {
309 cmptMin(cf[i], f[i]);
310 }
311}
312
313template<template<class> class Field, class Type>
315(
317)
318{
319 typedef typename FieldField<Field, Type>::cmptType cmptType;
320
321 auto tres
322 (
324 );
325
326 cmptMin(tres.ref(), f);
327 return tres;
328}
329
330template<template<class> class Field, class Type>
332(
334)
335{
336 typedef typename FieldField<Field, Type>::cmptType cmptType;
337
338 auto tres
339 (
341 );
342
343 cmptMin(tres.ref(), tf());
344 tf.clear();
345 return tres;
346}
347
348
349template<template<class> class Field, class Type>
351(
354)
355{
356 forAll(cf, i)
357 {
358 cmptAv(cf[i], f[i]);
359 }
360}
361
362template<template<class> class Field, class Type>
364(
366)
367{
368 typedef typename FieldField<Field, Type>::cmptType cmptType;
369
370 auto tres
371 (
373 );
374
375 cmptAv(tres.ref(), f);
376 return tres;
377}
378
379template<template<class> class Field, class Type>
381(
383)
384{
385 typedef typename FieldField<Field, Type>::cmptType cmptType;
386
387 auto tres
388 (
390 );
391
392 cmptAv(tres.ref(), tf());
393 tf.clear();
394 return tres;
395}
396
397
398template<template<class> class Field, class Type>
400(
403)
404{
405 forAll(cf, i)
406 {
407 cmptMag(cf[i], f[i]);
408 }
409}
410
411template<template<class> class Field, class Type>
413(
415)
416{
417 auto tres
418 (
420 );
421
422 cmptMag(tres.ref(), f);
423 return tres;
424}
425
426template<template<class> class Field, class Type>
428(
430)
431{
433 cmptMag(tres.ref(), tf());
434 tf.clear();
435 return tres;
436}
437
438
439#define TMP_UNARY_FUNCTION(returnType, func) \
440 \
441template<template<class> class Field, class Type> \
442returnType func(const tmp<FieldField<Field, Type>>& tf1) \
443{ \
444 returnType res = func(tf1()); \
445 tf1.clear(); \
446 return res; \
447}
448
449template<template<class> class Field, class Type>
451{
452 Type result = pTraits<Type>::min;
453
454 forAll(f, i)
455 {
456 if (f[i].size())
457 {
458 result = max(max(f[i]), result);
459 }
460
461 }
462
463 return result;
464}
465
467
468
469template<template<class> class Field, class Type>
471{
472 Type result = pTraits<Type>::max;
473
474 forAll(f, i)
475 {
476 if (f[i].size())
477 {
478 result = min(min(f[i]), result);
479 }
480 }
481
482 return result;
483}
484
486
487
488template<template<class> class Field, class Type>
490{
491 Type Sum = Zero;
492
493 forAll(f, i)
494 {
495 Sum += sum(f[i]);
496 }
497
498 return Sum;
499}
500
502
503template<template<class> class Field, class Type>
505{
506 typedef typename typeOfMag<Type>::type magType;
507
508 magType result = Zero;
509
510 forAll(f, i)
511 {
512 result += sumMag(f[i]);
513 }
514
515 return result;
516}
517
519
520template<template<class> class Field, class Type>
522{
523 if (f.size())
524 {
525 label n = 0;
526
527 forAll(f, i)
528 {
529 n += f[i].size();
530 }
531
532 if (n)
533 {
534 Type avrg = sum(f)/n;
535
536 return avrg;
537 }
538 }
539
541 << "empty fieldField, returning zero" << endl;
542
543 return Zero;
544}
545
547
548
549template<template<class> class Field, class Type>
551{
552 MinMax<Type> result;
553
554 forAll(f, i)
555 {
556 result += minMax(f[i]);
557 }
558
559 return result;
560}
561
563
564template<template<class> class Field, class Type>
566{
567 scalarMinMax result;
568
569 forAll(f, i)
570 {
571 result += minMaxMag(f[i]);
572 }
573
574 return result;
575}
576
578
579
580// With reduction on ReturnType
581#define G_UNARY_FUNCTION(ReturnType, gFunc, func, rFunc) \
582 \
583template<template<class> class Field, class Type> \
584ReturnType gFunc(const FieldField<Field, Type>& f) \
585{ \
586 ReturnType res = func(f); \
587 reduce(res, rFunc##Op<ReturnType>()); \
588 return res; \
589} \
590TMP_UNARY_FUNCTION(ReturnType, gFunc)
591
597
599
600#undef G_UNARY_FUNCTION
601
602
603template<template<class> class Field, class Type>
605{
606 label n = 0;
607
608 forAll(f, i)
609 {
610 n += f[i].size();
611 }
612
614
615 if (n)
616 {
617 Type avrg = gSum(f)/n;
618
619 return avrg;
620 }
621
623 << "empty fieldField, returning zero" << endl;
624
625 return Zero;
626}
627
629
630#undef TMP_UNARY_FUNCTION
631
632
633BINARY_FUNCTION(Type, Type, Type, max)
634BINARY_FUNCTION(Type, Type, Type, min)
635BINARY_FUNCTION(Type, Type, Type, cmptMultiply)
636BINARY_FUNCTION(Type, Type, Type, cmptDivide)
637
638BINARY_TYPE_FUNCTION(Type, Type, Type, max)
639BINARY_TYPE_FUNCTION(Type, Type, Type, min)
640BINARY_TYPE_FUNCTION(Type, Type, Type, cmptMultiply)
641BINARY_TYPE_FUNCTION(Type, Type, Type, cmptDivide)
642
644
645
646/* * * * * * * * * * * * * * * * Global operators * * * * * * * * * * * * * */
647
648UNARY_OPERATOR(Type, Type, -, negate)
649
650BINARY_OPERATOR(Type, Type, scalar, *, multiply)
651BINARY_OPERATOR(Type, scalar, Type, *, multiply)
652BINARY_OPERATOR(Type, Type, scalar, /, divide)
653
654BINARY_TYPE_OPERATOR_SF(Type, scalar, Type, *, multiply)
655BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, *, multiply)
656
657BINARY_TYPE_OPERATOR_FS(Type, Type, scalar, /, divide)
658
659
660// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
661
662#define PRODUCT_OPERATOR(product, op, opFunc) \
663 \
664template \
665< \
666 template<class> class Field1, \
667 template<class> class Field2, \
668 class Type1, \
669 class Type2 \
670> \
671void opFunc \
672( \
673 FieldField<Field1, typename product<Type1, Type2>::type>& f, \
674 const FieldField<Field1, Type1>& f1, \
675 const FieldField<Field2, Type2>& f2 \
676) \
677{ \
678 forAll(f, i) \
679 { \
680 opFunc(f[i], f1[i], f2[i]); \
681 } \
682} \
683 \
684template \
685< \
686 template<class> class Field1, \
687 template<class> class Field2, \
688 class Type1, \
689 class Type2 \
690> \
691tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
692operator op \
693( \
694 const FieldField<Field1, Type1>& f1, \
695 const FieldField<Field2, Type2>& f2 \
696) \
697{ \
698 typedef typename product<Type1, Type2>::type productType; \
699 tmp<FieldField<Field1, productType>> tres \
700 ( \
701 FieldField<Field1, productType>::NewCalculatedType(f1) \
702 ); \
703 opFunc(tres.ref(), f1, f2); \
704 return tres; \
705} \
706 \
707template<template<class> class Field, class Type1, class Type2> \
708tmp<FieldField<Field, typename product<Type1, Type2>::type>> \
709operator op \
710( \
711 const FieldField<Field, Type1>& f1, \
712 const tmp<FieldField<Field, Type2>>& tf2 \
713) \
714{ \
715 typedef typename product<Type1, Type2>::type productType; \
716 tmp<FieldField<Field, productType>> tres \
717 ( \
718 reuseTmpFieldField<Field, productType, Type2>::New(tf2) \
719 ); \
720 opFunc(tres.ref(), f1, tf2()); \
721 tf2.clear(); \
722 return tres; \
723} \
724 \
725template \
726< \
727 template<class> class Field1, \
728 template<class> class Field2, \
729 class Type1, \
730 class Type2 \
731> \
732tmp<FieldField<Field, typename product<Type1, Type2>::type>> \
733operator op \
734( \
735 const FieldField<Field1, Type1>& f1, \
736 const tmp<FieldField<Field2, Type2>>& tf2 \
737) \
738{ \
739 typedef typename product<Type1, Type2>::type productType; \
740 tmp<FieldField<Field1, productType>> tres \
741 ( \
742 FieldField<Field1, productType>::NewCalculatedType(f1) \
743 ); \
744 opFunc(tres.ref(), f1, tf2()); \
745 tf2.clear(); \
746 return tres; \
747} \
748 \
749template \
750< \
751 template<class> class Field1, \
752 template<class> class Field2, \
753 class Type1, \
754 class Type2 \
755> \
756tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
757operator op \
758( \
759 const tmp<FieldField<Field1, Type1>>& tf1, \
760 const FieldField<Field2, Type2>& f2 \
761) \
762{ \
763 typedef typename product<Type1, Type2>::type productType; \
764 tmp<FieldField<Field1, productType>> tres \
765 ( \
766 reuseTmpFieldField<Field1, productType, Type1>::New(tf1) \
767 ); \
768 opFunc(tres.ref(), tf1(), f2); \
769 tf1.clear(); \
770 return tres; \
771} \
772 \
773template \
774< \
775 template<class> class Field1, \
776 template<class> class Field2, \
777 class Type1, \
778 class Type2 \
779> \
780tmp<FieldField<Field1, typename product<Type1, Type2>::type>> \
781operator op \
782( \
783 const tmp<FieldField<Field1, Type1>>& tf1, \
784 const tmp<FieldField<Field2, Type2>>& tf2 \
785) \
786{ \
787 typedef typename product<Type1, Type2>::type productType; \
788 tmp<FieldField<Field1, productType>> tres \
789 ( \
790 reuseTmpTmpFieldField<Field1, productType, Type1, Type1, Type2>::New \
791 (tf1, tf2) \
792 ); \
793 opFunc(tres.ref(), tf1(), tf2()); \
794 tf1.clear(); \
795 tf2.clear(); \
796 return tres; \
797} \
798 \
799template \
800< \
801 template<class> class Field, \
802 class Type, \
803 class Form, \
804 class Cmpt, \
805 direction nCmpt \
806> \
807void opFunc \
808( \
809 FieldField<Field, typename product<Type, Form>::type>& f, \
810 const FieldField<Field, Type>& f1, \
811 const VectorSpace<Form,Cmpt,nCmpt>& vs \
812) \
813{ \
814 forAll(f, i) \
815 { \
816 opFunc(f[i], f1[i], vs); \
817 } \
818} \
819 \
820template \
821< \
822 template<class> class Field, \
823 class Type, \
824 class Form, \
825 class Cmpt, \
826 direction nCmpt \
827> \
828tmp<FieldField<Field, typename product<Type, Form>::type>> \
829operator op \
830( \
831 const FieldField<Field, Type>& f1, \
832 const VectorSpace<Form,Cmpt,nCmpt>& vs \
833) \
834{ \
835 typedef typename product<Type, Form>::type productType; \
836 tmp<FieldField<Field, productType>> tres \
837 ( \
838 FieldField<Field, productType>::NewCalculatedType(f1) \
839 ); \
840 opFunc(tres.ref(), f1, static_cast<const Form&>(vs)); \
841 return tres; \
842} \
843 \
844template \
845< \
846 template<class> class Field, \
847 class Type, \
848 class Form, \
849 class Cmpt, \
850 direction nCmpt \
851> \
852tmp<FieldField<Field, typename product<Type, Form>::type>> \
853operator op \
854( \
855 const tmp<FieldField<Field, Type>>& tf1, \
856 const VectorSpace<Form,Cmpt,nCmpt>& vs \
857) \
858{ \
859 typedef typename product<Type, Form>::type productType; \
860 tmp<FieldField<Field, productType>> tres \
861 ( \
862 reuseTmpFieldField<Field, productType, Type>::New(tf1) \
863 ); \
864 opFunc(tres.ref(), tf1(), static_cast<const Form&>(vs)); \
865 tf1.clear(); \
866 return tres; \
867} \
868 \
869template \
870< \
871 template<class> class Field, \
872 class Type, \
873 class Form, \
874 class Cmpt, \
875 direction nCmpt \
876> \
877void opFunc \
878( \
879 FieldField<Field, typename product<Form, Type>::type>& f, \
880 const VectorSpace<Form,Cmpt,nCmpt>& vs, \
881 const FieldField<Field, Type>& f1 \
882) \
883{ \
884 forAll(f, i) \
885 { \
886 opFunc(f[i], vs, f1[i]); \
887 } \
888} \
889 \
890template \
891< \
892 template<class> class Field, \
893 class Type, \
894 class Form, \
895 class Cmpt, \
896 direction nCmpt \
897> \
898tmp<FieldField<Field, typename product<Form, Type>::type>> \
899operator op \
900( \
901 const VectorSpace<Form,Cmpt,nCmpt>& vs, \
902 const FieldField<Field, Type>& f1 \
903) \
904{ \
905 typedef typename product<Form, Type>::type productType; \
906 tmp<FieldField<Field, productType>> tres \
907 ( \
908 FieldField<Field, productType>::NewCalculatedType(f1) \
909 ); \
910 opFunc(tres.ref(), static_cast<const Form&>(vs), f1); \
911 return tres; \
912} \
913 \
914template \
915< \
916 template<class> class Field, \
917 class Type, \
918 class Form, \
919 class Cmpt, \
920 direction nCmpt \
921> \
922tmp<FieldField<Field, typename product<Form, Type>::type>> \
923operator op \
924( \
925 const VectorSpace<Form,Cmpt,nCmpt>& vs, \
926 const tmp<FieldField<Field, Type>>& tf1 \
927) \
928{ \
929 typedef typename product<Form, Type>::type productType; \
930 tmp<FieldField<Field, productType>> tres \
931 ( \
932 reuseTmpFieldField<Field, productType, Type>::New(tf1) \
933 ); \
934 opFunc(tres.ref(), static_cast<const Form&>(vs), tf1()); \
935 tf1.clear(); \
936 return tres; \
937}
938
941
946
947#undef PRODUCT_OPERATOR
948
949
950// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
951
952} // End namespace Foam
953
954// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
955
956#include "undefFieldFunctionsM.H"
957
958// ************************************************************************* //
#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 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 G_UNARY_FUNCTION(ReturnType, gFunc, func, rFunc)
#define TMP_UNARY_FUNCTION(returnType, func)
Inter-processor communication reduction functions.
label n
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
pTraits< Type >::cmptType cmptType
Component type.
Definition: FieldField.H:85
Generic templated field type.
Definition: Field.H:82
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition: MinMax.H:128
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
The magnitude type for given argument.
Definition: products.H:89
pTraits< typenamepTraits< arg1 >::cmptType >::magType type
Definition: products.H:92
type
Volume classification types.
Definition: volumeType.H:66
const volScalarField & T
#define PRODUCT_OPERATOR(product, op, opFunc)
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionSet clip(const dimensionSet &ds1, const dimensionSet &ds2)
Definition: dimensionSet.C:278
Type gSum(const FieldField< Field, Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
scalarMinMax gMinMaxMag(const FieldField< Field, Type > &f)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:61
dimensioned< typename typeOfMag< Type >::type > sumMag(const DimensionedField< Type, GeoMesh > &df)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
void negate(FieldField< Field, Type > &res, const FieldField< Field, Type > &f)
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
uint8_t direction
Definition: direction.H:56
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Type gAverage(const FieldField< Field, Type > &f)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
typeOfMag< Type >::type gSumMag(const FieldField< Field, Type > &f)
void dotdot(FieldField< Field1, typename scalarProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Type gMin(const FieldField< Field, Type > &f)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
void cross(FieldField< Field1, typename crossProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh > > cmptAv(const DimensionedField< Type, GeoMesh > &df)
dimensioned< scalarMinMax > minMaxMag(const DimensionedField< Type, GeoMesh > &df)
Type gMax(const FieldField< Field, Type > &f)
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
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)
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333