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