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
ef95b084
Commit
ef95b084
authored
Nov 29, 2012
by
andy
Browse files
ENH: Added MRF and porosity functionality to pimpleFoam
parent
f5d95cd5
Changes
4
Hide whitespace changes
Inline
Side-by-side
applications/solvers/incompressible/pimpleFoam/UEqn.H
View file @
ef95b084
...
...
@@ -7,6 +7,10 @@ tmp<fvVectorMatrix> UEqn
+
turbulence
->
divDevReff
(
U
)
);
mrfZones
.
addCoriolis
(
UEqn
());
pZones
.
addResistance
(
UEqn
());
UEqn
().
relax
();
sources
.
constrain
(
UEqn
());
...
...
applications/solvers/incompressible/pimpleFoam/createZones.H
0 → 100644
View file @
ef95b084
IOMRFZoneList
mrfZones
(
mesh
);
mrfZones
.
correctBoundaryVelocity
(
U
);
IOporosityModelList
pZones
(
mesh
);
applications/solvers/incompressible/pimpleFoam/pEqn.H
View file @
ef95b084
...
...
@@ -15,6 +15,8 @@ surfaceScalarField phiHbyA
adjustPhi
(
phiHbyA
,
U
,
p
);
mrfZones
.
relativeFlux
(
phiHbyA
);
// Non-orthogonal pressure corrector loop
while
(
pimple
.
correctNonOrthogonal
())
{
...
...
applications/solvers/incompressible/pimpleFoam/pimpleFoam.C
View file @
ef95b084
...
...
@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011
-2012
OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
...
...
@@ -28,7 +28,10 @@ Description
Large time-step transient solver for incompressible, flow using the PIMPLE
(merged PISO-SIMPLE) algorithm.
Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
Sub-models include:
- turbulence modelling, i.e. laminar, RAS or LES
- porosity (explicit treatment)
- Multiple Reference Frame (MRF)
\*---------------------------------------------------------------------------*/
...
...
@@ -37,6 +40,8 @@ Description
#include
"turbulenceModel.H"
#include
"pimpleControl.H"
#include
"IObasicSourceList.H"
#include
"IOporosityModelList.H"
#include
"IOMRFZoneList.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
...
...
@@ -46,6 +51,7 @@ int main(int argc, char *argv[])
#include
"createTime.H"
#include
"createMesh.H"
#include
"createFields.H"
#include
"createZones.H"
#include
"initContinuityErrs.H"
pimpleControl
pimple
(
mesh
);
...
...
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