EulerDdtScheme.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::EulerDdtScheme
26 
27 Description
28  Basic first-order Euler implicit/explicit ddt using only the current and
29  previous time-step values.
30 
31 SourceFiles
32  EulerDdtScheme.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef EulerDdtScheme_H
37 #define EulerDdtScheme_H
38 
39 #include "ddtScheme.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace fv
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class EulerDdtScheme Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 template<class Type>
56 class EulerDdtScheme
57 :
58  public ddtScheme<Type>
59 {
60 public:
61 
62  //- Runtime type information
63  TypeName("Euler");
64 
65 
66  // Constructors
67 
68  //- Construct from mesh
69  EulerDdtScheme(const fvMesh& mesh)
70  :
71  ddtScheme<Type>(mesh)
72  {}
73 
74  //- Construct from mesh and Istream
75  EulerDdtScheme(const fvMesh& mesh, Istream& is)
76  :
77  ddtScheme<Type>(mesh, is)
78  {}
79 
80  //- Disallow default bitwise copy construction
81  EulerDdtScheme(const EulerDdtScheme&) = delete;
82 
83 
84  // Member Functions
85 
87 
89  (
90  const dimensioned<Type>&
91  );
92 
94  (
96  );
97 
99  (
100  const dimensionedScalar&,
102  );
103 
105  (
106  const volScalarField&,
108  );
109 
111  (
112  const volScalarField& alpha,
113  const volScalarField& rho,
115  );
116 
118  (
120  );
121 
122  virtual tmp<fvMatrix<Type>> fvmDdt
123  (
125  );
126 
127  virtual tmp<fvMatrix<Type>> fvmDdt
128  (
129  const dimensionedScalar&,
131  );
132 
133  virtual tmp<fvMatrix<Type>> fvmDdt
134  (
135  const volScalarField&,
137  );
138 
139  virtual tmp<fvMatrix<Type>> fvmDdt
140  (
141  const volScalarField& alpha,
142  const volScalarField& rho,
144  );
147 
149  (
152  );
153 
155  (
157  const fluxFieldType& phi
158  );
159 
161  (
162  const volScalarField& rho,
165  );
166 
168  (
169  const volScalarField& rho,
171  const fluxFieldType& phi
172  );
173 
175  (
177  );
178 
179  virtual tmp<scalarField> meshPhi
180  (
182  const label patchi
183  );
184 
185 
186  // Member Operators
187 
188  //- Disallow default bitwise assignment
189  void operator=(const EulerDdtScheme&) = delete;
190 };
191 
192 
193 template<>
195 (
198 );
199 
200 template<>
202 (
203  const volScalarField& U,
204  const surfaceScalarField& phi
205 );
206 
207 template<>
209 (
210  const volScalarField& rho,
211  const volScalarField& U,
212  const surfaceScalarField& Uf
213 );
214 
215 template<>
217 (
218  const volScalarField& rho,
219  const volScalarField& U,
220  const surfaceScalarField& phi
221 );
222 
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 } // End namespace fv
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace Foam
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 #ifdef NoRepository
235  #include "EulerDdtScheme.C"
236 #endif
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 #endif
241 
242 // ************************************************************************* //
EulerDdtScheme(const fvMesh &mesh)
Construct from mesh.
U
Definition: pEqn.H:72
TypeName("Euler")
Runtime type information.
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 > &)
Generic GeometricField class.
Generic dimensioned Type class.
void operator=(const EulerDdtScheme &)=delete
Disallow default bitwise assignment.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:64
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:129
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
labelList fv(nPoints)
const volScalarField & psi
phi
Definition: correctPhi.H:3
autoPtr< surfaceVectorField > Uf
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
A class for managing temporary objects.
Definition: PtrList.H:53
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values...
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
Namespace for OpenFOAM.