ENH: interfaceOxideRate: New mass transfer model

Based on:

    Cao, L., Sun, F., Chen, T., Tang, Y., & Liao, D. (2018).
    Quantitative prediction of oxide inclusion defects inside
    the casting and on the walls during cast-filling processes.
    International Journal of Heat and Mass Transfer, 119, 614-623.
    DOI:10.1016/j.ijheatmasstransfer.2017.11.127

Co-authored-by: Kutalmis Bercin <kutalmis.bercin@esi-group.com>
This commit is contained in:
sergio 2021-09-03 13:27:23 -07:00 committed by Andrew Heather
parent 6580f1bbef
commit ea9d424e24
3 changed files with 478 additions and 1 deletions

View File

@ -5,7 +5,7 @@
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2017 OpenCFD Ltd.
Copyright (C) 2017-2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -51,6 +51,7 @@ License
#include "kineticGasEvaporation.H"
#include "Lee.H"
#include "interfaceHeatResistance.H"
#include "interfaceOxideRate.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -573,6 +574,37 @@ namespace Foam
hPowerSolidThermoPhysics
);
// interfaceOxideRate model definitions
//From pure phase (tabulated) to solid phase (const)
makeInterfacePureType
(
interfaceOxideRate,
heRhoThermo,
rhoThermo,
pureMixture,
tabulatedThermoPhysics,
heSolidThermo,
solidThermo,
pureMixture,
hConstSolidThermoPhysics
);
// From pure phase (rho const) to phase (rho const)
makeInterfacePureType
(
interfaceOxideRate,
heRhoThermo,
rhoThermo,
pureMixture,
constRhoHThermoPhysics,
heRhoThermo,
rhoThermo,
pureMixture,
constRhoHThermoPhysics
);
}

View File

