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 | Copyright (C) 2011-2016 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 CombThermoType, class ThermoType>
36 infinitelyFastChemistry<CombThermoType, ThermoType>::infinitelyFastChemistry
37 (
38  const word& modelType,
39  const fvMesh& mesh,
40  const word& combustionProperties,
41  const word& phaseName
42 )
43 :
45  (
46  modelType,
47  mesh,
48  combustionProperties,
49  phaseName
50  ),
51  C_(readScalar(this->coeffs().lookup("C")))
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
56 
57 template<class CombThermoType, class ThermoType>
59 {}
60 
61 
62 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
63 
64 template<class CombThermoType, class ThermoType>
66 {
67  this->wFuel_ ==
69 
70  if (this->active())
71  {
72  this->singleMixturePtr_->fresCorrect();
73 
74  const label fuelI = this->singleMixturePtr_->fuelIndex();
75 
76  const volScalarField& YFuel =
77  this->thermoPtr_->composition().Y()[fuelI];
78 
79  const dimensionedScalar s = this->singleMixturePtr_->s();
80 
81  if (this->thermoPtr_->composition().contains("O2"))
82  {
83  const volScalarField& YO2 = this->thermoPtr_->composition().Y("O2");
84 
85  this->wFuel_ ==
86  this->rho()/(this->mesh().time().deltaT()*C_)
87  *min(YFuel, YO2/s.value());
88  }
89  }
90 }
91 
92 
93 template<class CombThermoType, class ThermoType>
95 {
97  {
98  this->coeffs().lookup("C") >> C_ ;
99  return true;
100  }
101  else
102  {
103  return false;
104  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 } // End namespace combustionModels
111 } // End namespace Foam
112 
113 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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
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.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
virtual void correct()
Correct combustion rate.
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.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
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.