regionModel.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-2018 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 
26 #include "regionModel.H"
27 #include "fvMesh.H"
28 #include "Time.H"
29 #include "mappedWallPolyPatch.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace regionModels
37 {
38  defineTypeNameAndDebug(regionModel, 0);
39 }
40 }
41 
42 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
43 
44 void Foam::regionModels::regionModel::constructMeshObjects()
45 {
46  // construct region mesh
47  if (!time_.foundObject<fvMesh>(regionName_))
48  {
49  regionMeshPtr_.reset
50  (
51  new fvMesh
52  (
53  IOobject
54  (
55  regionName_,
56  time_.timeName(),
57  time_,
59  )
60  )
61  );
62  }
63 }
64 
65 
66 void Foam::regionModels::regionModel::initialise()
67 {
68  if (debug)
69  {
70  Pout<< "regionModel::initialise()" << endl;
71  }
72 
73  label nBoundaryFaces = 0;
74  DynamicList<label> primaryPatchIDs;
75  DynamicList<label> intCoupledPatchIDs;
76  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
77 
78  forAll(rbm, patchi)
79  {
80  const polyPatch& regionPatch = rbm[patchi];
81  if (isA<mappedPatchBase>(regionPatch))
82  {
83  if (debug)
84  {
85  Pout<< "found " << mappedWallPolyPatch::typeName
86  << " " << regionPatch.name() << endl;
87  }
88 
89  intCoupledPatchIDs.append(patchi);
90 
91  nBoundaryFaces += regionPatch.faceCells().size();
92 
93  const mappedPatchBase& mapPatch =
94  refCast<const mappedPatchBase>(regionPatch);
95 
96  if
97  (
98  primaryMesh_.time().foundObject<polyMesh>
99  (
100  mapPatch.sampleRegion()
101  )
102  )
103  {
104 
105  const label primaryPatchi = mapPatch.samplePolyPatch().index();
106  primaryPatchIDs.append(primaryPatchi);
107  }
108  }
109  }
110 
111  primaryPatchIDs_.transfer(primaryPatchIDs);
112  intCoupledPatchIDs_.transfer(intCoupledPatchIDs);
113 
114  if (returnReduce(nBoundaryFaces, sumOp<label>()) == 0)
115  {
117  << "Region model has no mapped boundary conditions - transfer "
118  << "between regions will not be possible" << endl;
119  }
120 
121  if (!outputPropertiesPtr_.valid())
122  {
123  const fileName uniformPath(word("uniform")/"regionModels");
124 
125  outputPropertiesPtr_.reset
126  (
127  new IOdictionary
128  (
129  IOobject
130  (
131  regionName_ + "OutputProperties",
132  time_.timeName(),
133  uniformPath/regionName_,
134  primaryMesh_,
137  )
138  )
139  );
140  }
141 }
142 
143 
144 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
145 
147 {
148  if (regIOobject::read())
149  {
150  if (active_)
151  {
152  if (const dictionary* dictPtr = subDictPtr(modelName_ + "Coeffs"))
153  {
154  coeffs_ <<= *dictPtr;
155  }
156 
157  infoOutput_.readIfPresent("infoOutput", *this);
158  }
159 
160  return true;
161  }
162  else
163  {
164  return false;
165  }
166 }
167 
168 
170 {
171  if (active_)
172  {
173  if (const dictionary* dictPtr = dict.subDictPtr(modelName_ + "Coeffs"))
174  {
175  coeffs_ <<= *dictPtr;
176  }
177 
178  infoOutput_.readIfPresent("infoOutput", dict);
179 
180  return true;
181  }
182  else
183  {
184  return false;
185  }
186 }
187 
188 
191 (
192  const regionModel& nbrRegion,
193  const label regionPatchi,
194  const label nbrPatchi,
195  const bool flip
196 ) const
197 {
198  label nbrRegionID = findIndex(interRegionAMINames_, nbrRegion.name());
199 
200  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
201 
202  if (nbrRegionID != -1)
203  {
204  if (!interRegionAMI_[nbrRegionID].set(regionPatchi))
205  {
206  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
207  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
208 
209  int oldTag = UPstream::msgType();
210  UPstream::msgType() = oldTag + 1;
211 
212  interRegionAMI_[nbrRegionID].set
213  (
214  regionPatchi,
216  (
217  p,
218  nbrP,
220  true,
222  -1,
223  flip
224  )
225  );
226 
227  UPstream::msgType() = oldTag;
228  }
229 
230  return interRegionAMI_[nbrRegionID][regionPatchi];
231  }
232  else
233  {
234  label nbrRegionID = interRegionAMINames_.size();
235 
236  interRegionAMINames_.append(nbrRegion.name());
237 
238  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
239  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
240 
241  label nPatch = regionMesh().boundaryMesh().size();
242 
243 
244  interRegionAMI_.resize(nbrRegionID + 1);
245 
246  interRegionAMI_.set
247  (
248  nbrRegionID,
250  );
251 
252  int oldTag = UPstream::msgType();
253  UPstream::msgType() = oldTag + 1;
254 
255  interRegionAMI_[nbrRegionID].set
256  (
257  regionPatchi,
259  (
260  p,
261  nbrP,
263  true,
265  -1,
266  flip
267  )
268  );
269 
270  UPstream::msgType() = oldTag;
271 
272  return interRegionAMI_[nbrRegionID][regionPatchi];
273  }
274 }
275 
276 
278 (
279  const regionModel& nbrRegion,
280  const label regionPatchi
281 ) const
282 {
283  label nbrPatchi = -1;
284 
285  // region
286  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
287 
288  // boundary mesh
289  const polyBoundaryMesh& nbrPbm = nbrRegionMesh.boundaryMesh();
290 
291  const polyBoundaryMesh& pbm = regionMesh().boundaryMesh();
292 
293  if (regionPatchi > pbm.size() - 1)
294  {
296  << "region patch index out of bounds: "
297  << "region patch index = " << regionPatchi
298  << ", maximum index = " << pbm.size() - 1
299  << abort(FatalError);
300  }
301 
302  const polyPatch& pp = regionMesh().boundaryMesh()[regionPatchi];
303 
304  if (!isA<mappedPatchBase>(pp))
305  {
307  << "Expected a " << mappedPatchBase::typeName
308  << " patch, but found a " << pp.type() << abort(FatalError);
309  }
310 
311  const mappedPatchBase& mpb = refCast<const mappedPatchBase>(pp);
312 
313  // sample patch name on the primary region
314  const word& primaryPatchName = mpb.samplePatch();
315 
316  // find patch on nbr region that has the same sample patch name
317  forAll(nbrRegion.intCoupledPatchIDs(), j)
318  {
319  const label nbrRegionPatchi = nbrRegion.intCoupledPatchIDs()[j];
320 
321  const mappedPatchBase& mpb =
322  refCast<const mappedPatchBase>(nbrPbm[nbrRegionPatchi]);
323 
324  if (mpb.samplePatch() == primaryPatchName)
325  {
326  nbrPatchi = nbrRegionPatchi;
327  break;
328  }
329  }
330 
331  if (nbrPatchi == -1)
332  {
333  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
334 
336  << "Unable to find patch pair for local patch "
337  << p.name() << " and region " << nbrRegion.name()
338  << abort(FatalError);
339  }
340 
341  return nbrPatchi;
342 }
343 
344 
345 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
346 
347 Foam::regionModels::regionModel::regionModel
348 (
349  const fvMesh& mesh,
350  const word& regionType
351 )
352 :
354  (
355  IOobject
356  (
357  regionType + "Properties",
358  mesh.time().constant(),
359  mesh.time(),
362  )
363  ),
364  primaryMesh_(mesh),
365  time_(mesh.time()),
366  active_(false),
367  infoOutput_(false),
368  modelName_("none"),
369  regionMeshPtr_(nullptr),
370  coeffs_(dictionary::null),
371  outputPropertiesPtr_(nullptr),
372  primaryPatchIDs_(),
373  intCoupledPatchIDs_(),
374  regionName_("none"),
375  functions_(*this),
376  interRegionAMINames_(),
377  interRegionAMI_()
378 {}
379 
380 
381 Foam::regionModels::regionModel::regionModel
382 (
383  const fvMesh& mesh,
384  const word& regionType,
385  const word& modelName,
386  bool readFields
387 )
388 :
390  (
391  IOobject
392  (
393  regionType + "Properties",
394  mesh.time().constant(),
395  mesh.time(),
398  )
399  ),
400  primaryMesh_(mesh),
401  time_(mesh.time()),
402  active_(lookup("active")),
403  infoOutput_(true),
404  modelName_(modelName),
405  regionMeshPtr_(nullptr),
406  coeffs_(subOrEmptyDict(modelName + "Coeffs")),
407  outputPropertiesPtr_(nullptr),
408  primaryPatchIDs_(),
409  intCoupledPatchIDs_(),
410  regionName_(lookup("regionName")),
411  functions_(*this, subOrEmptyDict("functions"))
412 {
413  if (active_)
414  {
415  constructMeshObjects();
416  initialise();
417 
418  if (readFields)
419  {
420  read();
421  }
422  }
423 }
424 
425 
426 Foam::regionModels::regionModel::regionModel
427 (
428  const fvMesh& mesh,
429  const word& regionType,
430  const word& modelName,
431  const dictionary& dict,
432  bool readFields
433 )
434 :
436  (
437  IOobject
438  (
439  regionType + "Properties",
440  mesh.time().constant(),
441  mesh.time(),
444  true
445  ),
446  dict
447  ),
448  primaryMesh_(mesh),
449  time_(mesh.time()),
450  active_(dict.lookup("active")),
451  infoOutput_(false),
452  modelName_(modelName),
453  regionMeshPtr_(nullptr),
454  coeffs_(dict.subOrEmptyDict(modelName + "Coeffs")),
455  outputPropertiesPtr_(nullptr),
456  primaryPatchIDs_(),
457  intCoupledPatchIDs_(),
458  regionName_(dict.lookup("regionName")),
459  functions_(*this, subOrEmptyDict("functions"))
460 {
461  if (active_)
462  {
463  constructMeshObjects();
464  initialise();
465 
466  if (readFields)
467  {
468  read(dict);
469  }
470  }
471 }
472 
473 
474 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
475 
477 {}
478 
479 
480 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
481 
483 {
484  if (active_)
485  {
486  Info<< "\nEvolving " << modelName_ << " for region "
487  << regionMesh().name() << endl;
488 
489  // read();
490 
491  preEvolveRegion();
492 
493  evolveRegion();
494 
495  postEvolveRegion();
496 
497  // Provide some feedback
498  if (infoOutput_)
499  {
500  Info<< incrIndent;
501  info();
502  Info<< endl << decrIndent;
503  }
504 
505  if (time_.writeTime())
506  {
507  outputProperties().writeObject
508  (
511  time_.writeCompression(),
512  true
513  );
514  }
515  }
516 }
517 
518 
520 {
521  functions_.preEvolveRegion();
522 }
523 
524 
526 {}
527 
528 
530 {
531  functions_.postEvolveRegion();
532 }
533 
534 
536 {}
537 
538 
539 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const word & name() const
Return name.
virtual void evolveRegion()
Evolve the region.
Definition: regionModel.C:525
virtual void info()
Provide some feedback.
Definition: regionModel.C:535
virtual bool read()
Read object.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:278
static const dictionary null
Null dictionary.
Definition: dictionary.H:202
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:474
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
defineTypeNameAndDebug(regionModel, 0)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:519
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel.C:146
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
Definition: regionModel.C:278
virtual void postEvolveRegion()
Post-evolve region.
Definition: regionModel.C:529
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
const word & name() const
Name function is needed to disambiguate those inherited.
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:191
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dictionary * subDictPtr(const word &) const
Find and return a sub-dictionary pointer if present.
Definition: dictionary.C:662
Foam::polyBoundaryMesh.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:240
virtual ~regionModel()
Destructor.
Definition: regionModel.C:476
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
static const versionNumber currentVersion
Current version number.
Definition: IOstream.H:206
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
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:482
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
messageStream Info
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
Base class for region models.
Definition: regionModel.H:57
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
volScalarField & p
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:179
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:233
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:727
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576