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-2024 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, template<class> class PrimitiveField>
29 #define TEMPLATE2 \
30  template \
31  < \
32  class GeoMesh, \
33  template<class> class PrimitiveField1, \
34  template<class> class PrimitiveField2 \
35  >
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 
43 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
44 
45 template<class GeoMesh, template<class> class PrimitiveField>
47 (
49  const dimensioned<scalar>& ds
50 )
51 {
53  (
55  (
56  "stabilise(" + dsf.name() + ',' + ds.name() + ')',
57  dsf.mesh(),
58  dsf.dimensions() + ds.dimensions()
59  )
60  );
61 
62  stabilise(tRes.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
63 
64  return tRes;
65 }
66 
67 
68 template<class GeoMesh, template<class> class PrimitiveField>
70 (
72  const dimensioned<scalar>& ds
73 )
74 {
75  const DimensionedField<scalar, GeoMesh, Field>& dsf = tdsf();
76 
78  (
79  tdsf,
80  "stabilise(" + dsf.name() + ',' + ds.name() + ')',
81  dsf.dimensions() + ds.dimensions()
82  );
83 
84  stabilise(tRes.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
85 
86  tdsf.clear();
87 
88  return tRes;
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93 
94 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, '+', add)
95 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, '-', subtract)
96 
97 BINARY_OPERATOR(scalar, scalar, scalar, *, '*', multiply)
98 BINARY_OPERATOR(scalar, scalar, scalar, /, '|', divide)
99 
100 BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, '|', divide)
101 
102 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
103 
104 template<class GeoMesh, template<class> class PrimitiveField>
106 (
109 )
110 {
111  if (!dsf1.dimensions().dimensionless())
112  {
114  << "Base field is not dimensionless: " << dsf1.dimensions()
115  << exit(FatalError);
116  }
117 
118  if (!dsf2.dimensions().dimensionless())
119  {
121  << "Exponent field is not dimensionless: " << dsf2.dimensions()
122  << exit(FatalError);
123  }
124 
126  (
128  (
129  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
130  dsf1.mesh(),
131  dimless
132  )
133  );
134 
135  pow
136  (
137  tPow.ref().primitiveFieldRef(),
138  dsf1.primitiveField(),
139  dsf2.primitiveField()
140  );
141 
142  return tPow;
143 }
144 
145 
146 template
147 <
148  class GeoMesh,
149  template<class> class PrimitiveField1,
150  template<class> class PrimitiveField2
151 >
153 (
156 )
157 {
159 
160  if (!dsf1.dimensions().dimensionless())
161  {
163  << "Base field is not dimensionless: " << dsf1.dimensions()
164  << exit(FatalError);
165  }
166 
167  if (!dsf2.dimensions().dimensionless())
168  {
170  << "Exponent field is not dimensionless: " << dsf2.dimensions()
171  << exit(FatalError);
172  }
173 
175  (
176  tdsf1,
177  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
178  dimless
179  );
180 
181  pow
182  (
183  tPow.ref().primitiveFieldRef(),
184  dsf1.primitiveField(),
185  dsf2.primitiveField()
186  );
187 
188  tdsf1.clear();
189 
190  return tPow;
191 }
192 
193 
194 template
195 <
196  class GeoMesh,
197  template<class> class PrimitiveField1,
198  template<class> class PrimitiveField2
199 >
201 (
204 )
205 {
207 
208  if (!dsf1.dimensions().dimensionless())
209  {
211  << "Base field is not dimensionless: " << dsf1.dimensions()
212  << exit(FatalError);
213  }
214 
215  if (!dsf2.dimensions().dimensionless())
216  {
218  << "Exponent field is not dimensionless: " << dsf2.dimensions()
219  << exit(FatalError);
220  }
221 
223  (
224  tdsf2,
225  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
226  dimless
227  );
228 
229  pow
230  (
231  tPow.ref().primitiveFieldRef(),
232  dsf1.primitiveField(),
233  dsf2.primitiveField()
234  );
235 
236  tdsf2.clear();
237 
238  return tPow;
239 }
240 
241 
242 template
243 <
244  class GeoMesh,
245  template<class> class PrimitiveField1,
246  template<class> class PrimitiveField2
247 >
249 (
252 )
253 {
256 
257  if (!dsf1.dimensions().dimensionless())
258  {
260  << "Base field is not dimensionless: " << dsf1.dimensions()
261  << exit(FatalError);
262  }
263 
264  if (!dsf2.dimensions().dimensionless())
265  {
267  << "Exponent field is not dimensionless: " << dsf2.dimensions()
268  << exit(FatalError);
269  }
270 
273  <
274  scalar,
275  scalar,
276  scalar,
277  GeoMesh,
278  PrimitiveField1,
279  PrimitiveField2
280  >::New
281  (
282  tdsf1,
283  tdsf2,
284  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
285  dimless
286  );
287 
288  pow
289  (
290  tPow.ref().primitiveFieldRef(),
291  dsf1.primitiveField(),
292  dsf2.primitiveField()
293  );
294 
295  tdsf1.clear();
296  tdsf2.clear();
297 
298  return tPow;
299 }
300 
301 
302 template<class GeoMesh, template<class> class PrimitiveField>
304 (
306  const dimensionedScalar& ds
307 )
308 {
309  if (!ds.dimensions().dimensionless())
310  {
312  << "Exponent is not dimensionless: " << ds.dimensions()
313  << exit(FatalError);
314  }
315 
317  (
319  (
320  "pow(" + dsf.name() + ',' + ds.name() + ')',
321  dsf.mesh(),
322  pow(dsf.dimensions(), ds)
323  )
324  );
325 
326  pow(tPow.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
327 
328  return tPow;
329 }
330 
331 
332 template<class GeoMesh, template<class> class PrimitiveField>
334 (
336  const dimensionedScalar& ds
337 )
338 {
339  if (!ds.dimensions().dimensionless())
340  {
342  << "Exponent is not dimensionless: " << ds.dimensions()
343  << exit(FatalError);
344  }
345 
347 
349  (
350  tdsf,
351  "pow(" + dsf.name() + ',' + ds.name() + ')',
352  pow(dsf.dimensions(), ds)
353  );
354 
355  pow(tPow.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
356 
357  tdsf.clear();
358 
359  return tPow;
360 }
361 
362 
363 template<class GeoMesh, template<class> class PrimitiveField>
365 (
367  const scalar& s
368 )
369 {
370  return pow(dsf, dimensionedScalar(s));
371 }
372 
373 
374 template<class GeoMesh, template<class> class PrimitiveField>
376 (
378  const scalar& s
379 )
380 {
381  return pow(tdsf, dimensionedScalar(s));
382 }
383 
384 
385 template<class GeoMesh, template<class> class PrimitiveField>
387 (
388  const dimensionedScalar& ds,
390 )
391 {
392  if (!ds.dimensions().dimensionless())
393  {
395  << "Base scalar is not dimensionless: " << ds.dimensions()
396  << exit(FatalError);
397  }
398 
399  if (!dsf.dimensions().dimensionless())
400  {
402  << "Exponent field is not dimensionless: " << dsf.dimensions()
403  << exit(FatalError);
404  }
405 
407  (
409  (
410  "pow(" + ds.name() + ',' + dsf.name() + ')',
411  dsf.mesh(),
412  dimless
413  )
414  );
415 
416  pow(tPow.ref().primitiveFieldRef(), ds.value(), dsf.primitiveField());
417 
418  return tPow;
419 }
420 
421 
422 template<class GeoMesh, template<class> class PrimitiveField>
424 (
425  const dimensionedScalar& ds,
427 )
428 {
430 
431  if (!ds.dimensions().dimensionless())
432  {
434  << "Base scalar is not dimensionless: " << ds.dimensions()
435  << exit(FatalError);
436  }
437 
438  if (!dsf.dimensions().dimensionless())
439  {
441  << "Exponent field is not dimensionless: " << dsf.dimensions()
442  << exit(FatalError);
443  }
444 
446  (
447  tdsf,
448  "pow(" + ds.name() + ',' + dsf.name() + ')',
449  dimless
450  );
451 
452  pow(tPow.ref().primitiveFieldRef(), ds.value(), dsf.primitiveField());
453 
454  tdsf.clear();
455 
456  return tPow;
457 }
458 
459 template<class GeoMesh, template<class> class PrimitiveField>
461 (
462  const scalar& s,
464 )
465 {
466  return pow(dimensionedScalar(s), dsf);
467 }
468 
469 template<class GeoMesh, template<class> class PrimitiveField>
471 (
472  const scalar& s,
474 )
475 {
476  return pow(dimensionedScalar(s), tdsf);
477 }
478 
479 
480 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
481 
482 template
483 <
484  class GeoMesh,
485  template<class> class PrimitiveField1,
486  template<class> class PrimitiveField2
487 >
489 (
492 )
493 {
495  (
497  (
498  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
499  dsf1.mesh(),
500  atan2(dsf1.dimensions(), dsf2.dimensions())
501  )
502  );
503 
504  atan2
505  (
506  tAtan2.ref().primitiveFieldRef(),
507  dsf1.primitiveField(),
508  dsf2.primitiveField()
509  );
510 
511  return tAtan2;
512 }
513 
514 
515 template
516 <
517  class GeoMesh,
518  template<class> class PrimitiveField1,
519  template<class> class PrimitiveField2
520 >
522 (
525 )
526 {
528 
530  (
531  tdsf1,
532  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
533  atan2(dsf1.dimensions(), dsf2.dimensions())
534  );
535 
536  atan2
537  (
538  tAtan2.ref().primitiveFieldRef(),
539  dsf1.primitiveField(),
540  dsf2.primitiveField()
541  );
542 
543  tdsf1.clear();
544 
545  return tAtan2;
546 }
547 
548 
549 template
550 <
551  class GeoMesh,
552  template<class> class PrimitiveField1,
553  template<class> class PrimitiveField2
554 >
556 (
559 )
560 {
562 
564  (
565  tdsf2,
566  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
567  atan2(dsf1.dimensions(), dsf2.dimensions())
568  );
569 
570  atan2
571  (
572  tAtan2.ref().primitiveFieldRef(),
573  dsf1.primitiveField(),
574  dsf2.primitiveField()
575  );
576 
577  tdsf2.clear();
578 
579  return tAtan2;
580 }
581 
582 template
583 <
584  class GeoMesh,
585  template<class> class PrimitiveField1,
586  template<class> class PrimitiveField2
587 >
589 (
592 )
593 {
596 
599  <
600  scalar,
601  scalar,
602  scalar,
603  GeoMesh,
604  PrimitiveField1,
605  PrimitiveField2
606  >::New
607  (
608  tdsf1,
609  tdsf2,
610  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
611  atan2(dsf1.dimensions(), dsf2.dimensions())
612  );
613 
614  atan2
615  (
616  tAtan2.ref().primitiveFieldRef(),
617  dsf1.primitiveField(),
618  dsf2.primitiveField()
619  );
620 
621  tdsf1.clear();
622  tdsf2.clear();
623 
624  return tAtan2;
625 }
626 
627 
628 template<class GeoMesh, template<class> class PrimitiveField>
630 (
632  const dimensionedScalar& ds
633 )
634 {
636  (
638  (
639  "atan2(" + dsf.name() + ',' + ds.name() + ')',
640  dsf.mesh(),
641  atan2(dsf.dimensions(), ds)
642  )
643  );
644 
645  atan2(tAtan2.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
646 
647  return tAtan2;
648 }
649 
650 template<class GeoMesh, template<class> class PrimitiveField>
652 (
654  const dimensionedScalar& ds
655 )
656 {
658 
660  (
661  tdsf,
662  "atan2(" + dsf.name() + ',' + ds.name() + ')',
663  atan2(dsf.dimensions(), ds)
664  );
665 
666  atan2(tAtan2.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
667 
668  tdsf.clear();
669 
670  return tAtan2;
671 }
672 
673 template<class GeoMesh, template<class> class PrimitiveField>
675 (
677  const scalar& s
678 )
679 {
680  return atan2(dsf, dimensionedScalar(s));
681 }
682 
683 template<class GeoMesh, template<class> class PrimitiveField>
685 (
687  const scalar& s
688 )
689 {
690  return atan2(tdsf, dimensionedScalar(s));
691 }
692 
693 
694 template<class GeoMesh, template<class> class PrimitiveField>
696 (
697  const dimensionedScalar& ds,
699 )
700 {
702  (
704  (
705  "atan2(" + ds.name() + ',' + dsf.name() + ')',
706  dsf.mesh(),
707  atan2(ds, dsf.dimensions())
708  )
709  );
710 
711  atan2(tAtan2.ref().primitiveFieldRef(), ds.value(), dsf.primitiveField());
712 
713  return tAtan2;
714 }
715 
716 
717 template<class GeoMesh, template<class> class PrimitiveField>
719 (
720  const dimensionedScalar& ds,
722 )
723 {
725 
727  (
728  tdsf,
729  "atan2(" + ds.name() + ',' + dsf.name() + ')',
730  atan2(ds, dsf.dimensions())
731  );
732 
733  atan2(tAtan2.ref().primitiveFieldRef(), ds.value(), dsf.primitiveField());
734 
735  tdsf.clear();
736 
737  return tAtan2;
738 }
739 
740 template<class GeoMesh, template<class> class PrimitiveField>
742 (
743  const scalar& s,
745 )
746 {
747  return atan2(dimensionedScalar(s), dsf);
748 }
749 
750 template<class GeoMesh, template<class> class PrimitiveField>
752 (
753  const scalar& s,
755 )
756 {
757  return atan2(dimensionedScalar(s), tdsf);
758 }
759 
760 
761 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
762 
763 UNARY_FUNCTION(scalar, scalar, pow3, pow3)
764 UNARY_FUNCTION(scalar, scalar, pow4, pow4)
765 UNARY_FUNCTION(scalar, scalar, pow5, pow5)
766 UNARY_FUNCTION(scalar, scalar, pow6, pow6)
767 UNARY_FUNCTION(scalar, scalar, pow025, pow025)
768 UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
769 UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
770 UNARY_FUNCTION(scalar, scalar, sign, sign)
771 UNARY_FUNCTION(scalar, scalar, pos, pos)
772 UNARY_FUNCTION(scalar, scalar, pos0, pos0)
773 UNARY_FUNCTION(scalar, scalar, neg, neg)
774 UNARY_FUNCTION(scalar, scalar, neg0, neg0)
775 UNARY_FUNCTION(scalar, scalar, posPart, posPart)
776 UNARY_FUNCTION(scalar, scalar, negPart, negPart)
777 
778 UNARY_FUNCTION(scalar, scalar, exp, trans)
779 UNARY_FUNCTION(scalar, scalar, log, trans)
780 UNARY_FUNCTION(scalar, scalar, log10, trans)
781 UNARY_FUNCTION(scalar, scalar, sin, trans)
782 UNARY_FUNCTION(scalar, scalar, cos, trans)
783 UNARY_FUNCTION(scalar, scalar, tan, trans)
784 UNARY_FUNCTION(scalar, scalar, asin, trans)
785 UNARY_FUNCTION(scalar, scalar, acos, trans)
786 UNARY_FUNCTION(scalar, scalar, atan, trans)
787 UNARY_FUNCTION(scalar, scalar, sinh, trans)
788 UNARY_FUNCTION(scalar, scalar, cosh, trans)
789 UNARY_FUNCTION(scalar, scalar, tanh, trans)
790 UNARY_FUNCTION(scalar, scalar, asinh, trans)
791 UNARY_FUNCTION(scalar, scalar, acosh, trans)
792 UNARY_FUNCTION(scalar, scalar, atanh, trans)
793 UNARY_FUNCTION(scalar, scalar, erf, trans)
794 UNARY_FUNCTION(scalar, scalar, erfc, trans)
795 UNARY_FUNCTION(scalar, scalar, lgamma, trans)
796 UNARY_FUNCTION(scalar, scalar, j0, trans)
797 UNARY_FUNCTION(scalar, scalar, j1, trans)
798 UNARY_FUNCTION(scalar, scalar, y0, trans)
799 UNARY_FUNCTION(scalar, scalar, y1, trans)
800 
801 
802 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
803 
804 #define BesselFunc(func) \
805  \
806 template<class GeoMesh, template<class> class PrimitiveField> \
807 tmp<DimensionedField<scalar, GeoMesh, Field>> func \
808 ( \
809  const int n, \
810  const DimensionedField<scalar, GeoMesh, PrimitiveField>& dsf \
811 ) \
812 { \
813  if (!dsf.dimensions().dimensionless()) \
814  { \
815  FatalErrorInFunction \
816  << "dsf not dimensionless" \
817  << abort(FatalError); \
818  } \
819  \
820  tmp<DimensionedField<scalar, GeoMesh, Field>> tFunc \
821  ( \
822  DimensionedField<scalar, GeoMesh, Field>::New \
823  ( \
824  #func "(" + name(n) + ',' + dsf.name() + ')', \
825  dsf.mesh(), \
826  dimless \
827  ) \
828  ); \
829  \
830  func(tFunc.ref().primitiveFieldRef(), n, dsf.primitiveField()); \
831  \
832  return tFunc; \
833 } \
834  \
835 template<class GeoMesh, template<class> class PrimitiveField> \
836 tmp<DimensionedField<scalar, GeoMesh, Field>> func \
837 ( \
838  const int n, \
839  const tmp<DimensionedField<scalar, GeoMesh, PrimitiveField>>& tdsf \
840 ) \
841 { \
842  const DimensionedField<scalar, GeoMesh, PrimitiveField>& dsf = tdsf(); \
843  \
844  if (!dsf.dimensions().dimensionless()) \
845  { \
846  FatalErrorInFunction \
847  << " : dsf not dimensionless" \
848  << abort(FatalError); \
849  } \
850  \
851  tmp<DimensionedField<scalar, GeoMesh, Field>> tFunc \
852  ( \
853  New \
854  ( \
855  tdsf, \
856  #func "(" + name(n) + ',' + dsf.name() + ')', \
857  dimless \
858  ) \
859  ); \
860  \
861  func(tFunc.ref().primitiveFieldRef(), n, dsf.primitiveField()); \
862  \
863  tdsf.clear(); \
864  \
865  return tFunc; \
866 }
867 
870 
871 #undef BesselFunc
872 
873 
874 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
875 
876 } // End namespace Foam
877 
878 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
879 
880 #include "undefFieldFunctionsM.H"
881 
882 // ************************************************************************* //
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
#define BesselFunc(func)
Scalar specific part of the implementation of DimensionedField.
#define BINARY_OPERATOR(Template, Type, Type1, Type2, op, opFunc)
#define UNARY_FUNCTION(Template, Type, Type1, func)
#define BINARY_TYPE_OPERATOR_SF(TYPE, op, opFunc)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
const PrimitiveField< Type > & primitiveField() const
Return a const-reference to the primitive field.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
const word & name() const
Return name.
Definition: IOobject.H:307
bool dimensionless() const
Return true if it is dimensionless.
Definition: dimensionSet.C:107
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:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
void subtract(LagrangianPatchField< typename typeOfSum< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar tan(const dimensionedScalar &ds)
void divide(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1, const LagrangianPatchField< scalar > &f2)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
void pow025(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
const dimensionSet dimless
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)
tmp< DimensionedField< scalar, GeoMesh, Field > > stabilise(const DimensionedField< scalar, GeoMesh, PrimitiveField > &dsf, const dimensioned< scalar > &ds)
void pow4(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void pow6(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
void multiply(LagrangianPatchField< Type > &f, const LagrangianPatchField< scalar > &f1, const LagrangianPatchField< Type > &f2)
void pow5(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar neg(const dimensionedScalar &ds)
void cbrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar atanh(const dimensionedScalar &ds)
error FatalError
dimensionedScalar neg0(const dimensionedScalar &ds)
void add(LagrangianPatchField< typename typeOfSum< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)
dimensionedScalar atan(const dimensionedScalar &ds)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionSet trans(const dimensionSet &)
Definition: dimensionSet.C:474
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar j0(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)