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 OpenCFD Ltd
10 -------------------------------------------------------------------------------
11 License
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 
31 #include "emptyPolyPatch.H"
32 #include "mappedWallPolyPatch.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace compressible
39 {
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const fvPatch& p,
47 )
48 :
50  owner_(false),
51  internal_(true),
52  baffle_(),
53  dict_(dictionary::null),
54  extrudeMeshPtr_()
55 {}
56 
57 
59 (
61  const fvPatch& p,
63  const fvPatchFieldMapper& mapper
64 )
65 :
67  (
68  ptf,
69  p,
70  iF,
71  mapper
72  ),
73  owner_(ptf.owner_),
74  internal_(ptf.internal_),
75  baffle_(),
76  dict_(ptf.dict_),
77  extrudeMeshPtr_()
78 {}
79 
80 
82 (
83  const fvPatch& p,
85  const dictionary& dict
86 )
87 :
89  owner_(false),
90  internal_(true),
91  baffle_(),
92  dict_(dict),
93  extrudeMeshPtr_()
94 {
95 
96  const fvMesh& thisMesh = patch().boundaryMesh().mesh();
97 
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  !thisMesh.time().foundObject<fvMesh>(regionName)
110  && regionName != "none"
111  )
112  {
113  if (!extrudeMeshPtr_)
114  {
115  createPatchMesh();
116  }
117 
118  baffle_.reset(baffle::New(thisMesh, dict).ptr());
119  owner_ = true;
120  baffle_->rename(baffleName);
121  }
122 }
123 
124 
126 (
129 )
130 :
132  owner_(ptf.owner_),
133  internal_(ptf.internal_),
134  baffle_(),
135  dict_(ptf.dict_),
136  extrudeMeshPtr_()
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
142 
144 (
145  const fvPatchFieldMapper& m
146 )
147 {
148  mixedFvPatchScalarField::autoMap(m);
149 }
150 
151 
153 (
154  const fvPatchScalarField& ptf,
155  const labelList& addr
156 )
157 {
158  mixedFvPatchScalarField::rmap(ptf, addr);
159 }
160 
161 
162 void thermalBaffleFvPatchScalarField::createPatchMesh()
163 {
164  const fvMesh& thisMesh = patch().boundaryMesh().mesh();
165 
166  const word regionName(dict_.get<word>("region"));
167 
168  List<polyPatch*> regionPatches(3);
169  List<word> patchNames(regionPatches.size());
170  List<word> patchTypes(regionPatches.size());
171  List<dictionary> dicts(regionPatches.size());
172 
173  patchNames[bottomPatchID] = word("bottom");
174  patchNames[sidePatchID] = word("side");
175  patchNames[topPatchID] = word("top");
176 
177  patchTypes[bottomPatchID] = mappedWallPolyPatch::typeName;
178 
179  if (internal_)
180  {
181  patchTypes[topPatchID] = mappedWallPolyPatch::typeName;
182  }
183  else
184  {
185  patchTypes[topPatchID] = polyPatch::typeName;
186  }
187 
188  if (dict_.get<bool>("columnCells"))
189  {
190  patchTypes[sidePatchID] = emptyPolyPatch::typeName;
191  }
192  else
193  {
194  patchTypes[sidePatchID] = polyPatch::typeName;
195  }
196 
197  const mappedPatchBase& mpp =
198  refCast<const mappedPatchBase>(patch().patch(), dict_);
199 
200  const word coupleGroup(mpp.coupleGroup());
201 
202  wordList inGroups(1);
203  inGroups[0] = coupleGroup;
204 
205  // The bottomPatchID is coupled with this patch
206  dicts[bottomPatchID].add("coupleGroup", coupleGroup);
207  dicts[bottomPatchID].add("inGroups", inGroups);
208  dicts[bottomPatchID].add("sampleMode", mpp.sampleModeNames_[mpp.mode()]);
209  dicts[bottomPatchID].add("samplePatch", patch().name());
210  dicts[bottomPatchID].add("sampleRegion", thisMesh.name());
211 
212  // Internal baffle needs a coupled on the topPatchID
213  if (internal_)
214  {
215  const word coupleGroupSlave =
216  coupleGroup.substr(0, coupleGroup.find('_')) + "_slave";
217 
218  inGroups[0] = coupleGroupSlave;
219  dicts[topPatchID].add("coupleGroup", coupleGroupSlave);
220  dicts[topPatchID].add("inGroups", inGroups);
221  dicts[topPatchID].add("sampleMode", mpp.sampleModeNames_[mpp.mode()]);
222  }
223 
224 
225  forAll(regionPatches, patchi)
226  {
227  dictionary& patchDict = dicts[patchi];
228  patchDict.set("nFaces", 0);
229  patchDict.set("startFace", 0);
230 
231  regionPatches[patchi] = polyPatch::New
232  (
233  patchTypes[patchi],
234  patchNames[patchi],
235  dicts[patchi],
236  patchi,
237  thisMesh.boundaryMesh()
238  ).ptr();
239  }
240 
241  extrudeMeshPtr_.reset
242  (
243  new extrudePatchMesh
244  (
245  thisMesh,
246  patch(),
247  dict_,
248  regionName,
249  regionPatches
250  )
251  );
252 }
253 
254 
256 {
257  if (this->updated())
258  {
259  return;
260  }
261 
262  if (owner_)
263  {
264  baffle_->evolve();
265  }
266 
268 }
269 
270 
272 {
274 
275  if (owner_)
276  {
277  os.writeEntry("extrudeModel", dict_.get<word>("extrudeModel"));
278 
279  os.writeEntry("nLayers", dict_.get<label>("nLayers"));
280 
281  os.writeEntry("expansionRatio", dict_.get<scalar>("expansionRatio"));
282 
283  os.writeEntry("columnCells", dict_.get<Switch>("columnCells"));
284 
285  const word extrudeModel(dict_.get<word>("extrudeModel") + "Coeffs");
286 
288 
289  os.writeEntry("region", dict_.get<word>("region"));
290 
291  os.writeEntryIfDifferent<bool>("internal", true, internal_);
292 
293  os.writeEntry("active", dict_.get<Switch>("active"));
294 
295  dict_.subDict("thermoType").writeEntry("thermoType", os);
296  dict_.subDict("mixture").writeEntry("mixture", os);
297  dict_.subDict("radiation").writeEntry("radiation", os);
298  }
299 }
300 
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
305 (
308 );
309 
310 
311 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
312 
313 } // End namespace compressible
314 } // End namespace Foam
315 
316 
317 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:244
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
mappedWallPolyPatch.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::extrudeModel
Top level extrusion model class.
Definition: extrudeModel.H:76
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.C:179
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
Foam::compressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
patchTypes
wordList patchTypes(nPatches)
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:385
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:59
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::regionModels::thermalBaffleModels::thermalBaffleModel
Definition: thermalBaffleModel.H:60
Foam::dictionary::writeEntry
void writeEntry(Ostream &os) const
Write sub-dictionary with its dictName as its header.
Definition: dictionaryIO.C:164
Foam::compressible::thermalBaffleFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: thermalBaffleFvPatchScalarField.C:144
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:528
thermalBaffleFvPatchScalarField.H
patchNames
wordList patchNames(nPatches)
compressible
bool compressible
Definition: pEqn.H:2
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.C:349
Foam::polyPatch::New
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:35
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
emptyPolyPatch.H
Foam::compressible::thermalBaffleFvPatchScalarField
This boundary condition provides a coupled temperature condition between multiple mesh regions.
Definition: thermalBaffleFvPatchScalarField.H:310
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::compressible::thermalBaffleFvPatchScalarField::thermalBaffleFvPatchScalarField
thermalBaffleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: thermalBaffleFvPatchScalarField.C:44
Foam::compressible::turbulentTemperatureRadCoupledMixedFvPatchScalarField
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
Definition: turbulentTemperatureRadCoupledMixedFvPatchScalarField.H:141
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::compressible::thermalBaffleFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: thermalBaffleFvPatchScalarField.C:153
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::compressible::thermalBaffleFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: thermalBaffleFvPatchScalarField.C:255
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::fvMesh::name
const word & name() const
Return reference to name.
Definition: fvMesh.H:295
Foam::compressible::thermalBaffleFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: thermalBaffleFvPatchScalarField.C:271