CMULES.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2025 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>
64 void correct
65 (
66  const RdeltaTType& rDeltaT,
67  const RhoType& rho,
69  const surfaceScalarField& phiCorr,
70  const SpType& Sp
71 );
72 
73 template<class RhoType>
74 void correct
75 (
76  const RhoType& rho,
78  const surfaceScalarField& phiCorr
79 );
80 
81 template<class RhoType, class SpType>
82 void correct
83 (
84  const RhoType& rho,
86  const surfaceScalarField& phiCorr,
87  const SpType& Sp
88 );
89 
90 template<class RhoType, class PsiMaxType, class PsiMinType>
91 void correct
92 (
93  const control& controls,
94  const RhoType& rho,
96  const surfaceScalarField& phiBD,
97  surfaceScalarField& phiCorr,
98  const PsiMaxType& psiMax,
99  const PsiMinType& psiMin
100 );
101 
102 template
103 <
104  class RhoType,
105  class SpType,
106  class PsiMaxType,
107  class PsiMinType
108 >
109 void correct
110 (
111  const control& controls,
112  const RhoType& rho,
114  const surfaceScalarField& phiBD,
115  surfaceScalarField& phiCorr,
116  const SpType& Sp,
117  const PsiMaxType& psiMax,
118  const PsiMinType& psiMin
119 );
120 
121 template
122 <
123  class RdeltaTType,
124  class RhoType,
125  class SpType,
126  class PsiMaxType,
127  class PsiMinType
128 >
129 void limitCorr
130 (
131  const control& controls,
132  const RdeltaTType& rDeltaT,
133  const RhoType& rho,
134  const volScalarField& psi,
135  const surfaceScalarField& phiBD,
136  surfaceScalarField& phiCorr,
137  const SpType& Sp,
138  const PsiMaxType& psiMax,
139  const PsiMinType& psiMin
140 );
141 
142 template
143 <
144  class RhoType,
145  class SpType,
146  class PsiMaxType,
147  class PsiMinType
148 >
149 void limitCorr
150 (
151  const control& controls,
152  const RhoType& rho,
153  const volScalarField& psi,
154  const surfaceScalarField& phiBD,
155  surfaceScalarField& phiCorr,
156  const SpType& Sp,
157  const PsiMaxType& psiMax,
158  const PsiMinType& psiMin
159 );
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 } // End namespace MULES
165 } // End namespace Foam
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #ifdef NoRepository
170  #include "CMULESTemplates.C"
171 #endif
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 #endif
176 
177 // ************************************************************************* //
MULES: Multidimensional universal limiter for explicit solution.
Generic GeometricField class.
const volScalarField & psi
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
void limitCorr(const control &controls, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, surfaceScalarField &phiCorr, const SpType &Sp, const PsiMaxType &psiMax, const PsiMinType &psiMin)
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
Namespace for OpenFOAM.