12 #include "../config/config.hpp"
17 #include "../fem/fem.hpp"
27 #if MFEM_HYPRE_VERSION < 21400
29 #define mfem_hypre_TAlloc(type, size) hypre_TAlloc(type, size)
30 #define mfem_hypre_CTAlloc(type, size) hypre_CTAlloc(type, size)
31 #define mfem_hypre_TFree(ptr) hypre_TFree(ptr)
33 #else // MFEM_HYPRE_VERSION >= 21400
35 #define mfem_hypre_TAlloc(type, size) \
36 hypre_TAlloc(type, size, HYPRE_MEMORY_HOST)
37 #define mfem_hypre_CTAlloc(type, size) \
38 hypre_CTAlloc(type, size, HYPRE_MEMORY_HOST)
39 #define mfem_hypre_TFree(ptr) hypre_TFree(ptr, HYPRE_MEMORY_HOST)
71 #endif // #if MFEM_HYPRE_VERSION < 21400
78 template<
typename TargetT,
typename SourceT>
79 static TargetT *DuplicateAs(
const SourceT *array,
int size,
80 bool cplusplus =
true)
82 TargetT *target_array = cplusplus ?
new TargetT[size]
83 : mfem_hypre_TAlloc(TargetT, size);
84 for (
int i = 0; i < size; i++)
86 target_array[i] = array[i];
91 inline void HypreParVector::_SetDataAndSize_()
93 SetDataAndSize(hypre_VectorData(hypre_ParVectorLocalVector(x)),
95 hypre_VectorSize(hypre_ParVectorLocalVector(x))));
98 HypreParVector::HypreParVector(MPI_Comm comm, HYPRE_Int glob_size,
101 x = hypre_ParVectorCreate(comm,glob_size,col);
102 hypre_ParVectorInitialize(x);
103 hypre_ParVectorSetPartitioningOwner(x,0);
105 hypre_ParVectorSetDataOwner(x,1);
106 hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),1);
112 double *_data, HYPRE_Int *col) :
Vector()
114 x = hypre_ParVectorCreate(comm,glob_size,col);
115 hypre_ParVectorSetDataOwner(x,1);
116 hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),0);
117 hypre_ParVectorSetPartitioningOwner(x,0);
119 hypre_VectorData(hypre_ParVectorLocalVector(x)) = &tmp;
122 hypre_ParVectorInitialize(x);
124 hypre_VectorData(hypre_ParVectorLocalVector(x)) = _data;
131 x = hypre_ParVectorCreate(y.x -> comm, y.x -> global_size,
132 y.x -> partitioning);
133 hypre_ParVectorInitialize(x);
134 hypre_ParVectorSetPartitioningOwner(x,0);
135 hypre_ParVectorSetDataOwner(x,1);
136 hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),1);
146 x = hypre_ParVectorInDomainOf(const_cast<HypreParMatrix&>(A));
150 x = hypre_ParVectorInRangeOf(const_cast<HypreParMatrix&>(A));
158 x = (hypre_ParVector *) y;
167 hypre_ParVectorInitialize(x);
168 hypre_ParVectorSetPartitioningOwner(x,0);
170 hypre_ParVectorSetDataOwner(x,1);
171 hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),1);
178 hypre_Vector *hv = hypre_ParVectorToVectorAll(*
this);
181 hypre_SeqVectorSetDataOwner(hv,0);
182 hypre_SeqVectorDestroy(hv);
188 hypre_ParVectorSetConstantValues(x,d);
195 if (size != y.
Size())
201 for (
int i = 0; i <
size; i++)
210 Vector::data = hypre_VectorData(hypre_ParVectorLocalVector(x)) = _data;
215 return hypre_ParVectorSetRandomValues(x,seed);
220 hypre_ParVectorPrint(x,fname);
227 hypre_ParVectorDestroy(x);
231 #ifdef MFEM_USE_SUNDIALS
234 #define SUNFALSE FALSE
239 MFEM_ASSERT(nv && N_VGetVectorID(nv) == SUNDIALS_NVEC_PARHYP,
241 N_VectorContent_ParHyp nv_c = (N_VectorContent_ParHyp)(nv->content);
242 MFEM_ASSERT(nv_c->own_parvector == SUNFALSE,
"invalid N_Vector");
243 nv_c->local_length = x->local_vector->size;
244 nv_c->global_length = x->global_size;
245 nv_c->comm = x->comm;
249 #endif // MFEM_USE_SUNDIALS
254 return hypre_ParVectorInnerProd(*x, *y);
259 return hypre_ParVectorInnerProd(x, y);
268 double loc_norm = vec.
Norml1();
269 MPI_Allreduce(&loc_norm, &norm, 1, MPI_DOUBLE, MPI_SUM, comm);
273 double loc_norm = vec*vec;
274 MPI_Allreduce(&loc_norm, &norm, 1, MPI_DOUBLE, MPI_SUM, comm);
280 for (
int i = 0; i < vec.
Size(); i++)
282 sum += pow(fabs(vec(i)), p);
284 MPI_Allreduce(&sum, &norm, 1, MPI_DOUBLE, MPI_SUM, comm);
285 norm = pow(norm, 1.0/p);
290 MPI_Allreduce(&loc_norm, &norm, 1, MPI_DOUBLE, MPI_MAX, comm);
296 void HypreParMatrix::Init()
300 diagOwner = offdOwner = colMapOwner = -1;
310 char HypreParMatrix::CopyCSR(
SparseMatrix *csr, hypre_CSRMatrix *hypre_csr)
312 hypre_CSRMatrixData(hypre_csr) = csr->
GetData();
314 hypre_CSRMatrixI(hypre_csr) = csr->
GetI();
315 hypre_CSRMatrixJ(hypre_csr) = csr->
GetJ();
319 hypre_CSRMatrixI(hypre_csr) =
320 DuplicateAs<HYPRE_Int>(csr->
GetI(), csr->
Height()+1);
321 hypre_CSRMatrixJ(hypre_csr) =
328 char HypreParMatrix::CopyBoolCSR(Table *bool_csr, hypre_CSRMatrix *hypre_csr)
330 int nnz = bool_csr->Size_of_connections();
331 double *data =
new double[nnz];
332 for (
int i = 0; i < nnz; i++)
336 hypre_CSRMatrixData(hypre_csr) = data;
338 hypre_CSRMatrixI(hypre_csr) = bool_csr->GetI();
339 hypre_CSRMatrixJ(hypre_csr) = bool_csr->GetJ();
343 hypre_CSRMatrixI(hypre_csr) =
344 DuplicateAs<HYPRE_Int>(bool_csr->GetI(), bool_csr->Size()+1);
345 hypre_CSRMatrixJ(hypre_csr) =
346 DuplicateAs<HYPRE_Int>(bool_csr->GetJ(), nnz);
352 void HypreParMatrix::CopyCSR_J(hypre_CSRMatrix *hypre_csr,
int *J)
354 HYPRE_Int nnz = hypre_CSRMatrixNumNonzeros(hypre_csr);
355 for (HYPRE_Int j = 0; j < nnz; j++)
357 J[j] = int(hypre_CSRMatrixJ(hypre_csr)[j]);
364 :
Operator(diag->Height(), diag->Width())
367 A = hypre_ParCSRMatrixCreate(comm, glob_size, glob_size, row_starts,
369 hypre_ParCSRMatrixSetDataOwner(A,1);
370 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
371 hypre_ParCSRMatrixSetColStartsOwner(A,0);
373 hypre_CSRMatrixSetDataOwner(A->diag,0);
374 diagOwner = CopyCSR(diag, A->diag);
375 hypre_CSRMatrixSetRownnz(A->diag);
377 hypre_CSRMatrixSetDataOwner(A->offd,1);
378 hypre_CSRMatrixI(A->offd) = mfem_hypre_CTAlloc(HYPRE_Int, diag->
Height()+1);
385 hypre_ParCSRMatrixSetNumNonzeros(A);
388 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
390 CopyCSR_J(A->diag, diag->
GetJ());
393 hypre_MatvecCommPkgCreate(A);
398 HYPRE_Int global_num_rows,
399 HYPRE_Int global_num_cols,
400 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
402 :
Operator(diag->Height(), diag->Width())
405 A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
406 row_starts, col_starts,
408 hypre_ParCSRMatrixSetDataOwner(A,1);
409 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
410 hypre_ParCSRMatrixSetColStartsOwner(A,0);
412 hypre_CSRMatrixSetDataOwner(A->diag,0);
413 diagOwner = CopyCSR(diag, A->diag);
414 hypre_CSRMatrixSetRownnz(A->diag);
416 hypre_CSRMatrixSetDataOwner(A->offd,1);
417 hypre_CSRMatrixI(A->offd) = mfem_hypre_CTAlloc(HYPRE_Int, diag->
Height()+1);
419 hypre_ParCSRMatrixSetNumNonzeros(A);
422 if (row_starts == col_starts)
424 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
426 CopyCSR_J(A->diag, diag->
GetJ());
430 hypre_MatvecCommPkgCreate(A);
435 HYPRE_Int global_num_rows,
436 HYPRE_Int global_num_cols,
437 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
440 :
Operator(diag->Height(), diag->Width())
443 A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
444 row_starts, col_starts,
447 hypre_ParCSRMatrixSetDataOwner(A,1);
448 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
449 hypre_ParCSRMatrixSetColStartsOwner(A,0);
451 hypre_CSRMatrixSetDataOwner(A->diag,0);
452 diagOwner = CopyCSR(diag, A->diag);
453 hypre_CSRMatrixSetRownnz(A->diag);
455 hypre_CSRMatrixSetDataOwner(A->offd,0);
456 offdOwner = CopyCSR(offd, A->offd);
457 hypre_CSRMatrixSetRownnz(A->offd);
459 hypre_ParCSRMatrixColMapOffd(A) = cmap;
463 hypre_ParCSRMatrixSetNumNonzeros(A);
466 if (row_starts == col_starts)
468 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
470 CopyCSR_J(A->diag, diag->
GetJ());
474 hypre_MatvecCommPkgCreate(A);
480 HYPRE_Int global_num_rows, HYPRE_Int global_num_cols,
481 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
482 HYPRE_Int *diag_i, HYPRE_Int *diag_j,
double *diag_data,
483 HYPRE_Int *offd_i, HYPRE_Int *offd_j,
double *offd_data,
484 HYPRE_Int offd_num_cols, HYPRE_Int *offd_col_map)
487 A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
488 row_starts, col_starts, offd_num_cols, 0, 0);
489 hypre_ParCSRMatrixSetDataOwner(A,1);
490 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
491 hypre_ParCSRMatrixSetColStartsOwner(A,0);
493 HYPRE_Int local_num_rows = hypre_CSRMatrixNumRows(A->diag);
495 hypre_CSRMatrixSetDataOwner(A->diag,0);
496 hypre_CSRMatrixI(A->diag) = diag_i;
497 hypre_CSRMatrixJ(A->diag) = diag_j;
498 hypre_CSRMatrixData(A->diag) = diag_data;
499 hypre_CSRMatrixNumNonzeros(A->diag) = diag_i[local_num_rows];
500 hypre_CSRMatrixSetRownnz(A->diag);
504 hypre_CSRMatrixSetDataOwner(A->offd,0);
505 hypre_CSRMatrixI(A->offd) = offd_i;
506 hypre_CSRMatrixJ(A->offd) = offd_j;
507 hypre_CSRMatrixData(A->offd) = offd_data;
508 hypre_CSRMatrixNumNonzeros(A->offd) = offd_i[local_num_rows];
509 hypre_CSRMatrixSetRownnz(A->offd);
513 hypre_ParCSRMatrixColMapOffd(A) = offd_col_map;
517 hypre_ParCSRMatrixSetNumNonzeros(A);
520 if (row_starts == col_starts)
522 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
525 hypre_MatvecCommPkgCreate(A);
533 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
536 MFEM_ASSERT(sm_a != NULL,
"invalid input");
537 MFEM_VERIFY(!HYPRE_AssumedPartitionCheck(),
538 "this method can not be used with assumed partition");
542 hypre_CSRMatrix *csr_a;
543 csr_a = hypre_CSRMatrixCreate(sm_a ->
Height(), sm_a ->
Width(),
544 sm_a -> NumNonZeroElems());
546 hypre_CSRMatrixSetDataOwner(csr_a,0);
547 CopyCSR(sm_a, csr_a);
548 hypre_CSRMatrixSetRownnz(csr_a);
550 A = hypre_CSRMatrixToParCSRMatrix(comm, csr_a, row_starts, col_starts);
553 delete [] hypre_CSRMatrixI(csr_a);
554 delete [] hypre_CSRMatrixJ(csr_a);
556 hypre_CSRMatrixI(csr_a) = NULL;
557 hypre_CSRMatrixDestroy(csr_a);
563 if (row_starts == col_starts)
565 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
568 hypre_MatvecCommPkgCreate(A);
573 HYPRE_Int global_num_rows,
574 HYPRE_Int global_num_cols,
575 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
580 A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
581 row_starts, col_starts, 0, nnz, 0);
582 hypre_ParCSRMatrixSetDataOwner(A,1);
583 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
584 hypre_ParCSRMatrixSetColStartsOwner(A,0);
586 hypre_CSRMatrixSetDataOwner(A->diag,0);
587 diagOwner = CopyBoolCSR(diag, A->diag);
588 hypre_CSRMatrixSetRownnz(A->diag);
590 hypre_CSRMatrixSetDataOwner(A->offd,1);
591 hypre_CSRMatrixI(A->offd) = mfem_hypre_CTAlloc(HYPRE_Int, diag->
Size()+1);
593 hypre_ParCSRMatrixSetNumNonzeros(A);
596 if (row_starts == col_starts)
598 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
600 CopyCSR_J(A->diag, diag->
GetJ());
604 hypre_MatvecCommPkgCreate(A);
613 HYPRE_Int *row, HYPRE_Int *col,
614 HYPRE_Int *i_diag, HYPRE_Int *j_diag,
615 HYPRE_Int *i_offd, HYPRE_Int *j_offd,
616 HYPRE_Int *cmap, HYPRE_Int cmap_size)
618 HYPRE_Int diag_nnz, offd_nnz;
621 if (HYPRE_AssumedPartitionCheck())
623 diag_nnz = i_diag[row[1]-row[0]];
624 offd_nnz = i_offd[row[1]-row[0]];
626 A = hypre_ParCSRMatrixCreate(comm, row[2], col[2], row, col,
627 cmap_size, diag_nnz, offd_nnz);
631 diag_nnz = i_diag[row[
id+1]-row[id]];
632 offd_nnz = i_offd[row[
id+1]-row[id]];
634 A = hypre_ParCSRMatrixCreate(comm, row[np], col[np], row, col,
635 cmap_size, diag_nnz, offd_nnz);
638 hypre_ParCSRMatrixSetDataOwner(A,1);
639 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
640 hypre_ParCSRMatrixSetColStartsOwner(A,0);
644 double *a_diag =
new double[diag_nnz];
645 for (i = 0; i < diag_nnz; i++)
650 double *a_offd =
new double[offd_nnz];
651 for (i = 0; i < offd_nnz; i++)
656 hypre_CSRMatrixSetDataOwner(A->diag,0);
657 hypre_CSRMatrixI(A->diag) = i_diag;
658 hypre_CSRMatrixJ(A->diag) = j_diag;
659 hypre_CSRMatrixData(A->diag) = a_diag;
660 hypre_CSRMatrixSetRownnz(A->diag);
664 hypre_CSRMatrixSetDataOwner(A->offd,0);
665 hypre_CSRMatrixI(A->offd) = i_offd;
666 hypre_CSRMatrixJ(A->offd) = j_offd;
667 hypre_CSRMatrixData(A->offd) = a_offd;
668 hypre_CSRMatrixSetRownnz(A->offd);
672 hypre_ParCSRMatrixColMapOffd(A) = cmap;
676 hypre_ParCSRMatrixSetNumNonzeros(A);
681 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
684 hypre_MatvecCommPkgCreate(A);
694 HYPRE_Int glob_ncols,
int *I, HYPRE_Int *J,
695 double *data, HYPRE_Int *rows, HYPRE_Int *cols)
701 HYPRE_Int my_col_start, my_col_end;
702 if (HYPRE_AssumedPartitionCheck())
705 my_col_start = cols[0];
706 my_col_end = cols[1];
711 MPI_Comm_rank(comm, &myid);
712 MPI_Comm_size(comm, &part_size);
714 my_col_start = cols[myid];
715 my_col_end = cols[myid+1];
719 HYPRE_Int *row_starts, *col_starts;
722 row_starts = col_starts = mfem_hypre_TAlloc(HYPRE_Int, part_size);
723 for (
int i = 0; i < part_size; i++)
725 row_starts[i] = rows[i];
730 row_starts = mfem_hypre_TAlloc(HYPRE_Int, part_size);
731 col_starts = mfem_hypre_TAlloc(HYPRE_Int, part_size);
732 for (
int i = 0; i < part_size; i++)
734 row_starts[i] = rows[i];
735 col_starts[i] = cols[i];
741 HYPRE_Int diag_nnz = 0, offd_nnz = 0, offd_num_cols = 0;
742 map<HYPRE_Int, HYPRE_Int> offd_map;
743 for (HYPRE_Int j = 0, loc_nnz = I[nrows]; j < loc_nnz; j++)
745 HYPRE_Int glob_col = J[j];
746 if (my_col_start <= glob_col && glob_col < my_col_end)
752 offd_map.insert(pair<const HYPRE_Int, HYPRE_Int>(glob_col, -1));
757 for (map<HYPRE_Int, HYPRE_Int>::iterator it = offd_map.begin();
758 it != offd_map.end(); ++it)
760 it->second = offd_num_cols++;
764 A = hypre_ParCSRMatrixCreate(comm, glob_nrows, glob_ncols,
765 row_starts, col_starts, offd_num_cols,
767 hypre_ParCSRMatrixInitialize(A);
769 HYPRE_Int *diag_i, *diag_j, *offd_i, *offd_j, *offd_col_map;
770 double *diag_data, *offd_data;
773 diag_data = A->diag->data;
776 offd_data = A->offd->data;
777 offd_col_map = A->col_map_offd;
779 diag_nnz = offd_nnz = 0;
780 for (HYPRE_Int i = 0, j = 0; i < nrows; i++)
782 diag_i[i] = diag_nnz;
783 offd_i[i] = offd_nnz;
784 for (HYPRE_Int j_end = I[i+1]; j < j_end; j++)
786 HYPRE_Int glob_col = J[j];
787 if (my_col_start <= glob_col && glob_col < my_col_end)
789 diag_j[diag_nnz] = glob_col - my_col_start;
790 diag_data[diag_nnz] = data[j];
795 offd_j[offd_nnz] = offd_map[glob_col];
796 offd_data[offd_nnz] = data[j];
801 diag_i[nrows] = diag_nnz;
802 offd_i[nrows] = offd_nnz;
803 for (map<HYPRE_Int, HYPRE_Int>::iterator it = offd_map.begin();
804 it != offd_map.end(); ++it)
806 offd_col_map[it->second] = it->first;
809 hypre_ParCSRMatrixSetNumNonzeros(A);
811 if (row_starts == col_starts)
813 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
815 hypre_MatvecCommPkgCreate(A);
836 MFEM_ASSERT(diagOwner == -1 && offdOwner == -1 && colMapOwner == -1,
"");
837 MFEM_ASSERT(ParCSROwner,
"");
838 hypre_ParCSRMatrix *R = A;
847 if (!A || hypre_ParCSRMatrixOwnsRowStarts(A) ||
848 (hypre_ParCSRMatrixRowStarts(A) == hypre_ParCSRMatrixColStarts(A) &&
849 hypre_ParCSRMatrixOwnsColStarts(A)))
855 if (HYPRE_AssumedPartitionCheck())
861 MPI_Comm_size(hypre_ParCSRMatrixComm(A), &row_starts_size);
865 HYPRE_Int *old_row_starts = hypre_ParCSRMatrixRowStarts(A);
866 HYPRE_Int *new_row_starts = mfem_hypre_CTAlloc(HYPRE_Int, row_starts_size);
867 for (
int i = 0; i < row_starts_size; i++)
869 new_row_starts[i] = old_row_starts[i];
872 hypre_ParCSRMatrixRowStarts(A) = new_row_starts;
873 hypre_ParCSRMatrixOwnsRowStarts(A) = 1;
875 if (hypre_ParCSRMatrixColStarts(A) == old_row_starts)
877 hypre_ParCSRMatrixColStarts(A) = new_row_starts;
878 hypre_ParCSRMatrixOwnsColStarts(A) = 0;
884 if (!A || hypre_ParCSRMatrixOwnsColStarts(A) ||
885 (hypre_ParCSRMatrixRowStarts(A) == hypre_ParCSRMatrixColStarts(A) &&
886 hypre_ParCSRMatrixOwnsRowStarts(A)))
892 if (HYPRE_AssumedPartitionCheck())
898 MPI_Comm_size(hypre_ParCSRMatrixComm(A), &col_starts_size);
902 HYPRE_Int *old_col_starts = hypre_ParCSRMatrixColStarts(A);
903 HYPRE_Int *new_col_starts = mfem_hypre_CTAlloc(HYPRE_Int, col_starts_size);
904 for (
int i = 0; i < col_starts_size; i++)
906 new_col_starts[i] = old_col_starts[i];
909 hypre_ParCSRMatrixColStarts(A) = new_col_starts;
911 if (hypre_ParCSRMatrixRowStarts(A) == old_col_starts)
913 hypre_ParCSRMatrixRowStarts(A) = new_col_starts;
914 hypre_ParCSRMatrixOwnsRowStarts(A) = 1;
915 hypre_ParCSRMatrixOwnsColStarts(A) = 0;
919 hypre_ParCSRMatrixOwnsColStarts(A) = 1;
927 for (
int j = 0; j < size; j++)
929 diag(j) = A->diag->data[A->diag->i[j]];
930 MFEM_ASSERT(A->diag->j[A->diag->i[j]] == j,
931 "the first entry in each row must be the diagonal one");
935 static void MakeWrapper(
const hypre_CSRMatrix *mat,
SparseMatrix &wrapper)
937 HYPRE_Int nr = hypre_CSRMatrixNumRows(mat);
938 HYPRE_Int nc = hypre_CSRMatrixNumCols(mat);
941 hypre_CSRMatrixJ(mat),
942 hypre_CSRMatrixData(mat),
943 nr, nc,
false,
false,
false);
945 HYPRE_Int nnz = hypre_CSRMatrixNumNonzeros(mat);
946 SparseMatrix tmp(DuplicateAs<int>(hypre_CSRMatrixI(mat), nr+1),
947 DuplicateAs<int>(hypre_CSRMatrixJ(mat), nnz),
948 hypre_CSRMatrixData(mat),
949 nr, nc,
true,
false,
false);
956 MakeWrapper(A->diag, diag);
961 MakeWrapper(A->offd, offd);
962 cmap = A->col_map_offd;
966 bool interleaved_rows,
967 bool interleaved_cols)
const
972 hypre_ParCSRMatrix **hypre_blocks =
new hypre_ParCSRMatrix*[nr * nc];
973 internal::hypre_ParCSRMatrixSplit(A, nr, nc, hypre_blocks,
974 interleaved_rows, interleaved_cols);
976 for (
int i = 0; i < nr; i++)
978 for (
int j = 0; j < nc; j++)
984 delete [] hypre_blocks;
989 hypre_ParCSRMatrix * At;
990 hypre_ParCSRMatrixTranspose(A, &At, 1);
991 hypre_ParCSRMatrixSetNumNonzeros(At);
993 hypre_MatvecCommPkgCreate(At);
999 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(At));
1008 return hypre_ParCSRMatrixMatvec(a, A, x, b, y);
1013 MFEM_ASSERT(x.
Size() ==
Width(),
"invalid x.Size() = " << x.
Size()
1014 <<
", expected size = " <<
Width());
1015 MFEM_ASSERT(y.
Size() ==
Height(),
"invalid y.Size() = " << y.
Size()
1016 <<
", expected size = " <<
Height());
1035 hypre_ParCSRMatrixMatvec(a, A, *X, b, *Y);
1039 double b,
Vector &y)
const
1041 MFEM_ASSERT(x.
Size() ==
Height(),
"invalid x.Size() = " << x.
Size()
1042 <<
", expected size = " <<
Height());
1043 MFEM_ASSERT(y.
Size() ==
Width(),
"invalid y.Size() = " << y.
Size()
1044 <<
", expected size = " <<
Width());
1065 hypre_ParCSRMatrixMatvecT(a, A, *Y, b, *X);
1071 return hypre_ParCSRMatrixMatvec(a, A, (hypre_ParVector *) x, b,
1072 (hypre_ParVector *) y);
1078 return hypre_ParCSRMatrixMatvecT(a, A, x, b, y);
1082 HYPRE_Int* row_starts)
const
1084 const bool assumed_partition = HYPRE_AssumedPartitionCheck();
1085 const bool row_starts_given = (row_starts != NULL);
1086 if (!row_starts_given)
1088 row_starts = hypre_ParCSRMatrixRowStarts(A);
1089 MFEM_VERIFY(D.
Height() == hypre_CSRMatrixNumRows(A->diag),
1090 "the matrix D is NOT compatible with the row starts of"
1091 " this HypreParMatrix, row_starts must be given.");
1096 if (assumed_partition)
1102 MPI_Comm_rank(
GetComm(), &offset);
1104 int local_num_rows = row_starts[offset+1]-row_starts[offset];
1105 MFEM_VERIFY(local_num_rows == D.
Height(),
"the number of rows in D is "
1106 " not compatible with the given row_starts");
1112 HYPRE_Int global_num_rows;
1113 if (assumed_partition)
1116 if (row_starts_given)
1118 global_num_rows = row_starts[2];
1125 global_num_rows = hypre_ParCSRMatrixGlobalNumRows(A);
1130 MPI_Comm_size(
GetComm(), &part_size);
1131 global_num_rows = row_starts[part_size];
1135 HYPRE_Int *col_starts = hypre_ParCSRMatrixColStarts(A);
1136 HYPRE_Int *col_map_offd;
1141 GetOffd(A_offd, col_map_offd);
1149 global_num_rows, hypre_ParCSRMatrixGlobalNumCols(A),
1150 DuplicateAs<HYPRE_Int>(row_starts, part_size,
false),
1151 DuplicateAs<HYPRE_Int>(col_starts, part_size,
false),
1153 DuplicateAs<HYPRE_Int>(col_map_offd, A_offd.
Width()));
1158 #ifndef HYPRE_BIGINT
1169 hypre_ParCSRMatrixSetRowStartsOwner(DA->A, 1);
1170 hypre_ParCSRMatrixSetColStartsOwner(DA->A, 1);
1172 DA->diagOwner = DA->offdOwner = 3;
1173 DA->colMapOwner = 1;
1180 if (hypre_CSRMatrixNumRows(A->diag) != hypre_CSRMatrixNumRows(A->offd))
1185 if (hypre_CSRMatrixNumRows(A->diag) != diag.
Size())
1187 mfem_error(
"Note the Vector diag is not of compatible dimensions with A\n");
1191 double *Adiag_data = hypre_CSRMatrixData(A->diag);
1192 HYPRE_Int *Adiag_i = hypre_CSRMatrixI(A->diag);
1195 double *Aoffd_data = hypre_CSRMatrixData(A->offd);
1196 HYPRE_Int *Aoffd_i = hypre_CSRMatrixI(A->offd);
1199 for (
int i(0); i < size; ++i)
1202 for (jj = Adiag_i[i]; jj < Adiag_i[i+1]; ++jj)
1204 Adiag_data[jj] *= val;
1206 for (jj = Aoffd_i[i]; jj < Aoffd_i[i+1]; ++jj)
1208 Aoffd_data[jj] *= val;
1215 if (hypre_CSRMatrixNumRows(A->diag) != hypre_CSRMatrixNumRows(A->offd))
1220 if (hypre_CSRMatrixNumRows(A->diag) != diag.
Size())
1222 mfem_error(
"Note the Vector diag is not of compatible dimensions with A\n");
1226 double *Adiag_data = hypre_CSRMatrixData(A->diag);
1227 HYPRE_Int *Adiag_i = hypre_CSRMatrixI(A->diag);
1230 double *Aoffd_data = hypre_CSRMatrixData(A->offd);
1231 HYPRE_Int *Aoffd_i = hypre_CSRMatrixI(A->offd);
1234 for (
int i(0); i < size; ++i)
1239 mfem_error(
"HypreParMatrix::InvDiagScale : Division by 0");
1243 for (jj = Adiag_i[i]; jj < Adiag_i[i+1]; ++jj)
1245 Adiag_data[jj] *= val;
1247 for (jj = Aoffd_i[i]; jj < Aoffd_i[i+1]; ++jj)
1249 Aoffd_data[jj] *= val;
1256 if (hypre_CSRMatrixNumRows(A->diag) != hypre_CSRMatrixNumRows(A->offd))
1261 HYPRE_Int size=hypre_CSRMatrixNumRows(A->diag);
1264 double *Adiag_data = hypre_CSRMatrixData(A->diag);
1265 HYPRE_Int *Adiag_i = hypre_CSRMatrixI(A->diag);
1266 for (jj = 0; jj < Adiag_i[size]; ++jj)
1268 Adiag_data[jj] *= s;
1271 double *Aoffd_data = hypre_CSRMatrixData(A->offd);
1272 HYPRE_Int *Aoffd_i = hypre_CSRMatrixI(A->offd);
1273 for (jj = 0; jj < Aoffd_i[size]; ++jj)
1275 Aoffd_data[jj] *= s;
1279 static void get_sorted_rows_cols(
const Array<int> &rows_cols,
1284 for (
int i = 0; i < rows_cols.
Size(); i++)
1286 hypre_sorted[i] = rows_cols[i];
1287 if (i && rows_cols[i-1] > rows_cols[i]) { sorted =
false; }
1289 if (!sorted) { hypre_sorted.
Sort(); }
1297 hypre_CSRMatrix * csr_A;
1298 hypre_CSRMatrix * csr_A_wo_z;
1299 hypre_ParCSRMatrix * parcsr_A_ptr;
1300 HYPRE_Int * row_starts = NULL; HYPRE_Int * col_starts = NULL;
1301 HYPRE_Int row_start = -1; HYPRE_Int row_end = -1;
1302 HYPRE_Int col_start = -1; HYPRE_Int col_end = -1;
1304 comm = hypre_ParCSRMatrixComm(A);
1306 ierr += hypre_ParCSRMatrixGetLocalRange(A,
1307 &row_start,&row_end,
1308 &col_start,&col_end );
1310 row_starts = hypre_ParCSRMatrixRowStarts(A);
1311 col_starts = hypre_ParCSRMatrixColStarts(A);
1313 bool old_owns_row = hypre_ParCSRMatrixOwnsRowStarts(A);
1314 bool old_owns_col = hypre_ParCSRMatrixOwnsColStarts(A);
1315 HYPRE_Int global_num_rows = hypre_ParCSRMatrixGlobalNumRows(A);
1316 HYPRE_Int global_num_cols = hypre_ParCSRMatrixGlobalNumCols(A);
1317 parcsr_A_ptr = hypre_ParCSRMatrixCreate(comm, global_num_rows,
1319 row_starts, col_starts,
1321 hypre_ParCSRMatrixOwnsRowStarts(parcsr_A_ptr) = old_owns_row;
1322 hypre_ParCSRMatrixOwnsColStarts(parcsr_A_ptr) = old_owns_col;
1323 hypre_ParCSRMatrixOwnsRowStarts(A) = 0;
1324 hypre_ParCSRMatrixOwnsColStarts(A) = 0;
1326 csr_A = hypre_MergeDiagAndOffd(A);
1332 csr_A_wo_z = hypre_CSRMatrixDeleteZeros(csr_A,threshold);
1336 if (csr_A_wo_z == NULL)
1342 ierr += hypre_CSRMatrixDestroy(csr_A);
1348 ierr += GenerateDiagAndOffd(csr_A_wo_z,parcsr_A_ptr,
1351 ierr += hypre_CSRMatrixDestroy(csr_A_wo_z);
1353 MFEM_VERIFY(ierr == 0,
"");
1357 hypre_ParCSRMatrixSetNumNonzeros(A);
1359 if (row_starts == col_starts)
1361 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
1363 hypre_MatvecCommPkgCreate(A);
1373 get_sorted_rows_cols(rows_cols, rc_sorted);
1375 internal::hypre_ParCSRMatrixEliminateAXB(
1382 get_sorted_rows_cols(rows_cols, rc_sorted);
1384 hypre_ParCSRMatrix* Ae;
1385 internal::hypre_ParCSRMatrixEliminateAAe(
1393 hypre_ParCSRMatrixPrintIJ(A,offi,offj,fname);
1401 HYPRE_Int base_i, base_j;
1402 hypre_ParCSRMatrixReadIJ(comm, fname, &base_i, &base_j, &A);
1403 hypre_ParCSRMatrixSetNumNonzeros(A);
1405 hypre_MatvecCommPkgCreate(A);
1416 HYPRE_IJMatrix A_ij;
1417 HYPRE_IJMatrixRead(fname, comm, 5555, &A_ij);
1419 HYPRE_ParCSRMatrix A_parcsr;
1420 HYPRE_IJMatrixGetObject(A_ij, (
void**) &A_parcsr);
1422 A = (hypre_ParCSRMatrix*)A_parcsr;
1424 hypre_ParCSRMatrixSetNumNonzeros(A);
1426 hypre_MatvecCommPkgCreate(A);
1434 hypre_ParCSRCommPkg *comm_pkg = A->comm_pkg;
1435 MPI_Comm comm = A->comm;
1437 const int tag = 46801;
1439 MPI_Comm_rank(comm, &myid);
1440 MPI_Comm_size(comm, &nproc);
1444 MPI_Recv(&c, 1, MPI_CHAR, myid-1, tag, comm, MPI_STATUS_IGNORE);
1448 out <<
"\nHypreParMatrix: hypre_ParCSRCommPkg:\n";
1450 out <<
"Rank " << myid <<
":\n"
1451 " number of sends = " << comm_pkg->num_sends <<
1452 " (" <<
sizeof(double)*comm_pkg->send_map_starts[comm_pkg->num_sends] <<
1454 " number of recvs = " << comm_pkg->num_recvs <<
1455 " (" <<
sizeof(
double)*comm_pkg->recv_vec_starts[comm_pkg->num_recvs] <<
1457 if (myid != nproc-1)
1460 MPI_Send(&c, 1, MPI_CHAR, myid+1, tag, comm);
1469 void HypreParMatrix::Destroy()
1471 if ( X != NULL ) {
delete X; }
1472 if ( Y != NULL ) {
delete Y; }
1474 if (A == NULL) {
return; }
1480 delete [] hypre_CSRMatrixI(A->diag);
1481 delete [] hypre_CSRMatrixJ(A->diag);
1483 hypre_CSRMatrixI(A->diag) = NULL;
1484 hypre_CSRMatrixJ(A->diag) = NULL;
1487 delete [] hypre_CSRMatrixData(A->diag);
1489 hypre_CSRMatrixData(A->diag) = NULL;
1495 delete [] hypre_CSRMatrixI(A->offd);
1496 delete [] hypre_CSRMatrixJ(A->offd);
1498 hypre_CSRMatrixI(A->offd) = NULL;
1499 hypre_CSRMatrixJ(A->offd) = NULL;
1502 delete [] hypre_CSRMatrixData(A->offd);
1504 hypre_CSRMatrixData(A->offd) = NULL;
1506 if (colMapOwner >= 0)
1508 if (colMapOwner & 1)
1510 delete [] hypre_ParCSRMatrixColMapOffd(A);
1512 hypre_ParCSRMatrixColMapOffd(A) = NULL;
1517 hypre_ParCSRMatrixDestroy(A);
1524 hypre_ParCSRMatrix *C_hypre =
1525 internal::hypre_ParCSRMatrixAdd(const_cast<HypreParMatrix &>(A),
1526 const_cast<HypreParMatrix &>(B));
1527 MFEM_VERIFY(C_hypre,
"error in hypre_ParCSRMatrixAdd");
1529 hypre_MatvecCommPkgCreate(C_hypre);
1540 hypre_ParCSRMatrix * ab;
1541 ab = hypre_ParMatmul(*A,*B);
1542 hypre_ParCSRMatrixSetNumNonzeros(ab);
1544 hypre_MatvecCommPkgCreate(ab);
1551 hypre_ParCSRMatrix * C = internal::hypre_ParCSRMatrixAdd(*A,*B);
1553 hypre_MatvecCommPkgCreate(C);
1560 HYPRE_Int P_owns_its_col_starts =
1561 hypre_ParCSRMatrixOwnsColStarts((hypre_ParCSRMatrix*)(*P));
1563 hypre_ParCSRMatrix * rap;
1564 hypre_BoomerAMGBuildCoarseOperator(*P,*A,*P,&rap);
1565 hypre_ParCSRMatrixSetNumNonzeros(rap);
1570 hypre_ParCSRMatrixSetRowStartsOwner(rap,0);
1571 hypre_ParCSRMatrixSetColStartsOwner(rap,0);
1573 if (P_owns_its_col_starts)
1575 hypre_ParCSRMatrixSetColStartsOwner(*P, 1);
1584 HYPRE_Int P_owns_its_col_starts =
1585 hypre_ParCSRMatrixOwnsColStarts((hypre_ParCSRMatrix*)(*P));
1586 HYPRE_Int Rt_owns_its_col_starts =
1587 hypre_ParCSRMatrixOwnsColStarts((hypre_ParCSRMatrix*)(*Rt));
1589 hypre_ParCSRMatrix * rap;
1590 hypre_BoomerAMGBuildCoarseOperator(*Rt,*A,*P,&rap);
1592 hypre_ParCSRMatrixSetNumNonzeros(rap);
1597 hypre_ParCSRMatrixSetRowStartsOwner(rap,0);
1598 hypre_ParCSRMatrixSetColStartsOwner(rap,0);
1600 if (P_owns_its_col_starts)
1602 hypre_ParCSRMatrixSetColStartsOwner(*P, 1);
1604 if (Rt_owns_its_col_starts)
1606 hypre_ParCSRMatrixSetColStartsOwner(*Rt, 1);
1617 Ae.
Mult(-1.0, X, 1.0, B);
1619 hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag((hypre_ParCSRMatrix *)A);
1620 double *data = hypre_CSRMatrixData(A_diag);
1621 HYPRE_Int *I = hypre_CSRMatrixI(A_diag);
1623 HYPRE_Int *J = hypre_CSRMatrixJ(A_diag);
1624 hypre_CSRMatrix *A_offd = hypre_ParCSRMatrixOffd((hypre_ParCSRMatrix *)A);
1625 HYPRE_Int *I_offd = hypre_CSRMatrixI(A_offd);
1626 double *data_offd = hypre_CSRMatrixData(A_offd);
1629 for (
int i = 0; i < ess_dof_list.
Size(); i++)
1631 int r = ess_dof_list[i];
1632 B(r) = data[I[r]] * X(r);
1639 MFEM_ABORT(
"the diagonal entry must be the first entry in the row!");
1641 for (
int j = I[r]+1; j < I[r+1]; j++)
1645 MFEM_ABORT(
"all off-diagonal entries must be zero!");
1648 for (
int j = I_offd[r]; j < I_offd[r+1]; j++)
1650 if (data_offd[j] != 0.0)
1652 MFEM_ABORT(
"all off-diagonal entries must be zero!");
1672 hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag(A);
1673 HYPRE_Int num_rows = hypre_CSRMatrixNumRows(A_diag);
1675 double *u_data = hypre_VectorData(hypre_ParVectorLocalVector(u));
1676 double *r_data = hypre_VectorData(hypre_ParVectorLocalVector(r));
1678 for (
int i = 0; i < N; i++)
1681 hypre_ParVectorCopy(f, r);
1682 hypre_ParCSRMatrixMatvec(-1.0, A, u, 1.0, r);
1685 (0 == (i % 2)) ? coef = lambda : coef = mu;
1687 for (HYPRE_Int j = 0; j < num_rows; j++)
1689 u_data[j] += coef*r_data[j] / max_eig;
1705 hypre_ParVector *x0,
1706 hypre_ParVector *x1,
1707 hypre_ParVector *x2,
1708 hypre_ParVector *x3)
1711 hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag(A);
1712 HYPRE_Int num_rows = hypre_CSRMatrixNumRows(A_diag);
1714 double *u_data = hypre_VectorData(hypre_ParVectorLocalVector(u));
1716 double *x0_data = hypre_VectorData(hypre_ParVectorLocalVector(x0));
1717 double *x1_data = hypre_VectorData(hypre_ParVectorLocalVector(x1));
1718 double *x2_data = hypre_VectorData(hypre_ParVectorLocalVector(x2));
1719 double *x3_data = hypre_VectorData(hypre_ParVectorLocalVector(x3));
1721 hypre_ParVectorCopy(u, x0);
1724 hypre_ParVectorCopy(f, x1);
1725 hypre_ParCSRMatrixMatvec(-1.0, A, x0, 1.0, x1);
1727 for (HYPRE_Int i = 0; i < num_rows; i++)
1729 x1_data[i] /= -max_eig;
1733 for (HYPRE_Int i = 0; i < num_rows; i++)
1735 x1_data[i] = x0_data[i] -x1_data[i];
1739 for (HYPRE_Int i = 0; i < num_rows; i++)
1741 x3_data[i] = fir_coeffs[0]*x0_data[i] +fir_coeffs[1]*x1_data[i];
1744 for (
int n = 2; n <= poly_order; n++)
1747 hypre_ParVectorCopy(f, x2);
1748 hypre_ParCSRMatrixMatvec(-1.0, A, x1, 1.0, x2);
1750 for (HYPRE_Int i = 0; i < num_rows; i++)
1752 x2_data[i] /= -max_eig;
1760 for (HYPRE_Int i = 0; i < num_rows; i++)
1762 x2_data[i] = (x1_data[i]-x0_data[i]) +(x1_data[i]-2*x2_data[i]);
1763 x3_data[i] += fir_coeffs[n]*x2_data[i];
1764 x0_data[i] = x1_data[i];
1765 x1_data[i] = x2_data[i];
1769 for (HYPRE_Int i = 0; i < num_rows; i++)
1771 u_data[i] = x3_data[i];
1791 B =
X =
V =
Z = NULL;
1797 int _relax_times,
double _relax_weight,
double _omega,
1798 int _poly_order,
double _poly_fraction)
1809 B =
X =
V =
Z = NULL;
1818 type =
static_cast<int>(_type);
1844 double a = -1, b, c;
1845 if (!strcmp(name,
"Rectangular")) { a = 1.0, b = 0.0, c = 0.0; }
1846 if (!strcmp(name,
"Hanning")) { a = 0.5, b = 0.5, c = 0.0; }
1847 if (!strcmp(name,
"Hamming")) { a = 0.54, b = 0.46, c = 0.0; }
1848 if (!strcmp(name,
"Blackman")) { a = 0.42, b = 0.50, c = 0.08; }
1851 mfem_error(
"HypreSmoother::SetWindowByName : name not recognized!");
1869 mfem_error(
"HypreSmoother::SetOperator : not HypreParMatrix!");
1875 if (
B) {
delete B; }
1876 if (
X) {
delete X; }
1877 if (
V) {
delete V; }
1878 if (
Z) {
delete Z; }
1897 A->
Mult(ones, diag);
1906 for (
int i = 0; i <
height; i++)
1919 else if (
type == 1001 ||
type == 1002)
1950 double* window_coeffs =
new double[
poly_order+1];
1951 double* cheby_coeffs =
new double[
poly_order+1];
1959 window_coeffs[i] = a + b*cos(t) +c*cos(2*t);
1963 double theta_pb = acos(1.0 -0.5*k_pb);
1965 cheby_coeffs[0] = (theta_pb +
sigma)/M_PI;
1968 double t = i*(theta_pb+
sigma);
1969 cheby_coeffs[i] = 2.0*sin(t)/(i*M_PI);
1974 fir_coeffs[i] = window_coeffs[i]*cheby_coeffs[i];
1977 delete[] window_coeffs;
1978 delete[] cheby_coeffs;
1985 mfem_error(
"HypreSmoother::Mult (...) : HypreParMatrix A is missing");
1993 HYPRE_ParCSRDiagScale(NULL, *
A, b, x);
2017 else if (
type == 1002)
2032 hypre_ParCSRRelax(*
A, b,
type,
2037 hypre_ParCSRRelax(*
A, b,
type,
2048 mfem_error(
"HypreSmoother::Mult (...) : HypreParMatrix A is missing");
2054 A -> GetGlobalNumRows(),
2056 A -> GetRowStarts());
2058 A -> GetGlobalNumCols(),
2060 A -> GetColStarts());
2073 if (
B) {
delete B; }
2074 if (
X) {
delete X; }
2075 if (
V) {
delete V; }
2076 if (
Z) {
delete Z; }
2085 if (
X0) {
delete X0; }
2086 if (
X1) {
delete X1; }
2098 :
Solver(_A->Height(), _A->Width())
2109 mfem_error(
"HypreSolver::Mult (...) : HypreParMatrix A is missing");
2129 mfem_error(
"HypreSolver::Mult (...) : HypreParMatrix A is missing");
2135 A -> GetGlobalNumRows(),
2137 A -> GetRowStarts());
2139 A -> GetGlobalNumCols(),
2141 A -> GetColStarts());
2154 if (
B) {
delete B; }
2155 if (
X) {
delete X; }
2165 HYPRE_ParCSRMatrixGetComm(*
A, &comm);
2167 HYPRE_ParCSRPCGCreate(comm, &pcg_solver);
2172 HYPRE_PCGSetTol(pcg_solver, tol);
2177 HYPRE_PCGSetMaxIter(pcg_solver, max_iter);
2182 HYPRE_PCGSetLogging(pcg_solver, logging);
2187 HYPRE_ParCSRPCGSetPrintLevel(pcg_solver, print_lvl);
2192 HYPRE_ParCSRPCGSetPrecond(pcg_solver,
2200 HYPRE_PCGSetTwoNorm(pcg_solver, 1);
2201 if (res_frequency > 0)
2203 HYPRE_PCGSetRecomputeResidualP(pcg_solver, res_frequency);
2207 HYPRE_PCGSetResidualTol(pcg_solver, rtol);
2214 HYPRE_Int time_index = 0;
2215 HYPRE_Int num_iterations;
2216 double final_res_norm;
2218 HYPRE_Int print_level;
2220 HYPRE_PCGGetPrintLevel(pcg_solver, &print_level);
2221 HYPRE_ParCSRPCGSetPrintLevel(pcg_solver, print_level%3);
2223 HYPRE_ParCSRMatrixGetComm(*
A, &comm);
2227 if (print_level > 0 && print_level < 3)
2229 time_index = hypre_InitializeTiming(
"PCG Setup");
2230 hypre_BeginTiming(time_index);
2233 HYPRE_ParCSRPCGSetup(pcg_solver, *
A, b, x);
2236 if (print_level > 0 && print_level < 3)
2238 hypre_EndTiming(time_index);
2239 hypre_PrintTiming(
"Setup phase times", comm);
2240 hypre_FinalizeTiming(time_index);
2241 hypre_ClearTiming();
2245 if (print_level > 0 && print_level < 3)
2247 time_index = hypre_InitializeTiming(
"PCG Solve");
2248 hypre_BeginTiming(time_index);
2256 HYPRE_ParCSRPCGSolve(pcg_solver, *
A, b, x);
2258 if (print_level > 0)
2260 if (print_level < 3)
2262 hypre_EndTiming(time_index);
2263 hypre_PrintTiming(
"Solve phase times", comm);
2264 hypre_FinalizeTiming(time_index);
2265 hypre_ClearTiming();
2268 HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_iterations);
2269 HYPRE_ParCSRPCGGetFinalRelativeResidualNorm(pcg_solver,
2272 MPI_Comm_rank(comm, &myid);
2276 mfem::out <<
"PCG Iterations = " << num_iterations << endl
2277 <<
"Final PCG Relative Residual Norm = " << final_res_norm
2281 HYPRE_ParCSRPCGSetPrintLevel(pcg_solver, print_level);
2286 HYPRE_ParCSRPCGDestroy(pcg_solver);
2300 HYPRE_ParCSRMatrixGetComm(*
A, &comm);
2302 HYPRE_ParCSRGMRESCreate(comm, &gmres_solver);
2303 HYPRE_ParCSRGMRESSetKDim(gmres_solver, k_dim);
2304 HYPRE_ParCSRGMRESSetMaxIter(gmres_solver, max_iter);
2305 HYPRE_ParCSRGMRESSetTol(gmres_solver, tol);
2310 HYPRE_GMRESSetTol(gmres_solver, tol);
2315 HYPRE_GMRESSetMaxIter(gmres_solver, max_iter);
2320 HYPRE_GMRESSetKDim(gmres_solver, k_dim);
2325 HYPRE_GMRESSetLogging(gmres_solver, logging);
2330 HYPRE_GMRESSetPrintLevel(gmres_solver, print_lvl);
2335 HYPRE_ParCSRGMRESSetPrecond(gmres_solver,
2344 HYPRE_Int time_index = 0;
2345 HYPRE_Int num_iterations;
2346 double final_res_norm;
2348 HYPRE_Int print_level;
2350 HYPRE_GMRESGetPrintLevel(gmres_solver, &print_level);
2352 HYPRE_ParCSRMatrixGetComm(*
A, &comm);
2356 if (print_level > 0)
2358 time_index = hypre_InitializeTiming(
"GMRES Setup");
2359 hypre_BeginTiming(time_index);
2362 HYPRE_ParCSRGMRESSetup(gmres_solver, *
A, b, x);
2365 if (print_level > 0)
2367 hypre_EndTiming(time_index);
2368 hypre_PrintTiming(
"Setup phase times", comm);
2369 hypre_FinalizeTiming(time_index);
2370 hypre_ClearTiming();
2374 if (print_level > 0)
2376 time_index = hypre_InitializeTiming(
"GMRES Solve");
2377 hypre_BeginTiming(time_index);
2385 HYPRE_ParCSRGMRESSolve(gmres_solver, *
A, b, x);
2387 if (print_level > 0)
2389 hypre_EndTiming(time_index);
2390 hypre_PrintTiming(
"Solve phase times", comm);
2391 hypre_FinalizeTiming(time_index);
2392 hypre_ClearTiming();
2394 HYPRE_ParCSRGMRESGetNumIterations(gmres_solver, &num_iterations);
2395 HYPRE_ParCSRGMRESGetFinalRelativeResidualNorm(gmres_solver,
2398 MPI_Comm_rank(comm, &myid);
2402 mfem::out <<
"GMRES Iterations = " << num_iterations << endl
2403 <<
"Final GMRES Relative Residual Norm = " << final_res_norm
2411 HYPRE_ParCSRGMRESDestroy(gmres_solver);
2419 int sai_max_levels = 1;
2420 double sai_threshold = 0.1;
2421 double sai_filter = 0.1;
2423 double sai_loadbal = 0.0;
2425 int sai_logging = 1;
2427 HYPRE_ParCSRMatrixGetComm(A, &comm);
2429 HYPRE_ParaSailsCreate(comm, &sai_precond);
2430 HYPRE_ParaSailsSetParams(sai_precond, sai_threshold, sai_max_levels);
2431 HYPRE_ParaSailsSetFilter(sai_precond, sai_filter);
2432 HYPRE_ParaSailsSetSym(sai_precond, sai_sym);
2433 HYPRE_ParaSailsSetLoadbal(sai_precond, sai_loadbal);
2434 HYPRE_ParaSailsSetReuse(sai_precond, sai_reuse);
2435 HYPRE_ParaSailsSetLogging(sai_precond, sai_logging);
2440 HYPRE_ParaSailsSetSym(sai_precond, sym);
2445 HYPRE_ParaSailsDestroy(sai_precond);
2451 HYPRE_BoomerAMGCreate(&amg_precond);
2452 SetDefaultOptions();
2457 HYPRE_BoomerAMGCreate(&amg_precond);
2458 SetDefaultOptions();
2461 void HypreBoomerAMG::SetDefaultOptions()
2464 int coarsen_type = 10;
2466 double theta = 0.25;
2469 int interp_type = 6;
2474 int relax_sweeps = 1;
2477 int print_level = 1;
2478 int max_levels = 25;
2480 HYPRE_BoomerAMGSetCoarsenType(amg_precond, coarsen_type);
2481 HYPRE_BoomerAMGSetAggNumLevels(amg_precond, agg_levels);
2482 HYPRE_BoomerAMGSetRelaxType(amg_precond, relax_type);
2483 HYPRE_BoomerAMGSetNumSweeps(amg_precond, relax_sweeps);
2484 HYPRE_BoomerAMGSetStrongThreshold(amg_precond, theta);
2485 HYPRE_BoomerAMGSetInterpType(amg_precond, interp_type);
2486 HYPRE_BoomerAMGSetPMaxElmts(amg_precond, Pmax);
2487 HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level);
2488 HYPRE_BoomerAMGSetMaxLevels(amg_precond, max_levels);
2491 HYPRE_BoomerAMGSetMaxIter(amg_precond, 1);
2492 HYPRE_BoomerAMGSetTol(amg_precond, 0.0);
2495 void HypreBoomerAMG::ResetAMGPrecond()
2497 HYPRE_Int coarsen_type;
2498 HYPRE_Int agg_levels;
2499 HYPRE_Int relax_type;
2500 HYPRE_Int relax_sweeps;
2502 HYPRE_Int interp_type;
2504 HYPRE_Int print_level;
2506 HYPRE_Int nrbms = rbms.
Size();
2508 HYPRE_Int nodal_diag;
2509 HYPRE_Int relax_coarse;
2510 HYPRE_Int interp_vec_variant;
2512 HYPRE_Int smooth_interp_vectors;
2513 HYPRE_Int interp_refine;
2515 hypre_ParAMGData *amg_data = (hypre_ParAMGData *)amg_precond;
2518 HYPRE_BoomerAMGGetCoarsenType(amg_precond, &coarsen_type);
2519 agg_levels = hypre_ParAMGDataAggNumLevels(amg_data);
2520 relax_type = hypre_ParAMGDataUserRelaxType(amg_data);
2521 relax_sweeps = hypre_ParAMGDataUserNumSweeps(amg_data);
2522 HYPRE_BoomerAMGGetStrongThreshold(amg_precond, &theta);
2523 hypre_BoomerAMGGetInterpType(amg_precond, &interp_type);
2524 HYPRE_BoomerAMGGetPMaxElmts(amg_precond, &Pmax);
2525 HYPRE_BoomerAMGGetPrintLevel(amg_precond, &print_level);
2526 HYPRE_BoomerAMGGetNumFunctions(amg_precond, &dim);
2529 nodal = hypre_ParAMGDataNodal(amg_data);
2530 nodal_diag = hypre_ParAMGDataNodalDiag(amg_data);
2531 HYPRE_BoomerAMGGetCycleRelaxType(amg_precond, &relax_coarse, 3);
2532 interp_vec_variant = hypre_ParAMGInterpVecVariant(amg_data);
2533 q_max = hypre_ParAMGInterpVecQMax(amg_data);
2534 smooth_interp_vectors = hypre_ParAMGSmoothInterpVectors(amg_data);
2535 interp_refine = hypre_ParAMGInterpRefine(amg_data);
2538 HYPRE_BoomerAMGDestroy(amg_precond);
2539 HYPRE_BoomerAMGCreate(&amg_precond);
2541 HYPRE_BoomerAMGSetCoarsenType(amg_precond, coarsen_type);
2542 HYPRE_BoomerAMGSetAggNumLevels(amg_precond, agg_levels);
2543 HYPRE_BoomerAMGSetRelaxType(amg_precond, relax_type);
2544 HYPRE_BoomerAMGSetNumSweeps(amg_precond, relax_sweeps);
2545 HYPRE_BoomerAMGSetMaxLevels(amg_precond, 25);
2546 HYPRE_BoomerAMGSetTol(amg_precond, 0.0);
2547 HYPRE_BoomerAMGSetMaxIter(amg_precond, 1);
2548 HYPRE_BoomerAMGSetStrongThreshold(amg_precond, theta);
2549 HYPRE_BoomerAMGSetInterpType(amg_precond, interp_type);
2550 HYPRE_BoomerAMGSetPMaxElmts(amg_precond, Pmax);
2551 HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level);
2552 HYPRE_BoomerAMGSetNumFunctions(amg_precond, dim);
2555 HYPRE_BoomerAMGSetNodal(amg_precond, nodal);
2556 HYPRE_BoomerAMGSetNodalDiag(amg_precond, nodal_diag);
2557 HYPRE_BoomerAMGSetCycleRelaxType(amg_precond, relax_coarse, 3);
2558 HYPRE_BoomerAMGSetInterpVecVariant(amg_precond, interp_vec_variant);
2559 HYPRE_BoomerAMGSetInterpVecQMax(amg_precond, q_max);
2560 HYPRE_BoomerAMGSetSmoothInterpVectors(amg_precond, smooth_interp_vectors);
2561 HYPRE_BoomerAMGSetInterpRefine(amg_precond, interp_refine);
2563 HYPRE_BoomerAMGSetInterpVectors(amg_precond, rbms.
Size(), rbms.
GetData());
2570 MFEM_VERIFY(new_A,
"new Operator must be a HypreParMatrix!");
2572 if (
A) { ResetAMGPrecond(); }
2586 HYPRE_BoomerAMGSetNumFunctions(amg_precond, dim);
2589 HYPRE_BoomerAMGSetAggNumLevels(amg_precond, 0);
2590 HYPRE_BoomerAMGSetStrongThreshold(amg_precond, 0.5);
2596 y = 0.0; y(0) = x(1); y(1) = -x(0);
2598 static void func_ryz(
const Vector &x, Vector &y)
2600 y = 0.0; y(1) = x(2); y(2) = -x(1);
2602 static void func_rzx(
const Vector &x, Vector &y)
2604 y = 0.0; y(2) = x(0); y(0) = -x(2);
2607 void HypreBoomerAMG::RecomputeRBMs()
2610 Array<HypreParVector*> gf_rbms;
2613 for (
int i = 0; i < rbms.
Size(); i++)
2615 HYPRE_ParVectorDestroy(rbms[i]);
2622 VectorFunctionCoefficient coeff_rxy(2, func_rxy);
2624 ParGridFunction rbms_rxy(fespace);
2625 rbms_rxy.ProjectCoefficient(coeff_rxy);
2628 gf_rbms.SetSize(nrbms);
2629 gf_rbms[0] = rbms_rxy.ParallelAverage();
2635 VectorFunctionCoefficient coeff_rxy(3, func_rxy);
2636 VectorFunctionCoefficient coeff_ryz(3, func_ryz);
2637 VectorFunctionCoefficient coeff_rzx(3, func_rzx);
2639 ParGridFunction rbms_rxy(fespace);
2640 ParGridFunction rbms_ryz(fespace);
2641 ParGridFunction rbms_rzx(fespace);
2642 rbms_rxy.ProjectCoefficient(coeff_rxy);
2643 rbms_ryz.ProjectCoefficient(coeff_ryz);
2644 rbms_rzx.ProjectCoefficient(coeff_rzx);
2647 gf_rbms.SetSize(nrbms);
2648 gf_rbms[0] = rbms_rxy.ParallelAverage();
2649 gf_rbms[1] = rbms_ryz.ParallelAverage();
2650 gf_rbms[2] = rbms_rzx.ParallelAverage();
2659 for (
int i = 0; i < nrbms; i++)
2661 rbms[i] = gf_rbms[i]->StealParVector();
2669 this->fespace = fespace;
2679 int relax_coarse = 8;
2682 int interp_vec_variant = 2;
2684 int smooth_interp_vectors = 1;
2688 int interp_refine = 1;
2690 HYPRE_BoomerAMGSetNodal(amg_precond, nodal);
2691 HYPRE_BoomerAMGSetNodalDiag(amg_precond, nodal_diag);
2692 HYPRE_BoomerAMGSetCycleRelaxType(amg_precond, relax_coarse, 3);
2693 HYPRE_BoomerAMGSetInterpVecVariant(amg_precond, interp_vec_variant);
2694 HYPRE_BoomerAMGSetInterpVecQMax(amg_precond, q_max);
2695 HYPRE_BoomerAMGSetSmoothInterpVectors(amg_precond, smooth_interp_vectors);
2696 HYPRE_BoomerAMGSetInterpRefine(amg_precond, interp_refine);
2699 HYPRE_BoomerAMGSetInterpVectors(amg_precond, rbms.
Size(), rbms.
GetData());
2704 for (
int i = 0; i < rbms.
Size(); i++)
2706 HYPRE_ParVectorDestroy(rbms[i]);
2709 HYPRE_BoomerAMGDestroy(amg_precond);
2716 int cycle_type = 13;
2719 double rlx_weight = 1.0;
2720 double rlx_omega = 1.0;
2721 int amg_coarsen_type = 10;
2722 int amg_agg_levels = 1;
2723 int amg_rlx_type = 8;
2724 double theta = 0.25;
2725 int amg_interp_type = 6;
2732 bool trace_space, rt_trace_space;
2736 trace_space = trace_space || rt_trace_space;
2739 if (edge_fespace->
GetNE() > 0)
2744 if (dim == 2) { p++; }
2759 HYPRE_AMSCreate(&ams);
2761 HYPRE_AMSSetDimension(ams, sdim);
2762 HYPRE_AMSSetTol(ams, 0.0);
2763 HYPRE_AMSSetMaxIter(ams, 1);
2764 HYPRE_AMSSetCycleType(ams, cycle_type);
2765 HYPRE_AMSSetPrintLevel(ams, 1);
2787 for (
int i = 0; i < pmesh->
GetNV(); i++)
2789 coord = pmesh -> GetVertex(i);
2790 x_coord(i) = coord[0];
2791 y_coord(i) = coord[1];
2792 if (sdim == 3) { z_coord(i) = coord[2]; }
2799 HYPRE_AMSSetCoordinateVectors(ams, *x, *y, NULL);
2804 HYPRE_AMSSetCoordinateVectors(ams, *x, *y, *z);
2828 HYPRE_AMSSetDiscreteGradient(ams, *G);
2832 Pi = Pix = Piy = Piz = NULL;
2851 if (cycle_type < 10)
2859 Pix = Pi_blocks(0,0);
2860 Piy = Pi_blocks(0,1);
2861 if (sdim == 3) { Piz = Pi_blocks(0,2); }
2866 HYPRE_ParCSRMatrix HY_Pi = (Pi) ? (HYPRE_ParCSRMatrix) *Pi : NULL;
2867 HYPRE_ParCSRMatrix HY_Pix = (Pix) ? (HYPRE_ParCSRMatrix) *Pix : NULL;
2868 HYPRE_ParCSRMatrix HY_Piy = (Piy) ? (HYPRE_ParCSRMatrix) *Piy : NULL;
2869 HYPRE_ParCSRMatrix HY_Piz = (Piz) ? (HYPRE_ParCSRMatrix) *Piz : NULL;
2870 HYPRE_AMSSetInterpolations(ams, HY_Pi, HY_Pix, HY_Piy, HY_Piz);
2872 delete vert_fespace_d;
2875 delete vert_fespace;
2880 delete edge_fespace;
2885 HYPRE_AMSSetSmoothingOptions(ams, rlx_type, rlx_sweeps, rlx_weight, rlx_omega);
2886 HYPRE_AMSSetAlphaAMGOptions(ams, amg_coarsen_type, amg_agg_levels, amg_rlx_type,
2887 theta, amg_interp_type, amg_Pmax);
2888 HYPRE_AMSSetBetaAMGOptions(ams, amg_coarsen_type, amg_agg_levels, amg_rlx_type,
2889 theta, amg_interp_type, amg_Pmax);
2894 HYPRE_AMSDestroy(ams);
2909 HYPRE_AMSSetPrintLevel(ams, print_lvl);
2915 int cycle_type = 11;
2918 double rlx_weight = 1.0;
2919 double rlx_omega = 1.0;
2920 int amg_coarsen_type = 10;
2921 int amg_agg_levels = 1;
2922 int amg_rlx_type = 8;
2923 double theta = 0.25;
2924 int amg_interp_type = 6;
2926 int ams_cycle_type = 14;
2932 if (face_fespace->
GetNE() > 0)
2944 HYPRE_ADSCreate(&ads);
2946 HYPRE_ADSSetTol(ads, 0.0);
2947 HYPRE_ADSSetMaxIter(ads, 1);
2948 HYPRE_ADSSetCycleType(ads, cycle_type);
2949 HYPRE_ADSSetPrintLevel(ads, 1);
2977 for (
int i = 0; i < pmesh->
GetNV(); i++)
2979 coord = pmesh -> GetVertex(i);
2980 x_coord(i) = coord[0];
2981 y_coord(i) = coord[1];
2982 z_coord(i) = coord[2];
2987 HYPRE_ADSSetCoordinateVectors(ads, *x, *y, *z);
3011 HYPRE_ADSSetDiscreteCurl(ads, *C);
3030 HYPRE_ADSSetDiscreteGradient(ads, *G);
3034 RT_Pi = RT_Pix = RT_Piy = RT_Piz = NULL;
3035 ND_Pi = ND_Pix = ND_Piy = ND_Piz = NULL;
3054 if (ams_cycle_type < 10)
3064 ND_Pix = ND_Pi_blocks(0,0);
3065 ND_Piy = ND_Pi_blocks(0,1);
3066 ND_Piz = ND_Pi_blocks(0,2);
3084 if (cycle_type < 10)
3093 RT_Pix = RT_Pi_blocks(0,0);
3094 RT_Piy = RT_Pi_blocks(0,1);
3095 RT_Piz = RT_Pi_blocks(0,2);
3100 HYPRE_ParCSRMatrix HY_RT_Pi, HY_RT_Pix, HY_RT_Piy, HY_RT_Piz;
3101 HY_RT_Pi = (RT_Pi) ? (HYPRE_ParCSRMatrix) *RT_Pi : NULL;
3102 HY_RT_Pix = (RT_Pix) ? (HYPRE_ParCSRMatrix) *RT_Pix : NULL;
3103 HY_RT_Piy = (RT_Piy) ? (HYPRE_ParCSRMatrix) *RT_Piy : NULL;
3104 HY_RT_Piz = (RT_Piz) ? (HYPRE_ParCSRMatrix) *RT_Piz : NULL;
3105 HYPRE_ParCSRMatrix HY_ND_Pi, HY_ND_Pix, HY_ND_Piy, HY_ND_Piz;
3106 HY_ND_Pi = (ND_Pi) ? (HYPRE_ParCSRMatrix) *ND_Pi : NULL;
3107 HY_ND_Pix = (ND_Pix) ? (HYPRE_ParCSRMatrix) *ND_Pix : NULL;
3108 HY_ND_Piy = (ND_Piy) ? (HYPRE_ParCSRMatrix) *ND_Piy : NULL;
3109 HY_ND_Piz = (ND_Piz) ? (HYPRE_ParCSRMatrix) *ND_Piz : NULL;
3110 HYPRE_ADSSetInterpolations(ads,
3111 HY_RT_Pi, HY_RT_Pix, HY_RT_Piy, HY_RT_Piz,
3112 HY_ND_Pi, HY_ND_Pix, HY_ND_Piy, HY_ND_Piz);
3114 delete vert_fespace_d;
3118 delete vert_fespace;
3120 delete edge_fespace;
3123 HYPRE_ADSSetSmoothingOptions(ads, rlx_type, rlx_sweeps, rlx_weight, rlx_omega);
3124 HYPRE_ADSSetAMGOptions(ads, amg_coarsen_type, amg_agg_levels, amg_rlx_type,
3125 theta, amg_interp_type, amg_Pmax);
3126 HYPRE_ADSSetAMSOptions(ads, ams_cycle_type, amg_coarsen_type, amg_agg_levels,
3127 amg_rlx_type, theta, amg_interp_type, amg_Pmax);
3132 HYPRE_ADSDestroy(ads);
3154 HYPRE_ADSSetPrintLevel(ads, print_lvl);
3157 HypreLOBPCG::HypreMultiVector::HypreMultiVector(
int n,
HypreParVector & v,
3158 mv_InterfaceInterpreter & interpreter)
3162 mv_ptr = mv_MultiVectorCreateFromSampleVector(&interpreter, nv,
3163 (HYPRE_ParVector)v);
3165 HYPRE_ParVector* vecs = NULL;
3167 mv_TempMultiVector* tmp =
3168 (mv_TempMultiVector*)mv_MultiVectorGetData(mv_ptr);
3169 vecs = (HYPRE_ParVector*)(tmp -> vector);
3172 hpv =
new HypreParVector*[nv];
3173 for (
int i=0; i<nv; i++)
3175 hpv[i] =
new HypreParVector(vecs[i]);
3179 HypreLOBPCG::HypreMultiVector::~HypreMultiVector()
3183 for (
int i=0; i<nv; i++)
3190 mv_MultiVectorDestroy(mv_ptr);
3194 HypreLOBPCG::HypreMultiVector::Randomize(HYPRE_Int seed)
3196 mv_MultiVectorSetRandom(mv_ptr, seed);
3200 HypreLOBPCG::HypreMultiVector::GetVector(
unsigned int i)
3202 MFEM_ASSERT((
int)i < nv,
"index out of range");
3208 HypreLOBPCG::HypreMultiVector::StealVectors()
3210 HypreParVector ** hpv_ret = hpv;
3214 mv_TempMultiVector * mv_tmp =
3215 (mv_TempMultiVector*)mv_MultiVectorGetData(mv_ptr);
3217 mv_tmp->ownsVectors = 0;
3219 for (
int i=0; i<nv; i++)
3221 hpv_ret[i]->SetOwnership(1);
3239 MPI_Comm_size(comm,&numProcs);
3240 MPI_Comm_rank(comm,&myid);
3242 HYPRE_ParCSRSetupInterpreter(&interpreter);
3243 HYPRE_ParCSRSetupMatvec(&matvec_fn);
3244 HYPRE_LOBPCGCreate(&interpreter, &matvec_fn, &lobpcg_solver);
3253 HYPRE_LOBPCGDestroy(lobpcg_solver);
3259 HYPRE_LOBPCGSetTol(lobpcg_solver, tol);
3265 #if MFEM_HYPRE_VERSION >= 21101
3266 HYPRE_LOBPCGSetRTol(lobpcg_solver, rel_tol);
3268 MFEM_ABORT(
"This method requires HYPRE version >= 2.11.1");
3275 HYPRE_LOBPCGSetMaxIter(lobpcg_solver, max_iter);
3283 HYPRE_LOBPCGSetPrintLevel(lobpcg_solver, logging);
3290 HYPRE_LOBPCGSetPrecondUsageMode(lobpcg_solver, pcg_mode);
3296 HYPRE_LOBPCGSetPrecond(lobpcg_solver,
3297 (HYPRE_PtrToSolverFcn)this->PrecondSolve,
3298 (HYPRE_PtrToSolverFcn)this->PrecondSetup,
3299 (HYPRE_Solver)&precond);
3305 int locSize = A.
Width();
3307 if (HYPRE_AssumedPartitionCheck())
3309 part =
new HYPRE_Int[2];
3311 MPI_Scan(&locSize, &part[1], 1, HYPRE_MPI_INT, MPI_SUM, comm);
3313 part[0] = part[1] - locSize;
3315 MPI_Allreduce(&locSize, &glbSize, 1, HYPRE_MPI_INT, MPI_SUM, comm);
3319 part =
new HYPRE_Int[numProcs+1];
3321 MPI_Allgather(&locSize, 1, MPI_INT,
3322 &part[1], 1, HYPRE_MPI_INT, comm);
3325 for (
int i=0; i<numProcs; i++)
3327 part[i+1] += part[i];
3330 glbSize = part[numProcs];
3341 matvec_fn.MatvecCreate = this->OperatorMatvecCreate;
3342 matvec_fn.Matvec = this->OperatorMatvec;
3343 matvec_fn.MatvecDestroy = this->OperatorMatvecDestroy;
3345 HYPRE_LOBPCGSetup(lobpcg_solver,(HYPRE_Matrix)&A,NULL,NULL);
3351 matvec_fn.MatvecCreate = this->OperatorMatvecCreate;
3352 matvec_fn.Matvec = this->OperatorMatvec;
3353 matvec_fn.MatvecDestroy = this->OperatorMatvecDestroy;
3355 HYPRE_LOBPCGSetupB(lobpcg_solver,(HYPRE_Matrix)&M,NULL);
3364 for (
int i=0; i<nev; i++)
3366 eigs[i] = eigenvalues[i];
3373 return multi_vec->GetVector(i);
3380 if ( multi_vec == NULL )
3382 MFEM_ASSERT(x != NULL,
"In HypreLOBPCG::SetInitialVectors()");
3384 multi_vec =
new HypreMultiVector(nev, *x, interpreter);
3388 for (
int i=0; i < min(num_vecs,nev); i++)
3390 multi_vec->GetVector(i) = *vecs[i];
3394 for (
int i=min(num_vecs,nev); i < nev; i++)
3396 multi_vec->GetVector(i).Randomize(seed);
3400 if ( subSpaceProj != NULL )
3403 y = multi_vec->GetVector(0);
3405 for (
int i=1; i<nev; i++)
3407 subSpaceProj->
Mult(multi_vec->GetVector(i),
3408 multi_vec->GetVector(i-1));
3410 subSpaceProj->
Mult(y,
3411 multi_vec->GetVector(nev-1));
3419 if ( multi_vec == NULL )
3421 MFEM_ASSERT(x != NULL,
"In HypreLOBPCG::Solve()");
3423 multi_vec =
new HypreMultiVector(nev, *x, interpreter);
3424 multi_vec->Randomize(seed);
3426 if ( subSpaceProj != NULL )
3429 y = multi_vec->GetVector(0);
3431 for (
int i=1; i<nev; i++)
3433 subSpaceProj->
Mult(multi_vec->GetVector(i),
3434 multi_vec->GetVector(i-1));
3436 subSpaceProj->
Mult(y, multi_vec->GetVector(nev-1));
3447 HYPRE_LOBPCGSolve(lobpcg_solver, NULL, *multi_vec, eigenvalues);
3451 HypreLOBPCG::OperatorMatvecCreate(
void *A,
3458 return ( matvec_data );
3462 HypreLOBPCG::OperatorMatvec(
void *matvec_data,
3463 HYPRE_Complex
alpha,
3469 MFEM_VERIFY(alpha == 1.0 && beta == 0.0,
"values not supported");
3471 Operator *Aop = (Operator*)A;
3473 int width = Aop->Width();
3475 hypre_ParVector * xPar = (hypre_ParVector *)x;
3476 hypre_ParVector * yPar = (hypre_ParVector *)y;
3478 Vector xVec(xPar->local_vector->data, width);
3479 Vector yVec(yPar->local_vector->data, width);
3481 Aop->Mult( xVec, yVec );
3487 HypreLOBPCG::OperatorMatvecDestroy(
void *matvec_data )
3493 HypreLOBPCG::PrecondSolve(
void *solver,
3498 Solver *PC = (Solver*)solver;
3499 Operator *OP = (Operator*)A;
3501 int width = OP->Width();
3503 hypre_ParVector * bPar = (hypre_ParVector *)b;
3504 hypre_ParVector * xPar = (hypre_ParVector *)x;
3506 Vector bVec(bPar->local_vector->data, width);
3507 Vector xVec(xPar->local_vector->data, width);
3509 PC->Mult( bVec, xVec );
3515 HypreLOBPCG::PrecondSetup(
void *solver,
3533 MPI_Comm_size(comm,&numProcs);
3534 MPI_Comm_rank(comm,&myid);
3536 HYPRE_AMECreate(&ame_solver);
3537 HYPRE_AMESetPrintLevel(ame_solver, 0);
3544 mfem_hypre_TFree(multi_vec);
3549 for (
int i=0; i<nev; i++)
3551 delete eigenvectors[i];
3554 delete [] eigenvectors;
3558 mfem_hypre_TFree(eigenvalues);
3561 HYPRE_AMEDestroy(ame_solver);
3569 HYPRE_AMESetBlockSize(ame_solver, nev);
3575 HYPRE_AMESetTol(ame_solver, tol);
3581 #if MFEM_HYPRE_VERSION >= 21101
3582 HYPRE_AMESetRTol(ame_solver, rel_tol);
3584 MFEM_ABORT(
"This method requires HYPRE version >= 2.11.1");
3591 HYPRE_AMESetMaxIter(ame_solver, max_iter);
3599 HYPRE_AMESetPrintLevel(ame_solver, logging);
3606 ams_precond = &precond;
3614 HYPRE_Solver ams_precond_ptr = (HYPRE_Solver)*ams_precond;
3616 ams_precond->
SetupFcn()(*ams_precond,A,NULL,NULL);
3618 HYPRE_AMESetAMSSolver(ame_solver, ams_precond_ptr);
3621 HYPRE_AMESetup(ame_solver);
3627 HYPRE_ParCSRMatrix parcsr_M = M;
3628 HYPRE_AMESetMassMatrix(ame_solver,(HYPRE_ParCSRMatrix)parcsr_M);
3634 HYPRE_AMESolve(ame_solver);
3641 eigs.
SetSize(nev); eigs = -1.0;
3643 if ( eigenvalues == NULL )
3646 HYPRE_AMEGetEigenvalues(ame_solver,&eigenvalues);
3650 for (
int i=0; i<nev; i++)
3652 eigs[i] = eigenvalues[i];
3657 HypreAME::createDummyVectors()
3659 if ( multi_vec == NULL )
3661 HYPRE_AMEGetEigenvectors(ame_solver,&multi_vec);
3664 eigenvectors =
new HypreParVector*[nev];
3665 for (
int i=0; i<nev; i++)
3667 eigenvectors[i] =
new HypreParVector(multi_vec[i]);
3676 if ( eigenvectors == NULL )
3678 this->createDummyVectors();
3681 return *eigenvectors[i];
3687 if ( eigenvectors == NULL )
3689 this->createDummyVectors();
3694 eigenvectors = NULL;
virtual ~HypreBoomerAMG()
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
int Size() const
Logical size of the array.
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Vector()
Default constructor for Vector. Sets size = 0 and data = NULL.
MPI_Comm GetComm() const
MPI communicator.
Vector * GlobalVector() const
Returns the global vector in each processor.
HypreParVector * X0
FIR Filter Temporary Vectors.
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
double min_eig_est
Minimal eigenvalue estimate for polynomial smoothing.
void SetPrintLevel(int print_lvl)
void Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0)
Prints the locally owned rows in parallel.
HYPRE_Int N() const
Returns the global number of columns.
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
int setup_called
Was hypre's Setup function called already?
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to 'master'.
void Read_IJMatrix(MPI_Comm comm, const char *fname)
Read a matrix saved as a HYPRE_IJMatrix.
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_Int *row_starts=NULL) const
Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result a...
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A^t * x + beta * y.
void SetSize(int s)
Resize the vector to size s.
HypreParMatrix & Add(const double beta, const HypreParMatrix &B)
HypreParVector * B
Right-hand side and solution vector.
double window_params[3]
Parameters for windowing function of FIR filter.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
void SetInitialVectors(int num_vecs, HypreParVector **vecs)
HypreADS(HypreParMatrix &A, ParFiniteElementSpace *face_fespace)
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
int GetFaceOrder(int i) const
Returns the order of the i'th face finite element.
T * GetData()
Returns the data.
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
void SetPreconditioner(HypreSolver &precond)
int Size() const
Returns the size of the vector.
virtual N_Vector ToNVector()
Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARHYP.
void SetPreconditioner(Solver &precond)
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void SetMassMatrix(Operator &M)
Abstract parallel finite element space.
double Normlinf() const
Returns the l_infinity norm of the vector.
void SetPrintLevel(int logging)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void SetOperator(HypreParMatrix &A)
void LoseData()
Lose the ownership of the graph (I, J) and data (A) arrays.
bool pos_l1_norms
If set, take absolute values of the computed l1_norms.
double * GetData() const
Return a pointer to the beginning of the Vector data.
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
void SetPrintLevel(int print_lvl)
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
int Size_of_connections() const
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A * B.
void SetPrintLevel(int print_lvl)
int poly_scale
Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.
void AddDomainInterpolator(DiscreteInterpolator *di)
void SetWindowParameters(double a, double b, double c)
Set parameters for windowing function for FIR smoother.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's GMRES.
void Print(const char *fname) const
Prints the locally owned rows in parallel.
void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
void SetSystemsOptions(int dim)
HypreLOBPCG(MPI_Comm comm)
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
int GetNE() const
Returns number of elements in the mesh.
HYPRE_Int GetGlobalNumRows() const
void SetSymmetry(int sym)
virtual ~HypreParaSails()
HYPRE_Int GetGlobalNumCols() const
void SetMaxIter(int max_iter)
double * l1_norms
l1 norms of the rows of A
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void SetLogging(int logging)
Mesh * GetMesh() const
Returns the mesh.
HyprePCG(HypreParMatrix &_A)
virtual HYPRE_PtrToParSolverFcn SolveFcn() const =0
hypre's internal Solve function
void SetRelTol(double rel_tol)
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0)
int to_int(const std::string &str)
void GetBlocks(Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
HYPRE_Int GlobalTrueVSize() const
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
HypreParMatrix * Transpose() const
Returns the transpose of *this.
HYPRE_Int Randomize(HYPRE_Int seed)
Set random values.
double relax_weight
Damping coefficient (usually <= 1)
void SetElasticityOptions(ParFiniteElementSpace *fespace)
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
void Sort()
Sorts the array. This requires operator< to be defined for T.
int Size() const
Returns the number of TYPE I elements.
HypreParMatrix * A
The linear system matrix.
double * GetData() const
Return element data, i.e. array A.
HypreAMS(HypreParMatrix &A, ParFiniteElementSpace *edge_fespace)
int * GetI() const
Return the array I.
Dynamic 2D array using row-major layout.
virtual void SetOperator(const Operator &op)
virtual HYPRE_PtrToParSolverFcn SetupFcn() const =0
hypre's internal Setup function
void SetLogging(int logging)
HypreParMatrix()
An empty matrix to be used as a reference to an existing matrix.
void SetMaxIter(int max_iter)
int ParCSRRelax_FIR(hypre_ParCSRMatrix *A, hypre_ParVector *f, double max_eig, int poly_order, double *fir_coeffs, hypre_ParVector *u, hypre_ParVector *x0, hypre_ParVector *x1, hypre_ParVector *x2, hypre_ParVector *x3)
HypreParMatrix * ParallelAssemble() const
Returns the matrix "assembled" on the true dofs.
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
int SpaceDimension() const
Wrapper for hypre's parallel vector class.
double Norml1() const
Returns the l_1 norm of the vector.
int ParCSRRelax_Taubin(hypre_ParCSRMatrix *A, hypre_ParVector *f, double lambda, double mu, int N, double max_eig, hypre_ParVector *u, hypre_ParVector *r)
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
void SetPrintLevel(int print_lvl)
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
void Solve()
Solve the eigenproblem.
void AddTraceFaceInterpolator(DiscreteInterpolator *di)
int GetOrder(int i) const
Returns the order of the i'th finite element.
void SetNumModes(int num_eigs)
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A + B.
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
void SetRelTol(double rel_tol)
int height
Dimension of the output / number of rows in the matrix.
void SetMaxIter(int max_iter)
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
double * fir_coeffs
Combined coefficients for windowing and Chebyshev polynomials.
void Threshold(double threshold=0.0)
Remove values smaller in absolute value than some threshold.
double InnerProduct(HypreParVector *x, HypreParVector *y)
HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetData(double *_data)
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
void SetTaubinOptions(double lambda, double mu, int iter)
Set parameters for Taubin's lambda-mu method.
HypreParVector * B
Right-hand side and solution vectors.
void SetMassMatrix(HypreParMatrix &M)
HYPRE_Int M() const
Returns the global number of rows.
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
void SetOperator(Operator &A)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void SetPolyOptions(int poly_order, double poly_fraction)
Set parameters for polynomial smoothing.
int relax_times
Number of relaxation sweeps.
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
void SetPrintLevel(int logging)
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
double ParNormlp(const Vector &vec, double p, MPI_Comm comm)
Compute the l_p norm of the Vector which is split without overlap across the given communicator...
void SetOwnership(int own)
Sets ownership of the internal hypre_ParVector.
void SetDataOwner(bool owna)
Set the data ownership flag (A array).
Abstract class for hypre's solvers and preconditioners.
HypreParVector & operator=(double d)
Set constant values.
HYPRE_Int * GetTrueDofOffsets() const
Arbitrary order H(curl)-conforming Nedelec finite elements.
const FiniteElementCollection * FEColl() const
void PrintCommPkg(std::ostream &out=mfem::out) const
Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
HypreGMRES(HypreParMatrix &_A)
void GetParBlocks(Array2D< HypreParMatrix * > &blocks) const
hypre_ParCSRMatrix * StealData()
Changes the ownership of the the matrix.
void SetPrecondUsageMode(int pcg_mode)
HypreParMatrix * A
The linear system matrix.
Arbitrary order H1-conforming (continuous) finite elements.
void SetMaxIter(int max_iter)
void GetOffd(SparseMatrix &offd, HYPRE_Int *&cmap) const
Get the local off-diagonal block. NOTE: 'offd' will not own any data.
HYPRE_Int * GetColStarts() const
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
double omega
SOR parameter (usually in (0,2))
HYPRE_Int * GetRowStarts() const
Class for parallel grid function.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Wrapper for hypre's ParCSR matrix class.
void Swap(SparseMatrix &other)
~HypreParVector()
Calls hypre's destroy function.
virtual void Assemble(int skip_zeros=1)
int * GetJ() const
Return the array J.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
Class for parallel meshes.
void Solve()
Solve the eigenproblem.
void Read(MPI_Comm comm, const char *fname)
Reads the matrix from a file.
HypreParVector * V
Temporary vectors.
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
HypreParaSails(HypreParMatrix &A)
int width
Dimension of the input / number of columns in the matrix.
int poly_order
Order of the smoothing polynomial.
double lambda
Taubin's lambda-mu method parameters.
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
double sigma(const Vector &x)