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.
85 lines
2.9 KiB
C++
85 lines
2.9 KiB
C++
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration | Website: https://openfoam.org
|
|
\\ / A nd | Copyright (C) 2022-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 "compressibleVoF.H"
|
|
#include "fvcDiv.H"
|
|
|
|
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
|
|
|
void Foam::solvers::compressibleVoF::alphaSuSp
|
|
(
|
|
tmp<volScalarField::Internal>& tSu,
|
|
tmp<volScalarField::Internal>& tSp
|
|
)
|
|
{
|
|
const dimensionedScalar Szero(dgdt.dimensions(), 0);
|
|
|
|
tSp = volScalarField::Internal::New("Sp", mesh, Szero);
|
|
tSu = volScalarField::Internal::New("Su", mesh, Szero);
|
|
|
|
volScalarField::Internal& Sp = tSp.ref();
|
|
volScalarField::Internal& Su = tSu.ref();
|
|
|
|
if (fvModels().addsSupToField(mixture.rho1().name()))
|
|
{
|
|
const volScalarField::Internal alpha2ByRho1(alpha2()/mixture.rho1()());
|
|
const fvScalarMatrix alphaRho1Sup
|
|
(
|
|
fvModels().sourceProxy(alpha1, mixture.rho1(), alpha1)
|
|
);
|
|
|
|
Su += alpha2ByRho1*alphaRho1Sup.Su();
|
|
Sp += alpha2ByRho1*alphaRho1Sup.Sp();
|
|
}
|
|
|
|
if (fvModels().addsSupToField(mixture.rho2().name()))
|
|
{
|
|
const volScalarField::Internal alpha1ByRho2(alpha1()/mixture.rho2()());
|
|
const fvScalarMatrix alphaRho2Sup
|
|
(
|
|
fvModels().sourceProxy(alpha2, mixture.rho2(), alpha2)
|
|
);
|
|
|
|
Su -= alpha1ByRho2*(alphaRho2Sup.Su() + alphaRho2Sup.Sp());
|
|
Sp += alpha1ByRho2*alphaRho2Sup.Sp();
|
|
}
|
|
|
|
forAll(dgdt, celli)
|
|
{
|
|
if (dgdt[celli] > 0.0)
|
|
{
|
|
Sp[celli] -= dgdt[celli]/max(1.0 - alpha1[celli], 1e-4);
|
|
Su[celli] += dgdt[celli]/max(1.0 - alpha1[celli], 1e-4);
|
|
}
|
|
else if (dgdt[celli] < 0.0)
|
|
{
|
|
Sp[celli] += dgdt[celli]/max(alpha1[celli], 1e-4);
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|