mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: enthalpySorption: new temperature boundary condition
This commit is contained in:
@ -19,6 +19,7 @@ derivedFvPatchFields/fixedUnburntEnthalpy/fixedUnburntEnthalpyFvPatchScalarField
|
|||||||
derivedFvPatchFields/gradientUnburntEnthalpy/gradientUnburntEnthalpyFvPatchScalarField.C
|
derivedFvPatchFields/gradientUnburntEnthalpy/gradientUnburntEnthalpyFvPatchScalarField.C
|
||||||
derivedFvPatchFields/mixedUnburntEnthalpy/mixedUnburntEnthalpyFvPatchScalarField.C
|
derivedFvPatchFields/mixedUnburntEnthalpy/mixedUnburntEnthalpyFvPatchScalarField.C
|
||||||
derivedFvPatchFields/speciesSorption/speciesSorptionFvPatchScalarField.C
|
derivedFvPatchFields/speciesSorption/speciesSorptionFvPatchScalarField.C
|
||||||
|
derivedFvPatchFields/enthalpySorption/enthalpySorptionFvPatchScalarField.C
|
||||||
derivedFvPatchFields/boundarySourcePatch/boundarySourcePatch.C
|
derivedFvPatchFields/boundarySourcePatch/boundarySourcePatch.C
|
||||||
|
|
||||||
functionObjects/moleFractions/moleFractionsFunctionObjects.C
|
functionObjects/moleFractions/moleFractionsFunctionObjects.C
|
||||||
|
|||||||
@ -0,0 +1,350 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 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 "enthalpySorptionFvPatchScalarField.H"
|
||||||
|
#include "addToRunTimeSelectionTable.H"
|
||||||
|
#include "fvPatchFieldMapper.H"
|
||||||
|
#include "volFields.H"
|
||||||
|
#include "rhoReactionThermo.H"
|
||||||
|
#include "speciesSorptionFvPatchScalarField.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
const Foam::Enum
|
||||||
|
<
|
||||||
|
Foam::enthalpySorptionFvPatchScalarField::enthalpyModelType
|
||||||
|
>
|
||||||
|
Foam::enthalpySorptionFvPatchScalarField::enthalpyModelTypeNames
|
||||||
|
({
|
||||||
|
{ enthalpyModelType::estimated, "estimated" },
|
||||||
|
{ enthalpyModelType::calculated, "calculated" }
|
||||||
|
});
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::enthalpySorptionFvPatchScalarField::enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const fvPatch& p,
|
||||||
|
const DimensionedField<scalar, volMesh>& iF
|
||||||
|
)
|
||||||
|
:
|
||||||
|
zeroGradientFvPatchScalarField(p, iF),
|
||||||
|
enthalpyModel_(enthalpyModelType::estimated),
|
||||||
|
includeHs_(false),
|
||||||
|
enthalpyMassLoadPtr_(nullptr),
|
||||||
|
C_(0),
|
||||||
|
Hvap_(0),
|
||||||
|
speciesName_("none"),
|
||||||
|
pName_("p"),
|
||||||
|
TName_("T"),
|
||||||
|
dhdt_(p.size(), 0)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::enthalpySorptionFvPatchScalarField::enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const fvPatch& p,
|
||||||
|
const DimensionedField<scalar, volMesh>& iF,
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
:
|
||||||
|
zeroGradientFvPatchScalarField(p, iF, dict),
|
||||||
|
enthalpyModel_(enthalpyModelTypeNames.get("enthalpyModel", dict)),
|
||||||
|
includeHs_(dict.getOrDefault<bool>("includeHs", true)),
|
||||||
|
enthalpyMassLoadPtr_(nullptr),
|
||||||
|
C_(dict.getCheckOrDefault<scalar>("C", 0, scalarMinMax::ge(0))),
|
||||||
|
Hvap_(dict.getCheckOrDefault<scalar>("Hvap", 0, scalarMinMax::ge(0))),
|
||||||
|
speciesName_(dict.get<word>("species")),
|
||||||
|
pName_(dict.getOrDefault<word>("p", "p")),
|
||||||
|
TName_(dict.getOrDefault<word>("T", "T")),
|
||||||
|
dhdt_
|
||||||
|
(
|
||||||
|
dict.found("dhdt")
|
||||||
|
? scalarField("dhdt", dict, p.size())
|
||||||
|
: scalarField(p.size(), 0)
|
||||||
|
)
|
||||||
|
{
|
||||||
|
switch (enthalpyModel_)
|
||||||
|
{
|
||||||
|
case enthalpyModelType::calculated:
|
||||||
|
{
|
||||||
|
enthalpyMassLoadPtr_ =
|
||||||
|
Function1<scalar>::New("enthalpyTable", dict);
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
case enthalpyModelType::estimated:
|
||||||
|
{
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
if (dict.found("value"))
|
||||||
|
{
|
||||||
|
fvPatchScalarField::operator=
|
||||||
|
(
|
||||||
|
scalarField("value", dict, p.size())
|
||||||
|
);
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
fvPatchField<scalar>::operator=(Zero);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::enthalpySorptionFvPatchScalarField::enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const enthalpySorptionFvPatchScalarField& ptf,
|
||||||
|
const fvPatch& p,
|
||||||
|
const DimensionedField<scalar, volMesh>& iF,
|
||||||
|
const fvPatchFieldMapper& mapper
|
||||||
|
)
|
||||||
|
:
|
||||||
|
zeroGradientFvPatchScalarField(ptf, p, iF, mapper),
|
||||||
|
enthalpyModel_(ptf.enthalpyModel_),
|
||||||
|
includeHs_(ptf.includeHs_),
|
||||||
|
enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
|
||||||
|
C_(ptf.C_),
|
||||||
|
Hvap_(ptf.Hvap_),
|
||||||
|
speciesName_(ptf.speciesName_),
|
||||||
|
pName_(ptf.pName_),
|
||||||
|
TName_(ptf.TName_),
|
||||||
|
dhdt_(ptf.dhdt_, mapper)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::enthalpySorptionFvPatchScalarField::enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const enthalpySorptionFvPatchScalarField& ptf
|
||||||
|
)
|
||||||
|
:
|
||||||
|
zeroGradientFvPatchScalarField(ptf),
|
||||||
|
enthalpyModel_(ptf.enthalpyModel_),
|
||||||
|
includeHs_(ptf.includeHs_),
|
||||||
|
enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
|
||||||
|
C_(ptf.C_),
|
||||||
|
Hvap_(ptf.Hvap_),
|
||||||
|
speciesName_(ptf.speciesName_),
|
||||||
|
pName_(ptf.pName_),
|
||||||
|
TName_(ptf.TName_),
|
||||||
|
dhdt_(ptf.dhdt_)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::enthalpySorptionFvPatchScalarField::enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const enthalpySorptionFvPatchScalarField& ptf,
|
||||||
|
const DimensionedField<scalar, volMesh>& iF
|
||||||
|
)
|
||||||
|
:
|
||||||
|
zeroGradientFvPatchScalarField(ptf, iF),
|
||||||
|
enthalpyModel_(ptf.enthalpyModel_),
|
||||||
|
includeHs_(ptf.includeHs_),
|
||||||
|
enthalpyMassLoadPtr_(ptf.enthalpyMassLoadPtr_.clone()),
|
||||||
|
C_(ptf.C_),
|
||||||
|
Hvap_(ptf.Hvap_),
|
||||||
|
speciesName_(ptf.speciesName_),
|
||||||
|
pName_(ptf.pName_),
|
||||||
|
TName_(ptf.TName_),
|
||||||
|
dhdt_(ptf.dhdt_)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
void Foam::enthalpySorptionFvPatchScalarField::autoMap
|
||||||
|
(
|
||||||
|
const fvPatchFieldMapper& m
|
||||||
|
)
|
||||||
|
{
|
||||||
|
zeroGradientFvPatchScalarField::autoMap(m);
|
||||||
|
|
||||||
|
dhdt_.autoMap(m);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
void Foam::enthalpySorptionFvPatchScalarField::rmap
|
||||||
|
(
|
||||||
|
const fvPatchScalarField& ptf,
|
||||||
|
const labelList& addr
|
||||||
|
)
|
||||||
|
{
|
||||||
|
zeroGradientFvPatchScalarField::rmap(ptf, addr);
|
||||||
|
|
||||||
|
const auto& tiptf = refCast<const enthalpySorptionFvPatchScalarField>(ptf);
|
||||||
|
|
||||||
|
dhdt_.rmap(tiptf.dhdt_, addr);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::tmp<Foam::scalarField> Foam::enthalpySorptionFvPatchScalarField::
|
||||||
|
patchSource() const
|
||||||
|
{
|
||||||
|
const auto& Yp =
|
||||||
|
refCast<const speciesSorptionFvPatchScalarField>
|
||||||
|
(
|
||||||
|
patch().lookupPatchField<volScalarField, scalar>
|
||||||
|
(
|
||||||
|
speciesName_
|
||||||
|
)
|
||||||
|
);
|
||||||
|
|
||||||
|
// mass rate [kg/sec/m3]
|
||||||
|
tmp<scalarField> tmassb = Yp.patchSource();
|
||||||
|
const scalarField& massb = tmassb();
|
||||||
|
|
||||||
|
// The moles absorbed by the solid
|
||||||
|
// dhdt[J/kg] * kg/sec/m3 = [J/m3/s]
|
||||||
|
scalarField dhdt(dhdt_*massb);
|
||||||
|
|
||||||
|
if (includeHs_)
|
||||||
|
{
|
||||||
|
const fvPatchField<scalar>& pp =
|
||||||
|
patch().lookupPatchField<volScalarField, scalar>(pName_);
|
||||||
|
|
||||||
|
const fvPatchField<scalar>& Tp =
|
||||||
|
patch().lookupPatchField<volScalarField, scalar>(TName_);
|
||||||
|
|
||||||
|
const auto& thermo = db().lookupObject<rhoReactionThermo>
|
||||||
|
(
|
||||||
|
basicThermo::dictName
|
||||||
|
);
|
||||||
|
|
||||||
|
const basicSpecieMixture& composition = thermo.composition();
|
||||||
|
|
||||||
|
const label speicesId =
|
||||||
|
thermo.composition().species()[speciesName_];
|
||||||
|
|
||||||
|
scalarField hsp(this->patch().size(), 0);
|
||||||
|
|
||||||
|
forAll(pp, facei)
|
||||||
|
{
|
||||||
|
hsp[facei] = composition.Hs(speicesId, pp[facei], Tp[facei]);
|
||||||
|
}
|
||||||
|
|
||||||
|
dhdt += hsp*massb;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (debug)
|
||||||
|
{
|
||||||
|
Info<< " Patch enthalpy rate min/max [J/m3/sec]: "
|
||||||
|
<< gMin(dhdt) << " - " << gMax(dhdt) << endl;
|
||||||
|
}
|
||||||
|
|
||||||
|
return tmp<scalarField>::New(dhdt);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
void Foam::enthalpySorptionFvPatchScalarField::updateCoeffs()
|
||||||
|
{
|
||||||
|
if (updated())
|
||||||
|
{
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
|
||||||
|
const auto& Yp =
|
||||||
|
refCast<const speciesSorptionFvPatchScalarField>
|
||||||
|
(
|
||||||
|
patch().lookupPatchField<volScalarField, scalar>
|
||||||
|
(
|
||||||
|
speciesName_
|
||||||
|
)
|
||||||
|
);
|
||||||
|
|
||||||
|
switch (enthalpyModel_)
|
||||||
|
{
|
||||||
|
case enthalpyModelType::estimated:
|
||||||
|
{
|
||||||
|
dhdt_ = -C_*Hvap_;
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
case enthalpyModelType::calculated:
|
||||||
|
{
|
||||||
|
// mass [mol/kg]
|
||||||
|
tmp<scalarField> tmassb = Yp.mass();
|
||||||
|
const scalarField& massb = tmassb.ref();
|
||||||
|
|
||||||
|
forAll(massb, facei)
|
||||||
|
{
|
||||||
|
const scalar mfacei = massb[facei];
|
||||||
|
|
||||||
|
dhdt_[facei] = enthalpyMassLoadPtr_->value(mfacei);
|
||||||
|
}
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
default:
|
||||||
|
break;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (debug)
|
||||||
|
{
|
||||||
|
Info<< " Enthalpy change min/max [J/kg]: "
|
||||||
|
<< gMin(dhdt_) << " - " << gMax(dhdt_) << endl;
|
||||||
|
}
|
||||||
|
|
||||||
|
zeroGradientFvPatchScalarField::updateCoeffs();
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
void Foam::enthalpySorptionFvPatchScalarField::write(Ostream& os) const
|
||||||
|
{
|
||||||
|
fvPatchScalarField::write(os);
|
||||||
|
|
||||||
|
os.writeEntry("enthalpyModel", enthalpyModelTypeNames[enthalpyModel_]);
|
||||||
|
|
||||||
|
if (enthalpyMassLoadPtr_)
|
||||||
|
{
|
||||||
|
enthalpyMassLoadPtr_->writeData(os);
|
||||||
|
}
|
||||||
|
|
||||||
|
os.writeEntry("species", speciesName_);
|
||||||
|
|
||||||
|
os.writeEntryIfDifferent<bool>("includeHs", true, includeHs_);
|
||||||
|
os.writeEntryIfDifferent<scalar>("C", scalar(0), C_);
|
||||||
|
os.writeEntryIfDifferent<scalar>("Hvap", scalar(0), Hvap_);
|
||||||
|
os.writeEntryIfDifferent<word>("p", "p", pName_);
|
||||||
|
os.writeEntryIfDifferent<word>("T", "T", TName_);
|
||||||
|
|
||||||
|
dhdt_.writeEntry("dhdt", os);
|
||||||
|
|
||||||
|
writeEntry("value", os);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
makePatchTypeField
|
||||||
|
(
|
||||||
|
fvPatchScalarField,
|
||||||
|
enthalpySorptionFvPatchScalarField
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,284 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 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::enthalpySorptionFvPatchScalarField
|
||||||
|
|
||||||
|
Group
|
||||||
|
grpGenericBoundaryConditions
|
||||||
|
|
||||||
|
Description
|
||||||
|
This is a temperature boundary condition which works
|
||||||
|
in conjunction with the \c speciesSorption condition for species.
|
||||||
|
|
||||||
|
This boundary condition substracts or adds enthalpy associated with the
|
||||||
|
adsorption provided by the \c speciesSorption condition.
|
||||||
|
|
||||||
|
It can handle two enthalpy models:
|
||||||
|
|
||||||
|
1) Estimate
|
||||||
|
2) Calculated
|
||||||
|
|
||||||
|
On top of this, the enthalpy associated with the sensible enthalpy
|
||||||
|
corresponding with the species transport can be added using \c includeHs.
|
||||||
|
|
||||||
|
Usage
|
||||||
|
Example of the boundary condition specification:
|
||||||
|
\verbatim
|
||||||
|
<patchName>
|
||||||
|
{
|
||||||
|
// Mandatory entries
|
||||||
|
type enthalpySorption;
|
||||||
|
enthalpyModel <word>;
|
||||||
|
species <word>;
|
||||||
|
|
||||||
|
// Conditional mandatory entries
|
||||||
|
|
||||||
|
// when enthalpyModel == calculated
|
||||||
|
enthalpyTable <Function1<scalar>>
|
||||||
|
|
||||||
|
// enthalpyTable
|
||||||
|
// {
|
||||||
|
// type table;
|
||||||
|
// values ((0 0)(1 50));
|
||||||
|
// }
|
||||||
|
|
||||||
|
// Optional entries
|
||||||
|
includeHs <bool>;
|
||||||
|
C <scalar>;
|
||||||
|
Hvap <scalar>;
|
||||||
|
dhdt <scalarField>;
|
||||||
|
p <word>;
|
||||||
|
T <word>;
|
||||||
|
}
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
where the entries mean:
|
||||||
|
\table
|
||||||
|
Property | Description | Type | Reqd | Deflt
|
||||||
|
type | Type name: enthalpySorption | word | yes | -
|
||||||
|
enthalpyModel | Adsorption enthalpy model | word | yes | -
|
||||||
|
species | Name of associated species | word | yes | -
|
||||||
|
enthalpyTable | Calculated enthalpy model table <!--
|
||||||
|
--> | Function1\<scalar\> | no | -
|
||||||
|
includeHs | Include sensible enthalpy | bool | no | true
|
||||||
|
C | Estimate enthalpy model constant | scalar | no | 0
|
||||||
|
Hvap | Evaporation heat for species | scalar | no | 0
|
||||||
|
p | Name of operand pressure field | word | no | p
|
||||||
|
T | Name of operand temperature field | word | no | T
|
||||||
|
dhdt | Enthalpy change on cells next to patch [J/kg] <!--
|
||||||
|
--> | scalarField | no | Zero
|
||||||
|
\endtable
|
||||||
|
|
||||||
|
Options for the \c enthalpyModel entry:
|
||||||
|
\verbatim
|
||||||
|
estimated | Enthalpy is estimated
|
||||||
|
calculated | Enthalpy is calculated based on enthalpyTable
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
The inherited entries are elaborated in:
|
||||||
|
- \link zeroGradientFvPatchFields.H \endlink
|
||||||
|
- \link Function1.H \endlink
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
enthalpySorptionFvPatchScalarField.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef Foam_enthalpySorptionFvPatchScalarField_H
|
||||||
|
#define Foam_enthalpySorptionFvPatchScalarField_H
|
||||||
|
|
||||||
|
#include "boundarySourcePatch.H"
|
||||||
|
#include "zeroGradientFvPatchFields.H"
|
||||||
|
#include "Function1.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class enthalpySorptionFvPatchScalarField Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
class enthalpySorptionFvPatchScalarField
|
||||||
|
:
|
||||||
|
public zeroGradientFvPatchField<scalar>,
|
||||||
|
public boundarySourcePatch
|
||||||
|
{
|
||||||
|
// Private Enumeration
|
||||||
|
|
||||||
|
//- Options for the enthalpy model
|
||||||
|
enum enthalpyModelType : char
|
||||||
|
{
|
||||||
|
estimated,
|
||||||
|
calculated
|
||||||
|
};
|
||||||
|
|
||||||
|
//- Names for enthalpyModelType
|
||||||
|
static const Enum<enthalpyModelType> enthalpyModelTypeNames;
|
||||||
|
|
||||||
|
|
||||||
|
// Private Data
|
||||||
|
|
||||||
|
//- Enthalpy model
|
||||||
|
enum enthalpyModelType enthalpyModel_;
|
||||||
|
|
||||||
|
//- Include sensible enthalpy of the species
|
||||||
|
bool includeHs_;
|
||||||
|
|
||||||
|
//- Load enthalpy table for calculated model
|
||||||
|
autoPtr<Function1<scalar>> enthalpyMassLoadPtr_;
|
||||||
|
|
||||||
|
//- Estimated enthalpy model constant
|
||||||
|
scalar C_;
|
||||||
|
|
||||||
|
//- Heat of evaporation of species
|
||||||
|
scalar Hvap_;
|
||||||
|
|
||||||
|
//- Name of operand species field
|
||||||
|
word speciesName_;
|
||||||
|
|
||||||
|
//- Name of operand pressure field
|
||||||
|
word pName_;
|
||||||
|
|
||||||
|
//- Name of operand temperature field
|
||||||
|
word TName_;
|
||||||
|
|
||||||
|
//- Enthalpy change on cells next to patch [J/kg]
|
||||||
|
scalarField dhdt_;
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("enthalpySorption");
|
||||||
|
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from patch and internal field
|
||||||
|
enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const fvPatch&,
|
||||||
|
const DimensionedField<scalar, volMesh>&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct from patch, internal field and dictionary
|
||||||
|
enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const fvPatch&,
|
||||||
|
const DimensionedField<scalar, volMesh>&,
|
||||||
|
const dictionary&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct by mapping given
|
||||||
|
//- enthalpySorptionFvPatchScalarField onto a new patch
|
||||||
|
enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const enthalpySorptionFvPatchScalarField&,
|
||||||
|
const fvPatch&,
|
||||||
|
const DimensionedField<scalar, volMesh>&,
|
||||||
|
const fvPatchFieldMapper&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct as copy
|
||||||
|
enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const enthalpySorptionFvPatchScalarField&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct and return a clone
|
||||||
|
virtual tmp<fvPatchScalarField> clone() const
|
||||||
|
{
|
||||||
|
return tmp<fvPatchScalarField>
|
||||||
|
(
|
||||||
|
new enthalpySorptionFvPatchScalarField(*this)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
//- Construct as copy setting internal field reference
|
||||||
|
enthalpySorptionFvPatchScalarField
|
||||||
|
(
|
||||||
|
const enthalpySorptionFvPatchScalarField&,
|
||||||
|
const DimensionedField<scalar, volMesh>&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Construct and return a clone setting internal field reference
|
||||||
|
virtual tmp<fvPatchScalarField> clone
|
||||||
|
(
|
||||||
|
const DimensionedField<scalar, volMesh>& iF
|
||||||
|
) const
|
||||||
|
{
|
||||||
|
return tmp<fvPatchScalarField>
|
||||||
|
(
|
||||||
|
new enthalpySorptionFvPatchScalarField(*this, iF)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// Member Functions
|
||||||
|
|
||||||
|
// Mapping
|
||||||
|
|
||||||
|
//- Map (and resize as needed) from self given a mapping object
|
||||||
|
virtual void autoMap
|
||||||
|
(
|
||||||
|
const fvPatchFieldMapper&
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Reverse map the given fvPatchField onto this fvPatchField
|
||||||
|
virtual void rmap
|
||||||
|
(
|
||||||
|
const fvPatchScalarField&,
|
||||||
|
const labelList&
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
|
// Evaluation
|
||||||
|
|
||||||
|
//- Source of cells next to the patch
|
||||||
|
virtual tmp<scalarField> patchSource() const;
|
||||||
|
|
||||||
|
//- Update the coefficients associated with the patch field
|
||||||
|
virtual void updateCoeffs();
|
||||||
|
|
||||||
|
|
||||||
|
// I-O
|
||||||
|
|
||||||
|
//- Write
|
||||||
|
virtual void write(Ostream&) const;
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
Reference in New Issue
Block a user