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

fvScalarMatrix.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 "fvScalarMatrix.H"
00027 #include <finiteVolume/zeroGradientFvPatchFields.H>
00028 
00029 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00030 
00031 template<>
00032 void Foam::fvMatrix<Foam::scalar>::setComponentReference
00033 (
00034     const label patchi,
00035     const label facei,
00036     const direction,
00037     const scalar value
00038 )
00039 {
00040     if (psi_.needReference())
00041     {
00042         if (Pstream::master())
00043         {
00044             internalCoeffs_[patchi][facei] +=
00045                 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
00046 
00047             boundaryCoeffs_[patchi][facei] +=
00048                 diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]
00049                *value;
00050         }
00051     }
00052 }
00053 
00054 
00055 template<>
00056 Foam::autoPtr<Foam::fvMatrix<Foam::scalar>::fvSolver>
00057 Foam::fvMatrix<Foam::scalar>::solver
00058 (
00059     const dictionary& solverControls
00060 )
00061 {
00062     if (debug)
00063     {
00064         Info<< "fvMatrix<scalar>::solver(const dictionary& solverControls) : "
00065                "solver for fvMatrix<scalar>"
00066             << endl;
00067     }
00068 
00069     scalarField saveDiag = diag();
00070     addBoundaryDiag(diag(), 0);
00071 
00072     autoPtr<fvMatrix<scalar>::fvSolver> solverPtr
00073     (
00074         new fvMatrix<scalar>::fvSolver
00075         (
00076             *this,
00077             lduMatrix::solver::New
00078             (
00079                 psi_.name(),
00080                 *this,
00081                 boundaryCoeffs_,
00082                 internalCoeffs_,
00083                 psi_.boundaryField().interfaces(),
00084                 solverControls
00085             )
00086         )
00087     );
00088 
00089     diag() = saveDiag;
00090 
00091     return solverPtr;
00092 }
00093 
00094 
00095 template<>
00096 Foam::lduMatrix::solverPerformance Foam::fvMatrix<Foam::scalar>::fvSolver::solve
00097 (
00098     const dictionary& solverControls
00099 )
00100 {
00101     scalarField saveDiag = fvMat_.diag();
00102     fvMat_.addBoundaryDiag(fvMat_.diag(), 0);
00103 
00104     scalarField totalSource = fvMat_.source();
00105     fvMat_.addBoundarySource(totalSource, false);
00106 
00107     // assign new solver controls
00108     solver_->read(solverControls);
00109 
00110     lduMatrix::solverPerformance solverPerf =
00111         solver_->solve(fvMat_.psi().internalField(), totalSource);
00112 
00113     solverPerf.print();
00114 
00115     fvMat_.diag() = saveDiag;
00116 
00117     fvMat_.psi().correctBoundaryConditions();
00118 
00119     return solverPerf;
00120 }
00121 
00122 
00123 template<>
00124 Foam::lduMatrix::solverPerformance Foam::fvMatrix<Foam::scalar>::solve
00125 (
00126     const dictionary& solverControls
00127 )
00128 {
00129     if (debug)
00130     {
00131         Info<< "fvMatrix<scalar>::solve(const dictionary& solverControls) : "
00132                "solving fvMatrix<scalar>"
00133             << endl;
00134     }
00135 
00136     scalarField saveDiag = diag();
00137     addBoundaryDiag(diag(), 0);
00138 
00139     scalarField totalSource = source_;
00140     addBoundarySource(totalSource, false);
00141 
00142     // Solver call
00143     lduMatrix::solverPerformance solverPerf = lduMatrix::solver::New
00144     (
00145         psi_.name(),
00146         *this,
00147         boundaryCoeffs_,
00148         internalCoeffs_,
00149         psi_.boundaryField().interfaces(),
00150         solverControls
00151     )->solve(psi_.internalField(), totalSource);
00152 
00153     solverPerf.print();
00154 
00155     diag() = saveDiag;
00156 
00157     psi_.correctBoundaryConditions();
00158 
00159     return solverPerf;
00160 }
00161 
00162 
00163 template<>
00164 Foam::tmp<Foam::scalarField> Foam::fvMatrix<Foam::scalar>::residual() const
00165 {
00166     scalarField boundaryDiag(psi_.size(), 0.0);
00167     addBoundaryDiag(boundaryDiag, 0);
00168 
00169     tmp<scalarField> tres
00170     (
00171         lduMatrix::residual
00172         (
00173             psi_.internalField(),
00174             source_ - boundaryDiag*psi_.internalField(),
00175             boundaryCoeffs_,
00176             psi_.boundaryField().interfaces(),
00177             0
00178         )
00179     );
00180 
00181     addBoundarySource(tres());
00182 
00183     return tres;
00184 }
00185 
00186 
00187 template<>
00188 Foam::tmp<Foam::volScalarField> Foam::fvMatrix<Foam::scalar>::H() const
00189 {
00190     tmp<volScalarField> tHphi
00191     (
00192         new volScalarField
00193         (
00194             IOobject
00195             (
00196                 "H("+psi_.name()+')',
00197                 psi_.instance(),
00198                 psi_.mesh(),
00199                 IOobject::NO_READ,
00200                 IOobject::NO_WRITE
00201             ),
00202             psi_.mesh(),
00203             dimensions_/dimVol,
00204             zeroGradientFvPatchScalarField::typeName
00205         )
00206     );
00207     volScalarField& Hphi = tHphi();
00208 
00209     Hphi.internalField() = (lduMatrix::H(psi_.internalField()) + source_);
00210     addBoundarySource(Hphi.internalField());
00211 
00212     Hphi.internalField() /= psi_.mesh().V();
00213     Hphi.correctBoundaryConditions();
00214 
00215     return tHphi;
00216 }
00217 
00218 
00219 template<>
00220 Foam::tmp<Foam::volScalarField> Foam::fvMatrix<Foam::scalar>::H1() const
00221 {
00222     tmp<volScalarField> tH1
00223     (
00224         new volScalarField
00225         (
00226             IOobject
00227             (
00228                 "H(1)",
00229                 psi_.instance(),
00230                 psi_.mesh(),
00231                 IOobject::NO_READ,
00232                 IOobject::NO_WRITE
00233             ),
00234             psi_.mesh(),
00235             dimensions_/(dimVol*psi_.dimensions()),
00236             zeroGradientFvPatchScalarField::typeName
00237         )
00238     );
00239     volScalarField& H1_ = tH1();
00240 
00241     H1_.internalField() = lduMatrix::H1();
00242     //addBoundarySource(Hphi.internalField());
00243 
00244     H1_.internalField() /= psi_.mesh().V();
00245     H1_.correctBoundaryConditions();
00246 
00247     return tH1;
00248 }
00249 
00250 
00251 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines