SimplicialCholesky.h 33 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008-2012 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.

#ifndef EIGEN_SIMPLICIAL_CHOLESKY_H
#define EIGEN_SIMPLICIAL_CHOLESKY_H

// IWYU pragma: private
#include "./InternalHeaderCheck.h"

namespace Eigen {

enum SimplicialCholeskyMode { SimplicialCholeskyLLT, SimplicialCholeskyLDLT };

namespace internal {
template <typename CholMatrixType, typename InputMatrixType>
struct simplicial_cholesky_grab_input {
  typedef CholMatrixType const* ConstCholMatrixPtr;
  static void run(const InputMatrixType& input, ConstCholMatrixPtr& pmat, CholMatrixType& tmp) {
    tmp = input;
    pmat = &tmp;
  }
};

template <typename MatrixType>
struct simplicial_cholesky_grab_input<MatrixType, MatrixType> {
  typedef MatrixType const* ConstMatrixPtr;
  static void run(const MatrixType& input, ConstMatrixPtr& pmat, MatrixType& /*tmp*/) { pmat = &input; }
};
}  // end namespace internal

/** \ingroup SparseCholesky_Module
 * \brief A base class for direct sparse Cholesky factorizations
 *
 * This is a base class for LL^T and LDL^T Cholesky factorizations of sparse matrices that are
 * selfadjoint and positive definite. These factorizations allow for solving A.X = B where
 * X and B can be either dense or sparse.
 *
 * In order to reduce the fill-in, a symmetric permutation P is applied prior to the factorization
 * such that the factorized matrix is P A P^-1.
 *
 * \tparam Derived the type of the derived class, that is the actual factorization type.
 *
 */
template <typename Derived>
class SimplicialCholeskyBase : public SparseSolverBase<Derived> {
  typedef SparseSolverBase<Derived> Base;
  using Base::m_isInitialized;

 public:
  typedef typename internal::traits<Derived>::MatrixType MatrixType;
  typedef typename internal::traits<Derived>::OrderingType OrderingType;
  enum { UpLo = internal::traits<Derived>::UpLo };
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::RealScalar RealScalar;
  typedef typename internal::traits<Derived>::DiagonalScalar DiagonalScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
  typedef CholMatrixType const* ConstCholMatrixPtr;
  typedef Matrix<Scalar, Dynamic, 1> VectorType;
  typedef Matrix<StorageIndex, Dynamic, 1> VectorI;

  enum { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };

 public:
  using Base::derived;

  /** Default constructor */
  SimplicialCholeskyBase()
      : m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false), m_shiftOffset(0), m_shiftScale(1) {}

  explicit SimplicialCholeskyBase(const MatrixType& matrix)
      : m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false), m_shiftOffset(0), m_shiftScale(1) {
    derived().compute(matrix);
  }

  ~SimplicialCholeskyBase() {}

  Derived& derived() { return *static_cast<Derived*>(this); }
  const Derived& derived() const { return *static_cast<const Derived*>(this); }

  inline Index cols() const { return m_matrix.cols(); }
  inline Index rows() const { return m_matrix.rows(); }

  /** \brief Reports whether previous computation was successful.
   *
   * \returns \c Success if computation was successful,
   *          \c NumericalIssue if the matrix.appears to be negative.
   */
  ComputationInfo info() const {
    eigen_assert(m_isInitialized && "Decomposition is not initialized.");
    return m_info;
  }

  /** \returns the permutation P
   * \sa permutationPinv() */
  const PermutationMatrix<Dynamic, Dynamic, StorageIndex>& permutationP() const { return m_P; }

  /** \returns the inverse P^-1 of the permutation P
   * \sa permutationP() */
  const PermutationMatrix<Dynamic, Dynamic, StorageIndex>& permutationPinv() const { return m_Pinv; }

  /** Sets the shift parameters that will be used to adjust the diagonal coefficients during the numerical
   * factorization.
   *
   * During the numerical factorization, the diagonal coefficients are transformed by the following linear model:\n
   * \c d_ii = \a offset + \a scale * \c d_ii
   *
   * The default is the identity transformation with \a offset=0, and \a scale=1.
   *
   * \returns a reference to \c *this.
   */
  Derived& setShift(const DiagonalScalar& offset, const DiagonalScalar& scale = 1) {
    m_shiftOffset = offset;
    m_shiftScale = scale;
    return derived();
  }

