consumptionSpeed.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-2016 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 {
32  defineTypeNameAndDebug(consumptionSpeed, 0);
33 }
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 Foam::consumptionSpeed::consumptionSpeed
39 (
40  const dictionary& dict
41 )
42 : omega0_(readScalar(dict.lookup("omega0"))),
43  eta_(readScalar(dict.lookup("eta"))),
44  sigmaExt_(readScalar(dict.lookup("sigmaExt"))),
45  omegaMin_(readScalar(dict.lookup("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  (
85  new volScalarField
86  (
87  IOobject
88  (
89  "omega0",
90  sigma.time().timeName(),
91  sigma.db(),
94  ),
95  sigma.mesh(),
97  (
98  "omega0",
99  dimensionSet(1, -2, -1, 0, 0, 0, 0),
100  0
101  )
102  )
103  );
104 
105  volScalarField& omega0 = tomega0.ref();
106 
107  volScalarField::Internal& iomega0 = omega0;
108 
109  forAll(iomega0, celli)
110  {
111  iomega0[celli] = omega0Sigma(sigma[celli], 1.0);
112  }
113 
114  volScalarField::Boundary& bomega0 = omega0.boundaryFieldRef();
115 
116  forAll(bomega0, patchi)
117  {
118  forAll(bomega0[patchi], facei)
119  {
120  bomega0[patchi][facei] =
121  omega0Sigma
122  (
123  sigma.boundaryField()[patchi][facei],
124  1.0
125  );
126  }
127  }
128 
129  return tomega0;
130 }
131 
132 
134 {
135  dict.lookup("omega0") >> omega0_ ;
136  dict.lookup("eta") >> eta_ ;
137  dict.lookup("sigmaExt") >> sigmaExt_;
138  dict.lookup("omegaMin") >> omegaMin_;
139 }
140 
141 
142 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
void read(const dictionary &dict)
Update properties.
Dimension set for the base types.
Definition: dimensionSet.H:120
dimensionedScalar exp(const dimensionedScalar &ds)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
label patchi
virtual ~consumptionSpeed()
Destructor.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:337
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:331
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576