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
12ec31f0
Commit
12ec31f0
authored
Dec 05, 2012
by
andy
Browse files
ENH: Updated sources handling in pimpleFoam
parent
ef1d4a2a
Changes
4
Hide whitespace changes
Inline
Side-by-side
applications/solvers/incompressible/pimpleFoam/UEqn.H
View file @
12ec31f0
...
...
@@ -5,12 +5,10 @@ tmp<fvVectorMatrix> UEqn
fvm
::
ddt
(
U
)
+
fvm
::
div
(
phi
,
U
)
+
turbulence
->
divDevReff
(
U
)
==
sources
(
U
)
);
mrfZones
.
addCoriolis
(
UEqn
());
pZones
.
addResistance
(
UEqn
());
UEqn
().
relax
();
sources
.
constrain
(
UEqn
());
...
...
@@ -19,5 +17,5 @@ volScalarField rAU(1.0/UEqn().A());
if
(
pimple
.
momentumPredictor
())
{
solve
(
UEqn
()
==
-
fvc
::
grad
(
p
)
+
sources
(
U
)
);
solve
(
UEqn
()
==
-
fvc
::
grad
(
p
));
}
applications/solvers/incompressible/pimpleFoam/createZones.H
deleted
100644 → 0
View file @
ef1d4a2a
IOMRFZoneList
mrfZones
(
mesh
);
mrfZones
.
correctBoundaryVelocity
(
U
);
IOporosityModelList
pZones
(
mesh
);
applications/solvers/incompressible/pimpleFoam/pEqn.H
View file @
12ec31f0
volVectorField
HbyA
(
"HbyA"
,
U
);
HbyA
=
rAU
*
(
UEqn
()
==
sources
(
U
))()
.
H
();
HbyA
=
rAU
*
UEqn
().
H
();
if
(
pimple
.
nCorrPISO
()
<=
1
)
{
...
...
@@ -15,7 +15,7 @@ surfaceScalarField phiHbyA
adjustPhi
(
phiHbyA
,
U
,
p
);
mrfZon
es
.
relativeFlux
(
phiHbyA
);
sourc
es
.
relativeFlux
(
phiHbyA
);
// Non-orthogonal pressure corrector loop
while
(
pimple
.
correctNonOrthogonal
())
...
...
applications/solvers/incompressible/pimpleFoam/pimpleFoam.C
View file @
12ec31f0
...
...
@@ -30,8 +30,7 @@ Description
Sub-models include:
- turbulence modelling, i.e. laminar, RAS or LES
- porosity (explicit treatment)
- Multiple Reference Frame (MRF)
- run-time selectable sources, e.g. MRF, explicit porosity
\*---------------------------------------------------------------------------*/
...
...
@@ -51,7 +50,6 @@ 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