BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / combustionModels / infinitelyFastChemistry / infinitelyFastChemistry.C
blobc83cd2daa55a3e37f6ca74a35b82c216db1bf118
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
10     OpenFOAM is free software: you can redistribute it and/or modify it
11     under the terms of the GNU General Public License as published by
12     the Free Software Foundation, either version 3 of the License, or
13     (at your option) any later version.
15     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
16     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
17     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
18     for more details.
20     You should have received a copy of the GNU General Public License
21     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
23 \*---------------------------------------------------------------------------*/
25 #include "infinitelyFastChemistry.H"
26 #include "addToRunTimeSelectionTable.H"
27 #include "fvmSup.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 namespace Foam
33 namespace combustionModels
35     defineTypeNameAndDebug(infinitelyFastChemistry, 0);
36     addToRunTimeSelectionTable
37     (
38         combustionModel,
39         infinitelyFastChemistry,
40         dictionary
41     );
46 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
48 Foam::combustionModels::infinitelyFastChemistry::infinitelyFastChemistry
50     const dictionary& combustionProps,
51     hsCombustionThermo& thermo,
52     const compressible::turbulenceModel& turbulence,
53     const surfaceScalarField& phi,
54     const volScalarField& rho
57     combustionModel(typeName, combustionProps, thermo, turbulence, phi, rho),
58     C_(readScalar(coeffs_.lookup("C"))),
59     singleMixture_
60     (
61         dynamic_cast<singleStepReactingMixture<gasThermoPhysics>&>(thermo)
62     ),
63     wFuelNorm_
64     (
65         IOobject
66         (
67             "wFuelNorm",
68             mesh_.time().timeName(),
69             mesh_,
70             IOobject::NO_READ,
71             IOobject::NO_WRITE
72         ),
73         mesh_,
74         dimensionedScalar("zero", dimMass/pow3(dimLength)/dimTime, 0.0)
75     )
79 // * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
81 Foam::combustionModels::infinitelyFastChemistry::~infinitelyFastChemistry()
85 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
87 void Foam::combustionModels::infinitelyFastChemistry::correct()
89     singleMixture_.fresCorrect();
91     const label fuelI = singleMixture_.fuelIndex();
93     const volScalarField& YFuel = thermo_.composition().Y()[fuelI];
95     const dimensionedScalar s = singleMixture_.s();
97     if (thermo_.composition().contains("O2"))
98     {
99         const volScalarField& YO2 = thermo_.composition().Y("O2");
100         wFuelNorm_ == rho_/(mesh_.time().deltaT()*C_)*min(YFuel, YO2/s.value());
101     }
105 Foam::tmp<Foam::fvScalarMatrix>
106 Foam::combustionModels::infinitelyFastChemistry::R(volScalarField& Y) const
108     const label specieI = thermo_.composition().species()[Y.name()];
110     const label fNorm = singleMixture_.specieProd()[specieI];
112     const volScalarField fres(singleMixture_.fres(specieI));
114     const volScalarField wSpecie
115     (
116         wFuelNorm_*singleMixture_.specieStoichCoeffs()[specieI]
117       / max(fNorm*(Y - fres), scalar(0.001))
118     );
120     return -fNorm*wSpecie*fres + fNorm*fvm::Sp(wSpecie, Y);
124 Foam::tmp<Foam::volScalarField>
125 Foam::combustionModels::infinitelyFastChemistry::dQ() const
127     const label fuelI = singleMixture_.fuelIndex();
128     volScalarField& YFuel = thermo_.composition().Y(fuelI);
130     return -singleMixture_.qFuel()*(R(YFuel) & YFuel);
134 Foam::tmp<Foam::volScalarField>
135 Foam::combustionModels::infinitelyFastChemistry::wFuelNorm() const
137     return wFuelNorm_;
141 bool Foam::combustionModels::infinitelyFastChemistry::read
143     const dictionary& combustionProps
146     combustionModel::read(combustionProps);
147     coeffs_.lookup("C") >> C_ ;
149     return true;
153 // ************************************************************************* //