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>
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37 
38 template<class GeoMesh>
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().primitiveFieldRef(), dsf.primitiveField(), 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().primitiveFieldRef(), dsf.primitiveField(), 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
129  (
130  tPow.ref().primitiveFieldRef(),
131  dsf1.primitiveField(),
132  dsf2.primitiveField()
133  );
134 
135  return tPow;
136 }
137 
138 
139 template<class GeoMesh>
141 (
144 )
145 {
146  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
147 
148  if (!dsf1.dimensions().dimensionless())
149  {
151  << "Base field is not dimensionless: " << dsf1.dimensions()
152  << exit(FatalError);
153  }
154 
155  if (!dsf2.dimensions().dimensionless())
156  {
158  << "Exponent field is not dimensionless: " << dsf2.dimensions()
159  << exit(FatalError);
160  }
161 
163  (
164  tdsf1,
165  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
166  dimless
167  );
168 
169  pow
170  (
171  tPow.ref().primitiveFieldRef(),
172  dsf1.primitiveField(),
173  dsf2.primitiveField()
174  );
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
213  (
214  tPow.ref().primitiveFieldRef(),
215  dsf1.primitiveField(),
216  dsf2.primitiveField()
217  );
218 
219  tdsf2.clear();
220 
221  return tPow;
222 }
223 
224 
225 template<class GeoMesh>
227 (
230 )
231 {
232  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
233  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
234 
235  if (!dsf1.dimensions().dimensionless())
236  {
238  << "Base field is not dimensionless: " << dsf1.dimensions()
239  << exit(FatalError);
240  }
241 
242  if (!dsf2.dimensions().dimensionless())
243  {
245  << "Exponent field is not dimensionless: " << dsf2.dimensions()
246  << exit(FatalError);
247  }
248 
251  (
252  tdsf1,
253  tdsf2,
254  "pow(" + dsf1.name() + ',' + dsf2.name() + ')',
255  dimless
256  );
257 
258  pow
259  (
260  tPow.ref().primitiveFieldRef(),
261  dsf1.primitiveField(),
262  dsf2.primitiveField()
263  );
264 
265  tdsf1.clear();
266  tdsf2.clear();
267 
268  return tPow;
269 }
270 
271 
272 template<class GeoMesh>
274 (
276  const dimensionedScalar& ds
277 )
278 {
279  if (!ds.dimensions().dimensionless())
280  {
282  << "Exponent is not dimensionless: " << ds.dimensions()
283  << exit(FatalError);
284  }
285 
287  (
289  (
290  "pow(" + dsf.name() + ',' + ds.name() + ')',
291  dsf.mesh(),
292  pow(dsf.dimensions(), ds)
293  )
294  );
295 
296  pow(tPow.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
297 
298  return tPow;
299 }
300 
301 
302 template<class GeoMesh>
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 
316  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
317 
319  (
320  tdsf,
321  "pow(" + dsf.name() + ',' + ds.name() + ')',
322  pow(dsf.dimensions(), ds)
323  );
324 
325  pow(tPow.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
326 
327  tdsf.clear();
328 
329  return tPow;
330 }
331 
332 
333 template<class GeoMesh>
335 (
337  const scalar& s
338 )
339 {
340  return pow(dsf, dimensionedScalar(s));
341 }
342 
343 
344 template<class GeoMesh>
346 (
348  const scalar& s
349 )
350 {
351  return pow(tdsf, dimensionedScalar(s));
352 }
353 
354 
355 template<class GeoMesh>
357 (
358  const dimensionedScalar& ds,
360 )
361 {
362  if (!ds.dimensions().dimensionless())
363  {
365  << "Base scalar is not dimensionless: " << ds.dimensions()
366  << exit(FatalError);
367  }
368 
369  if (!dsf.dimensions().dimensionless())
370  {
372  << "Exponent field is not dimensionless: " << dsf.dimensions()
373  << exit(FatalError);
374  }
375 
377  (
379  (
380  "pow(" + ds.name() + ',' + dsf.name() + ')',
381  dsf.mesh(),
382  dimless
383  )
384  );
385 
386  pow(tPow.ref().primitiveFieldRef(), ds.value(), dsf.primitiveField());
387 
388  return tPow;
389 }
390 
391 
392 template<class GeoMesh>
394 (
395  const dimensionedScalar& ds,
397 )
398 {
399  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
400 
401  if (!ds.dimensions().dimensionless())
402  {
404  << "Base scalar is not dimensionless: " << ds.dimensions()
405  << exit(FatalError);
406  }
407 
408  if (!dsf.dimensions().dimensionless())
409  {
411  << "Exponent field is not dimensionless: " << dsf.dimensions()
412  << exit(FatalError);
413  }
414 
416  (
417  tdsf,
418  "pow(" + ds.name() + ',' + dsf.name() + ')',
419  dimless
420  );
421 
422  pow(tPow.ref().primitiveFieldRef(), ds.value(), dsf.primitiveField());
423 
424  tdsf.clear();
425 
426  return tPow;
427 }
428 
429 template<class GeoMesh>
431 (
432  const scalar& s,
434 )
435 {
436  return pow(dimensionedScalar(s), dsf);
437 }
438 
439 template<class GeoMesh>
441 (
442  const scalar& s,
444 )
445 {
446  return pow(dimensionedScalar(s), tdsf);
447 }
448 
449 
450 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
451 
452 template<class GeoMesh>
454 (
457 )
458 {
460  (
462  (
463  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
464  dsf1.mesh(),
465  atan2(dsf1.dimensions(), dsf2.dimensions())
466  )
467  );
468 
469  atan2
470  (
471  tAtan2.ref().primitiveFieldRef(),
472  dsf1.primitiveField(),
473  dsf2.primitiveField()
474  );
475 
476  return tAtan2;
477 }
478 
479 
480 template<class GeoMesh>
482 (
485 )
486 {
487  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
488 
490  (
491  tdsf1,
492  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
493  atan2(dsf1.dimensions(), dsf2.dimensions())
494  );
495 
496  atan2
497  (
498  tAtan2.ref().primitiveFieldRef(),
499  dsf1.primitiveField(),
500  dsf2.primitiveField()
501  );
502 
503  tdsf1.clear();
504 
505  return tAtan2;
506 }
507 
508 
509 template<class GeoMesh>
511 (
514 )
515 {
516  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
517 
519  (
520  tdsf2,
521  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
522  atan2(dsf1.dimensions(), dsf2.dimensions())
523  );
524 
525  atan2
526  (
527  tAtan2.ref().primitiveFieldRef(),
528  dsf1.primitiveField(),
529  dsf2.primitiveField()
530  );
531 
532  tdsf2.clear();
533 
534  return tAtan2;
535 }
536 
537 template<class GeoMesh>
539 (
542 )
543 {
544  const DimensionedField<scalar, GeoMesh>& dsf1 = tdsf1();
545  const DimensionedField<scalar, GeoMesh>& dsf2 = tdsf2();
546 
549  (
550  tdsf1,
551  tdsf2,
552  "atan2(" + dsf1.name() + ',' + dsf2.name() + ')',
553  atan2(dsf1.dimensions(), dsf2.dimensions())
554  );
555 
556  atan2
557  (
558  tAtan2.ref().primitiveFieldRef(),
559  dsf1.primitiveField(),
560  dsf2.primitiveField()
561  );
562 
563  tdsf1.clear();
564  tdsf2.clear();
565 
566  return tAtan2;
567 }
568 
569 
570 template<class GeoMesh>
572 (
574  const dimensionedScalar& ds
575 )
576 {
578  (
580  (
581  "atan2(" + dsf.name() + ',' + ds.name() + ')',
582  dsf.mesh(),
583  atan2(dsf.dimensions(), ds)
584  )
585  );
586 
587  atan2(tAtan2.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
588 
589  return tAtan2;
590 }
591 
592 template<class GeoMesh>
594 (
596  const dimensionedScalar& ds
597 )
598 {
599  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
600 
602  (
603  tdsf,
604  "atan2(" + dsf.name() + ',' + ds.name() + ')',
605  atan2(dsf.dimensions(), ds)
606  );
607 
608  atan2(tAtan2.ref().primitiveFieldRef(), dsf.primitiveField(), ds.value());
609 
610  tdsf.clear();
611 
612  return tAtan2;
613 }
614 
615 template<class GeoMesh>
617 (
619  const scalar& s
620 )
621 {
622  return atan2(dsf, dimensionedScalar(s));
623 }
624 
625 template<class GeoMesh>
627 (
629  const scalar& s
630 )
631 {
632  return atan2(tdsf, dimensionedScalar(s));
633 }
634 
635 
636 template<class GeoMesh>
638 (
639  const dimensionedScalar& ds,
641 )
642 {
644  (
646  (
647  "atan2(" + ds.name() + ',' + dsf.name() + ')',
648  dsf.mesh(),
649  atan2(ds, dsf.dimensions())
650  )
651  );
652 
653  atan2(tAtan2.ref().primitiveFieldRef(), ds.value(), dsf.primitiveField());
654 
655  return tAtan2;
656 }
657 
658 
659 template<class GeoMesh>
661 (
662  const dimensionedScalar& ds,
664 )
665 {
666  const DimensionedField<scalar, GeoMesh>& dsf = tdsf();
667 
669  (
670  tdsf,
671  "atan2(" + ds.name() + ',' + dsf.name() + ')',
672  atan2(ds, dsf.dimensions())
673  );
674 
675  atan2(tAtan2.ref().primitiveFieldRef(), ds.value(), dsf.primitiveField());
676 
677  tdsf.clear();
678 
679  return tAtan2;
680 }
681 
682 template<class GeoMesh>
684 (
685  const scalar& s,
687 )
688 {
689  return atan2(dimensionedScalar(s), dsf);
690 }
691 
692 template<class GeoMesh>
694 (
695  const scalar& s,
697 )
698 {
699  return atan2(dimensionedScalar(s), tdsf);
700 }
701 
702 
703 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
704 
705 UNARY_FUNCTION(scalar, scalar, pow3, pow3)
706 UNARY_FUNCTION(scalar, scalar, pow4, pow4)
707 UNARY_FUNCTION(scalar, scalar, pow5, pow5)
708 UNARY_FUNCTION(scalar, scalar, pow6, pow6)
709 UNARY_FUNCTION(scalar, scalar, pow025, pow025)
710 UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
711 UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
712 UNARY_FUNCTION(scalar, scalar, sign, sign)
713 UNARY_FUNCTION(scalar, scalar, pos, pos)
714 UNARY_FUNCTION(scalar, scalar, pos0, pos0)
715 UNARY_FUNCTION(scalar, scalar, neg, neg)
716 UNARY_FUNCTION(scalar, scalar, neg0, neg0)
717 UNARY_FUNCTION(scalar, scalar, posPart, posPart)
718 UNARY_FUNCTION(scalar, scalar, negPart, negPart)
719 
720 UNARY_FUNCTION(scalar, scalar, exp, trans)
721 UNARY_FUNCTION(scalar, scalar, log, trans)
722 UNARY_FUNCTION(scalar, scalar, log10, trans)
723 UNARY_FUNCTION(scalar, scalar, sin, trans)
724 UNARY_FUNCTION(scalar, scalar, cos, trans)
725 UNARY_FUNCTION(scalar, scalar, tan, trans)
726 UNARY_FUNCTION(scalar, scalar, asin, trans)
727 UNARY_FUNCTION(scalar, scalar, acos, trans)
728 UNARY_FUNCTION(scalar, scalar, atan, trans)
729 UNARY_FUNCTION(scalar, scalar, sinh, trans)
730 UNARY_FUNCTION(scalar, scalar, cosh, trans)
731 UNARY_FUNCTION(scalar, scalar, tanh, trans)
732 UNARY_FUNCTION(scalar, scalar, asinh, trans)
733 UNARY_FUNCTION(scalar, scalar, acosh, trans)
734 UNARY_FUNCTION(scalar, scalar, atanh, trans)
735 UNARY_FUNCTION(scalar, scalar, erf, trans)
736 UNARY_FUNCTION(scalar, scalar, erfc, trans)
737 UNARY_FUNCTION(scalar, scalar, lgamma, trans)
738 UNARY_FUNCTION(scalar, scalar, j0, trans)
739 UNARY_FUNCTION(scalar, scalar, j1, trans)
740 UNARY_FUNCTION(scalar, scalar, y0, trans)
741 UNARY_FUNCTION(scalar, scalar, y1, trans)
742 
743 
744 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
745 
746 #define BesselFunc(func) \
747  \
748 template<class GeoMesh> \
749 tmp<DimensionedField<scalar, GeoMesh>> func \
750 ( \
751  const int n, \
752  const DimensionedField<scalar, GeoMesh>& dsf \
753 ) \
754 { \
755  if (!dsf.dimensions().dimensionless()) \
756  { \
757  FatalErrorInFunction \
758  << "dsf not dimensionless" \
759  << abort(FatalError); \
760  } \
761  \
762  tmp<DimensionedField<scalar, GeoMesh>> tFunc \
763  ( \
764  DimensionedField<scalar, GeoMesh>::New \
765  ( \
766  #func "(" + name(n) + ',' + dsf.name() + ')', \
767  dsf.mesh(), \
768  dimless \
769  ) \
770  ); \
771  \
772  func(tFunc.ref().primitiveFieldRef(), n, dsf.primitiveField()); \
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().primitiveFieldRef(), n, dsf.primitiveField()); \
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 // ************************************************************************* //
#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 BesselFunc(func)
Scalar specific part of the implementation of DimensionedField.
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 Field< Type > & primitiveField() const
Return a const-reference to the primitive field.
const Mesh & mesh() const
Return mesh.
const word & name() const
Return name.
Definition: IOobject.H:310
bool dimensionless() const
Return true if it is dimensionless.
Definition: dimensionSet.C:110
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.
static tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< Type1, GeoMesh >> &tdf1, const tmp< DimensionedField< Type2, GeoMesh >> &tdf2, const word &name, const dimensionSet &dimensions)
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:181
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
dimensionedScalar pow6(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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)
dimensionedScalar tan(const dimensionedScalar &ds)
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &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)
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)
UNARY_FUNCTION(Type, Type, cmptMag, cmptMag)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
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)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
error FatalError
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar atan(const dimensionedScalar &ds)
dimensionSet trans(const dimensionSet &)
Definition: dimensionSet.C:477
void multiply(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)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar j0(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)