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-2018 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 
27 
28 namespace Foam
29 {
30 namespace combustionModels
31 {
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class ReactionThermo, class ThermoType>
37 (
38  const word& modelType,
39  ReactionThermo& thermo,
40  const compressibleTurbulenceModel& turb,
41  const word& combustionProperties
42 )
43 :
45  (
46  modelType,
47  thermo,
48  turb,
49  combustionProperties
50  ),
51  C_(readScalar(this->coeffs().lookup("C")))
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
56 
57 template<class ReactionThermo, class ThermoType>
59 {}
60 
61 
62 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
63 
64 template<class ReactionThermo, class ThermoType>
66 {
67  this->wFuel_ ==
69 
70  this->singleMixturePtr_->fresCorrect();
71 
72  const label fuelI = this->singleMixturePtr_->fuelIndex();
73 
74  const volScalarField& YFuel = this->thermo().composition().Y()[fuelI];
75 
76  const dimensionedScalar s = this->singleMixturePtr_->s();
77 
78  if (this->thermo().composition().contains("O2"))
79  {
80  const volScalarField& YO2 = this->thermo().composition().Y("O2");
81 
82  this->wFuel_ ==
83  this->rho()/(this->mesh().time().deltaT()*C_)
84  *min(YFuel, YO2/s.value());
85  }
86 }
87 
88 
89 template<class ReactionThermo, class ThermoType>
91 {
93  {
94  this->coeffs().lookup("C") >> C_ ;
95  return true;
96  }
97  else
98  {
99  return false;
100  }
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
105 
106 } // End namespace combustionModels
107 } // End namespace Foam
108 
109 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
basicSpecieMixture & composition
rhoReactionThermo & thermo
Definition: createFields.H:28
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))
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
const Type & value() const
Return const reference to value.
infinitelyFastChemistry(const word &modelType, ReactionThermo &thermo, const compressibleTurbulenceModel &turb, const word &combustionProperties)
Construct from components.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
virtual void correct()
Correct combustion rate.
Abstract base class for turbulence models (RAS, LES and laminar).
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Base class for combustion models using singleStepReactingMixture.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Namespace for OpenFOAM.