This BC solves a steady 1D thermal baffle. More...
#include <thermalBaffle1DFvPatchScalarField.H>


Public Member Functions | |
| TypeName ("compressible::thermalBaffle1D") | |
| Runtime type information. | |
| thermalBaffle1DFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &) | |
| Construct from patch and internal field. | |
| thermalBaffle1DFvPatchScalarField (const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &) | |
| Construct from patch, internal field and dictionary. | |
| thermalBaffle1DFvPatchScalarField (const thermalBaffle1DFvPatchScalarField &, const fvPatch &, const DimensionedField< scalar, volMesh > &, const fvPatchFieldMapper &) | |
| Construct by mapping given. | |
| thermalBaffle1DFvPatchScalarField (const thermalBaffle1DFvPatchScalarField &) | |
| Construct as copy. | |
| thermalBaffle1DFvPatchScalarField (const thermalBaffle1DFvPatchScalarField &, const DimensionedField< scalar, volMesh > &) | |
| Construct as copy setting internal field reference. | |
| virtual tmp< fvPatchField< scalar > > | clone () const |
| Return a clone. | |
| virtual tmp< fvPatchField< scalar > > | clone (const DimensionedField< scalar, volMesh > &iF) const |
| Clone with an internal field reference. | |
| virtual void | autoMap (const fvPatchFieldMapper &) |
| Map (and resize as needed) from self given a mapping object. | |
| virtual void | rmap (const fvPatchScalarField &, const labelList &) |
| Reverse map the given fvPatchField onto this fvPatchField. | |
| virtual void | updateCoeffs () |
| Update the coefficients associated with the patch field. | |
| virtual void | write (Ostream &) const |
| Write. | |
| Public Member Functions inherited from mappedPatchBase | |
| TypeName ("mappedPatchBase") | |
| Runtime type information. | |
| mappedPatchBase (const polyPatch &) | |
| Construct from patch. | |
| mappedPatchBase (const polyPatch &pp, const word &sampleRegion, const sampleMode sampleMode, const word &samplePatch, const vectorField &offsets) | |
| Construct with offsetMode=non-uniform. | |
| mappedPatchBase (const polyPatch &pp, const word &sampleRegion, const sampleMode sampleMode, const word &samplePatch, const vector &uniformOffset) | |
| Construct from offsetMode=uniform. | |
| mappedPatchBase (const polyPatch &pp, const word &sampleRegion, const sampleMode sampleMode, const word &samplePatch, const scalar normalDistance) | |
| Construct from offsetMode=normal and distance. | |
| mappedPatchBase (const polyPatch &, const dictionary &) | |
| Construct from dictionary. | |
| mappedPatchBase (const polyPatch &, const sampleMode, const dictionary &) | |
| Construct from dictionary and (collocated) sample mode. | |
| mappedPatchBase (const polyPatch &, const mappedPatchBase &) | |
| Construct as copy, resetting patch. | |
| mappedPatchBase (const polyPatch &, const mappedPatchBase &, const labelUList &mapAddressing) | |
| Construct as copy, resetting patch, map original data. | |
| virtual | ~mappedPatchBase () |
| Destructor. | |
| void | clearOut () |
| void | setOffset (const scalar normalDist) |
| Change to normal offset with given distance. | |
| void | setOffset (const vector &uniformOffset) |
| Change to uniform offset with value. | |
| void | setOffset (const vectorField &offsets) |
| Change to non-uniform offsets. | |
| sampleMode | mode () const noexcept |
| What to sample. | |
| const word & | sampleWorld () const noexcept |
| World to sample. | |
| const word & | sampleRegion () const |
| Region to sample. | |
| const word & | samplePatch () const |
| Patch (only if NEARESTPATCHFACE). | |
| const word & | coupleGroup () const |
| PatchGroup (only if NEARESTPATCHFACE). | |
| label | sampleSize () const |
| Return size of mapped mesh/patch/boundary. | |
| const vector & | offset () const noexcept |
| Offset vector (from patch faces to destination mesh objects). | |
| const vectorField & | offsets () const noexcept |
| Offset vectors (from patch faces to destination mesh objects). | |
| label | getCommunicator () const |
| Get the communicator (worldComm or world-to-world). | |
| label | comm () const |
| Identical to getCommunicator(). | |
| bool | sameWorld () const |
| Is sample world the local world? | |
| bool | masterWorld () const |
| Is my world ordered before the sampleWorld? | |
| bool | sameRegion () const noexcept |
| Cached sampleRegion != mesh.name(). | |
| uniformDimensionedScalarField & | updateSampleMeshTime () const |
| Local mesh update time. | |
| uniformDimensionedScalarField & | updateMeshTime () const |
| Sample mesh upate time. | |
| bool | upToDate () const |
| const mapDistribute & | map () const |
| Return reference to the parallel distribution map. | |
| const AMIPatchToPatchInterpolation & | AMI (const bool forceUpdate=false) const |
| Return reference to the AMI interpolator. | |
| bool | owner () const |
| Is it owner. | |
| const autoPtr< Foam::searchableSurface > & | surfPtr () const |
| Return a pointer to the AMI projection surface. | |
| const polyMesh & | sampleMesh () const |
| Get the region mesh. | |
| const polyPatch & | samplePolyPatch () const |
| Get the patch on the region. | |
| tmp< pointField > | samplePoints () const |
| Get the sample points. | |
| const fileName & | sampleDatabasePath () const |
| bool | sampleDatabase () const noexcept |
| 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. | |
| template<class Type, class CombineOp> | |
| void | distribute (List< Type > &lst, const CombineOp &cop) const |
| Wrapper around map/interpolate data distribution with operation. | |
| template<class Type> | |
| void | reverseDistribute (List< Type > &lst) const |
| Wrapper around map/interpolate data distribution. | |
| template<class Type, class CombineOp> | |
| void | reverseDistribute (List< Type > &lst, const CombineOp &cop) const |
| Wrapper around map/interpolate data distribution with operation. | |
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. | |
| typedef Tuple2< nearInfo, label > | nearInfoWorld |
| nearest + world | |
| 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: | |
| static fileName | sendPath (const fileName &root, const label proci) |
| Helper: return path to store data to be sent to processor i. | |
| static fileName | receivePath (const fileName &root, const label proci) |
| Helper: return path to store data to be received from processor i. | |
| static FOAM_NO_DANGLING_REFERENCE const objectRegistry & | subRegistry (const objectRegistry &obr, const fileName &path) |
| Lookup (sub)objectRegistry from '/' separated path (relative to objectRegistry). Creates non-existing intermediate ones. | |
| template<class Type> | |
| static void | storeField (objectRegistry &obr, const word &fieldName, const Field< Type > &values) |
| Store an IOField on the objectRegistry relative to obr. | |
| static void | writeDict (const objectRegistry &obr, dictionary &dict) |
| Convert objectRegistry contents into dictionary. | |
| static void | readDict (const dictionary &d, objectRegistry &obr) |
| (recursively) construct and register IOFields from dictionary | |
| Static Public Attributes inherited from mappedPatchBase | |
| static const Enum< sampleMode > | sampleModeNames_ |
| static const Enum< offsetMode > | offsetModeNames_ |
| Protected Member Functions inherited from mappedPatchBase | |
| virtual void | initGeometry (PstreamBuffers &) |
| Initialise the calculation of the patch geometry. | |
| virtual void | calcGeometry (PstreamBuffers &) |
| Calculate the patch geometry. | |
| virtual void | initMovePoints (PstreamBuffers &, const pointField &) |
| Initialise the patches for moving points. | |
| virtual void | movePoints (PstreamBuffers &, const pointField &) |
| Correct patches after moving points. | |
| virtual void | initUpdateMesh (PstreamBuffers &) |
| Initialise the update of the patch topology. | |
| virtual void | updateMesh (PstreamBuffers &) |
| Update of the patch topology. | |
| bool | addWorldConnection () |
| Add a world-world connection. | |
| label | getWorldCommunicator () const |
| Get the communicator for the world-world connection. | |
| const polyMesh & | lookupMesh (const word ®ion) const |
| Lookup mesh. | |
| const polyPatch & | lookupPatch (const word &sampleRegion, const word &samplePatch) const |
| Lookup patch. | |
| tmp< pointField > | facePoints (const polyPatch &) const |
| Get the points from face-centre-decomposition face centres and project them onto the face-diagonal-decomposition triangles. | |
| void | collectSamples (const label mySampleWorld, const pointField &facePoints, pointField &samples, labelList &patchFaceWorlds, labelList &patchFaceProcs, labelList &patchFaces, pointField &patchFc) const |
| Collect single list of samples and originating processor+face + wanted world. | |
| 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. | |
| 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. | |
| tmp< pointField > | samplePoints (const pointField &) const |
| Get the sample points given the face points. | |
| void | calcMapping () const |
| Calculate mapping. | |
| void | calcAMI () const |
| Calculate AMI interpolator. | |
| Static Protected Member Functions inherited from mappedPatchBase | |
| static autoPtr< fileName > | readDatabase (const dictionary &dict) |
| Read optional database name from dictionary. | |
| static FOAM_NO_DANGLING_REFERENCE const objectRegistry & | subRegistry (const objectRegistry &obr, const wordList &names, const label index) |
| Lookup (sub)objectRegistry by following names of sub registries. Creates non-existing intermediate ones. | |
| template<class Type> | |
| static bool | writeIOField (const regIOobject &obj, dictionary &dict) |
| Attempt to detect an IOField<Type> and write to dictionary. | |
| 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. | |
| Protected Attributes inherited from mappedPatchBase | |
| const polyPatch & | patch_ |
| Patch to sample. | |
| word | sampleWorld_ |
| World to sample. | |
| word | sampleRegion_ |
| Region to sample. | |
| const sampleMode | mode_ |
| What to sample. | |
| word | samplePatch_ |
| Patch (if in sampleMode NEARESTPATCH*). | |
| const coupleGroupIdentifier | coupleGroup_ |
| PatchGroup (if in sampleMode NEARESTPATCH*). | |
| const autoPtr< fileName > | sampleDatabasePtr_ |
| Empty or location of database. | |
| offsetMode | offsetMode_ |
| How to obtain samples. | |
| vector | offset_ |
| Offset vector (uniform). | |
| vectorField | offsets_ |
| Offset vector (nonuniform). | |
| scalar | distance_ |
| Offset distance (normal). | |
| label | communicator_ |
| Communicator. | |
| bool | sameRegion_ |
| Same region. | |
| autoPtr< mapDistribute > | mapPtr_ |
| Communication schedule: | |
| const bool | AMIReverse_ |
| Flag to indicate that slave patch should be reversed for AMI. | |
| autoPtr< AMIPatchToPatchInterpolation > | AMIPtr_ |
| Pointer to AMI interpolator. | |
| autoPtr< searchableSurface > | surfPtr_ |
| Pointer to projection surface employed by AMI interpolator. | |
| dictionary | surfDict_ |
| Dictionary storing projection surface description. | |
| autoPtr< uniformDimensionedScalarField > | updateMeshTimePtr_ |
| Local mesh update time. | |
| autoPtr< uniformDimensionedScalarField > | updateSampleMeshTimePtr_ |
| Sample mesh update time. | |
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.
<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;
}
Definition at line 107 of file thermalBaffle1DFvPatchScalarField.H.
| thermalBaffle1DFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF ) |
Construct from patch and internal field.
Definition at line 37 of file thermalBaffle1DFvPatchScalarField.C.
References mappedPatchBase::mappedPatchBase(), and p.
Referenced by rmap(), thermalBaffle1DFvPatchScalarField(), thermalBaffle1DFvPatchScalarField(), and thermalBaffle1DFvPatchScalarField().


