CoEulerDdtScheme.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::CoEulerDdtScheme
26 
27 Description
28  Courant number limited first-order Euler implicit/explicit ddt.
29 
30  The time-step is adjusted locally so that the local Courant number
31  does not exceed the specified value.
32 
33  This scheme should only be used for steady-state computations
34  using transient codes where local time-stepping is preferable to
35  under-relaxation for transport consistency reasons.
36 
37 SourceFiles
38  CoEulerDdtScheme.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef CoEulerDdtScheme_H
43 #define CoEulerDdtScheme_H
44 
45 #include "ddtScheme.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace fv
55 {
56 
57 /*---------------------------------------------------------------------------*\
58  Class CoEulerDdtScheme Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 template<class Type>
62 class CoEulerDdtScheme
63 :
64  public fv::ddtScheme<Type>
65 {
66  // Private Data
67 
68  //- Name of the flux field used to calculate the local time-step
69  word phiName_;
70 
71  //- Name of the density field used to obtain the volumetric flux
72  // from the mass flux if required
73  word rhoName_;
74 
75  //- Maximum local Courant number
76  scalar maxCo_;
77 
78 
79  // Private Member Functions
80 
81  //- Return the reciprocal of the Courant-number limited time-step
82  tmp<volScalarField> CorDeltaT() const;
83 
84  //- Return the reciprocal of the face-Courant-number limited time-step
85  tmp<surfaceScalarField> CofrDeltaT() const;
86 
87 
88 public:
89 
90  //- Runtime type information
91  TypeName("CoEuler");
92 
93 
94  // Constructors
95 
96  //- Construct from mesh and Istream
97  CoEulerDdtScheme(const fvMesh& mesh, Istream& is)
98  :
99  ddtScheme<Type>(mesh, is),
100  phiName_(is),
101  rhoName_(is),
102  maxCo_(readScalar(is))
103  {}
104 
105  //- Disallow default bitwise copy construction
106  CoEulerDdtScheme(const CoEulerDdtScheme&) = delete;
107 
108 
109  // Member Functions
110 
111  //- Return mesh reference
112  const fvMesh& mesh() const
113  {
114  return fv::ddtScheme<Type>::mesh();
115  }
116 
118  (
119  const dimensioned<Type>&
120  );
121 
123  (
125  );
126 
128  (
129  const dimensionedScalar&,
131  );
132 
134  (
135  const volScalarField&,
137  );
138 
140  (
141  const volScalarField& alpha,
142  const volScalarField& rho,
144  );
145 
146  virtual tmp<fvMatrix<Type>> fvmDdt
147  (
149  );
150 
151  virtual tmp<fvMatrix<Type>> fvmDdt
152  (
153  const dimensionedScalar&,
155  );
156 
157  virtual tmp<fvMatrix<Type>> fvmDdt
158  (
159  const volScalarField&,
161  );
162 
163  virtual tmp<fvMatrix<Type>> fvmDdt
164  (
165  const volScalarField& alpha,
166  const volScalarField& rho,
168  );
171 
173  (
176  );
177 
179  (
181  const fluxFieldType& phi
182  );
183 
185  (
186  const volScalarField& rho,
189  );
190 
192  (
193  const volScalarField& rho,
195  const fluxFieldType& phi
196  );
197 
199  (
201  );
202 
203  virtual tmp<scalarField> meshPhi
204  (
206  const label patchi
207  );
208 
209 
210  // Member Operators
211 
212  //- Disallow default bitwise assignment
213  void operator=(const CoEulerDdtScheme&) = delete;
214 };
215 
216 
217 template<>
219 (
222 );
223 
224 template<>
226 (
227  const volScalarField& U,
228  const surfaceScalarField& phi
229 );
230 
231 template<>
233 (
234  const volScalarField& rho,
235  const volScalarField& U,
236  const surfaceScalarField& Uf
237 );
238 
239 template<>
241 (
242  const volScalarField& rho,
243  const volScalarField& U,
244  const surfaceScalarField& phi
245 );
246 
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 } // End namespace fv
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 } // End namespace Foam
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #ifdef NoRepository
259  #include "CoEulerDdtScheme.C"
260 #endif
261 
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 
264 #endif
265 
266 // ************************************************************************* //
const fvMesh & mesh() const
Return mesh reference.
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDdt(const dimensioned< Type > &)
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))
Generic GeometricField class.
Generic dimensioned Type class.
TypeName("CoEuler")
Runtime type information.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:64
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:129
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
CoEulerDdtScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
const volScalarField & psi
Courant number limited first-order Euler implicit/explicit ddt.
phi
Definition: correctPhi.H:3
autoPtr< surfaceVectorField > Uf
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
label patchi
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
void operator=(const CoEulerDdtScheme &)=delete
Disallow default bitwise assignment.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
ddtScheme< Type >::fluxFieldType fluxFieldType
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.