Eigen  3.3.3
SparseLU_column_bmod.h
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 xcolumn_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_COLUMN_BMOD_H
00032 #define SPARSELU_COLUMN_BMOD_H
00033 
00034 namespace Eigen {
00035 
00036 namespace internal {
00052 template <typename Scalar, typename StorageIndex>
00053 Index SparseLUImpl<Scalar,StorageIndex>::column_bmod(const Index jcol, const Index nseg, BlockScalarVector dense, ScalarVector& tempv,
00054                                                      BlockIndexVector segrep, BlockIndexVector repfnz, Index fpanelc, GlobalLU_t& glu)
00055 {
00056   Index  jsupno, k, ksub, krep, ksupno; 
00057   Index lptr, nrow, isub, irow, nextlu, new_next, ufirst; 
00058   Index fsupc, nsupc, nsupr, luptr, kfnz, no_zeros; 
00059   /* krep = representative of current k-th supernode
00060     * fsupc =  first supernodal column
00061     * nsupc = number of columns in a supernode
00062     * nsupr = number of rows in a supernode
00063     * luptr = location of supernodal LU-block in storage
00064     * kfnz = first nonz in the k-th supernodal segment
00065     * no_zeros = no lf leading zeros in a supernodal U-segment
00066     */
00067   
00068   jsupno = glu.supno(jcol);
00069   // For each nonzero supernode segment of U[*,j] in topological order 
00070   k = nseg - 1; 
00071   Index d_fsupc; // distance between the first column of the current panel and the 
00072                // first column of the current snode
00073   Index fst_col; // First column within small LU update
00074   Index segsize; 
00075   for (ksub = 0; ksub < nseg; ksub++)
00076   {
00077     krep = segrep(k); k--; 
00078     ksupno = glu.supno(krep); 
00079     if (jsupno != ksupno )
00080     {
00081       // outside the rectangular supernode 
00082       fsupc = glu.xsup(ksupno); 
00083       fst_col = (std::max)(fsupc, fpanelc); 
00084       
00085       // Distance from the current supernode to the current panel; 
00086       // d_fsupc = 0 if fsupc > fpanelc
00087       d_fsupc = fst_col - fsupc; 
00088       
00089       luptr = glu.xlusup(fst_col) + d_fsupc; 
00090       lptr = glu.xlsub(fsupc) + d_fsupc; 
00091       
00092       kfnz = repfnz(krep); 
00093       kfnz = (std::max)(kfnz, fpanelc); 
00094       
00095       segsize = krep - kfnz + 1; 
00096       nsupc = krep - fst_col + 1; 
00097       nsupr = glu.xlsub(fsupc+1) - glu.xlsub(fsupc); 
00098       nrow = nsupr - d_fsupc - nsupc;
00099       Index lda = glu.xlusup(fst_col+1) - glu.xlusup(fst_col);
00100       
00101       
00102       // Perform a triangular solver and block update, 
00103       // then scatter the result of sup-col update to dense
00104       no_zeros = kfnz - fst_col; 
00105       if(segsize==1)
00106         LU_kernel_bmod<1>::run(segsize, dense, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
00107       else
00108         LU_kernel_bmod<Dynamic>::run(segsize, dense, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
00109     } // end if jsupno 
00110   } // end for each segment
00111   
00112   // Process the supernodal portion of  L\U[*,j]
00113   nextlu = glu.xlusup(jcol); 
00114   fsupc = glu.xsup(jsupno);
00115   
00116   // copy the SPA dense into L\U[*,j]
00117   Index mem; 
00118   new_next = nextlu + glu.xlsub(fsupc + 1) - glu.xlsub(fsupc); 
00119   Index offset = internal::first_multiple<Index>(new_next, internal::packet_traits<Scalar>::size) - new_next;
00120   if(offset)
00121     new_next += offset;
00122   while (new_next > glu.nzlumax )
00123   {
00124     mem = memXpand<ScalarVector>(glu.lusup, glu.nzlumax, nextlu, LUSUP, glu.num_expansions);  
00125     if (mem) return mem; 
00126   }
00127   
00128   for (isub = glu.xlsub(fsupc); isub < glu.xlsub(fsupc+1); isub++)
00129   {
00130     irow = glu.lsub(isub);
00131     glu.lusup(nextlu) = dense(irow);
00132     dense(irow) = Scalar(0.0); 
00133     ++nextlu; 
00134   }
00135   
00136   if(offset)
00137   {
00138     glu.lusup.segment(nextlu,offset).setZero();
00139     nextlu += offset;
00140   }
00141   glu.xlusup(jcol + 1) = StorageIndex(nextlu);  // close L\U(*,jcol); 
00142   
00143   /* For more updates within the panel (also within the current supernode),
00144    * should start from the first column of the panel, or the first column
00145    * of the supernode, whichever is bigger. There are two cases:
00146    *  1) fsupc < fpanelc, then fst_col <-- fpanelc
00147    *  2) fsupc >= fpanelc, then fst_col <-- fsupc
00148    */
00149   fst_col = (std::max)(fsupc, fpanelc); 
00150   
00151   if (fst_col  < jcol)
00152   {
00153     // Distance between the current supernode and the current panel
00154     // d_fsupc = 0 if fsupc >= fpanelc
00155     d_fsupc = fst_col - fsupc; 
00156     
00157     lptr = glu.xlsub(fsupc) + d_fsupc; 
00158     luptr = glu.xlusup(fst_col) + d_fsupc; 
00159     nsupr = glu.xlsub(fsupc+1) - glu.xlsub(fsupc); // leading dimension
00160     nsupc = jcol - fst_col; // excluding jcol 
00161     nrow = nsupr - d_fsupc - nsupc; 
00162     
00163     // points to the beginning of jcol in snode L\U(jsupno) 
00164     ufirst = glu.xlusup(jcol) + d_fsupc; 
00165     Index lda = glu.xlusup(jcol+1) - glu.xlusup(jcol);
00166     MappedMatrixBlock A( &(glu.lusup.data()[luptr]), nsupc, nsupc, OuterStride<>(lda) );
00167     VectorBlock<ScalarVector> u(glu.lusup, ufirst, nsupc); 
00168     u = A.template triangularView<UnitLower>().solve(u); 
00169     
00170     new (&A) MappedMatrixBlock ( &(glu.lusup.data()[luptr+nsupc]), nrow, nsupc, OuterStride<>(lda) );
00171     VectorBlock<ScalarVector> l(glu.lusup, ufirst+nsupc, nrow); 
00172     l.noalias() -= A * u;
00173     
00174   } // End if fst_col
00175   return 0; 
00176 }
00177 
00178 } // end namespace internal
00179 } // end namespace Eigen
00180 
00181 #endif // SPARSELU_COLUMN_BMOD_H
 All Classes Functions Variables Typedefs Enumerations Enumerator Friends