fvmDdt.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "volFields.H"
27 #include "surfaceFields.H"
28 #include "fvMatrix.H"
29 #include "ddtScheme.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace fvm
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 template<class Type>
44 tmp<fvMatrix<Type>>
45 ddt
46 (
48 )
49 {
51  (
52  vf.mesh(),
53  vf.mesh().schemes().ddt("ddt(" + vf.name() + ')')
54  ).ref().fvmDdt(vf);
55 }
56 
57 
58 template<class Type>
60 ddt
61 (
62  const one&,
64 )
65 {
66  return ddt(vf);
67 }
68 
69 
70 template<class Type>
72 ddt
73 (
74  const dimensionedScalar& rho,
76 )
77 {
79  (
80  vf.mesh(),
81  vf.mesh().schemes().ddt("ddt(" + rho.name() + ',' + vf.name() + ')')
82  ).ref().fvmDdt(rho, vf);
83 }
84 
85 
86 template<class Type>
88 ddt
89 (
90  const volScalarField& rho,
92 )
93 {
95  (
96  vf.mesh(),
97  vf.mesh().schemes().ddt("ddt(" + rho.name() + ',' + vf.name() + ')')
98  ).ref().fvmDdt(rho, vf);
99 }
100 
101 
102 template<class Type>
104 ddt
105 (
106  const volScalarField& alpha,
107  const volScalarField& rho,
109 )
110 {
112  (
113  vf.mesh(),
114  vf.mesh().schemes().ddt
115  (
116  "ddt("
117  + alpha.name() + ','
118  + rho.name() + ','
119  + vf.name() + ')'
120  )
121  ).ref().fvmDdt(alpha, rho, vf);
122 }
123 
124 
125 template<class Type>
127 ddt
128 (
129  const one&,
130  const one&,
132 )
133 {
134  return ddt(vf);
135 }
136 
137 
138 template<class Type>
140 ddt
141 (
142  const one&,
143  const volScalarField& rho,
145 )
146 {
147  return ddt(rho, vf);
148 }
149 
150 
151 template<class Type>
153 ddt
154 (
155  const volScalarField& alpha,
156  const one&,
158 )
159 {
160  return ddt(alpha, vf);
161 }
162 
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 } // End namespace fvm
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 } // End namespace Foam
171 
172 // ************************************************************************* //
Foam::surfaceFields.
const word & name() const
Return name.
Definition: IOobject.H:315
Generic GeometricField class.
Generic dimensioned Type class.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
const Mesh & mesh() const
Return mesh.
const word & name() const
Return const reference to name.
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:46
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:50