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 | namespace Eigen {
|
---|
35 |
|
---|
36 | namespace internal {
|
---|
37 | /**
|
---|
38 | * \brief Performs numeric block updates (sup-col) in topological order
|
---|
39 | *
|
---|
40 | * \param jcol current column to update
|
---|
41 | * \param nseg Number of segments in the U part
|
---|
42 | * \param dense Store the full representation of the column
|
---|
43 | * \param tempv working array
|
---|
44 | * \param segrep segment representative ...
|
---|
45 | * \param repfnz ??? First nonzero column in each row ??? ...
|
---|
46 | * \param fpanelc First column in the current panel
|
---|
47 | * \param glu Global LU data.
|
---|
48 | * \return 0 - successful return
|
---|
49 | * > 0 - number of bytes allocated when run out of space
|
---|
50 | *
|
---|
51 | */
|
---|
52 | template <typename Scalar, typename Index>
|
---|
53 | Index SparseLUImpl<Scalar,Index>::column_bmod(const Index jcol, const Index nseg, BlockScalarVector dense, ScalarVector& tempv, BlockIndexVector segrep, BlockIndexVector repfnz, Index fpanelc, GlobalLU_t& glu)
|
---|
54 | {
|
---|
55 | Index jsupno, k, ksub, krep, ksupno;
|
---|
56 | Index lptr, nrow, isub, irow, nextlu, new_next, ufirst;
|
---|
57 | Index fsupc, nsupc, nsupr, luptr, kfnz, no_zeros;
|
---|
58 | /* krep = representative of current k-th supernode
|
---|
59 | * fsupc = first supernodal column
|
---|
60 | * nsupc = number of columns in a supernode
|
---|
61 | * nsupr = number of rows in a supernode
|
---|
62 | * luptr = location of supernodal LU-block in storage
|
---|
63 | * kfnz = first nonz in the k-th supernodal segment
|
---|
64 | * no_zeros = no lf leading zeros in a supernodal U-segment
|
---|
65 | */
|
---|
66 |
|
---|
67 | jsupno = glu.supno(jcol);
|
---|
68 | // For each nonzero supernode segment of U[*,j] in topological order
|
---|
69 | k = nseg - 1;
|
---|
70 | Index d_fsupc; // distance between the first column of the current panel and the
|
---|
71 | // first column of the current snode
|
---|
72 | Index fst_col; // First column within small LU update
|
---|
73 | Index segsize;
|
---|
74 | for (ksub = 0; ksub < nseg; ksub++)
|
---|
75 | {
|
---|
76 | krep = segrep(k); k--;
|
---|
77 | ksupno = glu.supno(krep);
|
---|
78 | if (jsupno != ksupno )
|
---|
79 | {
|
---|
80 | // outside the rectangular supernode
|
---|
81 | fsupc = glu.xsup(ksupno);
|
---|
82 | fst_col = (std::max)(fsupc, fpanelc);
|
---|
83 |
|
---|
84 | // Distance from the current supernode to the current panel;
|
---|
85 | // d_fsupc = 0 if fsupc > fpanelc
|
---|
86 | d_fsupc = fst_col - fsupc;
|
---|
87 |
|
---|
88 | luptr = glu.xlusup(fst_col) + d_fsupc;
|
---|
89 | lptr = glu.xlsub(fsupc) + d_fsupc;
|
---|
90 |
|
---|
91 | kfnz = repfnz(krep);
|
---|
92 | kfnz = (std::max)(kfnz, fpanelc);
|
---|
93 |
|
---|
94 | segsize = krep - kfnz + 1;
|
---|
95 | nsupc = krep - fst_col + 1;
|
---|
96 | nsupr = glu.xlsub(fsupc+1) - glu.xlsub(fsupc);
|
---|
97 | nrow = nsupr - d_fsupc - nsupc;
|
---|
98 | Index lda = glu.xlusup(fst_col+1) - glu.xlusup(fst_col);
|
---|
99 |
|
---|
100 |
|
---|
101 | // Perform a triangular solver and block update,
|
---|
102 | // then scatter the result of sup-col update to dense
|
---|
103 | no_zeros = kfnz - fst_col;
|
---|
104 | if(segsize==1)
|
---|
105 | LU_kernel_bmod<1>::run(segsize, dense, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
|
---|
106 | else
|
---|
107 | LU_kernel_bmod<Dynamic>::run(segsize, dense, tempv, glu.lusup, luptr, lda, nrow, glu.lsub, lptr, no_zeros);
|
---|
108 | } // end if jsupno
|
---|
109 | } // end for each segment
|
---|
110 |
|
---|
111 | // Process the supernodal portion of L\U[*,j]
|
---|
112 | nextlu = glu.xlusup(jcol);
|
---|
113 | fsupc = glu.xsup(jsupno);
|
---|
114 |
|
---|
115 | // copy the SPA dense into L\U[*,j]
|
---|
116 | Index mem;
|
---|
117 | new_next = nextlu + glu.xlsub(fsupc + 1) - glu.xlsub(fsupc);
|
---|
118 | Index offset = internal::first_multiple<Index>(new_next, internal::packet_traits<Scalar>::size) - new_next;
|
---|
119 | if(offset)
|
---|
120 | new_next += offset;
|
---|
121 | while (new_next > glu.nzlumax )
|
---|
122 | {
|
---|
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 | {
|
---|
129 | irow = glu.lsub(isub);
|
---|
130 | glu.lusup(nextlu) = dense(irow);
|
---|
131 | dense(irow) = Scalar(0.0);
|
---|
132 | ++nextlu;
|
---|
133 | }
|
---|
134 |
|
---|
135 | if(offset)
|
---|
136 | {
|
---|
137 | glu.lusup.segment(nextlu,offset).setZero();
|
---|
138 | nextlu += offset;
|
---|
139 | }
|
---|
140 | glu.xlusup(jcol + 1) = nextlu; // close L\U(*,jcol);
|
---|
141 |
|
---|
142 | /* For more updates within the panel (also within the current supernode),
|
---|
143 | * should start from the first column of the panel, or the first column
|
---|
144 | * of the supernode, whichever is bigger. There are two cases:
|
---|
145 | * 1) fsupc < fpanelc, then fst_col <-- fpanelc
|
---|
146 | * 2) fsupc >= fpanelc, then fst_col <-- fsupc
|
---|
147 | */
|
---|
148 | fst_col = (std::max)(fsupc, fpanelc);
|
---|
149 |
|
---|
150 | if (fst_col < jcol)
|
---|
151 | {
|
---|
152 | // Distance between the current supernode and the current panel
|
---|
153 | // d_fsupc = 0 if fsupc >= fpanelc
|
---|
154 | d_fsupc = fst_col - fsupc;
|
---|
155 |
|
---|
156 | lptr = glu.xlsub(fsupc) + d_fsupc;
|
---|
157 | luptr = glu.xlusup(fst_col) + d_fsupc;
|
---|
158 | nsupr = glu.xlsub(fsupc+1) - glu.xlsub(fsupc); // leading dimension
|
---|
159 | nsupc = jcol - fst_col; // excluding jcol
|
---|
160 | nrow = nsupr - d_fsupc - nsupc;
|
---|
161 |
|
---|
162 | // points to the beginning of jcol in snode L\U(jsupno)
|
---|
163 | ufirst = glu.xlusup(jcol) + d_fsupc;
|
---|
164 | Index lda = glu.xlusup(jcol+1) - glu.xlusup(jcol);
|
---|
165 | MappedMatrixBlock A( &(glu.lusup.data()[luptr]), nsupc, nsupc, OuterStride<>(lda) );
|
---|
166 | VectorBlock<ScalarVector> u(glu.lusup, ufirst, nsupc);
|
---|
167 | u = A.template triangularView<UnitLower>().solve(u);
|
---|
168 |
|
---|
169 | new (&A) MappedMatrixBlock ( &(glu.lusup.data()[luptr+nsupc]), nrow, nsupc, OuterStride<>(lda) );
|
---|
170 | VectorBlock<ScalarVector> l(glu.lusup, ufirst+nsupc, nrow);
|
---|
171 | l.noalias() -= A * u;
|
---|
172 |
|
---|
173 | } // End if fst_col
|
---|
174 | return 0;
|
---|
175 | }
|
---|
176 |
|
---|
177 | } // end namespace internal
|
---|
178 | } // end namespace Eigen
|
---|
179 |
|
---|
180 | #endif // SPARSELU_COLUMN_BMOD_H
|
---|