#ifndef EIGEN_PARSED_BY_DOXYGEN
  /** \internal */
  template <typename Stream>
  void dumpMemory(Stream& s) {
    int total = 0;
    s << "  L:        "
      << ((total += (m_matrix.cols() + 1) * sizeof(int) + m_matrix.nonZeros() * (sizeof(int) + sizeof(Scalar))) >> 20)
      << "Mb"
      << "\n";
    s << "  diag:     " << ((total += m_diag.size() * sizeof(Scalar)) >> 20) << "Mb"
      << "\n";
    s << "  tree:     " << ((total += m_parent.size() * sizeof(int)) >> 20) << "Mb"
      << "\n";
    s << "  nonzeros: " << ((total += m_workSpace.size() * sizeof(int)) >> 20) << "Mb"
      << "\n";
    s << "  perm:     " << ((total += m_P.size() * sizeof(int)) >> 20) << "Mb"
      << "\n";
    s << "  perm^-1:  " << ((total += m_Pinv.size() * sizeof(int)) >> 20) << "Mb"
      << "\n";
    s << "  TOTAL:    " << (total >> 20) << "Mb"
      << "\n";
  }

  /** \internal */
  template <typename Rhs, typename Dest>
  void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const {
    eigen_assert(m_factorizationIsOk &&
                 "The decomposition is not in a valid state for solving, you must first call either compute() or "
                 "symbolic()/numeric()");
    eigen_assert(m_matrix.rows() == b.rows());

    if (m_info != Success) return;

    if (m_P.size() > 0)
      dest = m_P * b;
    else
      dest = b;

    if (m_matrix.nonZeros() > 0)  // otherwise L==I
      derived().matrixL().solveInPlace(dest);

    if (m_diag.size() > 0) dest = m_diag.asDiagonal().inverse() * dest;

    if (m_matrix.nonZeros() > 0)  // otherwise U==I
      derived().matrixU().solveInPlace(dest);

    if (m_P.size() > 0) dest = m_Pinv * dest;
  }

  template <typename Rhs, typename Dest>
  void _solve_impl(const SparseMatrixBase<Rhs>& b, SparseMatrixBase<Dest>& dest) const {
    internal::solve_sparse_through_dense_panels(derived(), b, dest);
  }

#endif  // EIGEN_PARSED_BY_DOXYGEN

 protected:
  /** Computes the sparse Cholesky decomposition of \a matrix */
  template <bool DoLDLT, bool NonHermitian>
  void compute(const MatrixType& matrix) {
    eigen_assert(matrix.rows() == matrix.cols());
    Index size = matrix.cols();
    CholMatrixType tmp(size, size);
    ConstCholMatrixPtr pmat;
    ordering<NonHermitian>(matrix, pmat, tmp);
    analyzePattern_preordered(*pmat, DoLDLT);
    factorize_preordered<DoLDLT, NonHermitian>(*pmat);
  }

  template <bool DoLDLT, bool NonHermitian>
  void factorize(const MatrixType& a) {
    eigen_assert(a.rows() == a.cols());
    Index size = a.cols();
    CholMatrixType tmp(size, size);
    ConstCholMatrixPtr pmat;

    if (m_P.size() == 0 && (int(UpLo) & int(Upper)) == Upper) {
      // If there is no ordering, try to directly use the input matrix without any copy
      internal::simplicial_cholesky_grab_input<CholMatrixType, MatrixType>::run(a, pmat, tmp);
    } else {
      internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, tmp, m_P.indices().data());
      pmat = &tmp;
    }

    factorize_preordered<DoLDLT, NonHermitian>(*pmat);
  }

  template <bool DoLDLT, bool NonHermitian>
  void factorize_preordered(const CholMatrixType& a);

  template <bool DoLDLT, bool NonHermitian>
  void analyzePattern(const MatrixType& a) {
    eigen_assert(a.rows() == a.cols());
    Index size = a.cols();
    CholMatrixType tmp(size, size);
    ConstCholMatrixPtr pmat;
    ordering<NonHermitian>(a, pmat, tmp);
    analyzePattern_preordered(*pmat, DoLDLT);
  }
  void analyzePattern_preordered(const CholMatrixType& a, bool doLDLT);

  template <bool NonHermitian>
  void ordering(const MatrixType& a, ConstCholMatrixPtr& pmat, CholMatrixType& ap);

  inline DiagonalScalar getDiag(Scalar x) { return internal::traits<Derived>::getDiag(x); }
  inline Scalar getSymm(Scalar x) { return internal::traits<Derived>::getSymm(x); }

  /** keeps off-diagonal entries; drops diagonal entries */
  struct keep_diag {
    inline bool operator()(const Index& row, const Index& col, const Scalar&) const { return row != col; }
  };

  mutable ComputationInfo m_info;
  bool m_factorizationIsOk;
  bool m_analysisIsOk;

  CholMatrixType m_matrix;
  VectorType m_diag;  // the diagonal coefficients (LDLT mode)
  VectorI m_parent;   // elimination tree
  VectorI m_workSpace;
  PermutationMatrix<Dynamic, Dynamic, StorageIndex> m_P;     // the permutation
  PermutationMatrix<Dynamic, Dynamic, StorageIndex> m_Pinv;  // the inverse permutation

  DiagonalScalar m_shiftOffset;
  DiagonalScalar m_shiftScale;
};

