SparseLU_column_bmod.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
5 // Copyright (C) 2012 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 /*
12 
13  * NOTE: This file is the modified version of xcolumn_bmod.c file in SuperLU
14 
15  * -- SuperLU routine (version 3.0) --
16  * Univ. of California Berkeley, Xerox Palo Alto Research Center,
17  * and Lawrence Berkeley National Lab.
18  * October 15, 2003
19  *
20  * Copyright (c) 1994 by Xerox Corporation. All rights reserved.
21  *
22  * THIS MATERIAL IS PROVIDED AS IS, WITH ABSOLUTELY NO WARRANTY
23  * EXPRESSED OR IMPLIED. ANY USE IS AT YOUR OWN RISK.
24  *
25  * Permission is hereby granted to use or copy this program for any
26  * purpose, provided the above notices are retained on all copies.
27  * Permission to modify the code and to distribute modified code is
28  * granted, provided the above notices are retained, and a notice that
29  * the code was modified is included with the above copyright notice.
30  */
31 #ifndef SPARSELU_COLUMN_BMOD_H
32 #define SPARSELU_COLUMN_BMOD_H
33 
34 // IWYU pragma: private
35 #include "./InternalHeaderCheck.h"
36 
37 namespace Eigen {
38 
39 namespace internal {
55 template <typename Scalar, typename StorageIndex>
57  ScalarVector& tempv, BlockIndexVector segrep,
58  BlockIndexVector repfnz, Index fpanelc, GlobalLU_t& glu) {
59  Index jsupno, k, ksub, krep, ksupno;
60  Index lptr, nrow, isub, irow, nextlu, new_next, ufirst;
61  Index fsupc, nsupc, nsupr, luptr, kfnz, no_zeros;
62  /* krep = representative of current k-th supernode
63  * fsupc = first supernodal column
64  * nsupc = number of columns in a supernode
65  * nsupr = number of rows in a supernode
66  * luptr = location of supernodal LU-block in storage
67  * kfnz = first nonz in the k-th supernodal segment
68  * no_zeros = no lf leading zeros in a supernodal U-segment
69  */
70 
71  jsupno = glu.supno(jcol);
72  // For each nonzero supernode segment of U[*,j] in topological order
73  k = nseg - 1;
74  Index d_fsupc; // distance between the first column of the current panel and the
75  // first column of the current snode
76  Index fst_col; // First column within small LU update
77  Index segsize;
78  for (ksub = 0; ksub < nseg; ksub++) {
79  krep = segrep(k);
80  k--;
81  ksupno = glu.supno(krep);
82  if (jsupno != ksupno) {
83  // outside the rectangular supernode
84  fsupc = glu.xsup(ksupno);
85  fst_col = (std::max)(fsupc, fpanelc);
86 
87  // Distance from the current supernode to the current panel;
88  // d_fsupc = 0 if fsupc > fpanelc
89  d_fsupc = fst_col - fsupc;
90 
91  luptr = glu.xlusup(fst_col) + d_fsupc;
92  lptr = glu.xlsub(fsupc) + d_fsupc;
93 
94  kfnz = repfnz(krep);
95  kfnz = (std::max)(kfnz, fpanelc);
96 
97  segsize = krep - kfnz + 1;
98  nsupc = krep - fst_col + 1;
99  nsupr = glu.xlsub(fsupc + 1) - glu.xlsub(fsupc);
100  nrow = nsupr - d_fsupc - nsupc;
101  Index lda = glu.xlusup(fst_col + 1) - glu.xlusup(fst_col);
102 
103  // Perform a triangular solver and block update,
104  // then scatter the result of sup-col update to dense
105  no_zeros = kfnz - fst_col;
106  if (segsize == 1)
107  LU_kernel_bmod<1>::run(segsize, dense, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
108  else
109  LU_kernel_bmod<Dynamic>::run(segsize, dense, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
110  } // end if jsupno
111  } // end for each segment
112 
113  // Process the supernodal portion of L\U[*,j]
114  nextlu = glu.xlusup(jcol);
115  fsupc = glu.xsup(jsupno);
116 
117  // copy the SPA dense into L\U[*,j]
118  Index mem;
119  new_next = nextlu + glu.xlsub(fsupc + 1) - glu.xlsub(fsupc);
120  Index offset = internal::first_multiple<Index>(new_next, internal::packet_traits<Scalar>::size) - new_next;
121  if (offset) new_next += offset;
122  while (new_next > glu.nzlumax) {
123  mem = memXpand<ScalarVector>(glu.lusup, glu.nzlumax, nextlu, LUSUP, glu.num_expansions);
124  if (mem) return mem;
125  }
126 
127  for (isub = glu.xlsub(fsupc); isub < glu.xlsub(fsupc + 1); isub++) {
128  irow = glu.lsub(isub);
129  glu.lusup(nextlu) = dense(irow);
130  dense(irow) = Scalar(0.0);
131  ++nextlu;
132  }
133 
134  if (offset) {
135  glu.lusup.segment(nextlu, offset).setZero();
136  nextlu += offset;
137  }
138  glu.xlusup(jcol + 1) = StorageIndex(nextlu); // close L\U(*,jcol);
139 
140  /* For more updates within the panel (also within the current supernode),
141  * should start from the first column of the panel, or the first column
142  * of the supernode, whichever is bigger. There are two cases:
143  * 1) fsupc < fpanelc, then fst_col <-- fpanelc
144  * 2) fsupc >= fpanelc, then fst_col <-- fsupc
145  */
146  fst_col = (std::max)(fsupc, fpanelc);
147 
148  if (fst_col < jcol) {
149  // Distance between the current supernode and the current panel
150  // d_fsupc = 0 if fsupc >= fpanelc
151  d_fsupc = fst_col - fsupc;
152 
153  lptr = glu.xlsub(fsupc) + d_fsupc;
154  luptr = glu.xlusup(fst_col) + d_fsupc;
155  nsupr = glu.xlsub(fsupc + 1) - glu.xlsub(fsupc); // leading dimension
156  nsupc = jcol - fst_col; // excluding jcol
157  nrow = nsupr - d_fsupc - nsupc;
158 
159  // points to the beginning of jcol in snode L\U(jsupno)
160  ufirst = glu.xlusup(jcol) + d_fsupc;
161  Index lda = glu.xlusup(jcol + 1) - glu.xlusup(jcol);
162  MappedMatrixBlock A(&(glu.lusup.data()[luptr]), nsupc, nsupc, OuterStride<>(lda));
163  VectorBlock<ScalarVector> u(glu.lusup, ufirst, nsupc);
164  u = A.template triangularView<UnitLower>().solve(u);
165 
166  new (&A) MappedMatrixBlock(&(glu.lusup.data()[luptr + nsupc]), nrow, nsupc, OuterStride<>(lda));
167  VectorBlock<ScalarVector> l(glu.lusup, ufirst + nsupc, nrow);
168  l.noalias() -= A * u;
169 
170  } // End if fst_col
171  return 0;
172 }
173 
174 } // end namespace internal
175 } // end namespace Eigen
176 
177 #endif // SPARSELU_COLUMN_BMOD_H
SCALAR Scalar
Definition: bench_gemm.cpp:45
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:264
Index column_bmod(const Index jcol, const Index nseg, BlockScalarVector dense, ScalarVector &tempv, BlockIndexVector segrep, BlockIndexVector repfnz, Index fpanelc, GlobalLU_t &glu)
Performs numeric block updates (sup-col) in topological order.
Definition: SparseLU_column_bmod.h:56
#define max(a, b)
Definition: datatypes.h:23
const char const int const RealScalar const RealScalar const int * lda
Definition: level2_cplx_impl.h:20
char char char int int * k
Definition: level2_impl.h:374
@ LUSUP
Definition: SparseLU_Structs.h:77
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
Definition: Eigen_Colamd.h:49
Definition: SparseLU_Structs.h:80
IndexVector xsup
Definition: SparseLU_Structs.h:82
Index num_expansions
Definition: SparseLU_Structs.h:95
IndexVector xlusup
Definition: SparseLU_Structs.h:86
IndexVector supno
Definition: SparseLU_Structs.h:83
IndexVector lsub
Definition: SparseLU_Structs.h:85
Index nzlumax
Definition: SparseLU_Structs.h:89
IndexVector xlsub
Definition: SparseLU_Structs.h:87
ScalarVector lusup
Definition: SparseLU_Structs.h:84
static EIGEN_DONT_INLINE void run(const Index segsize, BlockScalarVector &dense, ScalarVector &tempv, ScalarVector &lusup, Index &luptr, const Index lda, const Index nrow, IndexVector &lsub, const Index lptr, const Index no_zeros)
Performs numeric block updates from a given supernode to a single column.
Definition: SparseLU_kernel_bmod.h:43
Definition: GenericPacketMath.h:108