Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
src
finiteVolume
finiteVolume
ddtSchemes
EulerDdtScheme
EulerDdtScheme.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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
Class
25
Foam::fv::EulerDdtScheme
26
27
Description
28
Basic first-order Euler implicit/explicit ddt using only the current and
29
previous time-step values.
30
31
SourceFiles
32
EulerDdtScheme.C
33
34
\*---------------------------------------------------------------------------*/
35
36
#ifndef EulerDdtScheme_H
37
#define EulerDdtScheme_H
38
39
#include <
finiteVolume/ddtScheme.H
>
40
41
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43
namespace
Foam
44
{
45
46
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48
namespace
fv
49
{
50
51
/*---------------------------------------------------------------------------*\
52
Class EulerDdtScheme Declaration
53
\*---------------------------------------------------------------------------*/
54
55
template
<
class
Type>
56
class
EulerDdtScheme
57
:
58
public
ddtScheme
<Type>
59
{
60
// Private Member Functions
61
62
//- Disallow default bitwise copy construct
63
EulerDdtScheme
(
const
EulerDdtScheme
&);
64
65
//- Disallow default bitwise assignment
66
void
operator=(
const
EulerDdtScheme
&);
67
68
69
public
:
70
71
//- Runtime type information
72
TypeName
(
"Euler"
);
73
74
75
// Constructors
76
77
//- Construct from mesh
78
EulerDdtScheme
(
const
fvMesh
&
mesh
)
79
:
80
ddtScheme
<Type>(mesh)
81
{}
82
83
//- Construct from mesh and Istream
84
EulerDdtScheme
(
const
fvMesh
&
mesh
,
Istream
& is)
85
:
86
ddtScheme
<Type>(mesh, is)
87
{}
88
89
90
// Member Functions
91
92
//- Return mesh reference
93
const
fvMesh
&
mesh
()
const
94
{
95
return
fv::ddtScheme<Type>::mesh
();
96
}
97
98
tmp<GeometricField<Type, fvPatchField, volMesh>
>
fvcDdt
99
(
100
const
dimensioned<Type>
&
101
);
102
103
tmp<GeometricField<Type, fvPatchField, volMesh>
>
fvcDdt
104
(
105
const
GeometricField<Type, fvPatchField, volMesh>
&
106
);
107
108
tmp<GeometricField<Type, fvPatchField, volMesh>
>
fvcDdt
109
(
110
const
dimensionedScalar
&,
111
const
GeometricField<Type, fvPatchField, volMesh>
&
112
);
113
114
tmp<GeometricField<Type, fvPatchField, volMesh>
>
fvcDdt
115
(
116
const
volScalarField
&,
117
const
GeometricField<Type, fvPatchField, volMesh>
&
118
);
119
120
tmp<fvMatrix<Type>
>
fvmDdt
121
(
122
GeometricField<Type, fvPatchField, volMesh>
&
123
);
124
125
tmp<fvMatrix<Type>
>
fvmDdt
126
(
127
const
dimensionedScalar
&,
128
GeometricField<Type, fvPatchField, volMesh>
&
129
);
130
131
tmp<fvMatrix<Type>
>
fvmDdt
132
(
133
const
volScalarField
&,
134
GeometricField<Type, fvPatchField, volMesh>
&
135
);
136
137
typedef
typename
ddtScheme<Type>::fluxFieldType
fluxFieldType
;
138
139
tmp<fluxFieldType>
fvcDdtPhiCorr
140
(
141
const
volScalarField
& rA,
142
const
GeometricField<Type, fvPatchField, volMesh>
&
U
,
143
const
fluxFieldType
&
phi
144
);
145
146
tmp<fluxFieldType>
fvcDdtPhiCorr
147
(
148
const
volScalarField
& rA,
149
const
volScalarField
&
rho
,
150
const
GeometricField<Type, fvPatchField, volMesh>
&
U
,
151
const
fluxFieldType
&
phi
152
);
153
154
tmp<surfaceScalarField>
meshPhi
155
(
156
const
GeometricField<Type, fvPatchField, volMesh>
&
157
);
158
};
159
160
161
template
<>
162
tmp<surfaceScalarField>
EulerDdtScheme<scalar>::fvcDdtPhiCorr
163
(
164
const
volScalarField
& rA,
165
const
volScalarField
&
U
,
166
const
surfaceScalarField
&
phi
167
);
168
169
170
template
<>
171
tmp<surfaceScalarField>
EulerDdtScheme<scalar>::fvcDdtPhiCorr
172
(
173
const
volScalarField
& rA,
174
const
volScalarField
&
rho
,
175
const
volScalarField
&
U
,
176
const
surfaceScalarField
&
phi
177
);
178
179
180
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181
182
}
// End namespace fv
183
184
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185
186
}
// End namespace Foam
187
188
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189
190
#ifdef NoRepository
191
# include <
finiteVolume/EulerDdtScheme.C
>
192
#endif
193
194
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195
196
#endif
197
198
// ************************ vim: set sw=4 sts=4 et: ************************ //