faceAgglomerate.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-2024 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 Application
25  faceAgglomerate
26 
27 Description
28 
29  Agglomerate boundary faces using the pairPatchAgglomeration algorithm.
30  It writes a map from the fine to coarse grid.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "argList.H"
35 #include "fvMesh.H"
36 #include "Time.H"
37 #include "volFields.H"
38 #include "CompactListList.H"
39 #include "pairPatchAgglomeration.H"
40 #include "labelListIOList.H"
41 #include "syncTools.H"
42 #include "globalIndex.H"
43 #include "systemDict.H"
44 
45 using namespace Foam;
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 int main(int argc, char *argv[])
50 {
51  #include "addRegionOption.H"
52  #include "addDictOption.H"
53  #include "setRootCase.H"
54  #include "createTime.H"
56 
57  const word dictName("viewFactorsDict");
58 
59  // Read control dictionary
60  const IOdictionary agglomDict(systemDict(dictName, args, runTime));
61 
62  bool writeAgglom = readBool(agglomDict.lookup("writeFacesAgglomeration"));
63 
64 
65  const polyBoundaryMesh& boundary = mesh.boundaryMesh();
66 
67  labelListIOList finalAgglom
68  (
69  IOobject
70  (
71  "finalAgglom",
72  mesh.facesInstance(),
73  mesh,
76  false
77  ),
78  boundary.size()
79  );
80 
81  label nCoarseFaces = 0;
82 
83  forAllConstIter(dictionary, agglomDict, iter)
84  {
85  labelList patchids = boundary.findIndices(iter().keyword());
86  forAll(patchids, i)
87  {
88  label patchi = patchids[i];
89  const polyPatch& pp = boundary[patchi];
90  if (!pp.coupled())
91  {
92  Info << "\nAgglomerating patch : " << pp.name() << endl;
93  pairPatchAgglomeration agglomObject
94  (
95  pp,
96  agglomDict.subDict(pp.name())
97  );
98  agglomObject.agglomerate();
99  finalAgglom[patchi] =
100  agglomObject.restrictTopBottomAddressing();
101 
102  if (finalAgglom[patchi].size())
103  {
104  nCoarseFaces += max(finalAgglom[patchi] + 1);
105  }
106  }
107  }
108  }
109 
110 
111  // - All patches which are not agglomerated are identity for finalAgglom
112  forAll(boundary, patchid)
113  {
114  if (finalAgglom[patchid].size() == 0)
115  {
116  finalAgglom[patchid] = identityMap(boundary[patchid].size());
117  }
118  }
119 
120  // Sync agglomeration across coupled patches
121  labelList nbrAgglom(mesh.nFaces() - mesh.nInternalFaces(), -1);
122 
123  forAll(boundary, patchid)
124  {
125  const polyPatch& pp = boundary[patchid];
126  if (pp.coupled())
127  {
128  finalAgglom[patchid] = identityMap(pp.size());
129  forAll(pp, i)
130  {
131  nbrAgglom[pp.start() - mesh.nInternalFaces() + i] =
132  finalAgglom[patchid][i];
133  }
134  }
135  }
136 
137  syncTools::swapBoundaryFaceList(mesh, nbrAgglom);
138  forAll(boundary, patchid)
139  {
140  const polyPatch& pp = boundary[patchid];
141  if (pp.coupled() && !refCast<const coupledPolyPatch>(pp).owner())
142  {
143  forAll(pp, i)
144  {
145  finalAgglom[patchid][i] =
146  nbrAgglom[pp.start() - mesh.nInternalFaces() + i];
147  }
148  }
149  }
150 
151  finalAgglom.write();
152 
153  if (writeAgglom)
154  {
155  globalIndex index(nCoarseFaces);
156  volScalarField facesAgglomeration
157  (
158  IOobject
159  (
160  "facesAgglomeration",
161  mesh.time().name(),
162  mesh,
165  ),
166  mesh,
168  );
169 
170  volScalarField::Boundary& facesAgglomerationBf =
171  facesAgglomeration.boundaryFieldRef();
172 
173  label coarsePatchIndex = 0;
174  forAll(boundary, patchid)
175  {
176  const polyPatch& pp = boundary[patchid];
177  if (pp.size() > 0)
178  {
179  fvPatchScalarField& bFacesAgglomeration =
180  facesAgglomerationBf[patchid];
181 
182  forAll(bFacesAgglomeration, j)
183  {
184  bFacesAgglomeration[j] =
185  index.toGlobal
186  (
188  finalAgglom[patchid][j] + coarsePatchIndex
189  );
190  }
191 
192  coarsePatchIndex += max(finalAgglom[patchid]) + 1;
193  }
194  }
195 
196  Info<< "\nWriting facesAgglomeration" << endl;
197  facesAgglomeration.write();
198  }
199 
200  Info<< "End\n" << endl;
201  return 0;
202 }
203 
204 
205 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Generic GeometricBoundaryField class.
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual Ostream & write(const char)=0
Write character.
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
Primitive patch pair agglomerate method.
const word & name() const
Return name.
Foam::polyBoundaryMesh.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:290
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:436
A class for handling words, derived from string.
Definition: word.H:62
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
Namespace for OpenFOAM.
bool readBool(Istream &)
Definition: boolIO.C:66
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
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const dimensionSet dimless
messageStream Info
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
faceListList boundary(nPatches)
Foam::argList args(argc, argv)