initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / fvAgglomerationMethods / MGridGenGamgAgglomeration / MGridGenGAMGAgglomerate.C
blob4af895ca4cce9de0f0657cd230198c8285daa835
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Description
26     Agglomerate one level using the MGridGen algorithm.
28 \*---------------------------------------------------------------------------*/
30 #include "MGridGenGAMGAgglomeration.H"
31 #include "fvMesh.H"
32 #include "syncTools.H"
34 //extern "C"
35 //{
36 //#   include "mgridgen.h"
37 //}
39 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
41 void Foam::MGridGenGAMGAgglomeration::
42 makeCompactCellFaceAddressingAndFaceWeights
44     const lduAddressing& fineAddressing,
45     List<idxtype>& cellCells,
46     List<idxtype>& cellCellOffsets,
47     const vectorField& Si,
48     List<scalar>& faceWeights
51     const label nFineCells = fineAddressing.size();
52     const label nFineFaces = fineAddressing.upperAddr().size();
54     const unallocLabelList& upperAddr = fineAddressing.upperAddr();
55     const unallocLabelList& lowerAddr = fineAddressing.lowerAddr();
57     // Number of neighbours for each cell
58     labelList nNbrs(nFineCells, 0);
60     forAll (upperAddr, facei)
61     {
62         nNbrs[upperAddr[facei]]++;
63     }
65     forAll (lowerAddr, facei)
66     {
67         nNbrs[lowerAddr[facei]]++;
68     }
70     // Set the sizes of the addressing and faceWeights arrays
71     cellCellOffsets.setSize(nFineCells + 1);
72     cellCells.setSize(2*nFineFaces);
73     faceWeights.setSize(2*nFineFaces);
76     cellCellOffsets[0] = 0;
77     forAll (nNbrs, celli)
78     {
79         cellCellOffsets[celli+1] = cellCellOffsets[celli] + nNbrs[celli];
80     }
82     // reset the whole list to use as counter
83     nNbrs = 0;
85     forAll (upperAddr, facei)
86     {
87         label own = upperAddr[facei];
88         label nei = lowerAddr[facei];
90         label l1 = cellCellOffsets[own] + nNbrs[own]++;
91         label l2 = cellCellOffsets[nei] + nNbrs[nei]++;
93         cellCells[l1] = nei;
94         cellCells[l2] = own;
96         faceWeights[l1] = mag(Si[facei]);
97         faceWeights[l2] = mag(Si[facei]);
98     }
102 Foam::tmp<Foam::labelField> Foam::MGridGenGAMGAgglomeration::agglomerate
104     label& nCoarseCells,
105     const label minSize,
106     const label maxSize,
107     const lduAddressing& fineAddressing,
108     const scalarField& V,
109     const vectorField& Sf,
110     const scalarField& Sb
113     const label nFineCells = fineAddressing.size();
115     // Compact addressing for cellCells
116     List<idxtype> cellCells;
117     List<idxtype> cellCellOffsets;
119     // Face weights = face areas of the internal faces
120     List<scalar> faceWeights;
122     // Create the compact addressing for cellCells and faceWeights
123     makeCompactCellFaceAddressingAndFaceWeights
124     (
125         fineAddressing,
126         cellCells,
127         cellCellOffsets,
128         Sf,
129         faceWeights
130     );
132     // agglomeration options.
133     List<int> options(4, 0);
134     options[0] = 4;                   // globular agglom
135     options[1] = 6;                   // objective F3 and F2
136     options[2] = 128;                 // debugging output level
137     options[3] = fvMesh_.nGeometricD(); // Dimensionality of the grid
140     // output: cell -> processor addressing
141     List<int> finalAgglom(nFineCells);
142     int nMoves = -1;
143         
144     MGridGen
145     (
146         nFineCells,
147         cellCellOffsets.begin(),
148         const_cast<scalar*>(V.begin()),
149         const_cast<scalar*>(Sb.begin()),
150         cellCells.begin(),
151         faceWeights.begin(),
152         minSize,
153         maxSize,
154         options.begin(),
155         &nMoves,
156         &nCoarseCells,
157         finalAgglom.begin()
158     );
160     return tmp<labelField>(new labelField(finalAgglom));
164 // ************************************************************************* //