steadyStateD2dt2Scheme.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
29#include "fvcDiv.H"
30#include "fvMatrices.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace fv
40{
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44template<class Type>
45tmp<GeometricField<Type, fvPatchField, volMesh>>
47(
49)
50{
52 (
54 (
56 (
57 "d2dt2("+vf.name()+')',
58 mesh().time().timeName(),
59 mesh(),
62 ),
63 mesh(),
65 )
66 );
67}
68
69
70template<class Type>
73(
74 const volScalarField& rho,
76)
77{
79 (
81 (
82 "d2dt2("+rho.name()+','+vf.name()+')',
83 mesh().time().timeName(),
84 mesh(),
87 ),
88 mesh(),
90 (
91 rho.dimensions()*vf.dimensions()/dimTime/dimTime,
92 Zero
93 )
94 );
95}
96
97
98template<class Type>
101(
103)
104{
106 (
108 (
109 vf,
111 )
112 );
113
114 return tfvm;
115}
116
117
118template<class Type>
121(
122 const dimensionedScalar& rho,
124)
125{
127 (
129 (
130 vf,
131 rho.dimensions()*vf.dimensions()*dimVol/dimTime/dimTime
132 )
133 );
134
135 return tfvm;
136}
137
138
139template<class Type>
142(
143 const volScalarField& rho,
145)
146{
148 (
150 (
151 vf,
152 rho.dimensions()*vf.dimensions()*dimVol/dimTime/dimTime
153 )
154 );
155
156 return tfvm;
157}
158
159
160// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161
162} // End namespace fv
163
164// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165
166} // End namespace Foam
167
168// ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
Generic dimensioned Type class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
tmp< fvMatrix< Type > > fvmD2dt2(const GeometricField< Type, fvPatchField, volMesh > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcD2dt2(const GeometricField< Type, fvPatchField, volMesh > &)
A class for managing temporary objects.
Definition: tmp.H:65
dynamicFvMesh & mesh
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the divergence of the given field.
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:64
labelList fv(nPoints)