Public Member Functions | List of all members
semiPermeableBaffleMassFractionFvPatchScalarField Class Reference

This is a mass-fraction boundary condition for a semi-permeable baffle. More...

Inheritance diagram for semiPermeableBaffleMassFractionFvPatchScalarField:
Inheritance graph
[legend]
Collaboration diagram for semiPermeableBaffleMassFractionFvPatchScalarField:
Collaboration graph
[legend]

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 &)
 Construct as copy. More...
 
virtual tmp< fvPatchScalarFieldclone () const
 Construct and return a clone. More...
 
 semiPermeableBaffleMassFractionFvPatchScalarField (const semiPermeableBaffleMassFractionFvPatchScalarField &, const DimensionedField< scalar, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual tmp< fvPatchScalarFieldclone (const DimensionedField< scalar, volMesh > &iF) const
 Construct and return a clone setting internal field reference. More...
 
tmp< scalarFieldphiY () 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...
 
- Public Member Functions inherited from mappedPatchBase
 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 sampleModemode () const
 What to sample. More...
 
const wordsampleRegion () const
 Region to sample. More...
 
const wordsamplePatch () const
 Patch (only if NEARESTPATCHFACE) More...
 
const wordcoupleGroup () const
 PatchGroup (only if NEARESTPATCHFACE) More...
 
label sampleSize () const
 Return size of mapped mesh/patch/boundary. More...
 
const vectoroffset () const
 Offset vector (from patch faces to destination mesh objects) More...
 
const vectorFieldoffsets () const
 Offset vector (from patch faces to destination mesh objects) More...
 
bool sameRegion () const
 Cached sampleRegion != mesh.name() More...
 
const mapDistributemap () const
 Return reference to the parallel distribution map. More...
 
const AMIPatchToPatchInterpolationAMI (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 polyMeshsampleMesh () const
 Get the region mesh. More...
 
const polyPatchsamplePolyPatch () const
 Get the patch on the region. More...
 
tmp< pointFieldsamplePoints () 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...
 

Additional Inherited Members

- Public Types inherited from mappedPatchBase
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...
 
- Static Public Member Functions inherited from mappedPatchBase
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 inherited from mappedPatchBase
static const NamedEnum< sampleMode, 6 > sampleModeNames_
 
static const NamedEnum< offsetMode, 3 > offsetModeNames_
 
- Protected Member Functions inherited from mappedPatchBase
tmp< pointFieldfacePoints (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< pointFieldsamplePoints (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 Protected Member Functions inherited from mappedPatchBase
static tmp< pointFieldreadListOrField (const word &keyword, const dictionary &dict, const label size)
 Helper to read field or non-uniform list from dictionary. More...
 
- Protected Attributes inherited from mappedPatchBase
const polyPatchpatch_
 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< mapDistributemapPtr_
 Communication schedule: More...
 
autoPtr< AMIPatchToPatchInterpolationAMIPtr_
 Pointer to AMI interpolator. More...
 
const bool AMIReverse_
 Flag to indicate that slave patch should be reversed for AMI. More...
 
autoPtr< searchableSurfacesurfPtr_
 Pointer to projection surface employed by AMI interpolator. More...
 
dictionary surfDict_
 Dictionary storing projection surface description. More...
 

Detailed Description

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 mass fraction across the baffle.

\[ \phi_{Yi} = c A (Y_i - Y_{i,n}) \]

where

$ \phi_{Yi} $ = flux of the permeable species [kg/s]
$ c $ = transfer coefficient [kg/m2/s]
$ A $ = patch face area [m2]
$ Y_i $ = mass fraction on the patch []
$ Y_{i,n} $ = mass fraction 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.

Usage
Property Description Req'd? Default
c Transfer coefficient no 0
phi Name of the flux field no phi
See also
Foam::semiPermeableBaffleVelocityFvPatchVectorField
Source files

Definition at line 122 of file semiPermeableBaffleMassFractionFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ semiPermeableBaffleMassFractionFvPatchScalarField() [1/5]

Construct from patch and internal field.

Definition at line 37 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.

References Foam::Zero.

Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::clone(), and semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField().

Here is the caller graph for this function:

◆ semiPermeableBaffleMassFractionFvPatchScalarField() [2/5]

semiPermeableBaffleMassFractionFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const dictionary dict 
)

Construct from patch, internal field and dictionary.

Definition at line 55 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.

References scalarField(), semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField(), fvPatch::size(), and Foam::Zero.

Here is the call graph for this function:

◆ semiPermeableBaffleMassFractionFvPatchScalarField() [3/5]

Construct by mapping given fixedValueTypeFvPatchField.

onto a new patch

Definition at line 76 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.

References semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField().

Here is the call graph for this function:

◆ semiPermeableBaffleMassFractionFvPatchScalarField() [4/5]

Construct as copy.

Definition at line 92 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.

References semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField().

Here is the call graph for this function:

◆ semiPermeableBaffleMassFractionFvPatchScalarField() [5/5]

Construct as copy setting internal field reference.

Definition at line 105 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.

Member Function Documentation

◆ TypeName()

TypeName ( "semiPermeableBaffleMassFraction"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 176 of file semiPermeableBaffleMassFractionFvPatchScalarField.H.

References semiPermeableBaffleMassFractionFvPatchScalarField::semiPermeableBaffleMassFractionFvPatchScalarField().

Here is the call graph for this function:

◆ clone() [2/2]

virtual tmp<fvPatchScalarField> clone ( const DimensionedField< scalar, volMesh > &  iF) const
inlinevirtual

◆ phiY()

Foam::tmp< Foam::scalarField > phiY ( ) const

◆ updateCoeffs()

void updateCoeffs ( )
virtual

Update the coefficients associated with the patch field.

Definition at line 141 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.

References semiPermeableBaffleMassFractionFvPatchScalarField::phiY(), turbulenceModel::propertiesName, and semiPermeableBaffleMassFractionFvPatchScalarField::write().

Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::clone().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ write()

void write ( Ostream os) const
virtual

Write.

Reimplemented from mappedPatchBase.

Definition at line 167 of file semiPermeableBaffleMassFractionFvPatchScalarField.C.

References Foam::makePatchTypeField(), mappedPatchBase::write(), and fvPatchField< Type >::write().

Referenced by semiPermeableBaffleMassFractionFvPatchScalarField::clone(), and semiPermeableBaffleMassFractionFvPatchScalarField::updateCoeffs().

Here is the call graph for this function:
Here is the caller graph for this function:

The documentation for this class was generated from the following files: