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-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::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 
88  virtual tmp<VolField<Type>> fvcDdt
89  (
90  const dimensioned<Type>&
91  );
92 
93  virtual tmp<VolField<Type>> fvcDdt
94  (
95  const VolField<Type>&
96  );
97 
98  virtual tmp<VolField<Type>> fvcDdt
99  (
100  const dimensionedScalar&,
101  const VolField<Type>&
102  );
103 
104  virtual tmp<VolField<Type>> fvcDdt
105  (
106  const volScalarField&,
107  const VolField<Type>&
108  );
109 
110  virtual tmp<VolField<Type>> fvcDdt
111  (
112  const volScalarField& alpha,
113  const volScalarField& rho,
114  const VolField<Type>& psi
115  );
116 
118  (
119  const SurfaceField<Type>&
120  );
121 
122  virtual tmp<fvMatrix<Type>> fvmDdt
123  (
124  const VolField<Type>&
125  );
126 
127  virtual tmp<fvMatrix<Type>> fvmDdt
128  (
129  const dimensionedScalar&,
130  const VolField<Type>&
131  );
132 
133  virtual tmp<fvMatrix<Type>> fvmDdt
134  (
135  const volScalarField&,
136  const VolField<Type>&
137  );
138 
139  virtual tmp<fvMatrix<Type>> fvmDdt
140  (
141  const volScalarField& alpha,
142  const volScalarField& rho,
143  const VolField<Type>& psi
144  );
145 
147 
149  (
150  const VolField<Type>& U,
151  const SurfaceField<Type>& Uf
152  );
153 
155  (
156  const VolField<Type>& U,
157  const fluxFieldType& phi
158  );
159 
161  (
162  const volScalarField& rho,
163  const VolField<Type>& U,
164  const SurfaceField<Type>& rhoUf
165  );
166 
168  (
169  const volScalarField& rho,
170  const VolField<Type>& U,
171  const fluxFieldType& phi
172  );
173 
175  (
176  const volScalarField& alpha,
177  const volScalarField& rho,
178  const VolField<Type>& U,
179  const SurfaceField<Type>& Uf
180  );
181 
183  (
184  const volScalarField& alpha,
185  const volScalarField& rho,
186  const VolField<Type>& U,
187  const fluxFieldType& phi
188  );
189 
191  (
192  const VolField<Type>&
193  );
194 
195  virtual tmp<scalarField> meshPhi
196  (
197  const VolField<Type>&,
198  const label patchi
199  );
200 
201 
202  // Member Operators
203 
204  //- Disallow default bitwise assignment
205  void operator=(const EulerDdtScheme&) = delete;
206 };
207 
208 
209 template<>
211 (
212  const VolField<scalar>& U,
213  const SurfaceField<scalar>& Uf
214 );
215 
216 template<>
218 (
219  const volScalarField& U,
220  const surfaceScalarField& phi
221 );
222 
223 template<>
225 (
226  const volScalarField& rho,
227  const volScalarField& U,
228  const surfaceScalarField& Uf
229 );
230 
231 template<>
233 (
234  const volScalarField& rho,
235  const volScalarField& U,
236  const surfaceScalarField& phi
237 );
238 
239 template<>
241 (
242  const volScalarField& alpha,
243  const volScalarField& rho,
244  const volScalarField& U,
245  const surfaceScalarField& Uf
246 );
247 
248 template<>
250 (
251  const volScalarField& alpha,
252  const volScalarField& rho,
253  const volScalarField& U,
254  const surfaceScalarField& phi
255 );
256 
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 } // End namespace fv
261 
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 
264 } // End namespace Foam
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 #ifdef NoRepository
269  #include "EulerDdtScheme.C"
270 #endif
271 
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 
274 #endif
275 
276 // ************************************************************************* //
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:101
Basic first-order Euler implicit/explicit ddt using only the current and 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 EulerDdtScheme &)=delete
Disallow default bitwise assignment.
virtual tmp< surfaceScalarField > meshPhi(const VolField< Type > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
TypeName("Euler")
Runtime type information.
virtual tmp< fluxFieldType > fvcDdtUfCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)
EulerDdtScheme(const fvMesh &mesh)
Construct from mesh.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:67
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:129
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)