| thermalBaffle1DFvPatchScalarField | ( | const fvPatch & | p, |
| const DimensionedField< scalar, volMesh > & | iF, | ||
| const dictionary & | dict ) |
Construct from patch, internal field and dictionary.
Definition at line 83 of file thermalBaffle1DFvPatchScalarField.C.
References dict, mappedPatchBase::mappedPatchBase(), mappedPatchBase::NEARESTPATCHFACE, p, and Foam::Zero.

| 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 59 of file thermalBaffle1DFvPatchScalarField.C.
References mappedPatchBase::mappedPatchBase(), p, and thermalBaffle1DFvPatchScalarField().

| thermalBaffle1DFvPatchScalarField | ( | const thermalBaffle1DFvPatchScalarField< solidType > & | ptf | ) |
Construct as copy.
Definition at line 139 of file thermalBaffle1DFvPatchScalarField.C.
References mappedPatchBase::mappedPatchBase(), and thermalBaffle1DFvPatchScalarField().

| thermalBaffle1DFvPatchScalarField | ( | const thermalBaffle1DFvPatchScalarField< solidType > & | ptf, |
| const DimensionedField< scalar, volMesh > & | iF ) |
Construct as copy setting internal field reference.
Definition at line 160 of file thermalBaffle1DFvPatchScalarField.C.
References mappedPatchBase::mappedPatchBase(), and thermalBaffle1DFvPatchScalarField().

