fvmD2dt2.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 | Copyright (C) 2011-2016 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 "d2dt2Scheme.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 d2dt2
46 (
48 )
49 {
51  (
52  vf.mesh(),
53  vf.mesh().d2dt2Scheme("d2dt2(" + vf.name() + ')')
54  ).ref().fvmD2dt2(vf);
55 }
56 
57 
58 template<class Type>
60 d2dt2
61 (
62  const dimensionedScalar& rho,
64 )
65 {
67  (
68  vf.mesh(),
69  vf.mesh().d2dt2Scheme("d2dt2(" + rho.name() + ',' + vf.name() + ')')
70  ).ref().fvmD2dt2(rho, vf);
71 }
72 
73 
74 template<class Type>
76 d2dt2
77 (
78  const volScalarField& rho,
80 )
81 {
83  (
84  vf.mesh(),
85  vf.mesh().d2dt2Scheme("d2dt2(" + rho.name() + ',' + vf.name() + ')')
86  ).ref().fvmD2dt2(rho, vf);
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 
92 } // End namespace fvm
93 
94 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
95 
96 } // End namespace Foam
97 
98 // ************************************************************************* //
Foam::surfaceFields.
Generic GeometricField class.
const word & name() const
Return const reference to name.
tmp< fvMatrix< Type > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmD2dt2.C:46
static tmp< d2dt2Scheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new d2dt2Scheme created on freestore.
Definition: d2dt2Scheme.C:46
const Mesh & mesh() const
Return mesh.
A class for managing temporary objects.
Definition: PtrList.H:54
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.