fvmDdt.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "fvMatrix.H"
31 #include "ddtScheme.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace fvm
41 {
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 template<class Type>
46 tmp<fvMatrix<Type>>
47 ddt
48 (
50 )
51 {
53  (
54  vf.mesh(),
55  vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
56  ).ref().fvmDdt(vf);
57 }
58 
59 
60 template<class Type>
62 ddt
63 (
64  const one&,
66 )
67 {
68  return ddt(vf);
69 }
70 
71 
72 template<class Type>
74 ddt
75 (
76  const dimensionedScalar& rho,
78 )
79 {
81  (
82  vf.mesh(),
83  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
84  ).ref().fvmDdt(rho, vf);
85 }
86 
87 
88 template<class Type>
90 ddt
91 (
92  const volScalarField& rho,
94 )
95 {
97  (
98  vf.mesh(),
99  vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
100  ).ref().fvmDdt(rho, vf);
101 }
102 
103 
104 template<class Type>
106 ddt
107 (
108  const volScalarField& alpha,
109  const volScalarField& rho,
111 )
112 {
114  (
115  vf.mesh(),
116  vf.mesh().ddtScheme
117  (
118  "ddt("
119  + alpha.name() + ','
120  + rho.name() + ','
121  + vf.name() + ')'
122  )
123  ).ref().fvmDdt(alpha, rho, vf);
124 }
125 
126 
127 template<class Type>
129 ddt
130 (
131  const one&,
132  const one&,
134 )
135 {
136  return ddt(vf);
137 }
138 
139 
140 template<class Type>
142 ddt
143 (
144  const one&,
145  const volScalarField& rho,
147 )
148 {
149  return ddt(rho, vf);
150 }
151 
152 
153 template<class Type>
155 ddt
156 (
157  const volScalarField& alpha,
158  const one&,
160 )
161 {
162  return ddt(alpha, vf);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 } // End namespace fvm
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 } // End namespace Foam
173 
174 // ************************************************************************* //
volFields.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:406
Foam::one
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:61
fvMatrix.H
surfaceFields.H
Foam::surfaceFields.
rho
rho
Definition: readInitialConditions.H:88
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
ddtScheme.H
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::fv::ddtScheme::New
static tmp< ddtScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new ddtScheme created on freestore.
Definition: ddtScheme.C:50