ParGeMSLR
block_jacobi.hpp
Go to the documentation of this file.
1 #ifndef PARGEMSLR_BLOCK_JACOBI_H
2 #define PARGEMSLR_BLOCK_JACOBI_H
3 
8 #include <iostream>
9 
10 #include "../utils/utils.hpp"
11 #include "../vectors/int_vector.hpp"
12 #include "../vectors/sequential_vector.hpp"
13 #include "../vectors/parallel_vector.hpp"
14 #include "../matrices/csr_matrix.hpp"
15 #include "../matrices/parallel_csr_matrix.hpp"
16 #include "../solvers/solver.hpp"
17 
18 #ifdef PARGEMSLR_CUDA
19 #include "cusparse.h"
20 #endif
21 
22 namespace pargemslr
23 {
28  template <class ParallelMatrixType, class ParallelVectorType, class LocalMatrixType, class LocalVectorType, typename DataType>
29  class BlockJacobiClass: public SolverClass<ParallelMatrixType, ParallelVectorType, DataType>
30  {
31  private:
32 
38 
43  int _own_local_preconditioner;
44 
49  int _location;
50 
51  public:
57 
63  virtual int Clear();
64 
69  virtual ~BlockJacobiClass();
70 
78 
85 
94 
102 
110  virtual int Setup( ParallelVectorType &x, ParallelVectorType &rhs);
111 
119  virtual int Solve( ParallelVectorType &x, ParallelVectorType &rhs);
120 
126  virtual long int GetNumNonzeros();
127 
134 
141 
148 
154  int SetOwnLocalPreconditioner(bool own_local_preconditioner);
155 
162  virtual int SetSolveLocation( const int &location);
163 
170  virtual int MoveData( const int &location);
171 
172  };
173 
178 
179 }
180 
181 #endif
pargemslr::BlockJacobiClass::~BlockJacobiClass
virtual ~BlockJacobiClass()
The destructor of BlockJacobiClass.
Definition: block_jacobi.cpp:31
pargemslr::BlockJacobiClass::BlockJacobiClass
BlockJacobiClass()
The constructor of BlockJacobiClass.
Definition: block_jacobi.cpp:18
pargemslr::BlockJacobiClass::SetLocalPreconditioner
int SetLocalPreconditioner(SolverClass< LocalMatrixType, LocalVectorType, DataType > &local_precond)
Set the local preconditioenr.
Definition: block_jacobi.cpp:192
pargemslr::BlockJacobiClass::SetOwnLocalPreconditioner
int SetOwnLocalPreconditioner(bool own_local_preconditioner)
Set the own preconditioenr. If set to true, the preconditioner will be freed when destrooy.
Definition: block_jacobi.cpp:233
pargemslr::BlockJacobiClass::operator=
BlockJacobiClass< ParallelMatrixType, ParallelVectorType, LocalMatrixType, LocalVectorType, DataType > & operator=(const BlockJacobiClass< ParallelMatrixType, ParallelVectorType, LocalMatrixType, LocalVectorType, DataType > &solver)
The = operator of BlockJacobiClass. Note that this is not the true copy. We only copy the pointer to ...
Definition: block_jacobi.cpp:68
pargemslr::BlockJacobiClass::Solve
virtual int Solve(ParallelVectorType &x, ParallelVectorType &rhs)
Solve phase. Call this function after Setup. Solve with cusparse if unified memory/device memory is u...
Definition: block_jacobi.cpp:155
pargemslr::BlockJacobiClass::GetLocalPreconditionerP
SolverClass< LocalMatrixType, LocalVectorType, DataType > * GetLocalPreconditionerP()
Get the local preconditioenr.
Definition: block_jacobi.cpp:223
pargemslr::BlockJacobiClass::SetSolveLocation
virtual int SetSolveLocation(const int &location)
Set the data location that the preconditioner apply to.
Definition: block_jacobi.cpp:244
pargemslr::SolverClass
The base solver class.
Definition: solver.hpp:47
pargemslr::BlockJacobiClass::Setup
virtual int Setup(ParallelVectorType &x, ParallelVectorType &rhs)
Setup the precondioner phase. Will be called by the solver if not called directly.
Definition: block_jacobi.cpp:124
pargemslr::BlockJacobiClass::MoveData
virtual int MoveData(const int &location)
Move the preconditioner to another location. Only can be called after Setup.
Definition: block_jacobi.cpp:261
pargemslr::BlockJacobiClass
Block Jacobi preconditioner, only for parallel csr matrix.
Definition: block_jacobi.hpp:30
pargemslr::BlockJacobiClass::GetNumNonzeros
virtual long int GetNumNonzeros()
Get the total number of nonzeros.
Definition: block_jacobi.cpp:174
pargemslr::BlockJacobiClass::Clear
virtual int Clear()
Free the current precondioner.
Definition: block_jacobi.cpp:101
pargemslr::BlockJacobiClass::SetLocalPreconditionerP
int SetLocalPreconditionerP(SolverClass< LocalMatrixType, LocalVectorType, DataType > *local_precond)
Set the local preconditioenr.
Definition: block_jacobi.cpp:202