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 -------------------------------------------------------------------------------
11 License
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 
31 namespace Foam
32 {
33 namespace combustionModels
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class ReactionThermo, class ThermoType>
39 infinitelyFastChemistry<ReactionThermo, ThermoType>::infinitelyFastChemistry
40 (
41  const word& modelType,
42  ReactionThermo& thermo,
44  const word& combustionProperties
45 )
46 :
48  (
49  modelType,
50  thermo,
51  turb,
52  combustionProperties
53  ),
54  C_(this->coeffs().getScalar("C"))
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
60 template<class ReactionThermo, class ThermoType>
62 {}
63 
64 
65 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66 
67 template<class ReactionThermo, class ThermoType>
69 {
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 
95 template<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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
s
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))
Definition: gmvOutputSpray.H:25
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::combustionModels::infinitelyFastChemistry::~infinitelyFastChemistry
virtual ~infinitelyFastChemistry()
Destructor.
Definition: infinitelyFastChemistry.C:61
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
rho
rho
Definition: readInitialConditions.H:88
composition
basicSpecieMixture & composition
Definition: createFieldRefs.H:6
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::combustionModels::singleStepCombustion
Base class for combustion models using singleStepReactingMixture.
Definition: singleStepCombustion.H:58
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::combustionModels::infinitelyFastChemistry::read
virtual bool read()
Update properties.
Definition: infinitelyFastChemistry.C:96
Foam::combustionModels::infinitelyFastChemistry::correct
virtual void correct()
Correct combustion rate.
Definition: infinitelyFastChemistry.C:68
infinitelyFastChemistry.H
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
Foam::GeometricField< scalar, fvPatchField, volMesh >
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::compressibleTurbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: compressibleTurbulenceModel.H:54