thermalBaffle1DFvPatchScalarField< solidType > Class Template Reference

This BC solves a steady 1D thermal baffle. More...

Inheritance diagram for thermalBaffle1DFvPatchScalarField< solidType >:
[legend]
Collaboration diagram for thermalBaffle1DFvPatchScalarField< solidType >:
[legend]

Public Member Functions

 TypeName ("compressible::thermalBaffle1D")
 Runtime type information. More...
 
 thermalBaffle1DFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &)
 Construct from patch and internal field. More...
 
 thermalBaffle1DFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
 Construct from patch, internal field and dictionary. More...
 
 thermalBaffle1DFvPatchScalarField (const thermalBaffle1DFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &)
 Construct by mapping given. More...
 
 thermalBaffle1DFvPatchScalarField (const thermalBaffle1DFvPatchScalarField &)
 Construct as copy. More...
 
virtual tmp< fvPatchScalarFieldclone () const
 Construct and return a clone. More...
 
 thermalBaffle1DFvPatchScalarField (const thermalBaffle1DFvPatchScalarField &, 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...
 
virtual void autoMap (const fvPatchFieldMapper &)
 Map (and resize as needed) from self given a mapping object. More...
 
virtual void rmap (const fvPatchScalarField &, const labelList &)
 Reverse map the given fvPatchField onto this fvPatchField. 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 &uniformOffset)
 Construct from offsetMode=uniform. More...
 
 mappedPatchBase (const polyPatch &pp, const word &sampleRegion, const sampleMode sampleMode, const word &samplePatch, const scalar normalDistance)
 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 ()
 
void setOffset (const scalar normalDist)
 Change to normal offset with given distance. More...
 
void setOffset (const vector &uniformOffset)
 Change to uniform offset with value. More...
 
void setOffset (const vectorField &offsets)
 Change to non-uniform offsets. More...
 
sampleMode mode () const noexcept
 What to sample. More...
 
const wordsampleWorld () const noexcept
 World 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 noexcept
 Offset vector (from patch faces to destination mesh objects) More...
 
const vectorFieldoffsets () const noexcept
 Offset vectors (from patch faces to destination mesh objects) More...
 
label getCommunicator () const
 Get the communicator (worldComm or world-to-world) More...
 
label comm () const
 Identical to getCommunicator() More...
 
bool sameWorld () const
 Is sample world the local world? More...
 
bool masterWorld () const
 Is my world ordered before the sampleWorld? More...
 
bool sameRegion () const noexcept
 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...
 
bool owner () const
 Is it owner. 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...
 
const fileNamesampleDatabasePath () const
 
bool sampleDatabase () const
 
virtual fileName sendPath (const label proci) const
 
virtual fileName receivePath (const label proci) const
 
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...
 
typedef Tuple2< nearInfo, label > nearInfoWorld
 nearest + world 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 fileName sendPath (const fileName &root, const label proci)
 Helper: return path to store data to be sent to processor i. More...
 
static fileName receivePath (const fileName &root, const label proci)
 
static const objectRegistrysubRegistry (const objectRegistry &obr, const fileName &path)
 
template<class Type >
static void storeField (objectRegistry &obr, const word &fieldName, const Field< Type > &values)
 Store an IOField on the objectRegistry relative to obr. More...
 
static void writeDict (const objectRegistry &obr, dictionary &dict)
 Convert objectRegistry contents into dictionary. More...
 
static void readDict (const dictionary &d, objectRegistry &obr)
 (recursively) construct and register IOFields from dictionary More...
 
- Static Public Attributes inherited from mappedPatchBase
static const Enum< sampleModesampleModeNames_
 
static const Enum< offsetModeoffsetModeNames_
 
- Protected Member Functions inherited from mappedPatchBase
bool addWorldConnection ()
 Add a world-world connection. More...
 
label getWorldCommunicator () const
 Get the communicator for the world-world connection. More...
 
const polyMeshlookupMesh (const word &region) const
 Lookup mesh. More...
 
const polyPatchlookupPatch (const word &sampleRegion, const word &samplePatch) const
 Lookup patch. More...
 
tmp< pointFieldfacePoints (const polyPatch &) const
 
void collectSamples (const label mySampleWorld, const pointField &facePoints, pointField &samples, labelList &patchFaceWorlds, labelList &patchFaceProcs, labelList &patchFaces, pointField &patchFc) const
 
void findLocalSamples (const sampleMode mode, const label sampleWorld, const word &sampleRegion, const word &samplePatch, const pointField &samplePoints, List< nearInfoWorld > &nearest) const
 Find (local) cells/faces containing samples. More...
 
void findSamples (const sampleMode mode, const label myWorldIndex, const pointField &, const labelList &wantedWorlds, const labelList &origProcs, labelList &sampleProcs, labelList &sampleIndices, pointField &sampleLocations) const
 Find (global) 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 autoPtr< fileNamereadDatabase (const dictionary &dict)
 Read optional database name from dictionary. More...
 
static const objectRegistrysubRegistry (const objectRegistry &obr, const wordList &names, const label index)
 
template<class Type >
static bool writeIOField (const regIOobject &obj, dictionary &dict)
 Attempt to detect an IOField<Type> and write to dictionary. More...
 
template<class Type >
static bool constructIOField (const word &name, token &tok, Istream &is, objectRegistry &obr)
 Attempt to read an IOField<Type> and store on objectRegistry. More...
 
- Protected Attributes inherited from mappedPatchBase
const polyPatchpatch_
 Patch to sample. More...
 
word sampleWorld_
 World 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...
 
const autoPtr< fileNamesampleDatabasePtr_
 Empty or location of database. 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...
 
label communicator_
 Communicator. More...
 
