1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "transport.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 defineTypeNameAndDebug(transport, 0);
37 addToRunTimeSelectionTable(XiModel, transport, dictionary);
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 Foam::XiModels::transport::transport
46 const dictionary& XiProperties,
47 const hhuCombustionThermo& thermo,
48 const compressible::RASModel& turbulence,
49 const volScalarField& Su,
50 const volScalarField& rho,
51 const volScalarField& b,
52 const surfaceScalarField& phi
55 XiModel(XiProperties, thermo, turbulence, Su, rho, b, phi),
56 XiShapeCoef(readScalar(XiModelCoeffs_.lookup("XiShapeCoef"))),
57 XiEqModel_(XiEqModel::New(XiProperties, thermo, turbulence, Su)),
58 XiGModel_(XiGModel::New(XiProperties, thermo, turbulence, Su))
62 // * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
64 Foam::XiModels::transport::~transport()
68 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
70 Foam::tmp<Foam::volScalarField> Foam::XiModels::transport::Db() const
72 return XiGModel_->Db();
76 void Foam::XiModels::transport::correct
78 const fv::convectionScheme<scalar>& mvConvection
81 volScalarField XiEqEta = XiEqModel_->XiEq();
82 volScalarField GEta = XiGModel_->G();
84 volScalarField R = GEta*XiEqEta/(XiEqEta - 0.999);
86 volScalarField XiEqStar = R/(R - GEta);
89 1.0 + (1.0 + (2*XiShapeCoef)*(0.5 - b_))*(XiEqStar - 1.0);
91 volScalarField G = R*(XiEq - 1.0)/XiEq;
93 const objectRegistry& db = b_.db();
94 const volScalarField& betav = db.lookupObject<volScalarField>("betav");
95 const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
96 const surfaceScalarField& phiSt =
97 db.lookupObject<surfaceScalarField>("phiSt");
98 const volScalarField& Db = db.lookupObject<volScalarField>("Db");
99 const surfaceScalarField& nf = db.lookupObject<surfaceScalarField>("nf");
101 surfaceScalarField phiXi
106 - fvc::interpolate(fvc::laplacian(Db, b_)/mgb)*nf
107 + fvc::interpolate(rho_)*fvc::interpolate(Su_*(1.0/Xi_ - Xi_))*nf
113 betav*fvm::ddt(rho_, Xi_)
114 + mvConvection.fvmDiv(phi_, Xi_)
115 + fvm::div(phiXi, Xi_)
116 - fvm::Sp(fvc::div(phiXi), Xi_)
119 - fvm::Sp(betav*rho_*(R - G), Xi_)
122 // Correct boundedness of Xi
123 // ~~~~~~~~~~~~~~~~~~~~~~~~~
125 Xi_ = min(Xi_, 2.0*XiEq);
129 bool Foam::XiModels::transport::read(const dictionary& XiProperties)
131 XiModel::read(XiProperties);
133 XiModelCoeffs_.lookup("XiShapeCoef") >> XiShapeCoef;
139 // ************************************************************************* //