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-2019 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 {
69 
70  label nBoundaryFaces = 0;
71  DynamicList<label> primaryPatchIDs;
72  DynamicList<label> intCoupledPatchIDs;
73  const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
74 
75  forAll(rbm, patchi)
76  {
77  const polyPatch& regionPatch = rbm[patchi];
78  if (isA<mappedPatchBase>(regionPatch))
79  {
81  << "found " << mappedWallPolyPatch::typeName
82  << " " << regionPatch.name() << endl;
83 
84  intCoupledPatchIDs.append(patchi);
85 
86  nBoundaryFaces += regionPatch.faceCells().size();
87 
88  const mappedPatchBase& mapPatch =
89  refCast<const mappedPatchBase>(regionPatch);
90 
91  if
92  (
93  primaryMesh_.time().foundObject<polyMesh>
94  (
95  mapPatch.sampleRegion()
96  )
97  )
98  {
99 
100  const label primaryPatchi = mapPatch.samplePolyPatch().index();
101  primaryPatchIDs.append(primaryPatchi);
102  }
103  }
104  }
105 
106  primaryPatchIDs_.transfer(primaryPatchIDs);
107  intCoupledPatchIDs_.transfer(intCoupledPatchIDs);
108 
109  if (returnReduce(nBoundaryFaces, sumOp<label>()) == 0)
110  {
112  << "Region model has no mapped boundary conditions - transfer "
113  << "between regions will not be possible" << endl;
114  }
115 
116  if (!outputPropertiesPtr_.valid())
117  {
118  const fileName uniformPath(word("uniform")/"regionModels");
119 
120  outputPropertiesPtr_.reset
121  (
122  new IOdictionary
123  (
124  IOobject
125  (
126  regionName_ + "OutputProperties",
127  time_.timeName(),
128  uniformPath/regionName_,
129  primaryMesh_,
132  )
133  )
134  );
135  }
136 }
137 
138 
139 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
140 
142 {
143  if (regIOobject::read())
144  {
145  if (active_)
146  {
147  if (const dictionary* dictPtr = subDictPtr(modelName_ + "Coeffs"))
148  {
149  coeffs_ <<= *dictPtr;
150  }
151 
152  infoOutput_.readIfPresent("infoOutput", *this);
153  }
154 
155  return true;
156  }
157  else
158  {
159  return false;
160  }
161 }
162 
163 
165 {
166  if (active_)
167  {
168  if (const dictionary* dictPtr = dict.subDictPtr(modelName_ + "Coeffs"))
169  {
170  coeffs_ <<= *dictPtr;
171  }
172 
173  infoOutput_.readIfPresent("infoOutput", dict);
174 
175  return true;
176  }
177  else
178  {
179  return false;
180  }
181 }
182 
183 
186 (
187  const regionModel& nbrRegion,
188  const label regionPatchi,
189  const label nbrPatchi,
190  const bool flip
191 ) const
192 {
193  label nbrRegionID = findIndex(interRegionAMINames_, nbrRegion.name());
194 
195  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
196 
197  if (nbrRegionID != -1)
198  {
199  if (!interRegionAMI_[nbrRegionID].set(regionPatchi))
200  {
201  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
202  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
203 
204  int oldTag = UPstream::msgType();
205  UPstream::msgType() = oldTag + 1;
206 
207  interRegionAMI_[nbrRegionID].set
208  (
209  regionPatchi,
210  new AMIInterpolation
211  (
212  p,
213  nbrP,
215  true,
217  -1,
218  flip
219  )
220  );
221 
222  UPstream::msgType() = oldTag;
223  }
224 
225  return interRegionAMI_[nbrRegionID][regionPatchi];
226  }
227  else
228  {
229  label nbrRegionID = interRegionAMINames_.size();
230 
231  interRegionAMINames_.append(nbrRegion.name());
232 
233  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
234  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
235 
236  label nPatch = regionMesh().boundaryMesh().size();
237 
238 
239  interRegionAMI_.resize(nbrRegionID + 1);
240 
241  interRegionAMI_.set
242  (
243  nbrRegionID,
244  new PtrList<AMIInterpolation>(nPatch)
245  );
246 
247  int oldTag = UPstream::msgType();
248  UPstream::msgType() = oldTag + 1;
249 
250  interRegionAMI_[nbrRegionID].set
251  (
252  regionPatchi,
253  new AMIInterpolation
254  (
255  p,
256  nbrP,
258  true,
260  -1,
261  flip
262  )
263  );
264 
265  UPstream::msgType() = oldTag;
266 
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(),
357  )
358  ),
359  primaryMesh_(mesh),
360  time_(mesh.time()),
361  active_(false),
362  infoOutput_(false),
363  modelName_("none"),
364  regionMeshPtr_(nullptr),
365  coeffs_(dictionary::null),
366  outputPropertiesPtr_(nullptr),
367  primaryPatchIDs_(),
368  intCoupledPatchIDs_(),
369  regionName_("none"),
370  functions_(*this),
371  interRegionAMINames_(),
372  interRegionAMI_()
373 {}
374 
375 
376 Foam::regionModels::regionModel::regionModel
377 (
378  const fvMesh& mesh,
379  const word& regionType,
380  const word& modelName,
381  bool readFields
382 )
383 :
385  (
386  IOobject
387  (
388  regionType + "Properties",
389  mesh.time().constant(),
390  mesh.time(),
393  )
394  ),
395  primaryMesh_(mesh),
396  time_(mesh.time()),
397  active_(lookup("active")),
398  infoOutput_(true),
399  modelName_(modelName),
400  regionMeshPtr_(nullptr),
401  coeffs_(subOrEmptyDict(modelName + "Coeffs")),
402  outputPropertiesPtr_(nullptr),
403  primaryPatchIDs_(),
404  intCoupledPatchIDs_(),
405  regionName_(lookup("regionName")),
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 :
431  (
432  IOobject
433  (
434  regionType + "Properties",
435  mesh.time().constant(),
436  mesh.time(),
439  true
440  ),
441  dict
442  ),
443  primaryMesh_(mesh),
444  time_(mesh.time()),
445  active_(dict.lookup("active")),
446  infoOutput_(false),
447  modelName_(modelName),
448  regionMeshPtr_(nullptr),
449  coeffs_(dict.subOrEmptyDict(modelName + "Coeffs")),
450  outputPropertiesPtr_(nullptr),
451  primaryPatchIDs_(),
452  intCoupledPatchIDs_(),
453  regionName_(dict.lookup("regionName")),
454  functions_(*this, subOrEmptyDict("functions"))
455 {
456  if (active_)
457  {
458  constructMeshObjects();
459  initialise();
460 
461  if (readFields)
462  {
463  read(dict);
464  }
465  }
466 }
467 
468 
469 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
470 
472 {}
473 
474 
475 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
476 
478 {
479  if (active_)
480  {
481  Info<< "\nEvolving " << modelName_ << " for region "
482  << regionMesh().name() << endl;
483 
484  // read();
485 
486  preEvolveRegion();
487 
488  evolveRegion();
489 
490  postEvolveRegion();
491 
492  // Provide some feedback
493  if (infoOutput_)
494  {
495  Info<< incrIndent;
496  info();
497  Info<< endl << decrIndent;
498  }
499 
500  if (time_.writeTime())
501  {
502  outputProperties().writeObject
503  (
506  time_.writeCompression(),
507  true
508  );
509  }
510  }
511 }
512 
513 
515 {
516  functions_.preEvolveRegion();
517 }
518 
519 
521 {}
522 
523 
525 {
526  functions_.postEvolveRegion();
527 }
528 
529 
531 {}
532 
533 
534 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:520
virtual void info()
Provide some feedback.
Definition: regionModel.C:530
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:158
#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:241
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
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:514
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()
Read control parameters from dictionary.
Definition: regionModel.C:141
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
Definition: regionModel.C:273
virtual void postEvolveRegion()
Post-evolve region.
Definition: regionModel.C:524
#define DebugInFunction
Report an information message using Foam::Info.
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.
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:904
Foam::polyBoundaryMesh.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
virtual ~regionModel()
Destructor.
Definition: regionModel.C:471
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:70
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:477
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:55
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:228
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual const AMIInterpolation & 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:186
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:969
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812