DimensionedScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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  "stabilise(" + dsf.name() + ',' + ds.name() + ')',
50  dsf.mesh(),
51  dsf.dimensions() + ds.dimensions()
52  )
53  );
54 
55  stabilise(tRes.ref().field(), dsf.field(), ds.value());
56 
57  return tRes;
58 }
59 
60 
61 template<class GeoMesh>
63 (
65  const dimensioned<scalar>& ds
66 )
67 {
68  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
69 
71  (
72  tdsf,
73  "stabilise(" + dsf.name() + ',' + ds.name() + ')',
74  dsf.dimensions() + ds.dimensions()
75  );
76 
77  stabilise(tRes.ref().field(), dsf.field(), ds.value());
78 
79  tdsf.clear();
80 
81  return tRes;
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, '+', add)
88 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, '-', subtract)
89 
90 BINARY_OPERATOR(scalar, scalar, scalar, *, '*', multiply)
91 BINARY_OPERATOR(scalar, scalar, scalar, /, '|', divide)
92 
93 BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, '|', divide)
94 
95 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
96 
97 template<class GeoMesh>
99 (
102 )
103 {
104  if (!dsf1.dimensions().dimensionless())
105  {
107  << "Base field is not dimensionless: " << dsf1.dimensions()
108  << exit(FatalError);
109  }
110 
111  if (!dsf2.dimensions().dimensionless())
112  {
114  << "Exponent field is not dimensionless: " << dsf2.dimensions()
115  << exit(FatalError);
116  }
117 
119  (
121  (
122  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
123  dsf1.mesh(),
124  dimless
125  )
126  );
127 
128  pow(tPow.ref().field(), dsf1.field(), dsf2.field());
129 
130  return tPow;
131 }
132 
133 
134 template<class GeoMesh>
136 (
139 )
140 {
141  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
142 
143  if (!dsf1.dimensions().dimensionless())
144  {
146  << "Base field is not dimensionless: " << dsf1.dimensions()
147  << exit(FatalError);
148  }
149 
150  if (!dsf2.dimensions().dimensionless())
151  {
153  << "Exponent field is not dimensionless: " << dsf2.dimensions()
154  << exit(FatalError);
155  }
156 
158  (
159  tdsf1,
160  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
161  dimless
162  );
163 
164  pow(tPow.ref().field(), dsf1.field(), dsf2.field());
165 
166  tdsf1.clear();
167 
168  return tPow;
169 }
170 
171 
172 template<class GeoMesh>
174 (
177 )
178 {
179  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
180 
181  if (!dsf1.dimensions().dimensionless())
182  {
184  << "Base field is not dimensionless: " << dsf1.dimensions()
185  << exit(FatalError);
186  }
187 
188  if (!dsf2.dimensions().dimensionless())
189  {
191  << "Exponent field is not dimensionless: " << dsf2.dimensions()
192  << exit(FatalError);
193  }
194 
196  (
197  tdsf2,
198  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
199  dimless
200  );
201 
202  pow(tPow.ref().field(), dsf1.field(), dsf2.field());
203 
204  tdsf2.clear();
205 
206  return tPow;
207 }
208 
209 
210 template<class GeoMesh>
212 (
215 )
216 {
217  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
218  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
219 
220  if (!dsf1.dimensions().dimensionless())
221  {
223  << "Base field is not dimensionless: " << dsf1.dimensions()
224  << exit(FatalError);
225  }
226 
227  if (!dsf2.dimensions().dimensionless())
228  {
230  << "Exponent field is not dimensionless: " << dsf2.dimensions()
231  << exit(FatalError);
232  }
233 
236  New
237  (
238  tdsf1,
239  tdsf2,
240  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
241  dimless
242  );
243 
244  pow(tPow.ref().field(), dsf1.field(), dsf2.field());
245 
246  tdsf1.clear();
247  tdsf2.clear();
248 
249  return tPow;
250 }
251 
252 
253 template<class GeoMesh>
255 (
257  const dimensionedScalar& ds
258 )
259 {
260  if (!ds.dimensions().dimensionless())
261  {
263  << "Exponent is not dimensionless: " << ds.dimensions()
264  << exit(FatalError);
265  }
266 
268  (
270  (
271  "pow(" + dsf.name() + ',' + ds.name() + ')',
272  dsf.mesh(),
273  pow(dsf.dimensions(), ds)
274  )
275  );
276 
277  pow(tPow.ref().field(), dsf.field(), ds.value());
278 
279  return tPow;
280 }
281 
282 
283 template<class GeoMesh>
285 (
287  const dimensionedScalar& ds
288 )
289 {
290  if (!ds.dimensions().dimensionless())
291  {
293  << "Exponent is not dimensionless: " << ds.dimensions()
294  << exit(FatalError);
295  }
296 
297  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
298 
300  (
301  tdsf,
302  "pow(" + dsf.name() + ',' + ds.name() + ')',
303  pow(dsf.dimensions(), ds)
304  );
305 
306  pow(tPow.ref().field(), dsf.field(), ds.value());
307 
308  tdsf.clear();
309 
310  return tPow;
311 }
312 
313 
314 template<class GeoMesh>
316 (
318  const scalar& s
319 )
320 {
321  return pow(dsf, dimensionedScalar(s));
322 }
323 
324 
325 template<class GeoMesh>
327 (
329  const scalar& s
330 )
331 {
332  return pow(tdsf, dimensionedScalar(s));
333 }
334 
335 
336 template<class GeoMesh>
338 (
339  const dimensionedScalar& ds,
341 )
342 {
343  if (!ds.dimensions().dimensionless())
344  {
346  << "Base scalar is not dimensionless: " << ds.dimensions()
347  << exit(FatalError);
348  }
349 
350  if (!dsf.dimensions().dimensionless())
351  {
353  << "Exponent field is not dimensionless: " << dsf.dimensions()
354  << exit(FatalError);
355  }
356 
358  (
360  (
361  "pow(" + ds.name() + ',' + dsf.name() + ')',
362  dsf.mesh(),
363  dimless
364  )
365  );
366 
367  pow(tPow.ref().field(), ds.value(), dsf.field());
368 
369  return tPow;
370 }
371 
372 
373 template<class GeoMesh>
375 (
376  const dimensionedScalar& ds,
378 )
379 {
380  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
381 
382  if (!ds.dimensions().dimensionless())
383  {
385  << "Base scalar is not dimensionless: " << ds.dimensions()
386  << exit(FatalError);
387  }
388 
389  if (!dsf.dimensions().dimensionless())
390  {
392  << "Exponent field is not dimensionless: " << dsf.dimensions()
393  << exit(FatalError);
394  }
395 
397  (
398  tdsf,
399  "pow(" + ds.name() + ',' + dsf.name() + ')',
400  dimless
401  );
402 
403  pow(tPow.ref().field(), ds.value(), dsf.field());
404 
405  tdsf.clear();
406 
407  return tPow;
408 }
409 
410 template<class GeoMesh>
412 (
413  const scalar& s,
415 )
416 {
417  return pow(dimensionedScalar(s), dsf);
418 }
419 
420 template<class GeoMesh>
422 (
423  const scalar& s,
425 )
426 {
427  return pow(dimensionedScalar(s), tdsf);
428 }
429 
430 
431 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
432 
433 template<class GeoMesh>
435 (
438 )
439 {
441  (
443  (
444  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
445  dsf1.mesh(),
446  atan2(dsf1.dimensions(), dsf2.dimensions())
447  )
448  );
449 
450  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
451 
452  return tAtan2;
453 }
454 
455 
456 template<class GeoMesh>
458 (
461 )
462 {
463  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
464 
466  (
467  tdsf1,
468  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
469  atan2(dsf1.dimensions(), dsf2.dimensions())
470  );
471 
472  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
473 
474  tdsf1.clear();
475 
476  return tAtan2;
477 }
478 
479 
480 template<class GeoMesh>
482 (
485 )
486 {
487  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
488 
490  (
491  tdsf2,
492  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
493  atan2(dsf1.dimensions(), dsf2.dimensions())
494  );
495 
496  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
497 
498  tdsf2.clear();
499 
500  return tAtan2;
501 }
502 
503 template<class GeoMesh>
505 (
508 )
509 {
510  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
511  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
512 
515  New
516  (
517  tdsf1,
518  tdsf2,
519  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
520  atan2(dsf1.dimensions(), dsf2.dimensions())
521  );
522 
523  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
524 
525  tdsf1.clear();
526  tdsf2.clear();
527 
528  return tAtan2;
529 }
530 
531 
532 template<class GeoMesh>
534 (
536  const dimensionedScalar& ds
537 )
538 {
540  (
542  (
543  "atan2(" + dsf.name() + ',' + ds.name() + ')',
544  dsf.mesh(),
545  atan2(dsf.dimensions(), ds)
546  )
547  );
548 
549  atan2(tAtan2.ref().field(), dsf.field(), ds.value());
550 
551  return tAtan2;
552 }
553 
554 template<class GeoMesh>
556 (
558  const dimensionedScalar& ds
559 )
560 {
561  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
562 
564  (
565  tdsf,
566  "atan2(" + dsf.name() + ',' + ds.name() + ')',
567  atan2(dsf.dimensions(), ds)
568  );
569 
570  atan2(tAtan2.ref().field(), dsf.field(), ds.value());
571 
572  tdsf.clear();
573 
574  return tAtan2;
575 }
576 
577 template<class GeoMesh>
579 (
581  const scalar& s
582 )
583 {
584  return atan2(dsf, dimensionedScalar(s));
585 }
586 
587 template<class GeoMesh>
589 (
591  const scalar& s
592 )
593 {
594  return atan2(tdsf, dimensionedScalar(s));
595 }
596 
597 
598 template<class GeoMesh>
600 (
601  const dimensionedScalar& ds,
603 )
604 {
606  (
608  (
609  "atan2(" + ds.name() + ',' + dsf.name() + ')',
610  dsf.mesh(),
611  atan2(ds, dsf.dimensions())
612  )
613  );
614 
615  atan2(tAtan2.ref().field(), ds.value(), dsf.field());
616 
617  return tAtan2;
618 }
619 
620 
621 template<class GeoMesh>
623 (
624  const dimensionedScalar& ds,
626 )
627 {
628  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
629 
631  (
632  tdsf,
633  "atan2(" + ds.name() + ',' + dsf.name() + ')',
634  atan2(ds, dsf.dimensions())
635  );
636 
637  atan2(tAtan2.ref().field(), ds.value(), dsf.field());
638 
639  tdsf.clear();
640 
641  return tAtan2;
642 }
643 
644 template<class GeoMesh>
646 (
647  const scalar& s,
649 )
650 {
651  return atan2(dimensionedScalar(s), dsf);
652 }
653 
654 template<class GeoMesh>
656 (
657  const scalar& s,
659 )
660 {
661  return atan2(dimensionedScalar(s), tdsf);
662 }
663 
664 
665 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
666 
667 UNARY_FUNCTION(scalar, scalar, pow3, pow3)
668 UNARY_FUNCTION(scalar, scalar, pow4, pow4)
669 UNARY_FUNCTION(scalar, scalar, pow5, pow5)
670 UNARY_FUNCTION(scalar, scalar, pow6, pow6)
671 UNARY_FUNCTION(scalar, scalar, pow025, pow025)
672 UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
673 UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
674 UNARY_FUNCTION(scalar, scalar, sign, sign)
675 UNARY_FUNCTION(scalar, scalar, pos, pos)
676 UNARY_FUNCTION(scalar, scalar, pos0, pos0)
677 UNARY_FUNCTION(scalar, scalar, neg, neg)
678 UNARY_FUNCTION(scalar, scalar, neg0, neg0)
679 UNARY_FUNCTION(scalar, scalar, posPart, posPart)
680 UNARY_FUNCTION(scalar, scalar, negPart, negPart)
681 
682 UNARY_FUNCTION(scalar, scalar, exp, trans)
683 UNARY_FUNCTION(scalar, scalar, log, trans)
684 UNARY_FUNCTION(scalar, scalar, log10, trans)
685 UNARY_FUNCTION(scalar, scalar, sin, trans)
686 UNARY_FUNCTION(scalar, scalar, cos, trans)
687 UNARY_FUNCTION(scalar, scalar, tan, trans)
688 UNARY_FUNCTION(scalar, scalar, asin, trans)
689 UNARY_FUNCTION(scalar, scalar, acos, trans)
690 UNARY_FUNCTION(scalar, scalar, atan, trans)
691 UNARY_FUNCTION(scalar, scalar, sinh, trans)
692 UNARY_FUNCTION(scalar, scalar, cosh, trans)
693 UNARY_FUNCTION(scalar, scalar, tanh, trans)
694 UNARY_FUNCTION(scalar, scalar, asinh, trans)
695 UNARY_FUNCTION(scalar, scalar, acosh, trans)
696 UNARY_FUNCTION(scalar, scalar, atanh, trans)
697 UNARY_FUNCTION(scalar, scalar, erf, trans)
698 UNARY_FUNCTION(scalar, scalar, erfc, trans)
699 UNARY_FUNCTION(scalar, scalar, lgamma, trans)
700 UNARY_FUNCTION(scalar, scalar, j0, trans)
701 UNARY_FUNCTION(scalar, scalar, j1, trans)
702 UNARY_FUNCTION(scalar, scalar, y0, trans)
703 UNARY_FUNCTION(scalar, scalar, y1, trans)
704 
705 
706 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
707 
708 #define BesselFunc(func) \
709  \
710 template<class GeoMesh> \
711 tmp<DimensionedField<scalar, GeoMesh>> func \
712 ( \
713  const int n, \
714  const DimensionedField<scalar, GeoMesh>& dsf \
715 ) \
716 { \
717  if (!dsf.dimensions().dimensionless()) \
718  { \
719  FatalErrorInFunction \
720  << "dsf not dimensionless" \
721  << abort(FatalError); \
722  } \
723  \
724  tmp<DimensionedField<scalar, GeoMesh>> tFunc \
725  ( \
726  DimensionedField<scalar, GeoMesh>::New \
727  ( \
728  #func "(" + name(n) + ',' + dsf.name() + ')', \
729  dsf.mesh(), \
730  dimless \
731  ) \
732  ); \
733  \
734  func(tFunc.ref().field(), n, dsf.field()); \
735  \
736  return tFunc; \
737 } \
738  \
739 template<class GeoMesh> \
740 tmp<DimensionedField<scalar, GeoMesh>> func \
741 ( \
742  const int n, \
743  const tmp<DimensionedField<scalar, GeoMesh>>& tdsf \
744 ) \
745 { \
746  const DimensionedField<scalar, GeoMesh>& dsf = tdsf(); \
747  \
748  if (!dsf.dimensions().dimensionless()) \
749  { \
750  FatalErrorInFunction \
751  << " : dsf not dimensionless" \
752  << abort(FatalError); \
753  } \
754  \
755  tmp<DimensionedField<scalar, GeoMesh>> tFunc \
756  ( \
757  New \
758  ( \
759  tdsf, \
760  #func "(" + name(n) + ',' + dsf.name() + ')', \
761  dimless \
762  ) \
763  ); \
764  \
765  func(tFunc.ref().field(), n, dsf.field()); \
766  \
767  tdsf.clear(); \
768  \
769  return tFunc; \
770 }
771 
774 
775 #undef BesselFunc
776 
777 
778 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
779 
780 } // End namespace Foam
781 
782 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
783 
784 #include "undefFieldFunctionsM.H"
785 
786 // ************************************************************************* //
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:450
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)
const word & name() const
Return name.
Definition: IOobject.H:303
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
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)
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)
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet & dimensions() const
Return dimensions.
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)
dimensionedScalar neg0(const dimensionedScalar &ds)
#define BesselFunc(func)
const Type & value() const
Return const reference to value.
dimensionedScalar atanh(const dimensionedScalar &ds)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
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)
const Mesh & mesh() const
Return mesh.
const word & name() const
Return const reference to name.
const Field< Type > & field() const
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.
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)
const dimensionSet & dimensions() const
Return const reference to dimensions.
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:53
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
bool dimensionless() const
Return true if it is dimensionless.
Definition: dimensionSet.C:89
dimensionedScalar log10(const dimensionedScalar &ds)
Namespace for OpenFOAM.
dimensionedScalar negPart(const dimensionedScalar &ds)