| TypeName | ( | "compressible::thermalBaffle1D" | ) |
Runtime type information.
|
inlinevirtual |
Return a clone.
Definition at line 245 of file thermalBaffle1DFvPatchScalarField.H.
|
inlinevirtual |
Clone with an internal field reference.
Definition at line 253 of file thermalBaffle1DFvPatchScalarField.H.
|
virtual |
Map (and resize as needed) from self given a mapping object.
Definition at line 288 of file thermalBaffle1DFvPatchScalarField.C.
References mappedPatchBase::clearOut().

|
virtual |
Reverse map the given fvPatchField onto this fvPatchField.
Definition at line 306 of file thermalBaffle1DFvPatchScalarField.C.
References Foam::refCast(), and thermalBaffle1DFvPatchScalarField().

|
virtual |
Update the coefficients associated with the patch field.
Definition at line 326 of file thermalBaffle1DFvPatchScalarField.C.
References alpha, mapDistribute::distribute(), Foam::endl(), forAll, Foam::gAverage(), Foam::gMinMax(), UPstream::incrMsgType(), Foam::Info, ThermalDiffusivity< BasicTurbulenceModel >::kappaEff(), limits, mappedPatchBase::map(), UPstream::msgType(), fvPatch::name(), turbulenceModel::propertiesName, mappedPatchBase::samplePolyPatch(), UList< T >::size(), TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::transport(), and Foam::Zero.

|
virtual |
Write.
Reimplemented from mappedPatchBase.
Definition at line 420 of file thermalBaffle1DFvPatchScalarField.C.
References os(), mappedPatchBase::write(), and mixedFvPatchField< Type >::write().
