BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / applications / utilities / preProcessing / faceAgglomerate / faceAgglomerate.C
blob206c150e77e4396246d728e2492549ab9f28df91
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Application
25     faceAgglomerate
27 Description
29     Agglomerate boundary faces using the pairPatchAgglomeration algorithm.
30     It writes a map from the fine to coarse grid.
32 \*---------------------------------------------------------------------------*/
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"
44 using namespace Foam;
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 int main(int argc, char *argv[])
50     #include "addRegionOption.H"
51     argList::addOption
52     (
53         "dict",
54         "word",
55         "name of dictionary to provide patch agglomeration controls"
56     );
57     #include "setRootCase.H"
58     #include "createTime.H"
59     #include "createNamedMesh.H"
61     word agglomDictName
62     (
63         args.optionLookupOrDefault<word>("dict", "faceAgglomerateDict")
64     );
66     const polyBoundaryMesh& patches = mesh.boundaryMesh();
68     labelListIOList finalAgglom
69     (
70         IOobject
71         (
72             "finalAgglom",
73             mesh.facesInstance(),
74             mesh,
75             IOobject::NO_READ,
76             IOobject::NO_WRITE,
77             false
78         ),
79         patches.size()
80     );
83     // Read control dictionary
84     IOdictionary agglomDict
85     (
86         IOobject
87         (
88             agglomDictName,
89             runTime.constant(),
90             mesh,
91             IOobject::MUST_READ,
92             IOobject::NO_WRITE
93         )
94     );
96     bool writeAgglom = readBool(agglomDict.lookup("writeFacesAgglomeration"));
98     const polyBoundaryMesh& boundary = mesh.boundaryMesh();
100     forAll(boundary, patchId)
101     {
102         const polyPatch& pp = boundary[patchId];
104         label patchI = pp.index();
105         finalAgglom[patchI].setSize(pp.size(), 0);
107         if (!pp.coupled())
108         {
109             if (agglomDict.found(pp.name()))
110             {
111                 Info << "\nAgglomerating patch : " << pp.name() << endl;
112                 pairPatchAgglomeration agglomObject
113                 (
114                     pp,
115                     agglomDict.subDict(pp.name())
116                 );
118                 agglomObject.agglomerate();
120                 finalAgglom[patchI] =
121                     agglomObject.restrictTopBottomAddressing();
123             }
124             else
125             {
126                 FatalErrorIn(args.executable())
127                     << "Patch " << pp.name() << " not found in dictionary: "
128                     << agglomDict.name() << exit(FatalError);
129             }
130         }
131     }
133     // Sync agglomeration across coupled patches
134     labelList nbrAgglom(mesh.nFaces() - mesh.nInternalFaces(), -1);
136     forAll(boundary, patchId)
137     {
138         const polyPatch& pp = boundary[patchId];
139         if (pp.coupled())
140         {
141             finalAgglom[patchId] = identity(pp.size());
142             forAll(pp, i)
143             {
144                 nbrAgglom[pp.start() - mesh.nInternalFaces() + i] =
145                     finalAgglom[patchId][i];
146             }
147         }
148     }
150     syncTools::swapBoundaryFaceList(mesh, nbrAgglom);
151     forAll(boundary, patchId)
152     {
153         const polyPatch& pp = boundary[patchId];
154         if (pp.coupled() && !refCast<const coupledPolyPatch>(pp).owner())
155         {
156             forAll(pp, i)
157             {
158                 finalAgglom[patchId][i] =
159                     nbrAgglom[pp.start() - mesh.nInternalFaces() + i];
160             }
161         }
162     }
164     finalAgglom.write();
166     if (writeAgglom)
167     {
168         volScalarField facesAgglomeration
169         (
170             IOobject
171             (
172                 "facesAgglomeration",
173                 mesh.time().timeName(),
174                 mesh,
175                 IOobject::NO_READ,
176                 IOobject::NO_WRITE
177             ),
178             mesh,
179             dimensionedScalar("facesAgglomeration", dimless, 0)
180         );
182         forAll(boundary, patchId)
183         {
185             fvPatchScalarField& bFacesAgglomeration =
186                 facesAgglomeration.boundaryField()[patchId];
188             forAll(bFacesAgglomeration, j)
189             {
190                 bFacesAgglomeration[j] = finalAgglom[patchId][j];
191             }
192         }
194         Info << "\nWriting  facesAgglomeration" << endl;
195         facesAgglomeration.write();
196     }
198     Info<< "End\n" << endl;
199     return 0;
203 // ************************************************************************* //