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-2022 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 
74 public:
75 
76  //- Runtime type information
77  TypeName("backward");
78 
79 
80  // Constructors
81 
82  //- Construct from mesh
84  :
85  ddtScheme<Type>(mesh)
86  {
87  // Ensure the old-old-time cell volumes are available
88  // for moving meshes
89  if (mesh.moving())
90  {
91  mesh.V00();
92  }
93  }
94 
95  //- Construct from mesh and Istream
97  :
98  ddtScheme<Type>(mesh, is)
99  {
100  // Ensure the old-old-time cell volumes are available
101  // for moving meshes
102  if (mesh.moving())
103  {
104  mesh.V00();
105  }
106  }
107 
108  //- Disallow default bitwise copy construction
109  backwardDdtScheme(const backwardDdtScheme&) = delete;
110 
111 
112  // Member Functions
113 
114  //- Return mesh reference
115  const fvMesh& mesh() const
116  {
117  return fv::ddtScheme<Type>::mesh();
118  }
119 
121  (
122  const dimensioned<Type>&
123  );
124 
126  (
128  );
129 
131  (
132  const dimensionedScalar&,
134  );
135 
137  (
138  const volScalarField&,
140  );
141 
143  (
144  const volScalarField& alpha,
145  const volScalarField& rho,
147  );
148 
149  virtual tmp<fvMatrix<Type>> fvmDdt
150  (
152  );
153 
154  virtual tmp<fvMatrix<Type>> fvmDdt
155  (
156  const dimensionedScalar&,
158  );
159 
160  virtual tmp<fvMatrix<Type>> fvmDdt
161  (
162  const volScalarField&,
164  );
165 
166  virtual tmp<fvMatrix<Type>> fvmDdt
167  (
168  const volScalarField& alpha,
169  const volScalarField& rho,
171  );
174 
176  (
179  );
180 
182  (
184  const fluxFieldType& phi
185  );
186 
188  (
189  const volScalarField& rho,
192  );
193 
195  (
196  const volScalarField& rho,
198  const fluxFieldType& phi
199  );
200 
202  (
204  );
205 
206  virtual tmp<scalarField> meshPhi
207  (
209  const label patchi
210  );
211 
212 
213  // Member Operators
214 
215  //- Disallow default bitwise assignment
216  void operator=(const backwardDdtScheme&) = delete;
217 };
218 
219 
220 template<>
222 (
225 );
226 
227 template<>
229 (
230  const volScalarField& U,
231  const surfaceScalarField& phi
232 );
233 
234 template<>
236 (
237  const volScalarField& rho,
238  const volScalarField& U,
239  const surfaceScalarField& Uf
240 );
241 
242 template<>
244 (
245  const volScalarField& rho,
246  const volScalarField& U,
247  const surfaceScalarField& phi
248 );
249 
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 } // End namespace fv
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 } // End namespace Foam
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 #ifdef NoRepository
262  #include "backwardDdtScheme.C"
263 #endif
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 #endif
268 
269 // ************************************************************************* //
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
bool moving() const
Is mesh moving.
Definition: polyMesh.H:525
U
Definition: pEqn.H:72
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual 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.
virtual 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:129
labelList fv(nPoints)
const volScalarField & psi
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
phi
Definition: correctPhi.H:3
backwardDdtScheme(const fvMesh &mesh)
Construct from mesh.
ddtScheme< Type >::fluxFieldType fluxFieldType
autoPtr< surfaceVectorField > Uf
Second-order backward-differencing ddt using the current and two previous time-step values...
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
TypeName("backward")
Runtime type information.
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
void operator=(const backwardDdtScheme &)=delete
Disallow default bitwise assignment.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.