15 #include "../general/forall.hpp"
16 #include "../general/table.hpp"
17 #include "../general/sort_pairs.hpp"
18 #include "../general/backends.hpp"
27 #if defined(MFEM_USE_CUDA)
28 #define MFEM_cu_or_hip(stub) cu##stub
29 #define MFEM_Cu_or_Hip(stub) Cu##stub
30 #define MFEM_CU_or_HIP(stub) CU##stub
31 #define MFEM_CUDA_or_HIP(stub) CUDA##stub
33 #if CUSPARSE_VERSION >= 11400
34 #define MFEM_GPUSPARSE_ALG CUSPARSE_SPMV_CSR_ALG1
35 #else // CUSPARSE_VERSION >= 11400
36 #define MFEM_GPUSPARSE_ALG CUSPARSE_CSRMV_ALG1
37 #endif // CUSPARSE_VERSION >= 11400
39 #elif defined(MFEM_USE_HIP)
40 #define MFEM_cu_or_hip(stub) hip##stub
41 #define MFEM_Cu_or_Hip(stub) Hip##stub
42 #define MFEM_CU_or_HIP(stub) HIP##stub
43 #define MFEM_CUDA_or_HIP(stub) HIP##stub
46 #define MFEM_GPUSPARSE_ALG HIPSPARSE_CSRMV_ALG1
47 #endif // defined(MFEM_USE_CUDA)
54 #ifdef MFEM_USE_CUDA_OR_HIP
58 MFEM_cu_or_hip(sparseHandle_t) SparseMatrix::handle =
nullptr;
60 size_t SparseMatrix::bufferSize = 0;
61 void * SparseMatrix::dBuffer =
nullptr;
62 #endif // MFEM_USE_CUDA_OR_HIP
67 #ifdef MFEM_USE_CUDA_OR_HIP
70 if (!handle) { MFEM_cu_or_hip(sparseCreate)(&handle); }
78 #endif // MFEM_USE_CUDA_OR_HIP
83 #ifdef MFEM_USE_CUDA_OR_HIP
86 #if CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
87 MFEM_cu_or_hip(sparseDestroySpMat)(matA_descr);
88 MFEM_cu_or_hip(sparseDestroyDnVec)(vecX_descr);
89 MFEM_cu_or_hip(sparseDestroyDnVec)(vecY_descr);
91 cusparseDestroyMatDescr(matA_descr);
92 #endif // CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
95 #endif // MFEM_USE_CUDA_OR_HIP
100 Rows(new RowNode *[nrows]),
112 for (
int i = 0; i < nrows; i++)
117 #ifdef MFEM_USE_MEMALLOC
134 A.
Wrap(data,
I[height],
true);
136 #ifdef MFEM_USE_MEMALLOC
144 bool ownij,
bool owna,
bool issorted)
155 #ifdef MFEM_USE_MEMALLOC
161 A.
Wrap(data,
I[height], owna);
167 for (
int ii=0; ii<nnz; ++ii)
184 #ifdef MFEM_USE_MEMALLOC
188 J.
New(nrows * rowsize);
189 A.
New(nrows * rowsize);
191 for (
int i = 0; i <= nrows; i++)
225 #ifdef MFEM_USE_MEMALLOC
231 #ifdef MFEM_USE_MEMALLOC
235 for (
int i = 0; i <
height; i++)
237 RowNode **node_pp = &
Rows[i];
238 for (RowNode *node_p = mat.
Rows[i]; node_p; node_p = node_p->Prev)
240 #ifdef MFEM_USE_MEMALLOC
243 RowNode *new_node_p =
new RowNode;
245 new_node_p->Value = node_p->Value;
246 new_node_p->Column = node_p->Column;
247 *node_pp = new_node_p;
248 node_pp = &new_node_p->Prev;
276 #ifdef MFEM_USE_MEMALLOC
283 for (
int i = 0; i <=
height; i++)
288 for (
int r=0; r<
height; r++)
315 MFEM_ASSERT(master.
Finalized(),
"'master' must be finalized");
336 #ifdef MFEM_USE_MEMALLOC
354 return I[gi+1]-
I[gi];
358 RowNode *row =
Rows[gi];
359 for ( ; row != NULL; row = row->Prev)
360 if (row->Value != 0.0)
373 for (
int i=0; i <
height; ++i)
375 rowSize =
I[i+1]-
I[i];
376 max_row_size = (max_row_size > rowSize) ? max_row_size : rowSize;
381 for (
int i=0; i <
height; ++i)
384 max_row_size = (max_row_size > rowSize) ? max_row_size : rowSize;
393 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
400 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
407 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
414 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
421 if (newWidth ==
width)
426 else if (newWidth == -1)
432 else if (newWidth >
width)
443 ColPtrJ =
static_cast<int *
>(NULL);
451 "The new width needs to be bigger or equal to the actual width");
459 MFEM_VERIFY(
Finalized(),
"Matrix is not Finalized!");
466 #ifdef MFEM_USE_CUDA_OR_HIP
469 #if defined(MFEM_USE_CUDA)
470 size_t pBufferSizeInBytes = 0;
471 void *pBuffer = NULL;
474 const int m =
Width();
477 const int * d_ia =
ReadI();
479 csru2csrInfo_t sortInfoA;
482 cusparseCreateMatDescr( &matA_descr );
483 cusparseSetMatIndexBase( matA_descr, CUSPARSE_INDEX_BASE_ZERO );
484 cusparseSetMatType( matA_descr, CUSPARSE_MATRIX_TYPE_GENERAL );
486 cusparseCreateCsru2csrInfo( &sortInfoA );
488 cusparseDcsru2csr_bufferSizeExt(
handle, n, m, nnzA, d_a_sorted, d_ia,
489 d_ja_sorted, sortInfoA,
490 &pBufferSizeInBytes);
494 cusparseDcsru2csr(
handle, n, m, nnzA, matA_descr, d_a_sorted, d_ia,
495 d_ja_sorted, sortInfoA, pBuffer);
497 cusparseDestroyCsru2csrInfo( sortInfoA );
498 cusparseDestroyMatDescr( matA_descr );
505 #if defined(MFEM_USE_HIP)
506 size_t pBufferSizeInBytes = 0;
507 void *pBuffer = NULL;
511 const int m =
Width();
514 const int * d_ia =
ReadI();
517 hipsparseMatDescr_t descrA;
518 hipsparseCreateMatDescr( &descrA );
524 double *d_a_tmp = a_tmp.
Write();
526 hipsparseXcsrsort_bufferSizeExt(
handle, n, m, nnzA, d_ia, d_ja_sorted,
527 &pBufferSizeInBytes);
532 hipsparseCreateIdentityPermutation(
handle, nnzA, P);
533 hipsparseXcsrsort(
handle, n, m, nnzA, descrA, d_ia, d_ja_sorted, P, pBuffer);
535 hipsparseDgthr(
handle, nnzA, d_a_sorted, d_a_tmp, P,
536 HIPSPARSE_INDEX_BASE_ZERO);
539 hipsparseDestroyMatDescr( descrA );
546 #endif // MFEM_USE_CUDA_OR_HIP
553 for (
int j = 0, i = 0; i <
height; i++)
557 for (
int k = 0; k < row.
Size(); k++)
563 for (
int k = 0; k < row.
Size(); k++, j++)
575 MFEM_VERIFY(
Finalized(),
"Matrix is not Finalized!");
577 for (
int row = 0, end = 0; row <
height; row++)
581 for (j = start;
true; j++)
583 MFEM_VERIFY(j < end,
"diagonal entry not found in row = " << row);
584 if (
J[j] == row) {
break; }
586 const double diag =
A[j];
587 for ( ; j > start; j--)
609 MFEM_ASSERT(i < height && i >= 0 && j < width && j >= 0,
610 "Trying to access element outside of the matrix. "
611 <<
"height = " <<
height <<
", "
612 <<
"width = " <<
width <<
", "
613 <<
"i = " << i <<
", "
616 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
618 for (
int k =
I[i], end =
I[i+1]; k < end; k++)
626 MFEM_ABORT(
"Did not find i = " << i <<
", j = " << j <<
" in matrix.");
632 static const double zero = 0.0;
634 MFEM_ASSERT(i < height && i >= 0 && j < width && j >= 0,
635 "Trying to access element outside of the matrix. "
636 <<
"height = " <<
height <<
", "
637 <<
"width = " <<
width <<
", "
638 <<
"i = " << i <<
", "
643 for (
int k =
I[i], end =
I[i+1]; k < end; k++)
653 for (RowNode *node_p =
Rows[i]; node_p != NULL; node_p = node_p->Prev)
655 if (node_p->Column == j)
657 return node_p->Value;
667 MFEM_VERIFY(
height ==
width,
"Matrix must be square, not height = "
669 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
673 const auto II = this->
ReadI();
674 const auto JJ = this->
ReadJ();
680 const int begin = II[i];
681 const int end = II[i+1];
683 for (j = begin; j < end; j++)
701 int num_rows = this->
Height();
702 int num_cols = this->
Width();
717 for (
int r=0; r<
height; r++)
722 for (
int cj=0; cj<this->
RowSize(r); cj++)
724 B(r, col[cj]) = val[cj];
738 MFEM_ASSERT(
width == x.
Size(),
"Input vector size (" << x.
Size()
739 <<
") must match matrix width (" <<
width <<
")");
740 MFEM_ASSERT(
height == y.
Size(),
"Output vector size (" << y.
Size()
741 <<
") must match matrix height (" <<
height <<
")");
749 for (
int i = 0; i <
height; i++)
751 RowNode *row =
Rows[i];
753 for ( ; row != NULL; row = row->Prev)
755 b += row->Value * xp[row->Column];
763 #ifndef MFEM_USE_LEGACY_OPENMP
766 auto d_I =
Read(
I, height+1);
767 auto d_J =
Read(
J, nnz);
768 auto d_A =
Read(
A, nnz);
773 if (nnz == 0) {
return;}
776 #ifdef MFEM_USE_CUDA_OR_HIP
778 const double beta = 1.0;
783 #if CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
785 MFEM_cu_or_hip(sparseCreateCsr)(
789 const_cast<int *
>(d_I),
790 const_cast<int *>(d_J),
791 const_cast<double *
>(d_A),
792 MFEM_CU_or_HIP(SPARSE_INDEX_32I),
793 MFEM_CU_or_HIP(SPARSE_INDEX_32I),
794 MFEM_CU_or_HIP(SPARSE_INDEX_BASE_ZERO),
795 MFEM_CUDA_or_HIP(_R_64F));
798 MFEM_cu_or_hip(sparseCreateDnVec)(&
vecX_descr,
800 const_cast<double *
>(d_x),
801 MFEM_CUDA_or_HIP(_R_64F));
803 MFEM_CUDA_or_HIP(_R_64F));
806 cusparseSetMatIndexBase(
matA_descr, CUSPARSE_INDEX_BASE_ZERO);
807 cusparseSetMatType(
matA_descr, CUSPARSE_MATRIX_TYPE_GENERAL);
808 #endif // CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
812 size_t newBufferSize = 0;
814 MFEM_cu_or_hip(sparseSpMV_bufferSize)(
816 MFEM_CU_or_HIP(SPARSE_OPERATION_NON_TRANSPOSE),
822 MFEM_CUDA_or_HIP(_R_64F),
827 if (newBufferSize > bufferSize)
829 bufferSize = newBufferSize;
834 #if CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
836 MFEM_cu_or_hip(sparseDnVecSetValues)(
vecX_descr,
837 const_cast<double *
>(d_x));
838 MFEM_cu_or_hip(sparseDnVecSetValues)(
vecY_descr, d_y);
841 MFEM_cu_or_hip(sparseSpMV)(
843 MFEM_CU_or_HIP(SPARSE_OPERATION_NON_TRANSPOSE),
849 MFEM_CUDA_or_HIP(_R_64F),
854 CUSPARSE_OPERATION_NON_TRANSPOSE,
860 const_cast<double *
>(d_A),
861 const_cast<int *>(d_I),
862 const_cast<int *
>(d_J),
863 const_cast<double *>(d_x),
866 #endif // CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
867 #endif // MFEM_USE_CUDA_OR_HIP
872 MFEM_FORALL(i, height,
875 const int end = d_I[i+1];
876 for (
int j = d_I[i]; j < end; j++)
878 d += d_A[j] * d_x[d_J[j]];
885 #else // MFEM_USE_LEGACY_OPENMP
886 const double *Ap =
A, *xp = x.
GetData();
888 const int *Jp =
J, *Ip =
I;
890 #pragma omp parallel for
891 for (
int i = 0; i <
height; i++)
894 const int end = Ip[i+1];
895 for (
int j = Ip[i]; j < end; j++)
897 d += Ap[j] * xp[Jp[j]];
901 #endif // MFEM_USE_LEGACY_OPENMP
912 const double a)
const
915 <<
") must match matrix height (" <<
height <<
")");
916 MFEM_ASSERT(
width == y.
Size(),
"Output vector size (" << y.
Size()
917 <<
") must match matrix width (" <<
width <<
")");
924 for (
int i = 0; i <
height; i++)
926 RowNode *row =
Rows[i];
927 double b = a * xp[i];
928 for ( ; row != NULL; row = row->Prev)
930 yp[row->Column] += row->Value *
b;
943 for (
int i = 0; i <
height; i++)
945 const double xi = a * x[i];
946 const int end =
I[i+1];
947 for (
int j =
I[i]; j < end; j++)
981 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
983 const int n = rows.
Size();
985 auto d_rows = rows.
Read();
987 auto d_J =
Read(
J, nnz);
988 auto d_A =
Read(
A, nnz);
990 auto d_y = y.
Write();
993 const int r = d_rows[i];
994 const int end = d_I[r + 1];
996 for (
int j = d_I[r]; j < end; j++)
998 a += d_A[j] * d_x[d_J[j]];
1007 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
1009 for (
int i = 0; i < rows.
Size(); i++)
1014 for (
int j =
I[r]; j < end; j++)
1016 val +=
A[j] * x(
J[j]);
1024 MFEM_ASSERT(
Finalized(),
"Matrix must be finalized.");
1025 MFEM_ASSERT(x.
Size() ==
Width(),
"Input vector size (" << x.
Size()
1026 <<
") must match matrix width (" <<
Width() <<
")");
1032 auto d_I =
Read(
I, height+1);
1033 auto d_J =
Read(
J, nnz);
1036 MFEM_FORALL(i, height,
1039 const int end = d_I[i+1];
1040 for (
int j = d_I[i]; j < end; j++)
1055 MFEM_ASSERT(
Finalized(),
"Matrix must be finalized.");
1056 MFEM_ASSERT(x.
Size() ==
Height(),
"Input vector size (" << x.
Size()
1057 <<
") must match matrix height (" <<
Height() <<
")");
1062 for (
int i = 0; i <
Height(); i++)
1067 for (
int j =
I[i]; j < end; j++)
1077 MFEM_ASSERT(
width == x.
Size(),
"Input vector size (" << x.
Size()
1078 <<
") must match matrix width (" <<
width <<
")");
1079 MFEM_ASSERT(
height == y.
Size(),
"Output vector size (" << y.
Size()
1080 <<
") must match matrix height (" <<
height <<
")");
1091 for (
int i = 0; i <
height; i++)
1093 RowNode *row =
Rows[i];
1095 for ( ; row != NULL; row = row->Prev)
1097 b += std::abs(row->Value) * xp[row->Column];
1107 auto d_I =
Read(
I, height+1);
1108 auto d_J =
Read(
J, nnz);
1109 auto d_A =
Read(
A, nnz);
1110 auto d_x = x.
Read();
1112 MFEM_FORALL(i, height,
1115 const int end = d_I[i+1];
1116 for (
int j = d_I[i]; j < end; j++)
1118 d += std::abs(d_A[j]) * d_x[d_J[j]];
1126 MFEM_ASSERT(
height == x.
Size(),
"Input vector size (" << x.
Size()
1127 <<
") must match matrix height (" <<
height <<
")");
1128 MFEM_ASSERT(
width == y.
Size(),
"Output vector size (" << y.
Size()
1129 <<
") must match matrix width (" <<
width <<
")");
1137 for (
int i = 0; i <
height; i++)
1139 RowNode *row =
Rows[i];
1141 for ( ; row != NULL; row = row->Prev)
1143 yp[row->Column] += fabs(row->Value) *
b;
1156 for (
int i = 0; i <
height; i++)
1158 const double xi = x[i];
1159 const int end =
I[i+1];
1160 for (
int j =
I[i]; j < end; j++)
1162 const int Jj =
J[j];
1163 y[Jj] += std::abs(
A[j]) * xi;
1172 <<
" must be equal to Width() = " <<
Width());
1174 <<
" must be equal to Height() = " <<
Height());
1187 for (
int i = 0; i <
height; i++)
1192 for (
int j =
I[i], end =
I[i+1]; j < end; j++)
1194 a +=
A[j] * x(
J[j]);
1199 for (RowNode *np =
Rows[i]; np != NULL; np = np->Prev)
1201 a += np->Value * x(np->Column);
1212 for (
int i = 0; i <
height; i++)
1217 for (
int j =
I[i], end =
I[i+1]; j < end; j++)
1224 for (RowNode *np =
Rows[i]; np != NULL; np = np->Prev)
1235 MFEM_VERIFY(irow <
height,
1236 "row " << irow <<
" not in matrix with height " <<
height);
1241 for (
int j =
I[irow], end =
I[irow+1]; j < end; j++)
1248 for (RowNode *np =
Rows[irow]; np != NULL; np = np->Prev)
1250 a += fabs(np->Value);
1259 MFEM_ASSERT(
Finalized(),
"Matrix must be finalized.");
1261 atol = std::abs(tol);
1263 fix_empty_rows =
height ==
width ? fix_empty_rows :
false;
1271 for (i = 0, nz = 0; i <
height; i++)
1274 for (j =
I[i]; j <
I[i+1]; j++)
1275 if (std::abs(
A[j]) > atol)
1280 if (fix_empty_rows && !found) { nz++; }
1288 for (i = 0, nz = 0; i <
height; i++)
1292 for (j =
I[i]; j <
I[i+1]; j++)
1293 if (std::abs(
A[j]) > atol)
1298 if ( lastCol > newJ[nz] )
1305 if (fix_empty_rows && !found)
1313 I.
Wrap(newI, height+1,
true);
1314 J.
Wrap(newJ,
I[height],
true);
1315 A.
Wrap(newA,
I[height],
true);
1333 for (i = 1; i <=
height; i++)
1336 for (aux =
Rows[i-1]; aux != NULL; aux = aux->Prev)
1338 if (skip_zeros && aux->Value == 0.0)
1340 if (skip_zeros == 2) {
continue; }
1341 if ((i-1) != aux->Column) {
continue; }
1345 for (RowNode *other =
Rows[aux->Column]; other != NULL; other = other->Prev)
1347 if (other->Column == (i-1))
1350 found_val = other->Value;
1354 if (found && found_val == 0.0) {
continue; }
1359 if (fix_empty_rows && !nr) { nr = 1; }
1368 for (j = i = 0; i <
height; i++)
1372 for (aux =
Rows[i]; aux != NULL; aux = aux->Prev)
1374 if (skip_zeros && aux->Value == 0.0)
1376 if (skip_zeros == 2) {
continue; }
1377 if (i != aux->Column) {
continue; }
1381 for (RowNode *other =
Rows[aux->Column]; other != NULL; other = other->Prev)
1383 if (other->Column == i)
1386 found_val = other->Value;
1390 if (found && found_val == 0.0) {
continue; }
1396 if ( lastCol >
J[j] )
1405 if (fix_empty_rows && !nr)
1413 #ifdef MFEM_USE_MEMALLOC
1417 for (i = 0; i <
height; i++)
1419 RowNode *node_p =
Rows[i];
1420 while (node_p != NULL)
1423 node_p = node_p->Prev;
1436 int nr = (
height + br - 1)/br, nc = (
width + bc - 1)/bc;
1438 for (
int j = 0; j < bc; j++)
1440 for (
int i = 0; i < br; i++)
1443 for (
int k = 0; k <= nr; k++)
1447 blocks(i,j) =
new SparseMatrix(bI, NULL, NULL, nr, nc);
1451 for (
int gr = 0; gr <
height; gr++)
1453 int bi = gr/nr, i = gr%nr + 1;
1456 for (
int j =
I[gr]; j <
I[gr+1]; j++)
1460 blocks(bi,
J[j]/nc)->I[i]++;
1466 for (RowNode *n_p =
Rows[gr]; n_p != NULL; n_p = n_p->Prev)
1468 if (n_p->Value != 0.0)
1470 blocks(bi, n_p->Column/nc)->I[i]++;
1476 for (
int j = 0; j < bc; j++)
1478 for (
int i = 0; i < br; i++)
1482 for (
int k = 1; k <= nr; k++)
1484 rs = b.
I[k], b.
I[k] = nnz, nnz += rs;
1491 for (
int gr = 0; gr <
height; gr++)
1493 int bi = gr/nr, i = gr%nr + 1;
1496 for (
int j =
I[gr]; j <
I[gr+1]; j++)
1501 b.
J[b.
I[i]] =
J[j] % nc;
1509 for (RowNode *n_p =
Rows[gr]; n_p != NULL; n_p = n_p->Prev)
1511 if (n_p->Value != 0.0)
1514 b.
J[b.
I[i]] = n_p->Column % nc;
1515 b.
A[b.
I[i]] = n_p->Value;
1537 for (
int i = 1; i <
height; i++)
1539 for (
int j =
I[i]; j <
I[i+1]; j++)
1543 symm = std::max(symm, std::abs(
A[j]-(*
this)(
J[j],i)));
1550 for (
int i = 0; i <
height; i++)
1552 for (RowNode *node_p =
Rows[i]; node_p != NULL; node_p = node_p->Prev)
1554 int col = node_p->Column;
1557 symm = std::max(symm, std::abs(node_p->Value-(*
this)(col,i)));
1567 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
1570 for (i = 1; i <
height; i++)
1572 for (j =
I[i]; j <
I[i+1]; j++)
1576 A[j] += (*this)(
J[j],i);
1578 (*this)(J[j],i) = A[j];
1595 for (
int i = 0; i <
height; i++)
1597 for (RowNode *node_p =
Rows[i]; node_p != NULL; node_p = node_p->Prev)
1614 for (
int j = 0; j < nnz; j++)
1616 m = std::max(m, std::abs(
A[j]));
1621 for (
int i = 0; i <
height; i++)
1623 for (RowNode *n_p =
Rows[i]; n_p != NULL; n_p = n_p->Prev)
1625 m = std::max(m, std::abs(n_p->Value));
1639 const double *Ap =
A;
1641 for (
int i = 0; i < nz; i++)
1643 counter += (std::abs(Ap[i]) <= tol);
1648 for (
int i = 0; i <
height; i++)
1650 for (RowNode *aux =
Rows[i]; aux != NULL; aux = aux->Prev)
1652 counter += (std::abs(aux->Value) <= tol);
1673 for (
int i = 0; i <
height; i++)
1675 for (RowNode *aux =
Rows[i]; aux != NULL; aux = aux->Prev)
1693 MFEM_ASSERT(row < height && row >= 0,
1694 "Row " << row <<
" not in matrix of height " <<
height);
1696 MFEM_VERIFY(!
Finalized(),
"Matrix must NOT be finalized.");
1698 for (aux =
Rows[row]; aux != NULL; aux = aux->Prev)
1700 rhs(aux->Column) -= sol * aux->Value;
1709 MFEM_ASSERT(row < height && row >= 0,
1710 "Row " << row <<
" not in matrix of height " <<
height);
1711 MFEM_ASSERT(dpolicy !=
DIAG_KEEP,
"Diagonal policy must not be DIAG_KEEP");
1713 "if dpolicy == DIAG_ONE, matrix must be square, not height = "
1718 for (
int i=
I[row]; i <
I[row+1]; ++i)
1725 for (aux =
Rows[row]; aux != NULL; aux = aux->Prev)
1739 MFEM_ASSERT(col < width && col >= 0,
1740 "Col " << col <<
" not in matrix of width " <<
width);
1741 MFEM_ASSERT(dpolicy !=
DIAG_KEEP,
"Diagonal policy must not be DIAG_KEEP");
1743 "if dpolicy == DIAG_ONE, matrix must be square, not height = "
1749 for (
int jpos = 0; jpos != nnz; ++jpos)
1759 for (
int i = 0; i <
height; i++)
1761 for (RowNode *aux =
Rows[i]; aux != NULL; aux = aux->Prev)
1763 if (aux->Column == col)
1782 for (
int i = 0; i <
height; i++)
1784 for (
int jpos =
I[i]; jpos !=
I[i+1]; ++jpos)
1790 (*b)(i) -=
A[jpos] * (*x)(
J[jpos] );
1799 for (
int i = 0; i <
height; i++)
1801 for (RowNode *aux =
Rows[i]; aux != NULL; aux = aux->Prev)
1803 if (cols[aux -> Column])
1807 (*b)(i) -= aux -> Value * (*x)(aux -> Column);
1821 for (
int row = 0; row <
height; row++)
1823 for (nd =
Rows[row]; nd != NULL; nd = nd->Prev)
1825 if (col_marker[nd->Column])
1827 Ae.
Add(row, nd->Column, nd->Value);
1835 for (
int row = 0; row <
height; row++)
1837 for (
int j =
I[row]; j <
I[row+1]; j++)
1839 if (col_marker[
J[j]])
1841 Ae.
Add(row,
J[j],
A[j]);
1853 MFEM_ASSERT(rc < height && rc >= 0,
1854 "Row " << rc <<
" not in matrix of height " <<
height);
1861 for (
int j =
I[rc]; j <
I[rc+1]; j++)
1863 const int col =
J[j];
1869 rhs(rc) =
A[j] * sol;
1880 mfem_error(
"SparseMatrix::EliminateRowCol () #2");
1887 for (
int k = I[col]; 1; k++)
1891 mfem_error(
"SparseMatrix::EliminateRowCol () #3");
1893 else if (
J[k] == rc)
1895 rhs(col) -= sol *
A[k];
1905 for (RowNode *aux =
Rows[rc]; aux != NULL; aux = aux->Prev)
1907 const int col = aux->Column;
1913 rhs(rc) = aux->Value * sol;
1924 mfem_error(
"SparseMatrix::EliminateRowCol () #4");
1931 for (RowNode *node =
Rows[col]; 1; node = node->Prev)
1935 mfem_error(
"SparseMatrix::EliminateRowCol () #5");
1937 else if (node->Column == rc)
1939 rhs(col) -= sol * node->Value;
1953 MFEM_ASSERT(rc < height && rc >= 0,
1954 "Row " << rc <<
" not in matrix of height " <<
height);
1955 MFEM_ASSERT(sol.
Size() == rhs.
Width(),
"solution size (" << sol.
Size()
1956 <<
") must match rhs width (" << rhs.
Width() <<
")");
1958 const int num_rhs = rhs.
Width();
1961 for (
int j =
I[rc]; j <
I[rc+1]; j++)
1963 const int col =
J[j];
1969 for (
int r = 0; r < num_rhs; r++)
1971 rhs(rc,r) =
A[j] * sol(r);
1976 for (
int r = 0; r < num_rhs; r++)
1983 for (
int r = 0; r < num_rhs; r++)
1989 mfem_error(
"SparseMatrix::EliminateRowColMultipleRHS() #3");
1996 for (
int k = I[col]; 1; k++)
2000 mfem_error(
"SparseMatrix::EliminateRowColMultipleRHS() #4");
2002 else if (
J[k] == rc)
2004 for (
int r = 0; r < num_rhs; r++)
2006 rhs(col,r) -= sol(r) *
A[k];
2017 for (RowNode *aux =
Rows[rc]; aux != NULL; aux = aux->Prev)
2019 const int col = aux->Column;
2025 for (
int r = 0; r < num_rhs; r++)
2027 rhs(rc,r) = aux->Value * sol(r);
2032 for (
int r = 0; r < num_rhs; r++)
2039 for (
int r = 0; r < num_rhs; r++)
2045 mfem_error(
"SparseMatrix::EliminateRowColMultipleRHS() #5");
2052 for (RowNode *node =
Rows[col]; 1; node = node->Prev)
2056 mfem_error(
"SparseMatrix::EliminateRowColMultipleRHS() #6");
2058 else if (node->Column == rc)
2060 for (
int r = 0; r < num_rhs; r++)
2062 rhs(col,r) -= sol(r) * node->Value;
2075 MFEM_ASSERT(rc < height && rc >= 0,
2076 "Row " << rc <<
" not in matrix of height " <<
height);
2080 const auto &II = this->
I;
2081 const auto &JJ = this->
J;
2082 for (
int j = II[rc]; j < II[rc+1]; j++)
2084 const int col = JJ[j];
2099 for (
int k = II[col]; 1; k++)
2103 mfem_error(
"SparseMatrix::EliminateRowCol() #2");
2105 else if (JJ[k] == rc)
2116 RowNode *aux, *node;
2118 for (aux =
Rows[rc]; aux != NULL; aux = aux->Prev)
2120 const int col = aux->Column;
2135 for (node =
Rows[col]; 1; node = node->Prev)
2139 mfem_error(
"SparseMatrix::EliminateRowCol() #3");
2141 else if (node->Column == rc)
2156 MFEM_ASSERT(rc < height && rc >= 0,
2157 "Row " << rc <<
" not in matrix of height " <<
height);
2161 for (
int j =
I[rc]; j <
I[rc+1]; j++)
2163 const int col =
J[j];
2171 for (
int k = I[col]; 1; k++)
2175 mfem_error(
"SparseMatrix::EliminateRowCol() #2");
2177 else if (
J[k] == rc)
2188 RowNode *aux, *node;
2190 for (aux =
Rows[rc]; aux != NULL; aux = aux->Prev)
2192 const int col = aux->Column;
2200 for (node =
Rows[col]; 1; node = node->Prev)
2204 mfem_error(
"SparseMatrix::EliminateRowCol() #3");
2206 else if (node->Column == rc)
2223 for (nd =
Rows[rc]; nd != NULL; nd = nd->Prev)
2225 const int col = nd->Column;
2231 Ae.
Add(rc, rc, nd->Value - 1.0);
2235 Ae.
Add(rc, rc, nd->Value);
2241 mfem_error(
"SparseMatrix::EliminateRowCol #1");
2247 Ae.
Add(rc, col, nd->Value);
2249 for (nd2 =
Rows[col]; 1; nd2 = nd2->Prev)
2253 mfem_error(
"SparseMatrix::EliminateRowCol #2");
2255 else if (nd2->Column == rc)
2257 Ae.
Add(col, rc, nd2->Value);
2267 for (
int j =
I[rc]; j <
I[rc+1]; j++)
2269 const int col =
J[j];
2275 Ae.
Add(rc, rc,
A[j] - 1.0);
2279 Ae.
Add(rc, rc,
A[j]);
2285 mfem_error(
"SparseMatrix::EliminateRowCol #3");
2291 Ae.
Add(rc, col,
A[j]);
2293 for (
int k = I[col];
true; k++)
2297 mfem_error(
"SparseMatrix::EliminateRowCol #4");
2299 else if (
J[k] == rc)
2301 Ae.
Add(col, rc,
A[k]);
2314 const int n_ess_dofs = ess_dofs.
Size();
2315 const auto ess_dofs_d = ess_dofs.
Read();
2316 const auto dI =
ReadI();
2317 const auto dJ =
ReadJ();
2320 MFEM_FORALL(i, n_ess_dofs,
2322 const int idof = ess_dofs_d[i];
2323 for (
int j=dI[idof]; j<dI[idof+1]; ++j)
2325 const int jdof = dJ[j];
2329 for (
int k=dI[jdof]; k<dI[jdof+1]; ++k)
2340 if (diag_policy == DiagonalPolicy::DIAG_ONE)
2344 else if (diag_policy == DiagonalPolicy::DIAG_ZERO)
2356 for (
int i = 0; i <
height; i++)
2358 if (
I[i+1] ==
I[i]+1 && fabs(
A[
I[i]]) < 1e-16)
2367 for (
int i = 0; i <
height; i++)
2370 for (
int j =
I[i]; j <
I[i+1]; j++)
2374 if (zero <= threshold)
2376 for (
int j = I[i]; j < I[i+1]; j++)
2378 A[j] = (
J[j] == i) ? 1.0 : 0.0;
2389 const double *xp = x.
GetData();
2390 RowNode *diag_p, *n_p, **R =
Rows;
2393 for (
int i = 0; i <
s; i++)
2397 for (n_p = R[i]; n_p != NULL; n_p = n_p->Prev)
2399 const int c = n_p->Column;
2406 sum += n_p->Value * yp[c];
2410 if (diag_p != NULL && diag_p->Value != 0.0)
2412 yp[i] = (xp[i] - sum) / diag_p->Value;
2414 else if (xp[i] == sum)
2420 mfem_error(
"SparseMatrix::Gauss_Seidel_forw()");
2434 for (
int i = 0, j = Ip[0]; i <
s; i++)
2436 const int end = Ip[i+1];
2439 for ( ; j < end; j++)
2441 const int c = Jp[j];
2448 sum += Ap[j] * yp[c];
2452 if (d >= 0 && Ap[d] != 0.0)
2454 yp[i] = (xp[i] - sum) / Ap[d];
2456 else if (xp[i] == sum)
2462 mfem_error(
"SparseMatrix::Gauss_Seidel_forw(...) #2");
2473 const double *xp = x.
GetData();
2474 RowNode *diag_p, *n_p, **R =
Rows;
2476 for (
int i =
height-1; i >= 0; i--)
2480 for (n_p = R[i]; n_p != NULL; n_p = n_p->Prev)
2482 const int c = n_p->Column;
2489 sum += n_p->Value * yp[c];
2493 if (diag_p != NULL && diag_p->Value != 0.0)
2495 yp[i] = (xp[i] - sum) / diag_p->Value;
2497 else if (xp[i] == sum)
2503 mfem_error(
"SparseMatrix::Gauss_Seidel_back()");
2517 for (
int i = s-1, j = Ip[s]-1; i >= 0; i--)
2519 const int beg = Ip[i];
2522 for ( ; j >= beg; j--)
2524 const int c = Jp[j];
2531 sum += Ap[j] * yp[c];
2535 if (d >= 0 && Ap[d] != 0.0)
2537 yp[i] = (xp[i] - sum) / Ap[d];
2539 else if (xp[i] == sum)
2545 mfem_error(
"SparseMatrix::Gauss_Seidel_back(...) #2");
2553 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
2556 for (
int i = 0; i <
height; i++)
2560 for (
int j =
I[i]; j <
I[i+1]; j++)
2568 if (d >= 0 &&
A[d] != 0.0)
2570 double a = 1.8 * fabs(
A[d]) / norm;
2578 mfem_error(
"SparseMatrix::GetJacobiScaling() #2");
2585 double sc,
bool use_abs_diag)
const
2587 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
2589 for (
int i = 0; i <
height; i++)
2593 for (
int j =
I[i]; j <
I[i+1]; j++)
2601 sum -=
A[j] * x0(
J[j]);
2604 if (d >= 0 &&
A[d] != 0.0)
2606 const double diag = (use_abs_diag) ? fabs(
A[d]) :
A[d];
2607 x1(i) = sc * (sum / diag) + (1.0 - sc) * x0(i);
2617 double sc,
bool use_abs_diag)
const
2619 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
2625 const auto Ap =
Read(
A, nnz, use_dev);
2627 const auto Jp =
Read(
J, nnz, use_dev);
2629 const auto bp = b.
Read(use_dev);
2630 auto xp = x.
Write(use_dev);
2632 MFEM_FORALL_SWITCH(use_dev, i, H,
2634 const int end = Ip[i+1];
2635 for (
int j = Ip[i];
true; j++)
2639 MFEM_ABORT_KERNEL(
"Diagonal not found in SparseMatrix::DiagScale");
2643 const double diag = (use_abs_diag) ? fabs(Ap[j]) : Ap[j];
2646 MFEM_ABORT_KERNEL(
"Zero diagonal in SparseMatrix::DiagScale");
2648 xp[i] = sc * bp[i] / diag;
2655 template <
bool useFabs>
2663 const auto bp = b.
Read(useDevice);
2664 const auto x0p = x0.
Read(useDevice);
2665 auto x1p = x1.
Write(useDevice);
2667 const auto Ip =
Read(I, height+1, useDevice);
2671 MFEM_FORALL_SWITCH(useDevice, i, height,
2673 double resi = bp[i], norm = 0.0;
2674 for (
int j = Ip[i]; j < Ip[i+1]; j++)
2676 resi -= Ap[j] * x0p[Jp[j]];
2679 norm += fabs(Ap[j]);
2688 x1p[i] = x0p[i] + sc * resi / norm;
2694 MFEM_ABORT_KERNEL(
"L1 norm of row is zero.");
2698 MFEM_ABORT_KERNEL(
"sum of row is zero.");
2707 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
2708 JacobiDispatch<true>(
b,x0,x1,
I,
J,
A,
height,sc);
2714 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
2715 JacobiDispatch<false>(
b,x0,x1,
I,
J,
A,
height,sc);
2721 int i, j, gi, gj,
s,
t;
2731 for (i = 0; i < rows.
Size(); i++)
2733 if ((gi=rows[i]) < 0) { gi = -1-gi, s = -1; }
2735 MFEM_ASSERT(gi < height,
2736 "Trying to insert a row " << gi <<
" outside the matrix height "
2739 for (j = 0; j < cols.
Size(); j++)
2741 if ((gj=cols[j]) < 0) { gj = -1-gj, t = -
s; }
2743 MFEM_ASSERT(gj <
width,
2744 "Trying to insert a column " << gj <<
" outside the matrix width "
2747 if (skip_zeros && a == 0.0)
2752 if (skip_zeros == 2 || &rows != &cols || subm(j, i) == 0.0)
2757 if (t < 0) { a = -
a; }
2769 if ((gi=i) < 0) { gi = -1-gi, s = -1; }
2771 MFEM_ASSERT(gi < height,
2772 "Trying to set a row " << gi <<
" outside the matrix height "
2774 if ((gj=j) < 0) { gj = -1-gj, t = -
s; }
2776 MFEM_ASSERT(gj <
width,
2777 "Trying to set a column " << gj <<
" outside the matrix width "
2779 if (t < 0) { a = -
a; }
2788 if ((gi=i) < 0) { gi = -1-gi, s = -1; }
2790 MFEM_ASSERT(gi < height,
2791 "Trying to insert a row " << gi <<
" outside the matrix height "
2793 if ((gj=j) < 0) { gj = -1-gj, t = -
s; }
2795 MFEM_ASSERT(gj <
width,
2796 "Trying to insert a column " << gj <<
" outside the matrix width "
2798 if (t < 0) { a = -
a; }
2805 int i, j, gi, gj,
s,
t;
2808 for (i = 0; i < rows.
Size(); i++)
2810 if ((gi=rows[i]) < 0) { gi = -1-gi, s = -1; }
2812 MFEM_ASSERT(gi < height,
2813 "Trying to set a row " << gi <<
" outside the matrix height "
2816 for (j = 0; j < cols.
Size(); j++)
2819 if (skip_zeros && a == 0.0)
2824 if (skip_zeros == 2 || &rows != &cols || subm(j, i) == 0.0)
2829 if ((gj=cols[j]) < 0) { gj = -1-gj, t = -
s; }
2831 MFEM_ASSERT(gj <
width,
2832 "Trying to set a column " << gj <<
" outside the matrix width "
2834 if (t < 0) { a = -
a; }
2846 int i, j, gi, gj,
s,
t;
2849 for (i = 0; i < rows.
Size(); i++)
2851 if ((gi=rows[i]) < 0) { gi = -1-gi, s = -1; }
2853 MFEM_ASSERT(gi < height,
2854 "Trying to set a row " << gi <<
" outside the matrix height "
2857 for (j = 0; j < cols.
Size(); j++)
2860 if (skip_zeros && a == 0.0)
2865 if (skip_zeros == 2 || &rows != &cols || subm(j, i) == 0.0)
2870 if ((gj=cols[j]) < 0) { gj = -1-gj, t = -
s; }
2872 MFEM_ASSERT(gj <
width,
2873 "Trying to set a column " << gj <<
" outside the matrix width "
2875 if (t < 0) { a = -
a; }
2885 int i, j, gi, gj,
s,
t;
2888 for (i = 0; i < rows.
Size(); i++)
2890 if ((gi=rows[i]) < 0) { gi = -1-gi, s = -1; }
2892 MFEM_ASSERT(gi < height,
2893 "Trying to read a row " << gi <<
" outside the matrix height "
2896 for (j = 0; j < cols.
Size(); j++)
2898 if ((gj=cols[j]) < 0) { gj = -1-gj, t = -
s; }
2900 MFEM_ASSERT(gj <
width,
2901 "Trying to read a column " << gj <<
" outside the matrix width "
2904 subm(i, j) = (t < 0) ? (-a) : (
a);
2918 MFEM_ASSERT(gi < height,
2919 "Trying to query a row " << gi <<
" outside the matrix height "
2923 return (
Rows[gi] == NULL);
2927 return (I[gi] == I[gi+1]);
2936 if ((gi=row) < 0) { gi = -1-gi; }
2937 MFEM_ASSERT(gi < height,
2938 "Trying to read a row " << gi <<
" outside the matrix height "
2942 for (n =
Rows[gi], j = 0; n; n = n->Prev)
2948 for (n =
Rows[gi], j = 0; n; n = n->Prev, j++)
2950 cols[j] = n->Column;
2963 cols.
MakeRef(const_cast<int*>((
const int*)J) + j, I[gi+1]-j);
2965 const_cast<double*>((
const double*)A) + j, cols.
Size());
2966 MFEM_ASSERT(row >= 0,
"Row not valid: " << row <<
", height: " << height);
2977 if ((gi=row) < 0) { gi = -1-gi, s = -1; }
2979 MFEM_ASSERT(gi < height,
2980 "Trying to set a row " << gi <<
" outside the matrix height "
2986 for (
int j = 0; j < cols.
Size(); j++)
2988 if ((gj=cols[j]) < 0) { gj = -1-gj, t = -
s; }
2990 MFEM_ASSERT(gj <
width,
2991 "Trying to set a column " << gj <<
" outside the matrix"
2992 " width " <<
width);
2994 if (t < 0) { a = -
a; }
3002 MFEM_ASSERT(cols.
Size() == srow.
Size(),
"");
3004 for (
int i = I[gi], j = 0; j < cols.
Size(); j++, i++)
3006 if ((gj=cols[j]) < 0) { gj = -1-gj, t = -
s; }
3008 MFEM_ASSERT(gj <
width,
3009 "Trying to set a column " << gj <<
" outside the matrix"
3010 " width " <<
width);
3021 int j, gi, gj,
s,
t;
3024 MFEM_VERIFY(!
Finalized(),
"Matrix must NOT be finalized.");
3026 if ((gi=row) < 0) { gi = -1-gi, s = -1; }
3028 MFEM_ASSERT(gi < height,
3029 "Trying to insert a row " << gi <<
" outside the matrix height "
3032 for (j = 0; j < cols.
Size(); j++)
3034 if ((gj=cols[j]) < 0) { gj = -1-gj, t = -
s; }
3036 MFEM_ASSERT(gj <
width,
3037 "Trying to insert a column " << gj <<
" outside the matrix width "
3044 if (t < 0) { a = -
a; }
3062 for (aux =
Rows[i]; aux != NULL; aux = aux -> Prev)
3064 aux -> Value *= scale;
3069 int j, end = I[i+1];
3071 for (j = I[i]; j < end; j++)
3084 for (
int i=0; i <
height; ++i)
3087 for (aux =
Rows[i]; aux != NULL; aux = aux -> Prev)
3089 aux -> Value *= scale;
3097 for (
int i=0; i <
height; ++i)
3101 for (j = I[i]; j < end; j++)
3114 for (
int i=0; i <
height; ++i)
3116 for (aux =
Rows[i]; aux != NULL; aux = aux -> Prev)
3118 aux -> Value *= sr(aux->Column);
3126 for (
int i=0; i <
height; ++i)
3129 for (j = I[i]; j < end; j++)
3140 "Mismatch of this matrix size and rhs. This height = "
3141 << height <<
", width = " <<
width <<
", B.height = "
3144 for (
int i = 0; i <
height; i++)
3149 for (RowNode *aux = B.
Rows[i]; aux != NULL; aux = aux->Prev)
3151 _Add_(aux->Column, aux->Value);
3156 for (
int j = B.
I[i]; j < B.
I[i+1]; j++)
3169 for (
int i = 0; i <
height; i++)
3174 for (RowNode *np =
Rows[i]; np != NULL; np = np->Prev)
3176 np->Value += a * B.
_Get_(np->Column);
3181 for (
int j = I[i]; j < I[i+1]; j++)
3183 A[j] += a * B.
_Get_(J[j]);
3196 for (
int i = 0; i < nnz; i++)
3203 for (
int i = 0; i <
height; i++)
3205 for (RowNode *node_p =
Rows[i]; node_p != NULL;
3206 node_p = node_p -> Prev)
3208 node_p -> Value =
a;
3220 for (
int i = 0, nnz = I[height]; i < nnz; i++)
3227 for (
int i = 0; i <
height; i++)
3229 for (RowNode *node_p =
Rows[i]; node_p != NULL;
3230 node_p = node_p -> Prev)
3232 node_p -> Value *=
a;
3247 for (i = 0; i <
height; i++)
3249 os <<
"[row " << i <<
"]\n";
3250 for (nd =
Rows[i], j = 0; nd != NULL; nd = nd->Prev, j++)
3252 os <<
" (" << nd->Column <<
"," << nd->Value <<
")";
3253 if ( !((j+1) % width_) )
3270 for (i = 0; i <
height; i++)
3272 os <<
"[row " << i <<
"]\n";
3273 for (j = I[i]; j < I[i+1]; j++)
3275 os <<
" (" << J[j] <<
"," << A[j] <<
")";
3276 if ( !((j+1-I[i]) % width_) )
3281 if ((j-I[i]) % width_)
3290 os <<
"% size " << height <<
" " <<
width <<
"\n";
3293 ios::fmtflags old_fmt = os.flags();
3294 os.setf(ios::scientific);
3295 std::streamsize old_prec = os.precision(14);
3297 for (i = 0; i <
height; i++)
3299 for (j = I[i]; j < I[i+1]; j++)
3301 os << i+1 <<
" " << J[j]+1 <<
" " << A[j] <<
'\n';
3305 os << height <<
" " << width <<
" 0.0\n";
3306 os.precision(old_prec);
3313 ios::fmtflags old_fmt = os.flags();
3314 os.setf(ios::scientific);
3315 std::streamsize old_prec = os.precision(14);
3317 os <<
"%%MatrixMarket matrix coordinate real general" <<
'\n'
3318 <<
"% Generated by MFEM" <<
'\n';
3321 for (i = 0; i <
height; i++)
3323 for (j = I[i]; j < I[i+1]; j++)
3325 os << i+1 <<
" " << J[j]+1 <<
" " << A[j] <<
'\n';
3328 os.precision(old_prec);
3334 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
3338 os << height <<
'\n';
3340 for (i = 0; i <=
height; i++)
3342 os << I[i]+1 <<
'\n';
3345 for (i = 0; i < I[
height]; i++)
3347 os << J[i]+1 <<
'\n';
3350 for (i = 0; i < I[
height]; i++)
3358 MFEM_VERIFY(
Finalized(),
"Matrix must be finalized.");
3362 os << height <<
'\n';
3363 os <<
width <<
'\n';
3365 for (i = 0; i <=
height; i++)
3370 for (i = 0; i < I[
height]; i++)
3375 for (i = 0; i < I[
height]; i++)
3383 const double MiB = 1024.*1024;
3385 double pz = 100./nnz;
3395 "SparseMatrix statistics:\n"
3398 " Dimensions : " << height <<
" x " <<
width <<
"\n"
3399 " Number of entries (total) : " << nnz <<
"\n"
3400 " Number of entries (per row) : " << 1.*nnz/
Height() <<
"\n"
3401 " Number of stored zeros : " << nz*pz <<
"% (" << nz <<
")\n"
3402 " Number of Inf/Nan entries : " << nnf*pz <<
"% ("<< nnf <<
")\n"
3403 " Norm, max |a_ij| : " << max_norm <<
"\n"
3404 " Symmetry, max |a_ij-a_ji| : " << symm <<
"\n"
3405 " Number of small entries:\n"
3406 " |a_ij| <= 1e-12*Norm : " << ns12*pz <<
"% (" << ns12 <<
")\n"
3407 " |a_ij| <= 1e-15*Norm : " << ns15*pz <<
"% (" << ns15 <<
")\n"
3408 " |a_ij| <= 1e-18*Norm : " << ns18*pz <<
"% (" << ns18 <<
")\n";
3411 os <<
" Memory used by CSR : " <<
3412 (
sizeof(int)*(height+1+nnz)+
sizeof(double)*nnz)/MiB <<
" MiB\n";
3416 size_t used_mem =
sizeof(RowNode*)*height;
3417 #ifdef MFEM_USE_MEMALLOC
3420 for (
int i = 0; i <
height; i++)
3422 for (RowNode *aux =
Rows[i]; aux != NULL; aux = aux->Prev)
3424 used_mem +=
sizeof(RowNode);
3428 os <<
" Memory used by LIL : " << used_mem/MiB <<
" MiB\n";
3440 #if !defined(MFEM_USE_MEMALLOC)
3441 for (
int i = 0; i <
height; i++)
3443 RowNode *aux, *node_p =
Rows[i];
3444 while (node_p != NULL)
3447 node_p = node_p->Prev;
3457 #ifdef MFEM_USE_MEMALLOC
3470 const int *start_j =
J;
3471 const int *end_j = J + I[
height];
3472 for (
const int *jptr = start_j; jptr != end_j; ++jptr)
3474 awidth = std::max(awidth, *jptr + 1);
3480 for (
int i = 0; i <
height; i++)
3482 for (aux =
Rows[i]; aux != NULL; aux = aux->Prev)
3484 awidth = std::max(awidth, aux->Column + 1);
3496 for (
int i = 0; i < n; i++)
3506 "Finalize must be called before Transpose. Use TransposeRowMatrix instead");
3509 const int *A_i, *A_j;
3510 int m, n, nnz, *At_i, *At_j;
3511 const double *A_data;
3525 for (i = 0; i <= n; i++)
3529 for (i = 0; i < nnz; i++)
3533 for (i = 1; i < n; i++)
3535 At_i[i+1] += At_i[i];
3538 for (i = j = 0; i < m; i++)
3541 for ( ; j < end; j++)
3543 At_j[At_i[A_j[j]]] = i;
3544 At_data[At_i[A_j[j]]] = A_data[j];
3549 for (i = n; i > 0; i--)
3551 At_i[i] = At_i[i-1];
3562 int m, n, nnz, *At_i, *At_j;
3572 for (i = 0; i < m; i++)
3574 A.
GetRow(i, Acols, Avals);
3596 for (i = 0; i <= n; i++)
3601 for (i = 0; i < m; i++)
3603 A.
GetRow(i, Acols, Avals);
3604 for (j = 0; j<Acols.
Size(); ++j)
3609 for (i = 1; i < n; i++)
3611 At_i[i+1] += At_i[i];
3614 for (i = 0; i < m; i++)
3616 A.
GetRow(i, Acols, Avals);
3617 for (j = 0; j<Acols.
Size(); ++j)
3619 At_j[At_i[Acols[j]]] = i;
3620 At_data[At_i[Acols[j]]] = Avals[j];
3625 for (i = n; i > 0; i--)
3627 At_i[i] = At_i[i-1];
3638 int nrowsA, ncolsA, nrowsB, ncolsB;
3639 const int *A_i, *A_j, *B_i, *B_j;
3640 int *C_i, *C_j, *B_marker;
3641 const double *A_data, *B_data;
3643 int ia, ib, ic, ja, jb, num_nonzeros;
3644 int row_start, counter;
3645 double a_entry, b_entry;
3653 MFEM_VERIFY(ncolsA == nrowsB,
3654 "number of columns of A (" << ncolsA
3655 <<
") must equal number of rows of B (" << nrowsB <<
")");
3664 B_marker =
new int[ncolsB];
3666 for (ib = 0; ib < ncolsB; ib++)
3675 C_i[0] = num_nonzeros = 0;
3676 for (ic = 0; ic < nrowsA; ic++)
3678 for (ia = A_i[ic]; ia < A_i[ic+1]; ia++)
3681 for (ib = B_i[ja]; ib < B_i[ja+1]; ib++)
3684 if (B_marker[jb] != ic)
3691 C_i[ic+1] = num_nonzeros;
3697 C =
new SparseMatrix(C_i, C_j, C_data, nrowsA, ncolsB);
3699 for (ib = 0; ib < ncolsB; ib++)
3708 MFEM_VERIFY(nrowsA == C->
Height() && ncolsB == C->
Width(),
3709 "Input matrix sizes do not match output sizes"
3710 <<
" nrowsA = " << nrowsA
3711 <<
", C->Height() = " << C->
Height()
3712 <<
" ncolsB = " << ncolsB
3713 <<
", C->Width() = " << C->
Width());
3721 for (ic = 0; ic < nrowsA; ic++)
3724 row_start = counter;
3725 for (ia = A_i[ic]; ia < A_i[ic+1]; ia++)
3728 a_entry = A_data[ia];
3729 for (ib = B_i[ja]; ib < B_i[ja+1]; ib++)
3732 b_entry = B_data[ib];
3733 if (B_marker[jb] < row_start)
3735 B_marker[jb] = counter;
3740 C_data[counter] = a_entry*b_entry;
3745 C_data[B_marker[jb]] += a_entry*b_entry;
3753 "With pre-allocated output matrix, number of non-zeros ("
3755 <<
") did not match number of entries changed from matrix-matrix multiply, "
3774 int nrowsA, ncolsA, nrowsB, ncolsB;
3775 int *C_i, *C_j, *B_marker;
3777 int ia, ib, ic, ja, jb, num_nonzeros;
3778 int row_start, counter;
3779 double a_entry, b_entry;
3787 MFEM_VERIFY(ncolsA == nrowsB,
3788 "number of columns of A (" << ncolsA
3789 <<
") must equal number of rows of B (" << nrowsB <<
")");
3791 B_marker =
new int[ncolsB];
3793 for (ib = 0; ib < ncolsB; ib++)
3800 C_i[0] = num_nonzeros = 0;
3804 for (ic = 0; ic < nrowsA; ic++)
3806 A.
GetRow(ic, colsA, dataA);
3807 for (ia = 0; ia < colsA.
Size(); ia++)
3810 B.
GetRow(ja, colsB, dataB);
3811 for (ib = 0; ib < colsB.
Size(); ib++)
3814 if (B_marker[jb] != ic)
3821 C_i[ic+1] = num_nonzeros;
3827 C =
new SparseMatrix(C_i, C_j, C_data, nrowsA, ncolsB);
3829 for (ib = 0; ib < ncolsB; ib++)
3835 for (ic = 0; ic < nrowsA; ic++)
3837 row_start = counter;
3838 A.
GetRow(ic, colsA, dataA);
3839 for (ia = 0; ia < colsA.
Size(); ia++)
3842 a_entry = dataA[ia];
3843 B.
GetRow(ja, colsB, dataB);
3844 for (ib = 0; ib < colsB.
Size(); ib++)
3847 b_entry = dataB[ib];
3848 if (B_marker[jb] < row_start)
3850 B_marker[jb] = counter;
3852 C_data[counter] = a_entry*b_entry;
3857 C_data[B_marker[jb]] += a_entry*b_entry;
3872 for (
int j = 0; j < B.
Width(); ++j)
3876 A.
Mult(columnB, columnC);
3886 Mult (R, *AP, *RAP_);
3929 int i, At_nnz, *At_j;
3933 At_nnz = At -> NumNonZeroElems();
3934 At_j = At -> GetJ();
3935 At_data = At -> GetData();
3936 for (i = 0; i < At_nnz; i++)
3938 At_data[i] *= D(At_j[i]);
3949 int ncols = A.
Width();
3955 const int *A_i = A.
GetI();
3956 const int *A_j = A.
GetJ();
3957 const double *A_data = A.
GetData();
3959 const int *B_i = B.
GetI();
3960 const int *B_j = B.
GetJ();
3961 const double *B_data = B.
GetData();
3963 int * marker =
new int[ncols];
3964 std::fill(marker, marker+ncols, -1);
3966 int num_nonzeros = 0, jcol;
3968 for (
int ic = 0; ic < nrows; ic++)
3970 for (
int ia = A_i[ic]; ia < A_i[ic+1]; ia++)
3976 for (
int ib = B_i[ic]; ib < B_i[ic+1]; ib++)
3979 if (marker[jcol] != ic)
3985 C_i[ic+1] = num_nonzeros;
3991 for (
int ia = 0; ia < ncols; ia++)
3997 for (
int ic = 0; ic < nrows; ic++)
3999 for (
int ia = A_i[ic]; ia < A_i[ic+1]; ia++)
4003 C_data[pos] = a*A_data[ia];
4007 for (
int ib = B_i[ic]; ib < B_i[ic+1]; ib++)
4010 if (marker[jcol] < C_i[ic])
4013 C_data[pos] = b*B_data[ib];
4019 C_data[marker[jcol]] += b*B_data[ib];
4025 return new SparseMatrix(C_i, C_j, C_data, nrows, ncols);
4030 return Add(1.,A,1.,B);
4035 MFEM_ASSERT(Ai.
Size() > 0,
"invalid size Ai.Size() = " << Ai.
Size());
4040 for (
int i=1; i < Ai.
Size(); ++i)
4042 result =
Add(*accumulate, *Ai[i]);
4048 accumulate = result;
4058 for (
int r = 0; r < B.
Height(); r++)
4062 for (
int i=0; i<A.
RowSize(r); i++)
4064 B(r, colA[i]) += alpha * valA[i];
4077 for (
int i=0; i<mA; i++)
4079 for (
int j=0; j<nA; j++)
4081 C->
AddMatrix(A(i,j), B, i * mB, j * nB);
4095 for (
int i=0; i<mA; i++)
4097 for (
int j=0; j<nA; j++)
4099 for (
int r=0; r<mB; r++)
4104 for (
int cj=0; cj<B.
RowSize(r); cj++)
4106 C->
Set(i * mB + r, j * nB + colB[cj], A(i,j) * valB[cj]);
4124 for (
int r=0; r<mA; r++)
4129 for (
int aj=0; aj<A.
RowSize(r); aj++)
4131 for (
int i=0; i<mB; i++)
4133 for (
int j=0; j<nB; j++)
4135 C->
Set(r * mB + i, colA[aj] * nB + j, valA[aj] * B(i, j));
4153 for (
int ar=0; ar<mA; ar++)
4158 for (
int aj=0; aj<A.
RowSize(ar); aj++)
4160 for (
int br=0; br<mB; br++)
4165 for (
int bj=0; bj<B.
RowSize(br); bj++)
4167 C->
Set(ar * mB + br, colA[aj] * nB + colB[bj],
4168 valA[aj] * valB[bj]);
4191 #ifdef MFEM_USE_MEMALLOC
4201 #ifdef MFEM_USE_CUDA_OR_HIP
4208 MFEM_cu_or_hip(sparseDestroy)(
handle);
4213 MFEM_Cu_or_Hip(MemFree)(
dBuffer);
4220 #endif // MFEM_USE_CUDA_OR_HIP
Memory< int > I
Array with size (height+1) containing the row offsets.
double GetJacobiScaling() const
Determine appropriate scaling for Jacobi iteration.
int CheckFinite() const
Count the number of entries that are NOT finite, i.e. Inf or Nan.
int Size() const
Return the logical size of the array.
void DiagScale(const Vector &b, Vector &x, double sc=1.0, bool use_abs_diag=false) const
x = sc b / A_ii. When use_abs_diag = true, |A_ii| is used.
SparseMatrix * TransposeAbstractSparseMatrix(const AbstractSparseMatrix &A, int useActualWidth)
Transpose of a sparse matrix. A does not need to be a CSR matrix.
int RowSize(const int i) const
Returns the number of elements in row i.
int CheckFinite(const double *v, const int n)
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
void * CuMemFree(void *dptr)
Frees device memory and returns destination ptr.
void _Add_(const int col, const double a)
Add a value to an entry in the "current row". See SetColPtr().
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
void Clear()
Clear the contents of the SparseMatrix.
void EliminateCol(int col, DiagonalPolicy dpolicy=DIAG_ZERO)
Eliminates the column col from the matrix.
Memory< T > & GetMemory()
Return a reference to the Memory object used by the Array.
SparseMatrix * At
Transpose of A. Owned. Used to perform MultTranspose() on devices.
size_t MemoryUsage() const
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
void MakeRef(const SparseMatrix &master)
Clear the contents of the SparseMatrix and make it a reference to master.
void SetSize(int s)
Resize the vector to size s.
double & SearchRow(const int col)
Perform a fast search for an entry in the "current row". See SetColPtr().
void Delete()
Delete the owned pointers and reset the Memory object.
DenseMatrix * ToDenseMatrix() const
Produces a DenseMatrix from a SparseMatrix.
bool Empty() const
Check if the SparseMatrix is empty.
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
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().
const double * ReadData(bool on_dev=true) const
int * GetJ()
Return the array J.
Biwise-OR of all HIP backends.
Abstract data type for sparse matrices.
SparseMatrix & operator+=(const SparseMatrix &B)
Add the sparse matrix 'B' to '*this'. This operation will cause an error if '*this' is finalized and ...
T * Write(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for write access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true...
Data type dense matrix using column-major storage.
void AddMult(const Vector &x, Vector &y, const double a=1.0) const
y += A * x (default) or y += a * A * x
int Size() const
Returns the size of the vector.
int * GetI()
Return the array I.
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
y += At * x (default) or y += a * At * x
Abstract data type for matrix inverse.
void EliminateRowColMultipleRHS(int rc, const Vector &sol, DenseMatrix &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Similar to EliminateRowCol(int, const double, Vector &, DiagonalPolicy), but multiple values for elim...
virtual int NumNonZeroElems() const =0
Returns the number of non-zeros in a matrix.
void ClearGPUSparse()
Clear the cuSPARSE/hipSPARSE descriptors. This must be called after releasing the device memory of A...
void PrintInfo(std::ostream &out) const
Print various sparse matrix statistics.
void CopyFrom(const Memory &src, int size)
Copy size entries from src to *this.
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc, bool use_abs_diag=false) const
void AddRow(const int row, const Array< int > &cols, const Vector &srow)
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
double * GetData() const
Return a pointer to the beginning of the Vector data.
void SortColumnIndices()
Sort the column indices corresponding to each row.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
void PrintMM(std::ostream &out=mfem::out) const
Prints matrix in Matrix Market sparse format.
int Capacity() const
Return the size of the allocated memory.
cusparseDnVecDescr_t vecY_descr
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
double * HostReadWriteData()
void _Set_(const int col, const double a)
Set an entry in the "current row". See SetColPtr().
void BooleanMultTranspose(const Array< int > &x, Array< int > &y) const
y = At * x, treating all entries as booleans (zero=false, nonzero=true).
SparseMatrix * Mult_AtDA(const SparseMatrix &A, const Vector &D, SparseMatrix *OAtDA)
Matrix multiplication A^t D A. All matrices must be finalized.
void PartAddMult(const Array< int > &rows, const Vector &x, Vector &y, const double a=1.0) const
bool RowIsEmpty(const int row) const
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
unsigned flags
Bit flags defined from the FlagMask enum.
void GetBlocks(Array2D< SparseMatrix * > &blocks) const
const int * ReadJ(bool on_dev=true) const
void ScaleRow(const int row, const double scale)
double * GetData()
Return the element data, i.e. the array A.
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
bool IsFinite(const double &val)
void MoveDiagonalFirst()
Move the diagonal entry to the first position in each row, preserving the order of the rest of the co...
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
void Set(const int i, const int j, const double val)
void Wrap(T *ptr, int size, bool own)
Wrap an externally allocated host pointer, ptr with the current host memory type returned by MemoryMa...
void ClearColPtr() const
Reset the "current row" set by calling SetColPtr(). This method must be called between any two calls ...
void ScaleColumns(const Vector &sr)
this = this * diag(sr);
virtual ~SparseMatrix()
Destroys sparse matrix.
void PrintCSR2(std::ostream &out) const
Prints a sparse matrix to stream out in CSR format.
double f(const Vector &xvec)
virtual int GetRow(const int row, Array< int > &cols, Vector &srow) const
Extract all column indices and values from a given row.
bool isSorted
Are the columns sorted already.
Memory< double > A
Array with size I[height], containing the actual entries of the sparse matrix, as indexed by the I ar...
SparseMatrix()
Create an empty SparseMatrix.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void EnsureMultTranspose() const
Ensures that the matrix is capable of performing MultTranspose(), AddMultTranspose(), and AbsMultTranspose().
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void SetWidth(int width_=-1)
Change the width of a SparseMatrix.
const double * HostReadData() const
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
double IsSymmetric() const
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.
void * HipMemFree(void *dptr)
Frees device memory.
const int * HostReadJ() const
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
static int SparseMatrixCount
void SetDiagIdentity()
If a row contains only one diag entry of zero, set it to 1.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
void ScaleRows(const Vector &sl)
this = diag(sl) * this;
static cusparseHandle_t handle
DenseMatrix * OuterProduct(const DenseMatrix &A, const DenseMatrix &B)
Produces a block matrix with blocks A_{ij}*B.
void Sort()
Sorts the array in ascending order. This requires operator< to be defined for T.
virtual int GetRow(const int row, Array< int > &cols, Vector &srow) const =0
Gets the columns indexes and values for row row.
Set the diagonal value to one.
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
virtual void PrintMatlab(std::ostream &out=mfem::out) const
Prints matrix in matlab format.
const T * Read(const Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true...
void Reset()
Reset the memory to be empty, ensuring that Delete() will be a no-op.
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Biwise-OR of all CUDA backends.
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
void Add(const int i, const int j, const double val)
Dynamic 2D array using row-major layout.
void SetSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
bool Finalized() const
Returns whether or not CSR format has been finalized.
Biwise-OR of all CPU backends.
void * HipMemAlloc(void **dptr, size_t bytes)
Allocates device memory.
T * HostWrite(Memory< T > &mem, int size)
Shortcut to Write(const Memory<T> &mem, int size, false)
void SparseMatrixFunction(SparseMatrix &S, double(*f)(double))
Applies f() to each element of the matrix (after it is finalized).
virtual MatrixInverse * Inverse() const
This virtual method is not supported: it always returns NULL.
void Swap(Array< T > &, Array< T > &)
void ClearOwnerFlags() const
Clear the ownership flags for the host and device pointers, as well as any internal data allocated by...
void EliminateRowColDiag(int rc, double value)
Perform elimination and set the diagonal entry to the given value.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void SetColPtr(const int row) const
Initialize the SparseMatrix for fast access to the entries of the given row which becomes the "curren...
MemoryType
Memory types supported by MFEM.
void SetSubMatrixTranspose(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
void SetHostPtrOwner(bool own) const
Set/clear the ownership flag for the host pointer. Ownership indicates whether the pointer will be de...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void EliminateRowCol(int rc, const double sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate row rc and column rc and modify the rhs using sol.
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
void GetColumnReference(int c, Vector &col)
void Threshold(double tol, bool fix_empty_rows=false)
Remove entries smaller in absolute value than a given tolerance tol. If fix_empty_rows is true...
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0<=i<A.Height, 0<=j<A.Width.
void Gauss_Seidel_back(const Vector &x, Vector &y) const
const int * ReadI(bool on_dev=true) const
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
void Gauss_Seidel_forw(const Vector &x, Vector &y) const
Gauss-Seidel forward and backward iterations over a vector x.
RowNode ** Rows
Array of linked lists, one for every row. This array represents the linked list (LIL) storage format...
int height
Dimension of the output / number of rows in the matrix.
void EliminateRow(int row, const double sol, Vector &rhs)
Eliminates a column from the transpose matrix.
const T * HostRead(const Memory< T > &mem, int size)
Shortcut to Read(const Memory<T> &mem, int size, false)
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
int ActualWidth() const
Returns the actual Width of the matrix.
void GetDiag(Vector &d) const
Returns the Diagonal of A.
MemAlloc< RowNode, 1024 > RowNodeAlloc
double GetRowNorml1(int irow) const
For i = irow compute .
void ResetTranspose() const
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
virtual void EliminateZeroRows(const double threshold=1e-12)
If a row contains only zeros, set its diagonal to 1.
void AbsMultTranspose(const Vector &x, Vector &y) const
y = |At| * x, using entry-wise absolute values of the transpose of matrix A
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
double & operator()(int i, int j)
Returns reference to A[i][j].
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
void EliminateCols(const Array< int > &cols, const Vector *x=NULL, Vector *b=NULL)
Eliminate all columns i for which cols[i] != 0.
void OverrideSize(int height_, int width_)
Sets the height and width of the matrix.
SparseMatrix & operator=(const SparseMatrix &rhs)
Assignment operator: deep copy.
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
SparseMatrix & operator*=(double a)
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
cusparseSpMatDescr_t matA_descr
T * Write(bool on_dev=true)
Shortcut for mfem::Write(a.GetMemory(), a.Size(), on_dev).
void MakeRef(T *, int)
Make this Array a reference to a pointer.
double * ReadWriteData(bool on_dev=true)
int * ReadWriteJ(bool on_dev=true)
void Print(std::ostream &out=mfem::out, int width_=4) const
Prints matrix to stream out.
const int * HostReadI() const
void Jacobi2(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
cusparseDnVecDescr_t vecX_descr
int MaxRowSize() const
Returns the maximum number of elements among all rows.
void Jacobi3(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
SparseMatrix * MultAbstractSparseMatrix(const AbstractSparseMatrix &A, const AbstractSparseMatrix &B)
Matrix product of sparse matrices. A and B do not need to be CSR matrices.
void AbsMult(const Vector &x, Vector &y) const
y = |A| * x, using entry-wise absolute values of matrix A
double _Get_(const int col) const
Read the value of an entry in the "current row". See SetColPtr().
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Memory< int > J
Array with size I[height], containing the column indices for all matrix entries, as indexed by the I ...
void Swap(SparseMatrix &other)
void PrintCSR(std::ostream &out) const
Prints matrix to stream out in hypre_CSRMatrix format.
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
void * CuMemAlloc(void **dptr, size_t bytes)
Allocates device memory and returns destination ptr.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Set the diagonal value to zero.
int width
Dimension of the input / number of columns in the matrix.
void EliminateBC(const Array< int > &ess_dofs, DiagonalPolicy diag_policy)
Eliminate essential (Dirichlet) boundary conditions.
int CountSmallElems(double tol) const
Count the number of entries with |a_ij| <= tol.
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
void GetRowSums(Vector &x) const
For all i compute .
void Neg()
(*this) = -(*this)