MGridGenGAMGAgglomerate.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-2013 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Description
27  Agglomerate one level using the MGridGen algorithm.
28 
29 \*---------------------------------------------------------------------------*/
30 
32 #include "fvMesh.H"
33 #include "syncTools.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 void Foam::MGridGenGAMGAgglomeration::
38 makeCompactCellFaceAddressingAndFaceWeights
39 (
40  const lduAddressing& fineAddressing,
41  List<idxtype>& cellCells,
42  List<idxtype>& cellCellOffsets,
43  const scalarField& magSi,
44  List<scalar>& faceWeights
45 )
46 {
47  const label nFineCells = fineAddressing.size();
48  const label nFineFaces = fineAddressing.upperAddr().size();
49 
50  const labelUList& upperAddr = fineAddressing.upperAddr();
51  const labelUList& lowerAddr = fineAddressing.lowerAddr();
52 
53  // Number of neighbours for each cell
54  labelList nNbrs(nFineCells, Zero);
55 
56  forAll(upperAddr, facei)
57  {
58  nNbrs[upperAddr[facei]]++;
59  }
60 
61  forAll(lowerAddr, facei)
62  {
63  nNbrs[lowerAddr[facei]]++;
64  }
65 
66  // Set the sizes of the addressing and faceWeights arrays
67  cellCellOffsets.setSize(nFineCells + 1);
68  cellCells.setSize(2*nFineFaces);
69  faceWeights.setSize(2*nFineFaces);
70 
71 
72  cellCellOffsets[0] = 0;
73  forAll(nNbrs, celli)
74  {
75  cellCellOffsets[celli+1] = cellCellOffsets[celli] + nNbrs[celli];
76  }
77 
78  // reset the whole list to use as counter
79  nNbrs = 0;
80 
81  forAll(upperAddr, facei)
82  {
83  label own = upperAddr[facei];
84  label nei = lowerAddr[facei];
85 
86  label l1 = cellCellOffsets[own] + nNbrs[own]++;
87  label l2 = cellCellOffsets[nei] + nNbrs[nei]++;
88 
89  cellCells[l1] = nei;
90  cellCells[l2] = own;
91 
92  faceWeights[l1] = magSi[facei];
93  faceWeights[l2] = magSi[facei];
94  }
95 }
96 
97 
98 Foam::tmp<Foam::labelField> Foam::MGridGenGAMGAgglomeration::agglomerate
99 (
100  label& nCoarseCells,
101  const label minSize,
102  const label maxSize,
103  const lduAddressing& fineAddressing,
104  const scalarField& V,
105  const scalarField& magSf,
106  const scalarField& magSb
107 )
108 {
109  const label nFineCells = fineAddressing.size();
110 
111  // Compact addressing for cellCells
112  List<idxtype> cellCells;
113  List<idxtype> cellCellOffsets;
114 
115  // Face weights = face areas of the internal faces
116  List<scalar> faceWeights;
117 
118  // Create the compact addressing for cellCells and faceWeights
119  makeCompactCellFaceAddressingAndFaceWeights
120  (
121  fineAddressing,
122  cellCells,
123  cellCellOffsets,
124  magSf,
125  faceWeights
126  );
127 
128  // agglomeration options.
129  List<int> options(4, Zero);
130  options[0] = 4; // globular agglom
131  options[1] = 6; // objective F3 and F2
132  options[2] = 128; // debugging output level
133  options[3] = fvMesh_.nGeometricD(); // Dimensionality of the grid
134 
135 
136  // output: cell -> processor addressing
137  List<int> finalAgglom(nFineCells);
138  int nMoves = -1;
139 
140  MGridGen
141  (
142  nFineCells,
143  cellCellOffsets.begin(),
144  const_cast<scalar*>(V.begin()),
145  const_cast<scalar*>(magSb.begin()),
146  cellCells.begin(),
147  faceWeights.begin(),
148  minSize,
149  maxSize,
150  options.begin(),
151  &nMoves,
152  &nCoarseCells,
153  finalAgglom.begin()
154  );
155 
156  {
157  label nNewCoarseCells = 0;
158  labelList newRestrictAddr;
159  bool ok = checkRestriction
160  (
161  newRestrictAddr,
162  nNewCoarseCells,
163  fineAddressing,
164  finalAgglom,
165  nCoarseCells
166  );
167 
168  if (!ok)
169  {
170  nCoarseCells = nNewCoarseCells;
171  finalAgglom.transfer(newRestrictAddr);
172  }
173  }
174 
175  return tmp<labelField>::New(finalAgglom);
176 }
177 
178 
179 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
MGridGen
void MGridGen(int, int *, Foam::scalar *, Foam::scalar *, int *, Foam::scalar *, int, int, int *, int *, int *, int *)
Definition: dummyMGridGen.C:48
syncTools.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
fvMesh.H
MGridGenGAMGAgglomeration.H
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85