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-2020 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 fvMesh& mesh,
55  const combustionModel& combModel
56 )
57 :
58  reactionRateFlameArea(modelType, dict, mesh, combModel),
59  correlation_(dict.optionalSubDict(typeName + "Coeffs").subDict(fuel_)),
60  C_(dict.optionalSubDict(typeName + "Coeffs").lookup<scalar>("C")),
61  alpha_
62  (
63  dict.optionalSubDict(typeName + "Coeffs").lookup<scalar>("alpha")
64  )
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 
77 (
78  const volScalarField& sigma
79 )
80 {
81  dimensionedScalar omega0
82  (
83  "omega0",
84  dimensionSet(1, -2, -1, 0, 0, 0, 0),
85  correlation_.omega0()
86  );
87 
88  dimensionedScalar sigmaExt
89  (
90  "sigmaExt",
91  dimensionSet(0, 0, -1, 0, 0, 0, 0),
92  correlation_.sigmaExt()
93  );
94 
95  dimensionedScalar omegaMin
96  (
97  "omegaMin",
98  omega0.dimensions(),
99  1e-4
100  );
101 
102  dimensionedScalar kMin
103  (
104  "kMin",
105  sqr(dimVelocity),
106  small
107  );
108 
110  combModel_.turbulence();
111 
112  // Total strain
113  const volScalarField sigmaTotal
114  (
115  sigma + alpha_*turbulence.epsilon()/(turbulence.k() + kMin)
116  );
117 
118  const volScalarField omegaInf(correlation_.omega0Sigma(sigmaTotal));
119 
120  dimensionedScalar sigma0("sigma0", sigma.dimensions(), 0.0);
121 
122  const volScalarField tau(C_*mag(sigmaTotal));
123 
124  volScalarField Rc
125  (
126  (tau*omegaInf*(omega0 - omegaInf) + sqr(omegaMin)*sigmaExt)
127  /(sqr(omega0 - omegaInf) + sqr(omegaMin))
128  );
129 
130  const volScalarField& rho = combModel_.rho();
131  const tmp<surfaceScalarField> tphi = combModel_.phi();
132  const surfaceScalarField& phi = tphi();
133 
134  solve
135  (
136  fvm::ddt(rho, omega_)
137  + fvm::div(phi, omega_)
138  ==
139  rho*Rc*omega0
140  - fvm::SuSp(rho*(tau + Rc), omega_)
141  );
142 
143  omega_.min(omega0);
144  omega_.max(0.0);
145 }
146 
147 
149 (
150  const dictionary& dict
151 )
152 {
154  {
155  coeffDict_ = dict.optionalSubDict(typeName + "Coeffs");
156  coeffDict_.lookup("C") >> C_;
157  coeffDict_.lookup("alpha") >> alpha_;
158  correlation_.read
159  (
160  coeffDict_.subDict(fuel_)
161  );
162  return true;
163  }
164  else
165  {
166  return false;
167  }
168 }
169 
170 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Base class for combustion models.
Base class for single-phase compressible turbulence models.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Dimension set for the base types.
Definition: dimensionSet.H:122
const dimensionSet & dimensions() const
Return const reference to dimensions.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
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:77
relaxation(const word modelType, const dictionary &dictCoeffs, const fvMesh &mesh, const combustionModel &combModel)
Construct from dictionary, mesh and combustion model.
Definition: relaxation.C:51
virtual bool read(const dictionary &dictProperties)
Update properties from given dictionary.
Definition: relaxation.C:149
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
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
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:105
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< scalar > mag(const dimensioned< Type > &)
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))