1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2008 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
25 \*---------------------------------------------------------------------------*/
28 #include "HashTable.H"
29 #include "surfaceInterpolate.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
44 tmp<ddtScheme<Type> > ddtScheme<Type>::New
52 Info<< "ddtScheme<Type>::New(const fvMesh&, Istream&) : "
53 "constructing ddtScheme<Type>"
61 "ddtScheme<Type>::New(const fvMesh&, Istream&)",
63 ) << "Ddt scheme not specified" << endl << endl
64 << "Valid ddt schemes are :" << endl
65 << IstreamConstructorTablePtr_->toc()
66 << exit(FatalIOError);
69 word schemeName(schemeData);
71 typename IstreamConstructorTable::iterator cstrIter =
72 IstreamConstructorTablePtr_->find(schemeName);
74 if (cstrIter == IstreamConstructorTablePtr_->end())
78 "ddtScheme<Type>::New(const fvMesh&, Istream&)",
80 ) << "unknown ddt scheme " << schemeName << endl << endl
81 << "Valid ddt schemes are :" << endl
82 << IstreamConstructorTablePtr_->toc()
83 << exit(FatalIOError);
86 return cstrIter()(mesh, schemeData);
90 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 ddtScheme<Type>::~ddtScheme()
97 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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)
111 /(mag(phi) + dimensionedScalar("small", phi.dimensions(), SMALL)),
115 surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
117 forAll (U.boundaryField(), patchi)
119 if (U.boundaryField()[patchi].fixesValue())
121 ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
127 Info<< "ddtCouplingCoeff mean max min = "
128 << gAverage(ddtCouplingCoeff.internalField())
129 << " " << gMax(ddtCouplingCoeff.internalField())
130 << " " << gMin(ddtCouplingCoeff.internalField())
134 return tddtCouplingCoeff;
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)
150 mag(phi - (mesh().Sf() & fvc::interpolate(U)))
151 /(mag(phi) + dimensionedScalar("small", phi.dimensions(), VSMALL)),
152 //(rDeltaT*mesh().magSf()/mesh().deltaCoeffs()),
156 surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
158 forAll (U.boundaryField(), patchi)
160 if (U.boundaryField()[patchi].fixesValue())
162 ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
168 Info<< "ddtCouplingCoeff mean max min = "
169 << gAverage(ddtCouplingCoeff.internalField())
170 << " " << gMax(ddtCouplingCoeff.internalField())
171 << " " << gMin(ddtCouplingCoeff.internalField())
175 return tddtCouplingCoeff;
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)
192 mag(phi - (mesh().Sf() & fvc::interpolate(rhoU)))
194 mag(phi) + dimensionedScalar("small", phi.dimensions(), VSMALL)
195 //fvc::interpolate(rho)*rDeltaT
196 //*mesh().magSf()/mesh().deltaCoeffs()
201 surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
203 forAll (rhoU.boundaryField(), patchi)
205 if (rhoU.boundaryField()[patchi].fixesValue())
207 ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
213 Info<< "ddtCouplingCoeff mean max min = "
214 << gAverage(ddtCouplingCoeff.internalField())
215 << " " << gMax(ddtCouplingCoeff.internalField())
216 << " " << gMin(ddtCouplingCoeff.internalField())
220 return tddtCouplingCoeff;
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 } // End namespace fv
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 } // End namespace Foam
232 // ************************************************************************* //