faceAgglomerate.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-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2020,2022 OpenCFD Ltd.
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 Application
28  faceAgglomerate
29 
30 Group
31  grpPreProcessingUtilities
32 
33 Description
34  Agglomerate boundary faces using the pairPatchAgglomeration algorithm.
35 
36  It writes a map from the fine to coarse grid.
37 
38 SeeAlso
39  pairPatchAgglomeration.H
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "argList.H"
44 #include "fvMesh.H"
45 #include "Time.H"
46 #include "volFields.H"
47 #include "unitConversion.H"
48 #include "pairPatchAgglomeration.H"
49 #include "labelListIOList.H"
50 #include "syncTools.H"
51 #include "globalIndex.H"
52 
53 using namespace Foam;
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 int main(int argc, char *argv[])
58 {
60  (
61  "Agglomerate boundary faces using the pairPatchAgglomeration"
62  " algorithm. Writes a map of fine to coarse grid."
63  );
64 
65  argList::addOption("dict", "file", "Alternative viewFactorsDict");
66  #include "addRegionOption.H"
67 
68  #include "setRootCase.H"
69  #include "createTime.H"
70  #include "createNamedMesh.H"
71 
72  const word dictName("viewFactorsDict");
73 
75 
76  // Read control dictionary
77  const IOdictionary agglomDict(dictIO);
78 
79  const bool writeAgglom(agglomDict.get<bool>("writeFacesAgglomeration"));
80 
82 
83  labelListIOList finalAgglom
84  (
85  IOobject
86  (
87  "finalAgglom",
89  mesh,
92  false
93  ),
94  boundary.size()
95  );
96 
97  label nCoarseFaces = 0;
98 
99 
100  const auto& patchesDict =
101  agglomDict.optionalSubDict
102  (
103  "patchAgglomeration",
105  );
106 
107  for (const entry& dEntry : patchesDict)
108  {
109  labelList patchids = boundary.indices(dEntry.keyword());
110 
111  for (const label patchi : patchids)
112  {
113  const polyPatch& pp = boundary[patchi];
114 
115  if (!pp.coupled())
116  {
117  Info << "\nAgglomerating patch : " << pp.name() << endl;
118 
119  pairPatchAgglomeration agglomObject
120  (
121  pp.localFaces(),
122  pp.localPoints(),
123  dEntry.dict()
124  );
125 
126  agglomObject.agglomerate();
127 
128  finalAgglom[patchi] =
129  agglomObject.restrictTopBottomAddressing();
130 
131  if (finalAgglom[patchi].size())
132  {
133  nCoarseFaces += max(finalAgglom[patchi] + 1);
134  }
135  }
136  }
137  }
138 
139 
140  // All patches which are not agglomerated are identity for finalAgglom
141  forAll(boundary, patchi)
142  {
143  if (finalAgglom[patchi].size() == 0)
144  {
145  finalAgglom[patchi] = identity(boundary[patchi].size());
146  }
147  }
148 
149  // Sync agglomeration across coupled patches
150  labelList nbrAgglom(mesh.nBoundaryFaces(), -1);
151 
152  forAll(boundary, patchi)
153  {
154  const polyPatch& pp = boundary[patchi];
155  if (pp.coupled())
156  {
157  finalAgglom[patchi] = identity(pp.size());
158  forAll(pp, i)
159  {
160  const label agglomi = pp.start() - mesh.nInternalFaces() + i;
161  nbrAgglom[agglomi] = finalAgglom[patchi][i];
162  }
163  }
164  }
165 
167  forAll(boundary, patchi)
168  {
169  const polyPatch& pp = boundary[patchi];
170  if (pp.coupled() && !refCast<const coupledPolyPatch>(pp).owner())
171  {
172  forAll(pp, i)
173  {
174  const label agglomi = pp.start() - mesh.nInternalFaces() + i;
175  finalAgglom[patchi][i] = nbrAgglom[agglomi];
176  }
177  }
178  }
179 
180  finalAgglom.write();
181 
182  if (writeAgglom)
183  {
184  globalIndex index(nCoarseFaces);
185  volScalarField facesAgglomeration
186  (
187  IOobject
188  (
189  "facesAgglomeration",
190  mesh.time().timeName(),
191  mesh,
194  ),
195  mesh,
197  );
198 
199  volScalarField::Boundary& facesAgglomerationBf =
200  facesAgglomeration.boundaryFieldRef();
201 
202  label coarsePatchIndex = 0;
203  forAll(boundary, patchi)
204  {
205  const polyPatch& pp = boundary[patchi];
206  if (pp.size() > 0)
207  {
208  fvPatchScalarField& bFacesAgglomeration =
209  facesAgglomerationBf[patchi];
210 
211  forAll(bFacesAgglomeration, j)
212  {
213  bFacesAgglomeration[j] =
214  index.toGlobal
215  (
217  finalAgglom[patchi][j] + coarsePatchIndex
218  );
219  }
220 
221  coarsePatchIndex += max(finalAgglom[patchi]) + 1;
222  }
223  }
224 
225  Info<< "\nWriting facesAgglomeration" << endl;
226  facesAgglomeration.write();
227  }
228 
229  Info<< "End\n" << endl;
230  return 0;
231 }
232 
233 
234 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
faceListList boundary
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:853
const Field< point_type > & localPoints() const
Return pointField of points in patch.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Unit conversion functions.
const word dictName("faMeshDefinition")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
Required Variables.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:688
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
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:464
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:31
void agglomerate()
Agglomerate patch.
A class for handling words, derived from Foam::string.
Definition: word.H:63
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:376
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
label nInternalFaces() const noexcept
Number of internal faces.
String literal.
Definition: keyType.H:82
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const word & name() const noexcept
The patch name.
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
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:441
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
IOobject dictIO
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
Namespace for OpenFOAM.
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:485
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157