All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-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::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 
91  virtual tmp<VolField<Type>> fvcDdt
92  (
93  const dimensioned<Type>&
94  );
95 
96  virtual tmp<VolField<Type>> fvcDdt
97  (
98  const VolField<Type>&
99  );
100 
101  virtual tmp<VolField<Type>> fvcDdt
102  (
103  const dimensionedScalar&,
104  const VolField<Type>&
105  );
106 
107  virtual tmp<VolField<Type>> fvcDdt
108  (
109  const volScalarField&,
110  const VolField<Type>&
111  );
112 
113  virtual tmp<VolField<Type>> fvcDdt
114  (
115  const volScalarField& alpha,
116  const volScalarField& rho,
117  const VolField<Type>& psi
118  );
119 
120  virtual tmp<fvMatrix<Type>> fvmDdt
121  (
122  const VolField<Type>&
123  );
124 
125  virtual tmp<fvMatrix<Type>> fvmDdt
126  (
127  const dimensionedScalar&,
128  const VolField<Type>&
129  );
130 
131  virtual tmp<fvMatrix<Type>> fvmDdt
132  (
133  const volScalarField&,
134  const VolField<Type>&
135  );
136 
137  virtual tmp<fvMatrix<Type>> fvmDdt
138  (
139  const volScalarField& alpha,
140  const volScalarField& rho,
141  const VolField<Type>& psi
142  );
143 
145 
147  (
148  const VolField<Type>& U,
149  const SurfaceField<Type>& Uf
150  );
151 
153  (
154  const VolField<Type>& U,
155  const fluxFieldType& phi
156  );
157 
159  (
160  const volScalarField& rho,
161  const VolField<Type>& U,
162  const SurfaceField<Type>& rhoUf
163  );
164 
166  (
167  const volScalarField& rho,
168  const VolField<Type>& U,
169  const fluxFieldType& phi
170  );
171 
173  (
174  const volScalarField& alpha,
175  const volScalarField& rho,
176  const VolField<Type>& U,
177  const SurfaceField<Type>& Uf
178  )
179  {
181  return fluxFieldType::null();
182  }
183 
185  (
186  const volScalarField& alpha,
187  const volScalarField& rho,
188  const VolField<Type>& U,
189  const fluxFieldType& phi
190  )
191  {
193  return fluxFieldType::null();
194  }
195 
197  (
198  const VolField<Type>&
199  );
200 
201  virtual tmp<scalarField> meshPhi
202  (
203  const VolField<Type>&,
204  const label patchi
205  );
206 
207 
208  // Member Operators
209 
210  //- Disallow default bitwise assignment
211  void operator=(const steadyStateDdtScheme&) = delete;
212 };
213 
214 
215 template<>
217 (
218  const VolField<scalar>& U,
219  const SurfaceField<scalar>& Uf
220 );
221 
222 template<>
224 (
225  const volScalarField& U,
226  const surfaceScalarField& phi
227 );
228 
229 template<>
231 (
232  const volScalarField& rho,
233  const volScalarField& U,
234  const surfaceScalarField& Uf
235 );
236 
237 template<>
239 (
240  const volScalarField& rho,
241  const volScalarField& U,
242  const surfaceScalarField& phi
243 );
244 
245 template<>
247 (
248  const volScalarField& alpha,
249  const volScalarField& rho,
250  const volScalarField& U,
251  const surfaceScalarField& Uf
252 );
253 
254 template<>
256 (
257  const volScalarField& alpha,
258  const volScalarField& rho,
259  const volScalarField& U,
260  const surfaceScalarField& phi
261 );
262 
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
266 } // End namespace fv
267 
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
269 
270 } // End namespace Foam
271 
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273 
274 #ifdef NoRepository
275  #include "steadyStateDdtScheme.C"
276 #endif
277 
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 
280 #endif
281 
282 // ************************************************************************* //
Generic GeometricField class.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
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
Abstract base class for ddt schemes.
Definition: ddtScheme.H:68
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:130
SteadyState implicit/explicit ddt which returns 0.
virtual tmp< fluxFieldType > fvcDdtUfCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
virtual tmp< fvMatrix< Type > > fvmDdt(const VolField< Type > &)
const fvMesh & mesh() const
Return mesh reference.
steadyStateDdtScheme(const fvMesh &mesh)
Construct from mesh.
virtual tmp< surfaceScalarField > meshPhi(const VolField< Type > &)
ddtScheme< Type >::fluxFieldType fluxFieldType
TypeName("steadyState")
Runtime type information.
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const VolField< Type > &U, const fluxFieldType &phi)
virtual tmp< VolField< Type > > fvcDdt(const dimensioned< Type > &)
void operator=(const steadyStateDdtScheme &)=delete
Disallow default bitwise assignment.
A class for managing temporary objects.
Definition: tmp.H:55
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
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)