fvcDdt.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 InNamespace
25  Foam::fvc
26 
27 Description
28  Calculate the first temporal derivative.
29 
30 SourceFiles
31  fvcDdt.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 
36 #ifndef fvcDdt_H
37 #define fvcDdt_H
38 
39 #include "volFieldsFwd.H"
40 #include "surfaceFieldsFwd.H"
41 #include "dimensionedTypes.H"
42 #include "one.H"
43 #include "geometricZeroField.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Namespace fvc functions Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 namespace fvc
55 {
56  template<class Type>
57  tmp<VolField<Type>> ddt
58  (
59  const dimensioned<Type>,
60  const fvMesh&
61  );
62 
63  template<class Type>
64  tmp<VolField<Type>> ddt
65  (
66  const VolField<Type>&
67  );
68 
69  template<class Type>
70  tmp<VolField<Type>> ddt
71  (
72  const dimensionedScalar&,
73  const VolField<Type>&
74  );
75 
76  template<class Type>
77  tmp<VolField<Type>> ddt
78  (
79  const volScalarField&,
80  const VolField<Type>&
81  );
82 
83  template<class Type>
84  tmp<VolField<Type>> ddt
85  (
86  const one&,
87  const VolField<Type>&
88  );
89 
90  template<class Type>
91  tmp<VolField<Type>> ddt
92  (
93  const volScalarField&,
94  const volScalarField&,
95  const VolField<Type>&
96  );
97 
98  template<class Type>
99  tmp<VolField<Type>> ddt
100  (
101  const one&,
102  const one&,
103  const VolField<Type>&
104  );
105 
106  template<class Type>
107  tmp<VolField<Type>> ddt
108  (
109  const one&,
110  const volScalarField&,
111  const VolField<Type>&
112  );
113 
114  template<class Type>
115  tmp<VolField<Type>> ddt
116  (
117  const volScalarField&,
118  const one&,
119  const VolField<Type>&
120  );
121 
122  template<class Type>
123  tmp<SurfaceField<Type>> ddt
124  (
125  const SurfaceField<Type>&
126  );
127 
128  template<class Type>
129  tmp<SurfaceField<typename Foam::flux<Type>::type>> ddtCorr
130  (
131  const VolField<Type>& U,
132  const SurfaceField<Type>& Uf
133  );
134 
135  template<class Type>
136  tmp<SurfaceField<typename Foam::flux<Type>::type>> ddtCorr
137  (
138  const VolField<Type>& U,
139  const SurfaceField<typename Foam::flux<Type>::type>& phi
140  );
141 
142  template<class Type>
143  tmp<SurfaceField<typename Foam::flux<Type>::type>> ddtCorr
144  (
145  const VolField<Type>& U,
146  const SurfaceField<typename Foam::flux<Type>::type>& phi,
147  const autoPtr<SurfaceField<Type>>& Uf
148  );
149 
150  template<class Type>
151  tmp<SurfaceField<typename Foam::flux<Type>::type>> ddtCorr
152  (
153  const volScalarField& rho,
154  const VolField<Type>& U,
155  const SurfaceField<Type>& rhoUf
156  );
157 
158  template<class Type>
159  tmp<SurfaceField<typename Foam::flux<Type>::type>> ddtCorr
160  (
161  const volScalarField& rho,
162  const VolField<Type>& U,
163  const SurfaceField<typename Foam::flux<Type>::type>& phi
164  );
165 
166  template<class Type>
167  tmp<SurfaceField<typename Foam::flux<Type>::type>> ddtCorr
168  (
169  const volScalarField& rho,
170  const VolField<Type>& U,
171  const SurfaceField<typename Foam::flux<Type>::type>& phi,
172  const autoPtr<SurfaceField<Type>>& rhoUf
173  );
174 
175  template<class Type>
176  tmp<SurfaceField<typename Foam::flux<Type>::type>> ddtCorr
177  (
178  const volScalarField& alpha,
179  const volScalarField& rho,
180  const VolField<Type>& U,
181  const SurfaceField<Type>& Uf
182  );
183 
184  template<class Type>
185  tmp<SurfaceField<typename Foam::flux<Type>::type>> ddtCorr
186  (
187  const volScalarField& alpha,
188  const volScalarField& rho,
189  const VolField<Type>& U,
190  const SurfaceField<typename Foam::flux<Type>::type>& phi
191  );
192 
193  template<class Type>
194  tmp<SurfaceField<typename Foam::flux<Type>::type>> ddtCorr
195  (
196  const volScalarField& alpha,
197  const volScalarField& rho,
198  const VolField<Type>& U,
199  const SurfaceField<typename Foam::flux<Type>::type>& phi,
200  const autoPtr<SurfaceField<Type>>& Uf
201  );
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace Foam
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 #ifdef NoRepository
212  #include "fvcDdt.C"
213 #endif
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 #endif
218 
219 // ************************************************************************* //
typeOfRank< typename pTraits< vector >::cmptType, direction(pTraits< vector >::rank)+direction(pTraits< Type >::rank) - 2 >::type type
Definition: products.H:115
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< SurfaceField< typename Foam::flux< Type >::type > > ddtCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
Definition: fvcDdt.C:196
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
GeometricField< Type, fvsPatchField, surfaceMesh > SurfaceField