Loading...
Searching...
No Matches
SLTSDdtScheme.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-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::fv::SLTSDdtScheme
28
29Group
30 grpFvDdtSchemes
31
32Description
33 Stabilised local time-step first-order Euler implicit/explicit ddt.
34 The time-step is adjusted locally so that an advective equations remains
35 diagonally dominant.
36
37 This scheme should only be used for steady-state computations
38 using transient codes where local time-stepping is preferably to
39 under-relaxation for transport consistency reasons.
40
41See also
42 Foam::fv::CoEulerDdtScheme
43
44SourceFiles
45 SLTSDdtScheme.C
46
47\*---------------------------------------------------------------------------*/
48
49#ifndef SLTSDdtScheme_H
50#define SLTSDdtScheme_H
51
52#include "ddtScheme.H"
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56namespace Foam
57{
58
59// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61namespace fv
62{
63
64/*---------------------------------------------------------------------------*\
65 Class SLTSDdtScheme Declaration
66\*---------------------------------------------------------------------------*/
67
68template<class Type>
69class SLTSDdtScheme
70:
71 public fv::ddtScheme<Type>
72{
73 // Private Data
74
75 //- Name of the flux field used to calculate the local time-step
76 word phiName_;
77
78 //- Name of the density field used to obtain the volumetric flux
79 // from the mass flux if required
80 word rhoName_;
81
82 //- Under-relaxation factor
83 scalar alpha_;
84
85
86 // Private Member Functions
87
88 //- No copy construct
89 SLTSDdtScheme(const SLTSDdtScheme&) = delete;
90
91 //- No copy assignment
92 void operator=(const SLTSDdtScheme&) = delete;
93
94 //- Calculate a relaxed diagonal from the given flux field
95 void relaxedDiag(scalarField& rD, const surfaceScalarField& phi) const;
96
97 //- Return the reciprocal of the stabilised local time-step
98 tmp<volScalarField> SLrDeltaT() const;
99
100
101public:
102
103 //- Runtime type information
104 TypeName("SLTS");
105
106
107 // Constructors
108
109 //- Construct from mesh and Istream
110 SLTSDdtScheme(const fvMesh& mesh, Istream& is)
111 :
113 phiName_(is),
114 rhoName_(is),
115 alpha_(readScalar(is))
116 {}
117
118
119 // Member Functions
120
121 //- Return mesh reference
122 const fvMesh& mesh() const
123 {
131
133 (
135 );
136
142
144 (
145 const volScalarField&,
147 );
148
150 (
151 const volScalarField& alpha,
152 const volScalarField& rho,
154 );
155
157 (
159 );
160
162 (
163 const dimensionedScalar&,
165 );
166
168 (
169 const volScalarField&,
171 );
172
174 (
175 const volScalarField& alpha,
176 const volScalarField& rho,
178 );
179
181
183 (
186 );
187
189 (
191 const fluxFieldType& phi
192 );
193
195 (
199 );
200
202 (
203 const volScalarField& rho,
205 const fluxFieldType& phi
206 );
207
209 (
211 );
212};
213
214
215template<>
217(
220);
221
222template<>
224(
225 const volScalarField& U,
227);
228
229template<>
231(
233 const volScalarField& U,
235);
236
237template<>
240 const volScalarField& rho,
241 const volScalarField& U,
243);
244
245
246// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247
248} // End namespace fv
249
250// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251
252} // End namespace Foam
253
254// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255
256#ifdef NoRepository
257 #include "SLTSDdtScheme.C"
258#endif
259
260// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261
262#endif
263
264// ************************************************************************* //
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
SLTSDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
TypeName("SLTS")
Runtime type information.
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
const fvMesh & mesh() const
Return mesh reference.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Abstract base class for ddt schemes.
Definition ddtScheme.H:85
ddtScheme(const ddtScheme &)=delete
No copy construct.
const fvMesh & mesh() const
Return mesh reference.
Definition ddtScheme.H:179
GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > fluxFieldType
Definition ddtScheme.H:247
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
const volScalarField & psi
autoPtr< surfaceVectorField > Uf
Namespace for finite-volume.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68