template <typename MatrixType_, int UpLo_ = Lower,
          typename Ordering_ = AMDOrdering<typename MatrixType_::StorageIndex> >
class SimplicialLLT;
template <typename MatrixType_, int UpLo_ = Lower,
          typename Ordering_ = AMDOrdering<typename MatrixType_::StorageIndex> >
class SimplicialLDLT;
template <typename MatrixType_, int UpLo_ = Lower,
          typename Ordering_ = AMDOrdering<typename MatrixType_::StorageIndex> >
class SimplicialNonHermitianLLT;
template <typename MatrixType_, int UpLo_ = Lower,
          typename Ordering_ = AMDOrdering<typename MatrixType_::StorageIndex> >
class SimplicialNonHermitianLDLT;
template <typename MatrixType_, int UpLo_ = Lower,
          typename Ordering_ = AMDOrdering<typename MatrixType_::StorageIndex> >
class SimplicialCholesky;

namespace internal {

template <typename MatrixType_, int UpLo_, typename Ordering_>
struct traits<SimplicialLLT<MatrixType_, UpLo_, Ordering_> > {
  typedef MatrixType_ MatrixType;
  typedef Ordering_ OrderingType;
  enum { UpLo = UpLo_ };
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::RealScalar DiagonalScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
  typedef TriangularView<const CholMatrixType, Eigen::Lower> MatrixL;
  typedef TriangularView<const typename CholMatrixType::AdjointReturnType, Eigen::Upper> MatrixU;
  static inline MatrixL getL(const CholMatrixType& m) { return MatrixL(m); }
  static inline MatrixU getU(const CholMatrixType& m) { return MatrixU(m.adjoint()); }
  static inline DiagonalScalar getDiag(Scalar x) { return numext::real(x); }
  static inline Scalar getSymm(Scalar x) { return numext::conj(x); }
};

template <typename MatrixType_, int UpLo_, typename Ordering_>
struct traits<SimplicialLDLT<MatrixType_, UpLo_, Ordering_> > {
  typedef MatrixType_ MatrixType;
  typedef Ordering_ OrderingType;
  enum { UpLo = UpLo_ };
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::RealScalar DiagonalScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
  typedef TriangularView<const CholMatrixType, Eigen::UnitLower> MatrixL;
  typedef TriangularView<const typename CholMatrixType::AdjointReturnType, Eigen::UnitUpper> MatrixU;
  static inline MatrixL getL(const CholMatrixType& m) { return MatrixL(m); }
  static inline MatrixU getU(const CholMatrixType& m) { return MatrixU(m.adjoint()); }
  static inline DiagonalScalar getDiag(Scalar x) { return numext::real(x); }
  static inline Scalar getSymm(Scalar x) { return numext::conj(x); }
};

template <typename MatrixType_, int UpLo_, typename Ordering_>
struct traits<SimplicialNonHermitianLLT<MatrixType_, UpLo_, Ordering_> > {
  typedef MatrixType_ MatrixType;
  typedef Ordering_ OrderingType;
  enum { UpLo = UpLo_ };
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::Scalar DiagonalScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
  typedef TriangularView<const CholMatrixType, Eigen::Lower> MatrixL;
  typedef TriangularView<const typename CholMatrixType::ConstTransposeReturnType, Eigen::Upper> MatrixU;
  static inline MatrixL getL(const CholMatrixType& m) { return MatrixL(m); }
  static inline MatrixU getU(const CholMatrixType& m) { return MatrixU(m.transpose()); }
  static inline DiagonalScalar getDiag(Scalar x) { return x; }
  static inline Scalar getSymm(Scalar x) { return x; }
};

template <typename MatrixType_, int UpLo_, typename Ordering_>
struct traits<SimplicialNonHermitianLDLT<MatrixType_, UpLo_, Ordering_> > {
  typedef MatrixType_ MatrixType;
  typedef Ordering_ OrderingType;
  enum { UpLo = UpLo_ };
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::Scalar DiagonalScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
  typedef TriangularView<const CholMatrixType, Eigen::UnitLower> MatrixL;
  typedef TriangularView<const typename CholMatrixType::ConstTransposeReturnType, Eigen::UnitUpper> MatrixU;
  static inline MatrixL getL(const CholMatrixType& m) { return MatrixL(m); }
  static inline MatrixU getU(const CholMatrixType& m) { return MatrixU(m.transpose()); }
  static inline DiagonalScalar getDiag(Scalar x) { return x; }
  static inline Scalar getSymm(Scalar x) { return x; }
};

template <typename MatrixType_, int UpLo_, typename Ordering_>
struct traits<SimplicialCholesky<MatrixType_, UpLo_, Ordering_> > {
  typedef MatrixType_ MatrixType;
  typedef Ordering_ OrderingType;
  enum { UpLo = UpLo_ };
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::RealScalar DiagonalScalar;
  static inline DiagonalScalar getDiag(Scalar x) { return numext::real(x); }
  static inline Scalar getSymm(Scalar x) { return numext::conj(x); }
};

}  // namespace internal

