Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / dbns / multigrid / mgMeshLevel / mgMeshLevel.H
blob3c72e564709e0f223a8dbf62aa981d3ee810a1fa
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::mgMeshLevel
27 Description
28     Virtual base class for a level of multigrid hierarchy.
30     Cells and internal faces are coarsened by clustering
31     There is no clustering on the boundary: original faces are used and
32     face-cell patch addressing on coarse levels returns addressing into
33     coarse cell clusters.
35 SourceFiles
36     mgMeshLevel.C
38 \*---------------------------------------------------------------------------*/
40 #ifndef mgMeshLevel_H
41 #define mgMeshLevel_H
43 #include "typeInfo.H"
44 #include "fvMesh.H"
45 #include "primitiveFields.H"
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 namespace Foam
52 /*---------------------------------------------------------------------------*\
53                           Class mgMeshLevel Declaration
54 \*---------------------------------------------------------------------------*/
56 class mgMeshLevel
58     // Private data
60         //- Number of coarse cells
61         mutable label nCoarseCells_;
63         //- Child array: coarse cluster index for every fine cell
64         mutable labelList child_;
67     // Private Member Functions
69         //- Construct the CSR format addressing
70         void makeCompactAddressingAndWeights
71         (
72             labelList& cellCells,
73             labelList& cellCellOffsets,
74             scalarField& faceWeights,
75             scalarField& boundaryAreas
76         ) const;
78         //- Create child array
79         void makeChild() const;
82 public:
84     // Static Data
86         //- Runtime type information
87         TypeName("mgMeshLevel");
90         // Coarsening parameters
92             //- Minimum cluster size
93             static const debug::optimisationSwitch mgMinClusterSize_;
95             //- Maximum cluster size
96             static const debug::optimisationSwitch mgMaxClusterSize_;
99     // Constructor
101         //- Construct null
102         mgMeshLevel()
103         :
104             nCoarseCells_(-1),
105             child_()
106         {}
109     //- Destructor
111         virtual ~mgMeshLevel()
112         {}
115     // Member Functions
117         //- Is this the finest level?
118         virtual bool finest() const = 0;
120         //- Return reference to fine level
121         virtual const mgMeshLevel& fineLevel() const = 0;
123         //- Return reference to finest level
124         const mgMeshLevel& finestLevel() const;
126         //- Return access to mesh
127         virtual const fvMesh& mesh() const = 0;
130         // Sizes
132             //- Return number of cells
133             virtual label nCells() const = 0;
135             //- Return number of internal faces
136             virtual label nInternalFaces() const = 0;
138             //- Return number of patches
139             virtual label nPatches() const = 0;
142         // Addressing
144             //- Return owner addressing
145             virtual const unallocLabelList& owner() const = 0;
147             //- Return neighbour addressing
148             virtual const unallocLabelList& neighbour() const = 0;
150             //- Return face-cell addressing for all patches
151             virtual const unallocLabelList& faceCells
152             (
153                 const label patchNo
154             ) const = 0;
157         // Geometrical data
159             //- Return number of valid geometric dimensions in the mesh
160             virtual label nGeometricD() const = 0;
162             //- Return cell centres
163             virtual const vectorField& cellCentres() const = 0;
165             //- Return face centres
166             virtual const vectorField& faceCentres() const = 0;
168             //- Return cell volumes
169             virtual const scalarField& cellVolumes() const = 0;
171             //- Return face area vectors
172             virtual const vectorField& faceAreas() const = 0;
174             //- Return face area magnitudes
175             virtual const scalarField& magFaceAreas() const = 0;
177             //- Return face area vectors for patches
178             //  Note: no clustering of boundary faces allowed:
179             //        using fine-level geometry
180             virtual const vectorField& patchFaceAreas
181             (
182                 const label patchNo
183             ) const = 0;
185             //- Return face area magnitudes for patches
186             //  Note: no clustering of boundary faces allowed:
187             //        using fine-level geometry
188             virtual const scalarField& magPatchFaceAreas
189             (
190                 const label patchNo
191             ) const = 0;
193             //- Return cell to face vectors for patches
194             //  Note: no clustering of boundary faces allowed:
195             //        using fine-level geometry
196             virtual tmp<vectorField> patchDeltaR
197             (
198                 const label patchNo
199             ) const = 0;
201             //- Return neighbour cell to face vectors for coupled patches
202             //  Note: no clustering of boundary faces allowed:
203             //        using fine-level geometry
204             virtual tmp<vectorField> patchDeltaRNeighbour
205             (
206                 const label patchNo
207             ) const = 0;
210         // Coarsening
212             //- Return number of cells
213             label nCoarseCells() const;
215             //- Return child array
216             const labelList& child() const;
218             //- Create next level from current level
219             autoPtr<mgMeshLevel> makeNextLevel() const;
222          // Restriction and prolongation
224              //- Restrict (from fine to coarse).  Call from the level
225              //  where data originates
226              template<class Type>
227              void restrict
228              (
229                  UList<Type>& coarseData,
230                  const UList<Type>& fineData
231              ) const
232              {
233                  if
234                  (
235                      coarseData.size() != nCoarseCells()
236                   || fineData.size() != nCells()
237                  )
238                  {
239                      FatalErrorIn
240                      (
241                          "void mgMeshLevel::restrict\n"
242                          "(\n"
243                          "    UList<Type>& coarseData,\n"
244                          "    const UList<Type>& fineData\n"
245                          ") const\n"
246                      )   << "Incorrect data sizes ("
247                          << coarseData.size() << " " << fineData.size()
248                          << ").  Should be ("
249                          << nCoarseCells() << " " << nCells() << ")"
250                          << abort(FatalError);
251                  }
253                  // Get addressing
254                  const labelList& c = child();
256                  // Restriction = summation.  Set coarse data to zero
257                  coarseData = pTraits<Type>::zero;
259                  forAll (c, i)
260                  {
261                      coarseData[c[i]] += fineData[i];
262                  }
263              }
265              //- Prolong (from  coarse to fine).  Call from the level
266              //  where data originates
267              template<class Type>
268              void prolong
269              (
270                  const UList<Type>& coarseData,
271                  UList<Type>& fineData
272              ) const
273              {
274                  if (finest())
275                  {
276                      FatalErrorIn
277                      (
278                          "void mgMeshLevel::prolong\n"
279                          "(\n"
280                          "    UList<Type>& coarseData,\n"
281                          "    const UList<Type>& fineData\n"
282                          ") const\n"
283                      )   << "Requested prolong from finest level"
284                          << abort(FatalError);
285                  }
287                  if
288                  (
289                      coarseData.size() != nCells()
290                   || fineData.size() != fineLevel().nCells()
291                  )
292                  {
293                      FatalErrorIn
294                      (
295                          "void mgMeshLevel::prolong\n"
296                          "(\n"
297                          "    UList<Type>& coarseData,\n"
298                          "    const UList<Type>& fineData\n"
299                          ") const\n"
300                      )   << "Incorrect data sizes ("
301                          << coarseData.size() << " " << fineData.size()
302                          << ").  Should be ("
303                          << nCoarseCells() << " " << nCells() << ")"
304                          << abort(FatalError);
305                  }
307                  // Get addressing
308                  const labelList& c = fineLevel().child();
310                  // Prolongation = assignment.  Initialisation not needed
311                  forAll (c, i)
312                  {
313                      fineData[i] = coarseData[c[i]];
314                  }
315              }
317              //- Prolong to the finest level (from  coarse to fine)
318              template<class Type>
319              tmp<Field<Type> > prolongToFinest
320              (
321                  const UList<Type>& coarseData
322              ) const
323              {
324                  if (finest())
325                  {
326                      // Already on finest level
327                      return tmp<Field<Type> >(new Field<Type>(coarseData));
328                  }
329                  else if (fineLevel().finest())
330                  {
331                      tmp<Field<Type> > tresult
332                      (
333                          new Field<Type>(fineLevel().nCells())
334                      );
336                      prolong(coarseData, tresult());
338                      // Reached finest level, return
339                      return tresult;
340                  }
341                  else
342                  {
343                      Field<Type> intermediate(fineLevel().nCells());
345                      prolong(coarseData, intermediate);
347                      // Recursive return
348                      return fineLevel().prolongToFinest(intermediate);
349                  }
350              }
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
356 } // End namespace Foam
358 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
360 #endif
362 // ************************************************************************* //