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

DICPreconditioner.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 <OpenFOAM/DICPreconditioner.H>
00027 
00028 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00029 
00030 namespace Foam
00031 {
00032     defineTypeNameAndDebug(DICPreconditioner, 0);
00033 
00034     lduMatrix::preconditioner::
00035         addsymMatrixConstructorToTable<DICPreconditioner>
00036         addDICPreconditionerSymMatrixConstructorToTable_;
00037 }
00038 
00039 
00040 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00041 
00042 Foam::DICPreconditioner::DICPreconditioner
00043 (
00044     const lduMatrix::solver& sol,
00045     const dictionary&
00046 )
00047 :
00048     lduMatrix::preconditioner(sol),
00049     rD_(sol.matrix().diag())
00050 {
00051     calcReciprocalD(rD_, sol.matrix());
00052 }
00053 
00054 
00055 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00056 
00057 void Foam::DICPreconditioner::calcReciprocalD
00058 (
00059     scalarField& rD,
00060     const lduMatrix& matrix
00061 )
00062 {
00063     scalar* __restrict__ rDPtr = rD.begin();
00064 
00065     const label* const __restrict__ uPtr = matrix.lduAddr().upperAddr().begin();
00066     const label* const __restrict__ lPtr = matrix.lduAddr().lowerAddr().begin();
00067     const scalar* const __restrict__ upperPtr = matrix.upper().begin();
00068 
00069     // Calculate the DIC diagonal
00070     register const label nFaces = matrix.upper().size();
00071     for (register label face=0; face<nFaces; face++)
00072     {
00073         rDPtr[uPtr[face]] -= upperPtr[face]*upperPtr[face]/rDPtr[lPtr[face]];
00074     }
00075 
00076 
00077     // Calculate the reciprocal of the preconditioned diagonal
00078     register const label nCells = rD.size();
00079 
00080     for (register label cell=0; cell<nCells; cell++)
00081     {
00082         rDPtr[cell] = 1.0/rDPtr[cell];
00083     }
00084 }
00085 
00086 
00087 void Foam::DICPreconditioner::precondition
00088 (
00089     scalarField& wA,
00090     const scalarField& rA,
00091     const direction
00092 ) const
00093 {
00094     scalar* __restrict__ wAPtr = wA.begin();
00095     const scalar* __restrict__ rAPtr = rA.begin();
00096     const scalar* __restrict__ rDPtr = rD_.begin();
00097 
00098     const label* const __restrict__ uPtr =
00099         solver_.matrix().lduAddr().upperAddr().begin();
00100     const label* const __restrict__ lPtr =
00101         solver_.matrix().lduAddr().lowerAddr().begin();
00102     const scalar* const __restrict__ upperPtr = solver_.matrix().upper().begin();
00103 
00104     register label nCells = wA.size();
00105     register label nFaces = solver_.matrix().upper().size();
00106     register label nFacesM1 = nFaces - 1;
00107 
00108     for (register label cell=0; cell<nCells; cell++)
00109     {
00110         wAPtr[cell] = rDPtr[cell]*rAPtr[cell];
00111     }
00112 
00113     for (register label face=0; face<nFaces; face++)
00114     {
00115         wAPtr[uPtr[face]] -= rDPtr[uPtr[face]]*upperPtr[face]*wAPtr[lPtr[face]];
00116     }
00117 
00118     for (register label face=nFacesM1; face>=0; face--)
00119     {
00120         wAPtr[lPtr[face]] -= rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[uPtr[face]];
00121     }
00122 }
00123 
00124 
00125 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines