isothermalFluid::correctBuoyantPressure: Relax the net force momentum equation source

rather than p_rgh which introduces an imbalance between the pressure and
buoyancy forces.  The relaxation factor for p_rgh specified in fvSolution is
used to relax the net force as the intent is to relax the pressure and this
provides convenient usage and backwards-compatibility.

Optional relaxation for the thermodynamic pressure p is also available for case
this provides convergence benefit for steady cases by relaxing the pressure work
term is the energy equation.

Note these changes only relate to the operation of the isothermalFluid solver
module for buoyant cases.
This commit is contained in:
Henry Weller
2022-09-22 15:05:33 +01:00
parent f4ac5f8748
commit c22a5a7aa6
6 changed files with 123 additions and 38 deletions

View File

@ -102,8 +102,9 @@ void Foam::solvers::isothermalFluid::correctBuoyantPressure()
*fvc::relative(phiHbyA, rho, U)
);
const fvScalarMatrix divPhidp(fvm::div(phid, p));
phiHbyA -= divPhidp.flux();
// Subtract the compressible part
// The resulting flux will be zero for a perfect gas
phiHbyA -= fvc::interpolate(psi*p)*phiHbyA/fvc::interpolate(rho);
if (pimple.consistent())
{
@ -120,7 +121,7 @@ void Foam::solvers::isothermalFluid::correctBuoyantPressure()
fvScalarMatrix p_rghDDtEqn
(
fvc::ddt(rho) + psi*correction(fvm::ddt(p_rgh))
+ fvc::div(phiHbyA) + divPhidp
+ fvc::div(phiHbyA) + fvm::div(phid, p)
==
fvModels().source(psi, p_rgh, rho.name())
);
@ -183,6 +184,21 @@ void Foam::solvers::isothermalFluid::correctBuoyantPressure()
phi = phiHbyA + p_rghEqn.flux();
// Calculate and relax the net pressure-buoyancy force
netForce.ref().relax
(
fvc::reconstruct((ghGradRhof + p_rghEqn.flux()/rhorAAtUf)),
p_rgh.relaxationFactor()
);
// Correct the momentum source with the pressure gradient flux
// calculated from the relaxed pressure
U = HbyA + rAAtU*netForce();
U.correctBoundaryConditions();
fvConstraints().constrain(U);
K = 0.5*magSqr(U);
if (!mesh.schemes().steady())
{
p = p_rgh + rho*gh + pRef;
@ -202,19 +218,8 @@ void Foam::solvers::isothermalFluid::correctBuoyantPressure()
continuityErrors();
// Explicitly relax pressure for momentum corrector
p_rgh.relax();
p = p_rgh + rho*gh + pRef;
// Correct the momentum source with the pressure gradient flux
// calculated from the relaxed pressure
U = HbyA + rAAtU*fvc::reconstruct((ghGradRhof + p_rghEqn.flux()/rhorAAtUf));
U.correctBoundaryConditions();
fvConstraints().constrain(U);
K = 0.5*magSqr(U);
if (mesh.schemes().steady())
{
if (fvConstraints().constrain(p))
@ -234,6 +239,9 @@ void Foam::solvers::isothermalFluid::correctBuoyantPressure()
p_rgh.correctBoundaryConditions();
}
// Optionally relax pressure for the thermophysics
p.relax();
if (mesh.schemes().steady() || pimple.simpleRho() || adjustMass)
{
rho = thermo.rho();

View File

@ -170,6 +170,21 @@ Foam::solvers::isothermalFluid::isothermalFluid
thermo,
pimple.dict()
);
netForce = new volVectorField
(
IOobject
(
"netForce",
runTime.timeName(),
mesh
),
fvc::reconstruct
(
(-buoyancy->ghf*fvc::snGrad(rho) - fvc::snGrad(p_rgh))
*mesh.magSf()
)
);
}
if (mesh.dynamic())

View File

@ -147,6 +147,10 @@ protected:
// Cached temporary fields
//- Momentum equation net force source term
// Used for buoyant simulations only
tmp<volVectorField> netForce;
//- Pointer to the surface momentum field
// used to recreate the flux after mesh-change
autoPtr<surfaceVectorField> rhoUf;

View File

@ -51,13 +51,7 @@ void Foam::solvers::isothermalFluid::momentumPredictor()
(
UEqn
==
fvc::reconstruct
(
(
- buoyancy->ghf*fvc::snGrad(rho)
- fvc::snGrad(p_rgh)
)*mesh.magSf()
)
netForce()
);
}
else