/** \ingroup SparseCholesky_Module
 * \class SimplicialLLT
 * \brief A direct sparse LLT Cholesky factorizations
 *
 * This class provides a LL^T Cholesky factorizations of sparse matrices that are
 * selfadjoint and positive definite. The factorization allows for solving A.X = B where
 * X and B can be either dense or sparse.
 *
 * In order to reduce the fill-in, a symmetric permutation P is applied prior to the factorization
 * such that the factorized matrix is P A P^-1.
 *
 * \tparam MatrixType_ the type of the sparse matrix A, it must be a SparseMatrix<>
 * \tparam UpLo_ the triangular part that will be used for the computations. It can be Lower
 *               or Upper. Default is Lower.
 * \tparam Ordering_ The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>
 *
 * \implsparsesolverconcept
 *
 * \sa class SimplicialLDLT, class AMDOrdering, class NaturalOrdering
 */
template <typename MatrixType_, int UpLo_, typename Ordering_>
class SimplicialLLT : public SimplicialCholeskyBase<SimplicialLLT<MatrixType_, UpLo_, Ordering_> > {
 public:
  typedef MatrixType_ MatrixType;
  enum { UpLo = UpLo_ };
  typedef SimplicialCholeskyBase<SimplicialLLT> Base;
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::RealScalar RealScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, Index> CholMatrixType;
  typedef Matrix<Scalar, Dynamic, 1> VectorType;
  typedef internal::traits<SimplicialLLT> Traits;
  typedef typename Traits::MatrixL MatrixL;
  typedef typename Traits::MatrixU MatrixU;

 public:
  /** Default constructor */
  SimplicialLLT() : Base() {}
  /** Constructs and performs the LLT factorization of \a matrix */
  explicit SimplicialLLT(const MatrixType& matrix) : Base(matrix) {}

