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

lduMatrixOperations.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 Description
00025     lduMatrix member operations.
00026 
00027 \*---------------------------------------------------------------------------*/
00028 
00029 #include <OpenFOAM/lduMatrix.H>
00030 
00031 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00032 
00033 void Foam::lduMatrix::sumDiag()
00034 {
00035     const scalarField& Lower = const_cast<const lduMatrix&>(*this).lower();
00036     const scalarField& Upper = const_cast<const lduMatrix&>(*this).upper();
00037     scalarField& Diag = diag();
00038 
00039     const unallocLabelList& l = lduAddr().lowerAddr();
00040     const unallocLabelList& u = lduAddr().upperAddr();
00041 
00042     for (register label face=0; face<l.size(); face++)
00043     {
00044         Diag[l[face]] += Lower[face];
00045         Diag[u[face]] += Upper[face];
00046     }
00047 }
00048 
00049 
00050 void Foam::lduMatrix::negSumDiag()
00051 {
00052     const scalarField& Lower = const_cast<const lduMatrix&>(*this).lower();
00053     const scalarField& Upper = const_cast<const lduMatrix&>(*this).upper();
00054     scalarField& Diag = diag();
00055 
00056     const unallocLabelList& l = lduAddr().lowerAddr();
00057     const unallocLabelList& u = lduAddr().upperAddr();
00058 
00059     for (register label face=0; face<l.size(); face++)
00060     {
00061         Diag[l[face]] -= Lower[face];
00062         Diag[u[face]] -= Upper[face];
00063     }
00064 }
00065 
00066 
00067 void Foam::lduMatrix::sumMagOffDiag
00068 (
00069     scalarField& sumOff
00070 ) const
00071 {
00072     const scalarField& Lower = const_cast<const lduMatrix&>(*this).lower();
00073     const scalarField& Upper = const_cast<const lduMatrix&>(*this).upper();
00074 
00075     const unallocLabelList& l = lduAddr().lowerAddr();
00076     const unallocLabelList& u = lduAddr().upperAddr();
00077 
00078     for (register label face = 0; face < l.size(); face++)
00079     {
00080         sumOff[u[face]] += mag(Lower[face]);
00081         sumOff[l[face]] += mag(Upper[face]);
00082     }
00083 }
00084 
00085 
00086 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00087 
00088 void Foam::lduMatrix::operator=(const lduMatrix& A)
00089 {
00090     if (this == &A)
00091     {
00092         FatalError
00093             << "lduMatrix::operator=(const lduMatrix&) : "
00094             << "attempted assignment to self"
00095             << abort(FatalError);
00096     }
00097 
00098     if (A.lowerPtr_)
00099     {
00100         lower() = A.lower();
00101     }
00102     else if (lowerPtr_)
00103     {
00104         delete lowerPtr_;
00105         lowerPtr_ = NULL;
00106     }
00107 
00108     if (A.upperPtr_)
00109     {
00110         upper() = A.upper();
00111     }
00112     else if (upperPtr_)
00113     {
00114         delete upperPtr_;
00115         upperPtr_ = NULL;
00116     }
00117 
00118     if (A.diagPtr_)
00119     {
00120         diag() = A.diag();
00121     }
00122 }
00123 
00124 
00125 void Foam::lduMatrix::negate()
00126 {
00127     if (lowerPtr_)
00128     {
00129         lowerPtr_->negate();
00130     }
00131 
00132     if (upperPtr_)
00133     {
00134         upperPtr_->negate();
00135     }
00136 
00137     if (diagPtr_)
00138     {
00139         diagPtr_->negate();
00140     }
00141 }
00142 
00143 
00144 void Foam::lduMatrix::operator+=(const lduMatrix& A)
00145 {
00146     if (A.diagPtr_)
00147     {
00148         diag() += A.diag();
00149     }
00150 
00151     if (symmetric() && A.symmetric())
00152     {
00153         upper() += A.upper();
00154     }
00155     else if (symmetric() && A.asymmetric())
00156     {
00157         if (upperPtr_)
00158         {
00159             lower();
00160         }
00161         else
00162         {
00163             upper();
00164         }
00165 
00166         upper() += A.upper();
00167         lower() += A.lower();
00168     }
00169     else if (asymmetric() && A.symmetric())
00170     {
00171         if (A.upperPtr_)
00172         {
00173             lower() += A.upper();
00174             upper() += A.upper();
00175         }
00176         else
00177         {
00178             lower() += A.lower();
00179             upper() += A.lower();
00180         }
00181 
00182     }
00183     else if (asymmetric() && A.asymmetric())
00184     {
00185         lower() += A.lower();
00186         upper() += A.upper();
00187     }
00188     else if (diagonal())
00189     {
00190         if (A.upperPtr_)
00191         {
00192             upper() = A.upper();
00193         }
00194 
00195         if (A.lowerPtr_)
00196         {
00197             lower() = A.lower();
00198         }
00199     }
00200     else if (A.diagonal())
00201     {
00202     }
00203     else
00204     {
00205         FatalErrorIn("lduMatrix::operator+=(const lduMatrix& A)")
00206             << "Unknown matrix type combination"
00207             << abort(FatalError);
00208     }
00209 }
00210 
00211 
00212 void Foam::lduMatrix::operator-=(const lduMatrix& A)
00213 {
00214     if (A.diagPtr_)
00215     {
00216         diag() -= A.diag();
00217     }
00218 
00219     if (symmetric() && A.symmetric())
00220     {
00221         upper() -= A.upper();
00222     }
00223     else if (symmetric() && A.asymmetric())
00224     {
00225         if (upperPtr_)
00226         {
00227             lower();
00228         }
00229         else
00230         {
00231             upper();
00232         }
00233 
00234         upper() -= A.upper();
00235         lower() -= A.lower();
00236     }
00237     else if (asymmetric() && A.symmetric())
00238     {
00239         if (A.upperPtr_)
00240         {
00241             lower() -= A.upper();
00242             upper() -= A.upper();
00243         }
00244         else
00245         {
00246             lower() -= A.lower();
00247             upper() -= A.lower();
00248         }
00249 
00250     }
00251     else if (asymmetric() && A.asymmetric())
00252     {
00253         lower() -= A.lower();
00254         upper() -= A.upper();
00255     }
00256     else if (diagonal())
00257     {
00258         if (A.upperPtr_)
00259         {
00260             upper() = -A.upper();
00261         }
00262 
00263         if (A.lowerPtr_)
00264         {
00265             lower() = -A.lower();
00266         }
00267     }
00268     else if (A.diagonal())
00269     {
00270     }
00271     else
00272     {
00273         FatalErrorIn("lduMatrix::operator-=(const lduMatrix& A)")
00274             << "Unknown matrix type combination"
00275             << abort(FatalError);
00276     }
00277 }
00278 
00279 
00280 void Foam::lduMatrix::operator*=(const scalarField& sf)
00281 {
00282     if (diagPtr_)
00283     {
00284         *diagPtr_ *= sf;
00285     }
00286 
00287     if (upperPtr_)
00288     {
00289         scalarField& upper = *upperPtr_;
00290 
00291         const unallocLabelList& l = lduAddr().lowerAddr();
00292 
00293         for (register label face=0; face<upper.size(); face++)
00294         {
00295             upper[face] *= sf[l[face]];
00296         }
00297     }
00298 
00299     if (lowerPtr_)
00300     {
00301         scalarField& lower = *lowerPtr_;
00302 
00303         const unallocLabelList& u = lduAddr().upperAddr();
00304 
00305         for (register label face=0; face<lower.size(); face++)
00306         {
00307             lower[face] *= sf[u[face]];
00308         }
00309     }
00310 }
00311 
00312 
00313 void Foam::lduMatrix::operator*=(scalar s)
00314 {
00315     if (diagPtr_)
00316     {
00317         *diagPtr_ *= s;
00318     }
00319 
00320     if (upperPtr_)
00321     {
00322         *upperPtr_ *= s;
00323     }
00324 
00325     if (lowerPtr_)
00326     {
00327         *lowerPtr_ *= s;
00328     }
00329 }
00330 
00331 
00332 Foam::tmp<Foam::scalarField > Foam::lduMatrix::H1() const
00333 {
00334     tmp<scalarField > tH1
00335     (
00336         new scalarField(lduAddr().size(), 0.0)
00337     );
00338 
00339     if (lowerPtr_ || upperPtr_)
00340     {
00341         scalarField& H1_ = tH1();
00342 
00343         scalar* __restrict__ H1Ptr = H1_.begin();
00344 
00345         const label* __restrict__ uPtr = lduAddr().upperAddr().begin();
00346         const label* __restrict__ lPtr = lduAddr().lowerAddr().begin();
00347 
00348         const scalar* __restrict__ lowerPtr = lower().begin();
00349         const scalar* __restrict__ upperPtr = upper().begin();
00350 
00351         register const label nFaces = upper().size();
00352 
00353         for (register label face=0; face<nFaces; face++)
00354         {
00355             H1Ptr[uPtr[face]] -= lowerPtr[face];
00356             H1Ptr[lPtr[face]] -= upperPtr[face];
00357         }
00358     }
00359 
00360     return tH1;
00361 }
00362 
00363 
00364 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines