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-2013 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  (
219  patchTypes[patchI],
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  {
241  WarningIn
242  (
243  "thermalBaffleFvPatchScalarField::createPatchMesh()\n"
244  ) << "Specified IOobject::MUST_READ_IF_MODIFIED but class"
245  << " patchMeshPtr not set."
246  << endl;
247  }
248 }
249 
250 
252 {
253  if (this->updated())
254  {
255  return;
256  }
257 
258  const fvMesh& thisMesh = patch().boundaryMesh().mesh();
259 
260  if (owner_ && thisMesh.name() == polyMesh::defaultRegion)
261  {
262  baffle_->evolve();
263  }
264 
266 }
267 
268 
270 {
272 
273  const fvMesh& thisMesh = patch().boundaryMesh().mesh();
274 
275  if (thisMesh.name() == polyMesh::defaultRegion && owner_)
276  {
277 
278  os.writeKeyword("extrudeModel");
279  os << word(dict_.lookup("extrudeModel"))
280  << token::END_STATEMENT << nl;
281 
282  os.writeKeyword("nLayers");
283  os << readLabel(dict_.lookup("nLayers"))
284  << token::END_STATEMENT << nl;
285 
286  os.writeKeyword("expansionRatio");
287  os << readScalar(dict_.lookup("expansionRatio"))
288  << token::END_STATEMENT << nl;
289 
290  os.writeKeyword("columnCells");
291  os << readBool(dict_.lookup("columnCells"))
292  << token::END_STATEMENT << nl;
293 
294  word extrudeModel(word(dict_.lookup("extrudeModel")) + "Coeffs");
295  os.writeKeyword(extrudeModel);
296  os << dict_.subDict(extrudeModel) << nl;
297 
298  word regionName = dict_.lookup("regionName");
299  os.writeKeyword("regionName") << regionName
300  << token::END_STATEMENT << nl;
301 
302  bool active = readBool(dict_.lookup("active"));
303  os.writeKeyword("active") << active
304  << token::END_STATEMENT << nl;
305 
306  os.writeKeyword("thermoType");
307  os << dict_.subDict("thermoType") << nl;
308 
309  os.writeKeyword("mixture");
310  os << dict_.subDict("mixture") << nl;
311 
312  os.writeKeyword("radiation");
313  os << dict_.subDict("radiation") << nl;
314  }
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
321 (
324 );
325 
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
329 } // End namespace compressible
330 } // End namespace Foam
331 
332 
333 // ************************************************************************* //
Top level extrusion model class.
Definition: extrudeModel.H:51
bool readBool(Istream &)
Definition: boolIO.C:60
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const word & coupleGroup() const
PatchGroup (only if NEARESTPATCHFACE)
wordList patchNames(nPatches)
wordList patchTypes(nPatches)
bool foundObject(const word &name) const
Is the named Type found?
thermalBaffleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for handling words, derived from string.
Definition: word.H:59
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Foam::fvPatchFieldMapper.
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
Namespace for OpenFOAM.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
label readLabel(Istream &is)
Definition: label.H:64
static const dictionary null
Null dictionary.
Definition: dictionary.H:193
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
#define WarningIn(functionName)
Report a warning using Foam::Warning.
static const NamedEnum< sampleMode, 6 > sampleModeNames_
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Macros for easy insertion into run-time selection tables.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
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
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:866
const sampleMode & mode() const
What to sample.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
bool compressible
Definition: pEqn.H:40
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const polyMesh & mesh() const
Return the mesh reference.
Foam::word regionName
Mesh at a patch created on the fly. The following entry should be used on the field boundary dictiona...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
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.