@ -0,0 +1,207 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 "interfaceOxideRate.H"
#include "cutCellIso.H"
#include "volPointInterpolation.H"
#include "timeVaryingMassSorptionFvPatchScalarField.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class Thermo, class OtherThermo>
Foam::meltingEvaporationModels::interfaceOxideRate<Thermo, OtherThermo>
::interfaceOxideRate
(
const dictionary& dict,
const phasePair& pair
)
:
InterfaceCompositionModel<Thermo, OtherThermo>(dict, pair),
C_
(
dimensionedScalar
(
dimDensity/dimTime,
dict.getCheck<scalar>("C", scalarMinMax::ge(0))
)
),
Tliquidus_
(
dimensionedScalar
(
dimTemperature,
dict.getCheck<scalar>("Tliquidus", scalarMinMax::ge(0))
)
),
Tsolidus_
(
dimensionedScalar
(
dimTemperature,
dict.getCheck<scalar>("Tsolidus", scalarMinMax::ge(0))
)
),
oxideCrit_
(
dimensionedScalar
(
dimDensity,
dict.getCheck<scalar>("oxideCrit", scalarMinMax::ge(0))
)
),
mDotOxide_
(
IOobject
(
"mDotOxide",
this->mesh_.time().timeName(),
this->mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
this->mesh_,
dimensionedScalar(dimDensity/dimTime, Zero)
),
isoAlpha_(dict.getOrDefault<scalar>("isoAlpha", 0.5))
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
template<class Thermo, class OtherThermo>
Foam::tmp<Foam::volScalarField>
Foam::meltingEvaporationModels::interfaceOxideRate<Thermo, OtherThermo>::Kexp
(
const volScalarField& T
)
{
const volScalarField& from = this->pair().from();
const volScalarField& to = this->pair().to();
// (CSC:Eq. 2)
const fvMesh& mesh = this->mesh_;
scalarField ap
(
volPointInterpolation::New(mesh).interpolate(from)
);
cutCellIso cutCell(mesh, ap);
tmp<volScalarField> tSalpha = scalar(0)*from;
volScalarField& Salpha = tSalpha.ref();
forAll(Salpha, celli)
{
const label status = cutCell.calcSubCell(celli, isoAlpha_);
if (status == 0) // cell is cut
{
Salpha[celli] = scalar(1);
}
}
// (CSC:Eq. 5)
tmp<volScalarField> tSoxide =
max((oxideCrit_.value() - to)/oxideCrit_.value(), scalar(0));
// (CSC:Eq. 4)
tmp<volScalarField> tST =
Foam::exp
(
scalar(1)
- scalar(1)/max((T - Tsolidus_)/(Tliquidus_ - Tsolidus_),scalar(1e-6))
);
// (CSC:Eq. 6)
mDotOxide_ = C_*tSalpha*tSoxide*tST;
const volScalarField::Boundary& alphab = to.boundaryField();
forAll(alphab, patchi)
{
if (isA<timeVaryingMassSorptionFvPatchScalarField>(alphab[patchi]))
{
const auto& pp =
refCast<const timeVaryingMassSorptionFvPatchScalarField>
(
alphab[patchi]
);
const labelUList& fc = mesh.boundary()[patchi].faceCells();
tmp<scalarField> tsb = pp.source();
auto tRhoto = tmp<volScalarField>::New
(
IOobject
(
"tRhoto",
mesh.time().timeName(),
mesh
),
mesh,
dimensionedScalar(dimDensity, Zero)
);
volScalarField& rhoto = tRhoto.ref();
rhoto = this->pair().to().rho();
forAll(fc, faceI)
{
const label cellI = fc[faceI];
const scalar rhoI = rhoto[cellI];
mDotOxide_[cellI] += rhoI*tsb()[faceI];
}
}
}
return tmp<volScalarField>::New(mDotOxide_);
}
template<class Thermo, class OtherThermo>
Foam::tmp<Foam::volScalarField>
Foam::meltingEvaporationModels::interfaceOxideRate<Thermo, OtherThermo>::KSp
(
label variable,
const volScalarField& refValue
)
{
return nullptr;
}
template<class Thermo, class OtherThermo>
Foam::tmp<Foam::volScalarField>
Foam::meltingEvaporationModels::interfaceOxideRate<Thermo, OtherThermo>::KSu
(
label variable,
const volScalarField& refValue
)
{
return nullptr;
}
// ************************************************************************* //

View File

@ -0,0 +1,238 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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/>.
Class
Foam::meltingEvaporationModels::interfaceOxideRate
Description
The \c interfaceOxideRate is a simple model to calculate the formation
rate of oxide inclusions (\c mDotOxide), which is a local function of
the volume fraction of reducing agent (\c alpha), temperature (\c T)
and oxide-inclusion density (\c chi).
The oxide-inclusion formation rate is modelled as follows:
\f[
S_\chi = C_\chi S_\alpha S_T S_\rho
\f]
with
\f[
S_\alpha = \alpha (1 - \alpha)
\f]
\f[
S_T =
\exp{
1 - \frac{1}{max(T - T_{solidus}/(T_liquidus - T_solidus), 1e-6)}
}
\f]
\f[
S_\rho =
max \left(\frac{\chi_{crit} - \chi_{curr}}{\chi_{curr}}, 0\right)
\f]
where
\vartable
S_\chi | Oxide-inclusion formation rate [kg/(m^3 s)]
C_\chi | Oxide-inclusion formation rate constant [kg/(m^3 s)]
S_\alpha | Formation factor due to volume fraction of reducing agent [-]
S_T | Formation factor due to temperature [-]
S_\rho | Formation factor due to oxide-inclusion density [-]
\alpha | Volume fraction of reducing agent [-]
T | Local temperature [K]
T_{solidus} | Solidus temperature of reducing agent [K]
T_{liquidus} | Liquidus temperature of reducing agent [K]
\chi_{crit} | Critical oxide-inclusion density [kg/m^3]
\chi_{curr} | Current oxide-inclusion density [kg/m^3]
\endvartable
References:
\verbatim
Oxide-inclusion model (tag:CSC):
Cao, L., Sun, F., Chen, T., Tang, Y., & Liao, D. (2018).
Quantitative prediction of oxide inclusion defects inside
the casting and on the walls during cast-filling processes.
International Journal of Heat and Mass Transfer, 119, 614-623.
DOI:10.1016/j.ijheatmasstransfer.2017.11.127
\endverbatim
Usage
Minimal example by using \c constant/phaseProperties.massTransferModel:
\verbatim
massTransferModel
(
(liquid to oxide)
{
type interfaceOxideRate;
C <scalar>;
Tliquidus <scalar>;
Tsolidus <scalar>;
oxideCrit <scalar>;
isoAlpha <scalar>;
}
);
\endverbatim
where the entries mean:
\table
Property | Description | Type | Reqd | Deflt
type | Type name: interfaceOxideRate | word | yes | -
C | Oxide-inclusion formation rate constant | scalar | yes | -
Tliquidus | Liquidus temperature of reducing agent | scalar | yes | -
Tsolidus | Solidus temperature of reducing agent | scalar | yes | -
oxideCrit | Critical oxide-inclusion density | scalar | yes | -
isoAlpha | Location of the source | scalar | no | 0.5
\endtable
Note
- \c oxideCrit should be determined experimentally (CSC:p. 616).
- \c C should be determined by practical production (CSC:p. 616).
SourceFiles
interfaceOxideRate.C
\*---------------------------------------------------------------------------*/
#ifndef meltingEvaporationModels_interfaceOxideRate_H
#define meltingEvaporationModels_interfaceOxideRate_H
#include "InterfaceCompositionModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *//
namespace Foam
{
// Forward Declarations
class phasePair;
namespace meltingEvaporationModels
{
/*---------------------------------------------------------------------------*\
Class interfaceOxideRate Declaration
\*---------------------------------------------------------------------------*/
template<class Thermo, class OtherThermo>
class interfaceOxideRate
:
public InterfaceCompositionModel<Thermo, OtherThermo>
{
// Private Data
//- Oxide-inclusion formation rate constant
const dimensionedScalar C_;
//- Liquidus temperature of reducing agent (e.g. a casting metal)
const dimensionedScalar Tliquidus_;
//- Solidus temperature of reducing agent (e.g. a casting metal)
const dimensionedScalar Tsolidus_;
//- Critical oxide-inclusion density
const dimensionedScalar oxideCrit_;
//- Oxide-inclusion formation rate
volScalarField mDotOxide_;
//- Interface Iso-value
scalar isoAlpha_;
public:
//- Runtime type information
TypeName("interfaceOxideRate");
// Constructors
//- Construct from components
interfaceOxideRate
(
const dictionary& dict,
const phasePair& pair
);
//- Destructor
virtual ~interfaceOxideRate() = default;
// Member Functions
//- Explicit total mass transfer coefficient
virtual tmp<volScalarField> Kexp
(
const volScalarField& field
);
//- Implicit mass transfer coefficient
virtual tmp<volScalarField> KSp
(
label modelVariable,
const volScalarField& field
);
//- Explicit mass transfer coefficient
virtual tmp<volScalarField> KSu
(
label modelVariable,
const volScalarField& field
);
//- Return solidus temperature of reducing agent
virtual const dimensionedScalar& Tactivate() const noexcept
{
return Tsolidus_;
}
//- Add/subtract alpha*div(U) as a source term
//- for alpha, substituting div(U) = mDot(1/rho1 - 1/rho2)
virtual bool includeDivU() const noexcept
{
return true;
}
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace meltingEvaporationModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
# include "interfaceOxideRate.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //