/*---------------------------------------------------------------------------*\
========= |
\\ / 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 .
\*---------------------------------------------------------------------------*/
#include "interfaceOxideRate.H"
#include "cutCellIso.H"
#include "volPointInterpolation.H"
#include "timeVaryingMassSorptionFvPatchScalarField.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template
Foam::meltingEvaporationModels::interfaceOxideRate
::interfaceOxideRate
(
const dictionary& dict,
const phasePair& pair
)
:
InterfaceCompositionModel(dict, pair),
C_
(
dimensionedScalar
(
dimDensity/dimTime,
dict.getCheck("C", scalarMinMax::ge(0))
)
),
Tliquidus_
(
dimensionedScalar
(
dimTemperature,
dict.getCheck("Tliquidus", scalarMinMax::ge(0))
)
),
Tsolidus_
(
dimensionedScalar
(
dimTemperature,
dict.getCheck("Tsolidus", scalarMinMax::ge(0))
)
),
oxideCrit_
(
dimensionedScalar
(
dimDensity,
dict.getCheck("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("isoAlpha", 0.5))
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
template
Foam::tmp
Foam::meltingEvaporationModels::interfaceOxideRate::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 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 tSoxide =
max((oxideCrit_.value() - to)/oxideCrit_.value(), scalar(0));
// (CSC:Eq. 4)
tmp 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(alphab[patchi]))
{
const auto& pp =
refCast
(
alphab[patchi]
);
const labelUList& fc = mesh.boundary()[patchi].faceCells();
tmp tsb = pp.source();
auto tRhoto = volScalarField::New
(
"tRhoto",
IOobject::NO_REGISTER,
mesh,
dimensionedScalar(dimDensity, Zero)
);
auto& 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::New(mDotOxide_);
}
template
Foam::tmp
Foam::meltingEvaporationModels::interfaceOxideRate::KSp
(
label variable,
const volScalarField& refValue
)
{
return nullptr;
}
template
Foam::tmp
Foam::meltingEvaporationModels::interfaceOxideRate::KSu
(
label variable,
const volScalarField& refValue
)
{
return nullptr;
}
// ************************************************************************* //