diffusion.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) 2012-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 
26 #include "diffusion.H"
27 #include "fvcGrad.H"
28 
29 namespace Foam
30 {
31 namespace combustionModels
32 {
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class ReactionThermo, class ThermoType>
38 (
39  const word& modelType,
40  ReactionThermo& thermo,
41  const compressibleTurbulenceModel& turb,
42  const word& combustionProperties
43 )
44 :
46  (
47  modelType,
48  thermo,
49  turb,
50  combustionProperties
51  ),
52  C_(readScalar(this->coeffs().lookup("C"))),
53  oxidantName_(this->coeffs().template lookupOrDefault<word>("oxidant", "O2"))
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
59 template<class ReactionThermo, class ThermoType>
61 {}
62 
63 
64 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65 
66 template<class ReactionThermo, class ThermoType>
68 {
69  this->wFuel_ ==
71 
72  this->singleMixturePtr_->fresCorrect();
73 
74  const label fuelI = this->singleMixturePtr_->fuelIndex();
75 
76  const volScalarField& YFuel = this->thermo().composition().Y()[fuelI];
77 
78  if (this->thermo().composition().contains(oxidantName_))
79  {
80  const volScalarField& YO2 =
81  this->thermo().composition().Y(oxidantName_);
82 
83  this->wFuel_ ==
84  C_*this->turbulence().muEff()
85  *mag(fvc::grad(YFuel) & fvc::grad(YO2))
86  *pos0(YFuel)*pos0(YO2);
87  }
88 }
89 
90 
91 template<class ReactionThermo, class ThermoType>
93 {
95  {
96  this->coeffs().lookup("C") >> C_ ;
97  this->coeffs().readIfPresent("oxidant", oxidantName_);
98  return true;
99  }
100  else
101  {
102  return false;
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
108 
109 } // End namespace combustionModels
110 } // End namespace Foam
111 
112 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
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
diffusion(const word &modelType, ReactionThermo &thermo, const compressibleTurbulenceModel &turb, const word &combustionProperties)
Construct from components.
Definition: diffusion.C:38
virtual void correct()
Correct combustion rate.
Definition: diffusion.C:67
virtual bool read()
Update properties.
Definition: diffusion.C:92
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
rhoReactionThermo & thermo
Definition: createFields.H:28
stressControl lookup("compactNormalStress") >> compactNormalStress
Calculate the gradient of the given field.
A class for handling words, derived from string.
Definition: word.H:59
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
dimensionedScalar pos0(const dimensionedScalar &ds)
Abstract base class for turbulence models (RAS, LES and laminar).
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
virtual ~diffusion()
Destructor.
Definition: diffusion.C:60
dimensioned< scalar > mag(const dimensioned< Type > &)
Base class for combustion models using singleStepReactingMixture.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Namespace for OpenFOAM.