Forward compatibility: flex
[foam-extend-3.2.git] / src / dbns / multigrid / mgFieldLevel / mgFieldLevel.H
blobb45a70f1d7c46eca96aa102e5fe512990e2843f0
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::mgFieldLevel
27 Description
28     Virtual base class for a level of multigrid hierarchy.
30 SourceFiles
31     mgFieldLevel.C
33 \*---------------------------------------------------------------------------*/
35 #ifndef mgFieldLevel_H
36 #define mgFieldLevel_H
38 #include "typeInfo.H"
39 #include "fvMesh.H"
40 #include "primitiveFields.H"
41 #include "coarseMgMeshLevel.H"
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 namespace Foam
48 /*---------------------------------------------------------------------------*\
49                           Class mgFieldLevel Declaration
50 \*---------------------------------------------------------------------------*/
52 class mgFieldLevel
54 public:
56         //- Runtime type information
57         virtual const word& type() const = 0;
59         //- Access to child array
60         virtual const labelList& child() const = 0;
62         //- Access to fine level
63         virtual const mgFieldLevel& fineLevel() const = 0;
65         //- Is this the finest level?
66         virtual bool finest() const = 0;
69     //- Destructor
71         virtual ~mgFieldLevel()
72         {}
75     // Member Functions
78         //- Access
80             //- Return number of cells
81             virtual label nCells() const = 0;
83             //- Return number of internal faces
84             virtual label nInternalFaces() const = 0;
86             //- Return number of patches
87             virtual label nPatches() const = 0;
90         // Primitive variables
92             //- Access to p field
93             virtual const scalarField& pVar() const = 0;
95             //- Access to U field
96             virtual const vectorField& UVar() const = 0;
98             //- Access to T field
99             virtual const scalarField& TVar() const = 0;
101             //- Access to fine p field
102             virtual const volScalarField& p() const = 0;
104             //- Access to fine U field
105             virtual const volVectorField& U() const = 0;
107             //- Access to fine T field
108             virtual const volScalarField& T() const = 0;
110             //- Access to p boundary field
111             virtual const scalarField& patchP(const label patchNo) const = 0;
113             //- Access to U field
114             virtual const vectorField& patchU(const label patchNo) const = 0;
116             //- Access to T field
117             virtual const scalarField& patchT(const label patchNo) const = 0;
119             //- Access to Cv field
120             virtual const scalarField& CvVar() const = 0;
122             //- Access to fine Cv field
123             virtual const volScalarField& Cv() const = 0;
125             //- Access to Cv boundary field
126             virtual const scalarField& patchCv(const label patchNo) const = 0;
128             //- Access to R field
129             virtual const scalarField& RVar() const = 0;
131             //- Access to fine R field
132             virtual const volScalarField& R() const = 0;
134             //- Access to R boundary field
135             virtual const scalarField& patchR(const label patchNo) const = 0;
137         // Fine fluxes
139             //- Access to fine rhoFlux field
140             virtual surfaceScalarField& rhoFlux() const = 0;
142             //- Access to fine rhoUFlux field
143             virtual surfaceVectorField& rhoUFlux() const = 0;
145             //- Access to fine rhoEFlux field
146             virtual surfaceScalarField& rhoEFlux() const = 0;
149         // Residuals
151             //- Access to rhoRes field
152             virtual const scalarField& rhoResVar() const = 0;
154             //- Access to fine rhoRes field
155             virtual const volScalarField& rhoRes() const = 0;
157             //- Access to rhoRes boundary field
158             virtual const scalarField& patchRhoRes
159             (
160                 const label patchNo
161             ) const = 0;
163             //- Access to rhoURes field
164             virtual const vectorField& rhoUResVar() const = 0;
166             //- Access to fine rhoRes field
167             virtual const volVectorField& rhoURes() const = 0;
169             //- Access to rhoRes boundary field
170             virtual const vectorField& patchRhoURes
171             (
172                 const label patchNo
173             ) const = 0;
175             //- Access to rhoERes field
176             virtual const scalarField& rhoEResVar() const = 0;
178             //- Access to fine rhoERes field
179             virtual const volScalarField& rhoERes() const = 0;
181             //- Access to rhoERes boundary field
182             virtual const scalarField& patchRhoERes
183             (
184                 const label patchNo
185             ) const = 0;
188         // Conservative variables
190             //- Access to rho field
191             virtual const scalarField& rhoVar() const = 0;
193             //- Access to fine rho field
194             virtual const volScalarField& rho() const = 0;
196             //- Access to rho boundary field
197             virtual const scalarField& patchRho
198             (
199                 const label patchNo
200             ) const = 0;
202             //- Access to rhoURes field
203             virtual const vectorField& rhoUVar() const = 0;
205             //- Access to fine rhoRes field
206             virtual const volVectorField& rhoU() const = 0;
208             //- Access to rhoRes boundary field
209             virtual const vectorField& patchRhoU
210             (
211                 const label patchNo
212             ) const = 0;
214             //- Access to rhoERes field
215             virtual const scalarField& rhoEVar() const = 0;
217             //- Access to fine rhoERes field
218             virtual const volScalarField& rhoE() const = 0;
220             //- Access to rhoERes boundary field
221             virtual const scalarField& patchRhoE
222             (
223                 const label patchNo
224             ) const = 0;
226             //- Access to rhoERes boundary field
227             virtual label const& level() const = 0;
230         // Restriction and prolongation
232             //- Restrict (from fine to coarse).  Call from the level
233             //  where data originates
234             template<class T>
235             void restrict
236             (
237                 UList<T>& coarseData,
238                 const UList<T>& fineData
239             ) const
240             {
241                 // Get addressing
242                 const labelList& c = child();
244                 // Restriction = summation.  Set coarse data to zero
245                 coarseData = pTraits<T>::zero;
247                 forAll (c, i)
248                 {
249                     coarseData[c[i]] += fineData[i];
250                 }
251             }
253             //- Prolong (from  coarse to fine).  Call from the level
254             //  where data originates
255             template<class T>
256             void prolong
257             (
258                 const UList<T>& coarseData,
259                 UList<T>& fineData
260             ) const
261             {
262                 if (finest())
263                 {
264                     FatalErrorIn
265                     (
266                         "void mgMeshLevel::prolong\n"
267                         "(\n"
268                         "    UList<T>& coarseData,\n"
269                         "    const UList<T>& fineData\n"
270                         ") const\n"
271                     )   << "Requested prolong from finest level"
272                         << abort(FatalError);
273                 }
275                 if
276                 (
277                     coarseData.size() != nCells()
278                  || fineData.size() != fineLevel().nCells()
279                 )
280                 {
281                     FatalErrorIn
282                     (
283                         "void mgMeshLevel::prolong\n"
284                         "(\n"
285                         "    UList<T>& coarseData,\n"
286                         "    const UList<T>& fineData\n"
287                         ") const\n"
288                     )    << "Incorrect data sizes "
289                          << abort(FatalError);
290                 }
292                 // Get addressing
293                 const labelList& c = fineLevel().child();
295                 // Prolongation = assignment.  Initialisation not needed
296                 forAll (c, i)
297                 {
298                     fineData[i] = coarseData[c[i]];
299                 }
300             }
302             //- Prolong to the finest level (from  coarse to fine)
303             template<class T>
304             tmp<Field<T> > prolongToFinest
305             (
306                 const UList<T>& coarseData
307             ) const
308             {
309                 if (finest())
310                 {
311                     // Already on finest level
312                     return tmp<Field<T> >(new Field<T>(coarseData));
313                 }
314                 else if (fineLevel().finest())
315                 {
316                     tmp<Field<T> > tresult
317                     (
318                         new Field<T>(fineLevel().nCells())
319                     );
321                     prolong(coarseData, tresult());
323                     // Reached finest level, return
324                     return tresult;
325                 }
326                 else
327                 {
328                     Field<T> intermediate(fineLevel().nCells());
330                     prolong(coarseData, intermediate);
332                     // Recursive return
333                     return fineLevel().prolongToFinest(intermediate);
334                 }
335             }
337         // Coarsening
339             //- Create next level from current level
340             autoPtr<mgFieldLevel> makeNextLevel
341             (
342                 const mgMeshLevel& coarseMeshLevel
343             ) const;
348 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
350 } // End namespace Foam
352 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
354 #endif
356 // ************************************************************************* //