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

lduMatrix.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 "lduMatrix.H"
00027 #include <OpenFOAM/IOstreams.H>
00028 
00029 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033     defineTypeNameAndDebug(lduMatrix, 1);
00034 }
00035 
00036 const Foam::scalar Foam::lduMatrix::great_ = 1.0e+20;
00037 const Foam::scalar Foam::lduMatrix::small_ = 1.0e-20;
00038 
00039 
00040 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00041 
00042 Foam::lduMatrix::lduMatrix(const lduMesh& mesh)
00043 :
00044     lduMesh_(mesh),
00045     lowerPtr_(NULL),
00046     diagPtr_(NULL),
00047     upperPtr_(NULL)
00048 {}
00049 
00050 
00051 Foam::lduMatrix::lduMatrix(const lduMatrix& A)
00052 :
00053     lduMesh_(A.lduMesh_),
00054     lowerPtr_(NULL),
00055     diagPtr_(NULL),
00056     upperPtr_(NULL)
00057 {
00058     if (A.lowerPtr_)
00059     {
00060         lowerPtr_ = new scalarField(*(A.lowerPtr_));
00061     }
00062 
00063     if (A.diagPtr_)
00064     {
00065         diagPtr_ = new scalarField(*(A.diagPtr_));
00066     }
00067 
00068     if (A.upperPtr_)
00069     {
00070         upperPtr_ = new scalarField(*(A.upperPtr_));
00071     }
00072 }
00073 
00074 
00075 Foam::lduMatrix::lduMatrix(lduMatrix& A, bool reUse)
00076 :
00077     lduMesh_(A.lduMesh_),
00078     lowerPtr_(NULL),
00079     diagPtr_(NULL),
00080     upperPtr_(NULL)
00081 {
00082     if (reUse)
00083     {
00084         if (A.lowerPtr_)
00085         {
00086             lowerPtr_ = A.lowerPtr_;
00087             A.lowerPtr_ = NULL;
00088         }
00089 
00090         if (A.diagPtr_)
00091         {
00092             diagPtr_ = A.diagPtr_;
00093             A.diagPtr_ = NULL;
00094         }
00095 
00096         if (A.upperPtr_)
00097         {
00098             upperPtr_ = A.upperPtr_;
00099             A.upperPtr_ = NULL;
00100         }
00101     }
00102     else
00103     {
00104         if (A.lowerPtr_)
00105         {
00106             lowerPtr_ = new scalarField(*(A.lowerPtr_));
00107         }
00108 
00109         if (A.diagPtr_)
00110         {
00111             diagPtr_ = new scalarField(*(A.diagPtr_));
00112         }
00113 
00114         if (A.upperPtr_)
00115         {
00116             upperPtr_ = new scalarField(*(A.upperPtr_));
00117         }
00118     }
00119 }
00120 
00121 
00122 Foam::lduMatrix::lduMatrix
00123 (
00124     const lduMesh& mesh,
00125     Istream& is
00126 )
00127 :
00128     lduMesh_(mesh),
00129     lowerPtr_(new scalarField(is)),
00130     diagPtr_(new scalarField(is)),
00131     upperPtr_(new scalarField(is))
00132 {}
00133 
00134 
00135 Foam::lduMatrix::~lduMatrix()
00136 {
00137     if (lowerPtr_)
00138     {
00139         delete lowerPtr_;
00140     }
00141 
00142     if (diagPtr_)
00143     {
00144         delete diagPtr_;
00145     }
00146 
00147     if (upperPtr_)
00148     {
00149         delete upperPtr_;
00150     }
00151 }
00152 
00153 
00154 Foam::scalarField& Foam::lduMatrix::lower()
00155 {
00156     if (!lowerPtr_)
00157     {
00158         if (upperPtr_)
00159         {
00160             lowerPtr_ = new scalarField(*upperPtr_);
00161         }
00162         else
00163         {
00164             lowerPtr_ = new scalarField(lduAddr().lowerAddr().size(), 0.0);
00165         }
00166     }
00167 
00168     return *lowerPtr_;
00169 }
00170 
00171 
00172 Foam::scalarField& Foam::lduMatrix::diag()
00173 {
00174     if (!diagPtr_)
00175     {
00176         diagPtr_ = new scalarField(lduAddr().size(), 0.0);
00177     }
00178 
00179     return *diagPtr_;
00180 }
00181 
00182 
00183 Foam::scalarField& Foam::lduMatrix::upper()
00184 {
00185     if (!upperPtr_)
00186     {
00187         if (lowerPtr_)
00188         {
00189             upperPtr_ = new scalarField(*lowerPtr_);
00190         }
00191         else
00192         {
00193             upperPtr_ = new scalarField(lduAddr().lowerAddr().size(), 0.0);
00194         }
00195     }
00196 
00197     return *upperPtr_;
00198 }
00199 
00200 
00201 const Foam::scalarField& Foam::lduMatrix::lower() const
00202 {
00203     if (!lowerPtr_ && !upperPtr_)
00204     {
00205         FatalErrorIn("lduMatrix::lower() const")
00206             << "lowerPtr_ or upperPtr_ unallocated"
00207             << abort(FatalError);
00208     }
00209 
00210     if (lowerPtr_)
00211     {
00212         return *lowerPtr_;
00213     }
00214     else
00215     {
00216         return *upperPtr_;
00217     }
00218 }
00219 
00220 
00221 const Foam::scalarField& Foam::lduMatrix::diag() const
00222 {
00223     if (!diagPtr_)
00224     {
00225         FatalErrorIn("const scalarField& lduMatrix::diag() const")
00226             << "diagPtr_ unallocated"
00227             << abort(FatalError);
00228     }
00229 
00230     return *diagPtr_;
00231 }
00232 
00233 
00234 const Foam::scalarField& Foam::lduMatrix::upper() const
00235 {
00236     if (!lowerPtr_ && !upperPtr_)
00237     {
00238         FatalErrorIn("lduMatrix::upper() const")
00239             << "lowerPtr_ or upperPtr_ unallocated"
00240             << abort(FatalError);
00241     }
00242 
00243     if (upperPtr_)
00244     {
00245         return *upperPtr_;
00246     }
00247     else
00248     {
00249         return *lowerPtr_;
00250     }
00251 }
00252 
00253 
00254 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
00255 
00256 Foam::Ostream& Foam::operator<<(Ostream& os, const lduMatrix& ldum)
00257 {
00258     if (ldum.lowerPtr_)
00259     {
00260         os  << "Lower triangle = "
00261             << *ldum.lowerPtr_
00262             << endl << endl;
00263     }
00264 
00265     if (ldum.diagPtr_)
00266     {
00267         os  << "diagonal = "
00268             << *ldum.diagPtr_
00269             << endl << endl;
00270     }
00271 
00272     if (ldum.upperPtr_)
00273     {
00274         os  << "Upper triangle = "
00275             << *ldum.upperPtr_
00276             << endl << endl;
00277     }
00278 
00279     os.check("Ostream& operator<<(Ostream&, const lduMatrix&");
00280 
00281     return os;
00282 }
00283 
00284 
00285 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines