MULES.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) 2011-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  MULES: Multidimensional universal limiter for explicit solution.
29 
30  Solve a convective-only transport equation using an explicit universal
31  multi-dimensional limiter.
32 
33  Parameters are the variable to solve, the normal convective flux and the
34  actual explicit flux of the variable which is also used to return limited
35  flux used in the bounded-solution.
36 
37 SourceFiles
38  MULES.C
39  MULESTemplates.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef MULES_H
44 #define MULES_H
45 
46 #include "Switch.H"
47 #include "volFieldsFwd.H"
48 #include "surfaceFieldsFwd.H"
49 #include "primitiveFieldsFwd.H"
50 #include "geometricOneField.H"
51 #include "zero.H"
52 #include "zeroField.H"
53 #include "UPtrList.H"
54 #include "HashSet.H"
55 #include "UniformField.H"
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 namespace Foam
60 {
61 namespace MULES
62 {
63 
64 NamespaceName("MULES");
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 // MULES control structure
69 struct control
70 {
71  //- Optional maximum number of limiter iterations
72  // Defaults to 3
73  label nIter;
74 
75  //- Optional limiter convergence tolerance
76  // 1e-2 is usually sufficient
77  // but 1e-3 can be used for very tight convergence
78  // Defaults to 0, i.e. unused and nIter iterations are performed
79  scalar tol;
80 
81  //- Optional switch to select global bounds only
82  // rather than local and global bounds
84 
85  //- Optional coefficient to relax the local boundedness constraint
86  // by adding extremaCoeff*local_range to the allowed range
87  // allowing the solution to approach the global bounds in case temporary
88  // smearing has caused the solution to deviate from the entire range.
89  //
90  // A value of 0 enforces strict local boundedness,
91  // a value of 1 enforces the global bounds only.
92  // Defaults to 0
93  scalar extremaCoeff;
94 
95  //- Optional coefficient to relax the local boundedness constraint
96  // in the cells adjacent to boundaries only
97  // Defaults to 0
98  scalar boundaryExtremaCoeff;
99 
100  //- Optional coefficient to reduce the allowed range of the solution to
101  // smoothingCoeff*local_range.
102  //
103  // This option should only be used to improve the smoothness of the
104  // solution at low Courant number (explicit solution with sub-cycling)
105  // otherwise excessive smearing may result.
106  scalar smoothingCoeff;
107 
108  //- Null constructor
109  // Should be followed by a call to read(dict)
110  control()
111  {}
112 
113  //- Construct from dict and set the controls
114  control(const dictionary& dict);
115 
116  //- Read dict and set the controls
117  void read(const dictionary& dict);
118 };
119 
120 template<class RdeltaTType, class RhoType, class SpType, class SuType>
121 void explicitSolve
122 (
123  const RdeltaTType& rDeltaT,
124  const RhoType& rho,
126  const surfaceScalarField& psiPhi,
127  const SpType& Sp,
128  const SuType& Su
129 );
130 
131 template<class RhoType>
132 void explicitSolve
133 (
134  const RhoType& rho,
136  const surfaceScalarField& psiPhi
137 );
138 
139 template<class RhoType, class SpType, class SuType>
140 void explicitSolve
141 (
142  const RhoType& rho,
144  const surfaceScalarField& psiPhi,
145  const SpType& Sp,
146  const SuType& Su
147 );
148 
149 template<class RhoType, class PsiMaxType, class PsiMinType>
150 void explicitSolve
151 (
152  const control& controls,
153  const RhoType& rho,
155  const surfaceScalarField& phiBD,
156  surfaceScalarField& psiPhi,
157  const PsiMaxType& psiMax,
158  const PsiMinType& psiMin
159 );
160 
161 template
162 <
163  class RhoType,
164  class SpType,
165  class SuType,
166  class PsiMaxType,
167  class PsiMinType
168 >
169 void explicitSolve
170 (
171  const control& controls,
172  const RhoType& rho,
174  const surfaceScalarField& phiBD,
175  surfaceScalarField& psiPhi,
176  const SpType& Sp,
177  const SuType& Su,
178  const PsiMaxType& psiMax,
179  const PsiMinType& psiMin
180 );
181 
182 template
183 <
184  class RdeltaTType,
185  class RhoType,
186  class SpType,
187  class PsiMaxType,
188  class PsiMinType
189 >
190 void limiter
191 (
192  const control& controls,
194  const RdeltaTType& rDeltaT,
195  const RhoType& rho,
196  const volScalarField& psi,
197  const scalarField& SuCorr,
198  const surfaceScalarField& phi,
199  const surfaceScalarField& phiCorr,
200  const SpType& Sp,
201  const PsiMaxType& psiMax,
202  const PsiMinType& psiMin
203 );
204 
205 template
206 <
207  class RdeltaTType,
208  class RhoType,
209  class SpType,
210  class SuType,
211  class PsiMaxType,
212  class PsiMinType
213 >
214 void limit
215 (
216  const control& controls,
217  const RdeltaTType& rDeltaT,
218  const RhoType& rho,
219  const volScalarField& psi,
220  const surfaceScalarField& phi,
221  surfaceScalarField& psiPhi,
222  const SpType& Sp,
223  const SuType& Su,
224  const PsiMaxType& psiMax,
225  const PsiMinType& psiMin,
226  const bool returnCorr
227 );
228 
229 template
230 <
231  class RdeltaTType,
232  class RhoType,
233  class SpType,
234  class SuType,
235  class PsiMaxType,
236  class PsiMinType
237 >
238 void limit
239 (
240  const control& controls,
241  const RdeltaTType& rDeltaT,
242  const RhoType& rho,
243  const volScalarField& psi,
244  const surfaceScalarField& phi,
245  const surfaceScalarField& phiBD,
246  surfaceScalarField& psiPhi,
247  const SpType& Sp,
248  const SuType& Su,
249  const PsiMaxType& psiMax,
250  const PsiMinType& psiMin,
251  const bool returnCorr
252 );
253 
254 template
255 <
256  class RhoType,
257  class SpType,
258  class SuType,
259  class PsiMaxType,
260  class PsiMinType
261 >
262 void limit
263 (
264  const control& controls,
265  const RhoType& rho,
266  const volScalarField& psi,
267  const surfaceScalarField& phi,
268  surfaceScalarField& psiPhi,
269  const SpType& Sp,
270  const SuType& Su,
271  const PsiMaxType& psiMax,
272  const PsiMinType& psiMin,
273  const bool returnCorr
274 );
275 
276 template
277 <
278  class RhoType,
279  class SpType,
280  class SuType,
281  class PsiMaxType,
282  class PsiMinType
283 >
284 void limit
285 (
286  const control& controls,
287  const RhoType& rho,
288  const volScalarField& psi,
289  const surfaceScalarField& phi,
290  const surfaceScalarField& phiBD,
291  surfaceScalarField& psiPhi,
292  const SpType& Sp,
293  const SuType& Su,
294  const PsiMaxType& psiMax,
295  const PsiMinType& psiMin,
296  const bool returnCorr
297 );
298 
299 void limitSumCorr(UPtrList<scalarField>& psiPhiCorrs);
300 
301 void limitSumCorr
302 (
303  UPtrList<scalarField>& phiPsiCorrs,
304  const UPtrList<scalarField>& phiPsiFixedCorrs
305 );
306 
307 void limitSumCorr
308 (
309  const UPtrList<const volScalarField>& psis,
310  UPtrList<surfaceScalarField>& psiPhiCorrs,
311  const surfaceScalarField& phi
312 );
313 
314 void limitSum
315 (
316  const UPtrList<const volScalarField>& psis,
317  const PtrList<surfaceScalarField>& alphaPhiBDs,
319  const surfaceScalarField& phi
320 );
321 
322 void limitSum
323 (
324  const UPtrList<const volScalarField>& psis,
326  const surfaceScalarField& phi
327 );
328 
329 
330 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331 
332 } // End namespace MULES
333 } // End namespace Foam
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 
337 #ifdef NoRepository
338  #include "MULESlimiter.C"
339  #include "MULESTemplates.C"
340 #endif
341 
342 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
343 
344 #endif
345 
346 // ************************************************************************* //
Generic GeometricField class.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const volScalarField & psi
dimensionedScalar lambda(viscosity->lookup("lambda"))
void limitSum(const UPtrList< const volScalarField > &psis, const PtrList< surfaceScalarField > &alphaPhiBDs, UPtrList< surfaceScalarField > &psiPhis, const surfaceScalarField &phi)
Definition: MULES.C:272
void limitSumCorr(UPtrList< scalarField > &psiPhiCorrs)
Definition: MULES.C:78
void limiter(const control &controls, surfaceScalarField &lambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const scalarField &SuCorr, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const PsiMaxType &psiMax, const PsiMinType &psiMin)
Definition: MULESlimiter.C:42
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &psiPhi, const SpType &Sp, const SuType &Su)
void limit(const control &controls, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &psiPhi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
NamespaceName("MULES")
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
dictionary dict
Switch globalBounds
Optional switch to select global bounds only.
Definition: MULES.H:80
void read(const dictionary &dict)
Read dict and set the controls.
Definition: MULES.C:44
scalar smoothingCoeff
Optional coefficient to reduce the allowed range of the solution to.
Definition: MULES.H:103
scalar boundaryExtremaCoeff
Optional coefficient to relax the local boundedness constraint.
Definition: MULES.H:95
control()
Null constructor.
Definition: MULES.H:107
label nIter
Optional maximum number of limiter iterations.
Definition: MULES.H:70
scalar extremaCoeff
Optional coefficient to relax the local boundedness constraint.
Definition: MULES.H:90
scalar tol
Optional limiter convergence tolerance.
Definition: MULES.H:76