SparseLU_panel_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 [s,d,c,z]panel_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_PANEL_BMOD_H
32 #define SPARSELU_PANEL_BMOD_H
33 
34 // IWYU pragma: private
35 #include "./InternalHeaderCheck.h"
36 
37 namespace Eigen {
38 namespace internal {
39 
58 template <typename Scalar, typename StorageIndex>
59 void SparseLUImpl<Scalar, StorageIndex>::panel_bmod(const Index m, const Index w, const Index jcol, const Index nseg,
60  ScalarVector& dense, ScalarVector& tempv, IndexVector& segrep,
61  IndexVector& repfnz, GlobalLU_t& glu) {
62  Index ksub, jj, nextl_col;
63  Index fsupc, nsupc, nsupr, nrow;
64  Index krep, kfnz;
65  Index lptr; // points to the row subscripts of a supernode
66  Index luptr; // ...
67  Index segsize, no_zeros;
68  // For each nonz supernode segment of U[*,j] in topological order
69  Index k = nseg - 1;
71 
72  for (ksub = 0; ksub < nseg; ksub++) { // For each updating supernode
73  /* krep = representative of current k-th supernode
74  * fsupc = first supernodal column
75  * nsupc = number of columns in a supernode
76  * nsupr = number of rows in a supernode
77  */
78  krep = segrep(k);
79  k--;
80  fsupc = glu.xsup(glu.supno(krep));
81  nsupc = krep - fsupc + 1;
82  nsupr = glu.xlsub(fsupc + 1) - glu.xlsub(fsupc);
83  nrow = nsupr - nsupc;
84  lptr = glu.xlsub(fsupc);
85 
86  // loop over the panel columns to detect the actual number of columns and rows
87  Index u_rows = 0;
88  Index u_cols = 0;
89  for (jj = jcol; jj < jcol + w; jj++) {
90  nextl_col = (jj - jcol) * m;
91  VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
92 
93  kfnz = repfnz_col(krep);
94  if (kfnz == emptyIdxLU) continue; // skip any zero segment
95 
96  segsize = krep - kfnz + 1;
97  u_cols++;
98  u_rows = (std::max)(segsize, u_rows);
99  }
100 
101  if (nsupc >= 2) {
102  Index ldu = internal::first_multiple<Index>(u_rows, PacketSize);
103  Map<ScalarMatrix, Aligned, OuterStride<> > U(tempv.data(), u_rows, u_cols, OuterStride<>(ldu));
104 
105  // gather U
106  Index u_col = 0;
107  for (jj = jcol; jj < jcol + w; jj++) {
108  nextl_col = (jj - jcol) * m;
109  VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
110  VectorBlock<ScalarVector> dense_col(dense, nextl_col, m); // Scatter/gather entire matrix column from/to here
111 
112  kfnz = repfnz_col(krep);
113  if (kfnz == emptyIdxLU) continue; // skip any zero segment
114 
115  segsize = krep - kfnz + 1;
116  luptr = glu.xlusup(fsupc);
117  no_zeros = kfnz - fsupc;
118 
119  Index isub = lptr + no_zeros;
120  Index off = u_rows - segsize;
121  for (Index i = 0; i < off; i++) U(i, u_col) = 0;
122  for (Index i = 0; i < segsize; i++) {
123  Index irow = glu.lsub(isub);
124  U(i + off, u_col) = dense_col(irow);
125  ++isub;
126  }
127  u_col++;
128  }
129  // solve U = A^-1 U
130  luptr = glu.xlusup(fsupc);
131  Index lda = glu.xlusup(fsupc + 1) - glu.xlusup(fsupc);
132  no_zeros = (krep - u_rows + 1) - fsupc;
133  luptr += lda * no_zeros + no_zeros;
134  MappedMatrixBlock A(glu.lusup.data() + luptr, u_rows, u_rows, OuterStride<>(lda));
135  U = A.template triangularView<UnitLower>().solve(U);
136 
137  // update
138  luptr += u_rows;
139  MappedMatrixBlock B(glu.lusup.data() + luptr, nrow, u_rows, OuterStride<>(lda));
140  eigen_assert(tempv.size() > w * ldu + nrow * w + 1);
141 
142  Index ldl = internal::first_multiple<Index>(nrow, PacketSize);
143  Index offset = (PacketSize - internal::first_default_aligned(B.data(), PacketSize)) % PacketSize;
144  MappedMatrixBlock L(tempv.data() + w * ldu + offset, nrow, u_cols, OuterStride<>(ldl));
145 
146  L.noalias() = B * U;
147 
148  // scatter U and L
149  u_col = 0;
150  for (jj = jcol; jj < jcol + w; jj++) {
151  nextl_col = (jj - jcol) * m;
152  VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
153  VectorBlock<ScalarVector> dense_col(dense, nextl_col, m); // Scatter/gather entire matrix column from/to here
154 
155  kfnz = repfnz_col(krep);
156  if (kfnz == emptyIdxLU) continue; // skip any zero segment
157 
158  segsize = krep - kfnz + 1;
159  no_zeros = kfnz - fsupc;
160  Index isub = lptr + no_zeros;
161 
162  Index off = u_rows - segsize;
163  for (Index i = 0; i < segsize; i++) {
164  Index irow = glu.lsub(isub++);
165  dense_col(irow) = U.coeff(i + off, u_col);
166  U.coeffRef(i + off, u_col) = 0;
167  }
168 
169  // Scatter l into SPA dense[]
170  for (Index i = 0; i < nrow; i++) {
171  Index irow = glu.lsub(isub++);
172  dense_col(irow) -= L.coeff(i, u_col);
173  L.coeffRef(i, u_col) = 0;
174  }
175  u_col++;
176  }
177  } else // level 2 only
178  {
179  // Sequence through each column in the panel
180  for (jj = jcol; jj < jcol + w; jj++) {
181  nextl_col = (jj - jcol) * m;
182  VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
183  VectorBlock<ScalarVector> dense_col(dense, nextl_col, m); // Scatter/gather entire matrix column from/to here
184 
185  kfnz = repfnz_col(krep);
186  if (kfnz == emptyIdxLU) continue; // skip any zero segment
187 
188  segsize = krep - kfnz + 1;
189  luptr = glu.xlusup(fsupc);
190 
191  Index lda = glu.xlusup(fsupc + 1) - glu.xlusup(fsupc); // nsupr
192 
193  // Perform a trianglar solve and block update,
194  // then scatter the result of sup-col update to dense[]
195  no_zeros = kfnz - fsupc;
196  if (segsize == 1)
197  LU_kernel_bmod<1>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
198  else if (segsize == 2)
199  LU_kernel_bmod<2>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
200  else if (segsize == 3)
201  LU_kernel_bmod<3>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
202  else
203  LU_kernel_bmod<Dynamic>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr,
204  no_zeros);
205  } // End for each column in the panel
206  }
207 
208  } // End for each updating supernode
209 } // end panel bmod
210 
211 } // end namespace internal
212 
213 } // end namespace Eigen
214 
215 #endif // SPARSELU_PANEL_BMOD_H
int i
Definition: BiCGSTAB_step_by_step.cpp:9
MatrixXd L
Definition: LLT_example.cpp:6
#define eigen_assert(x)
Definition: Macros.h:910
RowVector3d w
Definition: Matrix_resize_int.cpp:3
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
Matrix< SCALARB, Dynamic, Dynamic, opt_B > B
Definition: bench_gemm.cpp:48
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Convenience specialization of Stride to specify only an outer stride See class Map for some examples.
Definition: Stride.h:104
constexpr EIGEN_DEVICE_FUNC const Scalar * data() const
Definition: PlainObjectBase.h:273
Expression of a fixed-size or dynamic-size sub-vector.
Definition: VectorBlock.h:58
void panel_bmod(const Index m, const Index w, const Index jcol, const Index nseg, ScalarVector &dense, ScalarVector &tempv, IndexVector &segrep, IndexVector &repfnz, GlobalLU_t &glu)
Performs numeric block updates (sup-panel) in topological order.
Definition: SparseLU_panel_bmod.h:59
Definition: matrices.h:74
#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
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
@ emptyIdxLU
Definition: SparseLU_Memory.h:41
static Index first_default_aligned(const DenseBase< Derived > &m)
Definition: DenseCoeffsBase.h:539
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
double U
Swimming speed.
Definition: two_d_variable_diff_adapt.cc:53
Definition: Eigen_Colamd.h:49
Definition: SparseLU_Structs.h:80
IndexVector xsup
Definition: SparseLU_Structs.h:82
IndexVector xlusup
Definition: SparseLU_Structs.h:86
IndexVector supno
Definition: SparseLU_Structs.h:83
IndexVector lsub
Definition: SparseLU_Structs.h:85
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