  /** \returns an expression of the factor L */
  inline const MatrixL matrixL() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
    return Traits::getL(Base::m_matrix);
  }

  /** \returns an expression of the factor U (= L^*) */
  inline const MatrixU matrixU() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
    return Traits::getU(Base::m_matrix);
  }

  /** Computes the sparse Cholesky decomposition of \a matrix */
  SimplicialLLT& compute(const MatrixType& matrix) {
    Base::template compute<false, false>(matrix);
    return *this;
  }

  /** Performs a symbolic decomposition on the sparcity of \a matrix.
   *
   * This function is particularly useful when solving for several problems having the same structure.
   *
   * \sa factorize()
   */
  void analyzePattern(const MatrixType& a) { Base::template analyzePattern<false, false>(a); }

  /** Performs a numeric decomposition of \a matrix
   *
   * The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.
   *
   * \sa analyzePattern()
   */
  void factorize(const MatrixType& a) { Base::template factorize<false, false>(a); }

  /** \returns the determinant of the underlying matrix from the current factorization */
  Scalar determinant() const {
    Scalar detL = Base::m_matrix.diagonal().prod();
    return numext::abs2(detL);
  }
};

/** \ingroup SparseCholesky_Module
 * \class SimplicialLDLT
 * \brief A direct sparse LDLT Cholesky factorizations without square root.
 *
 * This class provides a LDL^T Cholesky factorizations without square root of sparse matrices that are
 * selfadjoint and positive definite. The factorization allows for solving A.X = B where
 * X and B can be either dense or sparse.
 *
 * In order to reduce the fill-in, a symmetric permutation P is applied prior to the factorization
 * such that the factorized matrix is P A P^-1.
 *
 * \tparam MatrixType_ the type of the sparse matrix A, it must be a SparseMatrix<>
 * \tparam UpLo_ the triangular part that will be used for the computations. It can be Lower
 *               or Upper. Default is Lower.
 * \tparam Ordering_ The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>
 *
 * \implsparsesolverconcept
 *
 * \sa class SimplicialLLT, class AMDOrdering, class NaturalOrdering
 */
template <typename MatrixType_, int UpLo_, typename Ordering_>
class SimplicialLDLT : public SimplicialCholeskyBase<SimplicialLDLT<MatrixType_, UpLo_, Ordering_> > {
 public:
  typedef MatrixType_ MatrixType;
  enum { UpLo = UpLo_ };
  typedef SimplicialCholeskyBase<SimplicialLDLT> Base;
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::RealScalar RealScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
  typedef Matrix<Scalar, Dynamic, 1> VectorType;
  typedef internal::traits<SimplicialLDLT> Traits;
  typedef typename Traits::MatrixL MatrixL;
  typedef typename Traits::MatrixU MatrixU;

 public:
  /** Default constructor */
  SimplicialLDLT() : Base() {}

  /** Constructs and performs the LLT factorization of \a matrix */
  explicit SimplicialLDLT(const MatrixType& matrix) : Base(matrix) {}

  /** \returns a vector expression of the diagonal D */
  inline const VectorType vectorD() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
    return Base::m_diag;
  }
  /** \returns an expression of the factor L */
  inline const MatrixL matrixL() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
    return Traits::getL(Base::m_matrix);
  }

  /** \returns an expression of the factor U (= L^*) */
  inline const MatrixU matrixU() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
    return Traits::getU(Base::m_matrix);
  }

  /** Computes the sparse Cholesky decomposition of \a matrix */
  SimplicialLDLT& compute(const MatrixType& matrix) {
    Base::template compute<true, false>(matrix);
    return *this;
  }

  /** Performs a symbolic decomposition on the sparcity of \a matrix.
   *
   * This function is particularly useful when solving for several problems having the same structure.
   *
   * \sa factorize()
   */
  void analyzePattern(const MatrixType& a) { Base::template analyzePattern<true, false>(a); }

  /** Performs a numeric decomposition of \a matrix
   *
   * The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.
   *
   * \sa analyzePattern()
   */
  void factorize(const MatrixType& a) { Base::template factorize<true, false>(a); }

  /** \returns the determinant of the underlying matrix from the current factorization */
  Scalar determinant() const { return Base::m_diag.prod(); }
};

