CMULES.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2013-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Global
25  MULES
26 
27 Description
28  CMULES: Multidimensional universal limiter for
29  explicit corrected implicit solution.
30 
31  Solve a convective-only transport equation using an explicit universal
32  multi-dimensional limiter to correct an implicit conservative bounded
33  obtained using rigorously bounded schemes such as Euler-implicit in time
34  upwind in space.
35 
36  Parameters are the variable to solve, the normal convective flux and the
37  actual explicit flux of the variable which is also used to return limited
38  flux used in the bounded-solution.
39 
40 SourceFiles
41  CMULES.C
42  CMULESTemplates.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef CMULES_H
47 #define CMULES_H
48 
49 #include "MULES.H"
50 #include "EulerDdtScheme.H"
51 #include "localEulerDdtScheme.H"
52 #include "gaussConvectionScheme.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 namespace MULES
59 {
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 template<class RdeltaTType, class RhoType, class SpType, class SuType>
64 void correct
65 (
66  const RdeltaTType& rDeltaT,
67  const RhoType& rho,
69  const surfaceScalarField& phi,
70  const surfaceScalarField& phiCorr,
71  const SpType& Sp,
72  const SuType& Su
73 );
74 
75 template<class RhoType, class SpType, class SuType>
76 void correct
77 (
78  const RhoType& rho,
79  volScalarField& psi,
80  const surfaceScalarField& phi,
81  surfaceScalarField& phiCorr,
82  const SpType& Sp,
83  const SuType& Su,
84  const scalar psiMax,
85  const scalar psiMin
86 );
87 
88 void correct
89 (
90  volScalarField& psi,
91  const surfaceScalarField& phi,
92  surfaceScalarField& phiCorr,
93  const scalar psiMax,
94  const scalar psiMin
95 );
96 
97 template<class RdeltaTType, class RhoType, class SpType, class SuType>
98 void limiterCorr
99 (
100  scalarField& allLambda,
101  const RdeltaTType& rDeltaT,
102  const RhoType& rho,
103  const volScalarField& psi,
104  const surfaceScalarField& phi,
105  const surfaceScalarField& phiCorr,
106  const SpType& Sp,
107  const SuType& Su,
108  const scalar psiMax,
109  const scalar psiMin
110 );
111 
112 template<class RdeltaTType, class RhoType, class SpType, class SuType>
113 void limitCorr
114 (
115  const RdeltaTType& rDeltaT,
116  const RhoType& rho,
117  const volScalarField& psi,
118  const surfaceScalarField& phi,
119  surfaceScalarField& phiCorr,
120  const SpType& Sp,
121  const SuType& Su,
122  const scalar psiMax,
123  const scalar psiMin
124 );
125 
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 } // End namespace MULES
130 } // End namespace Foam
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 #ifdef NoRepository
135  #include "CMULESTemplates.C"
136 #endif
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 #endif
141 
142 // ************************************************************************* //
surfaceScalarField & phi
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
void limitCorr(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
MULES: Multidimensional universal limiter for explicit solution.
const volScalarField & psi
void limiterCorr(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.