All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
fvcDDt.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 "fvcDDt.H"
27 #include "fvcDiv.H"
28 #include "fvMesh.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 DDt
45 (
46  const surfaceScalarField& phi,
48 )
49 {
51  = fvc::ddt(psi) + fvc::div(phi, psi);
52 
53  if (phi.mesh().moving())
54  {
55  return ddtDivPhiPsi - fvc::div(phi + phi.mesh().phi())*psi;
56  }
57  else
58  {
59  return ddtDivPhiPsi - fvc::div(phi)*psi;
60  }
61 }
62 
63 
64 template<class Type>
66 DDt
67 (
68  const tmp<surfaceScalarField>& tphi,
70 )
71 {
73  (
74  fvc::DDt(tphi(), psi)
75  );
76  tphi.clear();
77  return DDtPsi;
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 } // End namespace fvc
84 
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87 } // End namespace Foam
88 
89 // ************************************************************************* //
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:45
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Generic GeometricField class.
Calculate the substantive (total) derivative.
phi
Definition: pEqn.H:104
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
Calculate the divergence of the given field.
const Mesh & mesh() const
Return mesh.
const volScalarField & psi
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.