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