/** \ingroup SparseCholesky_Module
 * \class SimplicialNonHermitianLLT
 * \brief A direct sparse LLT Cholesky factorizations, for symmetric non-hermitian matrices.
 *
 * This class provides a LL^T Cholesky factorizations of sparse matrices that are
 * symmetric but not hermitian. For real matrices, this is equivalent to the regular LLT factorization.
 * The factorization allows for solving A.X = B where X and B can be either dense or sparse.
 *
 * In order to reduce the fill-in, a symmetric permutation P is applied prior to the factorization
 * such that the factorized matrix is P A P^-1.
 *
 * \tparam MatrixType_ the type of the sparse matrix A, it must be a SparseMatrix<>
 * \tparam UpLo_ the triangular part that will be used for the computations. It can be Lower
 *               or Upper. Default is Lower.
 * \tparam Ordering_ The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>
 *
 * \implsparsesolverconcept
 *
 * \sa class SimplicialNonHermitianLDLT, SimplicialLLT, class AMDOrdering, class NaturalOrdering
 */
template <typename MatrixType_, int UpLo_, typename Ordering_>
class SimplicialNonHermitianLLT
    : public SimplicialCholeskyBase<SimplicialNonHermitianLLT<MatrixType_, UpLo_, Ordering_> > {
 public:
  typedef MatrixType_ MatrixType;
  enum { UpLo = UpLo_ };
  typedef SimplicialCholeskyBase<SimplicialNonHermitianLLT> Base;
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::RealScalar RealScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
  typedef Matrix<Scalar, Dynamic, 1> VectorType;
  typedef internal::traits<SimplicialNonHermitianLLT> Traits;
  typedef typename Traits::MatrixL MatrixL;
  typedef typename Traits::MatrixU MatrixU;

 public:
  /** Default constructor */
  SimplicialNonHermitianLLT() : Base() {}

  /** Constructs and performs the LLT factorization of \a matrix */
  explicit SimplicialNonHermitianLLT(const MatrixType& matrix) : Base(matrix) {}

  /** \returns an expression of the factor L */
  inline const MatrixL matrixL() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
    return Traits::getL(Base::m_matrix);
  }

  /** \returns an expression of the factor U (= L^*) */
  inline const MatrixU matrixU() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LLT not factorized");
    return Traits::getU(Base::m_matrix);
  }

  /** Computes the sparse Cholesky decomposition of \a matrix */
  SimplicialNonHermitianLLT& compute(const MatrixType& matrix) {
    Base::template compute<false, true>(matrix);
    return *this;
  }

  /** Performs a symbolic decomposition on the sparcity of \a matrix.
   *
   * This function is particularly useful when solving for several problems having the same structure.
   *
   * \sa factorize()
   */
  void analyzePattern(const MatrixType& a) { Base::template analyzePattern<false, true>(a); }

  /** Performs a numeric decomposition of \a matrix
   *
   * The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.
   *
   * \sa analyzePattern()
   */
  void factorize(const MatrixType& a) { Base::template factorize<false, true>(a); }

  /** \returns the determinant of the underlying matrix from the current factorization */
  Scalar determinant() const {
    Scalar detL = Base::m_matrix.diagonal().prod();
    return detL * detL;
  }
};

/** \ingroup SparseCholesky_Module
 * \class SimplicialNonHermitianLDLT
 * \brief A direct sparse LDLT Cholesky factorizations without square root, for symmetric non-hermitian matrices.
 *
 * This class provides a LDL^T Cholesky factorizations without square root of sparse matrices that are
 * symmetric but not hermitian. For real matrices, this is equivalent to the regular LDLT factorization.
 * The factorization allows for solving A.X = B where X and B can be either dense or sparse.
 *
 * In order to reduce the fill-in, a symmetric permutation P is applied prior to the factorization
 * such that the factorized matrix is P A P^-1.
 *
 * \tparam MatrixType_ the type of the sparse matrix A, it must be a SparseMatrix<>
 * \tparam UpLo_ the triangular part that will be used for the computations. It can be Lower
 *               or Upper. Default is Lower.
 * \tparam Ordering_ The ordering method to use, either AMDOrdering<> or NaturalOrdering<>. Default is AMDOrdering<>
 *
 * \implsparsesolverconcept
 *
 * \sa class SimplicialNonHermitianLLT, SimplicialLDLT, class AMDOrdering, class NaturalOrdering
 */
