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
91261b30
Commit
91261b30
authored
Nov 18, 2015
by
Henry Weller
Browse files
reactingFoam: Added support for PIMPLE-consistent and pressure relaxation
Pressure relaxation is useful with LTS to damp acoustic waves
parent
3a5cb136
Changes
6
Hide whitespace changes
Inline
Side-by-side
applications/solvers/combustion/reactingFoam/UEqn.H
View file @
91261b30
MRF
.
correctBoundaryVelocity
(
U
);
// Solve the Momentum equation
fvVectorMatrix
UEqn
(
fvm
::
ddt
(
rho
,
U
)
+
fvm
::
div
(
phi
,
U
)
+
MRF
.
DDt
(
rho
,
U
)
+
turbulence
->
divDevRhoReff
(
U
)
==
rho
*
g
+
fvOptions
(
rho
,
U
)
);
MRF
.
correctBoundaryVelocity
(
U
);
UEqn
.
relax
();
tmp
<
fvVectorMatrix
>
UEqn
(
fvm
::
ddt
(
rho
,
U
)
+
fvm
::
div
(
phi
,
U
)
+
MRF
.
DDt
(
rho
,
U
)
+
turbulence
->
divDevRhoReff
(
U
)
==
rho
*
g
+
fvOptions
(
rho
,
U
)
);
fvOptions
.
constrain
(
UEqn
);
UEqn
().
relax
(
);
if
(
pimple
.
momentumPredictor
())
{
solve
(
UEqn
==
-
fvc
::
grad
(
p
));
fvOptions
.
constrain
(
UEqn
());
fvOptions
.
correct
(
U
);
K
=
0
.
5
*
magSqr
(
U
);
}
if
(
pimple
.
momentumPredictor
())
{
solve
(
UEqn
()
==
-
fvc
::
grad
(
p
));
fvOptions
.
correct
(
U
);
K
=
0
.
5
*
magSqr
(
U
);
}
applications/solvers/combustion/reactingFoam/createFields.H
View file @
91261b30
...
...
@@ -45,6 +45,28 @@ const volScalarField& T = thermo.T();
#include
"compressibleCreatePhi.H"
dimensionedScalar
rhoMax
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMax"
,
pimple
.
dict
(),
dimDensity
,
GREAT
)
);
dimensionedScalar
rhoMin
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMin"
,
pimple
.
dict
(),
dimDensity
,
0
)
);
mesh
.
setFluxRequired
(
p
.
name
());
Info
<<
"Creating turbulence model.
\n
"
<<
nl
;
...
...
applications/solvers/combustion/reactingFoam/pEqn.H
View file @
91261b30
rho
=
thermo
.
rho
();
rho
=
max
(
rho
,
rhoMin
);
rho
=
min
(
rho
,
rhoMax
);
rho
.
relax
();
volScalarField
rAU
(
1
.
0
/
UEqn
.
A
());
volScalarField
rAU
(
1
.
0
/
UEqn
()
.
A
());
surfaceScalarField
rhorAUf
(
"rhorAUf"
,
fvc
::
interpolate
(
rho
*
rAU
));
volVectorField
HbyA
(
"HbyA"
,
U
);
HbyA
=
rAU
*
UEqn
.
H
();
HbyA
=
rAU
*
UEqn
().
H
();
if
(
pimple
.
nCorrPISO
()
<=
1
)
{
UEqn
.
clear
();
}
if
(
pimple
.
transonic
())
{
...
...
@@ -26,7 +34,7 @@ if (pimple.transonic())
(
fvm
::
ddt
(
psi
,
p
)
+
fvm
::
div
(
phid
,
p
)
-
fvm
::
laplacian
(
rho
*
rAU
,
p
)
-
fvm
::
laplacian
(
rhorAU
f
,
p
)
==
fvOptions
(
psi
,
p
,
rho
.
name
())
);
...
...
@@ -58,7 +66,7 @@ else
(
fvm
::
ddt
(
psi
,
p
)
+
fvc
::
div
(
phiHbyA
)
-
fvm
::
laplacian
(
rho
*
rAU
,
p
)
-
fvm
::
laplacian
(
rhorAU
f
,
p
)
==
fvOptions
(
psi
,
p
,
rho
.
name
())
);
...
...
@@ -75,6 +83,17 @@ else
#include
"rhoEqn.H"
#include
"compressibleContinuityErrs.H"
// Explicitly relax pressure for momentum corrector
p
.
relax
();
// Recalculate density from the relaxed pressure
rho
=
thermo
.
rho
();
rho
=
max
(
rho
,
rhoMin
);
rho
=
min
(
rho
,
rhoMax
);
rho
.
relax
();
Info
<<
"rho max/min : "
<<
max
(
rho
).
value
()
<<
" "
<<
min
(
rho
).
value
()
<<
endl
;
U
=
HbyA
-
rAU
*
fvc
::
grad
(
p
);
U
.
correctBoundaryConditions
();
fvOptions
.
correct
(
U
);
...
...
applications/solvers/combustion/reactingFoam/pcEqn.H
0 → 100644
View file @
91261b30
rho
=
thermo
.
rho
();
rho
=
max
(
rho
,
rhoMin
);
rho
=
min
(
rho
,
rhoMax
);
rho
.
relax
();
volScalarField
rAU
(
1
.
0
/
UEqn
().
A
());
volScalarField
rAtU
(
1
.
0
/
(
1
.
0
/
rAU
-
UEqn
().
H1
()));
volVectorField
HbyA
(
"HbyA"
,
U
);
HbyA
=
rAU
*
UEqn
().
H
();
if
(
pimple
.
nCorrPISO
()
<=
1
)
{
UEqn
.
clear
();
}
if
(
pimple
.
transonic
())
{
surfaceScalarField
phid
(
"phid"
,
fvc
::
interpolate
(
psi
)
*
(
(
fvc
::
interpolate
(
HbyA
)
&
mesh
.
Sf
())
+
fvc
::
interpolate
(
rho
*
rAU
)
*
fvc
::
ddtCorr
(
rho
,
U
,
phi
)
/
fvc
::
interpolate
(
rho
)
)
);
MRF
.
makeRelative
(
fvc
::
interpolate
(
psi
),
phid
);
surfaceScalarField
phic
(
"phic"
,
fvc
::
interpolate
(
rho
*
(
rAtU
-
rAU
))
*
fvc
::
snGrad
(
p
)
*
mesh
.
magSf
()
);
HbyA
-=
(
rAU
-
rAtU
)
*
fvc
::
grad
(
p
);
volScalarField
rhorAtU
(
"rhorAtU"
,
rho
*
rAtU
);
while
(
pimple
.
correctNonOrthogonal
())
{
fvScalarMatrix
pEqn
(
fvm
::
ddt
(
psi
,
p
)
+
fvm
::
div
(
phid
,
p
)
+
fvc
::
div
(
phic
)
-
fvm
::
laplacian
(
rhorAtU
,
p
)
==
fvOptions
(
psi
,
p
,
rho
.
name
())
);
pEqn
.
solve
(
mesh
.
solver
(
p
.
select
(
pimple
.
finalInnerIter
())));
if
(
pimple
.
finalNonOrthogonalIter
())
{
phi
==
phic
+
pEqn
.
flux
();
}
}
}
else
{
surfaceScalarField
phiHbyA
(
"phiHbyA"
,
(
(
fvc
::
interpolate
(
rho
*
HbyA
)
&
mesh
.
Sf
())
+
fvc
::
interpolate
(
rho
*
rAU
)
*
fvc
::
ddtCorr
(
rho
,
U
,
phi
)
)
);
MRF
.
makeRelative
(
fvc
::
interpolate
(
rho
),
phiHbyA
);
phiHbyA
+=
fvc
::
interpolate
(
rho
*
(
rAtU
-
rAU
))
*
fvc
::
snGrad
(
p
)
*
mesh
.
magSf
();
HbyA
-=
(
rAU
-
rAtU
)
*
fvc
::
grad
(
p
);
volScalarField
rhorAtU
(
"rhorAtU"
,
rho
*
rAtU
);
while
(
pimple
.
correctNonOrthogonal
())
{
fvScalarMatrix
pEqn
(
fvm
::
ddt
(
psi
,
p
)
+
fvc
::
div
(
phiHbyA
)
-
fvm
::
laplacian
(
rhorAtU
,
p
)
==
fvOptions
(
psi
,
p
,
rho
.
name
())
);
pEqn
.
solve
(
mesh
.
solver
(
p
.
select
(
pimple
.
finalInnerIter
())));
if
(
pimple
.
finalNonOrthogonalIter
())
{
phi
=
phiHbyA
+
pEqn
.
flux
();
}
}
}
#include
"rhoEqn.H"
#include
"compressibleContinuityErrs.H"
// Explicitly relax pressure for momentum corrector
p
.
relax
();
U
=
HbyA
-
rAtU
*
fvc
::
grad
(
p
);
U
.
correctBoundaryConditions
();
fvOptions
.
correct
(
U
);
K
=
0
.
5
*
magSqr
(
U
);
if
(
thermo
.
dpdt
())
{
dpdt
=
fvc
::
ddt
(
p
);
}
// Recalculate density from the relaxed pressure
rho
=
thermo
.
rho
();
rho
=
max
(
rho
,
rhoMin
);
rho
=
min
(
rho
,
rhoMax
);
if
(
!
pimple
.
transonic
())
{
rho
.
relax
();
}
Info
<<
"rho max/min : "
<<
max
(
rho
).
value
()
<<
" "
<<
min
(
rho
).
value
()
<<
endl
;
applications/solvers/combustion/reactingFoam/reactingFoam.C
View file @
91261b30
...
...
@@ -95,7 +95,14 @@ int main(int argc, char *argv[])
// --- Pressure corrector loop
while
(
pimple
.
correct
())
{
#include
"pEqn.H"
if
(
pimple
.
consistent
())
{
#include
"pcEqn.H"
}
else
{
#include
"pEqn.H"
}
}
if
(
pimple
.
turbCorr
())
...
...
applications/solvers/combustion/reactingFoam/setRDeltaT.H
View file @
91261b30
...
...
@@ -116,6 +116,9 @@ License
);
}
// Update tho boundary values of the reciprocal time-step
rDeltaT
.
correctBoundaryConditions
();
Info
<<
" Overall = "
<<
gMin
(
1
/
rDeltaT
.
internalField
())
<<
", "
<<
gMax
(
1
/
rDeltaT
.
internalField
())
<<
endl
;
...
...
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