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-2013 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 
79  dimensionedScalar omega0
80  (
81  "omega0",
82  dimensionSet(1, -2, -1, 0, 0, 0, 0),
83  correlation_.omega0()
84  );
85 
86  dimensionedScalar sigmaExt
87  (
88  "sigmaExt",
89  dimensionSet(0, 0, -1, 0, 0, 0, 0),
90  correlation_.sigmaExt()
91  );
92 
93  dimensionedScalar omegaMin
94  (
95  "omegaMin",
96  omega0.dimensions(),
97  1e-4
98  );
99 
100  const compressible::LESModel& lesModel =
101  omega_.db().lookupObject<compressible::LESModel>("LESProperties");
102 
103  // Total strain : resolved and sub-grid (just LES for now)
104  const volScalarField sigmaTotal
105  (
106  sigma + alpha_*lesModel.epsilon()/(lesModel.k() + lesModel.kMin())
107  );
108 
109  const volScalarField omegaInf(correlation_.omega0Sigma(sigmaTotal));
110 
111  dimensionedScalar sigma0("sigma0", sigma.dimensions(), 0.0);
112 
113  const volScalarField tau(C_*mag(sigmaTotal));
114 
115  volScalarField Rc
116  (
117  (tau*omegaInf*(omega0 - omegaInf) + sqr(omegaMin)*sigmaExt)
118  /(sqr(omega0 - omegaInf) + sqr(omegaMin))
119  );
120 
121  const volScalarField rho(combModel_.rho());
122  const surfaceScalarField phi(combModel_.phi());
123 
124  solve
125  (
126  fvm::ddt(rho, omega_)
127  + fvm::div(phi, omega_, "div(phi,omega)")
128  ==
129  rho*Rc*omega0
130  - fvm::SuSp(rho*(tau + Rc), omega_)
131  );
132 
133  omega_.min(omega0);
134  omega_.max(0.0);
135 }
136 
137 
139 (
140  const dictionary& dict
141 )
142 {
143  if (reactionRateFlameArea::read(dict))
144  {
145  coeffDict_ = dict.subDict(typeName + "Coeffs");
146  coeffDict_.lookup("C") >> C_;
147  coeffDict_.lookup("alpha") >> alpha_;
148  correlation_.read
149  (
150  coeffDict_.subDict(fuel_)
151  );
152  return true;
153  }
154  else
155  {
156  return false;
157  }
158 }
159 
160 // ************************************************************************* //
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
virtual void correct(const volScalarField &sigma)
Correct omega.
Definition: relaxation.C:74
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Abstract class for reaction rate per flame area unit.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual bool read(const dictionary &dictProperties)
Update from dictionary.
A class for handling words, derived from string.
Definition: word.H:59
Consumption rate per unit of flame area obtained from a relaxation equation.
Definition: relaxation.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
virtual Istream & read(token &)
Return next token from stream.
Definition: ITstream.C:56
Namespace for OpenFOAM.
const double e
Elementary charge.
Definition: doubleFloat.H:78
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Dimension set for the base types.
Definition: dimensionSet.H:116
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
virtual bool read(const dictionary &dictProperties)
Update properties from given dictionary.
Definition: relaxation.C:139
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
surfaceScalarField & phi
solverPerformance solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
addToRunTimeSelectionTable(reactionRateFlameArea, relaxation, dictionary)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
LESModel< EddyDiffusivity< turbulenceModel > > LESModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
Base class for combustion models.