Loading...
Searching...
No Matches
GeometricScalarField.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-2017 OpenFOAM Foundation
9 Copyright (C) 2023 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 template<template<class> class PatchField, class GeoMesh>
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
40
41template<template<class> class PatchField, class GeoMesh>
42void stabilise
43(
46 const dimensioned<scalar>& ds
47)
49 stabilise(result.primitiveFieldRef(), gsf.primitiveField(), ds.value());
50 stabilise(result.boundaryFieldRef(), gsf.boundaryField(), ds.value());
51}
52
53
54template<template<class> class PatchField, class GeoMesh>
55tmp<GeometricField<scalar, PatchField, GeoMesh>> stabilise
56(
57 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
58 const dimensioned<scalar>& ds
59)
60{
62 (
63 "stabilise(" + gsf.name() + ',' + ds.name() + ')',
64 gsf.mesh(),
65 ds.dimensions() + gsf.dimensions()
66 );
67
68 stabilise(tRes.ref(), gsf, ds);
69
70 return tRes;
71}
72
73
74template<template<class> class PatchField, class GeoMesh>
75tmp<GeometricField<scalar, PatchField, GeoMesh>> stabilise
76(
77 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
78 const dimensioned<scalar>& ds
79)
80{
81 const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf();
82
83 tmp<GeometricField<scalar, PatchField, GeoMesh>> tRes
84 (
85 New
86 (
87 tgsf,
88 "stabilise(" + gsf.name() + ',' + ds.name() + ')',
89 ds.dimensions() + gsf.dimensions()
90 )
91 );
92
93 stabilise(tRes.ref(), gsf, ds);
94
95 tgsf.clear();
97 return tRes;
98}
101// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
103BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, '+', add)
104BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, '-', subtract)
105
106BINARY_OPERATOR(scalar, scalar, scalar, *, '*', multiply)
107BINARY_OPERATOR(scalar, scalar, scalar, /, '|', divide)
109BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, '|', divide)
110
111
112// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113
114template<template<class> class PatchField, class GeoMesh>
115void pow
116(
117 GeometricField<scalar, PatchField, GeoMesh>& Pow,
118 const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
119 const GeometricField<scalar, PatchField, GeoMesh>& gsf2
120)
121{
122 pow(Pow.primitiveFieldRef(), gsf1.primitiveField(), gsf2.primitiveField());
123 pow(Pow.boundaryFieldRef(), gsf1.boundaryField(), gsf2.boundaryField());
124 Pow.correctLocalBoundaryConditions();
127 Pow.boundaryField().check();
128 }
129}
130
131
132template<template<class> class PatchField, class GeoMesh>
134(
137)
138{
139 if
140 (
142 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
143 )
144 {
146 << "pow() expects dimensionless parameters, but found" << nl;
147
148 if (!f1.dimensions().dimensionless())
149 {
151 << " Base field dimensions: " << f1.dimensions() << nl;
152 }
153 if (!f2.dimensions().dimensionless())
154 {
156 << " Exponent field dimensions: " << f2.dimensions() << nl;
157 }
159 }
160
162 (
163 "pow(" + f1.name() + ',' + f2.name() + ')',
164 f1.mesh(),
165 dimless
166 );
167
168 pow(tresult.ref(), f1, f2);
169
170 return tresult;
171}
172
173
174template<template<class> class PatchField, class GeoMesh>
175tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
176(
177 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf1,
178 const GeometricField<scalar, PatchField, GeoMesh>& f2
179)
180{
181 const auto& f1 = tf1();
182
183 if
184 (
186 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
187 )
188 {
190 << "pow() expects dimensionless parameters, but found" << nl;
191
192 if (!f1.dimensions().dimensionless())
193 {
195 << " Base field dimensions: " << f1.dimensions() << nl;
196 }
197 if (!f2.dimensions().dimensionless())
198 {
200 << " Exponent field dimensions: " << f2.dimensions() << nl;
201 }
203 }
204
206 (
207 New
208 (
209 tf1,
210 "pow(" + f1.name() + ',' + f2.name() + ')',
211 dimless
212 )
213 );
214
215 pow(tresult.ref(), f1, f2);
216
217 tf1.clear();
218
219 return tresult;
220}
221
222
223template<template<class> class PatchField, class GeoMesh>
224tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
225(
226 const GeometricField<scalar, PatchField, GeoMesh>& f1,
227 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf2
228)
229{
230 const auto& f2 = tf2();
231
232 if
233 (
235 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
236 )
237 {
239 << "pow() expects dimensionless parameters, but found" << nl;
240
241 if (!f1.dimensions().dimensionless())
242 {
244 << " Base field dimensions: " << f1.dimensions() << nl;
245 }
246 if (!f2.dimensions().dimensionless())
247 {
249 << " Exponent field dimensions: " << f2.dimensions() << nl;
250 }
252 }
253
255 (
256 New
257 (
258 tf2,
259 "pow(" + f1.name() + ',' + f2.name() + ')',
260 dimless
261 )
262 );
263
264 pow(tresult.ref(), f1, f2);
266 tf2.clear();
267
268 return tresult;
269}
270
271template<template<class> class PatchField, class GeoMesh>
272tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
273(
274 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf1,
275 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf2
276)
277{
278 const auto& f1 = tf1();
279 const auto& f2 = tf2();
280
281 if
282 (
284 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
285 )
286 {
288 << "pow() expects dimensionless parameters, but found" << nl;
289
290 if (!f1.dimensions().dimensionless())
291 {
293 << " Base field dimensions: " << f1.dimensions() << nl;
294 }
295 if (!f2.dimensions().dimensionless())
296 {
298 << " Exponent field dimensions: " << f2.dimensions() << nl;
299 }
301 }
302
304 (
306 <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::New
307 (
308 tf1,
309 tf2,
310 "pow(" + f1.name() + ',' + f2.name() + ')',
311 dimless
312 )
313 );
314
315 pow(tresult.ref(), f1, f2);
316
317 tf1.clear();
318 tf2.clear();
319
320 return tresult;
321}
322
323
324template<template<class> class PatchField, class GeoMesh>
325void pow
326(
327 GeometricField<scalar, PatchField, GeoMesh>& tPow,
328 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
329 const dimensioned<scalar>& ds
330)
331{
332 pow(tPow.primitiveFieldRef(), gsf.primitiveField(), ds.value());
333 pow(tPow.boundaryFieldRef(), gsf.boundaryField(), ds.value());
334 tPow.correctLocalBoundaryConditions();
337 tPow.boundaryField().check();
338 }
339}
340
341
342template<template<class> class PatchField, class GeoMesh>
344(
346 const dimensionedScalar& ds
347)
348{
349 if
350 (
352 && (!ds.dimensions().dimensionless())
353 )
354 {
356 << "pow() expects dimensionless parameters, but found" << nl
357 << " Exponent dimensions: " << ds.dimensions() << nl
358 << exit(FatalError);
359 }
360
362 (
363 "pow(" + f1.name() + ',' + ds.name() + ')',
364 f1.mesh(),
365 pow(f1.dimensions(), ds)
366 );
368 pow(tresult.ref(), f1, ds);
369
370 return tresult;
371}
372
373template<template<class> class PatchField, class GeoMesh>
374tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
375(
376 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf1,
377 const dimensionedScalar& ds
378)
379{
380 const auto& f1 = tf1();
381
382 if
383 (
385 && (!ds.dimensions().dimensionless())
386 )
387 {
389 << "pow() expects dimensionless parameters, but found" << nl
390 << " Exponent dimensions: " << ds.dimensions() << nl
391 << exit(FatalError);
392 }
393
395 (
396 New
397 (
398 tf1,
399 "pow(" + f1.name() + ',' + ds.name() + ')',
400 pow(f1.dimensions(), ds)
401 )
402 );
403
404 pow(tresult.ref(), f1, ds);
406 tf1.clear();
407
408 return tresult;
409}
410
411template<template<class> class PatchField, class GeoMesh>
412tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
413(
414 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
415 const scalar& s
416)
417{
418 return pow(gsf, dimensionedScalar(s));
419}
420
421template<template<class> class PatchField, class GeoMesh>
422tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
423(
424 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
425 const scalar& s
427{
428 return pow(tgsf, dimensionedScalar(s));
429}
430
431
432template<template<class> class PatchField, class GeoMesh>
433void pow
434(
435 GeometricField<scalar, PatchField, GeoMesh>& tPow,
436 const dimensioned<scalar>& ds,
437 const GeometricField<scalar, PatchField, GeoMesh>& gsf
438)
439{
440 pow(tPow.primitiveFieldRef(), ds.value(), gsf.primitiveField());
441 pow(tPow.boundaryFieldRef(), ds.value(), gsf.boundaryField());
442 tPow.correctLocalBoundaryConditions();
445 tPow.boundaryField().check();
446 }
447}
448
449
450template<template<class> class PatchField, class GeoMesh>
452(
453 const dimensionedScalar& ds,
455)
456{
457 if
458 (
460 && (!ds.dimensions().dimensionless() || !f2.dimensions().dimensionless())
461 )
462 {
464 << "pow() expects dimensionless parameters, but found" << nl;
465
466 if (!ds.dimensions().dimensionless())
467 {
469 << " Base scalar dimensions: " << ds.dimensions() << nl;
470 }
471 if (!f2.dimensions().dimensionless())
472 {
474 << " Exponent field dimensions: " << f2.dimensions() << nl;
475 }
477 }
478
480 (
481 "pow(" + ds.name() + ',' + f2.name() + ')',
482 f2.mesh(),
483 dimless
484 );
485
486 pow(tresult.ref(), ds, f2);
487
488 return tresult;
489}
490
491
492template<template<class> class PatchField, class GeoMesh>
493tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
494(
495 const dimensionedScalar& ds,
496 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tf2
497)
498{
499 const auto& f2 = tf2();
500
501 if
502 (
504 && (!ds.dimensions().dimensionless() || !f2.dimensions().dimensionless())
505 )
506 {
508 << "pow() expects dimensionless parameters, but found" << nl;
509
510 if (!ds.dimensions().dimensionless())
511 {
513 << " Base scalar dimensions: " << ds.dimensions() << nl;
514 }
515 if (!f2.dimensions().dimensionless())
516 {
518 << " Exponent field dimensions: " << f2.dimensions() << nl;
519 }
521 }
522
524 (
525 New
526 (
527 tf2,
528 "pow(" + ds.name() + ',' + f2.name() + ')',
529 dimless
530 )
531 );
532
533 pow(tresult.ref(), ds, f2);
535 tf2.clear();
536
537 return tresult;
538}
539
540template<template<class> class PatchField, class GeoMesh>
541tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
542(
543 const scalar& s,
545)
546{
547 return pow(dimensionedScalar(s), gsf);
548}
549
550template<template<class> class PatchField, class GeoMesh>
551tmp<GeometricField<scalar, PatchField, GeoMesh>> pow
552(
553 const scalar& s,
554 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf
555)
556{
557 return pow(dimensionedScalar(s), tgsf);
558}
559
560
561// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
562
563template<template<class> class PatchField, class GeoMesh>
564void atan2
565(
566 GeometricField<scalar, PatchField, GeoMesh>& Atan2,
567 const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
568 const GeometricField<scalar, PatchField, GeoMesh>& gsf2
569)
570{
571 atan2
572 (
573 Atan2.primitiveFieldRef(),
574 gsf1.primitiveField(),
575 gsf2.primitiveField()
576 );
577 atan2
578 (
579 Atan2.boundaryFieldRef(),
580 gsf1.boundaryField(),
581 gsf2.boundaryField()
582 );
583
584 Atan2.correctLocalBoundaryConditions();
585
588 Atan2.boundaryField().check();
589 }
590}
591
592
593template<template<class> class PatchField, class GeoMesh>
595(
598)
599{
601 (
602 "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
603 gsf1.mesh(),
604 atan2(gsf1.dimensions(), gsf2.dimensions())
605 );
606
607 atan2(tAtan2.ref(), gsf1, gsf2);
608
609 return tAtan2;
610}
611
612
613template<template<class> class PatchField, class GeoMesh>
614tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
615(
616 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf1,
617 const GeometricField<scalar, PatchField, GeoMesh>& gsf2
618)
619{
620 const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
621
622 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
623 (
624 New
625 (
626 tgsf1,
627 "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
628 atan2(gsf1.dimensions(), gsf2.dimensions())
629 )
630 );
631
632 atan2(tAtan2.ref(), gsf1, gsf2);
633
634 tgsf1.clear();
635
636 return tAtan2;
637}
638
639
640template<template<class> class PatchField, class GeoMesh>
641tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
642(
643 const GeometricField<scalar, PatchField, GeoMesh>& gsf1,
644 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf2
645)
646{
647 const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
648
649 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
650 (
651 New
652 (
653 tgsf2,
654 "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
655 atan2( gsf1.dimensions(), gsf2.dimensions())
656 )
657 );
658
659 atan2(tAtan2.ref(), gsf1, gsf2);
661 tgsf2.clear();
662
663 return tAtan2;
664}
665
666template<template<class> class PatchField, class GeoMesh>
667tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
668(
669 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf1,
670 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf2
671)
672{
673 const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
674 const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
675
676 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
677 (
678 reuseTmpTmpGeometricField
679 <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::New
680 (
681 tgsf1,
682 tgsf2,
683 "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
684 atan2(gsf1.dimensions(), gsf2.dimensions())
685 )
686 );
687
688 atan2(tAtan2.ref(), gsf1, gsf2);
689
690 tgsf1.clear();
691 tgsf2.clear();
692
693 return tAtan2;
694}
695
696
697template<template<class> class PatchField, class GeoMesh>
698void atan2
699(
700 GeometricField<scalar, PatchField, GeoMesh>& tAtan2,
701 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
702 const dimensioned<scalar>& ds
703)
704{
705 atan2(tAtan2.primitiveFieldRef(), gsf.primitiveField(), ds.value());
706 atan2(tAtan2.boundaryFieldRef(), gsf.boundaryField(), ds.value());
707 tAtan2.correctLocalBoundaryConditions();
710 tAtan2.boundaryField().check();
711 }
712}
713
714
715template<template<class> class PatchField, class GeoMesh>
717(
719 const dimensionedScalar& ds
720)
721{
723 (
724 "atan2(" + gsf.name() + ',' + ds.name() + ')',
725 gsf.mesh(),
726 atan2(gsf.dimensions(), ds)
727 );
729 atan2(tAtan2.ref(), gsf, ds);
730
731 return tAtan2;
732}
733
734template<template<class> class PatchField, class GeoMesh>
735tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
736(
737 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
738 const dimensionedScalar& ds
739)
740{
741 const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf();
742
743 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
744 (
745 New
746 (
747 tgsf,
748 "atan2(" + gsf.name() + ',' + ds.name() + ')',
749 atan2(gsf.dimensions(), ds)
750 )
751 );
752
753 atan2(tAtan2.ref(), gsf, ds);
755 tgsf.clear();
756
757 return tAtan2;
758}
759
760template<template<class> class PatchField, class GeoMesh>
761tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
762(
763 const GeometricField<scalar, PatchField, GeoMesh>& gsf,
764 const scalar& s
765)
766{
767 return atan2(gsf, dimensionedScalar(s));
768}
769
770template<template<class> class PatchField, class GeoMesh>
771tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
772(
773 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf,
774 const scalar& s
776{
777 return atan2(tgsf, dimensionedScalar(s));
778}
779
780
781template<template<class> class PatchField, class GeoMesh>
782void atan2
783(
784 GeometricField<scalar, PatchField, GeoMesh>& tAtan2,
785 const dimensioned<scalar>& ds,
786 const GeometricField<scalar, PatchField, GeoMesh>& gsf
787)
788{
789 atan2(tAtan2.primitiveFieldRef(), ds.value(), gsf.primitiveField());
790 atan2(tAtan2.boundaryFieldRef(), ds.value(), gsf.boundaryField());
791 tAtan2.correctLocalBoundaryConditions();
794 tAtan2.boundaryField().check();
795 }
796}
797
798
799template<template<class> class PatchField, class GeoMesh>
801(
802 const dimensionedScalar& ds,
804)
805{
807 (
808 "atan2(" + ds.name() + ',' + gsf.name() + ')',
809 gsf.mesh(),
810 atan2(ds, gsf.dimensions())
811 );
812
813 atan2(tAtan2.ref(), ds, gsf);
814
815 return tAtan2;
816}
817
818
819template<template<class> class PatchField, class GeoMesh>
820tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
821(
822 const dimensionedScalar& ds,
823 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf
824)
825{
826 const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf();
827
828 tmp<GeometricField<scalar, PatchField, GeoMesh>> tAtan2
829 (
830 New
831 (
832 tgsf,
833 "atan2(" + ds.name() + ',' + gsf.name() + ')',
834 atan2(ds, gsf.dimensions())
835 )
836 );
837
838 atan2(tAtan2.ref(), ds, gsf);
840 tgsf.clear();
841
842 return tAtan2;
843}
844
845template<template<class> class PatchField, class GeoMesh>
846tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
847(
848 const scalar& s,
850)
851{
852 return atan2(dimensionedScalar(s), gsf);
853}
854
855template<template<class> class PatchField, class GeoMesh>
856tmp<GeometricField<scalar, PatchField, GeoMesh>> atan2
857(
858 const scalar& s,
859 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf
860)
862 return atan2(dimensionedScalar(s), tgsf);
866// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
868UNARY_FUNCTION(scalar, scalar, pow3, pow3)
869UNARY_FUNCTION(scalar, scalar, pow4, pow4)
870UNARY_FUNCTION(scalar, scalar, pow5, pow5)
871UNARY_FUNCTION(scalar, scalar, pow6, pow6)
872UNARY_FUNCTION(scalar, scalar, pow025, pow025)
873UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
874UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
875UNARY_FUNCTION(scalar, scalar, sign, sign)
876UNARY_FUNCTION(scalar, scalar, pos, pos)
877UNARY_FUNCTION(scalar, scalar, pos0, pos0)
878UNARY_FUNCTION(scalar, scalar, neg, neg)
879UNARY_FUNCTION(scalar, scalar, neg0, neg0)
883UNARY_FUNCTION(scalar, scalar, exp, trans)
884UNARY_FUNCTION(scalar, scalar, log, trans)
885UNARY_FUNCTION(scalar, scalar, log10, trans)
886UNARY_FUNCTION(scalar, scalar, sin, trans)
887UNARY_FUNCTION(scalar, scalar, cos, trans)
888UNARY_FUNCTION(scalar, scalar, tan, trans)
889UNARY_FUNCTION(scalar, scalar, asin, trans)
890UNARY_FUNCTION(scalar, scalar, acos, trans)
891UNARY_FUNCTION(scalar, scalar, atan, trans)
892UNARY_FUNCTION(scalar, scalar, sinh, trans)
893UNARY_FUNCTION(scalar, scalar, cosh, trans)
894UNARY_FUNCTION(scalar, scalar, tanh, trans)
895UNARY_FUNCTION(scalar, scalar, asinh, trans)
896UNARY_FUNCTION(scalar, scalar, acosh, trans)
897UNARY_FUNCTION(scalar, scalar, atanh, trans)
898UNARY_FUNCTION(scalar, scalar, erf, trans)
899UNARY_FUNCTION(scalar, scalar, erfc, trans)
900UNARY_FUNCTION(scalar, scalar, lgamma, trans)
901UNARY_FUNCTION(scalar, scalar, j0, trans)
902UNARY_FUNCTION(scalar, scalar, j1, trans)
903UNARY_FUNCTION(scalar, scalar, y0, trans)
904UNARY_FUNCTION(scalar, scalar, y1, trans)
905
906
907// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
908
909#define BesselFunc(func) \
910 \
911template<template<class> class PatchField, class GeoMesh> \
912void func \
913( \
914 GeometricField<scalar, PatchField, GeoMesh>& gsf, \
915 const int n, \
916 const GeometricField<scalar, PatchField, GeoMesh>& gsf1 \
917) \
918{ \
919 func(gsf.primitiveFieldRef(), n, gsf1.primitiveField()); \
920 func(gsf.boundaryFieldRef(), n, gsf1.boundaryField()); \
921 gsf.correctLocalBoundaryConditions(); \
922 if (GeometricBoundaryField<scalar, PatchField, GeoMesh>::debug) \
923 { \
924 gsf.boundaryField().check(); \
925 } \
926} \
927 \
928template<template<class> class PatchField, class GeoMesh> \
929tmp<GeometricField<scalar, PatchField, GeoMesh>> func \
930( \
931 const int n, \
932 const GeometricField<scalar, PatchField, GeoMesh>& gsf \
933) \
934{ \
935 if (dimensionSet::checking() && !gsf.dimensions().dimensionless()) \
936 { \
937 FatalErrorInFunction \
938 << "Field is not dimensionless: " << gsf.dimensions() << nl \
939 << abort(FatalError); \
940 } \
941 \
942 auto tFunc = GeometricField<scalar, PatchField, GeoMesh>::New \
943 ( \
944 #func "(" + gsf.name() + ')', \
945 gsf.mesh(), \
946 dimless \
947 ); \
948 \
949 func(tFunc.ref(), n, gsf); \
950 \
951 return tFunc; \
952} \
953 \
954template<template<class> class PatchField, class GeoMesh> \
955tmp<GeometricField<scalar, PatchField, GeoMesh>> func \
956( \
957 const int n, \
958 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf \
959) \
960{ \
961 const auto& gsf = tgsf(); \
962 \
963 if (dimensionSet::checking() && !gsf.dimensions().dimensionless()) \
964 { \
965 FatalErrorInFunction \
966 << "Field is not dimensionless: " << gsf.dimensions() << nl \
967 << abort(FatalError); \
968 } \
969 \
970 tmp<GeometricField<scalar, PatchField, GeoMesh>> tFunc \
971 ( \
972 New \
973 ( \
974 tgsf, \
975 #func "(" + gsf.name() + ')', \
976 dimless \
977 ) \
978 ); \
979 \
980 func(tFunc.ref(), n, gsf); \
982 tgsf.clear(); \
983 \
984 return tFunc; \
985}
986
989
990#undef BesselFunc
991
992
993// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
994
995} // End namespace Foam
996
997// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
998
999#include "undefFieldFunctionsM.H"
1000
1001// ************************************************************************* //
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
#define BesselFunc(func)
Scalar specific part of the implementation of GeometricField.
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
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.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=PatchField< Type >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
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.
void correctLocalBoundaryConditions()
Correct boundary conditions after a purely local operation.
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
static bool checking(bool on) noexcept
Turn dimension checking on/off.
static bool checking() noexcept
True if dimension checking is enabled (the usual default).
bool dimensionless() const
Return true if it is dimensionless.
Generic dimensioned Type class.
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
const word & name() const noexcept
Return const reference to name.
const Type & value() const noexcept
Return const reference to value.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar lgamma(const dimensionedScalar &ds)
void subtract(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar log10(const dimensionedScalar &ds)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
dimensionSet trans(const dimensionSet &ds)
Check the argument is dimensionless (for transcendental functions).
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar atan(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dimensionedScalar asinh(const dimensionedScalar &ds)