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-2017 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.optionalSubDict(typeName + "Coeffs").subDict(fuel_)),
60  C_(readScalar(dict.optionalSubDict(typeName + "Coeffs").lookup("C"))),
61  alpha_
62  (
63  readScalar(dict.optionalSubDict(typeName + "Coeffs").lookup("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 
109  const compressibleTurbulenceModel& turbulence = combModel_.turbulence();
110 
111  // Total strain
112  const volScalarField sigmaTotal
113  (
114  sigma + alpha_*turbulence.epsilon()/(turbulence.k() + kMin)
115  );
116 
117  const volScalarField omegaInf(correlation_.omega0Sigma(sigmaTotal));
118 
119  dimensionedScalar sigma0("sigma0", sigma.dimensions(), 0.0);
120 
121  const volScalarField tau(C_*mag(sigmaTotal));
122 
123  volScalarField Rc
124  (
125  (tau*omegaInf*(omega0 - omegaInf) + sqr(omegaMin)*sigmaExt)
126  /(sqr(omega0 - omegaInf) + sqr(omegaMin))
127  );
128 
129  const volScalarField& rho = combModel_.rho();
130  const tmp<surfaceScalarField> tphi = combModel_.phi();
131  const surfaceScalarField& phi = tphi();
132 
133  solve
134  (
135  fvm::ddt(rho, omega_)
136  + fvm::div(phi, omega_)
137  ==
138  rho*Rc*omega0
139  - fvm::SuSp(rho*(tau + Rc), omega_)
140  );
141 
142  omega_.min(omega0);
143  omega_.max(0.0);
144 }
145 
146 
148 (
149  const dictionary& dict
150 )
151 {
152  if (reactionRateFlameArea::read(dict))
153  {
154  coeffDict_ = dict.optionalSubDict(typeName + "Coeffs");
155  coeffDict_.lookup("C") >> C_;
156  coeffDict_.lookup("alpha") >> alpha_;
157  correlation_.read
158  (
159  coeffDict_.subDict(fuel_)
160  );
161  return true;
162  }
163  else
164  {
165  return false;
166  }
167 }
168 
169 // ************************************************************************* //
surfaceScalarField & phi
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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:77
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
virtual bool read(const dictionary &dictProperties)
Update properties from given dictionary.
Definition: relaxation.C:148
const dimensionSet & dimensions() const
Return dimensions.
Dimension set for the base types.
Definition: dimensionSet.H:120
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
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)
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.
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
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 > &)
A class for managing temporary objects.
Definition: PtrList.H:53
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
const dimensionSet dimVelocity