FreeFOAM The Cross-Platform CFD Toolkit
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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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 "transport.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiModels
34 {
35  defineTypeNameAndDebug(transport, 0);
36  addToRunTimeSelectionTable(XiModel, transport, dictionary);
37 };
38 };
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::XiModels::transport::transport
44 (
45  const dictionary& XiProperties,
46  const hhuCombustionThermo& thermo,
48  const volScalarField& Su,
49  const volScalarField& rho,
50  const volScalarField& b,
51  const surfaceScalarField& phi
52 )
53 :
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))
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68 
70 {
71  return XiGModel_->Db();
72 }
73 
74 
76 (
77  const fv::convectionScheme<scalar>& mvConvection
78 )
79 {
80  volScalarField XiEqEta = XiEqModel_->XiEq();
81  volScalarField GEta = XiGModel_->G();
82 
83  volScalarField R = GEta*XiEqEta/(XiEqEta - 0.999);
84 
85  volScalarField XiEqStar = R/(R - GEta);
86 
87  volScalarField XiEq =
88  1.0 + (1.0 + (2*XiShapeCoef)*(0.5 - b_))*(XiEqStar - 1.0);
89 
90  volScalarField G = R*(XiEq - 1.0)/XiEq;
91 
92  const objectRegistry& db = b_.db();
93  const volScalarField& betav = db.lookupObject<volScalarField>("betav");
94  const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
95  const surfaceScalarField& phiSt =
96  db.lookupObject<surfaceScalarField>("phiSt");
97  const volScalarField& Db = db.lookupObject<volScalarField>("Db");
98  const surfaceScalarField& nf = db.lookupObject<surfaceScalarField>("nf");
99 
100  surfaceScalarField phiXi
101  (
102  "phiXi",
103  phiSt
104  + (
105  - fvc::interpolate(fvc::laplacian(Db, b_)/mgb)*nf
106  + fvc::interpolate(rho_)*fvc::interpolate(Su_*(1.0/Xi_ - Xi_))*nf
107  )
108  );
109 
110  solve
111  (
112  betav*fvm::ddt(rho_, Xi_)
113  + mvConvection.fvmDiv(phi_, Xi_)
114  + fvm::div(phiXi, Xi_)
115  - fvm::Sp(fvc::div(phiXi), Xi_)
116  ==
117  betav*rho_*R
118  - fvm::Sp(betav*rho_*(R - G), Xi_)
119  );
120 
121  // Correct boundedness of Xi
122  // ~~~~~~~~~~~~~~~~~~~~~~~~~
123  Xi_.max(1.0);
124  Xi_ = min(Xi_, 2.0*XiEq);
125 }
126 
127 
128 bool Foam::XiModels::transport::read(const dictionary& XiProperties)
129 {
130  XiModel::read(XiProperties);
131 
132  XiModelCoeffs_.lookup("XiShapeCoef") >> XiShapeCoef;
133 
134  return true;
135 }
136 
137 
138 // ************************ vim: set sw=4 sts=4 et: ************************ //