1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 namespace laminarFlameSpeedModels
35 defineTypeNameAndDebug(Gulders, 0);
37 addToRunTimeSelectionTable
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 Foam::laminarFlameSpeedModels::Gulders::Gulders
51 const dictionary& dict,
52 const hhuCombustionThermo& ct
55 laminarFlameSpeed(dict, ct),
57 coeffsDict_(dict.subDict(typeName + "Coeffs").subDict(fuel_)),
58 W_(readScalar(coeffsDict_.lookup("W"))),
59 eta_(readScalar(coeffsDict_.lookup("eta"))),
60 xi_(readScalar(coeffsDict_.lookup("xi"))),
61 f_(readScalar(coeffsDict_.lookup("f"))),
62 alpha_(readScalar(coeffsDict_.lookup("alpha"))),
63 beta_(readScalar(coeffsDict_.lookup("beta")))
67 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 Foam::laminarFlameSpeedModels::Gulders::~Gulders()
73 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::SuRef
82 return W_*pow(phi, eta_)*exp(-xi_*sqr(phi - 1.075));
91 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
99 static const scalar Tref = 300.0;
100 static const scalar pRef = 1.013e5;
102 return SuRef(phi)*pow((Tu/Tref), alpha_)*pow((p/pRef), beta_)*(1 - f_*Yres);
106 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
108 const volScalarField& p,
109 const volScalarField& Tu,
113 tmp<volScalarField> tSu0
126 dimensionedScalar("Su0", dimVelocity, 0.0)
130 volScalarField& Su0 = tSu0();
134 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi, 0.0);
137 forAll(Su0.boundaryField(), patchi)
139 forAll(Su0.boundaryField()[patchi], facei)
141 Su0.boundaryField()[patchi][facei] =
144 p.boundaryField()[patchi][facei],
145 Tu.boundaryField()[patchi][facei],
156 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
158 const volScalarField& p,
159 const volScalarField& Tu,
160 const volScalarField& phi
163 tmp<volScalarField> tSu0
176 dimensionedScalar("Su0", dimVelocity, 0.0)
180 volScalarField& Su0 = tSu0();
184 Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], 0.0);
187 forAll(Su0.boundaryField(), patchi)
189 forAll(Su0.boundaryField()[patchi], facei)
191 Su0.boundaryField()[patchi][facei] =
194 p.boundaryField()[patchi][facei],
195 Tu.boundaryField()[patchi][facei],
196 phi.boundaryField()[patchi][facei],
206 Foam::tmp<Foam::volScalarField>
207 Foam::laminarFlameSpeedModels::Gulders::operator()() const
209 if (hhuCombustionThermo_.composition().contains("ft"))
211 const volScalarField& ft = hhuCombustionThermo_.composition().Y("ft");
215 hhuCombustionThermo_.p(),
216 hhuCombustionThermo_.Tu(),
219 hhuCombustionThermo_.lookup("stoichiometricAirFuelMassRatio")
220 )*ft/((1 + SMALL) - ft)
227 hhuCombustionThermo_.p(),
228 hhuCombustionThermo_.Tu(),
235 // ************************************************************************* //