Loading...
Searching...
No Matches
fvcDDt.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 "fvcDDt.H"
29#include "fvcDiv.H"
30#include "fvMesh.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace fvc
40{
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44template<class Type>
46DDt
47(
50)
51{
54
55 if (phi.mesh().moving())
56 {
57 return ddtDivPhiPsi - fvc::div(phi + phi.mesh().phi())*psi;
58 }
59 else
60 {
61 return ddtDivPhiPsi - fvc::div(phi)*psi;
62 }
63}
64
65
66template<class Type>
67tmp<GeometricField<Type, fvPatchField, volMesh>>
68DDt
69(
70 const tmp<surfaceScalarField>& tphi,
71 const GeometricField<Type, fvPatchField, volMesh>& psi
72)
73{
74 tmp<GeometricField<Type, fvPatchField, volMesh>> DDtPsi
75 (
76 fvc::DDt(tphi(), psi)
77 );
78 tphi.clear();
79 return DDtPsi;
80}
81
82
83// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84
85} // End namespace fvc
86
87// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88
89} // End namespace Foam
90
91// ************************************************************************* //
Generic GeometricField class.
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
const volScalarField & psi
Calculate the substantive (total) derivative.
Calculate the divergence of the given field.
Namespace of functions to calculate explicit derivatives.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition fvcDDt.C:40
Namespace for OpenFOAM.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField