Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / dbns / multigrid / mgFieldLevel / coarseMgFieldLevel.H
blob25a40857d798345bcc2ed73fcbe782ba0be12584
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::coarseMgFieldLevel
27 Description
28     Coarse level of multigrid hierarchy, holding fields
30 SourceFiles
31     coarseMgFieldLevel.C
33 \*---------------------------------------------------------------------------*/
35 #ifndef coarseMgFieldLevel_H
36 #define coarseMgFieldLevel_H
38 #include "mgFieldLevel.H"
39 #include "mgMeshLevel.H"
40 #include "coarseMgMeshLevel.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 namespace Foam
47 /*---------------------------------------------------------------------------*\
48                        Class coarseMgFieldLevel Declaration
49 \*---------------------------------------------------------------------------*/
51 class coarseMgFieldLevel
53     public mgFieldLevel
55     // Private data
57         // Mesh levels
59             //- Reference to the current mesh level
60             //  corresponding to the current field level
61             const mgMeshLevel& meshLevel_;
63             //- Reference to fine data level
64            const mgFieldLevel& fineFieldLevel_;
67         // Geometry
69             //- Number of cells
70             label nCells_;
72             //- Number of internal face
73             label nInternalFaces_;
76         // Fields
78             //- Coarse p
79             scalarField p_;
81             //- Coarse U
82             vectorField U_;
84             //- Coarse T
85             scalarField T_;
88         // Thermodynamics fields
90             //- Coarse Cv
91             scalarField Cv_;
93             //- Coarse R
94             scalarField R_;
97         // Residuals
99             //- Coarse rhoRes
100             scalarField rhoRes_;
102             //- Coarse rhoURes
103             vectorField rhoURes_;
105             //- Coarse rhoERes
106             scalarField rhoERes_;
109         // Fluxes
111             //- Fine rhoFlux
112             surfaceScalarField& rhoFlux_;
114             //- Fine rhoUFlux
115             surfaceVectorField& rhoUFlux_;
117             //- Fine rhoEFlux
118             surfaceScalarField& rhoEFlux_;
120             //- Coarse rho
121             scalarField rho_;
123             //- Coarse rhoU
124             vectorField rhoU_;
126             //- Coarse rhoE
127             scalarField rhoE_;
130         //- Helper data
131         label level_;
134     // Private Member Functions
136         //- Disallow default bitwise copy construct
137         coarseMgFieldLevel(const coarseMgFieldLevel&);
139         //- Disallow default bitwise assignment
140         void operator=(const coarseMgFieldLevel&);
142         // Restriction operators used for interpolation
143        template <typename T>
144        void interpolate
145        (
146            Field<T>& c,
147            Field<T> f
148        )
149        {
150            c = pTraits<T>::zero;
151            f *= meshLevel_.fineLevel().cellVolumes();
152            fineFieldLevel_.restrict(c, f);
153            c /= meshLevel_.cellVolumes();
154        }
157 public:
159         //- Runtime type information
160         TypeName("coarseMgFieldLevel");
163     // Constructors
165         //- Construct from fine mgFieldLevel
166         coarseMgFieldLevel
167         (
168             const mgMeshLevel& coarseMeshLevel,
169             const mgFieldLevel& fineFieldLevel
170         )
171         :
172             meshLevel_(coarseMeshLevel),
173             fineFieldLevel_(fineFieldLevel),
174             nCells_(coarseMeshLevel.nCells()),
175             nInternalFaces_(coarseMeshLevel.nInternalFaces()),
176             p_(scalarField(coarseMeshLevel.nCells())),
177             U_(vectorField(coarseMeshLevel.nCells())),
178             T_(scalarField(coarseMeshLevel.nCells())),
179             Cv_(scalarField(coarseMeshLevel.nCells())),
180             R_(scalarField(coarseMeshLevel.nCells())),
181             rhoRes_(scalarField(coarseMeshLevel.nCells())),
182             rhoURes_(vectorField(coarseMeshLevel.nCells())),
183             rhoERes_(scalarField(coarseMeshLevel.nCells())),
184             rhoFlux_(fineFieldLevel.rhoFlux()),
185             rhoUFlux_(fineFieldLevel.rhoUFlux()),
186             rhoEFlux_(fineFieldLevel.rhoEFlux()),
187             rho_(scalarField(coarseMeshLevel.nCells())),
188             rhoU_(vectorField(coarseMeshLevel.nCells())),
189             rhoE_(scalarField(coarseMeshLevel.nCells())),
190             level_(fineFieldLevel.level())
191         {
192             ++level_;
194             interpolate(p_, fineFieldLevel_.pVar());
195             interpolate(U_, fineFieldLevel_.UVar());
196             interpolate(T_, fineFieldLevel_.TVar());
198             interpolate(Cv_, fineFieldLevel_.CvVar());
199             interpolate(R_, fineFieldLevel_.RVar());
201             interpolate(rho_, fineFieldLevel_.rhoVar());
202             interpolate(rhoU_, fineFieldLevel_.rhoUVar());
203             interpolate(rhoE_, fineFieldLevel_.rhoEVar());
205             if (level_ == 1)
206             {
207                 const scalarField& fineVol =
208                     meshLevel_.fineLevel().cellVolumes();
210                 rhoRes_ = 0;
211                 scalarField fineRhoRes = fineFieldLevel_.rhoResVar();
212                 fineRhoRes *= fineVol;
213                 fineFieldLevel_.restrict(rhoRes_, fineRhoRes);
215                 rhoURes_ = vector::zero;
216                 vectorField fineRhoURes = fineFieldLevel_.rhoUResVar();
217                 fineRhoURes *= fineVol;
218                 fineFieldLevel_.restrict(rhoURes_, fineRhoURes);
220                 rhoERes_ = 0;
221                 scalarField fineRhoERes = fineFieldLevel_.rhoEResVar();
222                 fineRhoERes *= fineVol;
223                 fineFieldLevel_.restrict(rhoERes_, fineRhoERes);
224             }
225             else
226             {
227                 rhoRes_ = 0;
228                 fineFieldLevel_.restrict
229                 (
230                     rhoRes_,
231                     fineFieldLevel_.rhoResVar()
232                 );
234                 rhoURes_ = vector::zero;
235                 fineFieldLevel_.restrict
236                 (
237                      rhoURes_,
238                      fineFieldLevel_.rhoUResVar()
239                 );
241                 rhoERes_ = 0;
242                 fineFieldLevel_.restrict
243                 (
244                     rhoERes_,
245                     fineFieldLevel_.rhoEResVar()
246                 );
247             }
248         }
251     //- Destructor
252     virtual ~coarseMgFieldLevel()
253     {}
256     // Member Functions
258         // Access
260             //- Is this the finest level?
261             virtual bool finest() const
262             {
263                 return meshLevel_.finest();
264             }
267         // Addressing
269             //- Return number of cells
270             virtual label nCells() const
271             {
272                 return nCells_;
273             }
275             //- Return number of internal faces
276             virtual label nInternalFaces() const
277             {
278                 return nInternalFaces_;
279             }
281             //- Return number of patches
282             virtual label nPatches() const
283             {
284                 return fineFieldLevel_.nPatches();
285             }
287             //- Access to child array
288             virtual const labelList& child() const
289             {
290                 return meshLevel_.child();
291             }
293             //- Access to fine level
294             virtual const mgFieldLevel& fineLevel() const
295             {
296                 return fineFieldLevel_;
297             }
300         // Field data
302             //- Return p
303             virtual const scalarField& pVar() const
304             {
305                 return p_;
306             }
308             //- Access to fine level pressure field
309             virtual const volScalarField& p() const
310             {
311                 return fineFieldLevel_.p();
312             }
314             //- Return p value for patches on fine level
315             virtual const scalarField& patchP
316             (
317                 const label patchNo
318             ) const
319             {
320                 return fineFieldLevel_.patchP(patchNo);
321             }
323              //- Return U
324             virtual const vectorField& UVar() const
325             {
326                 return U_;
327             }
329             virtual const volVectorField& U() const
330             {
331                 return fineFieldLevel_.U();
332             }
334             //- Return U value for patches
335             virtual const vectorField& patchU
336             (
337                 const label patchNo
338             ) const
339             {
340                 return fineFieldLevel_.patchU(patchNo);
341             }
343             //- Return T
344             virtual const scalarField& TVar() const
345             {
346                 return T_;
347             }
349             virtual const volScalarField& T() const
350             {
351                 return fineFieldLevel_.T();
352             }
354             //- Return T value for patches on fine level
355             virtual const scalarField& patchT
356             (
357                 const label patchNo
358             ) const
359             {
360                 return fineFieldLevel_.patchT(patchNo);
361             }
363             //- Return Cv
364             virtual const scalarField& CvVar() const
365             {
366                 return Cv_;
367             }
369             virtual const volScalarField& Cv() const
370             {
371                 return fineFieldLevel_.Cv();
372             }
374             //- Return Cv value for patches
375             virtual const scalarField& patchCv
376             (
377                 const label patchNo
378             ) const
379             {
380                 return fineFieldLevel_.patchCv(patchNo);
381             }
383             //- Return R
384             virtual const scalarField& RVar() const
385             {
386                 return R_;
387             }
389             virtual const volScalarField& R() const
390             {
391                 return fineFieldLevel_.R();
392             }
394             //- Return R value for patches on fine level
395             virtual const scalarField& patchR
396             (
397                 const label patchNo
398             ) const
399             {
400                 return fineFieldLevel_.patchR(patchNo);
401             }
404         // Fluxes on fine level
406             virtual surfaceScalarField& rhoFlux() const
407             {
408                 return fineFieldLevel_.rhoFlux();
409             }
411             virtual surfaceVectorField& rhoUFlux() const
412             {
413                 return fineFieldLevel_.rhoUFlux();
414             }
416             virtual surfaceScalarField& rhoEFlux() const
417             {
418                 return fineFieldLevel_.rhoEFlux();
419             }
421         // Residuals
423            //- Return rhoRes
424             virtual const scalarField& rhoResVar() const
425             {
426                 return rhoRes_;
427             }
429             virtual const volScalarField& rhoRes() const
430             {
431                 return fineFieldLevel_.rhoRes();
432             }
434             //- Return rhoRes value for patches on fine level
435             virtual const scalarField& patchRhoRes
436             (
437                 const label patchNo
438             ) const
439             {
440                 return fineFieldLevel_.patchRhoRes(patchNo);
441             }
443             //- Return rhoURes
444             virtual const vectorField& rhoUResVar() const
445             {
446                 return rhoURes_;
447             }
449             virtual const volVectorField& rhoURes() const
450             {
451                 return fineFieldLevel_.rhoURes();
452             }
454             //- Return rhoRes value for patches on fine level
455             virtual const vectorField& patchRhoURes
456             (
457                 const label patchNo
458             ) const
459             {
460                 return fineFieldLevel_.patchRhoURes(patchNo);
461             }
463             //- Return rhoERes
464             virtual const scalarField& rhoEResVar() const
465             {
466                 return rhoERes_;
467             }
469             virtual const volScalarField& rhoERes() const
470             {
471                 return fineFieldLevel_.rhoERes();
472             }
474             //- Return rhoERes value for patches on fine level
475             virtual const scalarField& patchRhoERes
476             (
477                 const label patchNo
478             ) const
479             {
480                 return fineFieldLevel_.patchRhoERes(patchNo);
481             }
484         // Conservative variables
486             //- Return rho
487             virtual const scalarField& rhoVar() const
488             {
489                 return rho_;
490             }
492             virtual const volScalarField& rho() const
493             {
494                 return fineFieldLevel_.rho();
495             }
497             //- Return rho value for patches on fine level
498             virtual const scalarField& patchRho
499             (
500                 const label patchNo
501             ) const
502             {
503                 return fineFieldLevel_.patchRho(patchNo);
504             }
506             //- Return rhoU
507             virtual const vectorField& rhoUVar() const
508             {
509                 return rhoU_;
510             }
512             virtual const volVectorField& rhoU() const
513             {
514                 return fineFieldLevel_.rhoU();
515             }
517             //- Return rhoU value for patches on fine level
518             virtual const vectorField& patchRhoU
519             (
520                 const label patchNo
521             ) const
522             {
523                 return fineFieldLevel_.patchRhoU(patchNo);
524             }
526             //- Return rhoE
527             virtual const scalarField& rhoEVar() const
528             {
529                 return rhoE_;
530             }
532             virtual const volScalarField& rhoE() const
533             {
534                 return fineFieldLevel_.rhoE();
535             }
537             //- Return rho value for patches on fine level
538             virtual const scalarField& patchRhoE
539             (
540                 const label patchNo
541             ) const
542             {
543                 return fineFieldLevel_.patchRhoE(patchNo);
544             }
546             //- Access to helper data
547             virtual const label& level() const
548             {
549                 return level_;
550             }
554 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
556 } // End namespace Foam
558 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
560 #endif
562 // ************************************************************************* //