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