regionModel.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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2022 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 
29 #include "regionModel.H"
30 #include "fvMesh.H"
31 #include "Time.H"
32 #include "mappedWallPolyPatch.H"
34 #include "faceAreaWeightAMI.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace regionModels
41 {
42  defineTypeNameAndDebug(regionModel, 0);
43 }
44 }
45 
46 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
47 
48 void Foam::regionModels::regionModel::constructMeshObjects()
49 {
50  if (!time_.foundObject<fvMesh>(regionName_))
51  {
52  fvMesh* regionMeshPtr =
53  new fvMesh
54  (
55  IOobject
56  (
58  time_.timeName(),
59  time_,
61  )
62  );
63 
64  regionMeshPtr->objectRegistry::store();
65  }
66 }
67 
68 
69 void Foam::regionModels::regionModel::initialise()
70 {
71  if (debug)
72  {
73  Pout<< "regionModel::initialise()" << endl;
74  }
75 
76  label nBoundaryFaces = 0;
77  DynamicList<label> primaryPatchIDs;
78  DynamicList<label> intCoupledPatchIDs;
79  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
80 
81  forAll(rbm, patchi)
82  {
83  const polyPatch& regionPatch = rbm[patchi];
84  if (isA<mappedPatchBase>(regionPatch))
85  {
86  if (debug)
87  {
88  Pout<< "found " << mappedWallPolyPatch::typeName
89  << " " << regionPatch.name() << endl;
90  }
91 
92  intCoupledPatchIDs.append(patchi);
93 
94  nBoundaryFaces += regionPatch.faceCells().size();
95 
96  const mappedPatchBase& mapPatch =
97  refCast<const mappedPatchBase>(regionPatch);
98 
99  if
100  (
101  primaryMesh_.time().foundObject<polyMesh>
102  (
103  mapPatch.sampleRegion()
104  )
105  )
106  {
107 
108  const label primaryPatchi = mapPatch.samplePolyPatch().index();
109  primaryPatchIDs.append(primaryPatchi);
110  }
111  }
112  }
113 
114  primaryPatchIDs_.transfer(primaryPatchIDs);
115  intCoupledPatchIDs_.transfer(intCoupledPatchIDs);
116 
117  if (!returnReduceOr(nBoundaryFaces))
118  {
120  << "Region model has no mapped boundary conditions - transfer "
121  << "between regions will not be possible" << endl;
122  }
123 
124  if (!outputPropertiesPtr_)
125  {
126  const fileName uniformPath(word("uniform")/"regionModels");
127 
128  outputPropertiesPtr_.reset
129  (
130  new IOdictionary
131  (
132  IOobject
133  (
134  regionName_ + "OutputProperties",
135  time_.timeName(),
136  uniformPath/regionName_,
137  primaryMesh_,
140  )
141  )
142  );
143  }
144 }
145 
146 
147 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
148 
150 {
151  if (regIOobject::read())
152  {
153  if (active_)
154  {
155  if (const dictionary* dictptr = findDict(modelName_ + "Coeffs"))
156  {
157  coeffs_ <<= *dictptr;
158  }
159 
160  infoOutput_.readIfPresent("infoOutput", *this);
161  }
162 
163  return true;
164  }
165 
166  return false;
167 }
168 
169 
171 {
172  if (active_)
173  {
174  if (const dictionary* dictptr = dict.findDict(modelName_ + "Coeffs"))
175  {
176  coeffs_ <<= *dictptr;
177  }
178 
179  infoOutput_.readIfPresent("infoOutput", dict);
180  return true;
181  }
183  return false;
184 }
185 
186 
189 (
190  const regionModel& nbrRegion,
191  const label regionPatchi,
192  const label nbrPatchi,
193  const bool flip
194 ) const
195 {
196  label nbrRegionID = interRegionAMINames_.find(nbrRegion.name());
197 
198  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
199 
200  if (nbrRegionID != -1)
201  {
202  if (!interRegionAMI_[nbrRegionID].set(regionPatchi))
203  {
204  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
205  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
206 
207  int oldTag = UPstream::msgType();
208  UPstream::msgType() = oldTag + 1;
209 
210  interRegionAMI_[nbrRegionID].set
211  (
212  regionPatchi,
214  (
215  faceAreaWeightAMI::typeName,
216  true, // requireMatch
217  flip
218  )
219  );
220 
221  interRegionAMI_[nbrRegionID][regionPatchi].calculate(p, nbrP);
222 
223  UPstream::msgType() = oldTag;
224  }
225 
226  return interRegionAMI_[nbrRegionID][regionPatchi];
227  }
228  else
229  {
230  label nbrRegionID = interRegionAMINames_.size();
231 
232  interRegionAMINames_.append(nbrRegion.name());
233 
234  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
235  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
236 
237  label nPatch = regionMesh().boundaryMesh().size();
238 
239 
240  interRegionAMI_.resize(nbrRegionID + 1);
241 
242  interRegionAMI_.set
243  (
244  nbrRegionID,
246  );
247 
248  int oldTag = UPstream::msgType();
249  UPstream::msgType() = oldTag + 1;
250 
251  interRegionAMI_[nbrRegionID].set
252  (
253  regionPatchi,
255  (
256  faceAreaWeightAMI::typeName,
257  true, // requireMatch
258  flip // reverse
259  )
260  );
261 
262  interRegionAMI_[nbrRegionID][regionPatchi].calculate(p, nbrP);
263 
264  UPstream::msgType() = oldTag;
266  return interRegionAMI_[nbrRegionID][regionPatchi];
267  }
268 }
269 
270 
272 (
273  const regionModel& nbrRegion,
274  const label regionPatchi
275 ) const
276 {
277  label nbrPatchi = -1;
278 
279  // region
280  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
281 
282  // boundary mesh
283  const polyBoundaryMesh& nbrPbm = nbrRegionMesh.boundaryMesh();
284 
285  const polyBoundaryMesh& pbm = regionMesh().boundaryMesh();
286 
287  if (regionPatchi > pbm.size() - 1)
288  {
290  << "region patch index out of bounds: "
291  << "region patch index = " << regionPatchi
292  << ", maximum index = " << pbm.size() - 1
293  << abort(FatalError);
294  }
295 
296  const polyPatch& pp = regionMesh().boundaryMesh()[regionPatchi];
297 
298  if (!isA<mappedPatchBase>(pp))
299  {
301  << "Expected a " << mappedPatchBase::typeName
302  << " patch, but found a " << pp.type() << abort(FatalError);
303  }
304 
305  const mappedPatchBase& mpb = refCast<const mappedPatchBase>(pp);
306 
307  // sample patch name on the primary region
308  const word& primaryPatchName = mpb.samplePatch();
309 
310  // find patch on nbr region that has the same sample patch name
311  forAll(nbrRegion.intCoupledPatchIDs(), j)
312  {
313  const label nbrRegionPatchi = nbrRegion.intCoupledPatchIDs()[j];
314 
315  const mappedPatchBase& mpb =
316  refCast<const mappedPatchBase>(nbrPbm[nbrRegionPatchi]);
317 
318  if (mpb.samplePatch() == primaryPatchName)
319  {
320  nbrPatchi = nbrRegionPatchi;
321  break;
322  }
323  }
324 
325  if (nbrPatchi == -1)
326  {
327  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
328 
330  << "Unable to find patch pair for local patch "
331  << p.name() << " and region " << nbrRegion.name()
332  << abort(FatalError);
333  }
334 
335  return nbrPatchi;
336 }
337 
338 
339 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
340 
341 Foam::regionModels::regionModel::regionModel
342 (
343  const fvMesh& mesh,
344  const word& regionType
345 )
346 :
348  (
349  IOobject
350  (
351  regionType + "Properties",
352  mesh.time().constant(),
353  mesh.time(),
354  IOobject::NO_READ,
355  IOobject::NO_WRITE
356  )
357  ),
358  primaryMesh_(mesh),
359  time_(mesh.time()),
360  active_(false),
361  infoOutput_(false),
362  modelName_("none"),
363  coeffs_(),
364  outputPropertiesPtr_(nullptr),
365  primaryPatchIDs_(),
366  intCoupledPatchIDs_(),
367  regionName_("none"),
368  functions_(*this),
369  interRegionAMINames_(),
370  interRegionAMI_()
371 {}
372 
373 
374 Foam::regionModels::regionModel::regionModel
375 (
376  const fvMesh& mesh,
377  const word& regionType,
378  const word& modelName,
379  bool readFields
380 )
381 :
383  (
384  IOobject
385  (
386  regionType + "Properties",
387  mesh.time().constant(),
388  mesh.time(),
389  IOobject::MUST_READ,
390  IOobject::NO_WRITE
391  )
392  ),
393  primaryMesh_(mesh),
394  time_(mesh.time()),
395  active_(get<Switch>("active")),
396  infoOutput_(true),
397  modelName_(modelName),
398  coeffs_(subOrEmptyDict(modelName + "Coeffs")),
399  outputPropertiesPtr_(nullptr),
400  primaryPatchIDs_(),
401  intCoupledPatchIDs_(),
402  regionName_(lookup("region")),
403  functions_(*this, subOrEmptyDict("functions"))
404 {
405  if (active_)
406  {
407  constructMeshObjects();
408  initialise();
409 
410  if (readFields)
411  {
412  read();
413  }
414  }
415 }
416 
417 
418 Foam::regionModels::regionModel::regionModel
419 (
420  const fvMesh& mesh,
421  const word& regionType,
422  const word& modelName,
423  const dictionary& dict,
424  bool readFields
425 )
426 :
427  IOdictionary
428  (
429  IOobject
430  (
431  regionType + "Properties",
432  mesh.time().constant(),
433  mesh.time(),
434  IOobject::NO_READ,
435  IOobject::NO_WRITE,
436  true
437  ),
438  dict
439  ),
440  primaryMesh_(mesh),
441  time_(mesh.time()),
442  active_(dict.get<Switch>("active")),
443  infoOutput_(false),
444  modelName_(modelName),
445  coeffs_(dict.subOrEmptyDict(modelName + "Coeffs")),
446  outputPropertiesPtr_(nullptr),
447  primaryPatchIDs_(),
448  intCoupledPatchIDs_(),
449  regionName_(dict.lookup("region")),
450  functions_(*this, subOrEmptyDict("functions"))
451 {
452  if (active_)
453  {
454  constructMeshObjects();
455  initialise();
456 
457  if (readFields)
458  {
459  read(dict);
460  }
461  }
462 }
463 
464 
465 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
466 
468 {
469  if (active_)
470  {
471  Info<< "\nEvolving " << modelName_ << " for region "
472  << regionMesh().name() << endl;
473 
474  preEvolveRegion();
475 
476  evolveRegion();
477 
478  postEvolveRegion();
479 
480  // Provide some feedback
481  if (infoOutput_)
482  {
483  Info<< incrIndent;
484  info();
485  Info<< endl << decrIndent;
486  }
487 
488  if (time_.writeTime())
489  {
490  outputProperties().writeObject
491  (
492  IOstreamOption(IOstreamOption::ASCII, time_.writeCompression()),
493  true
494  );
495  }
496  }
497 }
498 
501 {
502  functions_.preEvolveRegion();
503 }
504 
505 
507 {}
508 
511 {
512  functions_.postEvolveRegion();
513 }
514 
515 
517 {}
518 
519 
520 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
virtual const fileName & name() const
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
dictionary dict
virtual void evolveRegion()
Evolve the region.
Definition: regionModel.C:499
virtual void info()
Provide some feedback.
Definition: regionModel.C:509
virtual bool read()
Read object.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
"ascii" (normal default)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:806
Ignore writing from objectRegistry::writeObject()
Lookup type of boundary radiation properties.
Definition: lookup.H:57
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
dynamicFvMesh & mesh
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:493
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel.C:142
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
Definition: regionModel.C:265
virtual void postEvolveRegion()
Post-evolve region.
Definition: regionModel.C:503
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:99
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
virtual const AMIPatchToPatchInterpolation & interRegionAMI(const regionModel &nbrRegion, const label regionPatchi, const label nbrPatchi, const bool flip) const
Create or return a new inter-region AMI object.
Definition: regionModel.C:182
errorManip< error > abort(error &err)
Definition: errorManip.H:139
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
Switch active_
Active flag.
Definition: regionModel.H:102
int debug
Static debugging option.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:467
#define WarningInFunction
Report a warning using Foam::Warning.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
const Time & time_
Reference to the time database.
Definition: regionModel.H:97
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
IOobject(const IOobject &)=default
Copy construct.
virtual void evolve()
Main driver routing to evolve the region - calls other evolves.
Definition: regionModel.C:460
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
messageStream Info
Information stream (stdout output on master, null elsewhere)
Base class for region models.
Definition: regionModel.H:56
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type.
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
word regionName_
Region name.
Definition: regionModel.H:141
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:458
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
defineTypeNameAndDebug(KirchhoffShell, 0)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.