consumptionSpeed.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-2019 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 "consumptionSpeed.H"
27 
28 /* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
29 
30 namespace Foam
31 {
33 }
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const dictionary& dict
41 )
42 : omega0_(dict.lookup<scalar>("omega0")),
43  eta_(dict.lookup<scalar>("eta")),
44  sigmaExt_(dict.lookup<scalar>("sigmaExt")),
45  omegaMin_(dict.lookup<scalar>("omegaMin"))
46 {}
47 
48 
49 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
50 
52 {}
53 
54 
55 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
56 
57 Foam::scalar Foam::consumptionSpeed::omega0Sigma
58 (
59  scalar sigma,
60  scalar a
61 ) const
62 {
63  scalar omega0 = 0.0;
64 
65  if (sigma < sigmaExt_)
66  {
67  omega0 = max
68  (
69  a*omega0_*(1.0 - exp(eta_*max(sigma, 0.0))),
70  omegaMin_
71  ) ;
72  }
73 
74  return omega0;
75 }
76 
77 
78 Foam::tmp<Foam::volScalarField> Foam::consumptionSpeed::omega0Sigma
79 (
80  const volScalarField& sigma
81 )
82 {
83  tmp<volScalarField> tomega0
84  (
86  (
87  "omega0",
88  sigma.mesh(),
90  (
91  dimensionSet(1, -2, -1, 0, 0, 0, 0),
92  0
93  )
94  )
95  );
96 
97  volScalarField& omega0 = tomega0.ref();
98 
99  volScalarField::Internal& iomega0 = omega0;
100 
101  forAll(iomega0, celli)
102  {
103  iomega0[celli] = omega0Sigma(sigma[celli], 1.0);
104  }
105 
106  volScalarField::Boundary& bomega0 = omega0.boundaryFieldRef();
107 
108  forAll(bomega0, patchi)
109  {
110  forAll(bomega0[patchi], facei)
111  {
112  bomega0[patchi][facei] =
113  omega0Sigma
114  (
115  sigma.boundaryField()[patchi][facei],
116  1.0
117  );
118  }
119  }
120 
121  return tomega0;
122 }
123 
124 
126 {
127  dict.lookup("omega0") >> omega0_ ;
128  dict.lookup("eta") >> eta_ ;
129  dict.lookup("sigmaExt") >> sigmaExt_;
130  dict.lookup("omegaMin") >> omegaMin_;
131 }
132 
133 
134 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricBoundaryField class.
Generic GeometricField class.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Correlation function for laminar consumption speed obtained from flamelet solution at increasing stra...
consumptionSpeed(const dictionary &dict)
Construct from dictionary.
void read(const dictionary &dict)
Update properties.
virtual ~consumptionSpeed()
Destructor.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
Dimension set for the base types.
Definition: dimensionSet.H:125
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
label patchi
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict