123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299 |
- #ifndef CERES_INTERNAL_SUITESPARSE_H_
- #define CERES_INTERNAL_SUITESPARSE_H_
- #include "ceres/internal/config.h"
- #ifndef CERES_NO_SUITESPARSE
- #include <cstring>
- #include <memory>
- #include <string>
- #include <vector>
- #include "SuiteSparseQR.hpp"
- #include "ceres/block_structure.h"
- #include "ceres/internal/disable_warnings.h"
- #include "ceres/linear_solver.h"
- #include "ceres/sparse_cholesky.h"
- #include "cholmod.h"
- #include "glog/logging.h"
- namespace ceres::internal {
- class CompressedRowSparseMatrix;
- class TripletSparseMatrix;
- class CERES_NO_EXPORT SuiteSparse {
- public:
- SuiteSparse();
- ~SuiteSparse();
-
-
-
- cholmod_sparse* CreateSparseMatrix(TripletSparseMatrix* A);
-
-
- cholmod_sparse* CreateSparseMatrixTranspose(TripletSparseMatrix* A);
-
-
-
- cholmod_sparse CreateSparseMatrixTransposeView(CompressedRowSparseMatrix* A);
-
-
-
- cholmod_dense CreateDenseVectorView(const double* x, int size);
-
-
-
- cholmod_dense* CreateDenseVector(const double* x, int in_size, int out_size);
-
-
-
-
-
-
- void Scale(cholmod_dense* scale, int mode, cholmod_sparse* A) {
- cholmod_scale(scale, mode, A, &cc_);
- }
-
-
- cholmod_sparse* AATranspose(cholmod_sparse* A) {
- cholmod_sparse* m = cholmod_aat(A, nullptr, A->nrow, 1, &cc_);
- m->stype = 1;
- return m;
- }
-
- void SparseDenseMultiply(cholmod_sparse* A,
- double alpha,
- double beta,
- cholmod_dense* x,
- cholmod_dense* y) {
- double alpha_[2] = {alpha, 0};
- double beta_[2] = {beta, 0};
- cholmod_sdmult(A, 0, alpha_, beta_, x, y, &cc_);
- }
-
-
-
-
-
-
-
-
-
-
-
-
- cholmod_factor* AnalyzeCholesky(cholmod_sparse* A,
- OrderingType ordering_type,
- std::string* message);
-
- cholmod_factor* BlockAnalyzeCholesky(cholmod_sparse* A,
- OrderingType ordering_type,
- const std::vector<Block>& row_blocks,
- const std::vector<Block>& col_blocks,
- std::string* message);
-
-
-
-
-
-
-
-
-
-
-
- cholmod_factor* AnalyzeCholeskyWithGivenOrdering(
- cholmod_sparse* A,
- const std::vector<int>& ordering,
- std::string* message);
-
-
-
-
-
-
- LinearSolverTerminationType Cholesky(cholmod_sparse* A,
- cholmod_factor* L,
- std::string* message);
-
-
-
-
-
- cholmod_dense* Solve(cholmod_factor* L,
- cholmod_dense* b,
- std::string* message);
-
-
- bool Ordering(cholmod_sparse* matrix,
- OrderingType ordering_type,
- int* ordering);
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
-
- bool BlockOrdering(const cholmod_sparse* A,
- OrderingType ordering_type,
- const std::vector<Block>& row_blocks,
- const std::vector<Block>& col_blocks,
- std::vector<int>* ordering);
-
-
- static bool IsNestedDissectionAvailable();
-
-
-
-
-
-
-
-
-
-
- bool ConstrainedApproximateMinimumDegreeOrdering(cholmod_sparse* matrix,
- int* constraints,
- int* ordering);
- void Free(cholmod_sparse* m) { cholmod_free_sparse(&m, &cc_); }
- void Free(cholmod_dense* m) { cholmod_free_dense(&m, &cc_); }
- void Free(cholmod_factor* m) { cholmod_free_factor(&m, &cc_); }
- void Print(cholmod_sparse* m, const std::string& name) {
- cholmod_print_sparse(m, const_cast<char*>(name.c_str()), &cc_);
- }
- void Print(cholmod_dense* m, const std::string& name) {
- cholmod_print_dense(m, const_cast<char*>(name.c_str()), &cc_);
- }
- void Print(cholmod_triplet* m, const std::string& name) {
- cholmod_print_triplet(m, const_cast<char*>(name.c_str()), &cc_);
- }
- cholmod_common* mutable_cc() { return &cc_; }
- private:
- cholmod_common cc_;
- };
- class CERES_NO_EXPORT SuiteSparseCholesky final : public SparseCholesky {
- public:
- static std::unique_ptr<SparseCholesky> Create(OrderingType ordering_type);
-
- ~SuiteSparseCholesky() override;
- CompressedRowSparseMatrix::StorageType StorageType() const final;
- LinearSolverTerminationType Factorize(CompressedRowSparseMatrix* lhs,
- std::string* message) final;
- LinearSolverTerminationType Solve(const double* rhs,
- double* solution,
- std::string* message) final;
- private:
- explicit SuiteSparseCholesky(const OrderingType ordering_type);
- const OrderingType ordering_type_;
- SuiteSparse ss_;
- cholmod_factor* factor_;
- };
- }
- #include "ceres/internal/reenable_warnings.h"
- #else
- using cholmod_factor = void;
- #include "ceres/internal/disable_warnings.h"
- namespace ceres {
- namespace internal {
- class CERES_NO_EXPORT SuiteSparse {
- public:
-
-
- static bool IsNestedDissectionAvailable() { return false; }
- void Free(void* ) {}
- };
- }
- }
- #include "ceres/internal/reenable_warnings.h"
- #endif
- #endif
|