The interface for fvModels has been modified to improve its application
to "proxy" equations. That is, equations that are not straightforward
statements of conservation laws in OpenFOAM's usual convention.
A standard conservation law typically takes the following form:
fvMatrix<scalar> psiEqn
(
fvm::ddt(alpha, rho, psi)
+ <fluxes>
==
<sources>
);
A proxy equation, on the other hand, may be a derivation or
rearrangement of a law like this, and may be linearised in terms of a
different variable.
The pressure equation is the most common example of a proxy equation. It
represents a statement of the conservation of volume or mass, but it is
a rearrangement of the original continuity equation, and it has been
linearised in terms of a different variable; the pressure. Another
example is that in the pre-predictor of a VoF solver the
phase-continuity equation is constructed, but it is linearised in terms
of volume fraction rather than density.
In these situations, fvModels sources are now applied by calling:
fvModels().sourceProxy(<conserved-fields ...>, <equation-field>)
Where <conserved-fields ...> are (alpha, rho, psi), (rho, psi), just
(psi), or are omitted entirely (for volume continuity), and the
<equation-field> is the field associated with the proxy equation. This
produces a source term identical in value to the following call:
fvModels().source(<conserved-fields ...>)
It is only the linearisation in terms of <equation-field> that differs
between these two calls.
This change permits much greater flexibility in the handling of mass and
volume sources than the previous name-based system did. All the relevant
fields are available, dimensions can be used in the logic to determine
what sources are being constructed, and sources relating to a given
conservation law all share the same function.
This commit adds the functionality for injection-type sources in the
compressibleVoF solver. A following commit will add a volume source
model for use in incompressible solvers.
182 lines
5.5 KiB
C++
182 lines
5.5 KiB
C++
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration | Website: https://openfoam.org
|
|
\\ / A nd | Copyright (C) 2023 OpenFOAM Foundation
|
|
\\/ M anipulation |
|
|
-------------------------------------------------------------------------------
|
|
License
|
|
This file is part of OpenFOAM.
|
|
|
|
OpenFOAM is free software: you can redistribute it and/or modify it
|
|
under the terms of the GNU General Public License as published by
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
(at your option) any later version.
|
|
|
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
|
for more details.
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "compressibleMultiphaseVoF.H"
|
|
#include "constrainHbyA.H"
|
|
#include "constrainPressure.H"
|
|
#include "adjustPhi.H"
|
|
#include "findRefCell.H"
|
|
#include "fvcMeshPhi.H"
|
|
#include "fvcFlux.H"
|
|
#include "fvcDdt.H"
|
|
#include "fvcDiv.H"
|
|
#include "fvcSnGrad.H"
|
|
#include "fvcSup.H"
|
|
#include "fvcReconstruct.H"
|
|
#include "fvmLaplacian.H"
|
|
|
|
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
|
|
|
void Foam::solvers::compressibleMultiphaseVoF::pressureCorrector()
|
|
{
|
|
volVectorField& U = U_;
|
|
surfaceScalarField& phi(phi_);
|
|
|
|
fvVectorMatrix& UEqn = tUEqn.ref();
|
|
setrAU(UEqn);
|
|
|
|
const surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU()));
|
|
|
|
while (pimple.correct())
|
|
{
|
|
const volVectorField HbyA(constrainHbyA(rAU()*UEqn.H(), U, p_rgh));
|
|
surfaceScalarField phiHbyA
|
|
(
|
|
"phiHbyA",
|
|
fvc::flux(HbyA)
|
|
+ fvc::interpolate(rho*rAU())*fvc::ddtCorr(U, phi, Uf)
|
|
);
|
|
|
|
MRF.makeRelative(phiHbyA);
|
|
|
|
const surfaceScalarField phig
|
|
(
|
|
(
|
|
surfaceTensionForce()
|
|
- buoyancy.ghf*fvc::snGrad(rho)
|
|
)*rAUf*mesh.magSf()
|
|
);
|
|
|
|
phiHbyA += phig;
|
|
|
|
// Update the pressure BCs to ensure flux consistency
|
|
constrainPressure(p_rgh, U, phiHbyA, rAUf, MRF);
|
|
|
|
PtrList<fvScalarMatrix> p_rghEqnComps(phases.size());
|
|
|
|
forAll(phases, phasei)
|
|
{
|
|
const compressibleVoFphase& phase = phases[phasei];
|
|
const rhoFluidThermo& thermo = phase.thermo();
|
|
const volScalarField& rho = phases[phasei].thermo().rho();
|
|
|
|
p_rghEqnComps.set
|
|
(
|
|
phasei,
|
|
(
|
|
fvc::ddt(rho) + thermo.psi()*correction(fvm::ddt(p_rgh))
|
|
+ fvc::div(phi, rho) - fvc::Sp(fvc::div(phi), rho)
|
|
- fvModels().sourceProxy(phase, rho, p_rgh)
|
|
).ptr()
|
|
);
|
|
}
|
|
|
|
// Cache p_rgh prior to solve for density update
|
|
const volScalarField p_rgh_0(p_rgh);
|
|
|
|
while (pimple.correctNonOrthogonal())
|
|
{
|
|
fvScalarMatrix p_rghEqnIncomp
|
|
(
|
|
fvc::div(phiHbyA)
|
|
- fvm::laplacian(rAUf, p_rgh)
|
|
);
|
|
|
|
tmp<fvScalarMatrix> p_rghEqnComp;
|
|
|
|
forAll(phases, phasei)
|
|
{
|
|
const compressibleVoFphase& phase = phases[phasei];
|
|
|
|
tmp<fvScalarMatrix> p_rghEqnCompi
|
|
(
|
|
(max(phase, scalar(0))/phase.thermo().rho())
|
|
*p_rghEqnComps[phasei]
|
|
);
|
|
|
|
if (phasei == 0)
|
|
{
|
|
p_rghEqnComp = p_rghEqnCompi;
|
|
}
|
|
else
|
|
{
|
|
p_rghEqnComp.ref() += p_rghEqnCompi;
|
|
}
|
|
}
|
|
|
|
{
|
|
fvScalarMatrix p_rghEqn(p_rghEqnComp + p_rghEqnIncomp);
|
|
|
|
fvConstraints().constrain(p_rghEqn);
|
|
|
|
p_rghEqn.solve();
|
|
}
|
|
|
|
if (pimple.finalNonOrthogonalIter())
|
|
{
|
|
forAll(phases, phasei)
|
|
{
|
|
compressibleVoFphase& phase = phases[phasei];
|
|
|
|
phase.dgdt() =
|
|
pos0(phase)
|
|
*(p_rghEqnComps[phasei] & p_rgh)/phase.thermo().rho();
|
|
}
|
|
|
|
phi = phiHbyA + p_rghEqnIncomp.flux();
|
|
|
|
p = p_rgh + rho*buoyancy.gh;
|
|
fvConstraints().constrain(p);
|
|
p_rgh = p - rho*buoyancy.gh;
|
|
p_rgh.correctBoundaryConditions();
|
|
|
|
U = HbyA
|
|
+ rAU()*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
|
|
U.correctBoundaryConditions();
|
|
fvConstraints().constrain(U);
|
|
}
|
|
}
|
|
|
|
// Update densities from change in p_rgh
|
|
mixture.correctRho(p_rgh - p_rgh_0);
|
|
mixture.correct();
|
|
|
|
// Correct p_rgh for consistency with p and the updated densities
|
|
p_rgh = p - rho*buoyancy.gh;
|
|
p_rgh.correctBoundaryConditions();
|
|
}
|
|
|
|
// Correct Uf if the mesh is moving
|
|
fvc::correctUf(Uf, U, fvc::absolute(phi, U), MRF);
|
|
|
|
K = 0.5*magSqr(U);
|
|
|
|
clearrAU();
|
|
tUEqn.clear();
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|