relaxation.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 "relaxation.H"
28 #include "fvm.H"
29 #include "LESModel.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace reactionRateFlameAreaModels
36 {
39  (
41  relaxation,
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 Foam::reactionRateFlameAreaModels::relaxation::relaxation
51 (
52  const word modelType,
53  const dictionary& dict,
54  const fvMesh& mesh,
55  const combustionModel& combModel
56 )
57 :
58  reactionRateFlameArea(modelType, dict, mesh, combModel),
59  correlation_(dict.subDict(typeName + "Coeffs").subDict(fuel_)),
60  C_(readScalar(dict.subDict(typeName + "Coeffs").lookup("C"))),
61  alpha_(readScalar(dict.subDict(typeName + "Coeffs").lookup("alpha")))
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
66 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
72 
74 (
75  const volScalarField& sigma
76 )
77 {
78  dimensionedScalar omega0
79  (
80  "omega0",
81  dimensionSet(1, -2, -1, 0, 0, 0, 0),
82  correlation_.omega0()
83  );
84 
85  dimensionedScalar sigmaExt
86  (
87  "sigmaExt",
88  dimensionSet(0, 0, -1, 0, 0, 0, 0),
89  correlation_.sigmaExt()
90  );
91 
92  dimensionedScalar omegaMin
93  (
94  "omegaMin",
95  omega0.dimensions(),
96  1e-4
97  );
98 
100  (
101  "kMin",
102  sqr(dimVelocity),
103  SMALL
104  );
105 
106  const compressibleTurbulenceModel& turbulence = combModel_.turbulence();
107 
108  // Total strain
109  const volScalarField sigmaTotal
110  (
111  sigma + alpha_*turbulence.epsilon()/(turbulence.k() + kMin)
112  );
113 
114  const volScalarField omegaInf(correlation_.omega0Sigma(sigmaTotal));
115 
116  dimensionedScalar sigma0("sigma0", sigma.dimensions(), 0.0);
117 
118  const volScalarField tau(C_*mag(sigmaTotal));
119 
120  volScalarField Rc
121  (
122  (tau*omegaInf*(omega0 - omegaInf) + sqr(omegaMin)*sigmaExt)
123  /(sqr(omega0 - omegaInf) + sqr(omegaMin))
124  );
125 
126  const volScalarField rho(combModel_.rho());
127  const surfaceScalarField phi(combModel_.phi());
128 
129  solve
130  (
131  fvm::ddt(rho, omega_)
132  + fvm::div(phi, omega_, "div(phi,omega)")
133  ==
134  rho*Rc*omega0
135  - fvm::SuSp(rho*(tau + Rc), omega_)
136  );
137 
138  omega_.min(omega0);
139  omega_.max(0.0);
140 }
141 
142 
144 (
145  const dictionary& dict
146 )
147 {
148  if (reactionRateFlameArea::read(dict))
149  {
150  coeffDict_ = dict.subDict(typeName + "Coeffs");
151  coeffDict_.lookup("C") >> C_;
152  coeffDict_.lookup("alpha") >> alpha_;
153  correlation_.read
154  (
155  coeffDict_.subDict(fuel_)
156  );
157  return true;
158  }
159  else
160  {
161  return false;
162  }
163 }
164 
165 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
surfaceScalarField & phi
const double e
Elementary charge.
Definition: doubleFloat.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void correct(const volScalarField &sigma)
Correct omega.
Definition: relaxation.C:74
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
Abstract class for reaction rate per flame area unit.
Consumption rate per unit of flame area obtained from a relaxation equation.
Definition: relaxation.H:51
virtual Istream & read(token &)
Return next token from stream.
Definition: ITstream.C:56
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &dictProperties)
Update from dictionary.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual bool read(const dictionary &dictProperties)
Update properties from given dictionary.
Definition: relaxation.C:144
Dimension set for the base types.
Definition: dimensionSet.H:118
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
addToRunTimeSelectionTable(reactionRateFlameArea, relaxation, dictionary)
const dimensionSet & dimensions() const
Return dimensions.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Abstract base class for turbulence models (RAS, LES and laminar).
Base class for combustion models.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
dimensioned< scalar > mag(const dimensioned< Type > &)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
const dimensionSet dimVelocity