View File

@ -1274,18 +1274,23 @@ bool Foam::GeometricField<Type, PatchField, GeoMesh>::needReference() const
template<class Type, template<class> class PatchField, class GeoMesh>
void Foam::GeometricField<Type, PatchField, GeoMesh>::relax(const scalar alpha)
{
if (debug)
if (alpha < 1)
{
InfoInFunction
<< "Relaxing" << endl << this->info() << " by " << alpha << endl;
}
if (debug)
{
InfoInFunction
<< "Relaxing" << endl << this->info()
<< " by " << alpha << endl;
}
operator==(prevIter() + alpha*(*this - prevIter()));
operator==(prevIter() + alpha*(*this - prevIter()));
}
}
template<class Type, template<class> class PatchField, class GeoMesh>
void Foam::GeometricField<Type, PatchField, GeoMesh>::relax()
Foam::scalar
Foam::GeometricField<Type, PatchField, GeoMesh>::relaxationFactor() const
{
if
(
@ -1297,18 +1302,61 @@ void Foam::GeometricField<Type, PatchField, GeoMesh>::relax()
&& this->mesh().solution().relaxField(this->name() + "Final")
)
{
relax
return this->mesh().solution().fieldRelaxationFactor
(
this->mesh().solution().fieldRelaxationFactor
(
this->name() + "Final"
)
this->name() + "Final"
);
}
else if (this->mesh().solution().relaxField(this->name()))
{
relax(this->mesh().solution().fieldRelaxationFactor(this->name()));
return this->mesh().solution().fieldRelaxationFactor(this->name());
}
else
{
return 1;
}
}
template<class Type, template<class> class PatchField, class GeoMesh>
void Foam::GeometricField<Type, PatchField, GeoMesh>::relax()
{
relax(relaxationFactor());
}
template<class Type, template<class> class PatchField, class GeoMesh>
void Foam::GeometricField<Type, PatchField, GeoMesh>::relax
(
const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf,
const scalar alpha
)
{
if (alpha < 1)
{
if (debug)
{
InfoInFunction
<< "Relaxing" << endl << this->info()
<< " by " << alpha << endl;
}
operator==(*this + alpha*(tgf - *this));
}
else
{
operator==(tgf);
}
}
template<class Type, template<class> class PatchField, class GeoMesh>
void Foam::GeometricField<Type, PatchField, GeoMesh>::relax
(
const tmp<GeometricField<Type, PatchField, GeoMesh>>& tgf
)
{
relax(tgf, relaxationFactor());
}

View File

@ -462,15 +462,31 @@ public:
tmp<GeometricField<Type, PatchField, GeoMesh>> T() const;
//- Relax field (for steady-state solution).
// alpha = 1 : no relaxation
// alpha < 1 : relaxation
// alpha = 0 : do nothing
// alpha >= 1 : no relaxation
// alpha < 1 : relaxation
void relax(const scalar alpha);
//- Relax field (for steady-state solution).
// alpha is read from controlDict
//- Return the field relaxation factor read from fvSolution
// or 1 if not specified
scalar relaxationFactor() const;
//- Relax current field with respect to the cached previous iteration.
// Relaxation factor is read from fvSolution
void relax();
//- Relax given field with respect to the current field
// and reset the field to the result
void relax
(
const tmp<GeometricField<Type, PatchField, GeoMesh>>&,
const scalar alpha
);
//- Relax given field with respect to the current field
// and reset the field to the result
// Relaxation factor is read from fvSolution
void relax(const tmp<GeometricField<Type, PatchField, GeoMesh>>&);
//- Select the final iteration parameters if `final' is true
// by returning the field name + "Final"
// otherwise the standard parameters by returning the field name