1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "transport.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 defineTypeNameAndDebug(transport, 0);
36 addToRunTimeSelectionTable(XiModel, transport, dictionary);
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 Foam::XiModels::transport::transport
45 const dictionary& XiProperties,
46 const hhuCombustionThermo& thermo,
47 const compressible::RASModel& turbulence,
48 const volScalarField& Su,
49 const volScalarField& rho,
50 const volScalarField& b,
51 const surfaceScalarField& phi
54 XiModel(XiProperties, thermo, turbulence, Su, rho, b, phi),
55 XiShapeCoef(readScalar(XiModelCoeffs_.lookup("XiShapeCoef"))),
56 XiEqModel_(XiEqModel::New(XiProperties, thermo, turbulence, Su)),
57 XiGModel_(XiGModel::New(XiProperties, thermo, turbulence, Su))
61 // * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
63 Foam::XiModels::transport::~transport()
67 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
69 Foam::tmp<Foam::volScalarField> Foam::XiModels::transport::Db() const
71 return XiGModel_->Db();
75 void Foam::XiModels::transport::correct
77 const fv::convectionScheme<scalar>& mvConvection
80 volScalarField XiEqEta(XiEqModel_->XiEq());
81 volScalarField GEta(XiGModel_->G());
83 volScalarField R(GEta*XiEqEta/(XiEqEta - 0.999));
85 volScalarField XiEqStar(R/(R - GEta));
89 1.0 + (1.0 + (2*XiShapeCoef)*(0.5 - b_))*(XiEqStar - 1.0)
92 volScalarField G(R*(XiEq - 1.0)/XiEq);
94 const objectRegistry& db = b_.db();
95 const volScalarField& betav = db.lookupObject<volScalarField>("betav");
96 const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
97 const surfaceScalarField& phiSt =
98 db.lookupObject<surfaceScalarField>("phiSt");
99 const volScalarField& Db = db.lookupObject<volScalarField>("Db");
100 const surfaceScalarField& nf = db.lookupObject<surfaceScalarField>("nf");
102 surfaceScalarField phiXi
107 - fvc::interpolate(fvc::laplacian(Db, b_)/mgb)*nf
108 + fvc::interpolate(rho_)*fvc::interpolate(Su_*(1.0/Xi_ - Xi_))*nf
114 betav*fvm::ddt(rho_, Xi_)
115 + mvConvection.fvmDiv(phi_, Xi_)
116 + fvm::div(phiXi, Xi_)
117 - fvm::Sp(fvc::div(phiXi), Xi_)
120 - fvm::Sp(betav*rho_*(R - G), Xi_)
123 // Correct boundedness of Xi
124 // ~~~~~~~~~~~~~~~~~~~~~~~~~
126 Xi_ = min(Xi_, 2.0*XiEq);
130 bool Foam::XiModels::transport::read(const dictionary& XiProperties)
132 XiModel::read(XiProperties);
134 XiModelCoeffs_.lookup("XiShapeCoef") >> XiShapeCoef;
140 // ************************************************************************* //