Consistency update: simpleFoam
This commit is contained in:
parent
ed03625c8b
commit
5d7793cf06
3 changed files with 21 additions and 31 deletions
|
@ -1,18 +1,13 @@
|
|||
// Solve the momentum equation
|
||||
|
||||
tmp<fvVectorMatrix> HUEqn
|
||||
tmp<fvVectorMatrix> UEqn
|
||||
(
|
||||
fvm::div(phi, U)
|
||||
+ turbulence->divDevReff()
|
||||
);
|
||||
|
||||
// Get under-relaxation factor
|
||||
const scalar UUrf = mesh.solutionDict().equationRelaxationFactor(U.name());
|
||||
// Relax the equation
|
||||
UEqn().relax();
|
||||
|
||||
// Momentum solution
|
||||
solve
|
||||
(
|
||||
relax(HUEqn(), UUrf)
|
||||
==
|
||||
-fvc::grad(p)
|
||||
);
|
||||
solve(UEqn() == -fvc::grad(p));
|
||||
|
|
|
@ -1,20 +1,11 @@
|
|||
p.boundaryField().updateCoeffs();
|
||||
|
||||
// Prepare clean 1/Ap without contribution from under-relaxation
|
||||
// HJ, 26/Oct/2015
|
||||
volScalarField rUA
|
||||
(
|
||||
"(1|A(U))",
|
||||
1/HUEqn().A()
|
||||
);
|
||||
volScalarField rAU = 1.0/UEqn().A();
|
||||
U = rAU*UEqn().H();
|
||||
UEqn.clear();
|
||||
|
||||
// Store velocity under-relaxation point before using U for
|
||||
// the flux precursor
|
||||
U.storePrevIter();
|
||||
|
||||
U = rUA*HUEqn().H();
|
||||
HUEqn.clear();
|
||||
phi = fvc::interpolate(U) & mesh.Sf();
|
||||
// Calculate under-relaxation consistent flux
|
||||
simple.calcSteadyConsistentFlux(phi, U);
|
||||
adjustPhi(phi, U, p);
|
||||
|
||||
// Non-orthogonal pressure corrector loop
|
||||
|
@ -22,7 +13,14 @@
|
|||
{
|
||||
fvScalarMatrix pEqn
|
||||
(
|
||||
fvm::laplacian(rUA, p) == fvc::div(phi)
|
||||
fvm::laplacian
|
||||
(
|
||||
fvc::interpolate(rAU)/simple.aCoeff(),
|
||||
p,
|
||||
"laplacian(rAU" + p.name() + ')'
|
||||
)
|
||||
==
|
||||
fvc::div(phi)
|
||||
);
|
||||
|
||||
pEqn.setReference(pRefCell, pRefValue);
|
||||
|
@ -40,9 +38,5 @@
|
|||
// Explicitly relax pressure for momentum corrector
|
||||
p.relax();
|
||||
|
||||
// Momentum corrector
|
||||
// Note: since under-relaxation does not change aU, H/a in U can be
|
||||
// re-used. HJ, 22/Jan/2016
|
||||
U = UUrf*(U - rUA*fvc::grad(p)) + (1 - UUrf)*U.prevIter();
|
||||
U.correctBoundaryConditions();
|
||||
|
||||
// Reconstruct consistent velocity after pressure equation
|
||||
simple.reconstructSteadyVelocity(U, rAU, p);
|
||||
|
|
|
@ -25,8 +25,9 @@ Application
|
|||
simpleFoam
|
||||
|
||||
Description
|
||||
Steady-state solver for incompressible, turbulent flow
|
||||
Steady-state solver for incompressible, turbulent flow.
|
||||
Consistent formulation without time-step and relaxation dependence by Jasak
|
||||
and Tukovic.
|
||||
|
||||
Author
|
||||
Hrvoje Jasak, Wikki Ltd. All rights reserved
|
||||
|
|
Reference in a new issue