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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "DimensionedScalarField.H"
27 
28 #define TEMPLATE template<class GeoMesh>
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37 
38 template<class GeoMesh>
39 tmp<DimensionedField<scalar, GeoMesh>> stabilise
40 (
42  const dimensioned<scalar>& ds
43 )
44 {
46  (
48  (
49  IOobject
50  (
51  "stabilise(" + dsf.name() + ',' + ds.name() + ')',
52  dsf.instance(),
53  dsf.db()
54  ),
55  dsf.mesh(),
56  dsf.dimensions() + ds.dimensions()
57  )
58  );
59 
60  stabilise(tRes.ref().field(), dsf.field(), ds.value());
61 
62  return tRes;
63 }
64 
65 
66 template<class GeoMesh>
68 (
70  const dimensioned<scalar>& ds
71 )
72 {
73  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
74 
76  (
77  tdsf,
78  "stabilise(" + dsf.name() + ',' + ds.name() + ')',
79  dsf.dimensions() + ds.dimensions()
80  );
81 
82  stabilise(tRes.ref().field(), dsf.field(), ds.value());
83 
84  tdsf.clear();
85 
86  return tRes;
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 
92 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, '+', add)
93 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, '-', subtract)
94 
95 BINARY_OPERATOR(scalar, scalar, scalar, *, '*', multiply)
96 BINARY_OPERATOR(scalar, scalar, scalar, /, '|', divide)
97 
98 BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, '|', divide)
99 
100 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
101 
102 template<class GeoMesh>
104 (
107 )
108 {
109  if (!dsf1.dimensions().dimensionless())
110  {
112  << "Base field is not dimensionless: " << dsf1.dimensions()
113  << exit(FatalError);
114  }
115 
116  if (!dsf2.dimensions().dimensionless())
117  {
119  << "Exponent field is not dimensionless: " << dsf2.dimensions()
120  << exit(FatalError);
121  }
122 
124  (
126  (
127  IOobject
128  (
129  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
130  dsf1.instance(),
131  dsf1.db()
132  ),
133  dsf1.mesh(),
134  dimless
135  )
136  );
137 
138  pow(tPow.ref().field(), dsf1.field(), dsf2.field());
139 
140  return tPow;
141 }
142 
143 
144 template<class GeoMesh>
146 (
149 )
150 {
151  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
152 
153  if (!dsf1.dimensions().dimensionless())
154  {
156  << "Base field is not dimensionless: " << dsf1.dimensions()
157  << exit(FatalError);
158  }
159 
160  if (!dsf2.dimensions().dimensionless())
161  {
163  << "Exponent field is not dimensionless: " << dsf2.dimensions()
164  << exit(FatalError);
165  }
166 
168  (
169  tdsf1,
170  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
171  dimless
172  );
173 
174  pow(tPow.ref().field(), dsf1.field(), dsf2.field());
175 
176  tdsf1.clear();
177 
178  return tPow;
179 }
180 
181 
182 template<class GeoMesh>
184 (
187 )
188 {
189  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
190 
191  if (!dsf1.dimensions().dimensionless())
192  {
194  << "Base field is not dimensionless: " << dsf1.dimensions()
195  << exit(FatalError);
196  }
197 
198  if (!dsf2.dimensions().dimensionless())
199  {
201  << "Exponent field is not dimensionless: " << dsf2.dimensions()
202  << exit(FatalError);
203  }
204 
206  (
207  tdsf2,
208  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
209  dimless
210  );
211 
212  pow(tPow.ref().field(), dsf1.field(), dsf2.field());
213 
214  tdsf2.clear();
215 
216  return tPow;
217 }
218 
219 
220 template<class GeoMesh>
222 (
225 )
226 {
227  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
228  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
229 
230  if (!dsf1.dimensions().dimensionless())
231  {
233  << "Base field is not dimensionless: " << dsf1.dimensions()
234  << exit(FatalError);
235  }
236 
237  if (!dsf2.dimensions().dimensionless())
238  {
240  << "Exponent field is not dimensionless: " << dsf2.dimensions()
241  << exit(FatalError);
242  }
243 
246  New
247  (
248  tdsf1,
249  tdsf2,
250  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
251  dimless
252  );
253 
254  pow(tPow.ref().field(), dsf1.field(), dsf2.field());
255 
256  tdsf1.clear();
257  tdsf2.clear();
258 
259  return tPow;
260 }
261 
262 
263 template<class GeoMesh>
265 (
267  const dimensionedScalar& ds
268 )
269 {
270  if (!ds.dimensions().dimensionless())
271  {
273  << "Exponent is not dimensionless: " << ds.dimensions()
274  << exit(FatalError);
275  }
276 
278  (
280  (
281  IOobject
282  (
283  "pow(" + dsf.name() + ',' + ds.name() + ')',
284  dsf.instance(),
285  dsf.db()
286  ),
287  dsf.mesh(),
288  pow(dsf.dimensions(), ds)
289  )
290  );
291 
292  pow(tPow.ref().field(), dsf.field(), ds.value());
293 
294  return tPow;
295 }
296 
297 
298 template<class GeoMesh>
300 (
302  const dimensionedScalar& ds
303 )
304 {
305  if (!ds.dimensions().dimensionless())
306  {
308  << "Exponent is not dimensionless: " << ds.dimensions()
309  << exit(FatalError);
310  }
311 
312  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
313 
315  (
316  tdsf,
317  "pow(" + dsf.name() + ',' + ds.name() + ')',
318  pow(dsf.dimensions(), ds)
319  );
320 
321  pow(tPow.ref().field(), dsf.field(), ds.value());
322 
323  tdsf.clear();
324 
325  return tPow;
326 }
327 
328 
329 template<class GeoMesh>
331 (
333  const scalar& s
334 )
335 {
336  return pow(dsf, dimensionedScalar(s));
337 }
338 
339 
340 template<class GeoMesh>
342 (
344  const scalar& s
345 )
346 {
347  return pow(tdsf, dimensionedScalar(s));
348 }
349 
350 
351 template<class GeoMesh>
353 (
354  const dimensionedScalar& ds,
356 )
357 {
358  if (!ds.dimensions().dimensionless())
359  {
361  << "Base scalar is not dimensionless: " << ds.dimensions()
362  << exit(FatalError);
363  }
364 
365  if (!dsf.dimensions().dimensionless())
366  {
368  << "Exponent field is not dimensionless: " << dsf.dimensions()
369  << exit(FatalError);
370  }
371 
373  (
375  (
376  IOobject
377  (
378  "pow(" + ds.name() + ',' + dsf.name() + ')',
379  dsf.instance(),
380  dsf.db()
381  ),
382  dsf.mesh(),
383  dimless
384  )
385  );
386 
387  pow(tPow.ref().field(), ds.value(), dsf.field());
388 
389  return tPow;
390 }
391 
392 
393 template<class GeoMesh>
395 (
396  const dimensionedScalar& ds,
398 )
399 {
400  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
401 
402  if (!ds.dimensions().dimensionless())
403  {
405  << "Base scalar is not dimensionless: " << ds.dimensions()
406  << exit(FatalError);
407  }
408 
409  if (!dsf.dimensions().dimensionless())
410  {
412  << "Exponent field is not dimensionless: " << dsf.dimensions()
413  << exit(FatalError);
414  }
415 
417  (
418  tdsf,
419  "pow(" + ds.name() + ',' + dsf.name() + ')',
420  dimless
421  );
422 
423  pow(tPow.ref().field(), ds.value(), dsf.field());
424 
425  tdsf.clear();
426 
427  return tPow;
428 }
429 
430 template<class GeoMesh>
432 (
433  const scalar& s,
435 )
436 {
437  return pow(dimensionedScalar(s), dsf);
438 }
439 
440 template<class GeoMesh>
442 (
443  const scalar& s,
445 )
446 {
447  return pow(dimensionedScalar(s), tdsf);
448 }
449 
450 
451 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
452 
453 template<class GeoMesh>
455 (
458 )
459 {
461  (
463  (
464  IOobject
465  (
466  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
467  dsf1.instance(),
468  dsf1.db()
469  ),
470  dsf1.mesh(),
471  atan2(dsf1.dimensions(), dsf2.dimensions())
472  )
473  );
474 
475  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
476 
477  return tAtan2;
478 }
479 
480 
481 template<class GeoMesh>
483 (
486 )
487 {
488  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
489 
491  (
492  tdsf1,
493  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
494  atan2(dsf1.dimensions(), dsf2.dimensions())
495  );
496 
497  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
498 
499  tdsf1.clear();
500 
501  return tAtan2;
502 }
503 
504 
505 template<class GeoMesh>
507 (
510 )
511 {
512  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
513 
515  (
516  tdsf2,
517  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
518  atan2(dsf1.dimensions(), dsf2.dimensions())
519  );
520 
521  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
522 
523  tdsf2.clear();
524 
525  return tAtan2;
526 }
527 
528 template<class GeoMesh>
530 (
533 )
534 {
535  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
536  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
537 
540  New
541  (
542  tdsf1,
543  tdsf2,
544  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
545  atan2(dsf1.dimensions(), dsf2.dimensions())
546  );
547 
548  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
549 
550  tdsf1.clear();
551  tdsf2.clear();
552 
553  return tAtan2;
554 }
555 
556 
557 template<class GeoMesh>
559 (
561  const dimensionedScalar& ds
562 )
563 {
565  (
567  (
568  IOobject
569  (
570  "atan2(" + dsf.name() + ',' + ds.name() + ')',
571  dsf.instance(),
572  dsf.db()
573  ),
574  dsf.mesh(),
575  atan2(dsf.dimensions(), ds)
576  )
577  );
578 
579  atan2(tAtan2.ref().field(), dsf.field(), ds.value());
580 
581  return tAtan2;
582 }
583 
584 template<class GeoMesh>
586 (
588  const dimensionedScalar& ds
589 )
590 {
591  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
592 
594  (
595  tdsf,
596  "atan2(" + dsf.name() + ',' + ds.name() + ')',
597  atan2(dsf.dimensions(), ds)
598  );
599 
600  atan2(tAtan2.ref().field(), dsf.field(), ds.value());
601 
602  tdsf.clear();
603 
604  return tAtan2;
605 }
606 
607 template<class GeoMesh>
609 (
611  const scalar& s
612 )
613 {
614  return atan2(dsf, dimensionedScalar(s));
615 }
616 
617 template<class GeoMesh>
619 (
621  const scalar& s
622 )
623 {
624  return atan2(tdsf, dimensionedScalar(s));
625 }
626 
627 
628 template<class GeoMesh>
630 (
631  const dimensionedScalar& ds,
633 )
634 {
636  (
638  (
639  IOobject
640  (
641  "atan2(" + ds.name() + ',' + dsf.name() + ')',
642  dsf.instance(),
643  dsf.db()
644  ),
645  dsf.mesh(),
646  atan2(ds, dsf.dimensions())
647  )
648  );
649 
650  atan2(tAtan2.ref().field(), ds.value(), dsf.field());
651 
652  return tAtan2;
653 }
654 
655 
656 template<class GeoMesh>
658 (
659  const dimensionedScalar& ds,
661 )
662 {
663  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
664 
666  (
667  tdsf,
668  "atan2(" + ds.name() + ',' + dsf.name() + ')',
669  atan2(ds, dsf.dimensions())
670  );
671 
672  atan2(tAtan2.ref().field(), ds.value(), dsf.field());
673 
674  tdsf.clear();
675 
676  return tAtan2;
677 }
678 
679 template<class GeoMesh>
681 (
682  const scalar& s,
684 )
685 {
686  return atan2(dimensionedScalar(s), dsf);
687 }
688 
689 template<class GeoMesh>
691 (
692  const scalar& s,
694 )
695 {
696  return atan2(dimensionedScalar(s), tdsf);
697 }
698 
699 
700 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
701 
702 UNARY_FUNCTION(scalar, scalar, pow3, pow3)
703 UNARY_FUNCTION(scalar, scalar, pow4, pow4)
704 UNARY_FUNCTION(scalar, scalar, pow5, pow5)
705 UNARY_FUNCTION(scalar, scalar, pow6, pow6)
706 UNARY_FUNCTION(scalar, scalar, pow025, pow025)
707 UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
708 UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
709 UNARY_FUNCTION(scalar, scalar, sign, sign)
710 UNARY_FUNCTION(scalar, scalar, pos, pos)
711 UNARY_FUNCTION(scalar, scalar, neg, neg)
712 UNARY_FUNCTION(scalar, scalar, posPart, posPart)
713 UNARY_FUNCTION(scalar, scalar, negPart, negPart)
714 
715 UNARY_FUNCTION(scalar, scalar, exp, trans)
716 UNARY_FUNCTION(scalar, scalar, log, trans)
717 UNARY_FUNCTION(scalar, scalar, log10, trans)
718 UNARY_FUNCTION(scalar, scalar, sin, trans)
719 UNARY_FUNCTION(scalar, scalar, cos, trans)
720 UNARY_FUNCTION(scalar, scalar, tan, trans)
721 UNARY_FUNCTION(scalar, scalar, asin, trans)
722 UNARY_FUNCTION(scalar, scalar, acos, trans)
723 UNARY_FUNCTION(scalar, scalar, atan, trans)
724 UNARY_FUNCTION(scalar, scalar, sinh, trans)
725 UNARY_FUNCTION(scalar, scalar, cosh, trans)
726 UNARY_FUNCTION(scalar, scalar, tanh, trans)
727 UNARY_FUNCTION(scalar, scalar, asinh, trans)
728 UNARY_FUNCTION(scalar, scalar, acosh, trans)
729 UNARY_FUNCTION(scalar, scalar, atanh, trans)
730 UNARY_FUNCTION(scalar, scalar, erf, trans)
731 UNARY_FUNCTION(scalar, scalar, erfc, trans)
732 UNARY_FUNCTION(scalar, scalar, lgamma, trans)
733 UNARY_FUNCTION(scalar, scalar, j0, trans)
734 UNARY_FUNCTION(scalar, scalar, j1, trans)
735 UNARY_FUNCTION(scalar, scalar, y0, trans)
736 UNARY_FUNCTION(scalar, scalar, y1, trans)
737 
738 
739 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
740 
741 #define BesselFunc(func) \
742  \
743 template<class GeoMesh> \
744 tmp<DimensionedField<scalar, GeoMesh>> func \
745 ( \
746  const int n, \
747  const DimensionedField<scalar, GeoMesh>& dsf \
748 ) \
749 { \
750  if (!dsf.dimensions().dimensionless()) \
751  { \
752  FatalErrorInFunction \
753  << "dsf not dimensionless" \
754  << abort(FatalError); \
755  } \
756  \
757  tmp<DimensionedField<scalar, GeoMesh>> tFunc \
758  ( \
759  new DimensionedField<scalar, GeoMesh> \
760  ( \
761  IOobject \
762  ( \
763  #func "(" + name(n) + ',' + dsf.name() + ')', \
764  dsf.instance(), \
765  dsf.db() \
766  ), \
767  dsf.mesh(), \
768  dimless \
769  ) \
770  ); \
771  \
772  func(tFunc.ref().field(), n, dsf.field()); \
773  \
774  return tFunc; \
775 } \
776  \
777 template<class GeoMesh> \
778 tmp<DimensionedField<scalar, GeoMesh>> func \
779 ( \
780  const int n, \
781  const tmp<DimensionedField<scalar, GeoMesh>>& tdsf \
782 ) \
783 { \
784  const DimensionedField<scalar, GeoMesh>& dsf = tdsf(); \
785  \
786  if (!dsf.dimensions().dimensionless()) \
787  { \
788  FatalErrorInFunction \
789  << " : dsf not dimensionless" \
790  << abort(FatalError); \
791  } \
792  \
793  tmp<DimensionedField<scalar, GeoMesh>> tFunc \
794  ( \
795  New \
796  ( \
797  tdsf, \
798  #func "(" + name(n) + ',' + dsf.name() + ')', \
799  dimless \
800  ) \
801  ); \
802  \
803  func(tFunc.ref().field(), n, dsf.field()); \
804  \
805  tdsf.clear(); \
806  \
807  return tFunc; \
808 }
809 
812 
813 #undef BesselFunc
814 
815 
816 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
817 
818 } // End namespace Foam
819 
820 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
821 
822 #include "undefFieldFunctionsM.H"
823 
824 // ************************************************************************* //
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionSet trans(const dimensionSet &)
Definition: dimensionSet.C:438
dimensionedScalar acos(const dimensionedScalar &ds)
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
static tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< Type1, GeoMesh >> &tdf1, const tmp< DimensionedField< Type2, GeoMesh >> &tdf2, const word &name, const dimensionSet &dimensions)
void multiply(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionedScalar log(const dimensionedScalar &ds)
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const dimensionSet & dimensions() const
Return const reference to dimensions.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
const Type & value() const
Return const reference to value.
void subtract(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar j1(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const word & name() const
Return const reference to name.
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
#define BesselFunc(func)
dimensionedScalar atanh(const dimensionedScalar &ds)
dimensionedScalar y1(const dimensionedScalar &ds)
const dimensionSet & dimensions() const
Return dimensions.
const Field< Type > & field() const
dimensionedScalar sin(const dimensionedScalar &ds)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
bool dimensionless() const
Return true if it is dimensionless.
Definition: dimensionSet.C:89
const Mesh & mesh() const
Return mesh.
dimensionedScalar atan(const dimensionedScalar &ds)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar pow6(const dimensionedScalar &ds)
Scalar specific part of the implementation of DimensionedField.
dimensionedScalar cosh(const dimensionedScalar &ds)
A class for managing temporary objects.
Definition: PtrList.H:54
dimensionedScalar tan(const dimensionedScalar &ds)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedScalar j0(const dimensionedScalar &ds)
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const fileName & instance() const
Definition: IOobject.H:337
dimensionedScalar log10(const dimensionedScalar &ds)
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
dimensionedScalar negPart(const dimensionedScalar &ds)