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 "RASModel.H"
00027 #include <finiteVolume/wallFvPatch.H>
00028 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00029
00030
00031
00032 namespace Foam
00033 {
00034 namespace compressible
00035 {
00036
00037
00038
00039 defineTypeNameAndDebug(RASModel, 0);
00040 defineRunTimeSelectionTable(RASModel, dictionary);
00041 addToRunTimeSelectionTable(turbulenceModel, RASModel, turbulenceModel);
00042
00043
00044
00045 void RASModel::printCoeffs()
00046 {
00047 if (printCoeffs_)
00048 {
00049 Info<< type() << "Coeffs" << coeffDict_ << endl;
00050 }
00051 }
00052
00053
00054
00055
00056 RASModel::RASModel
00057 (
00058 const word& type,
00059 const volScalarField& rho,
00060 const volVectorField& U,
00061 const surfaceScalarField& phi,
00062 const basicThermo& thermophysicalModel
00063 )
00064 :
00065 turbulenceModel(rho, U, phi, thermophysicalModel),
00066
00067 IOdictionary
00068 (
00069 IOobject
00070 (
00071 "RASProperties",
00072 U.time().constant(),
00073 U.db(),
00074 IOobject::MUST_READ,
00075 IOobject::NO_WRITE
00076 )
00077 ),
00078
00079 turbulence_(lookup("turbulence")),
00080 printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
00081 coeffDict_(subOrEmptyDict(type + "Coeffs")),
00082
00083 k0_("k0", dimVelocity*dimVelocity, SMALL),
00084 epsilon0_("epsilon0", k0_.dimensions()/dimTime, SMALL),
00085 epsilonSmall_("epsilonSmall", epsilon0_.dimensions(), SMALL),
00086 omega0_("omega0", dimless/dimTime, SMALL),
00087 omegaSmall_("omegaSmall", omega0_.dimensions(), SMALL),
00088
00089 y_(mesh_)
00090 {
00091 k0_.readIfPresent(*this);
00092 epsilon0_.readIfPresent(*this);
00093 epsilonSmall_.readIfPresent(*this);
00094 omega0_.readIfPresent(*this);
00095 omegaSmall_.readIfPresent(*this);
00096
00097
00098
00099 mesh_.deltaCoeffs();
00100 }
00101
00102
00103
00104
00105 autoPtr<RASModel> RASModel::New
00106 (
00107 const volScalarField& rho,
00108 const volVectorField& U,
00109 const surfaceScalarField& phi,
00110 const basicThermo& thermophysicalModel
00111 )
00112 {
00113 word modelName;
00114
00115
00116
00117
00118 {
00119 IOdictionary dict
00120 (
00121 IOobject
00122 (
00123 "RASProperties",
00124 U.time().constant(),
00125 U.db(),
00126 IOobject::MUST_READ,
00127 IOobject::NO_WRITE
00128 )
00129 );
00130
00131 dict.lookup("RASModel") >> modelName;
00132 }
00133
00134 Info<< "Selecting RAS turbulence model " << modelName << endl;
00135
00136 dictionaryConstructorTable::iterator cstrIter =
00137 dictionaryConstructorTablePtr_->find(modelName);
00138
00139 if (cstrIter == dictionaryConstructorTablePtr_->end())
00140 {
00141 FatalErrorIn
00142 (
00143 "RASModel::New(const volScalarField&, "
00144 "const volVectorField&, const surfaceScalarField&, "
00145 "basicThermo&)"
00146 ) << "Unknown RASModel type " << modelName
00147 << endl << endl
00148 << "Valid RASModel types are :" << endl
00149 << dictionaryConstructorTablePtr_->sortedToc()
00150 << exit(FatalError);
00151 }
00152
00153 return autoPtr<RASModel>
00154 (
00155 cstrIter()(rho, U, phi, thermophysicalModel)
00156 );
00157 }
00158
00159
00160
00161
00162 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
00163 {
00164 scalar ypl = 11.0;
00165
00166 for (int i=0; i<10; i++)
00167 {
00168 ypl = log(max(E*ypl, 1))/kappa;
00169 }
00170
00171 return ypl;
00172 }
00173
00174
00175 tmp<scalarField> RASModel::yPlus(const label patchNo, const scalar Cmu) const
00176 {
00177 const fvPatch& curPatch = mesh_.boundary()[patchNo];
00178
00179 tmp<scalarField> tYp(new scalarField(curPatch.size()));
00180 scalarField& Yp = tYp();
00181
00182 if (isA<wallFvPatch>(curPatch))
00183 {
00184 Yp = pow(Cmu, 0.25)
00185 *y_[patchNo]
00186 *sqrt(k()().boundaryField()[patchNo].patchInternalField())
00187 /(
00188 mu().boundaryField()[patchNo].patchInternalField()
00189 /rho_.boundaryField()[patchNo]
00190 );
00191 }
00192 else
00193 {
00194 WarningIn
00195 (
00196 "tmp<scalarField> RASModel::yPlus(const label patchNo) const"
00197 ) << "Patch " << patchNo << " is not a wall. Returning null field"
00198 << nl << endl;
00199
00200 Yp.setSize(0);
00201 }
00202
00203 return tYp;
00204 }
00205
00206
00207 void RASModel::correct()
00208 {
00209 if (mesh_.changing())
00210 {
00211 y_.correct();
00212 }
00213 }
00214
00215
00216 bool RASModel::read()
00217 {
00218 if (regIOobject::read())
00219 {
00220 lookup("turbulence") >> turbulence_;
00221
00222 if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
00223 {
00224 coeffDict_ <<= *dictPtr;
00225 }
00226
00227 k0_.readIfPresent(*this);
00228 epsilon0_.readIfPresent(*this);
00229 epsilonSmall_.readIfPresent(*this);
00230 omega0_.readIfPresent(*this);
00231 omegaSmall_.readIfPresent(*this);
00232
00233 return true;
00234 }
00235 else
00236 {
00237 return false;
00238 }
00239 }
00240
00241
00242
00243
00244 }
00245 }
00246
00247