template <typename MatrixType_, int UpLo_, typename Ordering_>
class SimplicialNonHermitianLDLT
    : public SimplicialCholeskyBase<SimplicialNonHermitianLDLT<MatrixType_, UpLo_, Ordering_> > {
 public:
  typedef MatrixType_ MatrixType;
  enum { UpLo = UpLo_ };
  typedef SimplicialCholeskyBase<SimplicialNonHermitianLDLT> Base;
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::RealScalar RealScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
  typedef Matrix<Scalar, Dynamic, 1> VectorType;
  typedef internal::traits<SimplicialNonHermitianLDLT> Traits;
  typedef typename Traits::MatrixL MatrixL;
  typedef typename Traits::MatrixU MatrixU;

 public:
  /** Default constructor */
  SimplicialNonHermitianLDLT() : Base() {}

  /** Constructs and performs the LLT factorization of \a matrix */
  explicit SimplicialNonHermitianLDLT(const MatrixType& matrix) : Base(matrix) {}

  /** \returns a vector expression of the diagonal D */
  inline const VectorType vectorD() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
    return Base::m_diag;
  }
  /** \returns an expression of the factor L */
  inline const MatrixL matrixL() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
    return Traits::getL(Base::m_matrix);
  }

  /** \returns an expression of the factor U (= L^*) */
  inline const MatrixU matrixU() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial LDLT not factorized");
    return Traits::getU(Base::m_matrix);
  }

  /** Computes the sparse Cholesky decomposition of \a matrix */
  SimplicialNonHermitianLDLT& compute(const MatrixType& matrix) {
    Base::template compute<true, true>(matrix);
    return *this;
  }

  /** Performs a symbolic decomposition on the sparcity of \a matrix.
   *
   * This function is particularly useful when solving for several problems having the same structure.
   *
   * \sa factorize()
   */
  void analyzePattern(const MatrixType& a) { Base::template analyzePattern<true, true>(a); }

  /** Performs a numeric decomposition of \a matrix
   *
   * The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.
   *
   * \sa analyzePattern()
   */
  void factorize(const MatrixType& a) { Base::template factorize<true, true>(a); }

  /** \returns the determinant of the underlying matrix from the current factorization */
  Scalar determinant() const { return Base::m_diag.prod(); }
};

/** \deprecated use SimplicialLDLT or class SimplicialLLT
 * \ingroup SparseCholesky_Module
 * \class SimplicialCholesky
 *
 * \sa class SimplicialLDLT, class SimplicialLLT
 */
template <typename MatrixType_, int UpLo_, typename Ordering_>
class SimplicialCholesky : public SimplicialCholeskyBase<SimplicialCholesky<MatrixType_, UpLo_, Ordering_> > {
 public:
  typedef MatrixType_ MatrixType;
  enum { UpLo = UpLo_ };
  typedef SimplicialCholeskyBase<SimplicialCholesky> Base;
  typedef typename MatrixType::Scalar Scalar;
  typedef typename MatrixType::RealScalar RealScalar;
  typedef typename MatrixType::StorageIndex StorageIndex;
  typedef SparseMatrix<Scalar, ColMajor, StorageIndex> CholMatrixType;
  typedef Matrix<Scalar, Dynamic, 1> VectorType;
  typedef internal::traits<SimplicialLDLT<MatrixType, UpLo> > LDLTTraits;
  typedef internal::traits<SimplicialLLT<MatrixType, UpLo> > LLTTraits;

 public:
  SimplicialCholesky() : Base(), m_LDLT(true) {}

  explicit SimplicialCholesky(const MatrixType& matrix) : Base(), m_LDLT(true) { compute(matrix); }

  SimplicialCholesky& setMode(SimplicialCholeskyMode mode) {
    switch (mode) {
      case SimplicialCholeskyLLT:
        m_LDLT = false;
        break;
      case SimplicialCholeskyLDLT:
        m_LDLT = true;
        break;
      default:
        break;
    }

    return *this;
  }

