SparseCholeskyFactorisation.hxx
1 // Copyright (C) 2010 Marc DuruflĂ©
2 //
3 // This file is part of the linear-algebra library Seldon,
4 // http://seldon.sourceforge.net/.
5 //
6 // Seldon is free software; you can redistribute it and/or modify it under the
7 // terms of the GNU Lesser General Public License as published by the Free
8 // Software Foundation; either version 2.1 of the License, or (at your option)
9 // any later version.
10 //
11 // Seldon is distributed in the hope that it will be useful, but WITHOUT ANY
12 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
13 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
14 // more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public License
17 // along with Seldon. If not, see http://www.gnu.org/licenses/.
18 
19 
20 #ifndef SELDON_FILE_SPARSE_CHOLESKY_FACTORISATION_HXX
21 
22 
23 namespace Seldon
24 {
25 
27  template<class T>
29  {
30  protected :
40  int n;
47 
48  public :
49  // Available solvers.
50  enum {SELDON_SOLVER, CHOLMOD, PASTIX};
51 
53 
54  void HideMessages();
55  void ShowMessages();
56  void ShowFullHistory();
57 
58  void Clear();
59 
60  int GetM() const;
61  int GetN() const;
62  size_t GetMemorySize() const;
63 
64  int GetTypeOrdering() const;
65  void SetPermutation(const IVect&);
66  void SelectOrdering(int);
67 
68  void SelectDirectSolver(int);
69  int GetDirectSolver();
70 
71  void InitSolver();
72 
73  template<class Prop, class Storage, class Allocator>
75  bool keep_matrix = false);
76 
77  template<class T1>
78  void Solve(const SeldonTranspose& TransA, Vector<T1>& x, bool assemble = true);
79 
80  template<class T1>
81  void Mlt(const SeldonTranspose& TransA, Vector<T1>& x, bool assemble = true);
82 
83  };
84 
85 
86  template<class T, class Prop, class Allocator>
87  void GetCholesky(Matrix<T, Prop, ArrayRowSymSparse, Allocator>& A,
88  int print_level = 0);
89 
90  template<class T0, class Prop, class Allocator0,
91  class T1, class Storage, class Allocator1>
92  void SolveCholesky(const SeldonTranspose& TransA,
93  const Matrix<T0, Prop, ArrayRowSymSparse, Allocator0>& A,
94  Vector<T1, Storage, Allocator1>& x);
95 
96  template<class T0, class Prop, class Alloc0,
97  class T1, class Storage, class Allocator1>
98  void SolveCholesky(const SeldonTranspose& TransA,
99  const Matrix<T0, Prop, RowSymSparse, Alloc0>& A,
100  Vector<T1, Storage, Allocator1>& X);
101 
102  template<class T0, class Prop, class Allocator0,
103  class T1, class Storage, class Allocator1>
104  void MltCholesky(const SeldonTranspose& TransA,
105  const Matrix<T0, Prop, ArrayRowSymSparse, Allocator0>& A,
106  Vector<T1, Storage, Allocator1>& x);
107 
108  template<class T0, class Prop, class Allocator0,
109  class T1, class Storage, class Allocator1>
110  void MltCholesky(const SeldonTranspose& TransA,
111  const Matrix<T0, Prop, RowSymSparse, Allocator0>& A,
112  Vector<T1, Storage, Allocator1>& x);
113 
114 } // namespace Seldon.
115 
116 
117 #define SELDON_FILE_SPARSE_CHOLESKY_FACTORISATION_HXX
118 #endif
Seldon::SparseCholeskySolver::ShowMessages
void ShowMessages()
Displays only brief messages.
Definition: SparseCholeskyFactorisation.cxx:476
Seldon::SparseCholeskySolver::GetTypeOrdering
int GetTypeOrdering() const
Returns the type of ordering used.
Definition: SparseCholeskyFactorisation.cxx:540
Seldon::SparseCholeskySolver::mat_sym
Matrix< T, Symmetric, ArrayRowSymSparse > mat_sym
Cholesky factors.
Definition: SparseCholeskyFactorisation.hxx:42
Seldon::SeldonTranspose
Definition: MatrixFlag.hxx:32
Seldon::SparseCholeskySolver::GetM
int GetM() const
Returns the number of rows.
Definition: SparseCholeskyFactorisation.cxx:511
Seldon::SparseCholeskySolver::SetPermutation
void SetPermutation(const IVect &)
Modifies the ordering used.
Definition: SparseCholeskyFactorisation.cxx:548
Seldon::Vector< int, VectFull >
Seldon::SparseCholeskySolver::print_level
int print_level
Verbosity level.
Definition: SparseCholeskyFactorisation.hxx:32
Seldon::SparseCholeskySolver::Solve
void Solve(const SeldonTranspose &TransA, Vector< T1 > &x, bool assemble=true)
Solves L x = b or L^T x = b.
Definition: SparseCholeskyFactorisation.cxx:667
Seldon::Matrix< T, Symmetric, ArrayRowSymSparse >
Seldon::SparseCholeskySolver::SparseCholeskySolver
SparseCholeskySolver()
Default constructor.
Definition: SparseCholeskyFactorisation.cxx:445
Seldon::SparseCholeskySolver::permutation
IVect permutation
Permutation array.
Definition: SparseCholeskyFactorisation.hxx:36
Seldon::SparseCholeskySolver::solver
VirtualSparseDirectSolver< T > * solver
extern Cholesky solver
Definition: SparseCholeskyFactorisation.hxx:46
Seldon::SparseCholeskySolver::xtmp
Vector< T > xtmp
Temporary vector.
Definition: SparseCholeskyFactorisation.hxx:44
Seldon::SparseCholeskySolver
Class grouping different Cholesky solvers.
Definition: SparseCholeskyFactorisation.hxx:28
Seldon::SparseCholeskySolver::GetN
int GetN() const
Returns the number of rows.
Definition: SparseCholeskyFactorisation.cxx:519
Seldon::SparseCholeskySolver::Mlt
void Mlt(const SeldonTranspose &TransA, Vector< T1 > &x, bool assemble=true)
Computes L x or L^T.
Definition: SparseCholeskyFactorisation.cxx:718
Seldon::SparseCholeskySolver::SelectDirectSolver
void SelectDirectSolver(int)
Modifies the direct solver used.
Definition: SparseCholeskyFactorisation.cxx:565
Seldon::VirtualSparseDirectSolver
Base class for an interface with a direct solver.
Definition: SparseSolver.hxx:30
Seldon::SparseCholeskySolver::n
int n
Size of factorized linear system.
Definition: SparseCholeskyFactorisation.hxx:40
Seldon::SparseCholeskySolver::Clear
void Clear()
Clears Cholesky factors.
Definition: SparseCholeskyFactorisation.cxx:496
Seldon::SparseCholeskySolver::type_ordering
int type_ordering
Ordering to use.
Definition: SparseCholeskyFactorisation.hxx:34
Seldon::SparseCholeskySolver::type_solver
int type_solver
Solver to use.
Definition: SparseCholeskyFactorisation.hxx:38
Seldon::SparseCholeskySolver::HideMessages
void HideMessages()
Displays no messages.
Definition: SparseCholeskyFactorisation.cxx:466
Seldon
Seldon namespace.
Definition: Array.cxx:24
Seldon::SparseCholeskySolver::Factorize
void Factorize(Matrix< T, Prop, Storage, Allocator > &A, bool keep_matrix=false)
Performs Cholesky factorization.
Definition: SparseCholeskyFactorisation.cxx:616
Seldon::SparseCholeskySolver::GetMemorySize
size_t GetMemorySize() const
returns memory size used by the object in bytes
Definition: SparseCholeskyFactorisation.cxx:527
Seldon::SparseCholeskySolver::SelectOrdering
void SelectOrdering(int)
Modifies the type of ordering used.
Definition: SparseCholeskyFactorisation.cxx:557
Seldon::SparseCholeskySolver::GetDirectSolver
int GetDirectSolver()
returns the type of direct solver used
Definition: SparseCholeskyFactorisation.cxx:574
Seldon::SparseCholeskySolver::ShowFullHistory
void ShowFullHistory()
Displays a lot of messages.
Definition: SparseCholeskyFactorisation.cxx:486