Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "regionModel.H"
30#include "fvMesh.H"
31#include "Time.H"
32#include "mappedWallPolyPatch.H"
34#include "faceAreaWeightAMI.H"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace regionModels
41{
43}
44}
45
46// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
47
48void Foam::regionModels::regionModel::constructMeshObjects()
49{
50 fvMesh* regionMeshPtr = time_.getObjectPtr<fvMesh>(regionName_);
51
52 if (!regionMeshPtr)
53 {
54 regionMeshPtr = new fvMesh
55 (
57 (
60 time_,
64 )
65 );
66
67 regionMeshPtr->objectRegistry::store();
68 }
69}
70
71
72void Foam::regionModels::regionModel::initialise()
73{
74 if (debug)
75 {
76 Pout<< "regionModel::initialise()" << endl;
77 }
78
79 label nBoundaryFaces = 0;
80 DynamicList<label> primaryPatchIDs;
81 DynamicList<label> intCoupledPatchIDs;
82 const polyBoundaryMesh& rbm = regionMesh().boundaryMesh();
83
84 forAll(rbm, patchi)
85 {
86 const polyPatch& regionPatch = rbm[patchi];
87 if (isA<mappedPatchBase>(regionPatch))
88 {
89 if (debug)
90 {
91 Pout<< "found " << mappedWallPolyPatch::typeName
92 << " " << regionPatch.name() << endl;
93 }
94
95 intCoupledPatchIDs.append(patchi);
96
97 nBoundaryFaces += regionPatch.faceCells().size();
98
99 const mappedPatchBase& mapPatch =
101
102 if
103 (
104 primaryMesh_.time().foundObject<polyMesh>
105 (
106 mapPatch.sampleRegion()
107 )
108 )
109 {
110
111 const label primaryPatchi = mapPatch.samplePolyPatch().index();
112 primaryPatchIDs.append(primaryPatchi);
113 }
114 }
115 }
116
117 primaryPatchIDs_.transfer(primaryPatchIDs);
118 intCoupledPatchIDs_.transfer(intCoupledPatchIDs);
119
120 if (!returnReduceOr(nBoundaryFaces))
121 {
123 << "Region model has no mapped boundary conditions - transfer "
124 << "between regions will not be possible" << endl;
125 }
126
127 if (!outputPropertiesPtr_)
128 {
129 const fileName uniformPath(word("uniform")/"regionModels");
130
131 outputPropertiesPtr_.reset
132 (
133 new IOdictionary
134 (
135 IOobject
136 (
137 regionName_ + "OutputProperties",
138 time_.timeName(),
139 uniformPath/regionName_,
140 primaryMesh_,
143 )
144 )
145 );
146 }
147}
148
149
150// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
151
153{
154 if (regIOobject::read())
155 {
156 if (active_)
157 {
158 if (const dictionary* dictptr = findDict(modelName_ + "Coeffs"))
159 {
160 coeffs_ <<= *dictptr;
161 }
162
163 infoOutput_.readIfPresent("infoOutput", *this);
164 }
165
166 return true;
167 }
168
169 return false;
170}
171
172
174{
175 if (active_)
176 {
177 if (const dictionary* dictptr = dict.findDict(modelName_ + "Coeffs"))
178 {
179 coeffs_ <<= *dictptr;
180 }
181
182 infoOutput_.readIfPresent("infoOutput", dict);
183 return true;
185
186 return false;
187}
188
189
192(
193 const regionModel& nbrRegion,
194 const label regionPatchi,
195 const label nbrPatchi,
196 const bool flip
197) const
198{
199 label nbrRegionID = interRegionAMINames_.find(nbrRegion.name());
200
201 const fvMesh& nbrRegionMesh = nbrRegion.regionMesh();
202
203 if (nbrRegionID != -1)
204 {
205 if (!interRegionAMI_[nbrRegionID].set(regionPatchi))
206 {
207 const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
208 const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
209
210 const int oldTag = UPstream::incrMsgType();
211
212 interRegionAMI_[nbrRegionID].set
213 (
214 regionPatchi,
216 (
217 faceAreaWeightAMI::typeName,
218 true, // requireMatch
219 flip
220 )
221 );
222
223 interRegionAMI_[nbrRegionID][regionPatchi].calculate(p, nbrP);
224
225 UPstream::msgType(oldTag); // Restore tag
226 }
227
228 return interRegionAMI_[nbrRegionID][regionPatchi];
229 }
230 else
231 {
232 label nbrRegionID = interRegionAMINames_.size();
233
234 interRegionAMINames_.append(nbrRegion.name());
235
236 const polyPatch& p = regionMesh().boundaryMesh()[regionPatchi];
237 const polyPatch& nbrP = nbrRegionMesh.boundaryMesh()[nbrPatchi];
238
239 const label nPatch = regionMesh().boundaryMesh().size();
240
241
242 interRegionAMI_.resize(nbrRegionID + 1);
243
244 interRegionAMI_.set
245 (
246 nbrRegionID,
247 new PtrList<AMIPatchToPatchInterpolation>(nPatch)
248 );
249
250 const int oldTag = UPstream::incrMsgType();
251
252 interRegionAMI_[nbrRegionID].set
253 (
254 regionPatchi,
256 (
257 faceAreaWeightAMI::typeName,
258 true, // requireMatch
259 flip // reverse
260 )
261 );
262
263 interRegionAMI_[nbrRegionID][regionPatchi].calculate(p, nbrP);
264
265 UPstream::msgType(oldTag); // Restore tag
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
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 }
336 return nbrPatchi;
337}
338
339
340// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
341
342Foam::regionModels::regionModel::regionModel
343(
344 const fvMesh& mesh,
345 const word& regionType
346)
347:
349 (
351 (
352 regionType + "Properties",
353 mesh.time().constant(),
354 mesh.time(),
355 IOobject::NO_READ,
356 IOobject::NO_WRITE
357 )
358 ),
359 primaryMesh_(mesh),
360 time_(mesh.time()),
361 active_(false),
362 infoOutput_(false),
363 modelName_("none"),
364 coeffs_(),
365 outputPropertiesPtr_(nullptr),
366 primaryPatchIDs_(),
367 intCoupledPatchIDs_(),
368 regionName_("none"),
369 functions_(*this),
372{}
373
374
375Foam::regionModels::regionModel::regionModel
376(
377 const fvMesh& mesh,
378 const word& regionType,
379 const word& modelName,
380 bool readFields
381)
382:
384 (
386 (
387 regionType + "Properties",
388 mesh.time().constant(),
389 mesh.time(),
390 IOobject::MUST_READ,
391 IOobject::NO_WRITE
392 )
393 ),
394 primaryMesh_(mesh),
395 time_(mesh.time()),
396 active_(get<Switch>("active")),
397 infoOutput_(true),
398 modelName_(modelName),
399 coeffs_(subOrEmptyDict(modelName + "Coeffs")),
400 outputPropertiesPtr_(nullptr),
401 primaryPatchIDs_(),
402 intCoupledPatchIDs_(),
403 regionName_(lookup("region")),
404 functions_(*this, subOrEmptyDict("functions"))
405{
406 if (active_)
407 {
408 constructMeshObjects();
409 initialise();
410
411 if (readFields)
413 read();
414 }
415 }
416}
417
418
419Foam::regionModels::regionModel::regionModel
420(
421 const fvMesh& mesh,
422 const word& regionType,
423 const word& modelName,
424 const dictionary& dict,
425 bool readFields
426)
427:
428 IOdictionary
429 (
430 IOobject
431 (
432 regionType + "Properties",
433 mesh.time().constant(),
434 mesh.time(),
435 IOobject::NO_READ,
436 IOobject::NO_WRITE,
437 IOobject::REGISTER
438 ),
439 dict
440 ),
441 primaryMesh_(mesh),
442 time_(mesh.time()),
443 active_(dict.get<Switch>("active")),
444 infoOutput_(false),
445 modelName_(modelName),
446 coeffs_(dict.subOrEmptyDict(modelName + "Coeffs")),
447 outputPropertiesPtr_(nullptr),
448 primaryPatchIDs_(),
449 intCoupledPatchIDs_(),
450 regionName_(dict.lookup("region")),
451 functions_(*this, subOrEmptyDict("functions"))
452{
453 if (active_)
454 {
455 constructMeshObjects();
456 initialise();
457
458 if (readFields)
459 {
460 read(dict);
462 }
463}
464
465
466// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
467
469{
470 if (active_)
471 {
472 Info<< "\nEvolving " << modelName_ << " for region "
473 << regionMesh().name() << endl;
474
475 preEvolveRegion();
476
477 evolveRegion();
478
479 postEvolveRegion();
480
481 // Provide some feedback
482 if (infoOutput_)
483 {
485 info();
486 Info<< endl << decrIndent;
487 }
488
489 if (time_.writeTime())
490 {
491 outputProperties().writeObject
492 (
493 IOstreamOption(IOstreamOption::ASCII, time_.writeCompression()),
494 true
495 );
496 }
497 }
498}
499
502{
503 functions_.preEvolveRegion();
508{}
509
512{
513 functions_.postEvolveRegion();
514}
515
516
518{}
519
520
521// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A simple container for options an IOstream can normally have.
@ ASCII
"ascii" (normal default)
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name).
Definition OSstream.H:134
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary.
Definition dictionary.C:521
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream.
Definition dictionary.C:367
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE).
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
constant condensation/saturation model.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Lookup type of boundary radiation properties.
Definition lookup.H:60
virtual bool read()
Read object.
Base class for region models.
Definition regionModel.H:59
autoPtr< IOdictionary > outputPropertiesPtr_
Dictionary of output properties.
const Time & time_
Reference to the time database.
Definition regionModel.H:97
const word modelName_
Model name.
virtual void postEvolveRegion()
Post-evolve region.
Switch infoOutput_
Active information output.
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
labelList primaryPatchIDs_
List of patch IDs on the primary region coupled to this region.
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.
const Time & time() const noexcept
Return the reference to the time database.
dictionary coeffs_
Model coefficients dictionary.
const fvMesh & regionMesh() const
Return the region mesh database.
PtrList< PtrList< AMIPatchToPatchInterpolation > > interRegionAMI_
List of AMI objects per coupled region.
virtual void preEvolveRegion()
Pre-evolve region.
const labelList & intCoupledPatchIDs() const noexcept
List of patch IDs internally coupled with the primary region.
virtual void evolve()
Main driver routing to evolve the region - calls other evolves.
const fvMesh & primaryMesh_
Reference to the primary mesh database.
Definition regionModel.H:92
regionModelFunctionObjectList functions_
Region model function objects.
virtual void info()
Provide some feedback.
wordList interRegionAMINames_
List of region names this region is coupled to.
const IOdictionary & outputProperties() const
Return const access to the output properties dictionary.
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
const word & modelName() const noexcept
Return the model name.
virtual bool read()
Read control parameters from dictionary.
virtual void evolveRegion()
Evolve the region.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
#define WarningInFunction
Report a warning using Foam::Warning.
Different types of constants.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
AMIInterpolation AMIPatchToPatchInterpolation
Patch-to-patch interpolation == Foam::AMIInterpolation.
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299