1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
26 Foam::fv::boundedBackwardDdtScheme
29 Second-order bounded-backward-differencing ddt using the current and
30 two previous time-step values.
33 boundedBackwardDdtScheme.C
35 \*---------------------------------------------------------------------------*/
37 #ifndef boundedBackwardDdtScheme_H
38 #define boundedBackwardDdtScheme_H
40 #include "ddtScheme.H"
41 #include "fvMatrices.H"
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 /*---------------------------------------------------------------------------*\
54 Class boundedBackwardDdtScheme Declaration
55 \*---------------------------------------------------------------------------*/
57 class boundedBackwardDdtScheme
59 public fv::ddtScheme<scalar>
61 // Private Member Functions
63 //- Return the current time-step
64 scalar deltaT_() const;
66 //- Return the previous time-step
67 scalar deltaT0_() const;
69 //- Return the previous time-step or GREAT if the old timestep field
70 // wasn't available in which case Euler ddt is used
71 template<class GeoField>
72 scalar deltaT0_(const GeoField& vf) const
74 if (vf.oldTime().timeIndex() == vf.oldTime().oldTime().timeIndex())
85 //- Disallow default bitwise copy construct
86 boundedBackwardDdtScheme(const boundedBackwardDdtScheme&);
88 //- Disallow default bitwise assignment
89 void operator=(const boundedBackwardDdtScheme&);
94 //- Runtime type information
95 TypeName("boundedBackward");
100 //- Construct from mesh
101 boundedBackwardDdtScheme(const fvMesh& mesh)
103 ddtScheme<scalar>(mesh)
106 //- Construct from mesh and Istream
107 boundedBackwardDdtScheme(const fvMesh& mesh, Istream& is)
109 ddtScheme<scalar>(mesh, is)
115 //- Return mesh reference
116 const fvMesh& mesh() const
118 return fv::ddtScheme<scalar>::mesh();
121 tmp<volScalarField> fvcDdt
123 const dimensionedScalar&
126 tmp<volScalarField> fvcDdt
128 const volScalarField&
131 tmp<volScalarField> fvcDdt
133 const dimensionedScalar&,
134 const volScalarField&
137 tmp<volScalarField> fvcDdt
139 const volScalarField&,
140 const volScalarField&
143 tmp<fvScalarMatrix> fvmDdt
148 tmp<fvScalarMatrix> fvmDdt
150 const dimensionedScalar&,
154 tmp<fvScalarMatrix> fvmDdt
156 const volScalarField&,
160 tmp<surfaceScalarField> fvcDdtPhiCorr
162 const volScalarField& rA,
163 const volScalarField& U,
164 const surfaceScalarField& phi
167 tmp<surfaceScalarField> fvcDdtPhiCorr
169 const volScalarField& rA,
170 const volScalarField& rho,
171 const volScalarField& U,
172 const surfaceScalarField& phi
175 tmp<surfaceScalarField> meshPhi
177 const volScalarField&
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 } // End namespace fv
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 } // End namespace Foam
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 // ************************************************************************* //