This is a mass-fraction boundary condition for a semi-permeable baffle. More...
Public Types | |
enum | input { none, massFraction, moleFraction, partialPressure } |
Enumeration for the input variable driving the transfer. More... | |
![]() | |
enum | sampleMode { NEARESTCELL, NEARESTPATCHFACE, NEARESTPATCHFACEAMI, NEARESTPATCHPOINT, NEARESTFACE, NEARESTONLYCELL } |
Mesh items to sample. More... | |
enum | offsetMode { UNIFORM, NONUNIFORM, NORMAL } |
How to project face centres. More... | |
typedef Tuple2< pointIndexHit, Tuple2< scalar, label > > | nearInfo |
Helper class for finding nearest. More... | |
Public Member Functions | |
TypeName ("semiPermeableBaffleMassFraction") | |
Runtime type information. More... | |
semiPermeableBaffleMassFractionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
Construct from patch and internal field. More... | |
semiPermeableBaffleMassFractionFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
semiPermeableBaffleMassFractionFvPatchScalarField (const semiPermeableBaffleMassFractionFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
Construct by mapping given fixedValueTypeFvPatchField. More... | |
semiPermeableBaffleMassFractionFvPatchScalarField (const semiPermeableBaffleMassFractionFvPatchScalarField &) | |
Copy constructor. More... | |
virtual tmp< fvPatchScalarField > | clone () const |
Construct and return a clone. More... | |
semiPermeableBaffleMassFractionFvPatchScalarField (const semiPermeableBaffleMassFractionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
Copy constructor setting internal field reference. More... | |
virtual tmp< fvPatchScalarField > | clone (const DimensionedField< scalar, volMesh > &iF) const |
Construct and return a clone setting internal field reference. More... | |
tmp< scalarField > | phiY () const |
Return the flux of this species through the baffle. More... | |
virtual void | updateCoeffs () |
Update the coefficients associated with the patch field. More... | |
virtual void | write (Ostream &) const |
Write. More... | |
![]() | |
TypeName ("mappedPatchBase") | |
Runtime type information. More... | |
mappedPatchBase (const polyPatch &) | |
Construct from patch. More... | |
mappedPatchBase (const polyPatch &pp, const word &sampleRegion, const sampleMode sampleMode, const word &samplePatch, const vectorField &offsets) | |
Construct with offsetMode=non-uniform. More... | |
mappedPatchBase (const polyPatch &pp, const word &sampleRegion, const sampleMode sampleMode, const word &samplePatch, const vector &offset) | |
Construct from offsetMode=uniform. More... | |
mappedPatchBase (const polyPatch &pp, const word &sampleRegion, const sampleMode sampleMode, const word &samplePatch, const scalar distance) | |
Construct from offsetMode=normal and distance. More... | |
mappedPatchBase (const polyPatch &, const dictionary &) | |
Construct from dictionary. More... | |
mappedPatchBase (const polyPatch &, const sampleMode, const dictionary &) | |
Construct from dictionary and (collocated) sample mode. More... | |
mappedPatchBase (const polyPatch &, const mappedPatchBase &) | |
Construct as copy, resetting patch. More... | |
mappedPatchBase (const polyPatch &, const mappedPatchBase &, const labelUList &mapAddressing) | |
Construct as copy, resetting patch, map original data. More... | |
virtual | ~mappedPatchBase () |
Destructor. More... | |
void | clearOut () |
const sampleMode & | mode () const |
What to sample. More... | |
const word & | sampleRegion () const |
Region to sample. More... | |
const word & | samplePatch () const |
Patch (only if NEARESTPATCHFACE) More... | |
const word & | coupleGroup () const |
PatchGroup (only if NEARESTPATCHFACE) More... | |
label | sampleSize () const |
Return size of mapped mesh/patch/boundary. More... | |
const vector & | offset () const |
Offset vector (from patch faces to destination mesh objects) More... | |
const vectorField & | offsets () const |
Offset vector (from patch faces to destination mesh objects) More... | |
bool | sameRegion () const |
Cached sampleRegion != mesh.name() More... | |
const mapDistribute & | map () const |
Return reference to the parallel distribution map. More... | |
const AMIInterpolation & | AMI (const bool forceUpdate=false) const |
Return reference to the AMI interpolator. More... | |
const autoPtr< Foam::searchableSurface > & | surfPtr () const |
Return a pointer to the AMI projection surface. More... | |
const polyMesh & | sampleMesh () const |
Get the region mesh. More... | |
const polyPatch & | samplePolyPatch () const |
Get the patch on the region. More... | |
tmp< pointField > | samplePoints () const |
Get the sample points. More... | |
template<class Type > | |
void | distribute (List< Type > &lst) const |
Wrapper around map/interpolate data distribution. More... | |
template<class Type , class CombineOp > | |
void | distribute (List< Type > &lst, const CombineOp &cop) const |
Wrapper around map/interpolate data distribution with operation. More... | |
template<class Type > | |
void | reverseDistribute (List< Type > &lst) const |
Wrapper around map/interpolate data distribution. More... | |
template<class Type , class CombineOp > | |
void | reverseDistribute (List< Type > &lst, const CombineOp &cop) const |
Wrapper around map/interpolate data distribution with operation. More... | |
Static Public Member Functions | |
static const basicSpecieMixture & | composition (const objectRegistry &db) |
Access the composition for the given database. More... | |
![]() | |
static pointIndexHit | facePoint (const polyMesh &, const label facei, const polyMesh::cellDecomposition) |
Get a point on the face given a face decomposition method: More... | |
Static Public Attributes | |
static const NamedEnum< input, 4 > | inputNames_ |
Input variable type names. More... | |
![]() | |
static const NamedEnum< sampleMode, 6 > | sampleModeNames_ |
static const NamedEnum< offsetMode, 3 > | offsetModeNames_ |
Additional Inherited Members | |
![]() | |
tmp< pointField > | facePoints (const polyPatch &) const |
Get the points from face-centre-decomposition face centres. More... | |
void | collectSamples (const pointField &facePoints, pointField &, labelList &patchFaceProcs, labelList &patchFaces, pointField &patchFc) const |
Collect single list of samples and originating processor+face. More... | |
void | findSamples (const sampleMode mode, const pointField &, labelList &sampleProcs, labelList &sampleIndices, pointField &sampleLocations) const |
Find cells/faces containing samples. More... | |
tmp< pointField > | samplePoints (const pointField &) const |
Get the sample points given the face points. More... | |
void | calcMapping () const |
Calculate mapping. More... | |
void | calcAMI () const |
Calculate AMI interpolator. More... | |
![]() | |
static tmp< pointField > | readListOrField (const word &keyword, const dictionary &dict, const label size) |
Helper to read field or non-uniform list from dictionary. More... | |
![]() | |
const polyPatch & | patch_ |
Patch to sample. More... | |
word | sampleRegion_ |
Region to sample. More... | |
const sampleMode | mode_ |
What to sample. More... | |
word | samplePatch_ |
Patch (if in sampleMode NEARESTPATCH*) More... | |
const coupleGroupIdentifier | coupleGroup_ |
PatchGroup (if in sampleMode NEARESTPATCH*) More... | |
offsetMode | offsetMode_ |
How to obtain samples. More... | |
vector | offset_ |
Offset vector (uniform) More... | |
vectorField | offsets_ |
Offset vector (nonuniform) More... | |
scalar | distance_ |
Offset distance (normal) More... | |
bool | sameRegion_ |
Same region. More... | |
autoPtr< mapDistribute > | mapPtr_ |
Communication schedule: More... | |
autoPtr< AMIInterpolation > | AMIPtr_ |
Pointer to AMI interpolator. More... | |
const bool | AMIReverse_ |
Flag to indicate that slave patch should be reversed for AMI. More... | |
autoPtr< searchableSurface > | surfPtr_ |
Pointer to projection surface employed by AMI interpolator. More... | |
dictionary | surfDict_ |
Dictionary storing projection surface description. More... | |
This is a mass-fraction boundary condition for a semi-permeable baffle.
This condition models a baffle which is permeable to a some species and impermeable to others. It must be used in conjunction with the corresponding velocity condition, semiPermeableBaffleVelocityFvPatchVectorField.
The mass flux of a species is calculated as a coefficient multiplied by the difference in an input variable across the baffle.
where
![]() | = | flux of the permeable species [kg/s] |
![]() | = | transfer coefficient [kg/m^2/s/<input-dimensions>] |
![]() | = | patch face area [m^2] |
![]() | = | input variable on the patch [] |
![]() | = | input variable on the neighbour patch [] |
A species that the baffle is permeable to will, therefore, have a coefficient greater than zero, whilst a species that does not transfer will have a coefficient equal to zero.
This condition calculates the species flux. The fluxes are summed up by the velocity condition to generate the net mass transfer across the baffle. This mass-fraction condition then generates a corrective diffusive flux to ensure that the correct amounts of the permeable species are transferred.
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 |
Definition at line 136 of file semiPermeableBaffleMassFractionFvPatchScalarField.H.
enum input |
Enumeration for the input variable driving the transfer.
Enumerator | |
---|---|
none | |
massFraction | |
moleFraction | |
partialPressure |
Definition at line 144 of file semiPermeableBaffleMassFractionFvPatchScalarField.H.
semiPermeableBaffleMassFractionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 93 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.
References Foam::Zero.
Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::clone(), semiPermeableBaffleMassFractionFvPatchScalarField::composition(), and semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField().
semiPermeableBaffleMassFractionFvPatchScalarField | ( | const fvPatch & | p, |
const DimensionedField< scalar, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 113 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.
References fvPatchField< Type >::operator=(), semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField(), fvPatch::size(), and Foam::Zero.
semiPermeableBaffleMassFractionFvPatchScalarField | ( | const semiPermeableBaffleMassFractionFvPatchScalarField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< scalar, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given fixedValueTypeFvPatchField.
onto a new patch
Definition at line 141 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.
References semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField().
semiPermeableBaffleMassFractionFvPatchScalarField | ( | const semiPermeableBaffleMassFractionFvPatchScalarField & | ptf | ) |
Copy constructor.
Definition at line 159 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.
References semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField().
semiPermeableBaffleMassFractionFvPatchScalarField | ( | const semiPermeableBaffleMassFractionFvPatchScalarField & | ptf, |
const DimensionedField< scalar, volMesh > & | iF | ||
) |
Copy constructor setting internal field reference.
Definition at line 174 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.
TypeName | ( | "semiPermeableBaffleMassFraction" | ) |
Runtime type information.
|
static |
Access the composition for the given database.
Definition at line 64 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.
References composition, dictionaryName::dictName(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, objectRegistry::foundObject(), objectRegistry::lookupObject(), Foam::name(), and semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField().
Referenced by NamedEnum< directionType, 3 >::names().
|
inlinevirtual |
Construct and return a clone.
Definition at line 219 of file semiPermeableBaffleMassFractionFvPatchScalarField.H.
References semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField().
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Definition at line 236 of file semiPermeableBaffleMassFractionFvPatchScalarField.H.
References semiPermeableBaffleMassFractionFvPatchScalarField::phiY(), semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField(), semiPermeableBaffleMassFractionFvPatchScalarField::updateCoeffs(), and semiPermeableBaffleMassFractionFvPatchScalarField::write().
Foam::tmp< Foam::scalarField > phiY | ( | ) | const |
Return the flux of this species through the baffle.
Definition at line 191 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.
References fvPatch::boundaryMesh(), composition, dictionaryName::dictName(), mapDistribute::distribute(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, fvPatch::index(), objectRegistry::lookupObject(), fvPatch::lookupPatchField(), mappedPatchBase::map(), fvPatchField< Type >::patchInternalField(), basicMultiComponentMixture::species(), thermo, basicThermo::W(), basicSpecieMixture::Wi(), and Foam::Zero.
Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::clone().
|
virtual |
Update the coefficients associated with the patch field.
Definition at line 258 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.
References objectRegistry::lookupObject(), turbulenceModel::propertiesName, and semiPermeableBaffleMassFractionFvPatchScalarField::write().
Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::clone().
|
virtual |
Write.
Reimplemented from mappedPatchBase.
Definition at line 284 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.
References Foam::makePatchTypeField(), mappedPatchBase::write(), fvPatchField< Type >::write(), and Foam::writeEntry().
Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::clone(), and semiPermeableBaffleMassFractionFvPatchScalarField::updateCoeffs().
|
static |
Input variable type names.
Definition at line 153 of file semiPermeableBaffleMassFractionFvPatchScalarField.H.
Referenced by NamedEnum< directionType, 3 >::names().