fvcD2dt2.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-2018 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 "fvcD2dt2.H"
27 #include "fvMesh.H"
28 #include "d2dt2Scheme.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace fvc
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<class Type>
43 tmp<GeometricField<Type, fvPatchField, volMesh>>
44 d2dt2
45 (
47 )
48 {
50  (
51  vf.mesh(),
52  vf.mesh().ddtScheme("d2dt2(" + vf.name() + ')')
53  ).ref().fvcD2dt2(vf);
54 }
55 
56 
57 template<class Type>
59 d2dt2
60 (
61  const volScalarField& rho,
63 )
64 {
66  (
67  vf.mesh(),
68  vf.mesh().ddtScheme("d2dt2(" + rho.name() + ',' + vf.name() + ')')
69  ).ref().fvcD2dt2(rho, vf);
70 }
71 
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 } // End namespace fvc
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 } // End namespace Foam
80 
81 // ************************************************************************* //
tmp< GeometricField< Type, fvPatchField, volMesh > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcD2dt2.C:45
const word & name() const
Return name.
Definition: IOobject.H:303
Generic GeometricField class.
Calculate the second temporal derivative.
const Mesh & mesh() const
Return mesh.
static tmp< d2dt2Scheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new d2dt2Scheme created on freestore.
Definition: d2dt2Scheme.C:46
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.