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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "emptyPolyPatch.H"
29 #include "polyPatch.H"
30 #include "mappedWallPolyPatch.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace compressible
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
43 (
44  const fvPatch& p,
46 )
47 :
49  owner_(false),
50  baffle_(),
51  dict_(dictionary::null),
52  extrudeMeshPtr_()
53 {}
54 
55 
58 (
60  const fvPatch& p,
62  const fvPatchFieldMapper& mapper
63 )
64 :
66  (
67  ptf,
68  p,
69  iF,
70  mapper
71  ),
72  owner_(ptf.owner_),
73  baffle_(),
74  dict_(ptf.dict_),
75  extrudeMeshPtr_()
76 {}
77 
78 
81 (
82  const fvPatch& p,
84  const dictionary& dict
85 )
86 :
88  owner_(false),
89  baffle_(),
90  dict_(dict),
91  extrudeMeshPtr_()
92 {
93 
94  const fvMesh& thisMesh = patch().boundaryMesh().mesh();
95 
97 
98  if (thisMesh.name() == polyMesh::defaultRegion)
99  {
100  const word regionName =
101  dict_.lookupOrDefault<word>("regionName", "none");
102 
103  const word baffleName("3DBaffle" + regionName);
104 
105  if
106  (
107  !thisMesh.time().foundObject<fvMesh>(regionName)
108  && regionName != "none"
109  )
110  {
111  if (extrudeMeshPtr_.empty())
112  {
113  createPatchMesh();
114  }
115 
116  baffle_.reset(baffle::New(thisMesh, dict).ptr());
117  owner_ = true;
118  baffle_->rename(baffleName);
119  }
120  }
121 }
122 
123 
126 (
129 )
130 :
132  owner_(ptf.owner_),
133  baffle_(),
134  dict_(ptf.dict_),
135  extrudeMeshPtr_()
136 {}
137 
138 
139 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140 
141 
143 (
144  const fvPatchFieldMapper& m
145 )
146 {
147  mixedFvPatchScalarField::autoMap(m);
148 }
149 
150 
152 (
153  const fvPatchScalarField& ptf,
154  const labelList& addr
155 )
156 {
157  mixedFvPatchScalarField::rmap(ptf, addr);
158 }
159 
160 
161 void thermalBaffleFvPatchScalarField::createPatchMesh()
162 {
163 
164  const fvMesh& thisMesh = patch().boundaryMesh().mesh();
165 
166  word regionName = dict_.lookup("regionName");
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  patchTypes[topPatchID] = mappedWallPolyPatch::typeName;
179 
180  if (readBool(dict_.lookup("columnCells")))
181  {
182  patchTypes[sidePatchID] = emptyPolyPatch::typeName;
183  }
184  else
185  {
186  patchTypes[sidePatchID] = polyPatch::typeName;
187  }
188 
189  const mappedPatchBase& mpp =
190  refCast<const mappedPatchBase>(patch().patch());
191 
192  const word coupleGroup(mpp.coupleGroup());
193 
194  wordList inGroups(1);
195  inGroups[0] = coupleGroup;
196 
197  dicts[bottomPatchID].add("coupleGroup", coupleGroup);
198  dicts[bottomPatchID].add("inGroups", inGroups);
199  dicts[bottomPatchID].add("sampleMode", mpp.sampleModeNames_[mpp.mode()]);
200 
201  const label sepPos = coupleGroup.find('_');
202 
203  const word coupleGroupSlave = coupleGroup(0, sepPos) + "_slave";
204 
205  inGroups[0] = coupleGroupSlave;
206  dicts[topPatchID].add("coupleGroup", coupleGroupSlave);
207  dicts[topPatchID].add("inGroups", inGroups);
208  dicts[topPatchID].add("sampleMode", mpp.sampleModeNames_[mpp.mode()]);
209 
210 
211  forAll(regionPatches, patchi)
212  {
213  dictionary& patchDict = dicts[patchi];
214  patchDict.set("nFaces", 0);
215  patchDict.set("startFace", 0);
216 
217  regionPatches[patchi] = polyPatch::New
218  (
220  patchNames[patchi],
221  dicts[patchi],
222  patchi,
223  thisMesh.boundaryMesh()
224  ).ptr();
225  }
226 
227  extrudeMeshPtr_.reset
228  (
229  new extrudePatchMesh
230  (
231  thisMesh,
232  patch(),
233  dict_,
234  regionName,
235  regionPatches
236  )
237  );
238 
239  if (extrudeMeshPtr_.empty())
240  {
242  << "Specified IOobject::MUST_READ_IF_MODIFIED but class"
243  << " patchMeshPtr not set."
244  << endl;
245  }
246 }
247 
248 
250 {
251  if (this->updated())
252  {
253  return;
254  }
255 
256  const fvMesh& thisMesh = patch().boundaryMesh().mesh();
257 
258  if (owner_ && thisMesh.name() == polyMesh::defaultRegion)
259  {
260  baffle_->evolve();
261  }
262 
264 }
265 
266 
268 {
270 
271  const fvMesh& thisMesh = patch().boundaryMesh().mesh();
272 
273  if (thisMesh.name() == polyMesh::defaultRegion && owner_)
274  {
275 
276  os.writeKeyword("extrudeModel");
277  os << word(dict_.lookup("extrudeModel"))
278  << token::END_STATEMENT << nl;
279 
280  os.writeKeyword("nLayers");
281  os << readLabel(dict_.lookup("nLayers"))
282  << token::END_STATEMENT << nl;
283 
284  os.writeKeyword("expansionRatio");
285  os << readScalar(dict_.lookup("expansionRatio"))
286  << token::END_STATEMENT << nl;
287 
288  os.writeKeyword("columnCells");
289  os << readBool(dict_.lookup("columnCells"))
290  << token::END_STATEMENT << nl;
291 
292  word extrudeModel(word(dict_.lookup("extrudeModel")) + "Coeffs");
293  os.writeKeyword(extrudeModel);
294  os << dict_.subDict(extrudeModel) << nl;
295 
296  word regionName = dict_.lookup("regionName");
297  os.writeKeyword("regionName") << regionName
298  << token::END_STATEMENT << nl;
299 
300  bool active = readBool(dict_.lookup("active"));
301  os.writeKeyword("active") << active
302  << token::END_STATEMENT << nl;
303 
304  os.writeKeyword("thermoType");
305  os << dict_.subDict("thermoType") << nl;
306 
307  os.writeKeyword("mixture");
308  os << dict_.subDict("mixture") << nl;
309 
310  os.writeKeyword("radiation");
311  os << dict_.subDict("radiation") << nl;
312  }
313 }
314 
315 
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317 
319 (
322 );
323 
324 
325 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326 
327 } // End namespace compressible
328 } // End namespace Foam
329 
330 
331 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
bool foundObject(const word &name) const
Is the named Type found?
Foam::word regionName
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
static const dictionary null
Null dictionary.
Definition: dictionary.H:193
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
wordList patchTypes(nPatches)
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Mesh at a patch created on the fly. The following entry should be used on the field boundary dictiona...
bool readBool(Istream &)
Definition: boolIO.C:60
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Macros for easy insertion into run-time selection tables.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Top level extrusion model class.
Definition: extrudeModel.H:51
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
wordList patchNames(nPatches)
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
const polyMesh & mesh() const
Return the mesh reference.
label readLabel(Istream &is)
Definition: label.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
static const NamedEnum< sampleMode, 6 > sampleModeNames_
thermalBaffleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:864
const word & coupleGroup() const
PatchGroup (only if NEARESTPATCHFACE)
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:32
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
const sampleMode & mode() const
What to sample.
This boundary condition provides a coupled temperature condition between multiple mesh regions...
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
bool compressible
Definition: pEqn.H:30