Go to the documentation of this file.
74 internal_(ptf.internal_),
96 const fvMesh& thisMesh =
patch().boundaryMesh().mesh();
103 dict_.readIfPresent(
"internal", internal_);
113 if (!extrudeMeshPtr_)
120 baffle_->rename(baffleName);
133 internal_(ptf.internal_),
142 void thermalBaffleFvPatchScalarField::createPatchMesh()
144 const fvMesh& thisMesh =
patch().boundaryMesh().mesh();
157 patchTypes[bottomPatchID] = mappedWallPolyPatch::typeName;
161 patchTypes[topPatchID] = mappedWallPolyPatch::typeName;
168 if (dict_.
get<
bool>(
"columnCells"))
170 patchTypes[sidePatchID] = emptyPolyPatch::typeName;
174 patchTypes[sidePatchID] = polyPatch::typeName;
177 const mappedPatchBase& mpp =
178 refCast<const mappedPatchBase>(
patch().
patch(), dict_);
180 const word coupleGroup(mpp.coupleGroup());
183 inGroups[0] = coupleGroup;
186 dicts[bottomPatchID].add(
"coupleGroup", coupleGroup);
187 dicts[bottomPatchID].add(
"inGroups", inGroups);
188 dicts[bottomPatchID].add(
"sampleMode", mpp.sampleModeNames_[mpp.mode()]);
189 dicts[bottomPatchID].add(
"samplePatch",
patch().
name());
190 dicts[bottomPatchID].add(
"sampleRegion", thisMesh.
name());
195 const word coupleGroupSlave =
196 coupleGroup.substr(0, coupleGroup.find(
'_')) +
"_slave";
198 inGroups[0] = coupleGroupSlave;
199 dicts[topPatchID].add(
"coupleGroup", coupleGroupSlave);
200 dicts[topPatchID].add(
"inGroups", inGroups);
201 dicts[topPatchID].add(
"sampleMode", mpp.sampleModeNames_[mpp.mode()]);
205 forAll(regionPatches, patchi)
207 dictionary& patchDict = dicts[patchi];
208 patchDict.set(
"nFaces", 0);
209 patchDict.set(
"startFace", 0);
221 extrudeMeshPtr_.reset
257 os.writeEntry(
"extrudeModel", dict_.
get<
word>(
"extrudeModel"));
259 os.writeEntry(
"nLayers", dict_.
get<label>(
"nLayers"));
261 os.writeEntry(
"expansionRatio", dict_.
get<scalar>(
"expansionRatio"));
263 os.writeEntry(
"columnCells", dict_.
get<
Switch>(
"columnCells"));
269 os.writeEntry(
"region", dict_.
get<
word>(
"region"));
271 os.writeEntryIfDifferent<
bool>(
"internal",
true, internal_);
273 os.writeEntry(
"active", dict_.
get<
Switch>(
"active"));
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
A class for handling words, derived from Foam::string.
Top level extrusion model class.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
#define forAll(list, i)
Loop across all elements in list.
wordList patchTypes(nPatches)
List< word > wordList
A List of words.
void writeEntry(Ostream &os) const
Write sub-dictionary with its dictName as its header.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
wordList patchNames(nPatches)
virtual void write(Ostream &os) const
Write.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
This boundary condition provides a coupled temperature condition between multiple mesh regions.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
thermalBaffleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
const std::string patch
OpenFOAM patch number as a std::string.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
const Time & time() const
Return the top-level database.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const word & name() const
Return reference to name.
virtual void write(Ostream &) const
Write.