diffusion.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) 2012-2015 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 CombThermoType, class ThermoType>
38 (
39  const word& modelType,
40  const fvMesh& mesh,
41  const word& phaseName
42 )
43 :
45  (
46  modelType,
47  mesh,
48  phaseName
49  ),
50  C_(readScalar(this->coeffs().lookup("C"))),
51  oxidantName_(this->coeffs().template lookupOrDefault<word>("oxidant", "O2"))
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  if (this->thermoPtr_->composition().contains(oxidantName_))
80  {
81  const volScalarField& YO2 =
82  this->thermoPtr_->composition().Y(oxidantName_);
83 
84  this->wFuel_ ==
85  C_*this->turbulence().muEff()
86  *mag(fvc::grad(YFuel) & fvc::grad(YO2))
87  *pos(YFuel)*pos(YO2);
88  }
89  }
90 }
91 
92 
93 template<class CombThermoType, class ThermoType>
95 {
97  {
98  this->coeffs().lookup("C") >> C_ ;
99  this->coeffs().readIfPresent("oxidant", oxidantName_);
100  return true;
101  }
102  else
103  {
104  return false;
105  }
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110 
111 } // End namespace combustionModels
112 } // End namespace Foam
113 
114 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
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
virtual void correct()
Correct combustion rate.
Definition: diffusion.C:65
virtual bool read()
Update properties.
Definition: diffusion.C:94
dimensionedScalar pos(const dimensionedScalar &ds)
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
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 succesful.
Definition: doubleScalar.H:63
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
virtual ~diffusion()
Destructor.
Definition: diffusion.C:58
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
Simple diffusion-based combustion model based on the principle mixed is burnt. Additional parameter C...
Definition: diffusion.H:53
Namespace for OpenFOAM.