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 "basicThermo.H"
00027 #include <finiteVolume/fvMesh.H>
00028 #include <OpenFOAM/HashTable.H>
00029 #include <finiteVolume/zeroGradientFvPatchFields.H>
00030 #include <basicThermophysicalModels/fixedEnthalpyFvPatchScalarField.H>
00031 #include <basicThermophysicalModels/gradientEnthalpyFvPatchScalarField.H>
00032 #include <basicThermophysicalModels/mixedEnthalpyFvPatchScalarField.H>
00033 #include <basicThermophysicalModels/fixedInternalEnergyFvPatchScalarField.H>
00034 #include <basicThermophysicalModels/gradientInternalEnergyFvPatchScalarField.H>
00035 #include <basicThermophysicalModels/mixedInternalEnergyFvPatchScalarField.H>
00036
00037
00038
00039 namespace Foam
00040 {
00041 defineTypeNameAndDebug(basicThermo, 0);
00042 }
00043
00044
00045
00046 Foam::wordList Foam::basicThermo::hBoundaryTypes()
00047 {
00048 const volScalarField::GeometricBoundaryField& tbf = T_.boundaryField();
00049
00050 wordList hbt = tbf.types();
00051
00052 forAll(tbf, patchi)
00053 {
00054 if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
00055 {
00056 hbt[patchi] = fixedEnthalpyFvPatchScalarField::typeName;
00057 }
00058 else if
00059 (
00060 isA<zeroGradientFvPatchScalarField>(tbf[patchi])
00061 || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
00062 )
00063 {
00064 hbt[patchi] = gradientEnthalpyFvPatchScalarField::typeName;
00065 }
00066 else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
00067 {
00068 hbt[patchi] = mixedEnthalpyFvPatchScalarField::typeName;
00069 }
00070 }
00071
00072 return hbt;
00073 }
00074
00075
00076 void Foam::basicThermo::hBoundaryCorrection(volScalarField& h)
00077 {
00078 volScalarField::GeometricBoundaryField& hbf = h.boundaryField();
00079
00080 forAll(hbf, patchi)
00081 {
00082 if (isA<gradientEnthalpyFvPatchScalarField>(hbf[patchi]))
00083 {
00084 refCast<gradientEnthalpyFvPatchScalarField>(hbf[patchi]).gradient()
00085 = hbf[patchi].fvPatchField::snGrad();
00086 }
00087 else if (isA<mixedEnthalpyFvPatchScalarField>(hbf[patchi]))
00088 {
00089 refCast<mixedEnthalpyFvPatchScalarField>(hbf[patchi]).refGrad()
00090 = hbf[patchi].fvPatchField::snGrad();
00091 }
00092 }
00093 }
00094
00095
00096 Foam::wordList Foam::basicThermo::eBoundaryTypes()
00097 {
00098 const volScalarField::GeometricBoundaryField& tbf = T_.boundaryField();
00099
00100 wordList ebt = tbf.types();
00101
00102 forAll(tbf, patchi)
00103 {
00104 if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
00105 {
00106 ebt[patchi] = fixedInternalEnergyFvPatchScalarField::typeName;
00107 }
00108 else if
00109 (
00110 isA<zeroGradientFvPatchScalarField>(tbf[patchi])
00111 || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
00112 )
00113 {
00114 ebt[patchi] = gradientInternalEnergyFvPatchScalarField::typeName;
00115 }
00116 else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
00117 {
00118 ebt[patchi] = mixedInternalEnergyFvPatchScalarField::typeName;
00119 }
00120 }
00121
00122 return ebt;
00123 }
00124
00125
00126 void Foam::basicThermo::eBoundaryCorrection(volScalarField& e)
00127 {
00128 volScalarField::GeometricBoundaryField& ebf = e.boundaryField();
00129
00130 forAll(ebf, patchi)
00131 {
00132 if (isA<gradientInternalEnergyFvPatchScalarField>(ebf[patchi]))
00133 {
00134 refCast<gradientInternalEnergyFvPatchScalarField>(ebf[patchi])
00135 .gradient() = ebf[patchi].fvPatchField::snGrad();
00136 }
00137 else if (isA<mixedInternalEnergyFvPatchScalarField>(ebf[patchi]))
00138 {
00139 refCast<mixedInternalEnergyFvPatchScalarField>(ebf[patchi])
00140 .refGrad() = ebf[patchi].fvPatchField::snGrad();
00141 }
00142 }
00143 }
00144
00145
00146
00147 Foam::basicThermo::basicThermo(const fvMesh& mesh)
00148 :
00149 IOdictionary
00150 (
00151 IOobject
00152 (
00153 "thermophysicalProperties",
00154 mesh.time().constant(),
00155 mesh,
00156 IOobject::MUST_READ,
00157 IOobject::NO_WRITE
00158 )
00159 ),
00160
00161 p_
00162 (
00163 IOobject
00164 (
00165 "p",
00166 mesh.time().timeName(),
00167 mesh,
00168 IOobject::MUST_READ,
00169 IOobject::AUTO_WRITE
00170 ),
00171 mesh
00172 ),
00173
00174 psi_
00175 (
00176 IOobject
00177 (
00178 "psi",
00179 mesh.time().timeName(),
00180 mesh,
00181 IOobject::NO_READ,
00182 IOobject::NO_WRITE
00183 ),
00184 mesh,
00185 dimensionSet(0, -2, 2, 0, 0)
00186 ),
00187
00188 T_
00189 (
00190 IOobject
00191 (
00192 "T",
00193 mesh.time().timeName(),
00194 mesh,
00195 IOobject::MUST_READ,
00196 IOobject::AUTO_WRITE
00197 ),
00198 mesh
00199 ),
00200
00201 mu_
00202 (
00203 IOobject
00204 (
00205 "mu",
00206 mesh.time().timeName(),
00207 mesh,
00208 IOobject::NO_READ,
00209 IOobject::NO_WRITE
00210 ),
00211 mesh,
00212 dimensionSet(1, -1, -1, 0, 0)
00213 ),
00214
00215 alpha_
00216 (
00217 IOobject
00218 (
00219 "alpha",
00220 mesh.time().timeName(),
00221 mesh,
00222 IOobject::NO_READ,
00223 IOobject::NO_WRITE
00224 ),
00225 mesh,
00226 dimensionSet(1, -1, -1, 0, 0)
00227 )
00228 {}
00229
00230
00231
00232
00233 Foam::basicThermo::~basicThermo()
00234 {}
00235
00236
00237
00238
00239 Foam::volScalarField& Foam::basicThermo::p()
00240 {
00241 return p_;
00242 }
00243
00244
00245 const Foam::volScalarField& Foam::basicThermo::p() const
00246 {
00247 return p_;
00248 }
00249
00250
00251 const Foam::volScalarField& Foam::basicThermo::psi() const
00252 {
00253 return psi_;
00254 }
00255
00256
00257 Foam::volScalarField& Foam::basicThermo::h()
00258 {
00259 notImplemented("basicThermo::h()");
00260 return const_cast<volScalarField&>(volScalarField::null());
00261 }
00262
00263
00264 const Foam::volScalarField& Foam::basicThermo::h() const
00265 {
00266 notImplemented("basicThermo::h() const");
00267 return volScalarField::null();
00268 }
00269
00270
00271 Foam::tmp<Foam::scalarField> Foam::basicThermo::h
00272 (
00273 const scalarField& T,
00274 const labelList& cells
00275 ) const
00276 {
00277 notImplemented
00278 (
00279 "basicThermo::h"
00280 "(const scalarField& T, const labelList& cells) const"
00281 );
00282 return tmp<scalarField>(NULL);
00283 }
00284
00285
00286 Foam::tmp<Foam::scalarField> Foam::basicThermo::h
00287 (
00288 const scalarField& T,
00289 const label patchi
00290 ) const
00291 {
00292 notImplemented
00293 (
00294 "basicThermo::h"
00295 "(const scalarField& T, const label patchi) const"
00296 );
00297 return tmp<scalarField>(NULL);
00298 }
00299
00300
00301 Foam::volScalarField& Foam::basicThermo::hs()
00302 {
00303 notImplemented("basicThermo::hs()");
00304 return const_cast<volScalarField&>(volScalarField::null());
00305 }
00306
00307
00308 const Foam::volScalarField& Foam::basicThermo::hs() const
00309 {
00310 notImplemented("basicThermo::hs() const");
00311 return volScalarField::null();
00312 }
00313
00314
00315 Foam::tmp<Foam::scalarField> Foam::basicThermo::hs
00316 (
00317 const scalarField& T,
00318 const labelList& cells
00319 ) const
00320 {
00321 notImplemented
00322 (
00323 "basicThermo::hs"
00324 "(const scalarField& T, const labelList& cells) const"
00325 );
00326 return tmp<scalarField>(NULL);
00327 }
00328
00329
00330 Foam::tmp<Foam::scalarField> Foam::basicThermo::hs
00331 (
00332 const scalarField& T,
00333 const label patchi
00334 ) const
00335 {
00336 notImplemented
00337 (
00338 "basicThermo::hs"
00339 "(const scalarField& T, const label patchi) const"
00340 );
00341 return tmp<scalarField>(NULL);
00342 }
00343
00344
00345 Foam::tmp<Foam::volScalarField> Foam::basicThermo::hc() const
00346 {
00347 notImplemented("basicThermo::hc()");
00348 return volScalarField::null();
00349 }
00350
00351
00352 Foam::volScalarField& Foam::basicThermo::e()
00353 {
00354 notImplemented("basicThermo::e()");
00355 return const_cast<volScalarField&>(volScalarField::null());
00356 }
00357
00358
00359 const Foam::volScalarField& Foam::basicThermo::e() const
00360 {
00361 notImplemented("basicThermo::e()");
00362 return volScalarField::null();
00363 }
00364
00365
00366 Foam::tmp<Foam::scalarField> Foam::basicThermo::e
00367 (
00368 const scalarField& T,
00369 const labelList& cells
00370 ) const
00371 {
00372 notImplemented
00373 (
00374 "basicThermo::e"
00375 "(const scalarField& T, const labelList& cells) const"
00376 );
00377 return tmp<scalarField>(NULL);
00378 }
00379
00380
00381 Foam::tmp<Foam::scalarField> Foam::basicThermo::e
00382 (
00383 const scalarField& T,
00384 const label patchi
00385 ) const
00386 {
00387 notImplemented
00388 (
00389 "basicThermo::e"
00390 "(const scalarField& T, const label patchi) const"
00391 );
00392 return tmp<scalarField>(NULL);
00393 }
00394
00395
00396 const Foam::volScalarField& Foam::basicThermo::T() const
00397 {
00398 return T_;
00399 }
00400
00401
00402 Foam::tmp<Foam::scalarField> Foam::basicThermo::Cp
00403 (
00404 const scalarField& T,
00405 const label patchi
00406 ) const
00407 {
00408 notImplemented
00409 (
00410 "basicThermo::Cp"
00411 "(const scalarField& T, const label patchi) const"
00412 );
00413 return tmp<scalarField>(NULL);
00414 }
00415
00416
00417 Foam::tmp<Foam::volScalarField> Foam::basicThermo::Cp() const
00418 {
00419 notImplemented("basicThermo::Cp() const");
00420 return volScalarField::null();
00421 }
00422
00423
00424 Foam::tmp<Foam::scalarField> Foam::basicThermo::Cv
00425 (
00426 const scalarField& T,
00427 const label patchi
00428 ) const
00429 {
00430 notImplemented
00431 (
00432 "basicThermo::Cv"
00433 "(const scalarField& T, const label patchi) const"
00434 );
00435 return tmp<scalarField>(NULL);
00436 }
00437
00438
00439 Foam::tmp<Foam::volScalarField> Foam::basicThermo::Cv() const
00440 {
00441 notImplemented("basicThermo::Cv() const");
00442 return volScalarField::null();
00443 }
00444
00445
00446 const Foam::volScalarField& Foam::basicThermo::mu() const
00447 {
00448 return mu_;
00449 }
00450
00451
00452 const Foam::volScalarField& Foam::basicThermo::alpha() const
00453 {
00454 return alpha_;
00455 }
00456
00457
00458 bool Foam::basicThermo::read()
00459 {
00460 return regIOobject::read();
00461 }
00462
00463
00464