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

scalarMatricesTemplates.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 "scalarMatrices.H"
00027 #include <OpenFOAM/Swap.H>
00028 
00029 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00030 
00031 template<class Type>
00032 void Foam::solve
00033 (
00034     scalarSquareMatrix& tmpMatrix,
00035     Field<Type>& sourceSol
00036 )
00037 {
00038     label n = tmpMatrix.n();
00039 
00040     // Elimination
00041     for (register label i=0; i<n; i++)
00042     {
00043         label iMax = i;
00044         scalar largestCoeff = mag(tmpMatrix[iMax][i]);
00045 
00046         // Swap entries around to find a good pivot
00047         for (register label j=i+1; j<n; j++)
00048         {
00049             if (mag(tmpMatrix[j][i]) > largestCoeff)
00050             {
00051                 iMax = j;
00052                 largestCoeff = mag(tmpMatrix[iMax][i]);
00053             }
00054         }
00055 
00056         if (i != iMax)
00057         {
00058             //Info<< "Pivoted on " << i << " " << iMax << endl;
00059 
00060             for (register label k=i; k<n; k++)
00061             {
00062                 Swap(tmpMatrix[i][k], tmpMatrix[iMax][k]);
00063             }
00064             Swap(sourceSol[i], sourceSol[iMax]);
00065         }
00066 
00067         // Check that the system of equations isn't singular
00068         if (mag(tmpMatrix[i][i]) < 1e-20)
00069         {
00070             FatalErrorIn("solve(scalarSquareMatrix&, Field<Type>& sourceSol)")
00071                 << "Singular Matrix"
00072                 << exit(FatalError);
00073         }
00074 
00075         // Reduce to upper triangular form
00076         for (register label j=i+1; j<n; j++)
00077         {
00078             sourceSol[j] -= sourceSol[i]*(tmpMatrix[j][i]/tmpMatrix[i][i]);
00079 
00080             for (register label k=n-1; k>=i; k--)
00081             {
00082                 tmpMatrix[j][k] -=
00083                     tmpMatrix[i][k]*tmpMatrix[j][i]/tmpMatrix[i][i];
00084             }
00085         }
00086     }
00087 
00088     // Back-substitution
00089     for (register label j=n-1; j>=0; j--)
00090     {
00091         Type ntempvec = pTraits<Type>::zero;
00092 
00093         for (register label k=j+1; k<n; k++)
00094         {
00095             ntempvec += tmpMatrix[j][k]*sourceSol[k];
00096         }
00097 
00098         sourceSol[j] = (sourceSol[j] - ntempvec)/tmpMatrix[j][j];
00099     }
00100 }
00101 
00102 
00103 template<class Type>
00104 void Foam::solve
00105 (
00106     Field<Type>& psi,
00107     const scalarSquareMatrix& matrix,
00108     const Field<Type>& source
00109 )
00110 {
00111     scalarSquareMatrix tmpMatrix = matrix;
00112     psi = source;
00113     solve(tmpMatrix, psi);
00114 }
00115 
00116 
00117 template<class Type>
00118 void Foam::LUBacksubstitute
00119 (
00120     const scalarSquareMatrix& luMatrix,
00121     const labelList& pivotIndices,
00122     Field<Type>& sourceSol
00123 )
00124 {
00125     label n = luMatrix.n();
00126 
00127     label ii = 0;
00128 
00129     for (register label i=0; i<n; i++)
00130     {
00131         label ip = pivotIndices[i];
00132         Type sum = sourceSol[ip];
00133         sourceSol[ip] = sourceSol[i];
00134         const scalar* __restrict__ luMatrixi = luMatrix[i];
00135 
00136         if (ii != 0)
00137         {
00138             for (label j=ii-1; j<i; j++)
00139             {
00140                 sum -= luMatrixi[j]*sourceSol[j];
00141             }
00142         }
00143         else if (sum != pTraits<Type>::zero)
00144         {
00145             ii = i+1;
00146         }
00147 
00148         sourceSol[i] = sum;
00149     }
00150 
00151     for (register label i=n-1; i>=0; i--)
00152     {
00153         Type sum = sourceSol[i];
00154         const scalar* __restrict__ luMatrixi = luMatrix[i];
00155 
00156         for (register label j=i+1; j<n; j++)
00157         {
00158             sum -= luMatrixi[j]*sourceSol[j];
00159         }
00160 
00161         sourceSol[i] = sum/luMatrixi[i];
00162     }
00163 }
00164 
00165 
00166 template<class Type>
00167 void Foam::LUsolve
00168 (
00169     scalarSquareMatrix& matrix,
00170     Field<Type>& sourceSol
00171 )
00172 {
00173     labelList pivotIndices(matrix.n());
00174     LUDecompose(matrix, pivotIndices);
00175     LUBacksubstitute(matrix, pivotIndices, sourceSol);
00176 }
00177 
00178 
00179 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines