FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

ddtScheme.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include <finiteVolume/fv.H>
00027 #include <OpenFOAM/HashTable.H>
00028 #include <finiteVolume/surfaceInterpolate.H>
00029 
00030 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00031 
00032 namespace Foam
00033 {
00034 
00035 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00036 
00037 namespace fv
00038 {
00039 
00040 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
00041 
00042 template<class Type>
00043 tmp<ddtScheme<Type> > ddtScheme<Type>::New
00044 (
00045     const fvMesh& mesh,
00046     Istream& schemeData
00047 )
00048 {
00049     if (fv::debug)
00050     {
00051         Info<< "ddtScheme<Type>::New(const fvMesh&, Istream&) : "
00052                "constructing ddtScheme<Type>"
00053             << endl;
00054     }
00055 
00056     if (schemeData.eof())
00057     {
00058         FatalIOErrorIn
00059         (
00060             "ddtScheme<Type>::New(const fvMesh&, Istream&)",
00061             schemeData
00062         )   << "Ddt scheme not specified" << endl << endl
00063             << "Valid ddt schemes are :" << endl
00064             << IstreamConstructorTablePtr_->sortedToc()
00065             << exit(FatalIOError);
00066     }
00067 
00068     word schemeName(schemeData);
00069 
00070     typename IstreamConstructorTable::iterator cstrIter =
00071         IstreamConstructorTablePtr_->find(schemeName);
00072 
00073     if (cstrIter == IstreamConstructorTablePtr_->end())
00074     {
00075         FatalIOErrorIn
00076         (
00077             "ddtScheme<Type>::New(const fvMesh&, Istream&)",
00078             schemeData
00079         )   << "unknown ddt scheme " << schemeName << endl << endl
00080             << "Valid ddt schemes are :" << endl
00081             << IstreamConstructorTablePtr_->sortedToc()
00082             << exit(FatalIOError);
00083     }
00084 
00085     return cstrIter()(mesh, schemeData);
00086 }
00087 
00088 
00089 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00090 
00091 template<class Type>
00092 ddtScheme<Type>::~ddtScheme()
00093 {}
00094 
00095 
00096 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00097 
00098 template<class Type>
00099 tmp<surfaceScalarField> ddtScheme<Type>::fvcDdtPhiCoeff
00100 (
00101     const GeometricField<Type, fvPatchField, volMesh>& U,
00102     const fluxFieldType& phi,
00103     const fluxFieldType& phiCorr
00104 )
00105 {
00106     tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
00107       - min
00108         (
00109             mag(phiCorr)
00110            /(mag(phi) + dimensionedScalar("small", phi.dimensions(), SMALL)),
00111             scalar(1)
00112         );
00113 
00114     surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
00115 
00116     forAll (U.boundaryField(), patchi)
00117     {
00118         if (U.boundaryField()[patchi].fixesValue())
00119         {
00120             ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
00121         }
00122     }
00123 
00124     if (debug > 1)
00125     {
00126         Info<< "ddtCouplingCoeff mean max min = "
00127             << gAverage(ddtCouplingCoeff.internalField())
00128             << " " << gMax(ddtCouplingCoeff.internalField())
00129             << " " << gMin(ddtCouplingCoeff.internalField())
00130             << endl;
00131     }
00132 
00133     return tddtCouplingCoeff;
00134 }
00135 
00136 
00137 template<class Type>
00138 tmp<surfaceScalarField> ddtScheme<Type>::fvcDdtPhiCoeff
00139 (
00140     const GeometricField<Type, fvPatchField, volMesh>& U,
00141     const fluxFieldType& phi
00142 )
00143 {
00144     dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
00145 
00146     tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
00147       - min
00148         (
00149             mag(phi - (mesh().Sf() & fvc::interpolate(U)))
00150            /(mag(phi) + dimensionedScalar("small", phi.dimensions(), VSMALL)),
00151            //(rDeltaT*mesh().magSf()/mesh().deltaCoeffs()),
00152             scalar(1)
00153         );
00154 
00155     surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
00156 
00157     forAll (U.boundaryField(), patchi)
00158     {
00159         if (U.boundaryField()[patchi].fixesValue())
00160         {
00161             ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
00162         }
00163     }
00164 
00165     if (debug > 1)
00166     {
00167         Info<< "ddtCouplingCoeff mean max min = "
00168             << gAverage(ddtCouplingCoeff.internalField())
00169             << " " << gMax(ddtCouplingCoeff.internalField())
00170             << " " << gMin(ddtCouplingCoeff.internalField())
00171             << endl;
00172     }
00173 
00174     return tddtCouplingCoeff;
00175 }
00176 
00177 
00178 template<class Type>
00179 tmp<surfaceScalarField> ddtScheme<Type>::fvcDdtPhiCoeff
00180 (
00181     const volScalarField& rho,
00182     const GeometricField<Type, fvPatchField, volMesh>& rhoU,
00183     const fluxFieldType& phi
00184 )
00185 {
00186     dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
00187 
00188     tmp<surfaceScalarField> tddtCouplingCoeff = scalar(1)
00189       - min
00190         (
00191             mag(phi - (mesh().Sf() & fvc::interpolate(rhoU)))
00192            /(
00193                 mag(phi) + dimensionedScalar("small", phi.dimensions(), VSMALL)
00194                 //fvc::interpolate(rho)*rDeltaT
00195                 //*mesh().magSf()/mesh().deltaCoeffs()
00196             ),
00197             scalar(1)
00198         );
00199 
00200     surfaceScalarField& ddtCouplingCoeff = tddtCouplingCoeff();
00201 
00202     forAll (rhoU.boundaryField(), patchi)
00203     {
00204         if (rhoU.boundaryField()[patchi].fixesValue())
00205         {
00206             ddtCouplingCoeff.boundaryField()[patchi] = 0.0;
00207         }
00208     }
00209 
00210     if (debug > 1)
00211     {
00212         Info<< "ddtCouplingCoeff mean max min = "
00213             << gAverage(ddtCouplingCoeff.internalField())
00214             << " " << gMax(ddtCouplingCoeff.internalField())
00215             << " " << gMin(ddtCouplingCoeff.internalField())
00216             << endl;
00217     }
00218 
00219     return tddtCouplingCoeff;
00220 }
00221 
00222 
00223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00224 
00225 } // End namespace fv
00226 
00227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00228 
00229 } // End namespace Foam
00230 
00231 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines