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-2022 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 timeIOdictionary
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 
175 (
176  const regionModel& nbrRegion,
177  const label regionPatchi
178 ) const
179 {
180  label nbrPatchi = -1;
181 
182  // region
183  const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
184 
185  // boundary mesh
186  const polyBoundaryMesh& nbrPbm = nbrRegionMesh.boundaryMesh();
187 
188  const polyBoundaryMesh& pbm = regionMesh().boundaryMesh();
189 
190  if (regionPatchi > pbm.size() - 1)
191  {
193  << "region patch index out of bounds: "
194  << "region patch index = " << regionPatchi
195  << ", maximum index = " << pbm.size() - 1
196  << abort(FatalError);
197  }
198 
199  const polyPatch& pp = regionMesh().boundaryMesh()[regionPatchi];
200 
201  if (!isA<mappedPatchBase>(pp))
202  {
204  << "Expected a " << mappedPatchBase::typeName
205  << " patch, but found a " << pp.type() << abort(FatalError);
206  }
207 
208  const mappedPatchBase& mpb = refCast<const mappedPatchBase>(pp);
209 
210  // sample patch name on the primary region
211  const word& primaryPatchName = mpb.samplePatch();
212 
213  // find patch on nbr region that has the same sample patch name
214  forAll(nbrRegion.intCoupledPatchIDs(), j)
215  {
216  const label nbrRegionPatchi = nbrRegion.intCoupledPatchIDs()[j];
217 
218  const mappedPatchBase& mpb =
219  refCast<const mappedPatchBase>(nbrPbm[nbrRegionPatchi]);
220 
221  if (mpb.samplePatch() == primaryPatchName)
222  {
223  nbrPatchi = nbrRegionPatchi;
224  break;
225  }
226  }
227 
228  if (nbrPatchi == -1)
229  {
230  const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
231 
233  << "Unable to find patch pair for local patch "
234  << p.name() << " and region " << nbrRegion.name()
235  << abort(FatalError);
236  }
237 
238  return nbrPatchi;
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
243 
244 Foam::regionModels::regionModel::regionModel
245 (
246  const fvMesh& mesh,
247  const word& regionType
248 )
249 :
251  (
252  IOobject
253  (
254  regionType + "Properties",
255  mesh.time().constant(),
256  mesh.time(),
259  )
260  ),
261  primaryMesh_(mesh),
262  time_(mesh.time()),
263  infoOutput_(false),
264  modelName_("none"),
265  regionMeshPtr_(nullptr),
266  coeffs_(dictionary::null),
267  outputPropertiesPtr_(nullptr),
268  primaryPatchIDs_(),
269  intCoupledPatchIDs_(),
270  regionName_("none"),
271  functions_(*this)
272 {}
273 
274 
275 Foam::regionModels::regionModel::regionModel
276 (
277  const fvMesh& mesh,
278  const word& regionType,
279  const word& modelName,
280  bool readFields
281 )
282 :
284  (
285  IOobject
286  (
287  regionType + "Properties",
288  mesh.time().constant(),
289  mesh.time(),
292  )
293  ),
294  primaryMesh_(mesh),
295  time_(mesh.time()),
296  infoOutput_(true),
297  modelName_(modelName),
298  regionMeshPtr_(nullptr),
299  coeffs_(optionalSubDict(modelName + "Coeffs")),
300  outputPropertiesPtr_(nullptr),
301  primaryPatchIDs_(),
302  intCoupledPatchIDs_(),
303  regionName_(lookup("regionName")),
304  functions_(*this, subOrEmptyDict("functions"))
305 {
306  constructMeshObjects();
307  initialise();
308 
309  if (readFields)
310  {
311  read();
312  }
313 }
314 
315 
316 Foam::regionModels::regionModel::regionModel
317 (
318  const fvMesh& mesh,
319  const word& regionType,
320  const word& modelName,
321  const dictionary& dict,
322  bool readFields
323 )
324 :
326  (
327  IOobject
328  (
329  regionType + "Properties",
330  mesh.time().constant(),
331  mesh.time(),
334  true
335  ),
336  dict
337  ),
338  primaryMesh_(mesh),
339  time_(mesh.time()),
340  infoOutput_(false),
341  modelName_(modelName),
342  regionMeshPtr_(nullptr),
343  coeffs_(dict.optionalSubDict(modelName + "Coeffs")),
344  outputPropertiesPtr_(nullptr),
345  primaryPatchIDs_(),
346  intCoupledPatchIDs_(),
347  regionName_(dict.lookup("regionName")),
348  functions_(*this, subOrEmptyDict("functions"))
349 {
350  constructMeshObjects();
351  initialise();
352 
353  if (readFields)
354  {
355  read(dict);
356  }
357 }
358 
359 
360 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
361 
363 {}
364 
365 
366 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
367 
369 {
370  Info<< "\nEvolving " << modelName_ << " for region "
371  << regionMesh().name() << endl;
372 
373  preEvolveRegion();
374 
375  evolveRegion();
376 
377  postEvolveRegion();
378 
379  // Provide some feedback
380  if (infoOutput_)
381  {
382  Info<< incrIndent;
383  info();
384  Info<< endl << decrIndent;
385  }
386 
387  if (time_.writeTime())
388  {
389  outputProperties().writeObject
390  (
393  time_.writeCompression(),
394  true
395  );
396  }
397 }
398 
399 
401 {
402  functions_.preEvolveRegion();
403 }
404 
405 
407 {}
408 
409 
411 {
412  functions_.postEvolveRegion();
413 }
414 
415 
417 {}
418 
419 
420 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
const word & name() const
Return name.
Definition: IOobject.H:315
virtual void evolveRegion()
Evolve the region.
Definition: regionModel.C:406
virtual void info()
Provide some feedback.
Definition: regionModel.C:416
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:306
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:297
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
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:1060
virtual void preEvolveRegion()
Pre-evolve region.
Definition: regionModel.C:400
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:175
virtual void postEvolveRegion()
Post-evolve region.
Definition: regionModel.C:410
#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:972
Foam::polyBoundaryMesh.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
virtual ~regionModel()
Destructor.
Definition: regionModel.C:362
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual void evolve()
Main driver routing to evolve the region - calls other evolves.
Definition: regionModel.C:368
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:98
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864