12 #include "../config/config.hpp"
17 #include "../fem/fem.hpp"
29 template<
typename TargetT,
typename SourceT>
30 static TargetT *DuplicateAs(
const SourceT *array,
int size,
31 bool cplusplus =
true)
33 TargetT *target_array = cplusplus ?
new TargetT[size]
34 : hypre_TAlloc(TargetT, size);
35 for (
int i = 0; i < size; i++)
37 target_array[i] = array[i];
42 inline void HypreParVector::_SetDataAndSize_()
44 SetDataAndSize(hypre_VectorData(hypre_ParVectorLocalVector(x)),
46 hypre_VectorSize(hypre_ParVectorLocalVector(x))));
49 HypreParVector::HypreParVector(MPI_Comm comm, HYPRE_Int glob_size,
52 x = hypre_ParVectorCreate(comm,glob_size,col);
53 hypre_ParVectorInitialize(x);
54 hypre_ParVectorSetPartitioningOwner(x,0);
56 hypre_ParVectorSetDataOwner(x,1);
57 hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),1);
63 double *_data, HYPRE_Int *col) :
Vector()
65 x = hypre_ParVectorCreate(comm,glob_size,col);
66 hypre_ParVectorSetDataOwner(x,1);
67 hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),0);
68 hypre_ParVectorSetPartitioningOwner(x,0);
70 hypre_VectorData(hypre_ParVectorLocalVector(x)) = &tmp;
73 hypre_ParVectorInitialize(x);
75 hypre_VectorData(hypre_ParVectorLocalVector(x)) = _data;
82 x = hypre_ParVectorCreate(y.x -> comm, y.x -> global_size,
84 hypre_ParVectorInitialize(x);
85 hypre_ParVectorSetPartitioningOwner(x,0);
86 hypre_ParVectorSetDataOwner(x,1);
87 hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),1);
97 x = hypre_ParVectorInDomainOf(const_cast<HypreParMatrix&>(A));
101 x = hypre_ParVectorInRangeOf(const_cast<HypreParMatrix&>(A));
109 x = (hypre_ParVector *) y;
118 hypre_ParVectorInitialize(x);
119 hypre_ParVectorSetPartitioningOwner(x,0);
121 hypre_ParVectorSetDataOwner(x,1);
122 hypre_SeqVectorSetDataOwner(hypre_ParVectorLocalVector(x),1);
129 hypre_Vector *hv = hypre_ParVectorToVectorAll(*
this);
132 hypre_SeqVectorSetDataOwner(hv,0);
133 hypre_SeqVectorDestroy(hv);
139 hypre_ParVectorSetConstantValues(x,d);
146 if (size != y.
Size())
152 for (
int i = 0; i <
size; i++)
161 Vector::data = hypre_VectorData(hypre_ParVectorLocalVector(x)) = _data;
166 return hypre_ParVectorSetRandomValues(x,seed);
171 hypre_ParVectorPrint(x,fname);
178 hypre_ParVectorDestroy(x);
182 #ifdef MFEM_USE_SUNDIALS
186 MFEM_ASSERT(nv && N_VGetVectorID(nv) == SUNDIALS_NVEC_PARHYP,
188 N_VectorContent_ParHyp nv_c = (N_VectorContent_ParHyp)(nv->content);
189 MFEM_ASSERT(nv_c->own_parvector == FALSE,
"invalid N_Vector");
190 nv_c->local_length = x->local_vector->size;
191 nv_c->global_length = x->global_size;
192 nv_c->comm = x->comm;
196 #endif // MFEM_USE_SUNDIALS
201 return hypre_ParVectorInnerProd(*x, *y);
206 return hypre_ParVectorInnerProd(x, y);
215 double loc_norm = vec.
Norml1();
216 MPI_Allreduce(&loc_norm, &norm, 1, MPI_DOUBLE, MPI_SUM, comm);
220 double loc_norm = vec*vec;
221 MPI_Allreduce(&loc_norm, &norm, 1, MPI_DOUBLE, MPI_SUM, comm);
224 if (p < std::numeric_limits<double>::infinity())
227 for (
int i = 0; i < vec.
Size(); i++)
229 sum += pow(fabs(vec(i)), p);
231 MPI_Allreduce(&sum, &norm, 1, MPI_DOUBLE, MPI_SUM, comm);
232 norm = pow(norm, 1.0/p);
237 MPI_Allreduce(&loc_norm, &norm, 1, MPI_DOUBLE, MPI_MAX, comm);
243 void HypreParMatrix::Init()
247 diagOwner = offdOwner = colMapOwner = -1;
257 char HypreParMatrix::CopyCSR(
SparseMatrix *csr, hypre_CSRMatrix *hypre_csr)
259 hypre_CSRMatrixData(hypre_csr) = csr->
GetData();
261 hypre_CSRMatrixI(hypre_csr) = csr->
GetI();
262 hypre_CSRMatrixJ(hypre_csr) = csr->
GetJ();
266 hypre_CSRMatrixI(hypre_csr) =
267 DuplicateAs<HYPRE_Int>(csr->
GetI(), csr->
Height()+1);
268 hypre_CSRMatrixJ(hypre_csr) =
275 char HypreParMatrix::CopyBoolCSR(Table *bool_csr, hypre_CSRMatrix *hypre_csr)
277 int nnz = bool_csr->Size_of_connections();
278 double *data =
new double[nnz];
279 for (
int i = 0; i < nnz; i++)
283 hypre_CSRMatrixData(hypre_csr) = data;
285 hypre_CSRMatrixI(hypre_csr) = bool_csr->GetI();
286 hypre_CSRMatrixJ(hypre_csr) = bool_csr->GetJ();
290 hypre_CSRMatrixI(hypre_csr) =
291 DuplicateAs<HYPRE_Int>(bool_csr->GetI(), bool_csr->Size()+1);
292 hypre_CSRMatrixJ(hypre_csr) =
293 DuplicateAs<HYPRE_Int>(bool_csr->GetJ(), nnz);
299 void HypreParMatrix::CopyCSR_J(hypre_CSRMatrix *hypre_csr,
int *J)
301 HYPRE_Int nnz = hypre_CSRMatrixNumNonzeros(hypre_csr);
302 for (HYPRE_Int j = 0; j < nnz; j++)
304 J[j] = int(hypre_CSRMatrixJ(hypre_csr)[j]);
311 :
Operator(diag->Height(), diag->Width())
314 A = hypre_ParCSRMatrixCreate(comm, glob_size, glob_size, row_starts,
316 hypre_ParCSRMatrixSetDataOwner(A,1);
317 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
318 hypre_ParCSRMatrixSetColStartsOwner(A,0);
320 hypre_CSRMatrixSetDataOwner(A->diag,0);
321 diagOwner = CopyCSR(diag, A->diag);
322 hypre_CSRMatrixSetRownnz(A->diag);
324 hypre_CSRMatrixSetDataOwner(A->offd,1);
325 hypre_CSRMatrixI(A->offd) = hypre_CTAlloc(HYPRE_Int, diag->
Height()+1);
332 hypre_ParCSRMatrixSetNumNonzeros(A);
335 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
337 CopyCSR_J(A->diag, diag->
GetJ());
340 hypre_MatvecCommPkgCreate(A);
345 HYPRE_Int global_num_rows,
346 HYPRE_Int global_num_cols,
347 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
349 :
Operator(diag->Height(), diag->Width())
352 A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
353 row_starts, col_starts,
355 hypre_ParCSRMatrixSetDataOwner(A,1);
356 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
357 hypre_ParCSRMatrixSetColStartsOwner(A,0);
359 hypre_CSRMatrixSetDataOwner(A->diag,0);
360 diagOwner = CopyCSR(diag, A->diag);
361 hypre_CSRMatrixSetRownnz(A->diag);
363 hypre_CSRMatrixSetDataOwner(A->offd,1);
364 hypre_CSRMatrixI(A->offd) = hypre_CTAlloc(HYPRE_Int, diag->
Height()+1);
366 hypre_ParCSRMatrixSetNumNonzeros(A);
369 if (row_starts == col_starts)
371 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
373 CopyCSR_J(A->diag, diag->
GetJ());
377 hypre_MatvecCommPkgCreate(A);
382 HYPRE_Int global_num_rows,
383 HYPRE_Int global_num_cols,
384 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
387 :
Operator(diag->Height(), diag->Width())
390 A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
391 row_starts, col_starts,
394 hypre_ParCSRMatrixSetDataOwner(A,1);
395 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
396 hypre_ParCSRMatrixSetColStartsOwner(A,0);
398 hypre_CSRMatrixSetDataOwner(A->diag,0);
399 diagOwner = CopyCSR(diag, A->diag);
400 hypre_CSRMatrixSetRownnz(A->diag);
402 hypre_CSRMatrixSetDataOwner(A->offd,0);
403 offdOwner = CopyCSR(offd, A->offd);
404 hypre_CSRMatrixSetRownnz(A->offd);
406 hypre_ParCSRMatrixColMapOffd(A) = cmap;
410 hypre_ParCSRMatrixSetNumNonzeros(A);
413 if (row_starts == col_starts)
415 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
417 CopyCSR_J(A->diag, diag->
GetJ());
421 hypre_MatvecCommPkgCreate(A);
427 HYPRE_Int global_num_rows, HYPRE_Int global_num_cols,
428 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
429 HYPRE_Int *diag_i, HYPRE_Int *diag_j,
double *diag_data,
430 HYPRE_Int *offd_i, HYPRE_Int *offd_j,
double *offd_data,
431 HYPRE_Int offd_num_cols, HYPRE_Int *offd_col_map)
434 A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
435 row_starts, col_starts, offd_num_cols, 0, 0);
436 hypre_ParCSRMatrixSetDataOwner(A,1);
437 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
438 hypre_ParCSRMatrixSetColStartsOwner(A,0);
440 HYPRE_Int local_num_rows = hypre_CSRMatrixNumRows(A->diag);
442 hypre_CSRMatrixSetDataOwner(A->diag,0);
443 hypre_CSRMatrixI(A->diag) = diag_i;
444 hypre_CSRMatrixJ(A->diag) = diag_j;
445 hypre_CSRMatrixData(A->diag) = diag_data;
446 hypre_CSRMatrixNumNonzeros(A->diag) = diag_i[local_num_rows];
447 hypre_CSRMatrixSetRownnz(A->diag);
451 hypre_CSRMatrixSetDataOwner(A->offd,0);
452 hypre_CSRMatrixI(A->offd) = offd_i;
453 hypre_CSRMatrixJ(A->offd) = offd_j;
454 hypre_CSRMatrixData(A->offd) = offd_data;
455 hypre_CSRMatrixNumNonzeros(A->offd) = offd_i[local_num_rows];
456 hypre_CSRMatrixSetRownnz(A->offd);
460 hypre_ParCSRMatrixColMapOffd(A) = offd_col_map;
464 hypre_ParCSRMatrixSetNumNonzeros(A);
467 if (row_starts == col_starts)
469 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
472 hypre_MatvecCommPkgCreate(A);
480 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
483 MFEM_ASSERT(sm_a != NULL,
"invalid input");
484 MFEM_VERIFY(!HYPRE_AssumedPartitionCheck(),
485 "this method can not be used with assumed partition");
489 hypre_CSRMatrix *csr_a;
490 csr_a = hypre_CSRMatrixCreate(sm_a ->
Height(), sm_a ->
Width(),
491 sm_a -> NumNonZeroElems());
493 hypre_CSRMatrixSetDataOwner(csr_a,0);
494 CopyCSR(sm_a, csr_a);
495 hypre_CSRMatrixSetRownnz(csr_a);
497 A = hypre_CSRMatrixToParCSRMatrix(comm, csr_a, row_starts, col_starts);
500 delete [] hypre_CSRMatrixI(csr_a);
501 delete [] hypre_CSRMatrixJ(csr_a);
503 hypre_CSRMatrixI(csr_a) = NULL;
504 hypre_CSRMatrixDestroy(csr_a);
510 if (row_starts == col_starts)
512 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
515 hypre_MatvecCommPkgCreate(A);
520 HYPRE_Int global_num_rows,
521 HYPRE_Int global_num_cols,
522 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
527 A = hypre_ParCSRMatrixCreate(comm, global_num_rows, global_num_cols,
528 row_starts, col_starts, 0, nnz, 0);
529 hypre_ParCSRMatrixSetDataOwner(A,1);
530 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
531 hypre_ParCSRMatrixSetColStartsOwner(A,0);
533 hypre_CSRMatrixSetDataOwner(A->diag,0);
534 diagOwner = CopyBoolCSR(diag, A->diag);
535 hypre_CSRMatrixSetRownnz(A->diag);
537 hypre_CSRMatrixSetDataOwner(A->offd,1);
538 hypre_CSRMatrixI(A->offd) = hypre_CTAlloc(HYPRE_Int, diag->
Size()+1);
540 hypre_ParCSRMatrixSetNumNonzeros(A);
543 if (row_starts == col_starts)
545 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
547 CopyCSR_J(A->diag, diag->
GetJ());
551 hypre_MatvecCommPkgCreate(A);
559 HYPRE_Int *row, HYPRE_Int *col,
560 HYPRE_Int *i_diag, HYPRE_Int *j_diag,
561 HYPRE_Int *i_offd, HYPRE_Int *j_offd,
562 HYPRE_Int *cmap, HYPRE_Int cmap_size)
564 HYPRE_Int diag_nnz, offd_nnz;
567 if (HYPRE_AssumedPartitionCheck())
569 diag_nnz = i_diag[row[1]-row[0]];
570 offd_nnz = i_offd[row[1]-row[0]];
572 A = hypre_ParCSRMatrixCreate(comm, row[2], col[2], row, col,
573 cmap_size, diag_nnz, offd_nnz);
577 diag_nnz = i_diag[row[
id+1]-row[id]];
578 offd_nnz = i_offd[row[
id+1]-row[id]];
580 A = hypre_ParCSRMatrixCreate(comm, row[np], col[np], row, col,
581 cmap_size, diag_nnz, offd_nnz);
584 hypre_ParCSRMatrixSetDataOwner(A,1);
585 hypre_ParCSRMatrixSetRowStartsOwner(A,0);
586 hypre_ParCSRMatrixSetColStartsOwner(A,0);
590 double *a_diag =
new double[diag_nnz];
591 for (i = 0; i < diag_nnz; i++)
596 double *a_offd =
new double[offd_nnz];
597 for (i = 0; i < offd_nnz; i++)
602 hypre_CSRMatrixSetDataOwner(A->diag,0);
603 hypre_CSRMatrixI(A->diag) = i_diag;
604 hypre_CSRMatrixJ(A->diag) = j_diag;
605 hypre_CSRMatrixData(A->diag) = a_diag;
606 hypre_CSRMatrixSetRownnz(A->diag);
610 hypre_CSRMatrixSetDataOwner(A->offd,0);
611 hypre_CSRMatrixI(A->offd) = i_offd;
612 hypre_CSRMatrixJ(A->offd) = j_offd;
613 hypre_CSRMatrixData(A->offd) = a_offd;
614 hypre_CSRMatrixSetRownnz(A->offd);
618 hypre_ParCSRMatrixColMapOffd(A) = cmap;
622 hypre_ParCSRMatrixSetNumNonzeros(A);
627 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
630 hypre_MatvecCommPkgCreate(A);
639 HYPRE_Int glob_ncols,
int *I, HYPRE_Int *J,
640 double *data, HYPRE_Int *rows, HYPRE_Int *cols)
646 HYPRE_Int my_col_start, my_col_end;
647 if (HYPRE_AssumedPartitionCheck())
650 my_col_start = cols[0];
651 my_col_end = cols[1];
656 MPI_Comm_rank(comm, &myid);
657 MPI_Comm_size(comm, &part_size);
659 my_col_start = cols[myid];
660 my_col_end = cols[myid+1];
664 HYPRE_Int *row_starts, *col_starts;
667 row_starts = col_starts = hypre_TAlloc(HYPRE_Int, part_size);
668 for (
int i = 0; i < part_size; i++)
670 row_starts[i] = rows[i];
675 row_starts = hypre_TAlloc(HYPRE_Int, part_size);
676 col_starts = hypre_TAlloc(HYPRE_Int, part_size);
677 for (
int i = 0; i < part_size; i++)
679 row_starts[i] = rows[i];
680 col_starts[i] = cols[i];
686 HYPRE_Int diag_nnz = 0, offd_nnz = 0, offd_num_cols = 0;
687 map<HYPRE_Int, HYPRE_Int> offd_map;
688 for (HYPRE_Int j = 0, loc_nnz = I[nrows]; j < loc_nnz; j++)
690 HYPRE_Int glob_col = J[j];
691 if (my_col_start <= glob_col && glob_col < my_col_end)
697 offd_map.insert(pair<const HYPRE_Int, HYPRE_Int>(glob_col, -1));
702 for (map<HYPRE_Int, HYPRE_Int>::iterator it = offd_map.begin();
703 it != offd_map.end(); ++it)
705 it->second = offd_num_cols++;
709 A = hypre_ParCSRMatrixCreate(comm, glob_nrows, glob_ncols,
710 row_starts, col_starts, offd_num_cols,
712 hypre_ParCSRMatrixInitialize(A);
714 HYPRE_Int *diag_i, *diag_j, *offd_i, *offd_j, *offd_col_map;
715 double *diag_data, *offd_data;
718 diag_data = A->diag->data;
721 offd_data = A->offd->data;
722 offd_col_map = A->col_map_offd;
724 diag_nnz = offd_nnz = 0;
725 for (HYPRE_Int i = 0, j = 0; i < nrows; i++)
727 diag_i[i] = diag_nnz;
728 offd_i[i] = offd_nnz;
729 for (HYPRE_Int j_end = I[i+1]; j < j_end; j++)
731 HYPRE_Int glob_col = J[j];
732 if (my_col_start <= glob_col && glob_col < my_col_end)
734 diag_j[diag_nnz] = glob_col - my_col_start;
735 diag_data[diag_nnz] = data[j];
740 offd_j[offd_nnz] = offd_map[glob_col];
741 offd_data[offd_nnz] = data[j];
746 diag_i[nrows] = diag_nnz;
747 offd_i[nrows] = offd_nnz;
748 for (map<HYPRE_Int, HYPRE_Int>::iterator it = offd_map.begin();
749 it != offd_map.end(); ++it)
751 offd_col_map[it->second] = it->first;
754 hypre_ParCSRMatrixSetNumNonzeros(A);
756 if (row_starts == col_starts)
758 hypre_CSRMatrixReorder(hypre_ParCSRMatrixDiag(A));
760 hypre_MatvecCommPkgCreate(A);
781 MFEM_ASSERT(diagOwner == -1 && offdOwner == -1 && colMapOwner == -1,
"");
782 MFEM_ASSERT(ParCSROwner,
"");
783 hypre_ParCSRMatrix *R = A;
792 if (!A || hypre_ParCSRMatrixOwnsRowStarts(A) ||
793 (hypre_ParCSRMatrixRowStarts(A) == hypre_ParCSRMatrixColStarts(A) &&
794 hypre_ParCSRMatrixOwnsColStarts(A)))
800 if (HYPRE_AssumedPartitionCheck())
806 MPI_Comm_size(hypre_ParCSRMatrixComm(A), &row_starts_size);
810 HYPRE_Int *old_row_starts = hypre_ParCSRMatrixRowStarts(A);
811 HYPRE_Int *new_row_starts = hypre_CTAlloc(HYPRE_Int, row_starts_size);
812 for (
int i = 0; i < row_starts_size; i++)
814 new_row_starts[i] = old_row_starts[i];
817 hypre_ParCSRMatrixRowStarts(A) = new_row_starts;
818 hypre_ParCSRMatrixOwnsRowStarts(A) = 1;
820 if (hypre_ParCSRMatrixColStarts(A) == old_row_starts)
822 hypre_ParCSRMatrixColStarts(A) = new_row_starts;
823 hypre_ParCSRMatrixOwnsColStarts(A) = 0;
829 if (!A || hypre_ParCSRMatrixOwnsColStarts(A) ||
830 (hypre_ParCSRMatrixRowStarts(A) == hypre_ParCSRMatrixColStarts(A) &&
831 hypre_ParCSRMatrixOwnsRowStarts(A)))
837 if (HYPRE_AssumedPartitionCheck())
843 MPI_Comm_size(hypre_ParCSRMatrixComm(A), &col_starts_size);
847 HYPRE_Int *old_col_starts = hypre_ParCSRMatrixColStarts(A);
848 HYPRE_Int *new_col_starts = hypre_CTAlloc(HYPRE_Int, col_starts_size);
849 for (
int i = 0; i < col_starts_size; i++)
851 new_col_starts[i] = old_col_starts[i];
854 hypre_ParCSRMatrixColStarts(A) = new_col_starts;
856 if (hypre_ParCSRMatrixRowStarts(A) == old_col_starts)
858 hypre_ParCSRMatrixRowStarts(A) = new_col_starts;
859 hypre_ParCSRMatrixOwnsRowStarts(A) = 1;
860 hypre_ParCSRMatrixOwnsColStarts(A) = 0;
864 hypre_ParCSRMatrixOwnsColStarts(A) = 1;
872 for (
int j = 0; j < size; j++)
874 diag(j) = A->diag->data[A->diag->i[j]];
875 MFEM_ASSERT(A->diag->j[A->diag->i[j]] == j,
876 "the first entry in each row must be the diagonal one");
880 static void MakeWrapper(
const hypre_CSRMatrix *mat,
SparseMatrix &wrapper)
882 HYPRE_Int nr = hypre_CSRMatrixNumRows(mat);
883 HYPRE_Int nc = hypre_CSRMatrixNumCols(mat);
886 hypre_CSRMatrixJ(mat),
887 hypre_CSRMatrixData(mat),
888 nr, nc,
false,
false,
false);
890 HYPRE_Int nnz = hypre_CSRMatrixNumNonzeros(mat);
891 SparseMatrix tmp(DuplicateAs<int>(hypre_CSRMatrixI(mat), nr+1),
892 DuplicateAs<int>(hypre_CSRMatrixJ(mat), nnz),
893 hypre_CSRMatrixData(mat),
894 nr, nc,
true,
false,
false);
901 MakeWrapper(A->diag, diag);
906 MakeWrapper(A->offd, offd);
907 cmap = A->col_map_offd;
911 bool interleaved_rows,
912 bool interleaved_cols)
const
917 hypre_ParCSRMatrix **hypre_blocks =
new hypre_ParCSRMatrix*[nr * nc];
918 internal::hypre_ParCSRMatrixSplit(A, nr, nc, hypre_blocks,
919 interleaved_rows, interleaved_cols);
921 for (
int i = 0; i < nr; i++)
923 for (
int j = 0; j < nc; j++)
929 delete [] hypre_blocks;
934 hypre_ParCSRMatrix * At;
935 hypre_ParCSRMatrixTranspose(A, &At, 1);
936 hypre_ParCSRMatrixSetNumNonzeros(At);
938 hypre_MatvecCommPkgCreate(At);
946 return hypre_ParCSRMatrixMatvec(a, A, x, b, y);
951 MFEM_ASSERT(x.
Size() ==
Width(),
"invalid x.Size() = " << x.
Size()
952 <<
", expected size = " <<
Width());
953 MFEM_ASSERT(y.
Size() ==
Height(),
"invalid y.Size() = " << y.
Size()
954 <<
", expected size = " <<
Height());
973 hypre_ParCSRMatrixMatvec(a, A, *X, b, *Y);
977 double b,
Vector &y)
const
979 MFEM_ASSERT(x.
Size() ==
Height(),
"invalid x.Size() = " << x.
Size()
980 <<
", expected size = " <<
Height());
981 MFEM_ASSERT(y.
Size() ==
Width(),
"invalid y.Size() = " << y.
Size()
982 <<
", expected size = " <<
Width());
1003 hypre_ParCSRMatrixMatvecT(a, A, *Y, b, *X);
1009 return hypre_ParCSRMatrixMatvec(a, A, (hypre_ParVector *) x, b,
1010 (hypre_ParVector *) y);
1016 return hypre_ParCSRMatrixMatvecT(a, A, x, b, y);
1020 HYPRE_Int* row_starts)
const
1022 const bool assumed_partition = HYPRE_AssumedPartitionCheck();
1023 const bool same_rows = (D.
Height() == hypre_CSRMatrixNumRows(A->diag));
1026 if (assumed_partition)
1032 MPI_Comm_size(
GetComm(), &part_size);
1036 HYPRE_Int global_num_rows;
1039 row_starts = hypre_ParCSRMatrixRowStarts(A);
1040 global_num_rows = hypre_ParCSRMatrixGlobalNumRows(A);
1044 MFEM_VERIFY(row_starts != NULL,
"the number of rows in D and A is not "
1045 "the same; row_starts must be given (not NULL)");
1050 assumed_partition ? row_starts[2] : row_starts[part_size-1];
1053 HYPRE_Int *col_starts = hypre_ParCSRMatrixColStarts(A);
1054 HYPRE_Int *col_map_offd;
1059 GetOffd(A_offd, col_map_offd);
1067 global_num_rows, hypre_ParCSRMatrixGlobalNumCols(A),
1068 DuplicateAs<HYPRE_Int>(row_starts, part_size,
false),
1069 DuplicateAs<HYPRE_Int>(col_starts, part_size,
false),
1071 DuplicateAs<HYPRE_Int>(col_map_offd, A_offd.
Width()));
1076 #ifndef HYPRE_BIGINT
1087 hypre_ParCSRMatrixSetRowStartsOwner(DA->A, 1);
1088 hypre_ParCSRMatrixSetColStartsOwner(DA->A, 1);
1090 DA->diagOwner = DA->offdOwner = 3;
1091 DA->colMapOwner = 1;
1098 if (hypre_CSRMatrixNumRows(A->diag) != hypre_CSRMatrixNumRows(A->offd))
1103 if (hypre_CSRMatrixNumRows(A->diag) != diag.
Size())
1105 mfem_error(
"Note the Vector diag is not of compatible dimensions with A\n");
1109 double *Adiag_data = hypre_CSRMatrixData(A->diag);
1110 HYPRE_Int *Adiag_i = hypre_CSRMatrixI(A->diag);
1113 double *Aoffd_data = hypre_CSRMatrixData(A->offd);
1114 HYPRE_Int *Aoffd_i = hypre_CSRMatrixI(A->offd);
1117 for (
int i(0); i < size; ++i)
1120 for (jj = Adiag_i[i]; jj < Adiag_i[i+1]; ++jj)
1122 Adiag_data[jj] *= val;
1124 for (jj = Aoffd_i[i]; jj < Aoffd_i[i+1]; ++jj)
1126 Aoffd_data[jj] *= val;
1133 if (hypre_CSRMatrixNumRows(A->diag) != hypre_CSRMatrixNumRows(A->offd))
1138 if (hypre_CSRMatrixNumRows(A->diag) != diag.
Size())
1140 mfem_error(
"Note the Vector diag is not of compatible dimensions with A\n");
1144 double *Adiag_data = hypre_CSRMatrixData(A->diag);
1145 HYPRE_Int *Adiag_i = hypre_CSRMatrixI(A->diag);
1148 double *Aoffd_data = hypre_CSRMatrixData(A->offd);
1149 HYPRE_Int *Aoffd_i = hypre_CSRMatrixI(A->offd);
1152 for (
int i(0); i < size; ++i)
1157 mfem_error(
"HypreParMatrix::InvDiagScale : Division by 0");
1161 for (jj = Adiag_i[i]; jj < Adiag_i[i+1]; ++jj)
1163 Adiag_data[jj] *= val;
1165 for (jj = Aoffd_i[i]; jj < Aoffd_i[i+1]; ++jj)
1167 Aoffd_data[jj] *= val;
1174 if (hypre_CSRMatrixNumRows(A->diag) != hypre_CSRMatrixNumRows(A->offd))
1179 HYPRE_Int size=hypre_CSRMatrixNumRows(A->diag);
1182 double *Adiag_data = hypre_CSRMatrixData(A->diag);
1183 HYPRE_Int *Adiag_i = hypre_CSRMatrixI(A->diag);
1184 for (jj = 0; jj < Adiag_i[size]; ++jj)
1186 Adiag_data[jj] *= s;
1189 double *Aoffd_data = hypre_CSRMatrixData(A->offd);
1190 HYPRE_Int *Aoffd_i = hypre_CSRMatrixI(A->offd);
1191 for (jj = 0; jj < Aoffd_i[size]; ++jj)
1193 Aoffd_data[jj] *= s;
1197 static void get_sorted_rows_cols(
const Array<int> &rows_cols,
1202 for (
int i = 0; i < rows_cols.
Size(); i++)
1204 hypre_sorted[i] = rows_cols[i];
1205 if (i && rows_cols[i-1] > rows_cols[i]) { sorted =
false; }
1207 if (!sorted) { hypre_sorted.
Sort(); }
1216 hypre_CSRMatrix * csr_A;
1217 hypre_CSRMatrix * csr_A_wo_z;
1218 hypre_ParCSRMatrix * parcsr_A_ptr;
1219 HYPRE_Int * row_starts = NULL; HYPRE_Int * col_starts = NULL;
1220 HYPRE_Int row_start = -1; HYPRE_Int row_end = -1;
1221 HYPRE_Int col_start = -1; HYPRE_Int col_end = -1;
1223 comm = hypre_ParCSRMatrixComm(A);
1225 MPI_Comm_size(comm, &num_procs);
1227 ierr += hypre_ParCSRMatrixGetLocalRange(A,
1228 &row_start,&row_end,
1229 &col_start,&col_end );
1231 row_starts = hypre_ParCSRMatrixRowStarts(A);
1232 col_starts = hypre_ParCSRMatrixColStarts(A);
1234 parcsr_A_ptr = hypre_ParCSRMatrixCreate(comm,row_starts[num_procs],
1235 col_starts[num_procs],row_starts,
1238 csr_A = hypre_MergeDiagAndOffd(A);
1240 csr_A_wo_z = hypre_CSRMatrixDeleteZeros(csr_A,threshold);
1244 if (csr_A_wo_z == NULL)
1250 ierr += hypre_CSRMatrixDestroy(csr_A);
1253 ierr += GenerateDiagAndOffd(csr_A_wo_z,parcsr_A_ptr,
1256 ierr += hypre_CSRMatrixDestroy(csr_A_wo_z);
1258 ierr += hypre_ParCSRMatrixDestroy(A);
1268 get_sorted_rows_cols(rows_cols, rc_sorted);
1270 internal::hypre_ParCSRMatrixEliminateAXB(
1277 get_sorted_rows_cols(rows_cols, rc_sorted);
1279 hypre_ParCSRMatrix* Ae;
1280 internal::hypre_ParCSRMatrixEliminateAAe(
1288 hypre_ParCSRMatrixPrintIJ(A,offi,offj,fname);
1296 HYPRE_Int base_i, base_j;
1297 hypre_ParCSRMatrixReadIJ(comm, fname, &base_i, &base_j, &A);
1298 hypre_ParCSRMatrixSetNumNonzeros(A);
1300 hypre_MatvecCommPkgCreate(A);
1306 void HypreParMatrix::Destroy()
1308 if ( X != NULL ) {
delete X; }
1309 if ( Y != NULL ) {
delete Y; }
1311 if (A == NULL) {
return; }
1317 delete [] hypre_CSRMatrixI(A->diag);
1318 delete [] hypre_CSRMatrixJ(A->diag);
1320 hypre_CSRMatrixI(A->diag) = NULL;
1321 hypre_CSRMatrixJ(A->diag) = NULL;
1324 delete [] hypre_CSRMatrixData(A->diag);
1326 hypre_CSRMatrixData(A->diag) = NULL;
1332 delete [] hypre_CSRMatrixI(A->offd);
1333 delete [] hypre_CSRMatrixJ(A->offd);
1335 hypre_CSRMatrixI(A->offd) = NULL;
1336 hypre_CSRMatrixJ(A->offd) = NULL;
1339 delete [] hypre_CSRMatrixData(A->offd);
1341 hypre_CSRMatrixData(A->offd) = NULL;
1343 if (colMapOwner >= 0)
1345 if (colMapOwner & 1)
1347 delete [] hypre_ParCSRMatrixColMapOffd(A);
1349 hypre_ParCSRMatrixColMapOffd(A) = NULL;
1354 hypre_ParCSRMatrixDestroy(A);
1361 hypre_ParCSRMatrix *C_hypre =
1362 internal::hypre_ParCSRMatrixAdd(const_cast<HypreParMatrix &>(A),
1363 const_cast<HypreParMatrix &>(B));
1364 MFEM_VERIFY(C_hypre,
"error in hypre_ParCSRMatrixAdd");
1366 hypre_MatvecCommPkgCreate(C_hypre);
1377 hypre_ParCSRMatrix * ab;
1378 ab = hypre_ParMatmul(*A,*B);
1380 hypre_MatvecCommPkgCreate(ab);
1387 HYPRE_Int P_owns_its_col_starts =
1388 hypre_ParCSRMatrixOwnsColStarts((hypre_ParCSRMatrix*)(*P));
1390 hypre_ParCSRMatrix * rap;
1391 hypre_BoomerAMGBuildCoarseOperator(*P,*A,*P,&rap);
1392 hypre_ParCSRMatrixSetNumNonzeros(rap);
1397 hypre_ParCSRMatrixSetRowStartsOwner(rap,0);
1398 hypre_ParCSRMatrixSetColStartsOwner(rap,0);
1400 if (P_owns_its_col_starts)
1402 hypre_ParCSRMatrixSetColStartsOwner(*P, 1);
1410 HYPRE_Int P_owns_its_col_starts =
1411 hypre_ParCSRMatrixOwnsColStarts((hypre_ParCSRMatrix*)(*P));
1412 HYPRE_Int Rt_owns_its_col_starts =
1413 hypre_ParCSRMatrixOwnsColStarts((hypre_ParCSRMatrix*)(*Rt));
1415 hypre_ParCSRMatrix * rap;
1416 hypre_BoomerAMGBuildCoarseOperator(*Rt,*A,*P,&rap);
1418 hypre_ParCSRMatrixSetNumNonzeros(rap);
1420 if (!P_owns_its_col_starts)
1424 hypre_ParCSRMatrixSetColStartsOwner(rap,0);
1426 if (!Rt_owns_its_col_starts)
1430 hypre_ParCSRMatrixSetRowStartsOwner(rap,0);
1440 Ae.
Mult(-1.0, X, 1.0, B);
1442 hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag((hypre_ParCSRMatrix *)A);
1443 double *data = hypre_CSRMatrixData(A_diag);
1444 HYPRE_Int *I = hypre_CSRMatrixI(A_diag);
1446 HYPRE_Int *J = hypre_CSRMatrixJ(A_diag);
1447 hypre_CSRMatrix *A_offd = hypre_ParCSRMatrixOffd((hypre_ParCSRMatrix *)A);
1448 HYPRE_Int *I_offd = hypre_CSRMatrixI(A_offd);
1449 double *data_offd = hypre_CSRMatrixData(A_offd);
1452 for (
int i = 0; i < ess_dof_list.
Size(); i++)
1454 int r = ess_dof_list[i];
1455 B(r) = data[I[r]] * X(r);
1462 MFEM_ABORT(
"the diagonal entry must be the first entry in the row!");
1464 for (
int j = I[r]+1; j < I[r+1]; j++)
1468 MFEM_ABORT(
"all off-diagonal entries must be zero!");
1471 for (
int j = I_offd[r]; j < I_offd[r+1]; j++)
1473 if (data_offd[j] != 0.0)
1475 MFEM_ABORT(
"all off-diagonal entries must be zero!");
1495 hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag(A);
1496 HYPRE_Int num_rows = hypre_CSRMatrixNumRows(A_diag);
1498 double *u_data = hypre_VectorData(hypre_ParVectorLocalVector(u));
1499 double *r_data = hypre_VectorData(hypre_ParVectorLocalVector(r));
1501 for (
int i = 0; i < N; i++)
1504 hypre_ParVectorCopy(f, r);
1505 hypre_ParCSRMatrixMatvec(-1.0, A, u, 1.0, r);
1508 (0 == (i % 2)) ? coef = lambda : coef = mu;
1510 for (HYPRE_Int j = 0; j < num_rows; j++)
1512 u_data[j] += coef*r_data[j] / max_eig;
1528 hypre_ParVector *x0,
1529 hypre_ParVector *x1,
1530 hypre_ParVector *x2,
1531 hypre_ParVector *x3)
1534 hypre_CSRMatrix *A_diag = hypre_ParCSRMatrixDiag(A);
1535 HYPRE_Int num_rows = hypre_CSRMatrixNumRows(A_diag);
1537 double *u_data = hypre_VectorData(hypre_ParVectorLocalVector(u));
1539 double *x0_data = hypre_VectorData(hypre_ParVectorLocalVector(x0));
1540 double *x1_data = hypre_VectorData(hypre_ParVectorLocalVector(x1));
1541 double *x2_data = hypre_VectorData(hypre_ParVectorLocalVector(x2));
1542 double *x3_data = hypre_VectorData(hypre_ParVectorLocalVector(x3));
1544 hypre_ParVectorCopy(u, x0);
1547 hypre_ParVectorCopy(f, x1);
1548 hypre_ParCSRMatrixMatvec(-1.0, A, x0, 1.0, x1);
1550 for (HYPRE_Int i = 0; i < num_rows; i++)
1552 x1_data[i] /= -max_eig;
1556 for (HYPRE_Int i = 0; i < num_rows; i++)
1558 x1_data[i] = x0_data[i] -x1_data[i];
1562 for (HYPRE_Int i = 0; i < num_rows; i++)
1564 x3_data[i] = fir_coeffs[0]*x0_data[i] +fir_coeffs[1]*x1_data[i];
1567 for (
int n = 2; n <= poly_order; n++)
1570 hypre_ParVectorCopy(f, x2);
1571 hypre_ParCSRMatrixMatvec(-1.0, A, x1, 1.0, x2);
1573 for (HYPRE_Int i = 0; i < num_rows; i++)
1575 x2_data[i] /= -max_eig;
1583 for (HYPRE_Int i = 0; i < num_rows; i++)
1585 x2_data[i] = (x1_data[i]-x0_data[i]) +(x1_data[i]-2*x2_data[i]);
1586 x3_data[i] += fir_coeffs[n]*x2_data[i];
1587 x0_data[i] = x1_data[i];
1588 x1_data[i] = x2_data[i];
1592 for (HYPRE_Int i = 0; i < num_rows; i++)
1594 u_data[i] = x3_data[i];
1614 B =
X =
V =
Z = NULL;
1620 int _relax_times,
double _relax_weight,
double _omega,
1621 int _poly_order,
double _poly_fraction)
1632 B =
X =
V =
Z = NULL;
1641 type =
static_cast<int>(_type);
1667 double a = -1, b, c;
1668 if (!strcmp(name,
"Rectangular")) { a = 1.0, b = 0.0, c = 0.0; }
1669 if (!strcmp(name,
"Hanning")) { a = 0.5, b = 0.5, c = 0.0; }
1670 if (!strcmp(name,
"Hamming")) { a = 0.54, b = 0.46, c = 0.0; }
1671 if (!strcmp(name,
"Blackman")) { a = 0.42, b = 0.50, c = 0.08; }
1674 mfem_error(
"HypreSmoother::SetWindowByName : name not recognized!");
1692 mfem_error(
"HypreSmoother::SetOperator : not HypreParMatrix!");
1698 if (
B) {
delete B; }
1699 if (
X) {
delete X; }
1700 if (
V) {
delete V; }
1701 if (
Z) {
delete Z; }
1720 A->
Mult(ones, diag);
1729 for (
int i = 0; i <
height; i++)
1742 else if (
type == 1001 ||
type == 1002)
1773 double* window_coeffs =
new double[
poly_order+1];
1774 double* cheby_coeffs =
new double[
poly_order+1];
1782 window_coeffs[i] = a + b*cos(t) +c*cos(2*t);
1786 double theta_pb = acos(1.0 -0.5*k_pb);
1788 cheby_coeffs[0] = (theta_pb +sigma)/M_PI;
1791 double t = i*(theta_pb+sigma);
1792 cheby_coeffs[i] = 2.0*sin(t)/(i*M_PI);
1797 fir_coeffs[i] = window_coeffs[i]*cheby_coeffs[i];
1800 delete[] window_coeffs;
1801 delete[] cheby_coeffs;
1808 mfem_error(
"HypreSmoother::Mult (...) : HypreParMatrix A is missing");
1816 HYPRE_ParCSRDiagScale(NULL, *
A, b, x);
1840 else if (
type == 1002)
1855 hypre_ParCSRRelax(*
A, b,
type,
1860 hypre_ParCSRRelax(*
A, b,
type,
1871 mfem_error(
"HypreSmoother::Mult (...) : HypreParMatrix A is missing");
1877 A -> GetGlobalNumRows(),
1879 A -> GetRowStarts());
1881 A -> GetGlobalNumCols(),
1883 A -> GetColStarts());
1896 if (
B) {
delete B; }
1897 if (
X) {
delete X; }
1898 if (
V) {
delete V; }
1899 if (
Z) {
delete Z; }
1908 if (
X0) {
delete X0; }
1909 if (
X1) {
delete X1; }
1921 :
Solver(_A->Height(), _A->Width())
1932 mfem_error(
"HypreSolver::Mult (...) : HypreParMatrix A is missing");
1952 mfem_error(
"HypreSolver::Mult (...) : HypreParMatrix A is missing");
1958 A -> GetGlobalNumRows(),
1960 A -> GetRowStarts());
1962 A -> GetGlobalNumCols(),
1964 A -> GetColStarts());
1977 if (
B) {
delete B; }
1978 if (
X) {
delete X; }
1988 HYPRE_ParCSRMatrixGetComm(*
A, &comm);
1990 HYPRE_ParCSRPCGCreate(comm, &pcg_solver);
1995 HYPRE_PCGSetTol(pcg_solver, tol);
2000 HYPRE_PCGSetMaxIter(pcg_solver, max_iter);
2005 HYPRE_PCGSetLogging(pcg_solver, logging);
2010 HYPRE_ParCSRPCGSetPrintLevel(pcg_solver, print_lvl);
2015 HYPRE_ParCSRPCGSetPrecond(pcg_solver,
2023 HYPRE_PCGSetTwoNorm(pcg_solver, 1);
2024 if (res_frequency > 0)
2026 HYPRE_PCGSetRecomputeResidualP(pcg_solver, res_frequency);
2030 HYPRE_PCGSetResidualTol(pcg_solver, rtol);
2037 HYPRE_Int time_index = 0;
2038 HYPRE_Int num_iterations;
2039 double final_res_norm;
2041 HYPRE_Int print_level;
2043 HYPRE_PCGGetPrintLevel(pcg_solver, &print_level);
2044 HYPRE_ParCSRPCGSetPrintLevel(pcg_solver, print_level%3);
2046 HYPRE_ParCSRMatrixGetComm(*
A, &comm);
2050 if (print_level > 0 && print_level < 3)
2052 time_index = hypre_InitializeTiming(
"PCG Setup");
2053 hypre_BeginTiming(time_index);
2056 HYPRE_ParCSRPCGSetup(pcg_solver, *
A, b, x);
2059 if (print_level > 0 && print_level < 3)
2061 hypre_EndTiming(time_index);
2062 hypre_PrintTiming(
"Setup phase times", comm);
2063 hypre_FinalizeTiming(time_index);
2064 hypre_ClearTiming();
2068 if (print_level > 0 && print_level < 3)
2070 time_index = hypre_InitializeTiming(
"PCG Solve");
2071 hypre_BeginTiming(time_index);
2079 HYPRE_ParCSRPCGSolve(pcg_solver, *
A, b, x);
2081 if (print_level > 0)
2083 if (print_level < 3)
2085 hypre_EndTiming(time_index);
2086 hypre_PrintTiming(
"Solve phase times", comm);
2087 hypre_FinalizeTiming(time_index);
2088 hypre_ClearTiming();
2091 HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_iterations);
2092 HYPRE_ParCSRPCGGetFinalRelativeResidualNorm(pcg_solver,
2095 MPI_Comm_rank(comm, &myid);
2099 cout <<
"PCG Iterations = " << num_iterations << endl
2100 <<
"Final PCG Relative Residual Norm = " << final_res_norm
2104 HYPRE_ParCSRPCGSetPrintLevel(pcg_solver, print_level);
2109 HYPRE_ParCSRPCGDestroy(pcg_solver);
2123 HYPRE_ParCSRMatrixGetComm(*
A, &comm);
2125 HYPRE_ParCSRGMRESCreate(comm, &gmres_solver);
2126 HYPRE_ParCSRGMRESSetKDim(gmres_solver, k_dim);
2127 HYPRE_ParCSRGMRESSetMaxIter(gmres_solver, max_iter);
2128 HYPRE_ParCSRGMRESSetTol(gmres_solver, tol);
2133 HYPRE_GMRESSetTol(gmres_solver, tol);
2138 HYPRE_GMRESSetMaxIter(gmres_solver, max_iter);
2143 HYPRE_GMRESSetKDim(gmres_solver, k_dim);
2148 HYPRE_GMRESSetLogging(gmres_solver, logging);
2153 HYPRE_GMRESSetPrintLevel(gmres_solver, print_lvl);
2158 HYPRE_ParCSRGMRESSetPrecond(gmres_solver,
2167 HYPRE_Int time_index = 0;
2168 HYPRE_Int num_iterations;
2169 double final_res_norm;
2171 HYPRE_Int print_level;
2173 HYPRE_GMRESGetPrintLevel(gmres_solver, &print_level);
2175 HYPRE_ParCSRMatrixGetComm(*
A, &comm);
2179 if (print_level > 0)
2181 time_index = hypre_InitializeTiming(
"GMRES Setup");
2182 hypre_BeginTiming(time_index);
2185 HYPRE_ParCSRGMRESSetup(gmres_solver, *
A, b, x);
2188 if (print_level > 0)
2190 hypre_EndTiming(time_index);
2191 hypre_PrintTiming(
"Setup phase times", comm);
2192 hypre_FinalizeTiming(time_index);
2193 hypre_ClearTiming();
2197 if (print_level > 0)
2199 time_index = hypre_InitializeTiming(
"GMRES Solve");
2200 hypre_BeginTiming(time_index);
2208 HYPRE_ParCSRGMRESSolve(gmres_solver, *
A, b, x);
2210 if (print_level > 0)
2212 hypre_EndTiming(time_index);
2213 hypre_PrintTiming(
"Solve phase times", comm);
2214 hypre_FinalizeTiming(time_index);
2215 hypre_ClearTiming();
2217 HYPRE_ParCSRGMRESGetNumIterations(gmres_solver, &num_iterations);
2218 HYPRE_ParCSRGMRESGetFinalRelativeResidualNorm(gmres_solver,
2221 MPI_Comm_rank(comm, &myid);
2225 cout <<
"GMRES Iterations = " << num_iterations << endl
2226 <<
"Final GMRES Relative Residual Norm = " << final_res_norm
2234 HYPRE_ParCSRGMRESDestroy(gmres_solver);
2242 int sai_max_levels = 1;
2243 double sai_threshold = 0.1;
2244 double sai_filter = 0.1;
2246 double sai_loadbal = 0.0;
2248 int sai_logging = 1;
2250 HYPRE_ParCSRMatrixGetComm(A, &comm);
2252 HYPRE_ParaSailsCreate(comm, &sai_precond);
2253 HYPRE_ParaSailsSetParams(sai_precond, sai_threshold, sai_max_levels);
2254 HYPRE_ParaSailsSetFilter(sai_precond, sai_filter);
2255 HYPRE_ParaSailsSetSym(sai_precond, sai_sym);
2256 HYPRE_ParaSailsSetLoadbal(sai_precond, sai_loadbal);
2257 HYPRE_ParaSailsSetReuse(sai_precond, sai_reuse);
2258 HYPRE_ParaSailsSetLogging(sai_precond, sai_logging);
2263 HYPRE_ParaSailsSetSym(sai_precond, sym);
2268 HYPRE_ParaSailsDestroy(sai_precond);
2274 HYPRE_BoomerAMGCreate(&amg_precond);
2275 SetDefaultOptions();
2280 HYPRE_BoomerAMGCreate(&amg_precond);
2281 SetDefaultOptions();
2284 void HypreBoomerAMG::SetDefaultOptions()
2287 int coarsen_type = 10;
2289 double theta = 0.25;
2292 int interp_type = 6;
2297 int relax_sweeps = 1;
2300 int print_level = 1;
2301 int max_levels = 25;
2303 HYPRE_BoomerAMGSetCoarsenType(amg_precond, coarsen_type);
2304 HYPRE_BoomerAMGSetAggNumLevels(amg_precond, agg_levels);
2305 HYPRE_BoomerAMGSetRelaxType(amg_precond, relax_type);
2306 HYPRE_BoomerAMGSetNumSweeps(amg_precond, relax_sweeps);
2307 HYPRE_BoomerAMGSetStrongThreshold(amg_precond, theta);
2308 HYPRE_BoomerAMGSetInterpType(amg_precond, interp_type);
2309 HYPRE_BoomerAMGSetPMaxElmts(amg_precond, Pmax);
2310 HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level);
2311 HYPRE_BoomerAMGSetMaxLevels(amg_precond, max_levels);
2314 HYPRE_BoomerAMGSetMaxIter(amg_precond, 1);
2315 HYPRE_BoomerAMGSetTol(amg_precond, 0.0);
2318 void HypreBoomerAMG::ResetAMGPrecond()
2320 HYPRE_Int coarsen_type;
2321 HYPRE_Int agg_levels;
2322 HYPRE_Int relax_type;
2323 HYPRE_Int relax_sweeps;
2325 HYPRE_Int interp_type;
2327 HYPRE_Int print_level;
2329 HYPRE_Int nrbms = rbms.
Size();
2331 HYPRE_Int nodal_diag;
2332 HYPRE_Int relax_coarse;
2333 HYPRE_Int interp_vec_variant;
2335 HYPRE_Int smooth_interp_vectors;
2336 HYPRE_Int interp_refine;
2338 hypre_ParAMGData *amg_data = (hypre_ParAMGData *)amg_precond;
2341 HYPRE_BoomerAMGGetCoarsenType(amg_precond, &coarsen_type);
2342 agg_levels = hypre_ParAMGDataAggNumLevels(amg_data);
2343 relax_type = hypre_ParAMGDataUserRelaxType(amg_data);
2344 relax_sweeps = hypre_ParAMGDataUserNumSweeps(amg_data);
2345 HYPRE_BoomerAMGGetStrongThreshold(amg_precond, &theta);
2346 hypre_BoomerAMGGetInterpType(amg_precond, &interp_type);
2347 HYPRE_BoomerAMGGetPMaxElmts(amg_precond, &Pmax);
2348 HYPRE_BoomerAMGGetPrintLevel(amg_precond, &print_level);
2349 HYPRE_BoomerAMGGetNumFunctions(amg_precond, &dim);
2352 nodal = hypre_ParAMGDataNodal(amg_data);
2353 nodal_diag = hypre_ParAMGDataNodalDiag(amg_data);
2354 HYPRE_BoomerAMGGetCycleRelaxType(amg_precond, &relax_coarse, 3);
2355 interp_vec_variant = hypre_ParAMGInterpVecVariant(amg_data);
2356 q_max = hypre_ParAMGInterpVecQMax(amg_data);
2357 smooth_interp_vectors = hypre_ParAMGSmoothInterpVectors(amg_data);
2358 interp_refine = hypre_ParAMGInterpRefine(amg_data);
2361 HYPRE_BoomerAMGDestroy(amg_precond);
2362 HYPRE_BoomerAMGCreate(&amg_precond);
2364 HYPRE_BoomerAMGSetCoarsenType(amg_precond, coarsen_type);
2365 HYPRE_BoomerAMGSetAggNumLevels(amg_precond, agg_levels);
2366 HYPRE_BoomerAMGSetRelaxType(amg_precond, relax_type);
2367 HYPRE_BoomerAMGSetNumSweeps(amg_precond, relax_sweeps);
2368 HYPRE_BoomerAMGSetMaxLevels(amg_precond, 25);
2369 HYPRE_BoomerAMGSetTol(amg_precond, 0.0);
2370 HYPRE_BoomerAMGSetMaxIter(amg_precond, 1);
2371 HYPRE_BoomerAMGSetStrongThreshold(amg_precond, theta);
2372 HYPRE_BoomerAMGSetInterpType(amg_precond, interp_type);
2373 HYPRE_BoomerAMGSetPMaxElmts(amg_precond, Pmax);
2374 HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level);
2375 HYPRE_BoomerAMGSetNumFunctions(amg_precond, dim);
2378 HYPRE_BoomerAMGSetNodal(amg_precond, nodal);
2379 HYPRE_BoomerAMGSetNodalDiag(amg_precond, nodal_diag);
2380 HYPRE_BoomerAMGSetCycleRelaxType(amg_precond, relax_coarse, 3);
2381 HYPRE_BoomerAMGSetInterpVecVariant(amg_precond, interp_vec_variant);
2382 HYPRE_BoomerAMGSetInterpVecQMax(amg_precond, q_max);
2383 HYPRE_BoomerAMGSetSmoothInterpVectors(amg_precond, smooth_interp_vectors);
2384 HYPRE_BoomerAMGSetInterpRefine(amg_precond, interp_refine);
2386 HYPRE_BoomerAMGSetInterpVectors(amg_precond, rbms.
Size(), rbms.
GetData());
2393 MFEM_VERIFY(new_A,
"new Operator must be a HypreParMatrix!");
2395 if (
A) { ResetAMGPrecond(); }
2409 HYPRE_BoomerAMGSetNumFunctions(amg_precond, dim);
2412 HYPRE_BoomerAMGSetAggNumLevels(amg_precond, 0);
2413 HYPRE_BoomerAMGSetStrongThreshold(amg_precond, 0.5);
2419 y = 0.0; y(0) = x(1); y(1) = -x(0);
2421 static void func_ryz(
const Vector &x, Vector &y)
2423 y = 0.0; y(1) = x(2); y(2) = -x(1);
2425 static void func_rzx(
const Vector &x, Vector &y)
2427 y = 0.0; y(2) = x(0); y(0) = -x(2);
2430 void HypreBoomerAMG::RecomputeRBMs()
2433 Array<HypreParVector*> gf_rbms;
2436 for (
int i = 0; i < rbms.
Size(); i++)
2438 HYPRE_ParVectorDestroy(rbms[i]);
2445 VectorFunctionCoefficient coeff_rxy(2, func_rxy);
2447 ParGridFunction rbms_rxy(fespace);
2448 rbms_rxy.ProjectCoefficient(coeff_rxy);
2451 gf_rbms.SetSize(nrbms);
2452 gf_rbms[0] = rbms_rxy.ParallelAverage();
2458 VectorFunctionCoefficient coeff_rxy(3, func_rxy);
2459 VectorFunctionCoefficient coeff_ryz(3, func_ryz);
2460 VectorFunctionCoefficient coeff_rzx(3, func_rzx);
2462 ParGridFunction rbms_rxy(fespace);
2463 ParGridFunction rbms_ryz(fespace);
2464 ParGridFunction rbms_rzx(fespace);
2465 rbms_rxy.ProjectCoefficient(coeff_rxy);
2466 rbms_ryz.ProjectCoefficient(coeff_ryz);
2467 rbms_rzx.ProjectCoefficient(coeff_rzx);
2470 gf_rbms.SetSize(nrbms);
2471 gf_rbms[0] = rbms_rxy.ParallelAverage();
2472 gf_rbms[1] = rbms_ryz.ParallelAverage();
2473 gf_rbms[2] = rbms_rzx.ParallelAverage();
2482 for (
int i = 0; i < nrbms; i++)
2484 rbms[i] = gf_rbms[i]->StealParVector();
2492 this->fespace = fespace;
2502 int relax_coarse = 8;
2505 int interp_vec_variant = 2;
2507 int smooth_interp_vectors = 1;
2511 int interp_refine = 1;
2513 HYPRE_BoomerAMGSetNodal(amg_precond, nodal);
2514 HYPRE_BoomerAMGSetNodalDiag(amg_precond, nodal_diag);
2515 HYPRE_BoomerAMGSetCycleRelaxType(amg_precond, relax_coarse, 3);
2516 HYPRE_BoomerAMGSetInterpVecVariant(amg_precond, interp_vec_variant);
2517 HYPRE_BoomerAMGSetInterpVecQMax(amg_precond, q_max);
2518 HYPRE_BoomerAMGSetSmoothInterpVectors(amg_precond, smooth_interp_vectors);
2519 HYPRE_BoomerAMGSetInterpRefine(amg_precond, interp_refine);
2522 HYPRE_BoomerAMGSetInterpVectors(amg_precond, rbms.
Size(), rbms.
GetData());
2527 for (
int i = 0; i < rbms.
Size(); i++)
2529 HYPRE_ParVectorDestroy(rbms[i]);
2532 HYPRE_BoomerAMGDestroy(amg_precond);
2539 int cycle_type = 13;
2542 double rlx_weight = 1.0;
2543 double rlx_omega = 1.0;
2544 int amg_coarsen_type = 10;
2545 int amg_agg_levels = 1;
2546 int amg_rlx_type = 8;
2547 double theta = 0.25;
2548 int amg_interp_type = 6;
2555 bool trace_space, rt_trace_space;
2559 trace_space = trace_space || rt_trace_space;
2562 if (edge_fespace->
GetNE() > 0)
2567 if (dim == 2) { p++; }
2582 HYPRE_AMSCreate(&ams);
2584 HYPRE_AMSSetDimension(ams, sdim);
2585 HYPRE_AMSSetTol(ams, 0.0);
2586 HYPRE_AMSSetMaxIter(ams, 1);
2587 HYPRE_AMSSetCycleType(ams, cycle_type);
2588 HYPRE_AMSSetPrintLevel(ams, 1);
2610 for (
int i = 0; i < pmesh->
GetNV(); i++)
2612 coord = pmesh -> GetVertex(i);
2613 x_coord(i) = coord[0];
2614 y_coord(i) = coord[1];
2615 if (sdim == 3) { z_coord(i) = coord[2]; }
2622 HYPRE_AMSSetCoordinateVectors(ams, *x, *y, NULL);
2627 HYPRE_AMSSetCoordinateVectors(ams, *x, *y, *z);
2651 HYPRE_AMSSetDiscreteGradient(ams, *G);
2655 Pi = Pix = Piy = Piz = NULL;
2674 if (cycle_type < 10)
2682 Pix = Pi_blocks(0,0);
2683 Piy = Pi_blocks(0,1);
2684 if (sdim == 3) { Piz = Pi_blocks(0,2); }
2689 HYPRE_ParCSRMatrix HY_Pi = (Pi) ? (HYPRE_ParCSRMatrix) *Pi : NULL;
2690 HYPRE_ParCSRMatrix HY_Pix = (Pix) ? (HYPRE_ParCSRMatrix) *Pix : NULL;
2691 HYPRE_ParCSRMatrix HY_Piy = (Piy) ? (HYPRE_ParCSRMatrix) *Piy : NULL;
2692 HYPRE_ParCSRMatrix HY_Piz = (Piz) ? (HYPRE_ParCSRMatrix) *Piz : NULL;
2693 HYPRE_AMSSetInterpolations(ams, HY_Pi, HY_Pix, HY_Piy, HY_Piz);
2695 delete vert_fespace_d;
2698 delete vert_fespace;
2703 delete edge_fespace;
2708 HYPRE_AMSSetSmoothingOptions(ams, rlx_type, rlx_sweeps, rlx_weight, rlx_omega);
2709 HYPRE_AMSSetAlphaAMGOptions(ams, amg_coarsen_type, amg_agg_levels, amg_rlx_type,
2710 theta, amg_interp_type, amg_Pmax);
2711 HYPRE_AMSSetBetaAMGOptions(ams, amg_coarsen_type, amg_agg_levels, amg_rlx_type,
2712 theta, amg_interp_type, amg_Pmax);
2717 HYPRE_AMSDestroy(ams);
2732 HYPRE_AMSSetPrintLevel(ams, print_lvl);
2738 int cycle_type = 11;
2741 double rlx_weight = 1.0;
2742 double rlx_omega = 1.0;
2743 int amg_coarsen_type = 10;
2744 int amg_agg_levels = 1;
2745 int amg_rlx_type = 8;
2746 double theta = 0.25;
2747 int amg_interp_type = 6;
2749 int ams_cycle_type = 14;
2755 if (face_fespace->
GetNE() > 0)
2767 HYPRE_ADSCreate(&ads);
2769 HYPRE_ADSSetTol(ads, 0.0);
2770 HYPRE_ADSSetMaxIter(ads, 1);
2771 HYPRE_ADSSetCycleType(ads, cycle_type);
2772 HYPRE_ADSSetPrintLevel(ads, 1);
2800 for (
int i = 0; i < pmesh->
GetNV(); i++)
2802 coord = pmesh -> GetVertex(i);
2803 x_coord(i) = coord[0];
2804 y_coord(i) = coord[1];
2805 z_coord(i) = coord[2];
2810 HYPRE_ADSSetCoordinateVectors(ads, *x, *y, *z);
2834 HYPRE_ADSSetDiscreteCurl(ads, *C);
2853 HYPRE_ADSSetDiscreteGradient(ads, *G);
2857 RT_Pi = RT_Pix = RT_Piy = RT_Piz = NULL;
2858 ND_Pi = ND_Pix = ND_Piy = ND_Piz = NULL;
2877 if (ams_cycle_type < 10)
2887 ND_Pix = ND_Pi_blocks(0,0);
2888 ND_Piy = ND_Pi_blocks(0,1);
2889 ND_Piz = ND_Pi_blocks(0,2);
2907 if (cycle_type < 10)
2916 RT_Pix = RT_Pi_blocks(0,0);
2917 RT_Piy = RT_Pi_blocks(0,1);
2918 RT_Piz = RT_Pi_blocks(0,2);
2923 HYPRE_ParCSRMatrix HY_RT_Pi, HY_RT_Pix, HY_RT_Piy, HY_RT_Piz;
2924 HY_RT_Pi = (RT_Pi) ? (HYPRE_ParCSRMatrix) *RT_Pi : NULL;
2925 HY_RT_Pix = (RT_Pix) ? (HYPRE_ParCSRMatrix) *RT_Pix : NULL;
2926 HY_RT_Piy = (RT_Piy) ? (HYPRE_ParCSRMatrix) *RT_Piy : NULL;
2927 HY_RT_Piz = (RT_Piz) ? (HYPRE_ParCSRMatrix) *RT_Piz : NULL;
2928 HYPRE_ParCSRMatrix HY_ND_Pi, HY_ND_Pix, HY_ND_Piy, HY_ND_Piz;
2929 HY_ND_Pi = (ND_Pi) ? (HYPRE_ParCSRMatrix) *ND_Pi : NULL;
2930 HY_ND_Pix = (ND_Pix) ? (HYPRE_ParCSRMatrix) *ND_Pix : NULL;
2931 HY_ND_Piy = (ND_Piy) ? (HYPRE_ParCSRMatrix) *ND_Piy : NULL;
2932 HY_ND_Piz = (ND_Piz) ? (HYPRE_ParCSRMatrix) *ND_Piz : NULL;
2933 HYPRE_ADSSetInterpolations(ads,
2934 HY_RT_Pi, HY_RT_Pix, HY_RT_Piy, HY_RT_Piz,
2935 HY_ND_Pi, HY_ND_Pix, HY_ND_Piy, HY_ND_Piz);
2937 delete vert_fespace_d;
2941 delete vert_fespace;
2943 delete edge_fespace;
2946 HYPRE_ADSSetSmoothingOptions(ads, rlx_type, rlx_sweeps, rlx_weight, rlx_omega);
2947 HYPRE_ADSSetAMGOptions(ads, amg_coarsen_type, amg_agg_levels, amg_rlx_type,
2948 theta, amg_interp_type, amg_Pmax);
2949 HYPRE_ADSSetAMSOptions(ads, ams_cycle_type, amg_coarsen_type, amg_agg_levels,
2950 amg_rlx_type, theta, amg_interp_type, amg_Pmax);
2955 HYPRE_ADSDestroy(ads);
2977 HYPRE_ADSSetPrintLevel(ads, print_lvl);
2980 HypreLOBPCG::HypreMultiVector::HypreMultiVector(
int n,
HypreParVector & v,
2981 mv_InterfaceInterpreter & interpreter)
2985 mv_ptr = mv_MultiVectorCreateFromSampleVector(&interpreter, nv,
2986 (HYPRE_ParVector)v);
2988 HYPRE_ParVector* vecs = NULL;
2990 mv_TempMultiVector* tmp =
2991 (mv_TempMultiVector*)mv_MultiVectorGetData(mv_ptr);
2992 vecs = (HYPRE_ParVector*)(tmp -> vector);
2995 hpv =
new HypreParVector*[nv];
2996 for (
int i=0; i<nv; i++)
2998 hpv[i] =
new HypreParVector(vecs[i]);
3002 HypreLOBPCG::HypreMultiVector::~HypreMultiVector()
3006 for (
int i=0; i<nv; i++)
3013 mv_MultiVectorDestroy(mv_ptr);
3017 HypreLOBPCG::HypreMultiVector::Randomize(HYPRE_Int seed)
3019 mv_MultiVectorSetRandom(mv_ptr, seed);
3023 HypreLOBPCG::HypreMultiVector::GetVector(
unsigned int i)
3025 MFEM_ASSERT((
int)i < nv,
"index out of range");
3031 HypreLOBPCG::HypreMultiVector::StealVectors()
3033 HypreParVector ** hpv_ret = hpv;
3037 mv_TempMultiVector * mv_tmp =
3038 (mv_TempMultiVector*)mv_MultiVectorGetData(mv_ptr);
3040 mv_tmp->ownsVectors = 0;
3042 for (
int i=0; i<nv; i++)
3044 hpv_ret[i]->SetOwnership(1);
3061 MPI_Comm_size(comm,&numProcs);
3062 MPI_Comm_rank(comm,&myid);
3064 HYPRE_ParCSRSetupInterpreter(&interpreter);
3065 HYPRE_ParCSRSetupMatvec(&matvec_fn);
3066 HYPRE_LOBPCGCreate(&interpreter, &matvec_fn, &lobpcg_solver);
3075 HYPRE_LOBPCGDestroy(lobpcg_solver);
3081 HYPRE_LOBPCGSetTol(lobpcg_solver, tol);
3087 HYPRE_LOBPCGSetMaxIter(lobpcg_solver, max_iter);
3095 HYPRE_LOBPCGSetPrintLevel(lobpcg_solver, logging);
3102 HYPRE_LOBPCGSetPrecondUsageMode(lobpcg_solver, pcg_mode);
3108 HYPRE_LOBPCGSetPrecond(lobpcg_solver,
3109 (HYPRE_PtrToSolverFcn)this->PrecondSolve,
3110 (HYPRE_PtrToSolverFcn)this->PrecondSetup,
3111 (HYPRE_Solver)&precond);
3117 int locSize = A.
Width();
3119 if (HYPRE_AssumedPartitionCheck())
3121 part =
new HYPRE_Int[2];
3123 MPI_Scan(&locSize, &part[1], 1, HYPRE_MPI_INT, MPI_SUM, comm);
3125 part[0] = part[1] - locSize;
3127 MPI_Allreduce(&locSize, &glbSize, 1, HYPRE_MPI_INT, MPI_SUM, comm);
3131 part =
new HYPRE_Int[numProcs+1];
3133 MPI_Allgather(&locSize, 1, MPI_INT,
3134 &part[1], 1, HYPRE_MPI_INT, comm);
3137 for (
int i=0; i<numProcs; i++)
3139 part[i+1] += part[i];
3142 glbSize = part[numProcs];
3153 matvec_fn.MatvecCreate = this->OperatorMatvecCreate;
3154 matvec_fn.Matvec = this->OperatorMatvec;
3155 matvec_fn.MatvecDestroy = this->OperatorMatvecDestroy;
3157 HYPRE_LOBPCGSetup(lobpcg_solver,(HYPRE_Matrix)&A,NULL,NULL);
3163 matvec_fn.MatvecCreate = this->OperatorMatvecCreate;
3164 matvec_fn.Matvec = this->OperatorMatvec;
3165 matvec_fn.MatvecDestroy = this->OperatorMatvecDestroy;
3167 HYPRE_LOBPCGSetupB(lobpcg_solver,(HYPRE_Matrix)&M,NULL);
3176 for (
int i=0; i<nev; i++)
3178 eigs[i] = eigenvalues[i];
3185 return multi_vec->GetVector(i);
3192 if ( multi_vec == NULL )
3194 MFEM_ASSERT(x != NULL,
"In HypreLOBPCG::Solve()");
3196 multi_vec =
new HypreMultiVector(nev, *x, interpreter);
3197 multi_vec->Randomize(seed);
3207 HYPRE_LOBPCGSolve(lobpcg_solver, NULL, *multi_vec, eigenvalues);
3211 HypreLOBPCG::OperatorMatvecCreate(
void *A,
3218 return ( matvec_data );
3222 HypreLOBPCG::OperatorMatvec(
void *matvec_data,
3223 HYPRE_Complex
alpha,
3229 MFEM_VERIFY(alpha == 1.0 && beta == 0.0,
"values not supported");
3231 Operator *Aop = (Operator*)A;
3233 int width = Aop->Width();
3235 hypre_ParVector * xPar = (hypre_ParVector *)x;
3236 hypre_ParVector * yPar = (hypre_ParVector *)y;
3238 Vector xVec(xPar->local_vector->data, width);
3239 Vector yVec(yPar->local_vector->data, width);
3241 Aop->Mult( xVec, yVec );
3247 HypreLOBPCG::OperatorMatvecDestroy(
void *matvec_data )
3253 HypreLOBPCG::PrecondSolve(
void *solver,
3258 Solver *PC = (Solver*)solver;
3259 Operator *OP = (Operator*)A;
3261 int width = OP->Width();
3263 hypre_ParVector * bPar = (hypre_ParVector *)b;
3264 hypre_ParVector * xPar = (hypre_ParVector *)x;
3266 Vector bVec(bPar->local_vector->data, width);
3267 Vector xVec(xPar->local_vector->data, width);
3269 PC->Mult( bVec, xVec );
3275 HypreLOBPCG::PrecondSetup(
void *solver,
3293 MPI_Comm_size(comm,&numProcs);
3294 MPI_Comm_rank(comm,&myid);
3296 HYPRE_AMECreate(&ame_solver);
3297 HYPRE_AMESetPrintLevel(ame_solver, 0);
3304 hypre_TFree(multi_vec);
3309 for (
int i=0; i<nev; i++)
3311 delete eigenvectors[i];
3314 delete [] eigenvectors;
3318 hypre_TFree(eigenvalues);
3321 HYPRE_AMEDestroy(ame_solver);
3329 HYPRE_AMESetBlockSize(ame_solver, nev);
3335 HYPRE_AMESetTol(ame_solver, tol);
3341 HYPRE_AMESetMaxIter(ame_solver, max_iter);
3349 HYPRE_AMESetPrintLevel(ame_solver, logging);
3356 ams_precond = &precond;
3364 HYPRE_Solver ams_precond_ptr = (HYPRE_Solver)*ams_precond;
3366 ams_precond->
SetupFcn()(*ams_precond,A,NULL,NULL);
3368 HYPRE_AMESetAMSSolver(ame_solver, ams_precond_ptr);
3371 HYPRE_AMESetup(ame_solver);
3377 HYPRE_ParCSRMatrix parcsr_M = M;
3378 HYPRE_AMESetMassMatrix(ame_solver,(HYPRE_ParCSRMatrix)parcsr_M);
3384 HYPRE_AMESolve(ame_solver);
3391 eigs.
SetSize(nev); eigs = -1.0;
3393 if ( eigenvalues == NULL )
3396 HYPRE_AMEGetEigenvalues(ame_solver,&eigenvalues);
3400 for (
int i=0; i<nev; i++)
3402 eigs[i] = eigenvalues[i];
3407 HypreAME::createDummyVectors()
3409 if ( multi_vec == NULL )
3411 HYPRE_AMEGetEigenvectors(ame_solver,&multi_vec);
3414 eigenvectors =
new HypreParVector*[nev];
3415 for (
int i=0; i<nev; i++)
3417 eigenvectors[i] =
new HypreParVector(multi_vec[i]);
3426 if ( eigenvectors == NULL )
3428 this->createDummyVectors();
3431 return *eigenvectors[i];
3437 if ( eigenvectors == NULL )
3439 this->createDummyVectors();
3444 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.
int setup_called
Was hypre's Setup function called already?
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to 'master'.
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_Int *row_starts=NULL) const
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().
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.
HYPRE_Int GlobalTrueVSize()
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)
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.
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.
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)
HYPRE_Int * GetTrueDofOffsets()
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 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.
HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P)
Returns the matrix P^t * A * P.
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
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.
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 mfem_error(const char *msg)
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.
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
HypreParMatrix * Transpose()
Returns the transpose of *this.
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)
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 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.
Arbitrary order H(curl)-conforming Nedelec finite elements.
const FiniteElementCollection * FEColl() const
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.
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.
HypreParMatrix * ParMult(HypreParMatrix *A, HypreParMatrix *B)
Returns the matrix A * B.
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.