Loading...
Searching...
No Matches
infinitelyFastChemistry.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 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30
31namespace Foam
33namespace combustionModels
34{
35
36// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37
38template<class ReactionThermo, class ThermoType>
39infinitelyFastChemistry<ReactionThermo, ThermoType>::infinitelyFastChemistry
40(
41 const word& modelType,
42 ReactionThermo& thermo,
44 const word& combustionProperties
45)
46:
47 singleStepCombustion<ReactionThermo, ThermoType>
48 (
49 modelType,
50 thermo,
51 turb,
52 combustionProperties
53 ),
54 C_(this->coeffs().getScalar("C"))
55{}
56
57
58// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59
60template<class ReactionThermo, class ThermoType>
62{}
63
64
65// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66
67template<class ReactionThermo, class ThermoType>
69{
70 this->wFuel_ == Zero;
71
72 if (this->active())
73 {
74 this->singleMixturePtr_->fresCorrect();
75
76 const label fuelI = this->singleMixturePtr_->fuelIndex();
77
78 const volScalarField& YFuel =
79 this->thermo().composition().Y()[fuelI];
80
81 const dimensionedScalar s = this->singleMixturePtr_->s();
82
83 if (this->thermo().composition().contains("O2"))
84 {
85 const volScalarField& YO2 = this->thermo().composition().Y("O2");
86
87 this->wFuel_ ==
88 this->rho()/(this->mesh().time().deltaT()*C_)
89 *min(YFuel, YO2/s.value());
90 }
91 }
92}
93
94
95template<class ReactionThermo, class ThermoType>
97{
99 {
100 this->coeffs().readEntry("C", C_);
101 return true;
102 }
103
104 return false;
105}
106
107
108// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109
110} // End namespace combustionModels
111} // End namespace Foam
112
113// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
compressible::turbulenceModel & turb
virtual ReactionThermo & thermo()
Return access to the thermo package.
const dictionary & coeffs() const
Return const dictionary of the model.
Switch active() const noexcept
Is combustion active?
virtual void correct()
Correct combustion rate.
singleStepReactingMixture< ThermoType > * singleMixturePtr_
Pointer to singleStepReactingMixture mixture.
volScalarField wFuel_
Fuel consumption rate.
virtual bool read()
Update properties from given dictionary.
Abstract base class for turbulence models (RAS, LES and laminar).
scalar getScalar(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get<scalar>(const word&, keyType::option).
A class for handling words, derived from Foam::string.
Definition word.H:66
basicSpecieMixture & composition
dynamicFvMesh & mesh
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
psiReactionThermo & thermo