MULES.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Global
27  MULES
28 
29 Description
30  MULES: Multidimensional universal limiter for explicit solution.
31 
32  Solve a convective-only transport equation using an explicit universal
33  multi-dimensional limiter.
34 
35  Parameters are the variable to solve, the normal convective flux and the
36  actual explicit flux of the variable which is also used to return limited
37  flux used in the bounded-solution.
38 
39 SourceFiles
40  MULES.C
41  MULESTemplates.C
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #ifndef MULES_H
46 #define MULES_H
47 
48 #include "volFieldsFwd.H"
49 #include "surfaceFieldsFwd.H"
50 #include "primitiveFieldsFwd.H"
51 #include "geometricOneField.H"
52 #include "zero.H"
53 #include "zeroField.H"
54 #include "UPtrList.H"
55 #include "HashSet.H"
56 #include "UniformField.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
62 namespace MULES
63 {
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 template<class RdeltaTType, class RhoType, class SpType, class SuType>
68 void explicitSolve
69 (
70  const RdeltaTType& rDeltaT,
71  const RhoType& rho,
73  const surfaceScalarField& phiPsi,
74  const SpType& Sp,
75  const SuType& Su
76 );
77 
78 template<class RhoType>
79 void explicitSolve
80 (
81  const RhoType& rho,
83  const surfaceScalarField& phiPsi
84 );
85 
86 template<class RhoType, class SpType, class SuType>
87 void explicitSolve
88 (
89  const RhoType& rho,
91  const surfaceScalarField& phiPsi,
92  const SpType& Sp,
93  const SuType& Su
94 );
95 
96 template<class RhoType, class PsiMaxType, class PsiMinType>
97 void explicitSolve
98 (
99  const RhoType& rho,
101  const surfaceScalarField& phiBD,
102  surfaceScalarField& phiPsi,
103  const PsiMaxType& psiMax,
104  const PsiMinType& psiMin
105 );
106 
107 
108 template
109 <
110  class RhoType,
111  class SpType,
112  class SuType,
113  class PsiMaxType,
114  class PsiMinType
115 >
116 void explicitSolve
117 (
118  const RhoType& rho,
120  const surfaceScalarField& phiBD,
121  surfaceScalarField& phiPsi,
122  const SpType& Sp,
123  const SuType& Su,
124  const PsiMaxType& psiMax,
125  const PsiMinType& psiMin
126 );
127 
128 
129 template
130 <
131  class RdeltaTType,
132  class RhoType,
133  class SpType,
134  class SuType,
135  class PsiMaxType,
136  class PsiMinType
137 >
138 void limiter
139 (
140  scalarField& allLambda,
141  const RdeltaTType& rDeltaT,
142  const RhoType& rho,
143  const volScalarField& psi,
144  const surfaceScalarField& phiBD,
145  const surfaceScalarField& phiCorr,
146  const SpType& Sp,
147  const SuType& Su,
148  const PsiMaxType& psiMax,
149  const PsiMinType& psiMin
150 );
151 
152 
153 template
154 <
155  class RdeltaTType,
156  class RhoType,
157  class SpType,
158  class SuType,
159  class PsiMaxType,
160  class PsiMinType
161 >
162 void limit
163 (
164  const RdeltaTType& rDeltaT,
165  const RhoType& rho,
166  const volScalarField& psi,
167  const surfaceScalarField& phi,
168  surfaceScalarField& phiPsi,
169  const SpType& Sp,
170  const SuType& Su,
171  const PsiMaxType& psiMax,
172  const PsiMinType& psiMin,
173  const bool returnCorr
174 );
175 
176 
177 template
178 <
179  class RhoType,
180  class SpType,
181  class SuType,
182  class PsiMaxType,
183  class PsiMinType
184 >
185 void limit
186 (
187  const RhoType& rho,
188  const volScalarField& psi,
189  const surfaceScalarField& phi,
190  surfaceScalarField& phiPsi,
191  const SpType& Sp,
192  const SuType& Su,
193  const PsiMaxType& psiMax,
194  const PsiMinType& psiMin,
195  const bool returnCorr
196 );
197 
198 
199 void limitSum(UPtrList<scalarField>& phiPsiCorrs);
200 
201 template<class SurfaceScalarFieldList>
202 void limitSum(SurfaceScalarFieldList& phiPsiCorrs);
203 
204 void limitSum
205 (
206  const UPtrList<const scalarField>& alphas,
207  UPtrList<scalarField>& phiPsiCorrs,
208  const labelHashSet& fixed
209 );
210 
211 template<class SurfaceScalarFieldList>
212 void limitSum
213 (
214  const SurfaceScalarFieldList& alphas,
215  SurfaceScalarFieldList& phiPsiCorrs,
216  const labelHashSet& fixed
217 );
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace MULES
222 } // End namespace Foam
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 #ifdef NoRepository
227  #include "MULESTemplates.C"
228 #endif
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 #endif
233 
234 // ************************************************************************* //
volFieldsFwd.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
UPtrList.H
primitiveFieldsFwd.H
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
Foam::MULES::explicitSolve
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
Definition: MULESTemplates.C:41
Foam::MULES::limit
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)
Definition: MULESTemplates.C:581
Sp
zeroField Sp
Definition: alphaSuSp.H:2
rho
rho
Definition: readInitialConditions.H:88
zeroField.H
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::MULES::limitSum
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Definition: MULES.C:34
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::limiter
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)
Definition: MULESTemplates.C:195
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
zero.H
HashSet.H
Foam::fixed
IOstream & fixed(IOstream &io)
Definition: IOstream.H:449
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
MULESTemplates.C
geometricOneField.H
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
surfaceFieldsFwd.H
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:409
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
UniformField.H