infinitelyFastChemistry.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace combustionModels
35 {
36  defineTypeNameAndDebug(infinitelyFastChemistry, 0);
38  (
39  combustionModel,
40  infinitelyFastChemistry,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& modelType,
52  const fluidReactionThermo& thermo,
54  const word& combustionProperties
55 )
56 :
58  (
59  modelType,
60  thermo,
61  turb,
62  combustionProperties
63  ),
64  C_(this->coeffs().template lookup<scalar>("C"))
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 
77 {
78  this->wFuel_ ==
80 
81  this->fresCorrect();
82 
83  const label fuelI = this->fuelIndex();
84 
85  const volScalarField& YFuel = this->thermo().composition().Y()[fuelI];
86 
87  const dimensionedScalar s = this->s();
88 
89  if (this->thermo().composition().contains("O2"))
90  {
91  const volScalarField& YO2 = this->thermo().composition().Y("O2");
92 
93  this->wFuel_ ==
94  this->rho()/(this->mesh().time().deltaT()*C_)
95  *min(YFuel, YO2/s.value());
96  }
97 }
98 
99 
101 {
103  {
104  this->coeffs().lookup("C") >> C_ ;
105  return true;
106  }
107  else
108  {
109  return false;
110  }
111 }
112 
113 
114 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
basicSpecieMixture & composition
Base-class for multi-component fluid thermodynamic properties.
fvMesh & mesh
infinitelyFastChemistry(const word &modelType, const fluidReactionThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
const dimensionSet dimTime
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
A class for handling words, derived from string.
Definition: word.H:59
const Type & value() const
Return const reference to value.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
virtual void correct()
Correct combustion rate.
const dimensionSet dimMass
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual bool read()
Update properties from given dictionary.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Base class for combustion models using basicSpecieMixture.
defineTypeNameAndDebug(diffusion, 0)
Base class for single-phase compressible turbulence models.
addToRunTimeSelectionTable(combustionModel, diffusion, dictionary)
Namespace for OpenFOAM.