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
15bade14
Commit
15bade14
authored
Aug 26, 2008
by
henry
Browse files
Yet more updates to the transonic formulation.
parent
99917f38
Changes
1
Hide whitespace changes
Inline
Side-by-side
applications/solvers/compressible/rhoSimpleFoam/pEqn.H
View file @
15bade14
...
...
@@ -16,12 +16,13 @@ if (transonic)
for
(
int
nonOrth
=
0
;
nonOrth
<=
nNonOrthCorr
;
nonOrth
++
)
{
fvScalarMatrix
pEqn
0
fvScalarMatrix
pEqn
(
fvm
::
div
(
phid
,
p
)
-
fvm
::
laplacian
(
rho
*
rUA
,
p
)
);
fvScalarMatrix
pEqn
=
pEqn0
;
// Relax the pressure equation to ensure diagonal-dominance
pEqn
.
relax
(
mesh
.
relaxationFactor
(
"pEqn"
));
pEqn
.
setReference
(
pRefCell
,
pRefValue
);
...
...
@@ -39,14 +40,13 @@ if (transonic)
if
(
nonOrth
==
nNonOrthCorr
)
{
phi
==
pEqn
0
.
flux
();
phi
==
pEqn
.
flux
();
}
}
}
else
{
phi
=
fvc
::
interpolate
(
rho
)
*
(
fvc
::
interpolate
(
U
)
&
mesh
.
Sf
());
//phi = fvc::interpolate(rho*U) & mesh.Sf();
closedVolume
=
adjustPhi
(
phi
,
U
,
p
);
for
(
int
nonOrth
=
0
;
nonOrth
<=
nNonOrthCorr
;
nonOrth
++
)
...
...
@@ -58,7 +58,7 @@ else
pEqn
.
setReference
(
pRefCell
,
pRefValue
);
//
r
etain the residual from the first iteration
//
R
etain the residual from the first iteration
if
(
nonOrth
==
0
)
{
eqnResidual
=
pEqn
.
solve
().
initialResidual
();
...
...
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