All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2023 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 
120  virtual tmp<VolField<Type>> fvcDdt
121  (
122  const dimensioned<Type>&
123  );
124 
125  virtual tmp<VolField<Type>> fvcDdt
126  (
127  const VolField<Type>&
128  );
129 
130  virtual tmp<VolField<Type>> fvcDdt
131  (
132  const dimensionedScalar&,
133  const VolField<Type>&
134  );
135 
136  virtual tmp<VolField<Type>> fvcDdt
137  (
138  const volScalarField&,
139  const VolField<Type>&
140  );
141 
142  virtual tmp<VolField<Type>> fvcDdt
143  (
144  const volScalarField& alpha,
145  const volScalarField& rho,
146  const VolField<Type>& psi
147  );
148 
149  virtual tmp<fvMatrix<Type>> fvmDdt
150  (
151  const VolField<Type>&
152  );
153 
154  virtual tmp<fvMatrix<Type>> fvmDdt
155  (
156  const dimensionedScalar&,
157  const VolField<Type>&
158  );
159 
160  virtual tmp<fvMatrix<Type>> fvmDdt
161  (
162  const volScalarField&,
163  const VolField<Type>&
164  );
165 
166  virtual tmp<fvMatrix<Type>> fvmDdt
167  (
168  const volScalarField& alpha,
169  const volScalarField& rho,
170  const VolField<Type>& psi
171  );
172 
174 
176  (
177  const VolField<Type>& U,
178  const SurfaceField<Type>& Uf
179  );
180 
182  (
183  const VolField<Type>& U,
184  const fluxFieldType& phi
185  );
186 
188  (
189  const volScalarField& rho,
190  const VolField<Type>& U,
191  const SurfaceField<Type>& rhoUf
192  );
193 
195  (
196  const volScalarField& rho,
197  const VolField<Type>& U,
198  const fluxFieldType& phi
199  );
200 
202  (
203  const volScalarField& alpha,
204  const volScalarField& rho,
205  const VolField<Type>& U,
206  const SurfaceField<Type>& Uf
207  );
208 
210  (
211  const volScalarField& alpha,
212  const volScalarField& rho,
213  const VolField<Type>& U,
214  const fluxFieldType& phi
215  );
216 
218  (
219  const VolField<Type>&
220  );
221 
222  virtual tmp<scalarField> meshPhi
223  (
224  const VolField<Type>&,
225  const label patchi
226  );
227 
228 
229  // Member Operators
230 
231  //- Disallow default bitwise assignment
232  void operator=(const backwardDdtScheme&) = delete;
233 };
234 
235 
236 template<>
238 (
239  const VolField<scalar>& U,
240  const SurfaceField<scalar>& Uf
241 );
242 
243 template<>
245 (
246  const volScalarField& U,
247  const surfaceScalarField& phi
248 );
249 
250 template<>
252 (
253  const volScalarField& rho,
254  const volScalarField& U,
255  const surfaceScalarField& rhoUf
256 );
257 
258 template<>
260 (
261  const volScalarField& rho,
262  const volScalarField& U,
263  const surfaceScalarField& phi
264 );
265 
266 template<>
268 (
269  const volScalarField& alpha,
270  const volScalarField& rho,
271  const volScalarField& U,
272  const surfaceScalarField& Uf
273 );
274 
275 template<>
277 (
278  const volScalarField& alpha,
279  const volScalarField& rho,
280  const volScalarField& U,
281  const surfaceScalarField& phi
282 );
283 
284 template<>
286 (
287  const volScalarField& alpha,
288  const volScalarField& rho,
289  const volScalarField& U,
290  const surfaceScalarField& Uf
291 );
292 
293 template<>
295 (
296  const volScalarField& alpha,
297  const volScalarField& rho,
298  const volScalarField& U,
299  const surfaceScalarField& phi
300 );
301 
302 
303 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
304 
305 } // End namespace fv
306 
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 
309 } // End namespace Foam
310 
311 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
312 
313 #ifdef NoRepository
314  #include "backwardDdtScheme.C"
315 #endif
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 #endif
320 
321 // ************************************************************************* //
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
Second-order backward-differencing ddt using the current and two previous time-step values.
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const VolField< Type > &U, const fluxFieldType &phi)
virtual tmp< fvMatrix< Type > > fvmDdt(const VolField< Type > &)
void operator=(const backwardDdtScheme &)=delete
Disallow default bitwise assignment.
const fvMesh & mesh() const
Return mesh reference.
TypeName("backward")
Runtime type information.
virtual tmp< surfaceScalarField > meshPhi(const VolField< Type > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
backwardDdtScheme(const fvMesh &mesh)
Construct from mesh.
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)
virtual tmp< fluxFieldType > fvcDdtUfCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
Abstract base class for ddt schemes.
Definition: ddtScheme.H:68
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:130
bool moving() const
Is mesh moving.
Definition: polyMesh.H:476
A class for managing temporary objects.
Definition: tmp.H:55
label patchi
const volScalarField & psi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
labelList fv(nPoints)