relaxation.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-2024 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 
51 (
52  const word modelType,
53  const dictionary& dict,
54  const dictionary& coeffDict,
55  const fvMesh& mesh,
56  const combustionModel& combModel
57 )
58 :
59  reactionRateFlameArea(modelType, dict, mesh, combModel),
60  consumptionSpeed_(coeffDict.subDict(fuel_)),
61  C_(coeffDict.lookup<scalar>("C")),
62  alpha_(coeffDict.lookup<scalar>("alpha"))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
73 
75 (
76  const volScalarField& sigma
77 )
78 {
79  dimensionedScalar omega0
80  (
81  "omega0",
82  dimensionSet(1, -2, -1, 0, 0, 0, 0),
83  consumptionSpeed_.omega0()
84  );
85 
86  dimensionedScalar sigmaExt
87  (
88  "sigmaExt",
89  dimensionSet(0, 0, -1, 0, 0, 0, 0),
90  consumptionSpeed_.sigmaExt()
91  );
92 
93  dimensionedScalar omegaMin
94  (
95  "omegaMin",
96  omega0.dimensions(),
97  1e-4
98  );
99 
100  dimensionedScalar kMin
101  (
102  "kMin",
103  sqr(dimVelocity),
104  small
105  );
106 
108  combModel_.turbulence();
109 
110  // Total strain
111  const volScalarField sigmaTotal
112  (
113  sigma + alpha_*turbulence.epsilon()/(turbulence.k() + kMin)
114  );
115 
116  const volScalarField omegaInf(consumptionSpeed_.omega0Sigma(sigmaTotal));
117 
118  dimensionedScalar sigma0("sigma0", sigma.dimensions(), 0.0);
119 
120  const volScalarField tau(C_*mag(sigmaTotal));
121 
122  volScalarField Rc
123  (
124  (tau*omegaInf*(omega0 - omegaInf) + sqr(omegaMin)*sigmaExt)
125  /(sqr(omega0 - omegaInf) + sqr(omegaMin))
126  );
127 
128  const volScalarField& rho = combModel_.rho();
129  const tmp<surfaceScalarField> tphi = combModel_.phi();
130  const surfaceScalarField& phi = tphi();
131 
132  solve
133  (
134  fvm::ddt(rho, omega_)
135  + fvm::div(phi, omega_)
136  ==
137  rho*Rc*omega0
138  - fvm::SuSp(rho*(tau + Rc), omega_)
139  );
140 
141  omega_.min(omega0);
142  omega_.max(0.0);
143 }
144 
145 
147 (
148  const dictionary& dict
149 )
150 {
152  {
153  const dictionary& coeffDict = dict.optionalSubDict(typeName + "Coeffs");
154  coeffDict.lookup("C") >> C_;
155  coeffDict.lookup("alpha") >> alpha_;
156  consumptionSpeed_.read(coeffDict.subDict(fuel_));
157  return true;
158  }
159  else
160  {
161  return false;
162  }
163 }
164 
165 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
virtual Istream & read(token &)
Return next token from stream.
Definition: ITstream.C:56
Base class for combustion models.
Base class for single-phase compressible turbulence models.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:740
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
Dimension set for the base types.
Definition: dimensionSet.H:125
const dimensionSet & dimensions() const
Return const reference to dimensions.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Consumption rate per unit of flame area obtained from a relaxation equation.
Definition: relaxation.H:54
virtual void correct(const volScalarField &sigma)
Correct omega.
Definition: relaxation.C:75
virtual bool read(const dictionary &dict)
Update properties from given dictionary.
Definition: relaxation.C:147
relaxation(const word modelType, const dictionary &dict, const dictionary &coeffDict, const fvMesh &mesh, const combustionModel &combModel)
Construct from dictionary, mesh and combustion model.
Definition: relaxation.C:51
Abstract class for reaction rate per flame area unit.
virtual bool read(const dictionary &dictProperties)
Update from dictionary.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const vector tau
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
addToRunTimeSelectionTable(reactionRateFlameArea, relaxation, dictionary)
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVelocity
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))