initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / finiteVolume / ddtSchemes / ddtScheme / ddtScheme.C
blob84d620759ffc90b6223a2b59d3d06d8aefcb8ab4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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
25 \*---------------------------------------------------------------------------*/
27 #include "fv.H"
28 #include "HashTable.H"
29 #include "surfaceInterpolate.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace Foam
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 namespace fv
41 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
43 template<class Type>
44 tmp<ddtScheme<Type> > ddtScheme<Type>::New
46     const fvMesh& mesh,
47     Istream& schemeData
50     if (fv::debug)
51     {
52         Info<< "ddtScheme<Type>::New(const fvMesh&, Istream&) : "
53                "constructing ddtScheme<Type>"
54             << endl;
55     }
57     if (schemeData.eof())
58     {
59         FatalIOErrorIn
60         (
61             "ddtScheme<Type>::New(const fvMesh&, Istream&)",
62             schemeData
63         )   << "Ddt scheme not specified" << endl << endl
64             << "Valid ddt schemes are :" << endl
65             << IstreamConstructorTablePtr_->toc()
66             << exit(FatalIOError);
67     }
69     word schemeName(schemeData);
71     typename IstreamConstructorTable::iterator cstrIter =
72         IstreamConstructorTablePtr_->find(schemeName);
74     if (cstrIter == IstreamConstructorTablePtr_->end())
75     {
76         FatalIOErrorIn
77         (
78             "ddtScheme<Type>::New(const fvMesh&, Istream&)",
79             schemeData
80         )   << "unknown ddt scheme " << schemeName << endl << endl
81             << "Valid ddt schemes are :" << endl
82             << IstreamConstructorTablePtr_->toc()
83             << exit(FatalIOError);
84     }
86     return cstrIter()(mesh, schemeData);
90 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
92 template<class Type>
93 ddtScheme<Type>::~ddtScheme()
97 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99 template<class Type>
100 tmp<surfaceScalarField> ddtScheme<Type>::fvcDdtPhiCoeff
102     const GeometricField<Type, fvPatchField, volMesh>& U,
103     const fluxFieldType& phi,
104     const fluxFieldType& phiCorr
107     tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
108       - min
109         (
110             mag(phiCorr)
111            /(mag(phi) + dimensionedScalar("small", phi.dimensions(), SMALL)),
112             scalar(1)
113         );
115     surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
117     forAll (U.boundaryField(), patchi)
118     {
119         if (U.boundaryField()[patchi].fixesValue())
120         {
121             ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
122         }
123     }
125     if (debug > 1)
126     {
127         Info<< "ddtCouplingCoeff mean max min = "
128             << gAverage(ddtCouplingCoeff.internalField())
129             << " " << gMax(ddtCouplingCoeff.internalField())
130             << " " << gMin(ddtCouplingCoeff.internalField())
131             << endl;
132     }
134     return tddtCouplingCoeff;
138 template<class Type>
139 tmp<surfaceScalarField> ddtScheme<Type>::fvcDdtPhiCoeff
141     const GeometricField<Type, fvPatchField, volMesh>& U,
142     const fluxFieldType& phi
145     dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
147     tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
148       - min
149         (
150             mag(phi - (mesh().Sf() & fvc::interpolate(U)))
151            /(mag(phi) + dimensionedScalar("small", phi.dimensions(), VSMALL)),
152            //(rDeltaT*mesh().magSf()/mesh().deltaCoeffs()),
153             scalar(1)
154         );
156     surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
158     forAll (U.boundaryField(), patchi)
159     {
160         if (U.boundaryField()[patchi].fixesValue())
161         {
162             ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
163         }
164     }
166     if (debug > 1)
167     {
168         Info<< "ddtCouplingCoeff mean max min = "
169             << gAverage(ddtCouplingCoeff.internalField())
170             << " " << gMax(ddtCouplingCoeff.internalField())
171             << " " << gMin(ddtCouplingCoeff.internalField())
172             << endl;
173     }
175     return tddtCouplingCoeff;
179 template<class Type>
180 tmp<surfaceScalarField> ddtScheme<Type>::fvcDdtPhiCoeff
182     const volScalarField& rho,
183     const GeometricField<Type, fvPatchField, volMesh>& rhoU,
184     const fluxFieldType& phi
187     dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
189     tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
190       - min
191         (
192             mag(phi - (mesh().Sf() & fvc::interpolate(rhoU)))
193            /(
194                 mag(phi) + dimensionedScalar("small", phi.dimensions(), VSMALL)
195                 //fvc::interpolate(rho)*rDeltaT
196                 //*mesh().magSf()/mesh().deltaCoeffs()
197             ),
198             scalar(1)
199         );
201     surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
203     forAll (rhoU.boundaryField(), patchi)
204     {
205         if (rhoU.boundaryField()[patchi].fixesValue())
206         {
207             ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
208         }
209     }
211     if (debug > 1)
212     {
213         Info<< "ddtCouplingCoeff mean max min = "
214             << gAverage(ddtCouplingCoeff.internalField())
215             << " " << gMax(ddtCouplingCoeff.internalField())
216             << " " << gMin(ddtCouplingCoeff.internalField())
217             << endl;
218     }
220     return tddtCouplingCoeff;
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 } // End namespace fv
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 } // End namespace Foam
232 // ************************************************************************* //