buoyantBoussinesqSimpleFoam: Add support for radiative heat-transfer consistent with buoyantBoussinesqPimpleFoam

Patch provided by Daniel Jasinski: http://www.openfoam.org/mantisbt/view.php?id=1856
This commit is contained in:
Henry Weller
2015-10-30 13:58:17 +00:00
parent 7fdf0ff095
commit c6fe72c6ad
3 changed files with 8 additions and 1 deletions

View File

@ -3,6 +3,7 @@ EXE_INC = \
-I$(LIB_SRC)/TurbulenceModels/incompressible/lnInclude \
-I$(LIB_SRC)/transportModels \
-I$(LIB_SRC)/transportModels/incompressible/singlePhaseTransportModel \
-I$(LIB_SRC)/thermophysicalModels/radiation/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/sampling/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude \
@ -12,6 +13,7 @@ EXE_LIBS = \
-lturbulenceModels \
-lincompressibleTurbulenceModels \
-lincompressibleTransportModels \
-lradiationModels \
-lfiniteVolume \
-lsampling \
-lmeshTools \

View File

@ -9,7 +9,8 @@
fvm::div(phi, T)
- fvm::laplacian(alphaEff, T)
==
fvOptions(T)
radiation->ST(rhoCpRef, T)
+ fvOptions(T)
);
TEqn.relax();
@ -18,6 +19,8 @@
TEqn.solve();
radiation->correct();
fvOptions.correct(T);
rhok = 1.0 - beta*(T - TRef);

View File

@ -48,6 +48,7 @@ Description
#include "fvCFD.H"
#include "singlePhaseTransportModel.H"
#include "turbulentTransportModel.H"
#include "radiationModel.H"
#include "fvIOoptionList.H"
#include "simpleControl.H"
#include "fixedFluxPressureFvPatchScalarField.H"
@ -63,6 +64,7 @@ int main(int argc, char *argv[])
simpleControl simple(mesh);
#include "createFields.H"
#include "createIncompressibleRadiationModel.H"
#include "createMRF.H"
#include "createFvOptions.H"
#include "initContinuityErrs.H"