GeometricScalarField.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 | Copyright (C) 2011-2015 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 "GeometricScalarField.H"
27 
28 #define TEMPLATE template<template<class> class PatchField, class GeoMesh>
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37 
38 template<template<class> class PatchField, class GeoMesh>
39 void stabilise
40 (
43  const dimensioned<scalar>& ds
44 )
45 {
46  stabilise(result.internalField(), gsf.internalField(), ds.value());
47  stabilise(result.boundaryField(), gsf.boundaryField(), ds.value());
48 }
49 
50 
51 template<template<class> class PatchField, class GeoMesh>
52 tmp<GeometricField<scalar, PatchField, GeoMesh> > stabilise
53 (
55  const dimensioned<scalar>& ds
56 )
57 {
59  (
61  (
62  IOobject
63  (
64  "stabilise(" + gsf.name() + ',' + ds.name() + ')',
65  gsf.instance(),
66  gsf.db(),
69  ),
70  gsf.mesh(),
71  ds.dimensions() + gsf.dimensions()
72  )
73  );
74 
75  stabilise(tRes(), gsf, ds);
76 
77  return tRes;
78 }
79 
80 
81 template<template<class> class PatchField, class GeoMesh>
82 tmp<GeometricField<scalar, PatchField, GeoMesh> > stabilise
83 (
85  const dimensioned<scalar>& ds
86 )
87 {
89 
91  (
93  (
94  tgsf,
95  "stabilise(" + gsf.name() + ',' + ds.name() + ')',
96  ds.dimensions() + gsf.dimensions()
97  )
98  );
99 
100  stabilise(tRes(), gsf, ds);
101 
103 
104  return tRes;
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, +, '+', add)
111 BINARY_TYPE_OPERATOR(scalar, scalar, scalar, -, '-', subtract)
112 
113 BINARY_OPERATOR(scalar, scalar, scalar, *, '*', multiply)
114 BINARY_OPERATOR(scalar, scalar, scalar, /, '|', divide)
115 
116 BINARY_TYPE_OPERATOR_SF(scalar, scalar, scalar, /, '|', divide)
117 
118 
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120 
121 template<template<class> class PatchField, class GeoMesh>
122 void pow
123 (
127 )
128 {
129  pow(Pow.internalField(), gsf1.internalField(), gsf2.internalField());
130  pow(Pow.boundaryField(), gsf1.boundaryField(), gsf2.boundaryField());
131 }
132 
133 
134 template<template<class> class PatchField, class GeoMesh>
135 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
136 (
139 )
140 {
142  (
144  (
145  IOobject
146  (
147  "pow(" + gsf1.name() + ',' + gsf2.name() + ')',
148  gsf1.instance(),
149  gsf1.db(),
152  ),
153  gsf1.mesh(),
154  pow
155  (
156  gsf1.dimensions(),
157  dimensionedScalar("1", gsf2.dimensions(), 1.0)
158  )
159  )
160  );
161 
162  pow(tPow(), gsf1, gsf2);
163 
164  return tPow;
165 }
166 
167 
168 template<template<class> class PatchField, class GeoMesh>
169 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
170 (
173 )
174 {
175  const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
176 
178  (
180  (
181  tgsf1,
182  "pow(" + gsf1.name() + ',' + gsf2.name() + ')',
183  pow
184  (
185  gsf1.dimensions(),
186  dimensionedScalar("1", gsf2.dimensions(), 1.0)
187  )
188  )
189  );
190 
191  pow(tPow(), gsf1, gsf2);
192 
194 
195  return tPow;
196 }
197 
198 
199 template<template<class> class PatchField, class GeoMesh>
200 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
201 (
204 )
205 {
206  const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
207 
209  (
211  (
212  tgsf2,
213  "pow(" + gsf1.name() + ',' + gsf2.name() + ')',
214  pow
215  (
216  gsf1.dimensions(),
217  dimensionedScalar("1", gsf2.dimensions(), 1.0)
218  )
219  )
220  );
221 
222  pow(tPow(), gsf1, gsf2);
223 
225 
226  return tPow;
227 }
228 
229 template<template<class> class PatchField, class GeoMesh>
230 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
231 (
234 )
235 {
236  const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
237  const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
238 
240  (
242  <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::New
243  (
244  tgsf1,
245  tgsf2,
246  "pow(" + gsf1.name() + ',' + gsf2.name() + ')',
247  pow
248  (
249  gsf1.dimensions(),
250  dimensionedScalar("1", gsf2.dimensions(), 1.0)
251  )
252  )
253  );
254 
255  pow(tPow(), gsf1, gsf2);
256 
258  <scalar, scalar, scalar, scalar, PatchField, GeoMesh>
259  ::clear(tgsf1, tgsf2);
260 
261  return tPow;
262 }
263 
264 
265 template<template<class> class PatchField, class GeoMesh>
266 void pow
267 (
270  const dimensioned<scalar>& ds
271 )
272 {
273  pow(tPow.internalField(), gsf.internalField(), ds.value());
274  pow(tPow.boundaryField(), gsf.boundaryField(), ds.value());
275 }
276 
277 
278 template<template<class> class PatchField, class GeoMesh>
279 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
280 (
282  const dimensionedScalar& ds
283 )
284 {
286  (
288  (
289  IOobject
290  (
291  "pow(" + gsf.name() + ',' + ds.name() + ')',
292  gsf.instance(),
293  gsf.db(),
296  ),
297  gsf.mesh(),
298  pow(gsf.dimensions(), ds)
299  )
300  );
301 
302  pow(tPow(), gsf, ds);
303 
304  return tPow;
305 }
306 
307 template<template<class> class PatchField, class GeoMesh>
308 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
309 (
311  const dimensionedScalar& ds
312 )
313 {
315 
317  (
319  (
320  tgsf,
321  "pow(" + gsf.name() + ',' + ds.name() + ')',
322  pow(gsf.dimensions(), ds)
323  )
324  );
325 
326  pow(tPow(), gsf, ds);
327 
329 
330  return tPow;
331 }
332 
333 template<template<class> class PatchField, class GeoMesh>
334 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
335 (
337  const scalar& s
338 )
339 {
340  return pow(gsf, dimensionedScalar(s));
341 }
342 
343 template<template<class> class PatchField, class GeoMesh>
344 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
345 (
347  const scalar& s
348 )
349 {
350  return pow(tgsf, dimensionedScalar(s));
351 }
352 
353 
354 template<template<class> class PatchField, class GeoMesh>
355 void pow
356 (
358  const dimensioned<scalar>& ds,
360 )
361 {
362  pow(tPow.internalField(), ds.value(), gsf.internalField());
363  pow(tPow.boundaryField(), ds.value(), gsf.boundaryField());
364 }
365 
366 
367 template<template<class> class PatchField, class GeoMesh>
368 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
369 (
370  const dimensionedScalar& ds,
372 )
373 {
375  (
377  (
378  IOobject
379  (
380  "pow(" + ds.name() + ',' + gsf.name() + ')',
381  gsf.instance(),
382  gsf.db(),
385  ),
386  gsf.mesh(),
387  pow(ds, gsf.dimensions())
388  )
389  );
390 
391  pow(tPow(), ds, gsf);
392 
393  return tPow;
394 }
395 
396 
397 template<template<class> class PatchField, class GeoMesh>
398 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
399 (
400  const dimensionedScalar& ds,
402 )
403 {
405 
407  (
409  (
410  tgsf,
411  "pow(" + ds.name() + ',' + gsf.name() + ')',
412  pow(ds, gsf.dimensions())
413  )
414  );
415 
416  pow(tPow(), ds, gsf);
417 
419 
420  return tPow;
421 }
422 
423 template<template<class> class PatchField, class GeoMesh>
424 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
425 (
426  const scalar& s,
428 )
429 {
430  return pow(dimensionedScalar(s), gsf);
431 }
432 
433 template<template<class> class PatchField, class GeoMesh>
434 tmp<GeometricField<scalar, PatchField, GeoMesh> > pow
435 (
436  const scalar& s,
438 )
439 {
440  return pow(dimensionedScalar(s), tgsf);
441 }
442 
443 
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
445 
446 template<template<class> class PatchField, class GeoMesh>
447 void atan2
448 (
452 )
453 {
454  atan2(Atan2.internalField(), gsf1.internalField(), gsf2.internalField());
455  atan2(Atan2.boundaryField(), gsf1.boundaryField(), gsf2.boundaryField());
456 }
457 
458 
459 template<template<class> class PatchField, class GeoMesh>
460 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
461 (
464 )
465 {
467  (
469  (
470  IOobject
471  (
472  "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
473  gsf1.instance(),
474  gsf1.db(),
477  ),
478  gsf1.mesh(),
479  atan2(gsf1.dimensions(), gsf2.dimensions())
480  )
481  );
482 
483  atan2(tAtan2(), gsf1, gsf2);
484 
485  return tAtan2;
486 }
487 
488 
489 template<template<class> class PatchField, class GeoMesh>
490 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
491 (
494 )
495 {
496  const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
497 
499  (
501  (
502  tgsf1,
503  "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
504  atan2(gsf1.dimensions(), gsf2.dimensions())
505  )
506  );
507 
508  atan2(tAtan2(), gsf1, gsf2);
509 
511 
512  return tAtan2;
513 }
514 
515 
516 template<template<class> class PatchField, class GeoMesh>
517 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
518 (
521 )
522 {
523  const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
524 
526  (
528  (
529  tgsf2,
530  "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
531  atan2( gsf1.dimensions(), gsf2.dimensions())
532  )
533  );
534 
535  atan2(tAtan2(), gsf1, gsf2);
536 
538 
539  return tAtan2;
540 }
541 
542 template<template<class> class PatchField, class GeoMesh>
543 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
544 (
547 )
548 {
549  const GeometricField<scalar, PatchField, GeoMesh>& gsf1 = tgsf1();
550  const GeometricField<scalar, PatchField, GeoMesh>& gsf2 = tgsf2();
551 
553  (
555  <scalar, scalar, scalar, scalar, PatchField, GeoMesh>::New
556  (
557  tgsf1,
558  tgsf2,
559  "atan2(" + gsf1.name() + ',' + gsf2.name() + ')',
560  atan2(gsf1.dimensions(), gsf2.dimensions())
561  )
562  );
563 
564  atan2(tAtan2(), gsf1, gsf2);
565 
567  <scalar, scalar, scalar, scalar, PatchField, GeoMesh>
568  ::clear(tgsf1, tgsf2);
569 
570  return tAtan2;
571 }
572 
573 
574 template<template<class> class PatchField, class GeoMesh>
575 void atan2
576 (
579  const dimensioned<scalar>& ds
580 )
581 {
582  atan2(tAtan2.internalField(), gsf.internalField(), ds.value());
583  atan2(tAtan2.boundaryField(), gsf.boundaryField(), ds.value());
584 }
585 
586 
587 template<template<class> class PatchField, class GeoMesh>
588 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
589 (
591  const dimensionedScalar& ds
592 )
593 {
595  (
597  (
598  IOobject
599  (
600  "atan2(" + gsf.name() + ',' + ds.name() + ')',
601  gsf.instance(),
602  gsf.db(),
605  ),
606  gsf.mesh(),
607  atan2(gsf.dimensions(), ds)
608  )
609  );
610 
611  atan2(tAtan2(), gsf, ds);
612 
613  return tAtan2;
614 }
615 
616 template<template<class> class PatchField, class GeoMesh>
617 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
618 (
620  const dimensionedScalar& ds
621 )
622 {
624 
626  (
628  (
629  tgsf,
630  "atan2(" + gsf.name() + ',' + ds.name() + ')',
631  atan2(gsf.dimensions(), ds)
632  )
633  );
634 
635  atan2(tAtan2(), gsf, ds);
636 
638 
639  return tAtan2;
640 }
641 
642 template<template<class> class PatchField, class GeoMesh>
643 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
644 (
646  const scalar& s
647 )
648 {
649  return atan2(gsf, dimensionedScalar(s));
650 }
651 
652 template<template<class> class PatchField, class GeoMesh>
653 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
654 (
656  const scalar& s
657 )
658 {
659  return atan2(tgsf, dimensionedScalar(s));
660 }
661 
662 
663 template<template<class> class PatchField, class GeoMesh>
664 void atan2
665 (
667  const dimensioned<scalar>& ds,
669 )
670 {
671  atan2(tAtan2.internalField(), ds.value(), gsf.internalField());
672  atan2(tAtan2.boundaryField(), ds.value(), gsf.boundaryField());
673 }
674 
675 
676 template<template<class> class PatchField, class GeoMesh>
677 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
678 (
679  const dimensionedScalar& ds,
681 )
682 {
684  (
686  (
687  IOobject
688  (
689  "atan2(" + ds.name() + ',' + gsf.name() + ')',
690  gsf.instance(),
691  gsf.db(),
694  ),
695  gsf.mesh(),
696  atan2(ds, gsf.dimensions())
697  )
698  );
699 
700  atan2(tAtan2(), ds, gsf);
701 
702  return tAtan2;
703 }
704 
705 
706 template<template<class> class PatchField, class GeoMesh>
707 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
708 (
709  const dimensionedScalar& ds,
711 )
712 {
714 
716  (
718  (
719  tgsf,
720  "atan2(" + ds.name() + ',' + gsf.name() + ')',
721  atan2(ds, gsf.dimensions())
722  )
723  );
724 
725  atan2(tAtan2(), ds, gsf);
726 
728 
729  return tAtan2;
730 }
731 
732 template<template<class> class PatchField, class GeoMesh>
733 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
734 (
735  const scalar& s,
737 )
738 {
739  return atan2(dimensionedScalar(s), gsf);
740 }
741 
742 template<template<class> class PatchField, class GeoMesh>
743 tmp<GeometricField<scalar, PatchField, GeoMesh> > atan2
744 (
745  const scalar& s,
747 )
748 {
749  return atan2(dimensionedScalar(s), tgsf);
750 }
751 
752 
753 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
754 
755 UNARY_FUNCTION(scalar, scalar, pow3, pow3)
756 UNARY_FUNCTION(scalar, scalar, pow4, pow4)
757 UNARY_FUNCTION(scalar, scalar, pow5, pow5)
758 UNARY_FUNCTION(scalar, scalar, pow6, pow6)
759 UNARY_FUNCTION(scalar, scalar, pow025, pow025)
760 UNARY_FUNCTION(scalar, scalar, sqrt, sqrt)
761 UNARY_FUNCTION(scalar, scalar, cbrt, cbrt)
762 UNARY_FUNCTION(scalar, scalar, sign, sign)
763 UNARY_FUNCTION(scalar, scalar, pos, pos)
764 UNARY_FUNCTION(scalar, scalar, neg, neg)
765 UNARY_FUNCTION(scalar, scalar, posPart, posPart)
766 UNARY_FUNCTION(scalar, scalar, negPart, negPart)
767 
768 UNARY_FUNCTION(scalar, scalar, exp, trans)
769 UNARY_FUNCTION(scalar, scalar, log, trans)
770 UNARY_FUNCTION(scalar, scalar, log10, trans)
771 UNARY_FUNCTION(scalar, scalar, sin, trans)
772 UNARY_FUNCTION(scalar, scalar, cos, trans)
773 UNARY_FUNCTION(scalar, scalar, tan, trans)
774 UNARY_FUNCTION(scalar, scalar, asin, trans)
775 UNARY_FUNCTION(scalar, scalar, acos, trans)
776 UNARY_FUNCTION(scalar, scalar, atan, trans)
777 UNARY_FUNCTION(scalar, scalar, sinh, trans)
778 UNARY_FUNCTION(scalar, scalar, cosh, trans)
779 UNARY_FUNCTION(scalar, scalar, tanh, trans)
780 UNARY_FUNCTION(scalar, scalar, asinh, trans)
781 UNARY_FUNCTION(scalar, scalar, acosh, trans)
782 UNARY_FUNCTION(scalar, scalar, atanh, trans)
783 UNARY_FUNCTION(scalar, scalar, erf, trans)
784 UNARY_FUNCTION(scalar, scalar, erfc, trans)
785 UNARY_FUNCTION(scalar, scalar, lgamma, trans)
786 UNARY_FUNCTION(scalar, scalar, j0, trans)
787 UNARY_FUNCTION(scalar, scalar, j1, trans)
788 UNARY_FUNCTION(scalar, scalar, y0, trans)
789 UNARY_FUNCTION(scalar, scalar, y1, trans)
790 
791 
792 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
793 
794 #define BesselFunc(func) \
795  \
796 template<template<class> class PatchField, class GeoMesh> \
797 void func \
798 ( \
799  GeometricField<scalar, PatchField, GeoMesh>& gsf, \
800  const int n, \
801  const GeometricField<scalar, PatchField, GeoMesh>& gsf1 \
802 ) \
803 { \
804  func(gsf.internalField(), n, gsf1.internalField()); \
805  func(gsf.boundaryField(), n, gsf1.boundaryField()); \
806 } \
807  \
808 template<template<class> class PatchField, class GeoMesh> \
809 tmp<GeometricField<scalar, PatchField, GeoMesh> > func \
810 ( \
811  const int n, \
812  const GeometricField<scalar, PatchField, GeoMesh>& gsf \
813 ) \
814 { \
815  if (!gsf.dimensions().dimensionless()) \
816  { \
817  FatalErrorInFunction \
818  << "gsf not dimensionless" \
819  << abort(FatalError); \
820  } \
821  \
822  tmp<GeometricField<scalar, PatchField, GeoMesh> > tFunc \
823  ( \
824  new GeometricField<scalar, PatchField, GeoMesh> \
825  ( \
826  IOobject \
827  ( \
828  #func "(" + gsf.name() + ')', \
829  gsf.instance(), \
830  gsf.db(), \
831  IOobject::NO_READ, \
832  IOobject::NO_WRITE \
833  ), \
834  gsf.mesh(), \
835  dimless \
836  ) \
837  ); \
838  \
839  func(tFunc(), n, gsf); \
840  \
841  return tFunc; \
842 } \
843  \
844 template<template<class> class PatchField, class GeoMesh> \
845 tmp<GeometricField<scalar, PatchField, GeoMesh> > func \
846 ( \
847  const int n, \
848  const tmp<GeometricField<scalar, PatchField, GeoMesh> >& tgsf \
849 ) \
850 { \
851  const GeometricField<scalar, PatchField, GeoMesh>& gsf = tgsf(); \
852  \
853  if (!gsf.dimensions().dimensionless()) \
854  { \
855  FatalErrorInFunction \
856  << " : gsf not dimensionless" \
857  << abort(FatalError); \
858  } \
859  \
860  tmp<GeometricField<scalar, PatchField, GeoMesh> > tFunc \
861  ( \
862  reuseTmpGeometricField<scalar, scalar, PatchField, GeoMesh>::New \
863  ( \
864  tgsf, \
865  #func "(" + gsf.name() + ')', \
866  dimless \
867  ) \
868  ); \
869  \
870  func(tFunc(), n, gsf); \
871  \
872  reuseTmpGeometricField<scalar, scalar, PatchField, GeoMesh> \
873  ::clear(tgsf); \
874  \
875  return tFunc; \
876 }
877 
880 
881 #undef BesselFunc
882 
883 
884 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
885 
886 } // End namespace Foam
887 
888 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
889 
890 #include "undefFieldFunctionsM.H"
891 
892 // ************************************************************************* //
BINARY_TYPE_OPERATOR_SF
#define BINARY_TYPE_OPERATOR_SF(ReturnType, Type1, Type2, Op, OpName, OpFunc)
Definition: DimensionedFieldFunctionsM.C:531
BINARY_OPERATOR
#define BINARY_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
Definition: DimensionedFieldFunctionsM.C:417
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
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:871
Foam::tan
dimensionedScalar tan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:257
UNARY_FUNCTION
#define UNARY_FUNCTION(ReturnType, Type1, Func, Dfunc)
Definition: DimensionedFieldFunctionsM.C:30
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::cosh
dimensionedScalar cosh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
Foam::y1
dimensionedScalar y1(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:273
clear
UEqn clear()
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
Foam::jn
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
Definition: dimensionedScalar.C:296
Foam::posPart
dimensionedScalar posPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:212
Foam::atan2
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
Definition: dimensionedScalar.C:303
BINARY_TYPE_OPERATOR
#define BINARY_TYPE_OPERATOR(ReturnType, Type1, Type2, Op, OpName, OpFunc)
Definition: DimensionedFieldFunctionsM.C:684
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:235
undefFieldFunctionsM.H
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
BesselFunc
#define BesselFunc(func)
Definition: GeometricScalarField.C:794
GeometricScalarField.H
Scalar specific part of the implementation of GeometricField.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:179
Foam::erf
dimensionedScalar erf(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:267
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
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:131
Foam::atanh
dimensionedScalar atanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:266
Foam::lgamma
dimensionedScalar lgamma(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:269
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:98
Foam::pow6
dimensionedScalar pow6(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:120
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::stabilise
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
Definition: DimensionedScalarField.C:40
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::reuseTmpGeometricField
Definition: GeometricFieldReuseFunctions.H:43
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:254
GeometricFieldFunctionsM.C
Foam::y0
dimensionedScalar y0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:272
Foam::erfc
dimensionedScalar erfc(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::asinh
dimensionedScalar asinh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::reuseTmpGeometricField::clear
static void clear(const tmp< GeometricField< Type1, PatchField, GeoMesh > > &tdf1)
Definition: GeometricFieldReuseFunctions.H:73
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
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:870
Foam::pow5
dimensionedScalar pow5(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:109
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
Foam::dimensioned< scalar >
Foam::GeoMesh
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
s
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::yn
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
Definition: dimensionedScalar.C:297
Foam::acosh
dimensionedScalar acosh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::negPart
dimensionedScalar negPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:223
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::j0
dimensionedScalar j0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:270
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:259
Foam::atan
dimensionedScalar atan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:260
Foam::reuseTmpTmpGeometricField
Definition: GeometricFieldReuseFunctions.H:144
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:248
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:153
Foam::trans
dimensionSet trans(const dimensionSet &)
Definition: dimensionSet.C:438
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::j1
dimensionedScalar j1(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:271
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:201
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:258
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190
Foam::sinh
dimensionedScalar sinh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261