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  fvMesh* regionMeshPtr = time_.getObjectPtr<fvMesh>(regionName_);
51 
52  if (!regionMeshPtr)
53  {
54  regionMeshPtr = new fvMesh
55  (
56  IOobject
57  (
59  time_.timeName(),
60  time_,
64  )
65  );
66 
67  regionMeshPtr->objectRegistry::store();
68  }
69 }
70 
71 
72 void Foam::regionModels::regionModel::initialise()
73 {
74  if (debug)
75  {
76  Pout<< "regionModel::initialise()" << endl;
77  }
78 
79  label nBoundaryFaces = 0;
80  DynamicList<label> primaryPatchIDs;
81  DynamicList<label> intCoupledPatchIDs;
82  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
83 
84  forAll(rbm, patchi)
85  {
86  const polyPatch& regionPatch = rbm[patchi];
87  if (isA<mappedPatchBase>(regionPatch))
88  {
89  if (debug)
90  {
91  Pout<< "found " << mappedWallPolyPatch::typeName
92  << " " << regionPatch.name() << endl;
93  }
94 
95  intCoupledPatchIDs.append(patchi);
96 
97  nBoundaryFaces += regionPatch.faceCells().size();
98 
99  const mappedPatchBase& mapPatch =
100  refCast<const mappedPatchBase>(regionPatch);
101 
102  if
103  (
104  primaryMesh_.time().foundObject<polyMesh>
105  (
106  mapPatch.sampleRegion()
107  )
108  )
109  {
110 
111  const label primaryPatchi = mapPatch.samplePolyPatch().index();
112  primaryPatchIDs.append(primaryPatchi);
113  }
114  }
115  }
116 
117  primaryPatchIDs_.transfer(primaryPatchIDs);
118  intCoupledPatchIDs_.transfer(intCoupledPatchIDs);
119 
120  if (!returnReduceOr(nBoundaryFaces))
121  {
123  << "Region model has no mapped boundary conditions - transfer "
124  << "between regions will not be possible" << endl;
125  }
126 
127  if (!outputPropertiesPtr_)
128  {
129  const fileName uniformPath(word("uniform")/"regionModels");
130 
131  outputPropertiesPtr_.reset
132  (
133  new IOdictionary
134  (
135  IOobject
136  (
137  regionName_ + "OutputProperties",
138  time_.timeName(),
139  uniformPath/regionName_,
140  primaryMesh_,
143  )
144  )
145  );
146  }
147 }
148 
149 
150 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
151 
153 {
154  if (regIOobject::read())
155  {
156  if (active_)
157  {
158  if (const dictionary* dictptr = findDict(modelName_ + "Coeffs"))
159  {
160  coeffs_ <<= *dictptr;
161  }
162 
163  infoOutput_.readIfPresent("infoOutput", *this);
164  }
165 
166  return true;
167  }
168 
169  return false;
170 }
171 
172 
174 {
175  if (active_)
176  {
177  if (const dictionary* dictptr = dict.findDict(modelName_ + "Coeffs"))
178  {
179  coeffs_ <<= *dictptr;
180  }
181 
182  infoOutput_.readIfPresent("infoOutput", dict);
183  return true;
184  }
186  return false;
187 }
188 
189 
192 (
193  const regionModel& nbrRegion,
194  const label regionPatchi,
195  const label nbrPatchi,
196  const bool flip
197 ) const
198 {
199  label nbrRegionID = interRegionAMINames_.find(nbrRegion.name());
200 
201  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
202 
203  if (nbrRegionID != -1)
204  {
205  if (!interRegionAMI_[nbrRegionID].set(regionPatchi))
206  {
207  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
208  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
209 
210  const int oldTag = UPstream::incrMsgType();
211 
212  interRegionAMI_[nbrRegionID].set
213  (
214  regionPatchi,
216  (
217  faceAreaWeightAMI::typeName,
218  true, // requireMatch
219  flip
220  )
221  );
222 
223  interRegionAMI_[nbrRegionID][regionPatchi].calculate(p, nbrP);
224 
225  UPstream::msgType(oldTag); // Restore tag
226  }
227 
228  return interRegionAMI_[nbrRegionID][regionPatchi];
229  }
230  else
231  {
232  label nbrRegionID = interRegionAMINames_.size();
233 
234  interRegionAMINames_.append(nbrRegion.name());
235 
236  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
237  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
238 
239  const label nPatch = regionMesh().boundaryMesh().size();
240 
241 
242  interRegionAMI_.resize(nbrRegionID + 1);
243 
244  interRegionAMI_.set
245  (
246  nbrRegionID,
248  );
249 
250  const int oldTag = UPstream::incrMsgType();
251 
252  interRegionAMI_[nbrRegionID].set
253  (
254  regionPatchi,
256  (
257  faceAreaWeightAMI::typeName,
258  true, // requireMatch
259  flip // reverse
260  )
261  );
262 
263  interRegionAMI_[nbrRegionID][regionPatchi].calculate(p, nbrP);
264 
265  UPstream::msgType(oldTag); // Restore tag
267  return interRegionAMI_[nbrRegionID][regionPatchi];
268  }
269 }
270 
271 
273 (
274  const regionModel& nbrRegion,
275  const label regionPatchi
276 ) const
277 {
278  label nbrPatchi = -1;
279 
280  // region
281  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
282 
283  // boundary mesh
284  const polyBoundaryMesh& nbrPbm = nbrRegionMesh.boundaryMesh();
285 
286  const polyBoundaryMesh& pbm = regionMesh().boundaryMesh();
287 
288  if (regionPatchi > pbm.size() - 1)
289  {
291  << "region patch index out of bounds: "
292  << "region patch index = " << regionPatchi
293  << ", maximum index = " << pbm.size() - 1
294  << abort(FatalError);
295  }
296 
297  const polyPatch& pp = regionMesh().boundaryMesh()[regionPatchi];
298 
299  if (!isA<mappedPatchBase>(pp))
300  {
302  << "Expected a " << mappedPatchBase::typeName
303  << " patch, but found a " << pp.type() << abort(FatalError);
304  }
305 
306  const mappedPatchBase& mpb = refCast<const mappedPatchBase>(pp);
307 
308  // sample patch name on the primary region
309  const word& primaryPatchName = mpb.samplePatch();
310 
311  // find patch on nbr region that has the same sample patch name
312  forAll(nbrRegion.intCoupledPatchIDs(), j)
313  {
314  const label nbrRegionPatchi = nbrRegion.intCoupledPatchIDs()[j];
315 
316  const mappedPatchBase& mpb =
317  refCast<const mappedPatchBase>(nbrPbm[nbrRegionPatchi]);
318 
319  if (mpb.samplePatch() == primaryPatchName)
320  {
321  nbrPatchi = nbrRegionPatchi;
322  break;
323  }
324  }
325 
326  if (nbrPatchi == -1)
327  {
328  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
329 
331  << "Unable to find patch pair for local patch "
332  << p.name() << " and region " << nbrRegion.name()
333  << abort(FatalError);
334  }
335 
336  return nbrPatchi;
337 }
338 
339 
340 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
341 
342 Foam::regionModels::regionModel::regionModel
343 (
344  const fvMesh& mesh,
345  const word& regionType
346 )
347 :
349  (
350  IOobject
351  (
352  regionType + "Properties",
353  mesh.time().constant(),
354  mesh.time(),
355  IOobject::NO_READ,
356  IOobject::NO_WRITE
357  )
358  ),
359  primaryMesh_(mesh),
360  time_(mesh.time()),
361  active_(false),
362  infoOutput_(false),
363  modelName_("none"),
364  coeffs_(),
365  outputPropertiesPtr_(nullptr),
366  primaryPatchIDs_(),
367  intCoupledPatchIDs_(),
368  regionName_("none"),
369  functions_(*this),
370  interRegionAMINames_(),
371  interRegionAMI_()
372 {}
373 
374 
375 Foam::regionModels::regionModel::regionModel
376 (
377  const fvMesh& mesh,
378  const word& regionType,
379  const word& modelName,
380  bool readFields
381 )
382 :
384  (
385  IOobject
386  (
387  regionType + "Properties",
388  mesh.time().constant(),
389  mesh.time(),
390  IOobject::MUST_READ,
391  IOobject::NO_WRITE
392  )
393  ),
394  primaryMesh_(mesh),
395  time_(mesh.time()),
396  active_(get<Switch>("active")),
397  infoOutput_(true),
398  modelName_(modelName),
399  coeffs_(subOrEmptyDict(modelName + "Coeffs")),
400  outputPropertiesPtr_(nullptr),
401  primaryPatchIDs_(),
402  intCoupledPatchIDs_(),
403  regionName_(lookup("region")),
404  functions_(*this, subOrEmptyDict("functions"))
405 {
406  if (active_)
407  {
408  constructMeshObjects();
409  initialise();
410 
411  if (readFields)
412  {
413  read();
414  }
415  }
416 }
417 
418 
419 Foam::regionModels::regionModel::regionModel
420 (
421  const fvMesh& mesh,
422  const word& regionType,
423  const word& modelName,
424  const dictionary& dict,
425  bool readFields
426 )
427 :
428  IOdictionary
429  (
430  IOobject
431  (
432  regionType + "Properties",
433  mesh.time().constant(),
434  mesh.time(),
435  IOobject::NO_READ,
436  IOobject::NO_WRITE,
437  IOobject::REGISTER
438  ),
439  dict
440  ),
441  primaryMesh_(mesh),
442  time_(mesh.time()),
443  active_(dict.get<Switch>("active")),
444  infoOutput_(false),
445  modelName_(modelName),
446  coeffs_(dict.subOrEmptyDict(modelName + "Coeffs")),
447  outputPropertiesPtr_(nullptr),
448  primaryPatchIDs_(),
449  intCoupledPatchIDs_(),
450  regionName_(dict.lookup("region")),
451  functions_(*this, subOrEmptyDict("functions"))
452 {
453  if (active_)
454  {
455  constructMeshObjects();
456  initialise();
457 
458  if (readFields)
459  {
460  read(dict);
461  }
462  }
463 }
464 
465 
466 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
467 
469 {
470  if (active_)
471  {
472  Info<< "\nEvolving " << modelName_ << " for region "
473  << regionMesh().name() << endl;
474 
475  preEvolveRegion();
476 
477  evolveRegion();
478 
479  postEvolveRegion();
480 
481  // Provide some feedback
482  if (infoOutput_)
483  {
484  Info<< incrIndent;
485  info();
486  Info<< endl << decrIndent;
487  }
488 
489  if (time_.writeTime())
490  {
491  outputProperties().writeObject
492  (
493  IOstreamOption(IOstreamOption::ASCII, time_.writeCompression()),
494  true
495  );
496  }
497  }
498 }
499 
502 {
503  functions_.preEvolveRegion();
504 }
505 
506 
508 {}
509 
512 {
513  functions_.postEvolveRegion();
514 }
515 
516 
518 {}
519 
520 
521 // ************************************************************************* //
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 & pbm
dictionary dict
virtual void evolveRegion()
Evolve the region.
Definition: regionModel.C:500
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition: UPstream.H:1251
virtual void info()
Provide some feedback.
Definition: regionModel.C:510
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:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
"ascii" (normal default)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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:1229
Ignore writing from objectRegistry::writeObject()
Lookup type of boundary radiation properties.
Definition: lookup.H:57
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
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:421
dynamicFvMesh & mesh
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:494
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel.C:145
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
Definition: regionModel.C:266
virtual void postEvolveRegion()
Post-evolve region.
Definition: regionModel.C:504
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
Reading is optional [identical to LAZY_READ].
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:185
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 a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
Switch active_
Active flag.
Definition: regionModel.H:102
int debug
Static debugging option.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:511
#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:78
virtual void evolve()
Main driver routing to evolve the region - calls other evolves.
Definition: regionModel.C:461
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
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
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:502
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
defineTypeNameAndDebug(KirchhoffShell, 0)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.