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-2018 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 "volFieldsFwd.H"
47 #include "surfaceFieldsFwd.H"
48 #include "primitiveFieldsFwd.H"
49 #include "geometricOneField.H"
50 #include "zero.H"
51 #include "zeroField.H"
52 #include "UPtrList.H"
53 #include "HashSet.H"
54 #include "UniformField.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 namespace MULES
61 {
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 template<class RdeltaTType, class RhoType, class SpType, class SuType>
66 void explicitSolve
67 (
68  const RdeltaTType& rDeltaT,
69  const RhoType& rho,
71  const surfaceScalarField& phiPsi,
72  const SpType& Sp,
73  const SuType& Su
74 );
75 
76 template<class RhoType>
77 void explicitSolve
78 (
79  const RhoType& rho,
81  const surfaceScalarField& phiPsi
82 );
83 
84 template<class RhoType, class SpType, class SuType>
85 void explicitSolve
86 (
87  const RhoType& rho,
89  const surfaceScalarField& phiPsi,
90  const SpType& Sp,
91  const SuType& Su
92 );
93 
94 template<class RhoType, class PsiMaxType, class PsiMinType>
95 void explicitSolve
96 (
97  const RhoType& rho,
99  const surfaceScalarField& phiBD,
100  surfaceScalarField& phiPsi,
101  const PsiMaxType& psiMax,
102  const PsiMinType& psiMin
103 );
104 
105 template
106 <
107  class RhoType,
108  class SpType,
109  class SuType,
110  class PsiMaxType,
111  class PsiMinType
112 >
113 void explicitSolve
114 (
115  const RhoType& rho,
117  const surfaceScalarField& phiBD,
118  surfaceScalarField& phiPsi,
119  const SpType& Sp,
120  const SuType& Su,
121  const PsiMaxType& psiMax,
122  const PsiMinType& psiMin
123 );
124 
125 template
126 <
127  class RdeltaTType,
128  class RhoType,
129  class SpType,
130  class SuType,
131  class PsiMaxType,
132  class PsiMinType
133 >
134 void limiter
135 (
136  scalarField& allLambda,
137  const RdeltaTType& rDeltaT,
138  const RhoType& rho,
139  const volScalarField& psi,
140  const surfaceScalarField& phiBD,
141  const surfaceScalarField& phiCorr,
142  const SpType& Sp,
143  const SuType& Su,
144  const PsiMaxType& psiMax,
145  const PsiMinType& psiMin
146 );
147 
148 template
149 <
150  class RdeltaTType,
151  class RhoType,
152  class SpType,
153  class SuType,
154  class PsiMaxType,
155  class PsiMinType
156 >
157 void limit
158 (
159  const RdeltaTType& rDeltaT,
160  const RhoType& rho,
161  const volScalarField& psi,
162  const surfaceScalarField& phi,
163  surfaceScalarField& phiPsi,
164  const SpType& Sp,
165  const SuType& Su,
166  const PsiMaxType& psiMax,
167  const PsiMinType& psiMin,
168  const bool returnCorr
169 );
170 
171 template
172 <
173  class RhoType,
174  class SpType,
175  class SuType,
176  class PsiMaxType,
177  class PsiMinType
178 >
179 void limit
180 (
181  const RhoType& rho,
182  const volScalarField& psi,
183  const surfaceScalarField& phi,
184  surfaceScalarField& phiPsi,
185  const SpType& Sp,
186  const SuType& Su,
187  const PsiMaxType& psiMax,
188  const PsiMinType& psiMin,
189  const bool returnCorr
190 );
191 
192 void limitSum(UPtrList<scalarField>& phiPsiCorrs);
193 
194 template<class SurfaceScalarFieldList>
195 void limitSum(SurfaceScalarFieldList& phiPsiCorrs);
196 
197 void limitSum
198 (
199  const UPtrList<const scalarField>& alphas,
200  UPtrList<scalarField>& phiPsiCorrs,
201  const labelHashSet& fixed
202 );
203 
204 template<class SurfaceScalarFieldList>
205 void limitSum
206 (
207  const SurfaceScalarFieldList& alphas,
208  SurfaceScalarFieldList& phiPsiCorrs,
209  const labelHashSet& fixed
210 );
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 } // End namespace MULES
215 } // End namespace Foam
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #ifdef NoRepository
220  #include "MULESTemplates.C"
221 #endif
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 #endif
226 
227 // ************************************************************************* //
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
zeroField Su
Definition: alphaSuSp.H:1
surfaceScalarField & phi
IOstream & fixed(IOstream &io)
Definition: IOstream.H:576
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:208
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:30
void limiter(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const volScalarField & psi
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Namespace for OpenFOAM.
zeroField Sp
Definition: alphaSuSp.H:2