TMB Documentation  v1.9.11
simplicial_inverse_subset.hpp
1 namespace Eigen {
2 
3 template<class T>
4 struct SimplicialInverseSubset {
5  // typedefs
6  typedef SimplicialLLT<SparseMatrix<T> > Base;
7  typedef SimplicialLLT<SparseMatrix<double> > BaseD;
8  typedef typename Base::StorageIndex StorageIndex;
9  typedef typename Base::Scalar Scalar;
10  typedef SparseMatrix<T> SpMat;
11  enum Operation {
12  //IndexScatter,
13  Scatter,
14  Zero,
15  InnerProduct
16  };
17  // workspaces
18  std::shared_ptr<Base> factor;
19  std::vector<StorageIndex> idg;
20  std::vector<StorageIndex> index_gather(const SpMat &mat) {
21  SpMat L = factor->matrixL();
22  const StorageIndex* Lp = L.outerIndexPtr();
23  const StorageIndex* Li = L.innerIndexPtr();
24  std::vector<StorageIndex> ans;
25  const StorageIndex *Perm = factor->permutationP().indices().data();
26  std::vector<StorageIndex> wrk(mat.rows());
27  for (StorageIndex k=0; k<mat.outerSize(); ++k) {
28  // Index scatter
29  StorageIndex j = Perm[k];
30  for (StorageIndex l=Lp[j]; l<Lp[j+1]; l++) {
31  StorageIndex i = Li[l];
32  wrk[i] = l;
33  }
34  for (typename SpMat::InnerIterator it(mat, k); it; ++it) {
35  if (Perm[it.row()] < Perm[it.col()]) ans.push_back(-1); // Skip upper triangle
36  else ans.push_back( wrk[Perm[it.row()]] );
37  }
38  }
39  return ans;
40  }
41  void valuesGet(SpMat &mat, const SpMat &S) {
42  if (idg.size() == 0)
43  idg = index_gather(mat);
44  T* vptr = mat.valuePtr();
45  const T* ans = S.valuePtr();
46  for (size_t i=0; i<idg.size(); i++) {
47  if (idg[i] != -1) {
48  vptr[i] = ans[idg[i]];
49  }
50  }
51  }
52  template<Operation Op, class Scalar_, class Scalar>
53  Scalar column(SparseMatrix<Scalar_> &L,
54  StorageIndex j,
55  Scalar* x) {
56  const StorageIndex* Lp = L.outerIndexPtr();
57  const StorageIndex* Li = L.innerIndexPtr();
58  Scalar_* Lx = L.valuePtr();
59  Scalar s = 0;
60  for (StorageIndex k=Lp[j]; k<Lp[j+1]; k++) {
61  StorageIndex i = Li[k];
62  // if (Op == IndexScatter) {
63  // x[i] = k;
64  // }
65  if (Op == Scatter) {
66  x[i] = Lx[k];
67  }
68  if (Op == Zero) {
69  x[i] = 0.;
70  }
71  if (Op == InnerProduct) {
72  s += Lx[k] * x[i];
73  }
74  }
75  return s;
76  }
77  Eigen::SparseMatrix<StorageIndex> LT;
78  void init_transpose(SpMat L) {
79  if (LT.rows() > 0) return;
80  std::vector<StorageIndex> x(L.nonZeros());
81  for (size_t i=0; i<x.size(); i++) x[i] = i;
82  Eigen::SparseMatrix<StorageIndex> TMP =
83  Eigen::Map<const Eigen::SparseMatrix<StorageIndex> > (L.rows(),
84  L.cols(),
85  L.nonZeros(),
86  L.outerIndexPtr(),
87  L.innerIndexPtr(),
88  x.data(),
89  L.innerNonZeroPtr());
90  LT = TMP.transpose();
91  }
92  SpMat chol2inv() {
93  SpMat L = factor->matrixL();
94  init_transpose(L);
95  // Allocate result
96  SpMat S(L);
97  for (StorageIndex i=0; i<S.nonZeros(); i++)
98  S.valuePtr()[i] = 0;
99  const StorageIndex ncol = L.cols();
100  // L pointers
101  const StorageIndex* Lp = L.outerIndexPtr();
102  Scalar* Lx = L.valuePtr();
103  // S pointers (pattern same as L)
104  Scalar* Sx = S.valuePtr();
105  // LT pointers
106  const StorageIndex* LTp = LT.outerIndexPtr();
107  StorageIndex* LTi = LT.innerIndexPtr();
108  StorageIndex* LTx = LT.valuePtr();
109  // Workspace dense row
110  std::vector<T> wrk(ncol, 0);
111  T* S_row = wrk.data();
112  // Recursions
113  for (StorageIndex j = ncol; j > 0; ) {
114  j--;
115  // Scatter sub-column S(j:n,j)
116  column<Scatter> (S, j, S_row);
117  // Diagonal element S(j,j)
118  Scalar Sjj = 0;
119  for (StorageIndex k = Lp[j] + 1; k < Lp[j+1]; k++) {
120  Sjj += Lx[k] * Sx[k];
121  }
122  Scalar Ljj_inv = 1. / Lx[Lp[j]];
123  S_row[j] = -Ljj_inv * Sjj + Ljj_inv * Ljj_inv;
124  // Recursions for j'th row S(j,1:(j-1))
125  for (StorageIndex k = LTp[j+1]-1; k > LTp[j]; ) {
126  k--;
127  StorageIndex i = LTi[k];
128  S_row[i] = (-1. / Lx[Lp[i]]) * column<InnerProduct> (L, i, S_row);
129  }
130  for (StorageIndex k = LTp[j]; k < LTp[j+1]; k++) {
131  // Insert result in sparse output
132  StorageIndex i = LTi[k];
133  StorageIndex kT = LTx[k];
134  Sx[kT] = S_row[i];
135  }
136  // Clear workspace
137  column<Zero> (L, j, S_row);
138  column<Zero> (LT, j, S_row);
139  }
140  return S;
141  }
142  SpMat operator()(SpMat x) {
143  // Factor initialize
144  if (!factor) {
145  factor = std::make_shared<Base>(x);
146  }
147  // forward recursions
148  factor->factorize(x);
149  // Reverse
150  SpMat S = chol2inv();
151  // Get result
152  x = x * 0;
153  valuesGet (x, S);
154  return x;
155  }
156  // No simple way to change factor<double> to other numerical type
157  // SimplicialInverseSubset(std::shared_ptr<Base> factor) :
158  // factor(factor) { }
159  SimplicialInverseSubset(std::shared_ptr<BaseD> factor) { }
160 };
161 
162 } // Eigen
std::vector< I > factor(const std::vector< T > &x)
Replace each element in a vector by an integer code such that x[i] == x[j] implies f[i] == f[j] ...
Definition: radix.hpp:116
License: GPL v2