topOZones.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) 2007-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 FOSS GP
10 -------------------------------------------------------------------------------
11 License
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 "topOZones.H"
30 #include "cellSet.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 }
38 
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
43 (
44  const dictionary& dict,
45  const word& zoneGroup
46 )
47 {
48  wordList zoneNames(dict.getOrDefault<wordList>(zoneGroup, wordList(0)));
49  labelList IDs(zoneNames.size(), -1);
50  forAll(zoneNames, zI)
51  {
52  IDs[zI] = mesh_.cellZones().findZoneID(zoneNames[zI]);
53  }
54  return IDs;
55 }
56 
57 
59 {
60  // Gather IO cells
61  DynamicList<label> IOcells(0);
62  for (const fvPatch& patch : mesh_.boundary())
63  {
64  if (patch.type() == "patch")
65  {
66  IOcells.append(patch.faceCells());
67  }
68  }
69 
70  // Add zone to cellZoneMesh and populate it
71  cellZoneMesh& cellZones = const_cast<fvMesh&>(mesh_).cellZones();
72  cellZone& IOcellsZone = cellZones("IOcells", true);
73  IOcellsZone = IOcells;
74  IOcellsID_ = cellZones.size() - 1;
75 
76  // Write as set
77  cellSet IOcellList(mesh_, "IOcellList", IOcells);
78  IOcellList.write();
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
84 Foam::topOZones::topOZones
85 (
86  const fvMesh& mesh,
87  const dictionary& dict
88 )
89 :
90  mesh_(mesh),
91  dict_(dict),
92  fixedPorousIDs_(getZoneIDs(dict, "fixedPorousZones")),
93  fixedPorousValues_
94  (
95  dict.getOrDefault<scalarList>
96  (
97  "fixedPorousValues", scalarList(fixedPorousIDs_.size(), 1.)
98  )
99  ),
100  fixedZeroPorousIDs_(getZoneIDs(dict, "fixedZeroPorousZones")),
101  adjointPorousIDs_(getZoneIDs(dict, "adjointPorousZones")),
102  IOcellsID_(-1),
103  betaMaxPtr_(betaMax::New(mesh, dict))
104 {
105  addIOcellsZone();
107  {
109  << "Number of fixedPorousValues and fixedPorousZones don't agree!"
110  << "\nSize of fixedPorousIDs is " << fixedPorousIDs_.size()
111  << " and size of fixedPorousValues is " << fixedPorousValues_.size()
112  << endl << endl
113  << exit(FatalError);
114  }
115 }
116 
117 
118 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
labelList fixedPorousIDs_
IDs of cellZones holding cells with constant alpha values.
Definition: topOZones.H:69
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
labelList getZoneIDs(const dictionary &dict, const word &zoneGroup)
Get zone IDs corresponding to a wordList, read from a dict.
Definition: topOZones.C:36
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:629
A class for handling words, derived from Foam::string.
Definition: word.H:63
void addIOcellsZone()
Add a cellZone containing the cells next to IO patches.
Definition: topOZones.C:51
const fvMesh & mesh_
Mesh reference.
Definition: topOZones.H:59
Base class for selecting the betaMax value, i.e. the value multiplying the Brinkman penalisation term...
Definition: betaMax.H:49
defineTypeNameAndDebug(combustionModel, 0)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const std::string patch
OpenFOAM patch number as a std::string.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
scalarList fixedPorousValues_
The constant alpha values of fixedPorousIDs_.
Definition: topOZones.H:74
Namespace for OpenFOAM.