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-------------------------------------------------------------------------------
10License
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
35namespace Foam
36{
37
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40namespace fvm
41{
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45template<class Type>
46tmp<fvMatrix<Type>>
48(
50)
51{
53 (
54 vf.mesh(),
55 vf.mesh().ddtScheme("ddt(" + vf.name() + ')')
56 ).ref().fvmDdt(vf);
57}
58
59
60template<class Type>
63(
64 const one&,
66)
67{
68 return ddt(vf);
69}
70
71
72template<class Type>
75(
78)
79{
81 (
82 vf.mesh(),
83 vf.mesh().ddtScheme("ddt(" + rho.name() + ',' + vf.name() + ')')
84 ).ref().fvmDdt(rho, vf);
85}
86
87
88template<class Type>
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
104template<class Type>
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
127template<class Type>
130(
131 const one&,
132 const one&,
134)
135{
136 return ddt(vf);
137}
138
139
140template<class Type>
143(
144 const one&,
145 const volScalarField& rho,
147)
148{
149 return ddt(rho, vf);
150}
151
152
153template<class Type>
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// ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
A class for managing temporary objects.
Definition: tmp.H:65
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Namespace for OpenFOAM.
volScalarField & alpha
Foam::surfaceFields.