backwardDdtScheme.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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  // Ensure the old-old-time cell volumes are available
94  // for moving meshes
95  if (mesh.moving())
96  {
97  mesh.V00();
98  }
99  }
100 
101  //- Construct from mesh and Istream
102  backwardDdtScheme(const fvMesh& mesh, Istream& is)
103  :
104  ddtScheme<Type>(mesh, is)
105  {
106  // Ensure the old-old-time cell volumes are available
107  // for moving meshes
108  if (mesh.moving())
109  {
110  mesh.V00();
111  }
112  }
113 
114 
115  // Member Functions
116 
117  //- Return mesh reference
118  const fvMesh& mesh() const
119  {
120  return fv::ddtScheme<Type>::mesh();
121  }
122 
124  (
125  const dimensioned<Type>&
126  );
127 
129  (
131  );
132 
134  (
135  const dimensionedScalar&,
137  );
138 
140  (
141  const volScalarField&,
143  );
144 
146  (
147  const volScalarField& alpha,
148  const volScalarField& rho,
150  );
151 
153  (
155  );
156 
158  (
159  const dimensionedScalar&,
161  );
162 
164  (
165  const volScalarField&,
167  );
168 
170  (
171  const volScalarField& alpha,
172  const volScalarField& rho,
174  );
177 
179  (
182  );
183 
185  (
187  const fluxFieldType& phi
188  );
189 
191  (
192  const volScalarField& rho,
195  );
196 
198  (
199  const volScalarField& rho,
201  const fluxFieldType& phi
202  );
203 
205  (
207  );
208 };
209 
210 
211 template<>
213 (
216 );
217 
218 template<>
220 (
221  const volScalarField& U,
222  const surfaceScalarField& phi
223 );
224 
225 template<>
227 (
228  const volScalarField& rho,
229  const volScalarField& U,
230  const surfaceScalarField& Uf
231 );
232 
233 template<>
235 (
236  const volScalarField& rho,
237  const volScalarField& U,
238  const surfaceScalarField& phi
239 );
240 
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 } // End namespace fv
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 } // End namespace Foam
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 #ifdef NoRepository
253  #include "backwardDdtScheme.C"
254 #endif
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #endif
259 
260 // ************************************************************************* //
tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
surfaceScalarField & phi
bool moving() const
Is mesh moving.
Definition: polyMesh.H:502
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 > &)
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
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
const fvMesh & mesh() const
Return mesh reference.
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:135
labelList fv(nPoints)
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
ddtScheme< Type >::fluxFieldType fluxFieldType
autoPtr< surfaceVectorField > Uf
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.
U
Definition: pEqn.H:72
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:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.