steadyStateDdtScheme.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::steadyStateDdtScheme
26 
27 Description
28  SteadyState implicit/explicit ddt which returns 0.
29 
30 SourceFiles
31  steadyStateDdtScheme.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef steadyStateDdtScheme_H
36 #define steadyStateDdtScheme_H
37 
38 #include "ddtScheme.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace fv
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class steadyStateDdtScheme Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 template<class Type>
56 :
57  public fv::ddtScheme<Type>
58 {
59 public:
60 
61  //- Runtime type information
62  TypeName("steadyState");
63 
64 
65  // Constructors
66 
67  //- Construct from mesh
69  :
70  ddtScheme<Type>(mesh)
71  {}
72 
73  //- Construct from mesh and Istream
75  :
76  ddtScheme<Type>(mesh, is)
77  {}
78 
79  //- Disallow default bitwise copy construction
81 
82 
83  // Member Functions
84 
85  //- Return mesh reference
86  const fvMesh& mesh() const
87  {
89  }
90 
92  (
93  const dimensioned<Type>&
94  );
95 
97  (
99  );
100 
102  (
103  const dimensionedScalar&,
105  );
106 
108  (
109  const volScalarField&,
111  );
112 
114  (
115  const volScalarField& alpha,
116  const volScalarField& rho,
118  );
119 
120  virtual tmp<fvMatrix<Type>> fvmDdt
121  (
123  );
124 
125  virtual tmp<fvMatrix<Type>> fvmDdt
126  (
127  const dimensionedScalar&,
129  );
130 
131  virtual tmp<fvMatrix<Type>> fvmDdt
132  (
133  const volScalarField&,
135  );
136 
137  virtual tmp<fvMatrix<Type>> fvmDdt
138  (
139  const volScalarField& alpha,
140  const volScalarField& rho,
142  );
145 
147  (
150  );
151 
153  (
155  const fluxFieldType& phi
156  );
157 
159  (
160  const volScalarField& rho,
163  );
164 
166  (
167  const volScalarField& rho,
169  const fluxFieldType& phi
170  );
171 
173  (
175  );
176 
177  virtual tmp<scalarField> meshPhi
178  (
180  const label patchi
181  );
182 
183 
184  // Member Operators
185 
186  //- Disallow default bitwise assignment
187  void operator=(const steadyStateDdtScheme&) = delete;
188 };
189 
190 
191 template<>
193 (
196 );
197 
198 template<>
200 (
201  const volScalarField& U,
202  const surfaceScalarField& phi
203 );
204 
205 template<>
207 (
208  const volScalarField& rho,
209  const volScalarField& U,
210  const surfaceScalarField& Uf
211 );
212 
213 template<>
215 (
216  const volScalarField& rho,
217  const volScalarField& U,
218  const surfaceScalarField& phi
219 );
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace fv
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 } // End namespace Foam
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 #ifdef NoRepository
233  #include "steadyStateDdtScheme.C"
234 #endif
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 #endif
239 
240 // ************************************************************************* //
const fvMesh & mesh() const
Return mesh reference.
U
Definition: pEqn.H:72
ddtScheme< Type >::fluxFieldType fluxFieldType
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.
TypeName("steadyState")
Runtime type information.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:64
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:129
labelList fv(nPoints)
const volScalarField & psi
phi
Definition: correctPhi.H:3
SteadyState implicit/explicit ddt which returns 0.
autoPtr< surfaceVectorField > Uf
steadyStateDdtScheme(const fvMesh &mesh)
Construct from mesh.
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
label patchi
virtual tmp< surfaceScalarField > meshPhi(const GeometricField< Type, fvPatchField, volMesh > &)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual tmp< fvMatrix< Type > > fvmDdt(const GeometricField< Type, fvPatchField, volMesh > &)
void operator=(const steadyStateDdtScheme &)=delete
Disallow default bitwise assignment.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.