Loading...
Searching...
No Matches
thermalBaffleFvPatchScalarField.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020-2024 OpenCFD Ltd
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "dictionaryContent.H"
32#include "emptyPolyPatch.H"
33#include "mappedWallPolyPatch.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace compressible
40{
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
45(
46 const fvPatch& p,
48)
49:
50 parent_bctype(p, iF),
51 owner_(false),
52 internal_(true),
53 dict_()
54{}
55
56
58(
60 const fvPatch& p,
62 const fvPatchFieldMapper& mapper
63)
64:
65 parent_bctype(ptf, p, iF, mapper),
66 owner_(ptf.owner_),
67 internal_(ptf.internal_),
68 dict_(ptf.dict_)
69{}
70
71
73(
74 const fvPatch& p,
76 const dictionary& dict
77)
78:
79 parent_bctype(p, iF, dict),
80 owner_(false),
81 internal_(true),
82 dict_
83 (
84 // Copy dictionary, but without "heavy" data chunks
85 dictionaryContent::copyDict
86 (
87 dict,
88 wordList(), // allow
89 wordList // deny
90 ({
91 "type", // redundant
92 "value"
93 })
94 )
95 )
96{
97
98 const fvMesh& thisMesh = patch().boundaryMesh().mesh();
99
100 word regionName("none");
101 dict_.readIfPresent("region", regionName);
102
103 dict_.readIfPresent("internal", internal_);
104
105 const word baffleName("3DBaffle" + regionName);
106
107 if
108 (
109 (regionName != "none")
110 && !thisMesh.time().foundObject<fvMesh>(regionName)
111 )
112 {
113 if (!extrudeMeshPtr_)
114 {
115 createPatchMesh();
116 }
117
118 baffle_.reset(baffleType::New(thisMesh, dict));
119 owner_ = true;
120 baffle_->rename(baffleName);
121 }
122}
123
124
126(
127 const thermalBaffleFvPatchScalarField& ptf,
129)
130:
131 parent_bctype(ptf, iF),
132 owner_(ptf.owner_),
133 internal_(ptf.internal_),
134 dict_(ptf.dict_)
135{}
136
137
138// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139
140void thermalBaffleFvPatchScalarField::createPatchMesh()
141{
142 const fvMesh& thisMesh = patch().boundaryMesh().mesh();
143
144 const word regionName(dict_.get<word>("region"));
145
146 polyPatchList regionPatches(3);
147 List<dictionary> dicts(regionPatches.size());
148 List<word> patchNames(regionPatches.size());
149 List<word> patchTypes(regionPatches.size());
150
151 patchNames[bottomPatchID] = word("bottom");
152 patchNames[sidePatchID] = word("side");
153 patchNames[topPatchID] = word("top");
154
155 patchTypes[bottomPatchID] = mappedWallPolyPatch::typeName;
156
157 if (internal_)
158 {
159 patchTypes[topPatchID] = mappedWallPolyPatch::typeName;
160 }
161 else
162 {
163 patchTypes[topPatchID] = polyPatch::typeName;
164 }
165
166 if (dict_.get<bool>("columnCells"))
167 {
168 patchTypes[sidePatchID] = emptyPolyPatch::typeName;
169 }
170 else
171 {
172 patchTypes[sidePatchID] = polyPatch::typeName;
173 }
174
175 const auto& mpp = refCast<const mappedPatchBase>(patch().patch(), dict_);
176
177 const word coupleGroup(mpp.coupleGroup());
178
179 wordList inGroups(1);
180 inGroups[0] = coupleGroup;
181
182 // The bottomPatchID is coupled with this patch
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());
188
189 // Internal baffle needs a coupled on the topPatchID
190 if (internal_)
191 {
192 const word coupleGroupSlave =
193 coupleGroup.substr(0, coupleGroup.find('_')) + "_slave";
194
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()]);
199 }
200
201
202 forAll(regionPatches, patchi)
203 {
204 dictionary& patchDict = dicts[patchi];
205 patchDict.set("nFaces", 0);
206 patchDict.set("startFace", 0);
207
208 regionPatches.set
209 (
210 patchi,
212 (
213 patchTypes[patchi],
214 patchNames[patchi],
215 dicts[patchi],
216 patchi,
217 thisMesh.boundaryMesh()
218 )
219 );
220 }
221
222 extrudeMeshPtr_.reset
223 (
224 new extrudePatchMesh
225 (
226 thisMesh,
227 patch(),
228 dict_,
230 regionPatches
231 )
232 );
233
234 // Adjust top patch for the thickness - it needs to subtract the offset
235 // distance when trying to do the mapping.
236 const auto& extrPbm = extrudeMeshPtr_().boundaryMesh();
237 const auto* topPtr = isA<const mappedPatchBase>(extrPbm[topPatchID]);
238
239 if (topPtr)
240 {
241 const auto& top = extrPbm[topPatchID];
242 const auto& bottom = extrPbm[bottomPatchID];
243
244 if (top.size() != bottom.size())
245 {
246 WarningInFunction<< "Top patch " << top.name()
247 << " size " << top.size()
248 << " has different size from bottom patch " << bottom.name()
249 << " size " << bottom.size() << endl
250 << " Disabling mapping offset calculation." << endl;
251 }
252 else
253 {
254 // Adjust top patch offsets
255 const vectorField offsets(bottom.faceCentres()-top.faceCentres());
256 const_cast<mappedPatchBase&>(*topPtr).setOffset(offsets);
257
259 << "Adjusting patch " << top.name()
260 << " offsets to " << flatOutput(offsets)
261 << endl;
262
263 if (internal_)
264 {
265 // Find other side of the baffle using its group
266 const auto& thisPbm = thisMesh.boundaryMesh();
267 const auto& groupPatchLookup = thisPbm.groupPatchIDs();
268 const auto& group = topPtr->coupleGroup();
269 const labelList patchIDs(groupPatchLookup[group]);
270
271 if (patchIDs.size() != 1)
272 {
273 FatalErrorInFunction<< "Group " << group
274 << " on region " << thisMesh.name()
275 << " contains more than one patch : "
276 << patchIDs
277 << exit(FatalError);
278 }
279
280 const auto* thisPp =
282 if (thisPp)
283 {
284 const_cast<mappedPatchBase&>(*thisPp).setOffset(-offsets);
285
287 << "Adjusting patch " << thisPbm[patchIDs[0]].name()
288 << " offsets to " << thisPp->offsets()
289 << endl;
290 }
291 }
292
293 // Enforce re-writing so baffleType::New reads updated mesh
294 extrudeMeshPtr_->write();
295 }
296 }
297}
298
299
301{
302 if (this->updated())
303 {
304 return;
305 }
306
307 if (owner_)
308 {
309 baffle_->evolve();
310 }
311
313}
314
315
317{
319
320 if (owner_)
321 {
322 os.writeEntry("extrudeModel", dict_.get<word>("extrudeModel"));
323
324 os.writeEntry("nLayers", dict_.get<label>("nLayers"));
325
326 os.writeEntry("expansionRatio", dict_.get<scalar>("expansionRatio"));
327
328 os.writeEntry("columnCells", dict_.get<Switch>("columnCells"));
329
330 const word extrudeModel(dict_.get<word>("extrudeModel") + "Coeffs");
331
333
334 os.writeEntry("region", dict_.get<word>("region"));
335
336 os.writeEntryIfDifferent<bool>("internal", true, internal_);
337
338 os.writeEntry("active", dict_.get<Switch>("active"));
339
340 dict_.subDict("thermoType").writeEntry("thermoType", os);
341 dict_.subDict("mixture").writeEntry("mixture", os);
342 dict_.subDict("radiation").writeEntry("radiation", os);
343 }
344}
345
346
347// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
348
350(
352 thermalBaffleFvPatchScalarField
353);
354
355
356// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357
358} // End namespace compressible
359} // End namespace Foam
360
361
362// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
labelList patchIDs
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...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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,...
Definition Switch.H:81
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 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,...
Definition dictionary.H:133
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.
Definition dictionary.C:441
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.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const word & name() const
Return reference to name.
Definition fvMesh.H:387
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
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.
Definition polyMesh.H:609
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.
Definition word.H:66
volScalarField & p
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
auto & name
#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.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
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)
Definition errorManip.H:125
fvPatchField< scalar > fvPatchScalarField
wordList patchTypes(nPatches)
wordList patchNames(nPatches)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299