1 /*---------------------------------------------------------------------------*\
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 -------------------------------------------------------------------------------
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/>.
25 Foam::coarseMgFieldLevel
28 Coarse level of multigrid hierarchy, holding fields
33 \*---------------------------------------------------------------------------*/
35 #ifndef coarseMgFieldLevel_H
36 #define coarseMgFieldLevel_H
38 #include "mgFieldLevel.H"
39 #include "mgMeshLevel.H"
40 #include "coarseMgMeshLevel.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 /*---------------------------------------------------------------------------*\
48 Class coarseMgFieldLevel Declaration
49 \*---------------------------------------------------------------------------*/
51 class coarseMgFieldLevel
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_;
72 //- Number of internal face
73 label nInternalFaces_;
88 // Thermodynamics fields
103 vectorField rhoURes_;
106 scalarField rhoERes_;
112 surfaceScalarField& rhoFlux_;
115 surfaceVectorField& rhoUFlux_;
118 surfaceScalarField& rhoEFlux_;
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>
150 c = pTraits<T>::zero;
151 f *= meshLevel_.fineLevel().cellVolumes();
152 fineFieldLevel_.restrict(c, f);
153 c /= meshLevel_.cellVolumes();
159 //- Runtime type information
160 TypeName("coarseMgFieldLevel");
165 //- Construct from fine mgFieldLevel
168 const mgMeshLevel& coarseMeshLevel,
169 const mgFieldLevel& fineFieldLevel
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())
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());
207 const scalarField& fineVol =
208 meshLevel_.fineLevel().cellVolumes();
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);
221 scalarField fineRhoERes = fineFieldLevel_.rhoEResVar();
222 fineRhoERes *= fineVol;
223 fineFieldLevel_.restrict(rhoERes_, fineRhoERes);
228 fineFieldLevel_.restrict
231 fineFieldLevel_.rhoResVar()
234 rhoURes_ = vector::zero;
235 fineFieldLevel_.restrict
238 fineFieldLevel_.rhoUResVar()
242 fineFieldLevel_.restrict
245 fineFieldLevel_.rhoEResVar()
252 virtual ~coarseMgFieldLevel()
260 //- Is this the finest level?
261 virtual bool finest() const
263 return meshLevel_.finest();
269 //- Return number of cells
270 virtual label nCells() const
275 //- Return number of internal faces
276 virtual label nInternalFaces() const
278 return nInternalFaces_;
281 //- Return number of patches
282 virtual label nPatches() const
284 return fineFieldLevel_.nPatches();
287 //- Access to child array
288 virtual const labelList& child() const
290 return meshLevel_.child();
293 //- Access to fine level
294 virtual const mgFieldLevel& fineLevel() const
296 return fineFieldLevel_;
303 virtual const scalarField& pVar() const
308 //- Access to fine level pressure field
309 virtual const volScalarField& p() const
311 return fineFieldLevel_.p();
314 //- Return p value for patches on fine level
315 virtual const scalarField& patchP
320 return fineFieldLevel_.patchP(patchNo);
324 virtual const vectorField& UVar() const
329 virtual const volVectorField& U() const
331 return fineFieldLevel_.U();
334 //- Return U value for patches
335 virtual const vectorField& patchU
340 return fineFieldLevel_.patchU(patchNo);
344 virtual const scalarField& TVar() const
349 virtual const volScalarField& T() const
351 return fineFieldLevel_.T();
354 //- Return T value for patches on fine level
355 virtual const scalarField& patchT
360 return fineFieldLevel_.patchT(patchNo);
364 virtual const scalarField& CvVar() const
369 virtual const volScalarField& Cv() const
371 return fineFieldLevel_.Cv();
374 //- Return Cv value for patches
375 virtual const scalarField& patchCv
380 return fineFieldLevel_.patchCv(patchNo);
384 virtual const scalarField& RVar() const
389 virtual const volScalarField& R() const
391 return fineFieldLevel_.R();
394 //- Return R value for patches on fine level
395 virtual const scalarField& patchR
400 return fineFieldLevel_.patchR(patchNo);
404 // Fluxes on fine level
406 virtual surfaceScalarField& rhoFlux() const
408 return fineFieldLevel_.rhoFlux();
411 virtual surfaceVectorField& rhoUFlux() const
413 return fineFieldLevel_.rhoUFlux();
416 virtual surfaceScalarField& rhoEFlux() const
418 return fineFieldLevel_.rhoEFlux();
424 virtual const scalarField& rhoResVar() const
429 virtual const volScalarField& rhoRes() const
431 return fineFieldLevel_.rhoRes();
434 //- Return rhoRes value for patches on fine level
435 virtual const scalarField& patchRhoRes
440 return fineFieldLevel_.patchRhoRes(patchNo);
444 virtual const vectorField& rhoUResVar() const
449 virtual const volVectorField& rhoURes() const
451 return fineFieldLevel_.rhoURes();
454 //- Return rhoRes value for patches on fine level
455 virtual const vectorField& patchRhoURes
460 return fineFieldLevel_.patchRhoURes(patchNo);
464 virtual const scalarField& rhoEResVar() const
469 virtual const volScalarField& rhoERes() const
471 return fineFieldLevel_.rhoERes();
474 //- Return rhoERes value for patches on fine level
475 virtual const scalarField& patchRhoERes
480 return fineFieldLevel_.patchRhoERes(patchNo);
484 // Conservative variables
487 virtual const scalarField& rhoVar() const
492 virtual const volScalarField& rho() const
494 return fineFieldLevel_.rho();
497 //- Return rho value for patches on fine level
498 virtual const scalarField& patchRho
503 return fineFieldLevel_.patchRho(patchNo);
507 virtual const vectorField& rhoUVar() const
512 virtual const volVectorField& rhoU() const
514 return fineFieldLevel_.rhoU();
517 //- Return rhoU value for patches on fine level
518 virtual const vectorField& patchRhoU
523 return fineFieldLevel_.patchRhoU(patchNo);
527 virtual const scalarField& rhoEVar() const
532 virtual const volScalarField& rhoE() const
534 return fineFieldLevel_.rhoE();
537 //- Return rho value for patches on fine level
538 virtual const scalarField& patchRhoE
543 return fineFieldLevel_.patchRhoE(patchNo);
546 //- Access to helper data
547 virtual const label& level() const
554 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
556 } // End namespace Foam
558 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
562 // ************************************************************************* //