Transferred copyright to the OpenFOAM Foundation
[OpenFOAM-2.0.x.git] / src / thermophysicalModels / laminarFlameSpeed / Gulders / Gulders.C
blobe7d0fc9e65dd9f4b59fdde24963e71e4d7131d58
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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 \*---------------------------------------------------------------------------*/
26 #include "Gulders.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 namespace Foam
33 namespace laminarFlameSpeedModels
35     defineTypeNameAndDebug(Gulders, 0);
37     addToRunTimeSelectionTable
38     (
39         laminarFlameSpeed,
40         Gulders,
41         dictionary
42     );
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
77     scalar phi
78 ) const
80     if (phi > SMALL)
81     {
82         return W_*pow(phi, eta_)*exp(-xi_*sqr(phi - 1.075));
83     }
84     else
85     {
86         return 0.0;
87     }
91 inline Foam::scalar Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
93     scalar p,
94     scalar Tu,
95     scalar phi,
96     scalar Yres
97 ) const
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,
110     scalar phi
111 ) const
113     tmp<volScalarField> tSu0
114     (
115         new volScalarField
116         (
117             IOobject
118             (
119                 "Su0",
120                 p.time().timeName(),
121                 p.db(),
122                 IOobject::NO_READ,
123                 IOobject::NO_WRITE
124             ),
125             p.mesh(),
126             dimensionedScalar("Su0", dimVelocity, 0.0)
127         )
128     );
130     volScalarField& Su0 = tSu0();
132     forAll(Su0, celli)
133     {
134         Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi, 0.0);
135     }
137     forAll(Su0.boundaryField(), patchi)
138     {
139         forAll(Su0.boundaryField()[patchi], facei)
140         {
141             Su0.boundaryField()[patchi][facei] =
142                 Su0pTphi
143                 (
144                     p.boundaryField()[patchi][facei],
145                     Tu.boundaryField()[patchi][facei],
146                     phi,
147                     0.0
148                 );
149         }
150     }
152     return tSu0;
156 Foam::tmp<Foam::volScalarField> Foam::laminarFlameSpeedModels::Gulders::Su0pTphi
158     const volScalarField& p,
159     const volScalarField& Tu,
160     const volScalarField& phi
161 ) const
163     tmp<volScalarField> tSu0
164     (
165         new volScalarField
166         (
167             IOobject
168             (
169                 "Su0",
170                 p.time().timeName(),
171                 p.db(),
172                 IOobject::NO_READ,
173                 IOobject::NO_WRITE
174             ),
175             p.mesh(),
176             dimensionedScalar("Su0", dimVelocity, 0.0)
177         )
178     );
180     volScalarField& Su0 = tSu0();
182     forAll(Su0, celli)
183     {
184         Su0[celli] = Su0pTphi(p[celli], Tu[celli], phi[celli], 0.0);
185     }
187     forAll(Su0.boundaryField(), patchi)
188     {
189         forAll(Su0.boundaryField()[patchi], facei)
190         {
191             Su0.boundaryField()[patchi][facei] =
192                 Su0pTphi
193                 (
194                     p.boundaryField()[patchi][facei],
195                     Tu.boundaryField()[patchi][facei],
196                     phi.boundaryField()[patchi][facei],
197                     0.0
198                 );
199         }
200     }
202     return tSu0;
206 Foam::tmp<Foam::volScalarField>
207 Foam::laminarFlameSpeedModels::Gulders::operator()() const
209     if (hhuCombustionThermo_.composition().contains("ft"))
210     {
211         const volScalarField& ft = hhuCombustionThermo_.composition().Y("ft");
213         return Su0pTphi
214         (
215             hhuCombustionThermo_.p(),
216             hhuCombustionThermo_.Tu(),
217             dimensionedScalar
218             (
219                 hhuCombustionThermo_.lookup("stoichiometricAirFuelMassRatio")
220             )*ft/((1 + SMALL) - ft)
221         );
222     }
223     else
224     {
225         return Su0pTphi
226         (
227             hhuCombustionThermo_.p(),
228             hhuCombustionThermo_.Tu(),
229             equivalenceRatio_
230         );
231     }
235 // ************************************************************************* //