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-2020 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 (const dictionary* dictPtr = subDictPtr(modelName_ + "Coeffs"))
146  {
147  coeffs_ <<= *dictPtr;
148  }
149 
150  infoOutput_.readIfPresent("infoOutput", *this);
151 
152  return true;
153  }
154  else
155  {
156  return false;
157  }
158 }
159 
160 
162 {
163  if (const dictionary* dictPtr = dict.subDictPtr(modelName_ + "Coeffs"))
164  {
165  coeffs_ <<= *dictPtr;
166  }
167 
168  infoOutput_.readIfPresent("infoOutput", dict);
169 
170  return true;
171 }
172 
173 
176 (
177  const regionModel& nbrRegion,
178  const label regionPatchi,
179  const label nbrPatchi,
180  const bool flip
181 ) const
182 {
183  label nbrRegionID = findIndex(interRegionAMINames_, nbrRegion.name());
184 
185  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
186 
187  if (nbrRegionID != -1)
188  {
189  if (!interRegionAMI_[nbrRegionID].set(regionPatchi))
190  {
191  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
192  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
193 
194  int oldTag = UPstream::msgType();
195  UPstream::msgType() = oldTag + 1;
196 
197  interRegionAMI_[nbrRegionID].set
198  (
199  regionPatchi,
200  new AMIInterpolation
201  (
202  p,
203  nbrP,
205  true,
207  -1,
208  flip
209  )
210  );
211 
212  UPstream::msgType() = oldTag;
213  }
214 
215  return interRegionAMI_[nbrRegionID][regionPatchi];
216  }
217  else
218  {
219  label nbrRegionID = interRegionAMINames_.size();
220 
221  interRegionAMINames_.append(nbrRegion.name());
222 
223  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
224  const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
225 
226  label nPatch = regionMesh().boundaryMesh().size();
227 
228 
229  interRegionAMI_.resize(nbrRegionID + 1);
230 
231  interRegionAMI_.set
232  (
233  nbrRegionID,
234  new PtrList<AMIInterpolation>(nPatch)
235  );
236 
237  int oldTag = UPstream::msgType();
238  UPstream::msgType() = oldTag + 1;
239 
240  interRegionAMI_[nbrRegionID].set
241  (
242  regionPatchi,
243  new AMIInterpolation
244  (
245  p,
246  nbrP,
248  true,
250  -1,
251  flip
252  )
253  );
254 
255  UPstream::msgType() = oldTag;
256 
257  return interRegionAMI_[nbrRegionID][regionPatchi];
258  }
259 }
260 
261 
263 (
264  const regionModel& nbrRegion,
265  const label regionPatchi
266 ) const
267 {
268  label nbrPatchi = -1;
269 
270  // region
271  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
272 
273  // boundary mesh
274  const polyBoundaryMesh& nbrPbm = nbrRegionMesh.boundaryMesh();
275 
276  const polyBoundaryMesh& pbm = regionMesh().boundaryMesh();
277 
278  if (regionPatchi > pbm.size() - 1)
279  {
281  << "region patch index out of bounds: "
282  << "region patch index = " << regionPatchi
283  << ", maximum index = " << pbm.size() - 1
284  << abort(FatalError);
285  }
286 
287  const polyPatch& pp = regionMesh().boundaryMesh()[regionPatchi];
288 
289  if (!isA<mappedPatchBase>(pp))
290  {
292  << "Expected a " << mappedPatchBase::typeName
293  << " patch, but found a " << pp.type() << abort(FatalError);
294  }
295 
296  const mappedPatchBase& mpb = refCast<const mappedPatchBase>(pp);
297 
298  // sample patch name on the primary region
299  const word& primaryPatchName = mpb.samplePatch();
300 
301  // find patch on nbr region that has the same sample patch name
302  forAll(nbrRegion.intCoupledPatchIDs(), j)
303  {
304  const label nbrRegionPatchi = nbrRegion.intCoupledPatchIDs()[j];
305 
306  const mappedPatchBase& mpb =
307  refCast<const mappedPatchBase>(nbrPbm[nbrRegionPatchi]);
308 
309  if (mpb.samplePatch() == primaryPatchName)
310  {
311  nbrPatchi = nbrRegionPatchi;
312  break;
313  }
314  }
315 
316  if (nbrPatchi == -1)
317  {
318  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
319 
321  << "Unable to find patch pair for local patch "
322  << p.name() << " and region " << nbrRegion.name()
323  << abort(FatalError);
324  }
325 
326  return nbrPatchi;
327 }
328 
329 
330 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
331 
332 Foam::regionModels::regionModel::regionModel
333 (
334  const fvMesh& mesh,
335  const word& regionType
336 )
337 :
339  (
340  IOobject
341  (
342  regionType + "Properties",
343  mesh.time().constant(),
344  mesh.time(),
347  )
348  ),
349  primaryMesh_(mesh),
350  time_(mesh.time()),
351  infoOutput_(false),
352  modelName_("none"),
353  regionMeshPtr_(nullptr),
354  coeffs_(dictionary::null),
355  outputPropertiesPtr_(nullptr),
356  primaryPatchIDs_(),
357  intCoupledPatchIDs_(),
358  regionName_("none"),
359  functions_(*this),
360  interRegionAMINames_(),
361  interRegionAMI_()
362 {}
363 
364 
365 Foam::regionModels::regionModel::regionModel
366 (
367  const fvMesh& mesh,
368  const word& regionType,
369  const word& modelName,
370  bool readFields
371 )
372 :
374  (
375  IOobject
376  (
377  regionType + "Properties",
378  mesh.time().constant(),
379  mesh.time(),
382  )
383  ),
384  primaryMesh_(mesh),
385  time_(mesh.time()),
386  infoOutput_(true),
387  modelName_(modelName),
388  regionMeshPtr_(nullptr),
389  coeffs_(optionalSubDict(modelName + "Coeffs")),
390  outputPropertiesPtr_(nullptr),
391  primaryPatchIDs_(),
392  intCoupledPatchIDs_(),
393  regionName_(lookup("regionName")),
394  functions_(*this, subOrEmptyDict("functions"))
395 {
396  constructMeshObjects();
397  initialise();
398 
399  if (readFields)
400  {
401  read();
402  }
403 }
404 
405 
406 Foam::regionModels::regionModel::regionModel
407 (
408  const fvMesh& mesh,
409  const word& regionType,
410  const word& modelName,
411  const dictionary& dict,
412  bool readFields
413 )
414 :
416  (
417  IOobject
418  (
419  regionType + "Properties",
420  mesh.time().constant(),
421  mesh.time(),
424  true
425  ),
426  dict
427  ),
428  primaryMesh_(mesh),
429  time_(mesh.time()),
430  infoOutput_(false),
431  modelName_(modelName),
432  regionMeshPtr_(nullptr),
433  coeffs_(dict.optionalSubDict(modelName + "Coeffs")),
434  outputPropertiesPtr_(nullptr),
435  primaryPatchIDs_(),
436  intCoupledPatchIDs_(),
437  regionName_(dict.lookup("regionName")),
438  functions_(*this, subOrEmptyDict("functions"))
439 {
440  constructMeshObjects();
441  initialise();
442 
443  if (readFields)
444  {
445  read(dict);
446  }
447 }
448 
449 
450 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
451 
453 {}
454 
455 
456 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
457 
459 {
460  Info<< "\nEvolving " << modelName_ << " for region "
461  << regionMesh().name() << endl;
462 
463  preEvolveRegion();
464 
465  evolveRegion();
466 
467  postEvolveRegion();
468 
469  // Provide some feedback
470  if (infoOutput_)
471  {
472  Info<< incrIndent;
473  info();
474  Info<< endl << decrIndent;
475  }
476 
477  if (time_.writeTime())
478  {
479  outputProperties().writeObject
480  (
483  time_.writeCompression(),
484  true
485  );
486  }
487 }
488 
489 
491 {
492  functions_.preEvolveRegion();
493 }
494 
495 
497 {}
498 
499 
501 {
502  functions_.postEvolveRegion();
503 }
504 
505 
507 {}
508 
509 
510 // ************************************************************************* //
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.
const word & name() const
Return name.
Definition: IOobject.H:303
virtual void evolveRegion()
Evolve the region.
Definition: regionModel.C:496
virtual void info()
Provide some feedback.
Definition: regionModel.C:506
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:291
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:53
defineTypeNameAndDebug(regionModel, 0)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
stressControl lookup("compactNormalStress") >> compactNormalStress
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1040
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:490
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:263
virtual void postEvolveRegion()
Post-evolve region.
Definition: regionModel.C:500
#define DebugInFunction
Report an information message using Foam::Info.
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
static const dictionary null
Null dictionary.
Definition: dictionary.H:242
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:55
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
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:952
Foam::polyBoundaryMesh.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
virtual ~regionModel()
Destructor.
Definition: regionModel.C:452
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:203
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:458
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:173
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:176
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844