All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2021 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 "unitConversion.H"
40 #include "pairPatchAgglomeration.H"
41 #include "labelListIOList.H"
42 #include "syncTools.H"
43 #include "globalIndex.H"
44 #include "systemDict.H"
45 
46 using namespace Foam;
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 int main(int argc, char *argv[])
51 {
52  #include "addRegionOption.H"
53  #include "addDictOption.H"
54  #include "setRootCase.H"
55  #include "createTime.H"
56  #include "createNamedMesh.H"
57 
58  const word dictName("viewFactorsDict");
59 
60  // Read control dictionary
61  const IOdictionary agglomDict(systemDict(dictName, args, runTime));
62 
63  bool writeAgglom = readBool(agglomDict.lookup("writeFacesAgglomeration"));
64 
65 
66  const polyBoundaryMesh& boundary = mesh.boundaryMesh();
67 
68  labelListIOList finalAgglom
69  (
70  IOobject
71  (
72  "finalAgglom",
74  mesh,
77  false
78  ),
79  boundary.size()
80  );
81 
82  label nCoarseFaces = 0;
83 
84  forAllConstIter(dictionary, agglomDict, iter)
85  {
86  labelList patchids = boundary.findIndices(iter().keyword());
87  forAll(patchids, i)
88  {
89  label patchi = patchids[i];
90  const polyPatch& pp = boundary[patchi];
91  if (!pp.coupled())
92  {
93  Info << "\nAgglomerating patch : " << pp.name() << endl;
94  pairPatchAgglomeration agglomObject
95  (
96  pp,
97  agglomDict.subDict(pp.name())
98  );
99  agglomObject.agglomerate();
100  finalAgglom[patchi] =
101  agglomObject.restrictTopBottomAddressing();
102 
103  if (finalAgglom[patchi].size())
104  {
105  nCoarseFaces += max(finalAgglom[patchi] + 1);
106  }
107  }
108  }
109  }
110 
111 
112  // - All patches which are not agglomerated are identity for finalAgglom
113  forAll(boundary, patchid)
114  {
115  if (finalAgglom[patchid].size() == 0)
116  {
117  finalAgglom[patchid] = identity(boundary[patchid].size());
118  }
119  }
120 
121  // Sync agglomeration across coupled patches
122  labelList nbrAgglom(mesh.nFaces() - mesh.nInternalFaces(), -1);
123 
124  forAll(boundary, patchid)
125  {
126  const polyPatch& pp = boundary[patchid];
127  if (pp.coupled())
128  {
129  finalAgglom[patchid] = identity(pp.size());
130  forAll(pp, i)
131  {
132  nbrAgglom[pp.start() - mesh.nInternalFaces() + i] =
133  finalAgglom[patchid][i];
134  }
135  }
136  }
137 
139  forAll(boundary, patchid)
140  {
141  const polyPatch& pp = boundary[patchid];
142  if (pp.coupled() && !refCast<const coupledPolyPatch>(pp).owner())
143  {
144  forAll(pp, i)
145  {
146  finalAgglom[patchid][i] =
147  nbrAgglom[pp.start() - mesh.nInternalFaces() + i];
148  }
149  }
150  }
151 
152  finalAgglom.write();
153 
154  if (writeAgglom)
155  {
156  globalIndex index(nCoarseFaces);
157  volScalarField facesAgglomeration
158  (
159  IOobject
160  (
161  "facesAgglomeration",
162  mesh.time().timeName(),
163  mesh,
166  ),
167  mesh,
169  );
170 
171  volScalarField::Boundary& facesAgglomerationBf =
172  facesAgglomeration.boundaryFieldRef();
173 
174  label coarsePatchIndex = 0;
175  forAll(boundary, patchid)
176  {
177  const polyPatch& pp = boundary[patchid];
178  if (pp.size() > 0)
179  {
180  fvPatchScalarField& bFacesAgglomeration =
181  facesAgglomerationBf[patchid];
182 
183  forAll(bFacesAgglomeration, j)
184  {
185  bFacesAgglomeration[j] =
186  index.toGlobal
187  (
189  finalAgglom[patchid][j] + coarsePatchIndex
190  );
191  }
192 
193  coarsePatchIndex += max(finalAgglom[patchid]) + 1;
194  }
195  }
196 
197  Info<< "\nWriting facesAgglomeration" << endl;
198  facesAgglomeration.write();
199  }
200 
201  Info<< "End\n" << endl;
202  return 0;
203 }
204 
205 
206 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
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
const word & name() const
Return name.
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:808
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
label nInternalFaces() const
label nFaces() const
Unit conversion functions.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const dimensionSet dimless
bool readBool(Istream &)
Definition: boolIO.C:60
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
dynamicFvMesh & mesh
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:319
A class for handling words, derived from string.
Definition: word.H:59
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Foam::polyBoundaryMesh.
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Primitive patch pair agglomerate method.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:309
messageStream Info
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const word dictName("noiseDict")
Namespace for OpenFOAM.