Medial Code Documentation
Loading...
Searching...
No Matches
SparseLU_panel_bmod.h
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
34namespace Eigen {
35namespace internal {
36
55template <typename Scalar, typename StorageIndex>
59{
60
62 Index fsupc, nsupc, nsupr, nrow;
63 Index krep, kfnz;
64 Index lptr; // points to the row subscripts of a supernode
65 Index luptr; // ...
67 // For each nonz supernode segment of U[*,j] in topological order
68 Index k = nseg - 1;
70
71 for (ksub = 0; ksub < nseg; ksub++)
72 { // 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); k--;
79 fsupc = glu.xsup(glu.supno(krep));
80 nsupc = krep - fsupc + 1;
81 nsupr = glu.xlsub(fsupc+1) - glu.xlsub(fsupc);
82 nrow = nsupr - nsupc;
83 lptr = glu.xlsub(fsupc);
84
85 // loop over the panel columns to detect the actual number of columns and rows
86 Index u_rows = 0;
87 Index u_cols = 0;
88 for (jj = jcol; jj < jcol + w; jj++)
89 {
90 nextl_col = (jj-jcol) * m;
91 VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
92
94 if ( kfnz == emptyIdxLU )
95 continue; // skip any zero segment
96
97 segsize = krep - kfnz + 1;
98 u_cols++;
99 u_rows = (std::max)(segsize,u_rows);
100 }
101
102 if(nsupc >= 2)
103 {
104 Index ldu = internal::first_multiple<Index>(u_rows, PacketSize);
106
107 // gather U
108 Index u_col = 0;
109 for (jj = jcol; jj < jcol + w; jj++)
110 {
111 nextl_col = (jj-jcol) * m;
112 VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
113 VectorBlock<ScalarVector> dense_col(dense, nextl_col, m); // Scatter/gather entire matrix column from/to here
114
115 kfnz = repfnz_col(krep);
116 if ( kfnz == emptyIdxLU )
117 continue; // skip any zero segment
118
119 segsize = krep - kfnz + 1;
120 luptr = glu.xlusup(fsupc);
121 no_zeros = kfnz - fsupc;
122
125 for (Index i = 0; i < off; i++) U(i,u_col) = 0;
126 for (Index i = 0; i < segsize; i++)
127 {
128 Index irow = glu.lsub(isub);
129 U(i+off,u_col) = dense_col(irow);
130 ++isub;
131 }
132 u_col++;
133 }
134 // solve U = A^-1 U
135 luptr = glu.xlusup(fsupc);
136 Index lda = glu.xlusup(fsupc+1) - glu.xlusup(fsupc);
137 no_zeros = (krep - u_rows + 1) - fsupc;
138 luptr += lda * no_zeros + no_zeros;
139 MappedMatrixBlock A(glu.lusup.data()+luptr, u_rows, u_rows, OuterStride<>(lda) );
140 U = A.template triangularView<UnitLower>().solve(U);
141
142 // update
143 luptr += u_rows;
144 MappedMatrixBlock B(glu.lusup.data()+luptr, nrow, u_rows, OuterStride<>(lda) );
145 eigen_assert(tempv.size()>w*ldu + nrow*w + 1);
146
147 Index ldl = internal::first_multiple<Index>(nrow, PacketSize);
148 Index offset = (PacketSize-internal::first_default_aligned(B.data(), PacketSize)) % PacketSize;
149 MappedMatrixBlock L(tempv.data()+w*ldu+offset, nrow, u_cols, OuterStride<>(ldl));
150
151 L.noalias() = B * U;
152
153 // scatter U and L
154 u_col = 0;
155 for (jj = jcol; jj < jcol + w; jj++)
156 {
157 nextl_col = (jj-jcol) * m;
158 VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
159 VectorBlock<ScalarVector> dense_col(dense, nextl_col, m); // Scatter/gather entire matrix column from/to here
160
161 kfnz = repfnz_col(krep);
162 if ( kfnz == emptyIdxLU )
163 continue; // skip any zero segment
164
165 segsize = krep - kfnz + 1;
166 no_zeros = kfnz - fsupc;
168
170 for (Index i = 0; i < segsize; i++)
171 {
172 Index irow = glu.lsub(isub++);
173 dense_col(irow) = U.coeff(i+off,u_col);
174 U.coeffRef(i+off,u_col) = 0;
175 }
176
177 // Scatter l into SPA dense[]
178 for (Index i = 0; i < nrow; i++)
179 {
180 Index irow = glu.lsub(isub++);
181 dense_col(irow) -= L.coeff(i,u_col);
182 L.coeffRef(i,u_col) = 0;
183 }
184 u_col++;
185 }
186 }
187 else // level 2 only
188 {
189 // Sequence through each column in the panel
190 for (jj = jcol; jj < jcol + w; jj++)
191 {
192 nextl_col = (jj-jcol) * m;
193 VectorBlock<IndexVector> repfnz_col(repfnz, nextl_col, m); // First nonzero column index for each row
194 VectorBlock<ScalarVector> dense_col(dense, nextl_col, m); // Scatter/gather entire matrix column from/to here
195
196 kfnz = repfnz_col(krep);
197 if ( kfnz == emptyIdxLU )
198 continue; // skip any zero segment
199
200 segsize = krep - kfnz + 1;
201 luptr = glu.xlusup(fsupc);
202
203 Index lda = glu.xlusup(fsupc+1)-glu.xlusup(fsupc);// nsupr
204
205 // Perform a trianglar solve and block update,
206 // then scatter the result of sup-col update to dense[]
207 no_zeros = kfnz - fsupc;
208 if(segsize==1) LU_kernel_bmod<1>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
209 else if(segsize==2) LU_kernel_bmod<2>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
210 else if(segsize==3) LU_kernel_bmod<3>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
211 else LU_kernel_bmod<Dynamic>::run(segsize, dense_col, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
212 } // End for each column in the panel
213 }
214
215 } // End for each updating supernode
216} // end panel bmod
217
218} // end namespace internal
219
220} // end namespace Eigen
221
222#endif // SPARSELU_PANEL_BMOD_H
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
NoAlias< Derived, Eigen::MatrixBase > EIGEN_DEVICE_FUNC noalias()
Definition NoAlias.h:102
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:56
Namespace containing all symbols from the Eigen library.
Definition LDLT.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
Definition SparseLU_Structs.h:77
Definition SparseLU_kernel_bmod.h:18
Definition GenericPacketMath.h:107