bool sameRegion_
 Same region. More...
 
autoPtr< mapDistributemapPtr_
 Communication schedule: More...
 
const bool AMIReverse_
 Flag to indicate that slave patch should be reversed for AMI. More...
 
autoPtr< AMIPatchToPatchInterpolationAMIPtr_
 Pointer to AMI interpolator. More...
 
autoPtr< searchableSurfacesurfPtr_
 Pointer to projection surface employed by AMI interpolator. More...
 
dictionary surfDict_
 Dictionary storing projection surface description. More...
 

Detailed Description

template<class solidType>
class Foam::compressible::thermalBaffle1DFvPatchScalarField< solidType >

This BC solves a steady 1D thermal baffle.

The solid properties are specify as dictionary. Optionally radiative heat flux (qr) can be incorporated into the balance. Some under-relaxation might be needed on qr. Baffle and solid properties need to be specified on the master side of the baffle.

Usage
Example of the boundary condition specification using constant solid thermo :
<masterPatchName>
{
    type   compressible::thermalBaffle1D<hConstSolidThermoPhysics>;
    samplePatch     <slavePatchName>;

    thickness       uniform 0.005;  // thickness [m]
    qs              uniform 100;    // heat flux [W/m2]

    qr              none;
    relaxation      1;

    // Solid thermo
    specie
    {
        molWeight       20;
    }
    transport
    {
        kappa           1;
    }
    thermodynamics
    {
        Hf              0;
        Cp              10;
    }
    equationOfState
    {
        rho             10;
    }

    value               uniform 300;
}

<slavePatchName>
{
    type   compressible::thermalBaffle1D<hConstSolidThermoPhysics>;
    samplePatch     <masterPatchName>;

    qr              none;
    relaxation      1;
}
Source files

Definition at line 111 of file thermalBaffle1DFvPatchScalarField.H.

Constructor & Destructor Documentation

◆ thermalBaffle1DFvPatchScalarField() [1/5]

thermalBaffle1DFvPatchScalarField ( const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF 
)

Construct from patch and internal field.

Definition at line 46 of file thermalBaffle1DFvPatchScalarField.C.

Referenced by thermalBaffle1DFvPatchScalarField< solidType >::clone().

Here is the caller graph for this function:

◆ thermalBaffle1DFvPatchScalarField() [2/5]

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

Construct from patch, internal field and dictionary.

Definition at line 92 of file thermalBaffle1DFvPatchScalarField.C.

◆ thermalBaffle1DFvPatchScalarField() [3/5]

thermalBaffle1DFvPatchScalarField ( const thermalBaffle1DFvPatchScalarField< solidType > &  ptf,
const fvPatch p,
const DimensionedField< scalar, volMesh > &  iF,
const fvPatchFieldMapper mapper 
)

Construct by mapping given.

thermalBaffle1DFvPatchScalarField onto a new patch

Definition at line 68 of file thermalBaffle1DFvPatchScalarField.C.

◆ thermalBaffle1DFvPatchScalarField() [4/5]

Construct as copy.

Definition at line 151 of file thermalBaffle1DFvPatchScalarField.C.

◆ thermalBaffle1DFvPatchScalarField() [5/5]

thermalBaffle1DFvPatchScalarField ( const thermalBaffle1DFvPatchScalarField< solidType > &  ptf,
const DimensionedField< scalar, volMesh > &  iF 
)

Construct as copy setting internal field reference.

Definition at line 172 of file thermalBaffle1DFvPatchScalarField.C.

Member Function Documentation

◆ TypeName()

TypeName ( "compressible::thermalBaffle1D"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchScalarField> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 201 of file thermalBaffle1DFvPatchScalarField.H.

References thermalBaffle1DFvPatchScalarField< solidType >::thermalBaffle1DFvPatchScalarField().

Here is the call graph for this function:

◆ clone() [2/2]

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

Construct and return a clone setting internal field reference.

Definition at line 218 of file thermalBaffle1DFvPatchScalarField.H.

References thermalBaffle1DFvPatchScalarField< solidType >::thermalBaffle1DFvPatchScalarField().

Here is the call graph for this function:

◆ autoMap()

void autoMap ( const fvPatchFieldMapper m)
virtual

Map (and resize as needed) from self given a mapping object.

Definition at line 300 of file thermalBaffle1DFvPatchScalarField.C.

References mappedPatchBase::clearOut().

Here is the call graph for this function:

◆ rmap()

void rmap ( const fvPatchScalarField ptf,
const labelList addr 
)
virtual

Reverse map the given fvPatchField onto this fvPatchField.

Definition at line 318 of file thermalBaffle1DFvPatchScalarField.C.

◆ updateCoeffs()

void updateCoeffs ( )
virtual

Update the coefficients associated with the patch field.

Definition at line 337 of file thermalBaffle1DFvPatchScalarField.C.

References Foam::constant::atomic::alpha, Foam::expressions::patchExpr::debug, mapDistribute::distribute(), Foam::endl(), forAll, Foam::gAverage(), Foam::gMax(), Foam::gMin(), Foam::Info, ThermalDiffusivity< BasicTurbulenceModel >::kappaEff(), mappedPatchBase::map(), UPstream::msgType(), fvPatch::name(), Foam::foamVersion::patch, turbulenceModel::propertiesName, Foam::fvc::snGrad(), and Foam::Zero.

Here is the call graph for this function:

◆ write()

void write ( Ostream os) const
virtual

Write.

Reimplemented from mappedPatchBase.

Definition at line 430 of file thermalBaffle1DFvPatchScalarField.C.

References os(), Foam::vtk::write(), mappedPatchBase::write(), and Ostream::writeEntry().

Here is the call graph for this function:

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