67 internal_(ptf.internal_),
79 parent_bctype(
p, iF,
dict),
98 const fvMesh& thisMesh = patch().boundaryMesh().mesh();
103 dict_.readIfPresent(
"internal", internal_);
113 if (!extrudeMeshPtr_)
120 baffle_->rename(baffleName);
127 const thermalBaffleFvPatchScalarField& ptf,
131 parent_bctype(ptf, iF),
133 internal_(ptf.internal_),
140void thermalBaffleFvPatchScalarField::createPatchMesh()
142 const fvMesh& thisMesh = patch().boundaryMesh().mesh();
147 List<dictionary> dicts(regionPatches.size());
155 patchTypes[bottomPatchID] = mappedWallPolyPatch::typeName;
159 patchTypes[topPatchID] = mappedWallPolyPatch::typeName;
166 if (dict_.get<
bool>(
"columnCells"))
168 patchTypes[sidePatchID] = emptyPolyPatch::typeName;
172 patchTypes[sidePatchID] = polyPatch::typeName;
177 const word coupleGroup(mpp.coupleGroup());
180 inGroups[0] = coupleGroup;
183 dicts[bottomPatchID].add(
"coupleGroup", coupleGroup);
184 dicts[bottomPatchID].add(
"inGroups", inGroups);
185 dicts[bottomPatchID].add(
"sampleMode", mpp.sampleModeNames_[mpp.mode()]);
186 dicts[bottomPatchID].add(
"samplePatch",
patch().
name());
187 dicts[bottomPatchID].add(
"sampleRegion", thisMesh.name());
192 const word coupleGroupSlave =
193 coupleGroup.substr(0, coupleGroup.find(
'_')) +
"_slave";
195 inGroups[0] = coupleGroupSlave;
196 dicts[topPatchID].add(
"coupleGroup", coupleGroupSlave);
197 dicts[topPatchID].add(
"inGroups", inGroups);
198 dicts[topPatchID].add(
"sampleMode", mpp.sampleModeNames_[mpp.mode()]);
202 forAll(regionPatches, patchi)
205 patchDict.set(
"nFaces", 0);
206 patchDict.set(
"startFace", 0);
217 thisMesh.boundaryMesh()
222 extrudeMeshPtr_.reset
236 const auto& extrPbm = extrudeMeshPtr_().boundaryMesh();
241 const auto& top = extrPbm[topPatchID];
242 const auto& bottom = extrPbm[bottomPatchID];
244 if (top.size() != bottom.size())
247 <<
" size " << top.size()
248 <<
" has different size from bottom patch " << bottom.name()
249 <<
" size " << bottom.size() <<
endl
250 <<
" Disabling mapping offset calculation." <<
endl;
255 const vectorField offsets(bottom.faceCentres()-top.faceCentres());
256 const_cast<mappedPatchBase&
>(*topPtr).setOffset(offsets);
259 <<
"Adjusting patch " << top.name()
266 const auto& thisPbm = thisMesh.boundaryMesh();
267 const auto& groupPatchLookup = thisPbm.groupPatchIDs();
268 const auto&
group = topPtr->coupleGroup();
274 <<
" on region " << thisMesh.name()
275 <<
" contains more than one patch : "
284 const_cast<mappedPatchBase&
>(*thisPp).setOffset(-offsets);
287 <<
"Adjusting patch " << thisPbm[
patchIDs[0]].name()
288 <<
" offsets to " << thisPp->offsets()
294 extrudeMeshPtr_->
write();
322 os.writeEntry(
"extrudeModel", dict_.
get<
word>(
"extrudeModel"));
324 os.writeEntry(
"nLayers", dict_.
get<label>(
"nLayers"));
326 os.writeEntry(
"expansionRatio", dict_.
get<scalar>(
"expansionRatio"));
328 os.writeEntry(
"columnCells", dict_.
get<
Switch>(
"columnCells"));
334 os.writeEntry(
"region", dict_.
get<
word>(
"region"));
336 os.writeEntryIfDifferent<
bool>(
"internal",
true, internal_);
338 os.writeEntry(
"active", dict_.
get<
Switch>(
"active"));
352 thermalBaffleFvPatchScalarField
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
This boundary condition provides a coupled temperature condition between multiple mesh regions.
thermalBaffleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A wrapper for dictionary content, without operators that could affect inheritance patterns.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
void writeEntry(Ostream &os) const
Write sub-dictionary with its dictName as its header.
Top level extrusion model class.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
const word & name() const
Return reference to name.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF)
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
const HashTable< labelList > & groupPatchIDs() const
The patch indices per patch group.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
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.
static autoPtr< thermalBaffleModel > New(const fvMesh &mesh)
Return a reference to the selected model.
A class for handling words, derived from Foam::string.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
#define DebugPoutInFunction
Report an information message using Foam::Pout.
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr const char *const group
Group name for atomic constants.
const std::string patch
OpenFOAM patch number as a std::string.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
List< word > wordList
List of word.
List< label > labelList
A List of labels.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Field< vector > vectorField
Specialisation of Field<T> for vector.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
fvPatchField< scalar > fvPatchScalarField
wordList patchTypes(nPatches)
wordList patchNames(nPatches)
#define forAll(list, i)
Loop across all elements in list.