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-2020 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  (
237  tdsf1,
238  tdsf2,
239  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
240  dimless
241  );
242 
243  pow(tPow.ref().field(), dsf1.field(), dsf2.field());
244 
245  tdsf1.clear();
246  tdsf2.clear();
247 
248  return tPow;
249 }
250 
251 
252 template<class GeoMesh>
254 (
256  const dimensionedScalar& ds
257 )
258 {
259  if (!ds.dimensions().dimensionless())
260  {
262  << "Exponent is not dimensionless: " << ds.dimensions()
263  << exit(FatalError);
264  }
265 
267  (
269  (
270  "pow(" + dsf.name() + ',' + ds.name() + ')',
271  dsf.mesh(),
272  pow(dsf.dimensions(), ds)
273  )
274  );
275 
276  pow(tPow.ref().field(), dsf.field(), ds.value());
277 
278  return tPow;
279 }
280 
281 
282 template<class GeoMesh>
284 (
286  const dimensionedScalar& ds
287 )
288 {
289  if (!ds.dimensions().dimensionless())
290  {
292  << "Exponent is not dimensionless: " << ds.dimensions()
293  << exit(FatalError);
294  }
295 
296  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
297 
299  (
300  tdsf,
301  "pow(" + dsf.name() + ',' + ds.name() + ')',
302  pow(dsf.dimensions(), ds)
303  );
304 
305  pow(tPow.ref().field(), dsf.field(), ds.value());
306 
307  tdsf.clear();
308 
309  return tPow;
310 }
311 
312 
313 template<class GeoMesh>
315 (
317  const scalar& s
318 )
319 {
320  return pow(dsf, dimensionedScalar(s));
321 }
322 
323 
324 template<class GeoMesh>
326 (
328  const scalar& s
329 )
330 {
331  return pow(tdsf, dimensionedScalar(s));
332 }
333 
334 
335 template<class GeoMesh>
337 (
338  const dimensionedScalar& ds,
340 )
341 {
342  if (!ds.dimensions().dimensionless())
343  {
345  << "Base scalar is not dimensionless: " << ds.dimensions()
346  << exit(FatalError);
347  }
348 
349  if (!dsf.dimensions().dimensionless())
350  {
352  << "Exponent field is not dimensionless: " << dsf.dimensions()
353  << exit(FatalError);
354  }
355 
357  (
359  (
360  "pow(" + ds.name() + ',' + dsf.name() + ')',
361  dsf.mesh(),
362  dimless
363  )
364  );
365 
366  pow(tPow.ref().field(), ds.value(), dsf.field());
367 
368  return tPow;
369 }
370 
371 
372 template<class GeoMesh>
374 (
375  const dimensionedScalar& ds,
377 )
378 {
379  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
380 
381  if (!ds.dimensions().dimensionless())
382  {
384  << "Base scalar is not dimensionless: " << ds.dimensions()
385  << exit(FatalError);
386  }
387 
388  if (!dsf.dimensions().dimensionless())
389  {
391  << "Exponent field is not dimensionless: " << dsf.dimensions()
392  << exit(FatalError);
393  }
394 
396  (
397  tdsf,
398  "pow(" + ds.name() + ',' + dsf.name() + ')',
399  dimless
400  );
401 
402  pow(tPow.ref().field(), ds.value(), dsf.field());
403 
404  tdsf.clear();
405 
406  return tPow;
407 }
408 
409 template<class GeoMesh>
411 (
412  const scalar& s,
414 )
415 {
416  return pow(dimensionedScalar(s), dsf);
417 }
418 
419 template<class GeoMesh>
421 (
422  const scalar& s,
424 )
425 {
426  return pow(dimensionedScalar(s), tdsf);
427 }
428 
429 
430 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
431 
432 template<class GeoMesh>
434 (
437 )
438 {
440  (
442  (
443  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
444  dsf1.mesh(),
445  atan2(dsf1.dimensions(), dsf2.dimensions())
446  )
447  );
448 
449  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
450 
451  return tAtan2;
452 }
453 
454 
455 template<class GeoMesh>
457 (
460 )
461 {
462  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
463 
465  (
466  tdsf1,
467  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
468  atan2(dsf1.dimensions(), dsf2.dimensions())
469  );
470 
471  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
472 
473  tdsf1.clear();
474 
475  return tAtan2;
476 }
477 
478 
479 template<class GeoMesh>
481 (
484 )
485 {
486  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
487 
489  (
490  tdsf2,
491  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
492  atan2(dsf1.dimensions(), dsf2.dimensions())
493  );
494 
495  atan2(tAtan2.ref().field(), dsf1.field(), dsf2.field());
496 
497  tdsf2.clear();
498 
499  return tAtan2;
500 }
501 
502 template<class GeoMesh>
504 (
507 )
508 {
509  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
510  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
511 
514  (
515  tdsf1,
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  tdsf1.clear();
524  tdsf2.clear();
525 
526  return tAtan2;
527 }
528 
529 
530 template<class GeoMesh>
532 (
534  const dimensionedScalar& ds
535 )
536 {
538  (
540  (
541  "atan2(" + dsf.name() + ',' + ds.name() + ')',
542  dsf.mesh(),
543  atan2(dsf.dimensions(), ds)
544  )
545  );
546 
547  atan2(tAtan2.ref().field(), dsf.field(), ds.value());
548 
549  return tAtan2;
550 }
551 
552 template<class GeoMesh>
554 (
556  const dimensionedScalar& ds
557 )
558 {
559  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
560 
562  (
563  tdsf,
564  "atan2(" + dsf.name() + ',' + ds.name() + ')',
565  atan2(dsf.dimensions(), ds)
566  );
567 
568  atan2(tAtan2.ref().field(), dsf.field(), ds.value());
569 
570  tdsf.clear();
571 
572  return tAtan2;
573 }
574 
575 template<class GeoMesh>
577 (
579  const scalar& s
580 )
581 {
582  return atan2(dsf, dimensionedScalar(s));
583 }
584 
585 template<class GeoMesh>
587 (
589  const scalar& s
590 )
591 {
592  return atan2(tdsf, dimensionedScalar(s));
593 }
594 
595 
596 template<class GeoMesh>
598 (
599  const dimensionedScalar& ds,
601 )
602 {
604  (
606  (
607  "atan2(" + ds.name() + ',' + dsf.name() + ')',
608  dsf.mesh(),
609  atan2(ds, dsf.dimensions())
610  )
611  );
612 
613  atan2(tAtan2.ref().field(), ds.value(), dsf.field());
614 
615  return tAtan2;
616 }
617 
618 
619 template<class GeoMesh>
621 (
622  const dimensionedScalar& ds,
624 )
625 {
626  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
627 
629  (
630  tdsf,
631  "atan2(" + ds.name() + ',' + dsf.name() + ')',
632  atan2(ds, dsf.dimensions())
633  );
634 
635  atan2(tAtan2.ref().field(), ds.value(), dsf.field());
636 
637  tdsf.clear();
638 
639  return tAtan2;
640 }
641 
642 template<class GeoMesh>
644 (
645  const scalar& s,
647 )
648 {
649  return atan2(dimensionedScalar(s), dsf);
650 }
651 
652 template<class GeoMesh>
654 (
655  const scalar& s,
657 )
658 {
659  return atan2(dimensionedScalar(s), tdsf);
660 }
661 
662 
663 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
664 
665 UNARY_FUNCTION(scalar, scalar, pow3, pow3)
666 UNARY_FUNCTION(scalar, scalar, pow4, pow4)
667 UNARY_FUNCTION(scalar, scalar, pow5, pow5)
668 UNARY_FUNCTION(scalar, scalar, pow6, pow6)
669 UNARY_FUNCTION(scalar, scalar, pow025, pow025)
670 UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
671 UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
672 UNARY_FUNCTION(scalar, scalar, sign, sign)
673 UNARY_FUNCTION(scalar, scalar, pos, pos)
674 UNARY_FUNCTION(scalar, scalar, pos0, pos0)
675 UNARY_FUNCTION(scalar, scalar, neg, neg)
676 UNARY_FUNCTION(scalar, scalar, neg0, neg0)
677 UNARY_FUNCTION(scalar, scalar, posPart, posPart)
678 UNARY_FUNCTION(scalar, scalar, negPart, negPart)
679 
680 UNARY_FUNCTION(scalar, scalar, exp, trans)
681 UNARY_FUNCTION(scalar, scalar, log, trans)
682 UNARY_FUNCTION(scalar, scalar, log10, trans)
683 UNARY_FUNCTION(scalar, scalar, sin, trans)
684 UNARY_FUNCTION(scalar, scalar, cos, trans)
685 UNARY_FUNCTION(scalar, scalar, tan, trans)
686 UNARY_FUNCTION(scalar, scalar, asin, trans)
687 UNARY_FUNCTION(scalar, scalar, acos, trans)
688 UNARY_FUNCTION(scalar, scalar, atan, trans)
689 UNARY_FUNCTION(scalar, scalar, sinh, trans)
690 UNARY_FUNCTION(scalar, scalar, cosh, trans)
691 UNARY_FUNCTION(scalar, scalar, tanh, trans)
692 UNARY_FUNCTION(scalar, scalar, asinh, trans)
693 UNARY_FUNCTION(scalar, scalar, acosh, trans)
694 UNARY_FUNCTION(scalar, scalar, atanh, trans)
695 UNARY_FUNCTION(scalar, scalar, erf, trans)
696 UNARY_FUNCTION(scalar, scalar, erfc, trans)
697 UNARY_FUNCTION(scalar, scalar, lgamma, trans)
698 UNARY_FUNCTION(scalar, scalar, j0, trans)
699 UNARY_FUNCTION(scalar, scalar, j1, trans)
700 UNARY_FUNCTION(scalar, scalar, y0, trans)
701 UNARY_FUNCTION(scalar, scalar, y1, trans)
702 
703 
704 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
705 
706 #define BesselFunc(func) \
707  \
708 template<class GeoMesh> \
709 tmp<DimensionedField<scalar, GeoMesh>> func \
710 ( \
711  const int n, \
712  const DimensionedField<scalar, GeoMesh>& dsf \
713 ) \
714 { \
715  if (!dsf.dimensions().dimensionless()) \
716  { \
717  FatalErrorInFunction \
718  << "dsf not dimensionless" \
719  << abort(FatalError); \
720  } \
721  \
722  tmp<DimensionedField<scalar, GeoMesh>> tFunc \
723  ( \
724  DimensionedField<scalar, GeoMesh>::New \
725  ( \
726  #func "(" + name(n) + ',' + dsf.name() + ')', \
727  dsf.mesh(), \
728  dimless \
729  ) \
730  ); \
731  \
732  func(tFunc.ref().field(), n, dsf.field()); \
733  \
734  return tFunc; \
735 } \
736  \
737 template<class GeoMesh> \
738 tmp<DimensionedField<scalar, GeoMesh>> func \
739 ( \
740  const int n, \
741  const tmp<DimensionedField<scalar, GeoMesh>>& tdsf \
742 ) \
743 { \
744  const DimensionedField<scalar, GeoMesh>& dsf = tdsf(); \
745  \
746  if (!dsf.dimensions().dimensionless()) \
747  { \
748  FatalErrorInFunction \
749  << " : dsf not dimensionless" \
750  << abort(FatalError); \
751  } \
752  \
753  tmp<DimensionedField<scalar, GeoMesh>> tFunc \
754  ( \
755  New \
756  ( \
757  tdsf, \
758  #func "(" + name(n) + ',' + dsf.name() + ')', \
759  dimless \
760  ) \
761  ); \
762  \
763  func(tFunc.ref().field(), n, dsf.field()); \
764  \
765  tdsf.clear(); \
766  \
767  return tFunc; \
768 }
769 
772 
773 #undef BesselFunc
774 
775 
776 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
777 
778 } // End namespace Foam
779 
780 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
781 
782 #include "undefFieldFunctionsM.H"
783 
784 // ************************************************************************* //
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:323
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
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 dimensionSet dimless
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)
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)