Skip to content
GitLab
Menu
Projects
Groups
Snippets
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Development
openfoam
Commits
46c5a2b7
Commit
46c5a2b7
authored
Dec 06, 2012
by
andy
Browse files
ENH: Added radiation capability to buoyantPimpleFoam
parent
92261101
Changes
4
Hide whitespace changes
Inline
Side-by-side
applications/solvers/heatTransfer/buoyantPimpleFoam/EEqn.H
0 → 100644
View file @
46c5a2b7
{
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
();
}
applications/solvers/heatTransfer/buoyantPimpleFoam/Make/options
View file @
46c5a2b7
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
applications/solvers/heatTransfer/buoyantPimpleFoam/buoyantPimpleFoam.C
View file @
46c5a2b7
...
...
@@ -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"
...
...
applications/solvers/heatTransfer/buoyantPimpleFoam/createFields.H
View file @
46c5a2b7
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"
);
...
...
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment