Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / utilities / preProcessing / faceAgglomerate / faceAgglomerate.C
blob387b791592a7f036f8a33598b84b60eb1a065bf1
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011-2011 OpenCFD Ltd.
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.size() > 0 && !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                 );
117                 agglomObject.agglomerate();
118                 finalAgglom[patchI] =
119                     agglomObject.restrictTopBottomAddressing();
120             }
121             else
122             {
123                 FatalErrorIn(args.executable())
124                     << "Patch " << pp.name() << " not found in dictionary: "
125                     << agglomDict.name() << exit(FatalError);
126             }
127         }
128     }
130     // Sync agglomeration across coupled patches
131     labelList nbrAgglom(mesh.nFaces() - mesh.nInternalFaces(), -1);
133     forAll(boundary, patchId)
134     {
135         const polyPatch& pp = boundary[patchId];
136         if (pp.coupled())
137         {
138             finalAgglom[patchId] = identity(pp.size());
139             forAll(pp, i)
140             {
141                 nbrAgglom[pp.start() - mesh.nInternalFaces() + i] =
142                     finalAgglom[patchId][i];
143             }
144         }
145     }
147     syncTools::swapBoundaryFaceList(mesh, nbrAgglom);
148     forAll(boundary, patchId)
149     {
150         const polyPatch& pp = boundary[patchId];
151         if (pp.coupled() && !refCast<const coupledPolyPatch>(pp).owner())
152         {
153             forAll(pp, i)
154             {
155                 finalAgglom[patchId][i] =
156                     nbrAgglom[pp.start() - mesh.nInternalFaces() + i];
157             }
158         }
159     }
161     finalAgglom.write();
163     if (writeAgglom)
164     {
165         volScalarField facesAgglomeration
166         (
167             IOobject
168             (
169                 "facesAgglomeration",
170                 mesh.time().timeName(),
171                 mesh,
172                 IOobject::NO_READ,
173                 IOobject::NO_WRITE
174             ),
175             mesh,
176             dimensionedScalar("facesAgglomeration", dimless, 0)
177         );
179         forAll(boundary, patchId)
180         {
182             fvPatchScalarField& bFacesAgglomeration =
183                 facesAgglomeration.boundaryField()[patchId];
185             forAll(bFacesAgglomeration, j)
186             {
187                 bFacesAgglomeration[j] = finalAgglom[patchId][j];
188             }
189         }
191         Info << "\nWriting  facesAgglomeration" << endl;
192         facesAgglomeration.write();
193     }
195     Info<< "End\n" << endl;
196     return 0;
200 // ************************************************************************* //