thermalBaffleFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 
27 #include "mappedWallPolyPatch.H"
28 #include "symmetryPolyPatch.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace compressible
36 {
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * /
39 
40 bool thermalBaffleFvPatchScalarField::primary() const
41 {
42  return patch().boundaryMesh().mesh().name() == polyMesh::defaultRegion;
43 }
44 
45 
46 bool thermalBaffleFvPatchScalarField::owner() const
47 {
48  return
49  primary()
50  && patch().index() < patch().boundaryMesh()[nbrPatch_].index();
51 }
52 
53 
54 void thermalBaffleFvPatchScalarField::checkPatches() const
55 {
56  if (!primary()) return;
57 
58  const polyPatch& pp = patch().patch();
59  const polyPatch& nbrPp = patch().patch().boundaryMesh()[nbrPatch_];
60 
61  // The patches should be of mapped type
62  auto checkPatchIsMapped = [&](const polyPatch& pp)
63  {
64  if (!isA<mappedPatchBase>(pp))
65  {
67  << "Patch field of type \"" << typeName
68  << "\" specified for patch \"" << pp.name() << "\" of field \""
69  << internalField().name() << "\", but this patch is not of "
70  << "type \"" << mappedPatchBase::typeName << "\""
71  << exit(FatalError);
72  }
73  };
74  checkPatchIsMapped(pp);
75  checkPatchIsMapped(nbrPp);
76 
77  const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp);
78  const mappedPatchBase nbrMpp = refCast<const mappedPatchBase>(nbrPp);
79 
80  // The patches should sample a different region
81  auto checkPatchMapsDifferentRegion = [&](const mappedPatchBase& mpp)
82  {
83  if (mpp.sameRegion())
84  {
86  << "Patch field of type \"" << typeName
87  << "\" specified for patch \"" << pp.name() << "\" of field \""
88  << internalField().name() << "\", but this patch maps to "
89  << "another patch in the same region. It should map to a "
90  << "different region; i.e., that of the thermal baffle model."
91  << exit(FatalError);
92  }
93  };
94  checkPatchMapsDifferentRegion(mpp);
95  checkPatchMapsDifferentRegion(nbrMpp);
96 
97  // The sample region of this patch and it's neighbour should be the same,
98  // i.e., that of the thermal baffle model
99  if (mpp.sampleRegion() != nbrMpp.sampleRegion())
100  {
102  << "Patch fields of type \"" << typeName
103  << "\" specified for patches \"" << pp.name() << "\" and \""
104  << nbrPp.name() << "\" of field \"" << internalField().name()
105  << "\", but these patches map to different regions \""
106  << mpp.sampleRegion() << "\" and \"" << nbrMpp.sampleRegion()
107  << ". They should map to the same region; i.e., that of the "
108  << "thermal baffle model."
109  << exit(FatalError);
110  }
111 
112  // The sample patch of this patch and it's neighbour should be different,
113  // i.e., they should sample opposite ends of the thermal baffle mesh
114  if (mpp.samplePatch() == nbrMpp.samplePatch())
115  {
117  << "Patch fields of type \"" << typeName
118  << "\" specified for patches \"" << pp.name() << "\" and \""
119  << nbrPp.name() << "\" of field \"" << internalField().name()
120  << "\", but these patches map to the same patch; \""
121  << mpp.samplePatch() << "\" of region \"" << mpp.sampleRegion()
122  << ". They should map to different patches, as these will become "
123  << "the patches at opposite ends of the extruded baffle mesh."
124  << exit(FatalError);
125  }
126 }
127 
128 
129 void thermalBaffleFvPatchScalarField::checkPatchFields() const
130 {
131  if (!primary()) return;
132 
133  const fvPatch& fvp = patch();
134  const fvPatch& nbrFvp = patch().boundaryMesh()[nbrPatch_];
135 
136  const fvPatchScalarField& nbrTp =
137  nbrFvp.lookupPatchField<volScalarField, scalar>(internalField().name());
138 
139  // The neighbour patch field should be of the same type
140  if (!isA<thermalBaffleFvPatchScalarField>(nbrTp))
141  {
143  << "Patch field of type \"" << typeName
144  << "\" specified for patch \"" << fvp.name() << "\" of field \""
145  << internalField().name() << "\" but the field on the "
146  << "neighbouring patch \"" << nbrFvp.name()
147  << "\" is of a different type. Both should be of type \""
148  << typeName << "\"."
149  << exit(FatalError);
150  }
151 
152  // The neighbour patch field's neighbour should be this patch
153  const thermalBaffleFvPatchScalarField& nbrTBp =
154  refCast<const thermalBaffleFvPatchScalarField>(nbrTp);
155  if (nbrTBp.nbrPatch_ != patch().name())
156  {
158  << "Patch field of type \"" << typeName
159  << "\" on patch \"" << fvp.name() << "\" of field \""
160  << internalField().name() << "\" is specified to neighbour "
161  << "patch \"" << nbrPatch_ << "\", but this patch does not "
162  << "reciprocally neighbour patch \"" << fvp.name() << "\""
163  << exit(FatalError);
164  }
165 }
166 
167 
168 autoPtr<extrudePatchMesh>
169 thermalBaffleFvPatchScalarField::initBaffleMesh() const
170 {
171  if (!owner())
172  {
174  << "Baffle mesh is only available to the owner patch in the "
175  << "primary region" << exit(FatalError);
176  }
177 
178  checkPatches();
179 
180  const fvMesh& mesh = patch().boundaryMesh().mesh();
181 
182  const mappedPatchBase& mpp =
183  refCast<const mappedPatchBase>(patch().patch());
184 
185  const mappedPatchBase nbrMpp =
186  refCast<const mappedPatchBase>
187  (patch().patch().boundaryMesh()[nbrPatch_]);
188 
189  const List<word> patchNames
190  ({
191  mpp.samplePatch(),
192  nbrMpp.samplePatch(),
193  "sides"
194  });
195 
196  const List<word> patchTypes
197  ({
198  mappedWallPolyPatch::typeName,
199  mappedWallPolyPatch::typeName,
200  symmetryPolyPatch::typeName
201  });
202 
203  List<dictionary> patchDicts(3);
205  {
206  patchDicts[patchi].set("nFaces", 0);
207  patchDicts[patchi].set("startFace", 0);
208  }
209  patchDicts[0].add("sampleMode", mpp.sampleModeNames_[mpp.mode()]);
210  patchDicts[0].add("sampleRegion", mesh.name());
211  patchDicts[0].add("samplePatch", patch().name());
212  patchDicts[1].add("sampleMode", mpp.sampleModeNames_[nbrMpp.mode()]);
213  patchDicts[1].add("sampleRegion", mesh.name());
214  patchDicts[1].add("samplePatch", nbrPatch_);
215 
216  List<polyPatch*> patchPtrs(3);
217  forAll(patchPtrs, patchi)
218  {
219  patchPtrs[patchi] = polyPatch::New
220  (
222  patchNames[patchi],
223  patchDicts[patchi],
224  patchi,
225  mesh.boundaryMesh()
226  ).ptr();
227  }
228 
229  dictionary dict(dict_);
230  dict.add("columnCells", false);
231 
232  return
233  autoPtr<extrudePatchMesh>
234  (
235  new extrudePatchMesh
236  (
237  mesh,
238  patch(),
239  dict,
240  mpp.sampleRegion(),
241  patchPtrs
242  )
243  );
244 }
245 
246 
247 autoPtr<regionModels::thermalBaffleModel>
248 thermalBaffleFvPatchScalarField::initBaffle() const
249 {
250  if (!owner())
251  {
253  << "Baffle model is only available to the owner patch in the "
254  << "primary region" << exit(FatalError);
255  }
256 
257  checkPatches();
258 
259  const fvMesh& mesh = patch().boundaryMesh().mesh();
260 
261  const mappedPatchBase& mpp =
262  refCast<const mappedPatchBase>(patch().patch());
263 
264  dictionary dict(dict_);
265  dict.add("regionName", mpp.sampleRegion());
266 
267  return regionModels::thermalBaffleModel::New(mesh, dict);
268 }
269 
270 
271 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
272 
274 (
275  const fvPatch& p,
277 )
278 :
280  dict_(dictionary::null),
281  nbrPatch_(word::null),
282  baffleMeshPtr_(),
283  bafflePtr_()
284 {}
285 
286 
288 (
289  const fvPatch& p,
291  const dictionary& dict
292 )
293 :
295  dict_(dict),
296  nbrPatch_(primary() ? dict.lookup<word>("neighbourPatch") : word::null),
297  baffleMeshPtr_(owner() ? initBaffleMesh().ptr() : nullptr),
298  bafflePtr_(owner() ? initBaffle().ptr() : nullptr)
299 {}
300 
301 
303 (
305  const fvPatch& p,
307  const fvPatchFieldMapper& mapper
308 )
309 :
311  (
312  ptf,
313  p,
314  iF,
315  mapper
316  ),
317  dict_(ptf.dict_),
318  baffleMeshPtr_(),
319  bafflePtr_()
320 {}
321 
322 
324 (
327 )
328 :
330  dict_(ptf.dict_),
331  baffleMeshPtr_(),
332  bafflePtr_()
333 {}
334 
335 
336 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
337 
339 (
340  const fvPatchFieldMapper& m
341 )
342 {
343  turbulentTemperatureRadCoupledMixedFvPatchScalarField::autoMap(m);
344 }
345 
346 
348 (
349  const fvPatchScalarField& ptf,
350  const labelList& addr
351 )
352 {
353  turbulentTemperatureRadCoupledMixedFvPatchScalarField::rmap(ptf, addr);
354 }
355 
356 
358 (
359  const fvPatchScalarField& ptf
360 )
361 {
362  turbulentTemperatureRadCoupledMixedFvPatchScalarField::reset(ptf);
363 }
364 
365 
367 {
368  if (this->updated())
369  {
370  return;
371  }
372 
373  checkPatchFields();
374 
375  if (owner())
376  {
377  bafflePtr_->evolve();
378  }
379 
381 }
382 
383 
385 {
387 
388  if (owner())
389  {
390  forAllConstIter(dictionary, dict_, iter)
391  {
392  os << *iter;
393  }
394  }
395  else if (primary())
396  {
397  writeEntry(os, "neighbourPatch", nbrPatch_);
398  }
399 }
400 
401 
402 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 
405 (
408 );
409 
410 
411 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
412 
413 } // End namespace compressible
414 } // End namespace Foam
415 
416 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
wordList patchTypes(nPatches)
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
static autoPtr< thermalBaffleModel > New(const fvMesh &mesh)
Return a reference to the selected model.
Macros for easy insertion into run-time selection tables.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1153
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
wordList patchNames(nPatches)
static const dictionary null
Null dictionary.
Definition: dictionary.H:242
turbulentTemperatureRadCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static const word null
An empty word.
Definition: word.H:77
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
thermalBaffleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label patchi
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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.
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:537
This boundary condition provides a coupled temperature condition between multiple mesh regions...
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864