SparseLU_panel_bmod.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
00005 // Copyright (C) 2012 Gael Guennebaud <gael.guennebaud@inria.fr>
00006 //
00007 // This Source Code Form is subject to the terms of the Mozilla
00008 // Public License v. 2.0. If a copy of the MPL was not distributed
00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00010 
00011 /* 
00012  
00013  * NOTE: This file is the modified version of [s,d,c,z]panel_bmod.c file in SuperLU 
00014  
00015  * -- SuperLU routine (version 3.0) --
00016  * Univ. of California Berkeley, Xerox Palo Alto Research Center,
00017  * and Lawrence Berkeley National Lab.
00018  * October 15, 2003
00019  *
00020  * Copyright (c) 1994 by Xerox Corporation.  All rights reserved.
00021  *
00022  * THIS MATERIAL IS PROVIDED AS IS, WITH ABSOLUTELY NO WARRANTY
00023  * EXPRESSED OR IMPLIED.  ANY USE IS AT YOUR OWN RISK.
00024  *
00025  * Permission is hereby granted to use or copy this program for any
00026  * purpose, provided the above notices are retained on all copies.
00027  * Permission to modify the code and to distribute modified code is
00028  * granted, provided the above notices are retained, and a notice that
00029  * the code was modified is included with the above copyright notice.
00030  */
00031 #ifndef SPARSELU_PANEL_BMOD_H
00032 #define SPARSELU_PANEL_BMOD_H
00033 
00034 namespace Eigen {
00035 namespace internal {
00036 
00055 template <typename Scalar, typename Index>
00056 void SparseLUImpl<Scalar,Index>::panel_bmod(const Index m, const Index w, const Index jcol, 
00057                                             const Index nseg, ScalarVector& dense, ScalarVector& tempv,
00058                                             IndexVector& segrep, IndexVector& repfnz, GlobalLU_t& glu)
00059 {
00060   
00061   Index ksub,jj,nextl_col; 
00062   Index fsupc, nsupc, nsupr, nrow; 
00063   Index krep, kfnz; 
00064   Index lptr; // points to the row subscripts of a supernode 
00065   Index luptr; // ...
00066   Index segsize,no_zeros ; 
00067   // For each nonz supernode segment of U[*,j] in topological order
00068   Index k = nseg - 1; 
00069   const Index PacketSize = internal::packet_traits<Scalar>::size;
00070   
00071   for (ksub = 0; ksub < nseg; ksub++)
00072   { // For each updating supernode
00073     /* krep = representative of current k-th supernode
00074      * fsupc =  first supernodal column
00075      * nsupc = number of columns in a supernode
00076      * nsupr = number of rows in a supernode
00077      */
00078     krep = segrep(k); k--; 
00079     fsupc = glu.xsup(glu.supno(krep)); 
00080     nsupc = krep - fsupc + 1; 
00081     nsupr = glu.xlsub(fsupc+1) - glu.xlsub(fsupc); 
00082     nrow = nsupr - nsupc; 
00083     lptr = glu.xlsub(fsupc); 
00084     
00085     // loop over the panel columns to detect the actual number of columns and rows
00086     Index u_rows = 0;
00087     Index u_cols = 0;
00088     for (jj = jcol; jj < jcol + w; jj++)
00089     {
00090       nextl_col = (jj-jcol) * m; 
00091       VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
00092       
00093       kfnz = repfnz_col(krep); 
00094       if ( kfnz == emptyIdxLU ) 
00095         continue; // skip any zero segment
00096       
00097       segsize = krep - kfnz + 1;
00098       u_cols++;
00099       u_rows = (std::max)(segsize,u_rows);
00100     }
00101     
00102     if(nsupc >= 2)
00103     { 
00104       Index ldu = internal::first_multiple<Index>(u_rows, PacketSize);
00105       Map<Matrix<Scalar,Dynamic,Dynamic>, Aligned,  OuterStride<> > U(tempv.data(), u_rows, u_cols, OuterStride<>(ldu));
00106       
00107       // gather U
00108       Index u_col = 0;
00109       for (jj = jcol; jj < jcol + w; jj++)
00110       {
00111         nextl_col = (jj-jcol) * m; 
00112         VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
00113         VectorBlock<ScalarVector> dense_col(dense, nextl_col, m); // Scatter/gather entire matrix column from/to here
00114         
00115         kfnz = repfnz_col(krep); 
00116         if ( kfnz == emptyIdxLU ) 
00117           continue; // skip any zero segment
00118         
00119         segsize = krep - kfnz + 1;
00120         luptr = glu.xlusup(fsupc);    
00121         no_zeros = kfnz - fsupc; 
00122         
00123         Index isub = lptr + no_zeros;
00124         Index off = u_rows-segsize;
00125         for (Index i = 0; i < off; i++) U(i,u_col) = 0;
00126         for (Index i = 0; i < segsize; i++)
00127         {
00128           Index irow = glu.lsub(isub); 
00129           U(i+off,u_col) = dense_col(irow); 
00130           ++isub; 
00131         }
00132         u_col++;
00133       }
00134       // solve U = A^-1 U
00135       luptr = glu.xlusup(fsupc);
00136       Index lda = glu.xlusup(fsupc+1) - glu.xlusup(fsupc);
00137       no_zeros = (krep - u_rows + 1) - fsupc;
00138       luptr += lda * no_zeros + no_zeros;
00139       Map<Matrix<Scalar,Dynamic,Dynamic>, 0, OuterStride<> > A(glu.lusup.data()+luptr, u_rows, u_rows, OuterStride<>(lda) );
00140       U = A.template triangularView<UnitLower>().solve(U);
00141       
00142       // update
00143       luptr += u_rows;
00144       Map<Matrix<Scalar,Dynamic,Dynamic>, 0, OuterStride<> > B(glu.lusup.data()+luptr, nrow, u_rows, OuterStride<>(lda) );
00145       eigen_assert(tempv.size()>w*ldu + nrow*w + 1);
00146       
00147       Index ldl = internal::first_multiple<Index>(nrow, PacketSize);
00148       Index offset = (PacketSize-internal::first_aligned(B.data(), PacketSize)) % PacketSize;
00149       Map<Matrix<Scalar,Dynamic,Dynamic>, 0, OuterStride<> > L(tempv.data()+w*ldu+offset, nrow, u_cols, OuterStride<>(ldl));
00150       
00151       L.setZero();
00152       internal::sparselu_gemm<Scalar>(L.rows(), L.cols(), B.cols(), B.data(), B.outerStride(), U.data(), U.outerStride(), L.data(), L.outerStride());
00153       
00154       // scatter U and L
00155       u_col = 0;
00156       for (jj = jcol; jj < jcol + w; jj++)
00157       {
00158         nextl_col = (jj-jcol) * m; 
00159         VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
00160         VectorBlock<ScalarVector> dense_col(dense, nextl_col, m); // Scatter/gather entire matrix column from/to here
00161         
00162         kfnz = repfnz_col(krep); 
00163         if ( kfnz == emptyIdxLU ) 
00164           continue; // skip any zero segment
00165         
00166         segsize = krep - kfnz + 1;
00167         no_zeros = kfnz - fsupc; 
00168         Index isub = lptr + no_zeros;
00169         
00170         Index off = u_rows-segsize;
00171         for (Index i = 0; i < segsize; i++)
00172         {
00173           Index irow = glu.lsub(isub++); 
00174           dense_col(irow) = U.coeff(i+off,u_col);
00175           U.coeffRef(i+off,u_col) = 0;
00176         }
00177         
00178         // Scatter l into SPA dense[]
00179         for (Index i = 0; i < nrow; i++)
00180         {
00181           Index irow = glu.lsub(isub++); 
00182           dense_col(irow) -= L.coeff(i,u_col);
00183           L.coeffRef(i,u_col) = 0;
00184         }
00185         u_col++;
00186       }
00187     }
00188     else // level 2 only
00189     {
00190       // Sequence through each column in the panel
00191       for (jj = jcol; jj < jcol + w; jj++)
00192       {
00193         nextl_col = (jj-jcol) * m; 
00194         VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
00195         VectorBlock<ScalarVector> dense_col(dense, nextl_col, m); // Scatter/gather entire matrix column from/to here
00196         
00197         kfnz = repfnz_col(krep); 
00198         if ( kfnz == emptyIdxLU ) 
00199           continue; // skip any zero segment
00200         
00201         segsize = krep - kfnz + 1;
00202         luptr = glu.xlusup(fsupc);
00203         
00204         Index lda = glu.xlusup(fsupc+1)-glu.xlusup(fsupc);// nsupr
00205         
00206         // Perform a trianglar solve and block update, 
00207         // then scatter the result of sup-col update to dense[]
00208         no_zeros = kfnz - fsupc; 
00209               if(segsize==1)  LU_kernel_bmod<1>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
00210         else  if(segsize==2)  LU_kernel_bmod<2>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
00211         else  if(segsize==3)  LU_kernel_bmod<3>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
00212         else                  LU_kernel_bmod<Dynamic>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros); 
00213       } // End for each column in the panel 
00214     }
00215     
00216   } // End for each updating supernode
00217 } // end panel bmod
00218 
00219 } // end namespace internal
00220 
00221 } // end namespace Eigen
00222 
00223 #endif // SPARSELU_PANEL_BMOD_H


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Thu Aug 27 2015 12:00:34