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  int oldTag = UPstream::msgType();
211  UPstream::msgType() = oldTag + 1;
212 
213  interRegionAMI_[nbrRegionID].set
214  (
215  regionPatchi,
217  (
218  faceAreaWeightAMI::typeName,
219  true, // requireMatch
220  flip
221  )
222  );
223 
224  interRegionAMI_[nbrRegionID][regionPatchi].calculate(p, nbrP);
225 
226  UPstream::msgType() = oldTag;
227  }
228 
229  return interRegionAMI_[nbrRegionID][regionPatchi];
230  }
231  else
232  {
233  label nbrRegionID = interRegionAMINames_.size();
234 
235  interRegionAMINames_.append(nbrRegion.name());
236 
237  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
238  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
239 
240  label nPatch = regionMesh().boundaryMesh().size();
241 
242 
243  interRegionAMI_.resize(nbrRegionID + 1);
244 
245  interRegionAMI_.set
246  (
247  nbrRegionID,
249  );
250 
251  int oldTag = UPstream::msgType();
252  UPstream::msgType() = oldTag + 1;
253 
254  interRegionAMI_[nbrRegionID].set
255  (
256  regionPatchi,
258  (
259  faceAreaWeightAMI::typeName,
260  true, // requireMatch
261  flip // reverse
262  )
263  );
264 
265  interRegionAMI_[nbrRegionID][regionPatchi].calculate(p, nbrP);
266 
267  UPstream::msgType() = oldTag;
269  return interRegionAMI_[nbrRegionID][regionPatchi];
270  }
271 }
272 
273 
275 (
276  const regionModel& nbrRegion,
277  const label regionPatchi
278 ) const
279 {
280  label nbrPatchi = -1;
281 
282  // region
283  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
284 
285  // boundary mesh
286  const polyBoundaryMesh& nbrPbm = nbrRegionMesh.boundaryMesh();
287 
288  const polyBoundaryMesh& pbm = regionMesh().boundaryMesh();
289 
290  if (regionPatchi > pbm.size() - 1)
291  {
293  << "region patch index out of bounds: "
294  << "region patch index = " << regionPatchi
295  << ", maximum index = " << pbm.size() - 1
296  << abort(FatalError);
297  }
298 
299  const polyPatch& pp = regionMesh().boundaryMesh()[regionPatchi];
300 
301  if (!isA<mappedPatchBase>(pp))
302  {
304  << "Expected a " << mappedPatchBase::typeName
305  << " patch, but found a " << pp.type() << abort(FatalError);
306  }
307 
308  const mappedPatchBase& mpb = refCast<const mappedPatchBase>(pp);
309 
310  // sample patch name on the primary region
311  const word& primaryPatchName = mpb.samplePatch();
312 
313  // find patch on nbr region that has the same sample patch name
314  forAll(nbrRegion.intCoupledPatchIDs(), j)
315  {
316  const label nbrRegionPatchi = nbrRegion.intCoupledPatchIDs()[j];
317 
318  const mappedPatchBase& mpb =
319  refCast<const mappedPatchBase>(nbrPbm[nbrRegionPatchi]);
320 
321  if (mpb.samplePatch() == primaryPatchName)
322  {
323  nbrPatchi = nbrRegionPatchi;
324  break;
325  }
326  }
327 
328  if (nbrPatchi == -1)
329  {
330  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
331 
333  << "Unable to find patch pair for local patch "
334  << p.name() << " and region " << nbrRegion.name()
335  << abort(FatalError);
336  }
337 
338  return nbrPatchi;
339 }
340 
341 
342 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
343 
344 Foam::regionModels::regionModel::regionModel
345 (
346  const fvMesh& mesh,
347  const word& regionType
348 )
349 :
351  (
352  IOobject
353  (
354  regionType + "Properties",
355  mesh.time().constant(),
356  mesh.time(),
357  IOobject::NO_READ,
358  IOobject::NO_WRITE
359  )
360  ),
361  primaryMesh_(mesh),
362  time_(mesh.time()),
363  active_(false),
364  infoOutput_(false),
365  modelName_("none"),
366  coeffs_(),
367  outputPropertiesPtr_(nullptr),
368  primaryPatchIDs_(),
369  intCoupledPatchIDs_(),
370  regionName_("none"),
371  functions_(*this),
372  interRegionAMINames_(),
373  interRegionAMI_()
374 {}
375 
376 
377 Foam::regionModels::regionModel::regionModel
378 (
379  const fvMesh& mesh,
380  const word& regionType,
381  const word& modelName,
382  bool readFields
383 )
384 :
386  (
387  IOobject
388  (
389  regionType + "Properties",
390  mesh.time().constant(),
391  mesh.time(),
392  IOobject::MUST_READ,
393  IOobject::NO_WRITE
394  )
395  ),
396  primaryMesh_(mesh),
397  time_(mesh.time()),
398  active_(get<Switch>("active")),
399  infoOutput_(true),
400  modelName_(modelName),
401  coeffs_(subOrEmptyDict(modelName + "Coeffs")),
402  outputPropertiesPtr_(nullptr),
403  primaryPatchIDs_(),
404  intCoupledPatchIDs_(),
405  regionName_(lookup("region")),
406  functions_(*this, subOrEmptyDict("functions"))
407 {
408  if (active_)
409  {
410  constructMeshObjects();
411  initialise();
412 
413  if (readFields)
414  {
415  read();
416  }
417  }
418 }
419 
420 
421 Foam::regionModels::regionModel::regionModel
422 (
423  const fvMesh& mesh,
424  const word& regionType,
425  const word& modelName,
426  const dictionary& dict,
427  bool readFields
428 )
429 :
430  IOdictionary
431  (
432  IOobject
433  (
434  regionType + "Properties",
435  mesh.time().constant(),
436  mesh.time(),
437  IOobject::NO_READ,
438  IOobject::NO_WRITE,
439  IOobject::REGISTER
440  ),
441  dict
442  ),
443  primaryMesh_(mesh),
444  time_(mesh.time()),
445  active_(dict.get<Switch>("active")),
446  infoOutput_(false),
447  modelName_(modelName),
448  coeffs_(dict.subOrEmptyDict(modelName + "Coeffs")),
449  outputPropertiesPtr_(nullptr),
450  primaryPatchIDs_(),
451  intCoupledPatchIDs_(),
452  regionName_(dict.lookup("region")),
453  functions_(*this, subOrEmptyDict("functions"))
454 {
455  if (active_)
456  {
457  constructMeshObjects();
458  initialise();
459 
460  if (readFields)
461  {
462  read(dict);
463  }
464  }
465 }
466 
467 
468 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
469 
471 {
472  if (active_)
473  {
474  Info<< "\nEvolving " << modelName_ << " for region "
475  << regionMesh().name() << endl;
476 
477  preEvolveRegion();
478 
479  evolveRegion();
480 
481  postEvolveRegion();
482 
483  // Provide some feedback
484  if (infoOutput_)
485  {
486  Info<< incrIndent;
487  info();
488  Info<< endl << decrIndent;
489  }
490 
491  if (time_.writeTime())
492  {
493  outputProperties().writeObject
494  (
495  IOstreamOption(IOstreamOption::ASCII, time_.writeCompression()),
496  true
497  );
498  }
499  }
500 }
501 
504 {
505  functions_.preEvolveRegion();
506 }
507 
508 
510 {}
511 
514 {
515  functions_.postEvolveRegion();
516 }
517 
518 
520 {}
521 
522 
523 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
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
const polyBoundaryMesh & pbm
dictionary dict
virtual void evolveRegion()
Evolve the region.
Definition: regionModel.C:502
virtual void info()
Provide some feedback.
Definition: regionModel.C:512
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:1184
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:414
dynamicFvMesh & mesh
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:496
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:584
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:268
virtual void postEvolveRegion()
Post-evolve region.
Definition: regionModel.C:506
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:113
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 time name of given scalar time formatted with the given precision.
Definition: Time.C:770
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
virtual void evolve()
Main driver routing to evolve the region - calls other evolves.
Definition: regionModel.C:463
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.
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:171
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.