Loading...
Searching...
No Matches
DimensionedScalarField.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) 2019-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<class GeoMesh>
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
40
41template<class GeoMesh>
43(
45 const dimensioned<scalar>& ds
46)
47{
48 auto tres =
50 (
51 "stabilise(" + dsf.name() + ',' + ds.name() + ')',
52 dsf.mesh(),
53 dsf.dimensions() + ds.dimensions()
54 );
55
56 stabilise(tres.ref().field(), dsf.field(), ds.value());
57
58 return tres;
59}
60
61
62template<class GeoMesh>
63tmp<DimensionedField<scalar, GeoMesh>> stabilise
64(
65 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf,
66 const dimensioned<scalar>& ds
67)
68{
69 const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
70
71 tmp<DimensionedField<scalar, GeoMesh>> tres = New
72 (
73 tdsf,
74 "stabilise(" + dsf.name() + ',' + ds.name() + ')',
75 dsf.dimensions() + ds.dimensions()
76 );
77
78 stabilise(tres.ref().field(), dsf.field(), ds.value());
79
80 tdsf.clear();
82 return tres;
83}
86// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, '+', add)
89BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, '-', subtract)
90
91BINARY_OPERATOR(scalar, scalar, scalar, *, '*', multiply)
92BINARY_OPERATOR(scalar, scalar, scalar, /, '|', divide)
93
94BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, '|', divide)
95
96// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97
98template<class GeoMesh>
100(
101 const DimensionedField<scalar, GeoMesh>& f1,
102 const DimensionedField<scalar, GeoMesh>& f2
103)
104{
105 if
106 (
108 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
109 )
110 {
111 FatalErrorInFunction
112 << "pow() expects dimensionless parameters, but found" << nl;
113
114 if (!f1.dimensions().dimensionless())
115 {
116 FatalError
117 << " Base field dimensions: " << f1.dimensions() << nl;
118 }
119 if (!f2.dimensions().dimensionless())
120 {
121 FatalError
122 << " Exponent field dimensions: " << f2.dimensions() << nl;
123 }
125 }
126
127 auto tresult =
129 (
130 "pow(" + f1.name() + ',' + f2.name() + ')',
131 f1.mesh(),
132 dimless
133 );
134
135 pow(tresult.ref().field(), f1.field(), f2.field());
136
137 return tresult;
138}
139
140
141template<class GeoMesh>
142tmp<DimensionedField<scalar, GeoMesh>> pow
143(
144 const tmp<DimensionedField<scalar, GeoMesh>>& tf1,
145 const DimensionedField<scalar, GeoMesh>& f2
146)
147{
148 const auto& f1 = tf1();
149
150 if
151 (
152 dimensionSet::checking()
153 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
154 )
155 {
157 << "pow() expects dimensionless parameters, but found" << nl;
158
159 if (!f1.dimensions().dimensionless())
160 {
161 FatalError
162 << " Base field dimensions: " << f1.dimensions() << nl;
163 }
164 if (!f2.dimensions().dimensionless())
165 {
167 << " Exponent field dimensions: " << f2.dimensions() << nl;
168 }
169 FatalError << exit(FatalError);
170 }
171
172 tmp<DimensionedField<scalar, GeoMesh>> tresult = New
173 (
174 tf1,
175 "pow(" + f1.name() + ',' + f2.name() + ')',
176 dimless
177 );
178
179 pow(tresult.ref().field(), f1.field(), f2.field());
180
181 tf1.clear();
182
183 return tresult;
184}
185
186
187template<class GeoMesh>
188tmp<DimensionedField<scalar, GeoMesh>> pow
189(
190 const DimensionedField<scalar, GeoMesh>& f1,
191 const tmp<DimensionedField<scalar, GeoMesh>>& tf2
192)
193{
194 const auto& f2 = tf2();
195
196 if
197 (
198 dimensionSet::checking()
199 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
200 )
201 {
203 << "pow() expects dimensionless parameters, but found" << nl;
204
205 if (!f1.dimensions().dimensionless())
206 {
207 FatalError
208 << " Base field dimensions: " << f1.dimensions() << nl;
209 }
210 if (!f2.dimensions().dimensionless())
211 {
212 FatalError
213 << " Exponent field dimensions: " << f2.dimensions() << nl;
214 }
215 FatalError << exit(FatalError);
216 }
217
218 tmp<DimensionedField<scalar, GeoMesh>> tresult = New
219 (
220 tf2,
221 "pow(" + f1.name() + ',' + f2.name() + ')',
222 dimless
223 );
224
225 pow(tresult.ref().field(), f1.field(), f2.field());
226
227 tf2.clear();
228
229 return tresult;
230}
231
232
233template<class GeoMesh>
234tmp<DimensionedField<scalar, GeoMesh>> pow
235(
236 const tmp<DimensionedField<scalar, GeoMesh>>& tf1,
237 const tmp<DimensionedField<scalar, GeoMesh>>& tf2
238)
239{
240 const auto& f1 = tf1();
241 const auto& f2 = tf2();
242
243 if
244 (
245 dimensionSet::checking()
246 && (!f1.dimensions().dimensionless() || !f2.dimensions().dimensionless())
247 )
248 {
250 << "pow() expects dimensionless parameters, but found" << nl;
251
252 if (!f1.dimensions().dimensionless())
253 {
254 FatalError
255 << " Base field dimensions: " << f1.dimensions() << nl;
256 }
257 if (!f2.dimensions().dimensionless())
258 {
259 FatalError
260 << " Exponent field dimensions: " << f2.dimensions() << nl;
261 }
262 FatalError << exit(FatalError);
263 }
264
265 auto tresult =
266 reuseTmpTmpDimensionedField<scalar, scalar, scalar, scalar, GeoMesh>::
267 New
268 (
269 tf1,
270 tf2,
271 "pow(" + f1.name() + ',' + f2.name() + ')',
272 dimless
273 );
274
275 pow(tresult.ref().field(), f1.field(), f2.field());
276
277 tf1.clear();
278 tf2.clear();
279
280 return tresult;
281}
282
283
284template<class GeoMesh>
285tmp<DimensionedField<scalar, GeoMesh>> pow
286(
287 const DimensionedField<scalar, GeoMesh>& f1,
288 const dimensionedScalar& ds
289)
290{
291 if
292 (
293 dimensionSet::checking()
294 && (!ds.dimensions().dimensionless())
295 )
296 {
298 << "pow() expects dimensionless parameters, but found" << nl
299 << " Exponent dimensions: " << ds.dimensions() << nl
300 << exit(FatalError);
301 }
302
303 auto tresult =
304 DimensionedField<scalar, GeoMesh>::New
305 (
306 "pow(" + f1.name() + ',' + ds.name() + ')',
307 f1.mesh(),
308 pow(f1.dimensions(), ds)
309 );
310
311 pow(tresult.ref().field(), f1.field(), ds.value());
312
313 return tresult;
314}
315
316
317template<class GeoMesh>
318tmp<DimensionedField<scalar, GeoMesh>> pow
319(
320 const tmp<DimensionedField<scalar, GeoMesh>>& tf1,
321 const dimensionedScalar& ds
322)
323{
324 const auto& f1 = tf1();
325
326 if
327 (
328 dimensionSet::checking()
329 && (!ds.dimensions().dimensionless())
330 )
331 {
333 << "pow() expects dimensionless parameters, but found" << nl
334 << " Exponent dimensions: " << ds.dimensions() << nl
335 << exit(FatalError);
336 }
337
338 tmp<DimensionedField<scalar, GeoMesh>> tresult = New
339 (
340 tf1,
341 "pow(" + f1.name() + ',' + ds.name() + ')',
342 pow(f1.dimensions(), ds)
343 );
344
345 pow(tresult.ref().field(), f1.field(), ds.value());
346
347 tf1.clear();
348
349 return tresult;
350}
351
352
353template<class GeoMesh>
354tmp<DimensionedField<scalar, GeoMesh>> pow
355(
356 const DimensionedField<scalar, GeoMesh>& dsf,
357 const scalar& s
359{
360 return pow(dsf, dimensionedScalar(s));
361}
362
363
364template<class GeoMesh>
365tmp<DimensionedField<scalar, GeoMesh>> pow
366(
367 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf,
368 const scalar& s
370{
371 return pow(tdsf, dimensionedScalar(s));
372}
373
374
375template<class GeoMesh>
376tmp<DimensionedField<scalar, GeoMesh>> pow
377(
378 const dimensionedScalar& ds,
379 const DimensionedField<scalar, GeoMesh>& f2
380)
381{
382 if
383 (
384 dimensionSet::checking()
385 && (!ds.dimensions().dimensionless() || !f2.dimensions().dimensionless())
386 )
387 {
389 << "pow() expects dimensionless parameters, but found" << nl;
390
391 if (!ds.dimensions().dimensionless())
392 {
393 FatalError
394 << " Base scalar dimensions: " << ds.dimensions() << nl;
395 }
396 if (!f2.dimensions().dimensionless())
397 {
399 << " Exponent field dimensions: " << f2.dimensions() << nl;
400 }
401 FatalError << exit(FatalError);
402 }
403
404 auto tresult =
405 DimensionedField<scalar, GeoMesh>::New
406 (
407 "pow(" + ds.name() + ',' + f2.name() + ')',
408 f2.mesh(),
409 dimless
410 );
411
412 pow(tresult.ref().field(), ds.value(), f2.field());
413
414 return tresult;
415}
416
417
418template<class GeoMesh>
419tmp<DimensionedField<scalar, GeoMesh>> pow
420(
421 const dimensionedScalar& ds,
422 const tmp<DimensionedField<scalar, GeoMesh>>& tf2
423)
424{
425 const auto& f2 = tf2();
426
427 if
428 (
429 dimensionSet::checking()
430 && (!ds.dimensions().dimensionless() || !f2.dimensions().dimensionless())
431 )
432 {
434 << "pow() expects dimensionless parameters, but found" << nl;
435
436 if (!ds.dimensions().dimensionless())
437 {
438 FatalError
439 << " Base scalar dimensions: " << ds.dimensions() << nl;
440 }
441 if (!f2.dimensions().dimensionless())
442 {
444 << " Exponent field dimensions: " << f2.dimensions() << nl;
445 }
446
447 FatalError << exit(FatalError);
448 }
449
450 tmp<DimensionedField<scalar, GeoMesh>> tresult = New
451 (
452 tf2,
453 "pow(" + ds.name() + ',' + f2.name() + ')',
454 dimless
455 );
456
457 pow(tresult.ref().field(), ds.value(), f2.field());
459 tf2.clear();
460
461 return tresult;
462}
463
464template<class GeoMesh>
465tmp<DimensionedField<scalar, GeoMesh>> pow
466(
467 const scalar& s,
469)
470{
471 return pow(dimensionedScalar(s), dsf);
472}
473
474template<class GeoMesh>
475tmp<DimensionedField<scalar, GeoMesh>> pow
476(
477 const scalar& s,
478 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf
479)
480{
481 return pow(dimensionedScalar(s), tdsf);
482}
483
484
485// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
486
487template<class GeoMesh>
488tmp<DimensionedField<scalar, GeoMesh>> atan2
489(
490 const DimensionedField<scalar, GeoMesh>& dsf1,
491 const DimensionedField<scalar, GeoMesh>& dsf2
492)
493{
494 auto tres =
495 DimensionedField<scalar, GeoMesh>::New
496 (
497 "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
498 dsf1.mesh(),
499 atan2(dsf1.dimensions(), dsf2.dimensions())
500 );
501
502 atan2(tres.ref().field(), dsf1.field(), dsf2.field());
503
504 return tres;
505}
506
507
508template<class GeoMesh>
509tmp<DimensionedField<scalar, GeoMesh>> atan2
510(
511 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf1,
512 const DimensionedField<scalar, GeoMesh>& dsf2
513)
514{
515 const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
516
517 tmp<DimensionedField<scalar, GeoMesh>> tres = New
518 (
519 tdsf1,
520 "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
521 atan2(dsf1.dimensions(), dsf2.dimensions())
522 );
523
524 atan2(tres.ref().field(), dsf1.field(), dsf2.field());
525
526 tdsf1.clear();
527
528 return tres;
529}
530
531
532template<class GeoMesh>
533tmp<DimensionedField<scalar, GeoMesh>> atan2
534(
535 const DimensionedField<scalar, GeoMesh>& dsf1,
536 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf2
537)
538{
539 const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
540
541 tmp<DimensionedField<scalar, GeoMesh>> tres = New
542 (
543 tdsf2,
544 "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
545 atan2(dsf1.dimensions(), dsf2.dimensions())
546 );
547
548 atan2(tres.ref().field(), dsf1.field(), dsf2.field());
550 tdsf2.clear();
551
552 return tres;
553}
554
555template<class GeoMesh>
556tmp<DimensionedField<scalar, GeoMesh>> atan2
557(
558 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf1,
559 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf2
560)
561{
562 const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
563 const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
564
565 auto tres =
566 reuseTmpTmpDimensionedField<scalar, scalar, scalar, scalar, GeoMesh>::
567 New
568 (
569 tdsf1,
570 tdsf2,
571 "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
572 atan2(dsf1.dimensions(), dsf2.dimensions())
573 );
574
575 atan2(tres.ref().field(), dsf1.field(), dsf2.field());
576
577 tdsf1.clear();
578 tdsf2.clear();
579
580 return tres;
581}
582
583
584template<class GeoMesh>
585tmp<DimensionedField<scalar, GeoMesh>> atan2
586(
587 const DimensionedField<scalar, GeoMesh>& dsf,
588 const dimensionedScalar& ds
589)
590{
591 auto tres =
592 DimensionedField<scalar, GeoMesh>::New
593 (
594 "atan2(" + dsf.name() + ',' + ds.name() + ')',
595 dsf.mesh(),
596 atan2(dsf.dimensions(), ds)
597 );
599 atan2(tres.ref().field(), dsf.field(), ds.value());
600
601 return tres;
602}
603
604template<class GeoMesh>
605tmp<DimensionedField<scalar, GeoMesh>> atan2
606(
607 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf,
608 const dimensionedScalar& ds
609)
610{
611 const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
612
613 tmp<DimensionedField<scalar, GeoMesh>> tres = New
614 (
615 tdsf,
616 "atan2(" + dsf.name() + ',' + ds.name() + ')',
617 atan2(dsf.dimensions(), ds)
618 );
619
620 atan2(tres.ref().field(), dsf.field(), ds.value());
622 tdsf.clear();
623
624 return tres;
625}
626
627template<class GeoMesh>
628tmp<DimensionedField<scalar, GeoMesh>> atan2
629(
630 const DimensionedField<scalar, GeoMesh>& dsf,
631 const scalar& s
632)
633{
634 return atan2(dsf, dimensionedScalar(s));
635}
636
637template<class GeoMesh>
638tmp<DimensionedField<scalar, GeoMesh>> atan2
639(
640 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf,
641 const scalar& s
643{
644 return atan2(tdsf, dimensionedScalar(s));
645}
646
647
648template<class GeoMesh>
649tmp<DimensionedField<scalar, GeoMesh>> atan2
650(
651 const dimensionedScalar& ds,
652 const DimensionedField<scalar, GeoMesh>& dsf
653)
654{
655 auto tres =
656 DimensionedField<scalar, GeoMesh>::New
657 (
658 "atan2(" + ds.name() + ',' + dsf.name() + ')',
659 dsf.mesh(),
660 atan2(ds, dsf.dimensions())
661 );
662
663 atan2(tres.ref().field(), ds.value(), dsf.field());
664
665 return tres;
666}
667
668
669template<class GeoMesh>
670tmp<DimensionedField<scalar, GeoMesh>> atan2
671(
672 const dimensionedScalar& ds,
673 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf
674)
675{
676 const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
677
678 tmp<DimensionedField<scalar, GeoMesh>> tres = New
679 (
680 tdsf,
681 "atan2(" + ds.name() + ',' + dsf.name() + ')',
682 atan2(ds, dsf.dimensions())
683 );
684
685 atan2(tres.ref().field(), ds.value(), dsf.field());
687 tdsf.clear();
688
689 return tres;
690}
691
692template<class GeoMesh>
693tmp<DimensionedField<scalar, GeoMesh>> atan2
694(
695 const scalar& s,
697)
698{
699 return atan2(dimensionedScalar(s), dsf);
700}
701
702template<class GeoMesh>
703tmp<DimensionedField<scalar, GeoMesh>> atan2
704(
705 const scalar& s,
706 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf
707)
709 return atan2(dimensionedScalar(s), tdsf);
713// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
715UNARY_FUNCTION(scalar, scalar, pow3, pow3)
716UNARY_FUNCTION(scalar, scalar, pow4, pow4)
717UNARY_FUNCTION(scalar, scalar, pow5, pow5)
718UNARY_FUNCTION(scalar, scalar, pow6, pow6)
719UNARY_FUNCTION(scalar, scalar, pow025, pow025)
720UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
721UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
722UNARY_FUNCTION(scalar, scalar, sign, sign)
723UNARY_FUNCTION(scalar, scalar, pos, pos)
724UNARY_FUNCTION(scalar, scalar, pos0, pos0)
725UNARY_FUNCTION(scalar, scalar, neg, neg)
726UNARY_FUNCTION(scalar, scalar, neg0, neg0)
730UNARY_FUNCTION(scalar, scalar, exp, trans)
731UNARY_FUNCTION(scalar, scalar, log, trans)
732UNARY_FUNCTION(scalar, scalar, log10, trans)
733UNARY_FUNCTION(scalar, scalar, sin, trans)
734UNARY_FUNCTION(scalar, scalar, cos, trans)
735UNARY_FUNCTION(scalar, scalar, tan, trans)
736UNARY_FUNCTION(scalar, scalar, asin, trans)
737UNARY_FUNCTION(scalar, scalar, acos, trans)
738UNARY_FUNCTION(scalar, scalar, atan, trans)
739UNARY_FUNCTION(scalar, scalar, sinh, trans)
740UNARY_FUNCTION(scalar, scalar, cosh, trans)
741UNARY_FUNCTION(scalar, scalar, tanh, trans)
742UNARY_FUNCTION(scalar, scalar, asinh, trans)
743UNARY_FUNCTION(scalar, scalar, acosh, trans)
744UNARY_FUNCTION(scalar, scalar, atanh, trans)
745UNARY_FUNCTION(scalar, scalar, erf, trans)
746UNARY_FUNCTION(scalar, scalar, erfc, trans)
747UNARY_FUNCTION(scalar, scalar, lgamma, trans)
748UNARY_FUNCTION(scalar, scalar, j0, trans)
749UNARY_FUNCTION(scalar, scalar, j1, trans)
750UNARY_FUNCTION(scalar, scalar, y0, trans)
751UNARY_FUNCTION(scalar, scalar, y1, trans)
752
753
754// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
755
756#define BesselFunc(func) \
757 \
758template<class GeoMesh> \
759tmp<DimensionedField<scalar, GeoMesh>> func \
760( \
761 const int n, \
762 const DimensionedField<scalar, GeoMesh>& dsf \
763) \
764{ \
765 if (dimensionSet::checking() && !dsf.dimensions().dimensionless()) \
766 { \
767 FatalErrorInFunction \
768 << "Field is not dimensionless: " << dsf.dimensions() << nl \
769 << abort(FatalError); \
770 } \
771 \
772 auto tres = \
773 DimensionedField<scalar, GeoMesh>::New \
774 ( \
775 #func "(" + name(n) + ',' + dsf.name() + ')', \
776 dsf.mesh(), \
777 dimless \
778 ); \
779 \
780 func(tres.ref().field(), n, dsf.field()); \
781 \
782 return tres; \
783} \
784 \
785 \
786template<class GeoMesh> \
787tmp<DimensionedField<scalar, GeoMesh>> func \
788( \
789 const int n, \
790 const tmp<DimensionedField<scalar, GeoMesh>>& tdsf \
791) \
792{ \
793 const auto& dsf = tdsf(); \
794 \
795 if (dimensionSet::checking() && !dsf.dimensions().dimensionless()) \
796 { \
797 FatalErrorInFunction \
798 << "Field is not dimensionless: " << dsf.dimensions() << nl \
799 << abort(FatalError); \
800 } \
801 \
802 tmp<DimensionedField<scalar, GeoMesh>> tres \
803 ( \
804 New \
805 ( \
806 tdsf, \
807 #func "(" + name(n) + ',' + dsf.name() + ')', \
808 dimless \
809 ) \
810 ); \
811 \
812 func(tres.ref().field(), n, dsf.field()); \
814 tdsf.clear(); \
815 \
816 return tres; \
817}
818
821
822#undef BesselFunc
823
824
825// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
826
827} // End namespace Foam
828
829// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
830
831#include "undefFieldFunctionsM.H"
832
833// ************************************************************************* //
#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 DimensionedField.
if(patchID !=-1)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const DynamicField< Type > & field() const noexcept
Return const-reference to the primitive field values.
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field....
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition GeoMesh.H:46
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static bool checking(bool on) noexcept
Turn dimension checking on/off.
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)