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) 2011-2016 OpenFOAM Foundation
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 "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>
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
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
90  :
91  ddtScheme<Type>(mesh)
92  {}
93 
94  //- Construct from mesh and Istream
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 
110  (
111  const dimensioned<Type>&
112  );
113 
115  (
117  );
118 
120  (
121  const dimensionedScalar&,
123  );
124 
126  (
127  const volScalarField&,
129  );
130 
132  (
133  const volScalarField& alpha,
134  const volScalarField& rho,
136  );
137 
139  (
141  );
142 
144  (
145  const dimensionedScalar&,
147  );
148 
150  (
151  const volScalarField&,
153  );
154 
156  (
157  const volScalarField& alpha,
158  const volScalarField& rho,
160  );
163 
165  (
168  );
169 
171  (
173  const fluxFieldType& phi
174  );
175 
177  (
178  const volScalarField& rho,
181  );
182 
184  (
185  const volScalarField& rho,
187  const fluxFieldType& phi
188  );
189 
191  (
193  );
194 };
195 
196 
197 template<>
199 (
202 );
203 
204 template<>
206 (
207  const volScalarField& U,
208  const surfaceScalarField& phi
209 );
210 
211 template<>
213 (
214  const volScalarField& rho,
215  const volScalarField& U,
216  const surfaceScalarField& Uf
217 );
218 
219 template<>
221 (
222  const volScalarField& rho,
223  const volScalarField& U,
224  const surfaceScalarField& phi
225 );
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace fv
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace Foam
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 #ifdef NoRepository
239  #include "backwardDdtScheme.C"
240 #endif
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 #endif
245 
246 // ************************************************************************* //
surfaceScalarField & phi
const fvMesh & mesh() const
Return mesh reference.
U
Definition: pEqn.H:83
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:135
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
Uf
Definition: pEqn.H:67
Generic GeometricField class.
Generic dimensioned Type class.
tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Abstract base class for ddt schemes.
Definition: ddtScheme.H:64
labelList fv(nPoints)
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
ddtScheme< Type >::fluxFieldType fluxFieldType
Second-order backward-differencing ddt using the current and two previous time-step values...
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
TypeName("backward")
Runtime type information.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const volScalarField & psi
A class for managing temporary objects.
Definition: PtrList.H:54
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.