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-2020 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"
29 
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace combustionModels
36 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const word& modelType,
48  const fluidReactionThermo& thermo,
50  const word& combustionProperties
51 )
52 :
54  (
55  modelType,
56  thermo,
57  turb,
58  combustionProperties
59  ),
60  C_(this->coeffs().template lookup<scalar>("C")),
61  oxidantName_(this->coeffs().template lookupOrDefault<word>("oxidant", "O2"))
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
66 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
72 
74 {
75  this->wFuel_ ==
77 
78  this->fresCorrect();
79 
80  const label fuelI = this->fuelIndex();
81 
82  const volScalarField& YFuel = this->thermo().composition().Y()[fuelI];
83 
84  if (this->thermo().composition().contains(oxidantName_))
85  {
86  const volScalarField& YO2 =
87  this->thermo().composition().Y(oxidantName_);
88 
89  this->wFuel_ ==
90  C_*this->turbulence().muEff()
91  *mag(fvc::grad(YFuel) & fvc::grad(YO2))
92  *pos0(YFuel)*pos0(YO2);
93  }
94 }
95 
96 
98 {
100  {
101  this->coeffs().lookup("C") >> C_ ;
102  this->coeffs().readIfPresent("oxidant", oxidantName_);
103  return true;
104  }
105  else
106  {
107  return false;
108  }
109 }
110 
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
fluidReactionThermo & thermo
Definition: createFields.H:28
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
basicSpecieMixture & composition
virtual void correct()
Correct combustion rate.
Definition: diffusion.C:73
virtual bool read()
Update properties.
Definition: diffusion.C:97
Base-class for multi-component fluid thermodynamic properties.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
diffusion(const word &modelType, const fluidReactionThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
Definition: diffusion.C:46
const dimensionSet dimTime
Calculate the gradient of the given field.
A class for handling words, derived from string.
Definition: word.H:59
dimensionedScalar pos0(const dimensionedScalar &ds)
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::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Base class for combustion models.
const dimensionSet dimMass
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual bool read()
Update properties from given dictionary.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual ~diffusion()
Destructor.
Definition: diffusion.C:67
dimensioned< scalar > mag(const dimensioned< Type > &)
Base class for combustion models using basicSpecieMixture.
Simple diffusion-based combustion model based on the principle mixed is burnt. Additional parameter C...
Definition: diffusion.H:52
defineTypeNameAndDebug(diffusion, 0)
Abstract base class for turbulence models (RAS, LES and laminar).
addToRunTimeSelectionTable(combustionModel, diffusion, dictionary)
Namespace for OpenFOAM.