diff options
author | 2012-07-27 11:36:58 +0200 | |
---|---|---|
committer | 2012-07-27 11:36:58 +0200 | |
commit | c0fa5811ec233a5a3065cce78b1bca155a9b4fc8 (patch) | |
tree | cdad1af9ac4060906676c39ff303d71bdf0b7a71 /Eigen | |
parent | 925ace196c182759026d3eb3edc06565ab5f01ee (diff) |
Refactoring codes for numeric updates
Diffstat (limited to 'Eigen')
-rw-r--r-- | Eigen/src/SparseLU/SparseLU.h | 1 | ||||
-rw-r--r-- | Eigen/src/SparseLU/SparseLU_column_bmod.h | 43 | ||||
-rw-r--r-- | Eigen/src/SparseLU/SparseLU_kernel_bmod.h | 92 | ||||
-rw-r--r-- | Eigen/src/SparseLU/SparseLU_panel_bmod.h | 51 |
4 files changed, 98 insertions, 89 deletions
diff --git a/Eigen/src/SparseLU/SparseLU.h b/Eigen/src/SparseLU/SparseLU.h index 25fad0f29..474dfdedc 100644 --- a/Eigen/src/SparseLU/SparseLU.h +++ b/Eigen/src/SparseLU/SparseLU.h @@ -388,6 +388,7 @@ void SparseLU<MatrixType, OrderingType>::analyzePattern(const MatrixType& mat) #include "SparseLU_snode_bmod.h" #include "SparseLU_pivotL.h" #include "SparseLU_panel_dfs.h" +#include "SparseLU_kernel_bmod.h" #include "SparseLU_panel_bmod.h" #include "SparseLU_column_dfs.h" #include "SparseLU_column_bmod.h" diff --git a/Eigen/src/SparseLU/SparseLU_column_bmod.h b/Eigen/src/SparseLU/SparseLU_column_bmod.h index 00787721b..1457b6f35 100644 --- a/Eigen/src/SparseLU/SparseLU_column_bmod.h +++ b/Eigen/src/SparseLU/SparseLU_column_bmod.h @@ -66,7 +66,7 @@ int LU_column_bmod(const int jcol, const int nseg, BlockScalarVector& dense, Sca typedef typename IndexVector::Scalar Index; typedef typename ScalarVector::Scalar Scalar; int jsupno, k, ksub, krep, ksupno; - int lptr, nrow, isub, i, irow, nextlu, new_next, ufirst; + int lptr, nrow, isub, irow, nextlu, new_next, ufirst; int fsupc, nsupc, nsupr, luptr, kfnz, no_zeros; /* krep = representative of current k-th supernode * fsupc = first supernodal column @@ -122,46 +122,7 @@ int LU_column_bmod(const int jcol, const int nseg, BlockScalarVector& dense, Sca // Perform a triangular solver and block update, // then scatter the result of sup-col update to dense no_zeros = kfnz - fst_col; - // First, copy U[*,j] segment from dense(*) to tempv(*) - isub = lptr + no_zeros; - for (i = 0; i < segsize; i++) - { - irow = lsub(isub); - tempv(i) = dense(irow); - ++isub; - } - // Dense triangular solve -- start effective triangle - luptr += nsupr * no_zeros + no_zeros; - // Form Eigen matrix and vector - Map<Matrix<Scalar,Dynamic,Dynamic>, 0, OuterStride<> > A( &(lusup.data()[luptr]), segsize, segsize, OuterStride<>(nsupr) ); - VectorBlock<ScalarVector> u(tempv, 0, segsize); - - u = A.template triangularView<UnitLower>().solve(u); - - // Dense matrix-vector product y <-- A*x - luptr += segsize; - new (&A) Map<Matrix<Scalar,Dynamic, Dynamic>, 0, OuterStride<> > ( &(lusup.data()[luptr]), nrow, segsize, OuterStride<>(nsupr) ); - VectorBlock<ScalarVector> l(tempv, segsize, nrow); - l= A * u; - - // Scatter tempv[] into SPA dense[] as a temporary storage - isub = lptr + no_zeros; - for (i = 0; i < segsize; i++) - { - irow = lsub(isub); - dense(irow) = tempv(i); - tempv(i) = Scalar(0.0); - ++isub; - } - - // Scatter l into SPA dense[] - for (i = 0; i < nrow; i++) - { - irow = lsub(isub); - dense(irow) -= l(i); - l(i) = Scalar(0.0); - ++isub; - } + LU_kernel_bmod(segsize, dense, tempv, lusup, luptr, nsupr, nrow, lsub, lptr, no_zeros); } // end if jsupno } // end for each segment diff --git a/Eigen/src/SparseLU/SparseLU_kernel_bmod.h b/Eigen/src/SparseLU/SparseLU_kernel_bmod.h new file mode 100644 index 000000000..d5df70fd2 --- /dev/null +++ b/Eigen/src/SparseLU/SparseLU_kernel_bmod.h @@ -0,0 +1,92 @@ +// This file is part of Eigen, a lightweight C++ template library +// for linear algebra. +// +// Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr> +// +// Eigen is free software; you can redistribute it and/or +// modify it under the terms of the GNU Lesser General Public +// License as published by the Free Software Foundation; either +// version 3 of the License, or (at your option) any later version. +// +// Alternatively, you can redistribute it and/or +// modify it under the terms of the GNU General Public License as +// published by the Free Software Foundation; either version 2 of +// the License, or (at your option) any later version. +// +// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY +// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS +// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the +// GNU General Public License for more details. +// +// You should have received a copy of the GNU Lesser General Public +// License and a copy of the GNU General Public License along with +// Eigen. If not, see <http://www.gnu.org/licenses/>. + +#ifndef SPARSELU_KERNEL_BMOD_H +#define SPARSELU_KERNEL_BMOD_H + +/** + * \brief Performs numeric block updates from a given supernode to a single column + * + * \param segsize Size of the segment (and blocks ) to use for updates + * \param [in,out]dense Packed values of the original matrix + * \param tempv temporary vector to use for updates + * \param lusup array containing the supernodes + * \param nsupr Number of rows in the supernode + * \param nrow Number of rows in the rectangular part of the supernode + * \param lsub compressed row subscripts of supernodes + * \param lptr pointer to the first column of the current supernode in lsub + * \param no_zeros Number of nonzeros elements before the diagonal part of the supernode + * \return 0 on success + */ +template <typename BlockScalarVector, typename ScalarVector, typename IndexVector> +int LU_kernel_bmod(const int segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, int& luptr, const int nsupr, const int nrow, IndexVector& lsub, const int lptr, const int no_zeros) +{ + typedef typename ScalarVector::Scalar Scalar; + // First, copy U[*,j] segment from dense(*) to tempv(*) + // The result of triangular solve is in tempv[*]; + // The result of matric-vector update is in dense[*] + int isub = lptr + no_zeros; + int i, irow; + for (i = 0; i < segsize; i++) + { + irow = lsub(isub); + tempv(i) = dense(irow); + ++isub; + } + // Dense triangular solve -- start effective triangle + luptr += nsupr * no_zeros + no_zeros; + // Form Eigen matrix and vector + Map<Matrix<Scalar,Dynamic,Dynamic>, 0, OuterStride<> > A( &(lusup.data()[luptr]), segsize, segsize, OuterStride<>(nsupr) ); + VectorBlock<ScalarVector> u(tempv, 0, segsize); + + u = A.template triangularView<UnitLower>().solve(u); + + // Dense matrix-vector product y <-- A*x + luptr += segsize; + new (&A) Map<Matrix<Scalar,Dynamic, Dynamic>, 0, OuterStride<> > ( &(lusup.data()[luptr]), nrow, segsize, OuterStride<>(nsupr) ); + VectorBlock<ScalarVector> l(tempv, segsize, nrow); + l= A * u; + + // Scatter tempv[] into SPA dense[] as a temporary storage + isub = lptr + no_zeros; + for (i = 0; i < segsize; i++) + { + irow = lsub(isub); + dense(irow) = tempv(i); + tempv(i) = Scalar(0.0); + ++isub; + } + + // Scatter l into SPA dense[] + for (i = 0; i < nrow; i++) + { + irow = lsub(isub); + dense(irow) -= l(i); + l(i) = Scalar(0.0); + ++isub; + } + + return 0; +} +#endif
\ No newline at end of file diff --git a/Eigen/src/SparseLU/SparseLU_panel_bmod.h b/Eigen/src/SparseLU/SparseLU_panel_bmod.h index 59ec69ec8..ebff787ee 100644 --- a/Eigen/src/SparseLU/SparseLU_panel_bmod.h +++ b/Eigen/src/SparseLU/SparseLU_panel_bmod.h @@ -73,12 +73,12 @@ void LU_panel_bmod(const int m, const int w, const int jcol, const int nseg, Sca IndexVector& xlusup = glu.xlusup; ScalarVector& lusup = glu.lusup; - int i,ksub,jj,nextl_col,irow; + int ksub,jj,nextl_col; int fsupc, nsupc, nsupr, nrow; int krep, kfnz; int lptr; // points to the row subscripts of a supernode int luptr; // ... - int segsize,no_zeros,isub ; + int segsize,no_zeros ; // For each nonz supernode segment of U[*,j] in topological order int k = nseg - 1; for (ksub = 0; ksub < nseg; ksub++) @@ -118,52 +118,7 @@ void LU_panel_bmod(const int m, const int w, const int jcol, const int nseg, Sca // Perform a trianglar solve and block update, // then scatter the result of sup-col update to dense[] no_zeros = kfnz - fsupc; - // First Copy U[*,j] segment from dense[*] to tempv[*] : - // The result of triangular solve is in tempv[*]; - // The result of matric-vector update is in dense_col[*] - isub = lptr + no_zeros; - for (i = 0; i < segsize; ++i) - { - irow = lsub(isub); - tempv(i) = dense_col(irow); // Gather to a compact vector - ++isub; - } - // Start effective triangle - luptr += nsupr * no_zeros + no_zeros; - // triangular solve with Eigen - Map<Matrix<Scalar,Dynamic, Dynamic>, 0, OuterStride<> > A( &(lusup.data()[luptr]), segsize, segsize, OuterStride<>(nsupr) ); - VectorBlock<ScalarVector> u(tempv, 0, segsize); - u = A.template triangularView<UnitLower>().solve(u); - - luptr += segsize; - // Dense Matrix vector product y <-- A*x; - new (&A) Map<Matrix<Scalar,Dynamic, Dynamic>, 0, OuterStride<> > ( &(lusup.data()[luptr]), nrow, segsize, OuterStride<>(nsupr) ); - VectorBlock<ScalarVector> l(tempv, segsize, nrow); - l= A * u; - - // Scatter tempv(*) into SPA dense(*) such that tempv(*) - // can be used for the triangular solve of the next - // column of the panel. The y will be copied into ucol(*) - // after the whole panel has been finished... after column_dfs() and column_bmod() - - isub = lptr + no_zeros; - for (i = 0; i < segsize; i++) - { - irow = lsub(isub); - dense_col(irow) = tempv(i); - tempv(i) = Scalar(0.0); - isub++; - } - - // Scatter the update from &tempv[segsize] into SPA dense(*) - // Start dense rectangular L - for (i = 0; i < nrow; i++) - { - irow = lsub(isub); - dense_col(irow) -= l(i); - l(i) = Scalar(0); - ++isub; - } + LU_kernel_bmod(segsize, dense_col, tempv, lusup, luptr, nsupr, nrow, lsub, lptr, no_zeros); } // End for each column in the panel } // End for each updating supernode |