semiPermeableBaffle: Added mole-fraction and partial-pressure input options

The semiPermeableBaffleMassFraction boundary condition can now calculate
the mass flux as proportional to the difference in mole fraction or
partial pressure. A mass fraction difference driven transfer is also
still possible. An additional keyword, "input" has been added which is
used to select the variable used to calculate the transfer. An example
specification is as follows:

    baffle
    {
        type            semiPermeableBaffleMassFraction;
        samplePatch     membranePipe;
        c               0.1;
        input           massFraction;
        value           uniform 0;
    }

In order to facilitate this, a "W" method to get the molar mass on a
patch has been added to the thermodynamics. To avoid name-clashes,
methods that generate per-species molar masses have been renamed "Wi".

This work was supported by Georg Skillas, at Evonik
This commit is contained in:
Will Bainbridge
2018-07-06 13:54:31 +01:00
parent e0cf9706ca
commit 20653ee01e
32 changed files with 272 additions and 86 deletions

View File

@ -313,6 +313,17 @@ Foam::tmp<Foam::volScalarField> Foam::twoPhaseMixtureThermo::W() const
}
Foam::tmp<Foam::scalarField> Foam::twoPhaseMixtureThermo::W
(
const label patchi
) const
{
return
alpha1().boundaryField()[patchi]*thermo1_->W(patchi)
+ alpha2().boundaryField()[patchi]*thermo1_->W(patchi);
}
Foam::tmp<Foam::volScalarField> Foam::twoPhaseMixtureThermo::nu() const
{
return mu()/(alpha1()*thermo1_->rho() + alpha2()*thermo2_->rho());

View File

@ -248,6 +248,9 @@ public:
//- Molecular weight [kg/kmol]
virtual tmp<volScalarField> W() const;
//- Molecular weight for patch [kg/kmol]
virtual tmp<scalarField> W(const label patchi) const;
// Fields derived from transport state variables

View File

@ -551,6 +551,28 @@ Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::W() const
}
Foam::tmp<Foam::scalarField> Foam::multiphaseMixtureThermo::W
(
const label patchi
) const
{
PtrDictionary<phaseModel>::const_iterator phasei = phases_.begin();
tmp<scalarField> tW
(
phasei().boundaryField()[patchi]*phasei().thermo().W(patchi)
);
for (++phasei; phasei != phases_.end(); ++phasei)
{
tW.ref() +=
phasei().boundaryField()[patchi]*phasei().thermo().W(patchi);
}
return tW;
}
Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::nu() const
{
return mu()/rho();

View File

@ -379,6 +379,9 @@ public:
//- Molecular weight [kg/kmol]
virtual tmp<volScalarField> W() const;
//- Molecular weight for patch [kg/kmol]
virtual tmp<scalarField> W(const label patchi) const;
// Fields derived from transport state variables

View File

@ -165,7 +165,7 @@ update
(
"W",
dimMass/dimMoles,
this->thermo_.composition().W(species1Index_)
this->thermo_.composition().Wi(species1Index_)
)
);
@ -177,7 +177,7 @@ update
(
"W",
dimMass/dimMoles,
this->thermo_.composition().W(species2Index_)
this->thermo_.composition().Wi(species2Index_)
)
);

View File

@ -36,7 +36,7 @@ wRatioByP() const
(
"W",
dimMass/dimMoles,
this->thermo_.composition().W(saturatedIndex_)
this->thermo_.composition().Wi(saturatedIndex_)
);
return Wi/this->thermo_.W()/this->thermo_.p();

View File

@ -52,8 +52,8 @@ Foam::COxidationDiffusionLimitedRate<CloudType>::COxidationDiffusionLimitedRate
CsLocalId_ = owner.composition().localId(idSolid, "C");
// Set local copies of thermo properties
WO2_ = owner.thermo().carrier().W(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().W(CO2GlobalId_);
WO2_ = owner.thermo().carrier().Wi(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().Wi(CO2GlobalId_);
WC_ = WCO2 - WO2_;
HcCO2_ = owner.thermo().carrier().Hc(CO2GlobalId_);

View File

@ -52,8 +52,8 @@ Foam::COxidationHurtMitchell<CloudType>::COxidationHurtMitchell
ashLocalId_ = owner.composition().localId(idSolid, "ash", true);
// Set local copies of thermo properties
WO2_ = owner.thermo().carrier().W(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().W(CO2GlobalId_);
WO2_ = owner.thermo().carrier().Wi(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().Wi(CO2GlobalId_);
WC_ = WCO2 - WO2_;
HcCO2_ = owner.thermo().carrier().Hc(CO2GlobalId_);

View File

@ -58,8 +58,8 @@ Foam::COxidationIntrinsicRate<CloudType>::COxidationIntrinsicRate
CsLocalId_ = owner.composition().localId(idSolid, "C");
// Set local copies of thermo properties
WO2_ = owner.thermo().carrier().W(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().W(CO2GlobalId_);
WO2_ = owner.thermo().carrier().Wi(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().Wi(CO2GlobalId_);
WC_ = WCO2 - WO2_;
HcCO2_ = owner.thermo().carrier().Hc(CO2GlobalId_);

View File

@ -53,8 +53,8 @@ COxidationKineticDiffusionLimitedRate
CsLocalId_ = owner.composition().localId(idSolid, "C");
// Set local copies of thermo properties
WO2_ = owner.thermo().carrier().W(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().W(CO2GlobalId_);
WO2_ = owner.thermo().carrier().Wi(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().Wi(CO2GlobalId_);
WC_ = WCO2 - WO2_;
HcCO2_ = owner.thermo().carrier().Hc(CO2GlobalId_);

View File

@ -65,8 +65,8 @@ Foam::COxidationMurphyShaddix<CloudType>::COxidationMurphyShaddix
CsLocalId_ = owner.composition().localId(idSolid, "C");
// Set local copies of thermo properties
WO2_ = owner.thermo().carrier().W(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().W(CO2GlobalId_);
WO2_ = owner.thermo().carrier().Wi(O2GlobalId_);
const scalar WCO2 = owner.thermo().carrier().Wi(CO2GlobalId_);
WC_ = WCO2 - WO2_;
HcCO2_ = owner.thermo().carrier().Hc(CO2GlobalId_);

View File

@ -581,7 +581,7 @@ void Foam::ReactingMultiphaseParcel<ParcelType>::calcDevolatilisation
{
const label id = composition.localToCarrierId(GAS, i);
const scalar Cp = composition.carrier().Cp(id, td.pc(), Ts);
const scalar W = composition.carrier().W(id);
const scalar W = composition.carrier().Wi(id);
const scalar Ni = dMassDV[i]/(this->areaS(d)*dt*W);
// Dab calc'd using API vapour mass diffusivity function

View File

@ -127,7 +127,7 @@ void Foam::ReactingParcel<ParcelType>::calcPhaseChange
const label cid = composition.localToCarrierId(idPhase, i);
const scalar Cp = composition.carrier().Cp(cid, td.pc(), Tsdash);
const scalar W = composition.carrier().W(cid);
const scalar W = composition.carrier().Wi(cid);
const scalar Ni = dMassPC[i]/(this->areaS(d)*dt*W);
const scalar Dab =
@ -311,7 +311,7 @@ void Foam::ReactingParcel<ParcelType>::correctSurfaceValues
forAll(Xinf, i)
{
Xinf[i] = thermo.carrier().Y(i)[this->cell()]/thermo.carrier().W(i);
Xinf[i] = thermo.carrier().Y(i)[this->cell()]/thermo.carrier().Wi(i);
}
Xinf /= sum(Xinf);
@ -333,7 +333,7 @@ void Foam::ReactingParcel<ParcelType>::correctSurfaceValues
const scalar Csi = Cs[i] + Xsff*Xinf[i]*CsTot;
Xs[i] = (2.0*Csi + Xinf[i]*CsTot)/3.0;
Ys[i] = Xs[i]*thermo.carrier().W(i);
Ys[i] = Xs[i]*thermo.carrier().Wi(i);
}
Xs /= sum(Xs);
Ys /= sum(Ys);
@ -348,7 +348,7 @@ void Foam::ReactingParcel<ParcelType>::correctSurfaceValues
forAll(Ys, i)
{
const scalar W = thermo.carrier().W(i);
const scalar W = thermo.carrier().Wi(i);
const scalar sqrtW = sqrt(W);
const scalar cbrtW = cbrt(W);

View File

@ -245,7 +245,7 @@ Foam::scalarField Foam::CompositionModel<CloudType>::X
forAll(Y, i)
{
label cid = props.carrierIds()[i];
X[i] = Y[i]/thermo_.carrier().W(cid);
X[i] = Y[i]/thermo_.carrier().Wi(cid);
WInv += X[i];
}
break;

View File

@ -43,7 +43,7 @@ Foam::tmp<Foam::scalarField> Foam::LiquidEvaporation<CloudType>::calcXc
{
Xc[i] =
this->owner().thermo().carrier().Y()[i][celli]
/this->owner().thermo().carrier().W(i);
/this->owner().thermo().carrier().Wi(i);
}
return Xc/sum(Xc);

View File

@ -43,7 +43,7 @@ Foam::tmp<Foam::scalarField> Foam::LiquidEvaporationBoil<CloudType>::calcXc
{
Xc[i] =
this->owner().thermo().carrier().Y()[i][celli]
/this->owner().thermo().carrier().W(i);
/this->owner().thermo().carrier().Wi(i);
}
return Xc/sum(Xc);

View File

@ -124,7 +124,7 @@ void standardPhaseChange::correctModel
);
// Molecular weight of vapour [kg/kmol]
const scalar Wvap = thermo.carrier().W(vapId);
const scalar Wvap = thermo.carrier().Wi(vapId);
// Molecular weight of liquid [kg/kmol]
const scalar Wliq = filmThermo.W();

View File

@ -179,7 +179,7 @@ void waxSolventEvaporation::correctModel
);
// Molecular weight of vapour [kg/kmol]
const scalar Wvap = thermo.carrier().W(vapId);
const scalar Wvap = thermo.carrier().Wi(vapId);
const scalar Wwax = Wwax_.value();
const scalar Wsolvent = Wsolvent_.value();

View File

@ -29,6 +29,62 @@ License
#include "volFields.H"
#include "surfaceFields.H"
#include "turbulenceModel.H"
#include "psiReactionThermo.H"
#include "rhoReactionThermo.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
template<>
const char* NamedEnum
<
semiPermeableBaffleMassFractionFvPatchScalarField::input,
4
>::names[] =
{
"none",
"massFraction",
"moleFraction",
"partialPressure",
};
}
const Foam::NamedEnum
<
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::input,
4
> Foam::semiPermeableBaffleMassFractionFvPatchScalarField::inputNames_;
// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
const Foam::basicSpecieMixture&
Foam::semiPermeableBaffleMassFractionFvPatchScalarField::composition
(
const objectRegistry& db
)
{
const word& name = basicThermo::dictName;
if (db.foundObject<psiReactionThermo>(name))
{
return db.lookupObject<psiReactionThermo>(name).composition();
}
else if (db.foundObject<rhoReactionThermo>(name))
{
return db.lookupObject<rhoReactionThermo>(name).composition();
}
else
{
FatalErrorInFunction
<< "Could not find a multi-component thermodynamic model."
<< exit(FatalError);
return NullObjectRef<basicSpecieMixture>();
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
@ -42,7 +98,9 @@ semiPermeableBaffleMassFractionFvPatchScalarField
mappedPatchBase(p.patch()),
mixedFvPatchScalarField(p, iF),
c_(0),
phiName_("phi")
input_(none),
phiName_("phi"),
pName_("p")
{
refValue() = Zero;
refGrad() = Zero;
@ -61,7 +119,14 @@ semiPermeableBaffleMassFractionFvPatchScalarField
mappedPatchBase(p.patch(), NEARESTPATCHFACE, dict),
mixedFvPatchScalarField(p, iF),
c_(dict.lookupOrDefault<scalar>("c", scalar(0))),
phiName_(dict.lookupOrDefault<word>("phi", "phi"))
input_
(
c_ == scalar(0)
? none
: inputNames_.read(dict.lookup("input"))
),
phiName_(dict.lookupOrDefault<word>("phi", "phi")),
pName_(dict.lookupOrDefault<word>("p", "p"))
{
fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
@ -83,7 +148,9 @@ semiPermeableBaffleMassFractionFvPatchScalarField
mappedPatchBase(p.patch(), ptf),
mixedFvPatchScalarField(ptf, p, iF, mapper),
c_(ptf.c_),
phiName_(ptf.phiName_)
input_(ptf.input_),
phiName_(ptf.phiName_),
pName_(ptf.pName_)
{}
@ -96,7 +163,9 @@ semiPermeableBaffleMassFractionFvPatchScalarField
mappedPatchBase(ptf.patch().patch(), ptf),
mixedFvPatchScalarField(ptf),
c_(ptf.c_),
phiName_(ptf.phiName_)
input_(ptf.input_),
phiName_(ptf.phiName_),
pName_(ptf.pName_)
{}
@ -110,7 +179,9 @@ semiPermeableBaffleMassFractionFvPatchScalarField
mappedPatchBase(ptf.patch().patch(), ptf),
mixedFvPatchScalarField(ptf, iF),
c_(ptf.c_),
phiName_(ptf.phiName_)
input_(ptf.input_),
phiName_(ptf.phiName_),
pName_(ptf.pName_)
{}
@ -126,15 +197,61 @@ Foam::semiPermeableBaffleMassFractionFvPatchScalarField::phiY() const
const word& YName = internalField().name();
// Initialise the input variables to the mass fractions
scalarField psic(patchInternalField());
const label nbrPatchi = samplePolyPatch().index();
const fvPatch& nbrPatch = patch().boundaryMesh()[nbrPatchi];
const fvPatchScalarField& nbrYp =
nbrPatch.lookupPatchField<volScalarField, scalar>(YName);
scalarField nbrYc(nbrYp.patchInternalField());
mappedPatchBase::map().distribute(nbrYc);
scalarField nbrPsic(nbrYp.patchInternalField());
mappedPatchBase::map().distribute(nbrPsic);
return c_*patch().magSf()*(patchInternalField() - nbrYc);
switch (input_)
{
case none:
FatalErrorInFunction
<< "A none input cannot be used with a non-zero transfer "
<< "coefficient" << exit(FatalError);
case massFraction:
// Do nothing
break;
case partialPressure:
// Multiply by pressure
{
psic *=
patch().lookupPatchField<volScalarField, scalar>(pName_);
fvPatchScalarField nbrP
(
nbrPatch.lookupPatchField<volScalarField, scalar>(pName_)
);
mappedPatchBase::map().distribute(nbrP);
nbrPsic *= nbrP;
}
// Falls through ...
case moleFraction:
// Convert to mole fraction
{
const basicSpecieMixture& mixture = composition(db());
const scalar Wi(mixture.Wi(mixture.species()[YName]));
const basicThermo& thermo =
db().lookupObject<basicThermo>(basicThermo::dictName);
psic *= thermo.W(patch().index())/Wi;
scalarField nbrW(thermo.W(nbrPatch.index()));
mappedPatchBase::map().distribute(nbrW);
nbrPsic *= nbrW/Wi;
}
break;
}
return c_*patch().magSf()*(psic - nbrPsic);
}
@ -171,7 +288,10 @@ void Foam::semiPermeableBaffleMassFractionFvPatchScalarField::write
fvPatchScalarField::write(os);
mappedPatchBase::write(os);
writeEntryIfDifferent<scalar>(os, "c", scalar(0), c_);
os.writeKeyword("input") << inputNames_[input_]
<< token::END_STATEMENT << nl;
writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
writeEntryIfDifferent<word>(os, "p", "p", pName_);
writeEntry("value", os);
}

View File

@ -33,17 +33,17 @@ Description
semiPermeableBaffleVelocityFvPatchVectorField.
The mass flux of a species is calculated as a coefficient multiplied by the
difference in mass fraction across the baffle.
difference in an input variable across the baffle.
\f[
\phi_{Yi} = c A (Y_i - Y_{i,n})
\phi_{Yi} = c A (\psi_i - \psi_{i,n})
\f]
where
\vartable
\phi_{Yi} | flux of the permeable species [kg/s]
c | transfer coefficient [kg/m2/s]
c | transfer coefficient [kg/m2/s/<input-dimensions>]
A | patch face area [m2]
Y_i | mass fraction on the patch []
Y_{i,n} | mass fraction on the neighbour patch []
\psi_i | input variable on the patch []
\psi_{i,n} | input variable on the neighbour patch []
\endvartable
A species that the baffle is permeable to will, therefore, have a
@ -59,7 +59,10 @@ Usage
\table
Property | Description | Req'd? | Default
c | Transfer coefficient | no | 0
input | Input variable used to drive the transfer; massFraction, \\
moleFraction or partialPressure | if c is non-zero | none
phi | Name of the flux field | no | phi
p | Name of the pressure field | no | p
\endtable
See also
@ -81,6 +84,8 @@ SourceFiles
namespace Foam
{
class basicSpecieMixture;
/*---------------------------------------------------------------------------*\
Class semiPermeableBaffleMassFractionFvPatchScalarField Declaration
\*---------------------------------------------------------------------------*/
@ -90,14 +95,37 @@ class semiPermeableBaffleMassFractionFvPatchScalarField
public mappedPatchBase,
public mixedFvPatchScalarField
{
public:
//- Enumeration for the input variable driving the transfer
enum input
{
none,
massFraction,
moleFraction,
partialPressure,
};
//- Input variable type names
static const NamedEnum<input, 4> inputNames_;
private:
// Private data
//- Transfer coefficient
const scalar c_;
//- Input variable driving the transfer
const input input_;
//- Name of the flux field
const word phiName_;
//- Name of the pressure field. Only needed if mode is pressure.
const word pName_;
public:
@ -105,6 +133,12 @@ public:
TypeName("semiPermeableBaffleMassFraction");
// Static member functions
//- Access the composition for the given database
static const basicSpecieMixture& composition(const objectRegistry& db);
// Constructors
//- Construct from patch and internal field

View File

@ -32,32 +32,6 @@ License
#include "psiReactionThermo.H"
#include "rhoReactionThermo.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
const Foam::basicSpecieMixture&
Foam::semiPermeableBaffleVelocityFvPatchVectorField::composition() const
{
const word& name = basicThermo::dictName;
if (db().foundObject<psiReactionThermo>(name))
{
return db().lookupObject<psiReactionThermo>(name).composition();
}
else if (db().foundObject<rhoReactionThermo>(name))
{
return db().lookupObject<rhoReactionThermo>(name).composition();
}
else
{
FatalErrorInFunction
<< "Could not find a multi-component thermodynamic model."
<< exit(FatalError);
return NullObjectRef<basicSpecieMixture>();
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::semiPermeableBaffleVelocityFvPatchVectorField::
@ -138,7 +112,7 @@ void Foam::semiPermeableBaffleVelocityFvPatchVectorField::updateCoeffs()
const scalarField& rhop =
patch().lookupPatchField<volScalarField, scalar>(rhoName_);
const PtrList<volScalarField>& Y = composition().Y();
const PtrList<volScalarField>& Y = YBCType::composition(db()).Y();
scalarField phip(patch().size(), Zero);
forAll(Y, i)

View File

@ -60,8 +60,6 @@ SourceFiles
namespace Foam
{
class basicSpecieMixture;
/*---------------------------------------------------------------------------*\
Class semiPermeableBaffleVelocityFvPatchVectorField Declaration
\*---------------------------------------------------------------------------*/
@ -76,12 +74,6 @@ class semiPermeableBaffleVelocityFvPatchVectorField
const word rhoName_;
// Private Member Functions
//- Return the composition
const basicSpecieMixture& composition() const;
public:
//- Runtime type information

View File

@ -404,6 +404,9 @@ public:
//- Molecular weight [kg/kmol]
virtual tmp<volScalarField> W() const = 0;
//- Molecular weight for patch [kg/kmol]
virtual tmp<scalarField> W(const label patchi) const = 0;
// Access to transport state variables

View File

@ -768,6 +768,25 @@ Foam::tmp<Foam::volScalarField> Foam::heThermo<BasicThermo, MixtureType>::W
}
template<class BasicThermo, class MixtureType>
Foam::tmp<Foam::scalarField> Foam::heThermo<BasicThermo, MixtureType>::W
(
const label patchi
) const
{
const fvMesh& mesh = this->T_.mesh();
tmp<scalarField> tW(new scalarField(mesh.boundaryMesh()[patchi].size()));
scalarField& W = tW.ref();
forAll(W, facei)
{
W[facei] = this->patchFaceMixture(patchi, facei).W();
}
return tW;
}
template<class BasicThermo, class MixtureType>
Foam::tmp<Foam::volScalarField>
Foam::heThermo<BasicThermo, MixtureType>::kappa() const

View File

@ -262,6 +262,9 @@ public:
//- Molecular weight [kg/kmol]
virtual tmp<volScalarField> W() const;
//- Molecular weight for patch [kg/kmol]
virtual tmp<scalarField> W(const label patchi) const;
// Fields derived from transport state variables

View File

@ -225,11 +225,11 @@ Foam::radiation::greyMeanAbsorptionEmission::aCont(const label bandI) const
scalar invWt = 0.0;
forAll(mixture.Y(), s)
{
invWt += mixture.Y(s)[celli]/mixture.W(s);
invWt += mixture.Y(s)[celli]/mixture.Wi(s);
}
label index = mixture.species()[iter.key()];
scalar Xk = mixture.Y(index)[celli]/(mixture.W(index)*invWt);
scalar Xk = mixture.Y(index)[celli]/(mixture.Wi(index)*invWt);
Xipi = Xk*paToAtm(p[celli]);
}

View File

@ -239,13 +239,13 @@ Foam::radiation::wideBandAbsorptionEmission::aCont(const label bandi) const
scalar invWt = 0;
forAll(mixture.Y(), s)
{
invWt += mixture.Y(s)[celli]/mixture.W(s);
invWt += mixture.Y(s)[celli]/mixture.Wi(s);
}
const label index = mixture.species()[iter.key()];
const scalar Xk =
mixture.Y(index)[celli]/(mixture.W(index)*invWt);
mixture.Y(index)[celli]/(mixture.Wi(index)*invWt);
Xipi = Xk*paToAtm(p[celli]);
}

View File

@ -42,9 +42,9 @@ void Foam::moleFractions<ThermoType>::calculateMoleFractions()
{
const dimensionedScalar Wi
(
"W",
"Wi",
dimMass/dimMoles,
thermo.composition().W(i)
thermo.composition().Wi(i)
);
X_[i] = W*Y[i]/Wi;

View File

@ -49,7 +49,7 @@ Foam::SpecieMixture<MixtureType>::SpecieMixture
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class MixtureType>
Foam::scalar Foam::SpecieMixture<MixtureType>::W
Foam::scalar Foam::SpecieMixture<MixtureType>::Wi
(
const label speciei
) const

View File

@ -76,7 +76,7 @@ public:
// Per specie properties
//- Molecular weight of the given specie [kg/kmol]
virtual scalar W(const label speciei) const;
virtual scalar Wi(const label speciei) const;
// Per specie thermo properties

View File

@ -83,7 +83,7 @@ public:
// Per specie properties
//- Molecular weight of the given specie [kg/kmol]
virtual scalar W(const label speciei) const = 0;
virtual scalar Wi(const label speciei) const = 0;
// Per specie thermo properties

View File

@ -54,6 +54,7 @@ boundaryField
type semiPermeableBaffleMassFraction;
samplePatch membranePipe;
c 0.1;
input massFraction;
value uniform $:sleeve.CH4;
}
membranePipe
@ -61,6 +62,7 @@ boundaryField
type semiPermeableBaffleMassFraction;
samplePatch membraneSleeve;
c 0.1;
input massFraction;
value uniform $:pipe.CH4;
}
}