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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29
30#define TEMPLATE template<template<class> class PatchField, class GeoMesh>
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39
40template<template<class> class PatchField, class GeoMesh>
42(
45 const dimensioned<scalar>& ds
46)
47{
48 stabilise(result.primitiveFieldRef(), gsf.primitiveField(), ds.value());
49 stabilise(result.boundaryFieldRef(), gsf.boundaryField(), ds.value());
50}
51
52
53template<template<class> class PatchField, class GeoMesh>
55(
57 const dimensioned<scalar>& ds
58)
59{
61 (
63 (
65 (
66 "stabilise(" + gsf.name() + ',' + ds.name() + ')',
67 gsf.instance(),
68 gsf.db(),
71 ),
72 gsf.mesh(),
73 ds.dimensions() + gsf.dimensions()
74 )
75 );
76
77 stabilise(tRes.ref(), gsf, ds);
78
79 return tRes;
80}
81
82
83template<template<class> class PatchField, class GeoMesh>
85(
87 const dimensioned<scalar>& ds
88)
89{
91
93 (
94 New
95 (
96 tgsf,
97 "stabilise(" + gsf.name() + ',' + ds.name() + ')',
98 ds.dimensions() + gsf.dimensions()
99 )
100 );
101
102 stabilise(tRes.ref(), gsf, ds);
103
104 tgsf.clear();
105
106 return tRes;
107}
108
109
110// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
111
112BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, '+', add)
113BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, '-', subtract)
114
115BINARY_OPERATOR(scalar, scalar, scalar, *, '*', multiply)
116BINARY_OPERATOR(scalar, scalar, scalar, /, '|', divide)
117
118BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, '|', divide)
119
120
121// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122
123template<template<class> class PatchField, class GeoMesh>
124void pow
125(
129)
130{
131 pow(Pow.primitiveFieldRef(), gsf1.primitiveField(), gsf2.primitiveField());
132 pow(Pow.boundaryFieldRef(), gsf1.boundaryField(), gsf2.boundaryField());
133}
134
135
136template<template<class> class PatchField, class GeoMesh>
138(
141)
142{
143 if (!gsf1.dimensions().dimensionless())
144 {
146 << "Base field is not dimensionless: " << gsf1.dimensions()
147 << exit(FatalError);
148 }
149
150 if (!gsf2.dimensions().dimensionless())
151 {
153 << "Exponent field is not dimensionless: " << gsf2.dimensions()
154 << exit(FatalError);
155 }
156
158 (
160 (
162 (
163 "pow(" + gsf1.name() + ',' + gsf2.name() + ')',
164 gsf1.instance(),
165 gsf1.db(),
168 ),
169 gsf1.mesh(),
170 dimless
171 )
172 );
173
174 pow(tPow.ref(), gsf1, gsf2);
175
176 return tPow;
177}
178
179
180template<template<class> class PatchField, class GeoMesh>
182(
185)
186{
188
189 if (!gsf1.dimensions().dimensionless())
190 {
192 << "Base field is not dimensionless: " << gsf1.dimensions()
193 << exit(FatalError);
194 }
195
196 if (!gsf2.dimensions().dimensionless())
197 {
199 << "Exponent field is not dimensionless: " << gsf2.dimensions()
200 << exit(FatalError);
201 }
202
204 (
205 New
206 (
207 tgsf1,
208 "pow(" + gsf1.name() + ',' + gsf2.name() + ')',
209 dimless
210 )
211 );
212
213 pow(tPow.ref(), gsf1, gsf2);
214
215 tgsf1.clear();
216
217 return tPow;
218}
219
220
221template<template<class> class PatchField, class GeoMesh>
223(
226)
227{
229
230 if (!gsf1.dimensions().dimensionless())
231 {
233 << "Base field is not dimensionless: " << gsf1.dimensions()
234 << exit(FatalError);
235 }
236
237 if (!gsf2.dimensions().dimensionless())
238 {
240 << "Exponent field is not dimensionless: " << gsf2.dimensions()
241 << exit(FatalError);
242 }
243
245 (
246 New
247 (
248 tgsf2,
249 "pow(" + gsf1.name() + ',' + gsf2.name() + ')',
250 dimless
251 )
252 );
253
254 pow(tPow.ref(), gsf1, gsf2);
255
256 tgsf2.clear();
257
258 return tPow;
259}
260
261template<template<class> class PatchField, class GeoMesh>
263(
266)
267{
270
271 if (!gsf1.dimensions().dimensionless())
272 {
274 << "Base field is not dimensionless: " << gsf1.dimensions()
275 << exit(FatalError);
276 }
277
278 if (!gsf2.dimensions().dimensionless())
279 {
281 << "Exponent field is not dimensionless: " << gsf2.dimensions()
282 << exit(FatalError);
283 }
284
286 (
288 <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::New
289 (
290 tgsf1,
291 tgsf2,
292 "pow(" + gsf1.name() + ',' + gsf2.name() + ')',
293 dimless
294 )
295 );
296
297 pow(tPow.ref(), gsf1, gsf2);
298
299 tgsf1.clear();
300 tgsf2.clear();
301
302 return tPow;
303}
304
305
306template<template<class> class PatchField, class GeoMesh>
307void pow
308(
311 const dimensioned<scalar>& ds
312)
313{
314 pow(tPow.primitiveFieldRef(), gsf.primitiveField(), ds.value());
315 pow(tPow.boundaryFieldRef(), gsf.boundaryField(), ds.value());
316}
317
318
319template<template<class> class PatchField, class GeoMesh>
321(
323 const dimensionedScalar& ds
324)
325{
326 if (!ds.dimensions().dimensionless())
327 {
329 << "Exponent is not dimensionless: " << ds.dimensions()
330 << exit(FatalError);
331 }
332
334 (
336 (
338 (
339 "pow(" + gsf.name() + ',' + ds.name() + ')',
340 gsf.instance(),
341 gsf.db(),
344 ),
345 gsf.mesh(),
346 pow(gsf.dimensions(), ds)
347 )
348 );
349
350 pow(tPow.ref(), gsf, ds);
351
352 return tPow;
353}
354
355template<template<class> class PatchField, class GeoMesh>
357(
359 const dimensionedScalar& ds
360)
361{
362 if (!ds.dimensions().dimensionless())
363 {
365 << "Exponent is not dimensionless: " << ds.dimensions()
366 << exit(FatalError);
367 }
368
370
372 (
373 New
374 (
375 tgsf,
376 "pow(" + gsf.name() + ',' + ds.name() + ')',
377 pow(gsf.dimensions(), ds)
378 )
379 );
380
381 pow(tPow.ref(), gsf, ds);
382
383 tgsf.clear();
384
385 return tPow;
386}
387
388template<template<class> class PatchField, class GeoMesh>
390(
392 const scalar& s
393)
394{
395 return pow(gsf, dimensionedScalar(s));
396}
397
398template<template<class> class PatchField, class GeoMesh>
400(
402 const scalar& s
403)
404{
405 return pow(tgsf, dimensionedScalar(s));
406}
407
408
409template<template<class> class PatchField, class GeoMesh>
410void pow
411(
413 const dimensioned<scalar>& ds,
415)
416{
417 pow(tPow.primitiveFieldRef(), ds.value(), gsf.primitiveField());
418 pow(tPow.boundaryFieldRef(), ds.value(), gsf.boundaryField());
419}
420
421
422template<template<class> class PatchField, class GeoMesh>
424(
425 const dimensionedScalar& ds,
427)
428{
429 if (!ds.dimensions().dimensionless())
430 {
432 << "Base scalar is not dimensionless: " << ds.dimensions()
433 << exit(FatalError);
434 }
435
436 if (!gsf.dimensions().dimensionless())
437 {
439 << "Exponent field is not dimensionless: " << gsf.dimensions()
440 << exit(FatalError);
441 }
442
444 (
446 (
448 (
449 "pow(" + ds.name() + ',' + gsf.name() + ')',
450 gsf.instance(),
451 gsf.db(),
454 ),
455 gsf.mesh(),
456 dimless
457 )
458 );
459
460 pow(tPow.ref(), ds, gsf);
461
462 return tPow;
463}
464
465
466template<template<class> class PatchField, class GeoMesh>
468(
469 const dimensionedScalar& ds,
471)
472{
474
475 if (!ds.dimensions().dimensionless())
476 {
478 << "Base scalar is not dimensionless: " << ds.dimensions()
479 << exit(FatalError);
480 }
481
482 if (!gsf.dimensions().dimensionless())
483 {
485 << "Exponent field is not dimensionless: " << gsf.dimensions()
486 << exit(FatalError);
487 }
488
490 (
491 New
492 (
493 tgsf,
494 "pow(" + ds.name() + ',' + gsf.name() + ')',
495 dimless
496 )
497 );
498
499 pow(tPow.ref(), ds, gsf);
500
501 tgsf.clear();
502
503 return tPow;
504}
505
506template<template<class> class PatchField, class GeoMesh>
508(
509 const scalar& s,
511)
512{
513 return pow(dimensionedScalar(s), gsf);
514}
515
516template<template<class> class PatchField, class GeoMesh>
518(
519 const scalar& s,
521)
522{
523 return pow(dimensionedScalar(s), tgsf);
524}
525
526
527// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
528
529template<template<class> class PatchField, class GeoMesh>
531(
535)
536{
537 atan2
538 (
539 Atan2.primitiveFieldRef(),
540 gsf1.primitiveField(),
541 gsf2.primitiveField()
542 );
543 atan2
544 (
545 Atan2.boundaryFieldRef(),
546 gsf1.boundaryField(),
547 gsf2.boundaryField()
548 );
549}
550
551
552template<template<class> class PatchField, class GeoMesh>
554(
557)
558{
560 (
562 (
564 (
565 "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
566 gsf1.instance(),
567 gsf1.db(),
570 ),
571 gsf1.mesh(),
572 atan2(gsf1.dimensions(), gsf2.dimensions())
573 )
574 );
575
576 atan2(tAtan2.ref(), gsf1, gsf2);
577
578 return tAtan2;
579}
580
581
582template<template<class> class PatchField, class GeoMesh>
584(
587)
588{
590
592 (
593 New
594 (
595 tgsf1,
596 "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
597 atan2(gsf1.dimensions(), gsf2.dimensions())
598 )
599 );
600
601 atan2(tAtan2.ref(), gsf1, gsf2);
602
603 tgsf1.clear();
604
605 return tAtan2;
606}
607
608
609template<template<class> class PatchField, class GeoMesh>
611(
614)
615{
617
619 (
620 New
621 (
622 tgsf2,
623 "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
624 atan2( gsf1.dimensions(), gsf2.dimensions())
625 )
626 );
627
628 atan2(tAtan2.ref(), gsf1, gsf2);
629
630 tgsf2.clear();
631
632 return tAtan2;
633}
634
635template<template<class> class PatchField, class GeoMesh>
637(
640)
641{
644
646 (
648 <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::New
649 (
650 tgsf1,
651 tgsf2,
652 "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
653 atan2(gsf1.dimensions(), gsf2.dimensions())
654 )
655 );
656
657 atan2(tAtan2.ref(), gsf1, gsf2);
658
659 tgsf1.clear();
660 tgsf2.clear();
661
662 return tAtan2;
663}
664
665
666template<template<class> class PatchField, class GeoMesh>
668(
671 const dimensioned<scalar>& ds
672)
673{
674 atan2(tAtan2.primitiveFieldRef(), gsf.primitiveField(), ds.value());
675 atan2(tAtan2.boundaryFieldRef(), gsf.boundaryField(), ds.value());
676}
677
678
679template<template<class> class PatchField, class GeoMesh>
681(
683 const dimensionedScalar& ds
684)
685{
687 (
689 (
691 (
692 "atan2(" + gsf.name() + ',' + ds.name() + ')',
693 gsf.instance(),
694 gsf.db(),
697 ),
698 gsf.mesh(),
699 atan2(gsf.dimensions(), ds)
700 )
701 );
702
703 atan2(tAtan2.ref(), gsf, ds);
704
705 return tAtan2;
706}
707
708template<template<class> class PatchField, class GeoMesh>
710(
712 const dimensionedScalar& ds
713)
714{
716
718 (
719 New
720 (
721 tgsf,
722 "atan2(" + gsf.name() + ',' + ds.name() + ')',
723 atan2(gsf.dimensions(), ds)
724 )
725 );
726
727 atan2(tAtan2.ref(), gsf, ds);
728
729 tgsf.clear();
730
731 return tAtan2;
732}
733
734template<template<class> class PatchField, class GeoMesh>
736(
738 const scalar& s
739)
740{
741 return atan2(gsf, dimensionedScalar(s));
742}
743
744template<template<class> class PatchField, class GeoMesh>
746(
748 const scalar& s
749)
750{
751 return atan2(tgsf, dimensionedScalar(s));
752}
753
754
755template<template<class> class PatchField, class GeoMesh>
757(
759 const dimensioned<scalar>& ds,
761)
762{
763 atan2(tAtan2.primitiveFieldRef(), ds.value(), gsf.primitiveField());
764 atan2(tAtan2.boundaryFieldRef(), ds.value(), gsf.boundaryField());
765}
766
767
768template<template<class> class PatchField, class GeoMesh>
770(
771 const dimensionedScalar& ds,
773)
774{
776 (
778 (
780 (
781 "atan2(" + ds.name() + ',' + gsf.name() + ')',
782 gsf.instance(),
783 gsf.db(),
786 ),
787 gsf.mesh(),
788 atan2(ds, gsf.dimensions())
789 )
790 );
791
792 atan2(tAtan2.ref(), ds, gsf);
793
794 return tAtan2;
795}
796
797
798template<template<class> class PatchField, class GeoMesh>
800(
801 const dimensionedScalar& ds,
803)
804{
806
808 (
809 New
810 (
811 tgsf,
812 "atan2(" + ds.name() + ',' + gsf.name() + ')',
813 atan2(ds, gsf.dimensions())
814 )
815 );
816
817 atan2(tAtan2.ref(), ds, gsf);
818
819 tgsf.clear();
820
821 return tAtan2;
822}
823
824template<template<class> class PatchField, class GeoMesh>
826(
827 const scalar& s,
829)
830{
831 return atan2(dimensionedScalar(s), gsf);
832}
833
834template<template<class> class PatchField, class GeoMesh>
836(
837 const scalar& s,
839)
840{
841 return atan2(dimensionedScalar(s), tgsf);
842}
843
844
845// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
846
847UNARY_FUNCTION(scalar, scalar, pow3, pow3)
848UNARY_FUNCTION(scalar, scalar, pow4, pow4)
849UNARY_FUNCTION(scalar, scalar, pow5, pow5)
850UNARY_FUNCTION(scalar, scalar, pow6, pow6)
851UNARY_FUNCTION(scalar, scalar, pow025, pow025)
852UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
853UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
854UNARY_FUNCTION(scalar, scalar, sign, sign)
855UNARY_FUNCTION(scalar, scalar, pos, pos)
856UNARY_FUNCTION(scalar, scalar, pos0, pos0)
857UNARY_FUNCTION(scalar, scalar, neg, neg)
858UNARY_FUNCTION(scalar, scalar, neg0, neg0)
859UNARY_FUNCTION(scalar, scalar, posPart, posPart)
860UNARY_FUNCTION(scalar, scalar, negPart, negPart)
861
862UNARY_FUNCTION(scalar, scalar, exp, trans)
863UNARY_FUNCTION(scalar, scalar, log, trans)
864UNARY_FUNCTION(scalar, scalar, log10, trans)
865UNARY_FUNCTION(scalar, scalar, sin, trans)
866UNARY_FUNCTION(scalar, scalar, cos, trans)
867UNARY_FUNCTION(scalar, scalar, tan, trans)
868UNARY_FUNCTION(scalar, scalar, asin, trans)
869UNARY_FUNCTION(scalar, scalar, acos, trans)
870UNARY_FUNCTION(scalar, scalar, atan, trans)
871UNARY_FUNCTION(scalar, scalar, sinh, trans)
872UNARY_FUNCTION(scalar, scalar, cosh, trans)
873UNARY_FUNCTION(scalar, scalar, tanh, trans)
874UNARY_FUNCTION(scalar, scalar, asinh, trans)
875UNARY_FUNCTION(scalar, scalar, acosh, trans)
876UNARY_FUNCTION(scalar, scalar, atanh, trans)
877UNARY_FUNCTION(scalar, scalar, erf, trans)
878UNARY_FUNCTION(scalar, scalar, erfc, trans)
879UNARY_FUNCTION(scalar, scalar, lgamma, trans)
880UNARY_FUNCTION(scalar, scalar, j0, trans)
881UNARY_FUNCTION(scalar, scalar, j1, trans)
882UNARY_FUNCTION(scalar, scalar, y0, trans)
883UNARY_FUNCTION(scalar, scalar, y1, trans)
884
885
886// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
887
888#define BesselFunc(func) \
889 \
890template<template<class> class PatchField, class GeoMesh> \
891void func \
892( \
893 GeometricField<scalar, PatchField, GeoMesh>& gsf, \
894 const int n, \
895 const GeometricField<scalar, PatchField, GeoMesh>& gsf1 \
896) \
897{ \
898 func(gsf.primitiveFieldRef(), n, gsf1.primitiveField()); \
899 func(gsf.boundaryFieldRef(), n, gsf1.boundaryField()); \
900} \
901 \
902template<template<class> class PatchField, class GeoMesh> \
903tmp<GeometricField<scalar, PatchField, GeoMesh>> func \
904( \
905 const int n, \
906 const GeometricField<scalar, PatchField, GeoMesh>& gsf \
907) \
908{ \
909 if (!gsf.dimensions().dimensionless()) \
910 { \
911 FatalErrorInFunction \
912 << "gsf not dimensionless" \
913 << abort(FatalError); \
914 } \
915 \
916 tmp<GeometricField<scalar, PatchField, GeoMesh>> tFunc \
917 ( \
918 new GeometricField<scalar, PatchField, GeoMesh> \
919 ( \
920 IOobject \
921 ( \
922 #func "(" + gsf.name() + ')', \
923 gsf.instance(), \
924 gsf.db(), \
925 IOobject::NO_READ, \
926 IOobject::NO_WRITE \
927 ), \
928 gsf.mesh(), \
929 dimless \
930 ) \
931 ); \
932 \
933 func(tFunc.ref(), n, gsf); \
934 \
935 return tFunc; \
936} \
937 \
938template<template<class> class PatchField, class GeoMesh> \
939tmp<GeometricField<scalar, PatchField, GeoMesh>> func \
940( \
941 const int n, \
942 const tmp<GeometricField<scalar, PatchField, GeoMesh>>& tgsf \
943) \
944{ \
945 const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf(); \
946 \
947 if (!gsf.dimensions().dimensionless()) \
948 { \
949 FatalErrorInFunction \
950 << " : gsf not dimensionless" \
951 << abort(FatalError); \
952 } \
953 \
954 tmp<GeometricField<scalar, PatchField, GeoMesh>> tFunc \
955 ( \
956 New \
957 ( \
958 tgsf, \
959 #func "(" + gsf.name() + ')', \
960 dimless \
961 ) \
962 ); \
963 \
964 func(tFunc.ref(), n, gsf); \
965 \
966 tgsf.clear(); \
967 \
968 return tFunc; \
969}
970
973
974#undef BesselFunc
975
976
977// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
978
979} // End namespace Foam
980
981// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
982
983#include "undefFieldFunctionsM.H"
984
985// ************************************************************************* //
#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 dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:49
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
bool dimensionless() const
Return true if it is dimensionless.
Definition: dimensionSet.C:114
Generic dimensioned Type class.
const dimensionSet & dimensions() const
Return const reference to dimensions.
const Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
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)
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
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)
dimensionedScalar sqrt(const dimensionedScalar &ds)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionSet trans(const dimensionSet &ds)
Check the argument is dimensionless (for transcendental functions)
Definition: dimensionSet.C:486
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar atan(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
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)
dimensionedScalar asinh(const dimensionedScalar &ds)