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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "relaxation.H"
30 #include "fvm.H"
31 #include "LESModel.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace reactionRateFlameAreaModels
38 {
41  (
43  relaxation,
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 Foam::reactionRateFlameAreaModels::relaxation::relaxation
53 (
54  const word modelType,
55  const dictionary& dict,
56  const fvMesh& mesh,
57  const combustionModel& combModel
58 )
59 :
60  reactionRateFlameArea(modelType, dict, mesh, combModel),
61  correlation_(dict.optionalSubDict(typeName + "Coeffs").subDict(fuel_)),
62  C_(dict.optionalSubDict(typeName + "Coeffs").get<scalar>("C")),
63  alpha_(dict.optionalSubDict(typeName + "Coeffs").get<scalar>("alpha"))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68 
70 {}
71 
72 
73 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
74 
76 (
77  const volScalarField& sigma
78 )
79 {
80  dimensionedScalar omega0
81  (
82  "omega0",
83  dimensionSet(1, -2, -1, 0, 0, 0, 0),
84  correlation_.omega0()
85  );
86 
87  dimensionedScalar sigmaExt
88  (
89  "sigmaExt",
90  dimensionSet(0, 0, -1, 0, 0, 0, 0),
91  correlation_.sigmaExt()
92  );
93 
94  dimensionedScalar omegaMin
95  (
96  "omegaMin",
97  omega0.dimensions(),
98  1e-4
99  );
100 
101  dimensionedScalar kMin
102  (
103  "kMin",
104  sqr(dimVelocity),
105  SMALL
106  );
107 
108  const compressibleTurbulenceModel& turbulence = 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(correlation_.omega0Sigma(sigmaTotal));
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 tmp<surfaceScalarField> tphi = combModel_.phi();
128  const surfaceScalarField& phi = tphi();
129 
130  solve
131  (
132  fvm::ddt(rho, omega_)
133  + fvm::div(phi, omega_)
134  ==
135  rho*Rc*omega0
136  - fvm::SuSp(rho*(tau + Rc), omega_)
137  );
138 
139  omega_.clamp_range(0, omega0.value());
140 }
141 
142 
144 (
145  const dictionary& dict
146 )
147 {
149  {
150  coeffDict_ = dict.optionalSubDict(typeName + "Coeffs");
151  coeffDict_.readEntry("C", C_);
152  coeffDict_.readEntry("alpha", alpha_);
153  correlation_.read
154  (
155  coeffDict_.subDict(fuel_)
156  );
157  return true;
158  }
159 
160  return false;
161 }
162 
163 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void correct(const volScalarField &sigma)
Correct omega.
Definition: relaxation.C:69
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:49
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &dictProperties)
Update from dictionary.
virtual bool read(const dictionary &dictProperties)
Update properties from given dictionary.
Definition: relaxation.C:137
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
addToRunTimeSelectionTable(reactionRateFlameArea, relaxation, dictionary)
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.
const dimensionSet dimVelocity