Skip to content
Snippets Groups Projects
Commit 46c5a2b7 authored by andy's avatar andy
Browse files

ENH: Added radiation capability to buoyantPimpleFoam

parent 92261101
Branches
Tags
No related merge requests found
{
volScalarField& he = thermo.he();
fvScalarMatrix EEqn
(
fvm::ddt(rho, he) + fvm::div(phi, he)
+ fvc::ddt(rho, K) + fvc::div(phi, K)
+ (
he.name() == "e"
? fvc::div
(
fvc::absolute(phi/fvc::interpolate(rho), U),
p,
"div(phiv,p)"
)
: -dpdt
)
- fvm::laplacian(turbulence->alphaEff(), he)
==
sources(rho, he)
);
EEqn.relax();
sources.constrain(EEqn);
EEqn.solve();
thermo.correct();
}
EXE_INC = \
-I../../compressible/rhoPimpleFoam \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/turbulenceModels/compressible/turbulenceModel \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/sampling/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude \
-I$(LIB_SRC)/fieldSources/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/radiationModels/lnInclude \
-I$(LIB_SRC)/turbulenceModels/compressible/turbulenceModel
EXE_LIBS = \
-lfiniteVolume \
-lsampling \
-lmeshTools \
-lfieldSources \
-lfluidThermophysicalModels \
-lradiationModels \
-lspecie \
-lcompressibleTurbulenceModel \
-lcompressibleRASModels \
-lcompressibleLESModels \
-lfiniteVolume \
-lsampling \
-lmeshTools \
-lfieldSources
-lcompressibleLESModels
......@@ -36,6 +36,7 @@ Description
#include "fvCFD.H"
#include "rhoThermo.H"
#include "turbulenceModel.H"
#include "radiationModel.H"
#include "IObasicSourceList.H"
#include "pimpleControl.H"
......@@ -48,6 +49,7 @@ int main(int argc, char *argv[])
#include "createMesh.H"
#include "readGravitationalAcceleration.H"
#include "createFields.H"
#include "createRadiationModel.H"
#include "initContinuityErrs.H"
#include "readTimeControls.H"
#include "compressibleCourantNo.H"
......
Info<< "Reading thermophysical properties\n" << endl;
autoPtr<rhoThermo> pThermo
(
rhoThermo::New(mesh)
);
autoPtr<rhoThermo> pThermo(rhoThermo::New(mesh));
rhoThermo& thermo = pThermo();
thermo.validate(args.executable(), "h", "e");
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment