16 #ifndef MOAB_HYPREPARMATRIX
17 #define MOAB_HYPREPARMATRIX
22 #ifdef MOAB_HAVE_EIGEN3
24 #include <Eigen/Sparse>
26 #error Configure with Eigen3 enabled
39 #include "_hypre_parcsr_mv.h"
40 #include "_hypre_parcsr_ls.h"
41 #include "interpreter.h"
42 #include "HYPRE_MatvecFunctions.h"
44 #include "temp_multivector.h"
47 #error "MOAB does not work with HYPRE's complex numbers support"
59 inline int to_int( HYPRE_Int i )
62 MFEM_ASSERT( HYPRE_Int(
int( i ) ) == i,
"overflow converting HYPRE_Int to int" );
75 hypre_ParCSRMatrix* A_parcsr;
77 mv_InterfaceInterpreter* interpreter;
95 friend class HypreSolver;
98 typedef Eigen::VectorXd Vector;
99 typedef Eigen::SparseMatrix< double, Eigen::RowMajor > MOABSparseMatrix;
107 HypreParMatrix( HYPRE_IJMatrix a )
111 height = GetNumRows();
112 width = GetNumCols();
118 HypreParMatrix(
moab::ParallelComm* p_comm, HYPRE_Int glob_size, HYPRE_Int* row_starts, HYPRE_Int nnz_pr_diag = 0 );
124 HYPRE_Int global_num_rows,
125 HYPRE_Int global_num_cols,
126 HYPRE_Int* row_starts,
127 HYPRE_Int* col_starts,
128 HYPRE_Int nnz_pr_diag = 0,
129 HYPRE_Int onz_pr_diag = 0,
130 HYPRE_Int nnz_pr_offdiag = 0 );
141 void MakeRef(
const HypreParMatrix& master );
150 operator HYPRE_IJMatrix()
155 operator HYPRE_ParCSRMatrix()
157 return (HYPRE_ParCSRMatrix)A_parcsr;
161 HYPRE_IJMatrix StealData();
163 void StealData( HypreParMatrix& A );
175 inline HYPRE_Int NNZ()
177 return A_parcsr->num_nonzeros;
180 inline HYPRE_Int* RowPart()
182 return A->row_partitioning;
185 inline HYPRE_Int* ColPart()
187 return A->col_partitioning;
192 return A->global_num_rows;
197 return A->global_num_cols;
200 inline HYPRE_Int FirstRow()
202 return A->global_first_row;
205 inline HYPRE_Int FirstCol()
207 return A->global_first_col;
211 void GetDiag( Vector& diag )
const;
213 void GetDiag( MOABSparseMatrix& diag )
const;
215 void GetOffd( MOABSparseMatrix& offd, HYPRE_Int*& cmap )
const;
229 void resize( HYPRE_Int glob_size, HYPRE_Int* row_starts, HYPRE_Int nnz_pr_diag = 0, HYPRE_Int nnz_pr_offdiag = 0 );
234 void resize( HYPRE_Int global_num_rows,
235 HYPRE_Int global_num_cols,
236 HYPRE_Int* row_starts,
237 HYPRE_Int* col_starts,
238 HYPRE_Int* nnz_pr_diag = NULL,
239 HYPRE_Int* onz_pr_diag = NULL,
240 HYPRE_Int nnz_pr_offdiag = 0 );
243 int GetNumRows()
const
245 return internal::to_int( hypre_CSRMatrixNumRows( hypre_ParCSRMatrixDiag( A_parcsr ) ) );
249 int GetNumCols()
const
251 return internal::to_int( hypre_CSRMatrixNumCols( hypre_ParCSRMatrixDiag( A_parcsr ) ) );
255 HYPRE_Int Mult( HypreParVector& x, HypreParVector& y,
double alpha = 1.0,
double beta = 0.0 );
257 HYPRE_Int Mult( HYPRE_ParVector x, HYPRE_ParVector y,
double alpha = 1.0,
double beta = 0.0 );
259 HYPRE_Int MultTranspose( HypreParVector& x, HypreParVector& y,
double alpha = 1.0,
double beta = 0.0 );
261 void Mult(
double a,
const HypreParVector& x,
double b, HypreParVector& y )
const;
262 void MultTranspose(
double a,
const HypreParVector& x,
double b, HypreParVector& y )
const;
287 void ScaleRows(
const Vector& s );
289 void InvScaleRows(
const Vector& s );
291 void operator*=(
double s );
294 void Threshold(
double threshold = 0.0 );
297 void EliminateZeroRows()
299 hypre_ParCSRMatrixFixZeroRows( A_parcsr );
304 void EliminateRowsCols(
const std::vector< HYPRE_Int >& rows_cols,
const HypreParVector& X, HypreParVector& B );
309 HypreParMatrix* EliminateRowsCols(
const std::vector< HYPRE_Int >& rows_cols );
311 HYPRE_Int GetValues(
const HYPRE_Int nrows,
315 HYPRE_Complex* values );
316 HYPRE_Int SetValues(
const HYPRE_Int nrows,
318 const HYPRE_Int* rows,
319 const HYPRE_Int* cols,
320 const HYPRE_Complex* values );
321 HYPRE_Int AddToValues(
const HYPRE_Int nrows,
323 const HYPRE_Int* rows,
324 const HYPRE_Int* cols,
325 const HYPRE_Complex* values );
327 HYPRE_Int FinalizeAssembly();
329 HYPRE_Int verbosity(
const HYPRE_Int level );
332 void Print(
const char* fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0 );
334 void Read(
const char* fname );
337 virtual ~HypreParMatrix()
343 friend HypreParMatrix* ParMult( HypreParMatrix* A, HypreParMatrix* B );
346 friend HypreParMatrix* RAP( HypreParMatrix* A, HypreParMatrix* P );
348 friend HypreParMatrix* RAP( HypreParMatrix* Rt, HypreParMatrix* A, HypreParMatrix* P );
353 friend void EliminateBC( HypreParMatrix& A,
355 const std::vector< int >& ess_dof_list,
356 const HypreParVector& X,
359 friend class HypreSolver;