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
backwardDdtScheme
backwardDdtScheme.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::backwardDdtScheme
26
27
Description
28
Second-order backward-differencing ddt using the current and
29
two previous time-step values.
30
31
SourceFiles
32
backwardDdtScheme.C
33
34
\*---------------------------------------------------------------------------*/
35
36
#ifndef backwardDdtScheme_H
37
#define backwardDdtScheme_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 backwardDdtScheme Declaration
53
\*---------------------------------------------------------------------------*/
54
55
template
<
class
Type>
56
class
backwardDdtScheme
57
:
58
public
fv::ddtScheme
<Type>
59
{
60
// Private Member Functions
61
62
//- Return the current time-step
63
scalar deltaT_()
const
;
64
65
//- Return the previous time-step
66
scalar deltaT0_()
const
;
67
68
//- Return the previous time-step or GREAT if the old timestep field
69
// wasn't available in which case Euler ddt is used
70
template
<
class
GeoField>
71
scalar deltaT0_(
const
GeoField&)
const
;
72
73
//- Disallow default bitwise copy construct
74
backwardDdtScheme
(
const
backwardDdtScheme
&);
75
76
//- Disallow default bitwise assignment
77
void
operator=(
const
backwardDdtScheme
&);
78
79
80
public
:
81
82
//- Runtime type information
83
TypeName
(
"backward"
);
84
85
86
// Constructors
87
88
//- Construct from mesh
89
backwardDdtScheme
(
const
fvMesh
&
mesh
)
90
:
91
ddtScheme
<Type>(mesh)
92
{}
93
94
//- Construct from mesh and Istream
95
backwardDdtScheme
(
const
fvMesh
&
mesh
,
Istream
& is)
96
:
97
ddtScheme
<Type>(mesh, is)
98
{}
99
100
101
// Member Functions
102
103
//- Return mesh reference
104
const
fvMesh
&
mesh
()
const
105
{
106
return
fv::ddtScheme<Type>::mesh
();
107
}
108
109
tmp<GeometricField<Type, fvPatchField, volMesh>
>
fvcDdt
110
(
111
const
dimensioned<Type>
&
112
);
113
114
tmp<GeometricField<Type, fvPatchField, volMesh>
>
fvcDdt
115
(
116
const
GeometricField<Type, fvPatchField, volMesh>
&
117
);
118
119
tmp<GeometricField<Type, fvPatchField, volMesh>
>
fvcDdt
120
(
121
const
dimensionedScalar
&,
122
const
GeometricField<Type, fvPatchField, volMesh>
&
123
);
124
125
tmp<GeometricField<Type, fvPatchField, volMesh>
>
fvcDdt
126
(
127
const
volScalarField
&,
128
const
GeometricField<Type, fvPatchField, volMesh>
&
129
);
130
131
tmp<fvMatrix<Type>
>
fvmDdt
132
(
133
GeometricField<Type, fvPatchField, volMesh>
&
134
);
135
136
tmp<fvMatrix<Type>
>
fvmDdt
137
(
138
const
dimensionedScalar
&,
139
GeometricField<Type, fvPatchField, volMesh>
&
140
);
141
142
tmp<fvMatrix<Type>
>
fvmDdt
143
(
144
const
volScalarField
&,
145
GeometricField<Type, fvPatchField, volMesh>
&
146
);
147
148
typedef
typename
ddtScheme<Type>::fluxFieldType
fluxFieldType
;
149
150
tmp<fluxFieldType>
fvcDdtPhiCorr
151
(
152
const
volScalarField
& rA,
153
const
GeometricField<Type, fvPatchField, volMesh>
&
U
,
154
const
fluxFieldType
&
phi
155
);
156
157
tmp<fluxFieldType>
fvcDdtPhiCorr
158
(
159
const
volScalarField
& rA,
160
const
volScalarField
&
rho
,
161
const
GeometricField<Type, fvPatchField, volMesh>
&
U
,
162
const
fluxFieldType
&
phi
163
);
164
165
166
tmp<surfaceScalarField>
meshPhi
167
(
168
const
GeometricField<Type, fvPatchField, volMesh>
&
169
);
170
};
171
172
173
template
<>
174
tmp<surfaceScalarField>
backwardDdtScheme<scalar>::fvcDdtPhiCorr
175
(
176
const
volScalarField
& rA,
177
const
volScalarField
&
U
,
178
const
surfaceScalarField
&
phi
179
);
180
181
182
template
<>
183
tmp<surfaceScalarField>
backwardDdtScheme<scalar>::fvcDdtPhiCorr
184
(
185
const
volScalarField
& rA,
186
const
volScalarField
&
rho
,
187
const
volScalarField
&
U
,
188
const
surfaceScalarField
&
phi
189
);
190
191
192
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193
194
}
// End namespace fv
195
196
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197
198
}
// End namespace Foam
199
200
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201
202
#ifdef NoRepository
203
# include <
finiteVolume/backwardDdtScheme.C
>
204
#endif
205
206
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207
208
#endif
209
210
// ************************ vim: set sw=4 sts=4 et: ************************ //