  inline const VectorType vectorD() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial Cholesky not factorized");
    return Base::m_diag;
  }
  inline const CholMatrixType rawMatrix() const {
    eigen_assert(Base::m_factorizationIsOk && "Simplicial Cholesky not factorized");
    return Base::m_matrix;
  }

  /** Computes the sparse Cholesky decomposition of \a matrix */
  SimplicialCholesky& compute(const MatrixType& matrix) {
    if (m_LDLT)
      Base::template compute<true, false>(matrix);
    else
      Base::template compute<false, false>(matrix);
    return *this;
  }

  /** Performs a symbolic decomposition on the sparcity of \a matrix.
   *
   * This function is particularly useful when solving for several problems having the same structure.
   *
   * \sa factorize()
   */
  void analyzePattern(const MatrixType& a) {
    if (m_LDLT)
      Base::template analyzePattern<true, false>(a);
    else
      Base::template analyzePattern<false, false>(a);
  }

  /** Performs a numeric decomposition of \a matrix
   *
   * The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.
   *
   * \sa analyzePattern()
   */
  void factorize(const MatrixType& a) {
    if (m_LDLT)
      Base::template factorize<true, false>(a);
    else
      Base::template factorize<false, false>(a);
  }

  /** \internal */
  template <typename Rhs, typename Dest>
  void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const {
    eigen_assert(Base::m_factorizationIsOk &&
                 "The decomposition is not in a valid state for solving, you must first call either compute() or "
                 "symbolic()/numeric()");
    eigen_assert(Base::m_matrix.rows() == b.rows());

    if (Base::m_info != Success) return;

    if (Base::m_P.size() > 0)
      dest = Base::m_P * b;
    else
      dest = b;

    if (Base::m_matrix.nonZeros() > 0)  // otherwise L==I
    {
      if (m_LDLT)
        LDLTTraits::getL(Base::m_matrix).solveInPlace(dest);
      else
        LLTTraits::getL(Base::m_matrix).solveInPlace(dest);
    }

    if (Base::m_diag.size() > 0) dest = Base::m_diag.real().asDiagonal().inverse() * dest;

    if (Base::m_matrix.nonZeros() > 0)  // otherwise I==I
    {
      if (m_LDLT)
        LDLTTraits::getU(Base::m_matrix).solveInPlace(dest);
      else
        LLTTraits::getU(Base::m_matrix).solveInPlace(dest);
    }

    if (Base::m_P.size() > 0) dest = Base::m_Pinv * dest;
  }

  /** \internal */
  template <typename Rhs, typename Dest>
  void _solve_impl(const SparseMatrixBase<Rhs>& b, SparseMatrixBase<Dest>& dest) const {
    internal::solve_sparse_through_dense_panels(*this, b, dest);
  }

  Scalar determinant() const {
    if (m_LDLT) {
      return Base::m_diag.prod();
    } else {
      Scalar detL = Diagonal<const CholMatrixType>(Base::m_matrix).prod();
      return numext::abs2(detL);
    }
  }

 protected:
  bool m_LDLT;
};

template <typename Derived>
template <bool NonHermitian>
void SimplicialCholeskyBase<Derived>::ordering(const MatrixType& a, ConstCholMatrixPtr& pmat, CholMatrixType& ap) {
  eigen_assert(a.rows() == a.cols());
  const Index size = a.rows();
  pmat = &ap;
  // Note that ordering methods compute the inverse permutation
  if (!internal::is_same<OrderingType, NaturalOrdering<StorageIndex> >::value) {
    {
      CholMatrixType C;
      internal::permute_symm_to_fullsymm<UpLo, NonHermitian>(a, C, NULL);

      OrderingType ordering;
      ordering(C, m_Pinv);
    }

    if (m_Pinv.size() > 0)
      m_P = m_Pinv.inverse();
    else
      m_P.resize(0);

    ap.resize(size, size);
    internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, ap, m_P.indices().data());
  } else {
    m_Pinv.resize(0);
    m_P.resize(0);
    if (int(UpLo) == int(Lower) || MatrixType::IsRowMajor) {
      // we have to transpose the lower part to to the upper one
      ap.resize(size, size);
      internal::permute_symm_to_symm<UpLo, Upper, NonHermitian>(a, ap, NULL);
    } else
      internal::simplicial_cholesky_grab_input<CholMatrixType, MatrixType>::run(a, pmat, ap);
  }
}

}  // end namespace Eigen

#endif  // EIGEN_SIMPLICIAL_CHOLESKY_H