steadyStateD2dt2Scheme.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-2024 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 "steadyStateD2dt2Scheme.H"
27 #include "fvcDiv.H"
28 #include "fvMatrices.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace fv
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 template<class Type>
43 tmp<VolField<Type>>
45 (
46  const VolField<Type>& vf
47 )
48 {
49  return VolField<Type>::New
50  (
51  "d2dt2("+vf.name()+')',
52  mesh(),
54  (
55  "0",
57  Zero
58  )
59  );
60 }
61 
62 
63 template<class Type>
66 (
67  const volScalarField& rho,
68  const VolField<Type>& vf
69 )
70 {
71  return VolField<Type>::New
72  (
73  "d2dt2("+rho.name()+','+vf.name()+')',
74  mesh(),
76  (
77  "0",
78  rho.dimensions()*vf.dimensions()/dimTime/dimTime,
79  Zero
80  )
81  );
82 }
83 
84 
85 template<class Type>
88 (
89  const VolField<Type>& vf
90 )
91 {
92  tmp<fvMatrix<Type>> tfvm
93  (
94  new fvMatrix<Type>
95  (
96  vf,
98  )
99  );
100 
101  return tfvm;
102 }
103 
104 
105 template<class Type>
108 (
109  const dimensionedScalar& rho,
110  const VolField<Type>& vf
111 )
112 {
113  tmp<fvMatrix<Type>> tfvm
114  (
115  new fvMatrix<Type>
116  (
117  vf,
118  rho.dimensions()*vf.dimensions()*dimVolume/dimTime/dimTime
119  )
120  );
121 
122  return tfvm;
123 }
124 
125 
126 template<class Type>
129 (
130  const volScalarField& rho,
131  const VolField<Type>& vf
132 )
133 {
134  tmp<fvMatrix<Type>> tfvm
135  (
136  new fvMatrix<Type>
137  (
138  vf,
139  rho.dimensions()*vf.dimensions()*dimVolume/dimTime/dimTime
140  )
141  );
142 
143  return tfvm;
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 
149 } // End namespace fv
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 } // End namespace Foam
154 
155 // ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
const word & name() const
Return name.
Definition: IOobject.H:310
Generic dimensioned Type class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
tmp< fvMatrix< Type > > fvmD2dt2(const VolField< Type > &)
tmp< VolField< Type > > fvcD2dt2(const VolField< Type > &)
A class for managing temporary objects.
Definition: tmp.H:55
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the divergence of the given field.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const dimensionSet dimTime
const dimensionSet dimVolume
labelList fv(nPoints)