00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "basicXiSubXiEq.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028
00029
00030
00031 namespace Foam
00032 {
00033 namespace XiEqModels
00034 {
00035 defineTypeNameAndDebug(basicSubGrid, 0);
00036 addToRunTimeSelectionTable(XiEqModel, basicSubGrid, dictionary);
00037 };
00038 };
00039
00040
00041
00042
00043 Foam::XiEqModels::basicSubGrid::basicSubGrid
00044 (
00045 const dictionary& XiEqProperties,
00046 const hhuCombustionThermo& thermo,
00047 const compressible::RASModel& turbulence,
00048 const volScalarField& Su
00049 )
00050 :
00051 XiEqModel(XiEqProperties, thermo, turbulence, Su),
00052
00053 N_
00054 (
00055 IOobject
00056 (
00057 "N",
00058 Su.mesh().time().findInstance(polyMesh::meshSubDir, "N"),
00059 polyMesh::meshSubDir,
00060 Su.mesh(),
00061 IOobject::MUST_READ,
00062 IOobject::NO_WRITE
00063 ),
00064 Su.mesh()
00065 ),
00066
00067 ns_
00068 (
00069 IOobject
00070 (
00071 "ns",
00072 Su.mesh().time().findInstance(polyMesh::meshSubDir, "ns"),
00073 polyMesh::meshSubDir,
00074 Su.mesh(),
00075 IOobject::MUST_READ,
00076 IOobject::NO_WRITE
00077 ),
00078 Su.mesh()
00079 ),
00080
00081 B_
00082 (
00083 IOobject
00084 (
00085 "B",
00086 Su.mesh().time().findInstance(polyMesh::meshSubDir, "B"),
00087 polyMesh::meshSubDir,
00088 Su.mesh(),
00089 IOobject::MUST_READ,
00090 IOobject::NO_WRITE
00091 ),
00092 Su.mesh()
00093 ),
00094
00095 Lobs_
00096 (
00097 IOobject
00098 (
00099 "Lobs",
00100 Su.mesh().time().findInstance(polyMesh::meshSubDir, "Lobs"),
00101 polyMesh::meshSubDir,
00102 Su.mesh(),
00103 IOobject::MUST_READ,
00104 IOobject::NO_WRITE
00105 ),
00106 Su.mesh()
00107 ),
00108
00109 XiEqModel_(XiEqModel::New(XiEqModelCoeffs_, thermo, turbulence, Su))
00110 {}
00111
00112
00113
00114
00115 Foam::XiEqModels::basicSubGrid::~basicSubGrid()
00116 {}
00117
00118
00119
00120
00121 Foam::tmp<Foam::volScalarField> Foam::XiEqModels::basicSubGrid::XiEq() const
00122 {
00123 const objectRegistry& db = Su_.db();
00124 const volVectorField& U = db.lookupObject<volVectorField>("U");
00125
00126 volScalarField magU = mag(U);
00127 volVectorField Uhat =
00128 U/(mag(U) + dimensionedScalar("Usmall", U.dimensions(), 1e-4));
00129
00130 volScalarField n = max(N_ - (Uhat & ns_ & Uhat), scalar(1e-4));
00131
00132 volScalarField b = (Uhat & B_ & Uhat)/n;
00133
00134 volScalarField up = sqrt((2.0/3.0)*turbulence_.k());
00135
00136 volScalarField XiSubEq =
00137 scalar(1)
00138 + max(2.2*sqrt(b), min(0.34*magU/up, scalar(1.6)))
00139 *min(0.25*n, scalar(1));
00140
00141 return XiSubEq*XiEqModel_->XiEq();
00142 }
00143
00144
00145 bool Foam::XiEqModels::basicSubGrid::read(const dictionary& XiEqProperties)
00146 {
00147 XiEqModel::read(XiEqProperties);
00148
00149 return XiEqModel_->read(XiEqModelCoeffs_);
00150 }
00151
00152
00153