transport.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-2015 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "transport.H"
29 #include "surfaceInterpolate.H"
30 #include "fvmDdt.H"
31 #include "fvcLaplacian.H"
32 #include "fvmDiv.H"
33 #include "fvmSup.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace XiModels
41 {
42  defineTypeNameAndDebug(transport, 0);
43  addToRunTimeSelectionTable(XiModel, transport, dictionary);
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 Foam::XiModels::transport::transport
51 (
52  const dictionary& XiProperties,
53  const psiuReactionThermo& thermo,
55  const volScalarField& Su,
56  const volScalarField& rho,
57  const volScalarField& b,
58  const surfaceScalarField& phi
59 )
60 :
61  XiModel(XiProperties, thermo, turbulence, Su, rho, b, phi),
62  XiShapeCoef(XiModelCoeffs_.get<scalar>("XiShapeCoef")),
63  XiEqModel_(XiEqModel::New(XiProperties, thermo, turbulence, Su)),
64  XiGModel_(XiGModel::New(XiProperties, thermo, turbulence, Su))
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 
77 {
78  return XiGModel_->Db();
79 }
80 
81 
83 (
84  const fv::convectionScheme<scalar>& mvConvection
85 )
86 {
87  volScalarField XiEqEta(XiEqModel_->XiEq());
88  volScalarField GEta(XiGModel_->G());
89 
90  volScalarField R(GEta*XiEqEta/(XiEqEta - 0.999));
91 
92  volScalarField XiEqStar(R/(R - GEta));
93 
94  volScalarField XiEq
95  (
96  1.0 + (1.0 + (2*XiShapeCoef)*(0.5 - b_))*(XiEqStar - 1.0)
97  );
98 
99  volScalarField G(R*(XiEq - 1.0)/XiEq);
100 
101  const objectRegistry& db = b_.db();
102  const volScalarField& betav = db.lookupObject<volScalarField>("betav");
103  const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
104  const surfaceScalarField& phiSt =
105  db.lookupObject<surfaceScalarField>("phiSt");
106  const volScalarField& Db = db.lookupObject<volScalarField>("Db");
107  const surfaceScalarField& nf = db.lookupObject<surfaceScalarField>("nf");
108 
109  surfaceScalarField phiXi
110  (
111  "phiXi",
112  phiSt
113  + (
114  - fvc::interpolate(fvc::laplacian(Db, b_)/mgb)*nf
115  + fvc::interpolate(rho_)*fvc::interpolate(Su_*(1.0/Xi_ - Xi_))*nf
116  )
117  );
118 
119  solve
120  (
121  betav*fvm::ddt(rho_, Xi_)
122  + mvConvection.fvmDiv(phi_, Xi_)
123  + fvm::div(phiXi, Xi_)
124  - fvm::Sp(fvc::div(phiXi), Xi_)
125  ==
126  betav*rho_*R
127  - fvm::Sp(betav*rho_*(R - G), Xi_)
128  );
129 
130  // Correct boundedness of Xi
131  // ~~~~~~~~~~~~~~~~~~~~~~~~~
132  Xi_.max(1.0);
133  Xi_ = min(Xi_, 2.0*XiEq);
134 }
135 
136 
137 bool Foam::XiModels::transport::read(const dictionary& XiProperties)
138 {
139  XiModel::read(XiProperties);
140 
141  XiModelCoeffs_.readEntry("XiShapeCoef", XiShapeCoef);
142 
143  return true;
144 }
145 
146 
147 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
zeroField Su
Definition: alphaSuSp.H:1
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:40
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
virtual bool read(const dictionary &XiProperties)=0
Update properties from given dictionary.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Calculate the laplacian of the given field.
const volScalarField & betav
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
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("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", 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
defineTypeNameAndDebug(combustionModel, 0)
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
volScalarField Db("Db", turbulence->muEff())
Calculate the matrix for the divergence of the given field and flux.
virtual void correct()
Correct the flame-wrinkling Xi.
Definition: transport.H:130
#define R(A, B, C, D, E, F, K, M)
virtual bool read(const dictionary &XiProperties)
Update properties from given dictionary.
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.divScheme("div(phi,Yi_h)")))
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual ~transport()
Destructor.
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.