MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
sparsemat.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 // Implementation of sparse matrix
13 
14 #include "linalg.hpp"
15 #include "../general/forall.hpp"
16 #include "../general/table.hpp"
17 #include "../general/sort_pairs.hpp"
18 #include "../general/backends.hpp"
19 
20 #include <iostream>
21 #include <iomanip>
22 #include <cmath>
23 #include <algorithm>
24 #include <limits>
25 #include <cstring>
26 
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
32 
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
38 
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
44 
45 // https://hipsparse.readthedocs.io/en/latest/usermanual.html#hipsparsespmvalg-t
46 #define MFEM_GPUSPARSE_ALG HIPSPARSE_CSRMV_ALG1
47 #endif // defined(MFEM_USE_CUDA)
48 
49 namespace mfem
50 {
51 
52 using namespace std;
53 
54 #ifdef MFEM_USE_CUDA_OR_HIP
56 // doxygen doesn't like the macro-assisted typename so let's skip parsing it:
57 // \cond false
58 MFEM_cu_or_hip(sparseHandle_t) SparseMatrix::handle = nullptr;
59 // \endcond
60 size_t SparseMatrix::bufferSize = 0;
61 void * SparseMatrix::dBuffer = nullptr;
62 #endif // MFEM_USE_CUDA_OR_HIP
63 
65 {
66  // Initialize cuSPARSE/hipSPARSE library
67 #ifdef MFEM_USE_CUDA_OR_HIP
69  {
70  if (!handle) { MFEM_cu_or_hip(sparseCreate)(&handle); }
71  useGPUSparse=true;
72  SparseMatrixCount++;
73  }
74  else
75  {
76  useGPUSparse=false;
77  }
78 #endif // MFEM_USE_CUDA_OR_HIP
79 }
80 
82 {
83 #ifdef MFEM_USE_CUDA_OR_HIP
84  if (initBuffers)
85  {
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);
90 #else
91  cusparseDestroyMatDescr(matA_descr);
92 #endif // CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
93  initBuffers = false;
94  }
95 #endif // MFEM_USE_CUDA_OR_HIP
96 }
97 
98 SparseMatrix::SparseMatrix(int nrows, int ncols)
99  : AbstractSparseMatrix(nrows, (ncols >= 0) ? ncols : nrows),
100  Rows(new RowNode *[nrows]),
101  current_row(-1),
102  ColPtrJ(NULL),
103  ColPtrNode(NULL),
104  At(NULL),
105  isSorted(false)
106 {
107  // We probably do not need to set the ownership flags here.
108  I.SetHostPtrOwner(true);
109  J.SetHostPtrOwner(true);
110  A.SetHostPtrOwner(true);
111 
112  for (int i = 0; i < nrows; i++)
113  {
114  Rows[i] = NULL;
115  }
116 
117 #ifdef MFEM_USE_MEMALLOC
118  NodesMem = new RowNodeAlloc;
119 #endif
120 
121  InitGPUSparse();
122 }
123 
124 SparseMatrix::SparseMatrix(int *i, int *j, double *data, int m, int n)
125  : AbstractSparseMatrix(m, n),
126  Rows(NULL),
127  ColPtrJ(NULL),
128  ColPtrNode(NULL),
129  At(NULL),
130  isSorted(false)
131 {
132  I.Wrap(i, height+1, true);
133  J.Wrap(j, I[height], true);
134  A.Wrap(data, I[height], true);
135 
136 #ifdef MFEM_USE_MEMALLOC
137  NodesMem = NULL;
138 #endif
139 
140  InitGPUSparse();
141 }
142 
143 SparseMatrix::SparseMatrix(int *i, int *j, double *data, int m, int n,
144  bool ownij, bool owna, bool issorted)
145  : AbstractSparseMatrix(m, n),
146  Rows(NULL),
147  ColPtrJ(NULL),
148  ColPtrNode(NULL),
149  At(NULL),
150  isSorted(issorted)
151 {
152  I.Wrap(i, height+1, ownij);
153  J.Wrap(j, I[height], ownij);
154 
155 #ifdef MFEM_USE_MEMALLOC
156  NodesMem = NULL;
157 #endif
158 
159  if (data)
160  {
161  A.Wrap(data, I[height], owna);
162  }
163  else
164  {
165  const int nnz = I[height];
166  A.New(nnz);
167  for (int ii=0; ii<nnz; ++ii)
168  {
169  A[ii] = 0.0;
170  }
171  }
172 
173  InitGPUSparse();
174 }
175 
176 SparseMatrix::SparseMatrix(int nrows, int ncols, int rowsize)
177  : AbstractSparseMatrix(nrows, ncols)
178  , Rows(NULL)
179  , ColPtrJ(NULL)
180  , ColPtrNode(NULL)
181  , At(NULL)
182  , isSorted(false)
183 {
184 #ifdef MFEM_USE_MEMALLOC
185  NodesMem = NULL;
186 #endif
187  I.New(nrows + 1);
188  J.New(nrows * rowsize);
189  A.New(nrows * rowsize);
190 
191  for (int i = 0; i <= nrows; i++)
192  {
193  I[i] = i * rowsize;
194  }
195 
196  InitGPUSparse();
197 }
198 
199 SparseMatrix::SparseMatrix(const SparseMatrix &mat, bool copy_graph,
200  MemoryType mt)
201  : AbstractSparseMatrix(mat.Height(), mat.Width())
202 {
203  if (mat.Finalized())
204  {
205  mat.HostReadI();
206  const int nnz = mat.I[height];
207  if (copy_graph)
208  {
209  I.New(height+1, mt == MemoryType::PRESERVE ? mat.I.GetMemoryType() : mt);
210  J.New(nnz, mt == MemoryType::PRESERVE ? mat.J.GetMemoryType() : mt);
211  I.CopyFrom(mat.I, height+1);
212  J.CopyFrom(mat.J, nnz);
213  }
214  else
215  {
216  I = mat.I;
217  J = mat.J;
218  I.ClearOwnerFlags();
219  J.ClearOwnerFlags();
220  }
221  A.New(nnz, mt == MemoryType::PRESERVE ? mat.A.GetMemoryType() : mt);
222  A.CopyFrom(mat.A, nnz);
223 
224  Rows = NULL;
225 #ifdef MFEM_USE_MEMALLOC
226  NodesMem = NULL;
227 #endif
228  }
229  else
230  {
231 #ifdef MFEM_USE_MEMALLOC
232  NodesMem = new RowNodeAlloc;
233 #endif
234  Rows = new RowNode *[height];
235  for (int i = 0; i < height; i++)
236  {
237  RowNode **node_pp = &Rows[i];
238  for (RowNode *node_p = mat.Rows[i]; node_p; node_p = node_p->Prev)
239  {
240 #ifdef MFEM_USE_MEMALLOC
241  RowNode *new_node_p = NodesMem->Alloc();
242 #else
243  RowNode *new_node_p = new RowNode;
244 #endif
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;
249  }
250  *node_pp = NULL;
251  }
252 
253  // We probably do not need to set the ownership flags here.
254  I.SetHostPtrOwner(true);
255  J.SetHostPtrOwner(true);
256  A.SetHostPtrOwner(true);
257  }
258 
259  current_row = -1;
260  ColPtrJ = NULL;
261  ColPtrNode = NULL;
262  At = NULL;
263  isSorted = mat.isSorted;
264 
265  InitGPUSparse();
266 }
267 
269  : AbstractSparseMatrix(v.Size(), v.Size())
270  , Rows(NULL)
271  , ColPtrJ(NULL)
272  , ColPtrNode(NULL)
273  , At(NULL)
274  , isSorted(true)
275 {
276 #ifdef MFEM_USE_MEMALLOC
277  NodesMem = NULL;
278 #endif
279  I.New(height + 1);
280  J.New(height);
281  A.New(height);
282 
283  for (int i = 0; i <= height; i++)
284  {
285  I[i] = i;
286  }
287 
288  for (int r=0; r<height; r++)
289  {
290  J[r] = r;
291  A[r] = v[r];
292  }
293 
294  InitGPUSparse();
295 }
296 
297 void SparseMatrix::OverrideSize(int height_, int width_)
298 {
299  height = height_;
300  width = width_;
301 }
302 
304 {
305  Clear();
306 
307  SparseMatrix copy(rhs);
308  Swap(copy);
309 
310  return *this;
311 }
312 
314 {
315  MFEM_ASSERT(master.Finalized(), "'master' must be finalized");
316  Clear();
317  height = master.Height();
318  width = master.Width();
319  I = master.I; I.ClearOwnerFlags();
320  J = master.J; J.ClearOwnerFlags();
321  A = master.A; A.ClearOwnerFlags();
322  isSorted = master.isSorted;
323 }
324 
326 {
327  height = width = 0;
328  I.Reset();
329  J.Reset();
330  A.Reset();
331  Rows = NULL;
332  current_row = -1;
333  ColPtrJ = NULL;
334  ColPtrNode = NULL;
335  At = NULL;
336 #ifdef MFEM_USE_MEMALLOC
337  NodesMem = NULL;
338 #endif
339  isSorted = false;
340 
341  ClearGPUSparse();
342 }
343 
344 int SparseMatrix::RowSize(const int i) const
345 {
346  int gi = i;
347  if (gi < 0)
348  {
349  gi = -1-gi;
350  }
351 
352  if (I)
353  {
354  return I[gi+1]-I[gi];
355  }
356 
357  int s = 0;
358  RowNode *row = Rows[gi];
359  for ( ; row != NULL; row = row->Prev)
360  if (row->Value != 0.0)
361  {
362  s++;
363  }
364  return s;
365 }
366 
368 {
369  int max_row_size=0;
370  int rowSize=0;
371  if (I)
372  {
373  for (int i=0; i < height; ++i)
374  {
375  rowSize = I[i+1]-I[i];
376  max_row_size = (max_row_size > rowSize) ? max_row_size : rowSize;
377  }
378  }
379  else
380  {
381  for (int i=0; i < height; ++i)
382  {
383  rowSize = RowSize(i);
384  max_row_size = (max_row_size > rowSize) ? max_row_size : rowSize;
385  }
386  }
387 
388  return max_row_size;
389 }
390 
391 int *SparseMatrix::GetRowColumns(const int row)
392 {
393  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
394 
395  return J + I[row];
396 }
397 
398 const int *SparseMatrix::GetRowColumns(const int row) const
399 {
400  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
401 
402  return J + I[row];
403 }
404 
405 double *SparseMatrix::GetRowEntries(const int row)
406 {
407  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
408 
409  return A + I[row];
410 }
411 
412 const double *SparseMatrix::GetRowEntries(const int row) const
413 {
414  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
415 
416  return A + I[row];
417 }
418 
419 void SparseMatrix::SetWidth(int newWidth)
420 {
421  if (newWidth == width)
422  {
423  // Nothing to be done here
424  return;
425  }
426  else if (newWidth == -1)
427  {
428  // Compute the actual width
429  width = ActualWidth();
430  // No need to reset the ColPtr, since the new ColPtr will be shorter.
431  }
432  else if (newWidth > width)
433  {
434  // We need to reset ColPtr, since now we may have additional columns.
435  if (Rows != NULL)
436  {
437  delete [] ColPtrNode;
438  ColPtrNode = static_cast<RowNode **>(NULL);
439  }
440  else
441  {
442  delete [] ColPtrJ;
443  ColPtrJ = static_cast<int *>(NULL);
444  }
445  width = newWidth;
446  }
447  else
448  {
449  // Check that the new width is bigger or equal to the actual width.
450  MFEM_ASSERT(newWidth >= ActualWidth(),
451  "The new width needs to be bigger or equal to the actual width");
452  width = newWidth;
453  }
454 }
455 
456 
458 {
459  MFEM_VERIFY(Finalized(), "Matrix is not Finalized!");
460 
461  if (isSorted)
462  {
463  return;
464  }
465 
466 #ifdef MFEM_USE_CUDA_OR_HIP
468  {
469 #if defined(MFEM_USE_CUDA)
470  size_t pBufferSizeInBytes = 0;
471  void *pBuffer = NULL;
472 
473  const int n = Height();
474  const int m = Width();
475  const int nnzA = J.Capacity();
476  double * d_a_sorted = ReadWriteData();
477  const int * d_ia = ReadI();
478  int * d_ja_sorted = ReadWriteJ();
479  csru2csrInfo_t sortInfoA;
480 
481  cusparseMatDescr_t matA_descr;
482  cusparseCreateMatDescr( &matA_descr );
483  cusparseSetMatIndexBase( matA_descr, CUSPARSE_INDEX_BASE_ZERO );
484  cusparseSetMatType( matA_descr, CUSPARSE_MATRIX_TYPE_GENERAL );
485 
486  cusparseCreateCsru2csrInfo( &sortInfoA );
487 
488  cusparseDcsru2csr_bufferSizeExt( handle, n, m, nnzA, d_a_sorted, d_ia,
489  d_ja_sorted, sortInfoA,
490  &pBufferSizeInBytes);
491 
492  CuMemAlloc( &pBuffer, pBufferSizeInBytes );
493 
494  cusparseDcsru2csr( handle, n, m, nnzA, matA_descr, d_a_sorted, d_ia,
495  d_ja_sorted, sortInfoA, pBuffer);
496 
497  cusparseDestroyCsru2csrInfo( sortInfoA );
498  cusparseDestroyMatDescr( matA_descr );
499 
500  CuMemFree( pBuffer );
501 #endif
502  }
503  else if ( Device::Allows( Backend::HIP_MASK ))
504  {
505 #if defined(MFEM_USE_HIP)
506  size_t pBufferSizeInBytes = 0;
507  void *pBuffer = NULL;
508  int *P = NULL;
509 
510  const int n = Height();
511  const int m = Width();
512  const int nnzA = J.Capacity();
513  double * d_a_sorted = ReadWriteData();
514  const int * d_ia = ReadI();
515  int * d_ja_sorted = ReadWriteJ();
516 
517  hipsparseMatDescr_t descrA;
518  hipsparseCreateMatDescr( &descrA );
519  // FIXME: There is not in-place version of csr sort in hipSPARSE currently, so we make
520  // a temporary copy of the data for gthr, sort that, and then copy the sorted values
521  // back to the array being returned. Where there is an in-place version available,
522  // we should use it.
523  Array< double > a_tmp( nnzA );
524  double *d_a_tmp = a_tmp.Write();
525 
526  hipsparseXcsrsort_bufferSizeExt(handle, n, m, nnzA, d_ia, d_ja_sorted,
527  &pBufferSizeInBytes);
528 
529  HipMemAlloc( &pBuffer, pBufferSizeInBytes );
530  HipMemAlloc( &P, nnzA * sizeof(int) );
531 
532  hipsparseCreateIdentityPermutation(handle, nnzA, P);
533  hipsparseXcsrsort(handle, n, m, nnzA, descrA, d_ia, d_ja_sorted, P, pBuffer);
534 
535  hipsparseDgthr(handle, nnzA, d_a_sorted, d_a_tmp, P,
536  HIPSPARSE_INDEX_BASE_ZERO);
537 
538  A.CopyFrom( a_tmp.GetMemory(), nnzA );
539  hipsparseDestroyMatDescr( descrA );
540 
541  HipMemFree( pBuffer );
542  HipMemFree( P );
543 #endif
544  }
545  else
546 #endif // MFEM_USE_CUDA_OR_HIP
547  {
548  const int * Ip=HostReadI();
549  HostReadWriteJ();
551 
553  for (int j = 0, i = 0; i < height; i++)
554  {
555  int end = Ip[i+1];
556  row.SetSize(end - j);
557  for (int k = 0; k < row.Size(); k++)
558  {
559  row[k].one = J[j+k];
560  row[k].two = A[j+k];
561  }
562  row.Sort();
563  for (int k = 0; k < row.Size(); k++, j++)
564  {
565  J[j] = row[k].one;
566  A[j] = row[k].two;
567  }
568  }
569  }
570  isSorted = true;
571 }
572 
574 {
575  MFEM_VERIFY(Finalized(), "Matrix is not Finalized!");
576 
577  for (int row = 0, end = 0; row < height; row++)
578  {
579  int start = end, j;
580  end = I[row+1];
581  for (j = start; true; j++)
582  {
583  MFEM_VERIFY(j < end, "diagonal entry not found in row = " << row);
584  if (J[j] == row) { break; }
585  }
586  const double diag = A[j];
587  for ( ; j > start; j--)
588  {
589  J[j] = J[j-1];
590  A[j] = A[j-1];
591  }
592  J[start] = row;
593  A[start] = diag;
594  }
595 }
596 
597 double &SparseMatrix::Elem(int i, int j)
598 {
599  return operator()(i,j);
600 }
601 
602 const double &SparseMatrix::Elem(int i, int j) const
603 {
604  return operator()(i,j);
605 }
606 
607 double &SparseMatrix::operator()(int i, int j)
608 {
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 << ", "
614  << "j = " << j);
615 
616  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
617 
618  for (int k = I[i], end = I[i+1]; k < end; k++)
619  {
620  if (J[k] == j)
621  {
622  return A[k];
623  }
624  }
625 
626  MFEM_ABORT("Did not find i = " << i << ", j = " << j << " in matrix.");
627  return A[0];
628 }
629 
630 const double &SparseMatrix::operator()(int i, int j) const
631 {
632  static const double zero = 0.0;
633 
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 << ", "
639  << "j = " << j);
640 
641  if (Finalized())
642  {
643  for (int k = I[i], end = I[i+1]; k < end; k++)
644  {
645  if (J[k] == j)
646  {
647  return A[k];
648  }
649  }
650  }
651  else
652  {
653  for (RowNode *node_p = Rows[i]; node_p != NULL; node_p = node_p->Prev)
654  {
655  if (node_p->Column == j)
656  {
657  return node_p->Value;
658  }
659  }
660  }
661 
662  return zero;
663 }
664 
666 {
667  MFEM_VERIFY(height == width, "Matrix must be square, not height = "
668  << height << ", width = " << width);
669  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
670 
671  d.SetSize(height);
672 
673  const auto II = this->ReadI();
674  const auto JJ = this->ReadJ();
675  const auto AA = this->ReadData();
676  auto dd = d.Write();
677 
678  MFEM_FORALL(i, height,
679  {
680  const int begin = II[i];
681  const int end = II[i+1];
682  int j;
683  for (j = begin; j < end; j++)
684  {
685  if (JJ[j] == i)
686  {
687  dd[i] = AA[j];
688  break;
689  }
690  }
691  if (j == end)
692  {
693  dd[i] = 0.;
694  }
695  });
696 }
697 
698 /// Produces a DenseMatrix from a SparseMatrix
700 {
701  int num_rows = this->Height();
702  int num_cols = this->Width();
703 
704  DenseMatrix * B = new DenseMatrix(num_rows, num_cols);
705 
706  this->ToDenseMatrix(*B);
707 
708  return B;
709 }
710 
711 /// Produces a DenseMatrix from a SparseMatrix
713 {
714  B.SetSize(height, width);
715  B = 0.0;
716 
717  for (int r=0; r<height; r++)
718  {
719  const int * col = this->GetRowColumns(r);
720  const double * val = this->GetRowEntries(r);
721 
722  for (int cj=0; cj<this->RowSize(r); cj++)
723  {
724  B(r, col[cj]) = val[cj];
725  }
726  }
727 }
728 
729 void SparseMatrix::Mult(const Vector &x, Vector &y) const
730 {
731  if (Finalized()) { y.UseDevice(true); }
732  y = 0.0;
733  AddMult(x, y);
734 }
735 
736 void SparseMatrix::AddMult(const Vector &x, Vector &y, const double a) const
737 {
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 << ")");
742 
743  if (!Finalized())
744  {
745  const double *xp = x.HostRead();
746  double *yp = y.HostReadWrite();
747 
748  // The matrix is not finalized, but multiplication is still possible
749  for (int i = 0; i < height; i++)
750  {
751  RowNode *row = Rows[i];
752  double b = 0.0;
753  for ( ; row != NULL; row = row->Prev)
754  {
755  b += row->Value * xp[row->Column];
756  }
757  *yp += a * b;
758  yp++;
759  }
760  return;
761  }
762 
763 #ifndef MFEM_USE_LEGACY_OPENMP
764  const int height = this->height;
765  const int nnz = J.Capacity();
766  auto d_I = Read(I, height+1);
767  auto d_J = Read(J, nnz);
768  auto d_A = Read(A, nnz);
769  auto d_x = x.Read();
770  auto d_y = y.ReadWrite();
771 
772  // Skip if matrix has no non-zeros
773  if (nnz == 0) {return;}
775  {
776 #ifdef MFEM_USE_CUDA_OR_HIP
777  const double alpha = a;
778  const double beta = 1.0;
779 
780  // Setup descriptors
781  if (!initBuffers)
782  {
783 #if CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
784  // Setup matrix descriptor
785  MFEM_cu_or_hip(sparseCreateCsr)(
786  &matA_descr,Height(),
787  Width(),
788  J.Capacity(),
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));
796 
797  // Create handles for input/output vectors
798  MFEM_cu_or_hip(sparseCreateDnVec)(&vecX_descr,
799  x.Size(),
800  const_cast<double *>(d_x),
801  MFEM_CUDA_or_HIP(_R_64F));
802  MFEM_cu_or_hip(sparseCreateDnVec)(&vecY_descr, y.Size(), d_y,
803  MFEM_CUDA_or_HIP(_R_64F));
804 #else
805  cusparseCreateMatDescr(&matA_descr);
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)
809  initBuffers = true;
810  }
811  // Allocate kernel space. Buffer is shared between different sparsemats
812  size_t newBufferSize = 0;
813 
814  MFEM_cu_or_hip(sparseSpMV_bufferSize)(
815  handle,
816  MFEM_CU_or_HIP(SPARSE_OPERATION_NON_TRANSPOSE),
817  &alpha,
818  matA_descr,
819  vecX_descr,
820  &beta,
821  vecY_descr,
822  MFEM_CUDA_or_HIP(_R_64F),
823  MFEM_GPUSPARSE_ALG,
824  &newBufferSize);
825 
826  // Check if we need to resize
827  if (newBufferSize > bufferSize)
828  {
829  bufferSize = newBufferSize;
830  if (dBuffer != nullptr) { MFEM_Cu_or_Hip(MemFree)(dBuffer); }
831  MFEM_Cu_or_Hip(MemAlloc)(&dBuffer, bufferSize);
832  }
833 
834 #if CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
835  // Update input/output vectors
836  MFEM_cu_or_hip(sparseDnVecSetValues)(vecX_descr,
837  const_cast<double *>(d_x));
838  MFEM_cu_or_hip(sparseDnVecSetValues)(vecY_descr, d_y);
839 
840  // Y = alpha A * X + beta * Y
841  MFEM_cu_or_hip(sparseSpMV)(
842  handle,
843  MFEM_CU_or_HIP(SPARSE_OPERATION_NON_TRANSPOSE),
844  &alpha,
845  matA_descr,
846  vecX_descr,
847  &beta,
848  vecY_descr,
849  MFEM_CUDA_or_HIP(_R_64F),
850  MFEM_GPUSPARSE_ALG,
851  dBuffer);
852 #else
853  cusparseDcsrmv(handle,
854  CUSPARSE_OPERATION_NON_TRANSPOSE,
855  Height(),
856  Width(),
857  J.Capacity(),
858  &alpha,
859  matA_descr,
860  const_cast<double *>(d_A),
861  const_cast<int *>(d_I),
862  const_cast<int *>(d_J),
863  const_cast<double *>(d_x),
864  &beta,
865  d_y);
866 #endif // CUDA_VERSION >= 10010 || defined(MFEM_USE_HIP)
867 #endif // MFEM_USE_CUDA_OR_HIP
868  }
869  else
870  {
871  // Native version
872  MFEM_FORALL(i, height,
873  {
874  double d = 0.0;
875  const int end = d_I[i+1];
876  for (int j = d_I[i]; j < end; j++)
877  {
878  d += d_A[j] * d_x[d_J[j]];
879  }
880  d_y[i] += a * d;
881  });
882 
883  }
884 
885 #else // MFEM_USE_LEGACY_OPENMP
886  const double *Ap = A, *xp = x.GetData();
887  double *yp = y.GetData();
888  const int *Jp = J, *Ip = I;
889 
890  #pragma omp parallel for
891  for (int i = 0; i < height; i++)
892  {
893  double d = 0.0;
894  const int end = Ip[i+1];
895  for (int j = Ip[i]; j < end; j++)
896  {
897  d += Ap[j] * xp[Jp[j]];
898  }
899  yp[i] += a * d;
900  }
901 #endif // MFEM_USE_LEGACY_OPENMP
902 }
903 
904 void SparseMatrix::MultTranspose(const Vector &x, Vector &y) const
905 {
906  if (Finalized()) { y.UseDevice(true); }
907  y = 0.0;
908  AddMultTranspose(x, y);
909 }
910 
912  const double a) const
913 {
914  MFEM_ASSERT(height == x.Size(), "Input vector size (" << x.Size()
915  << ") must match matrix height (" << height << ")");
916  MFEM_ASSERT(width == y.Size(), "Output vector size (" << y.Size()
917  << ") must match matrix width (" << width << ")");
918 
919  if (!Finalized())
920  {
921  double *yp = y.HostReadWrite();
922  const double *xp = x.HostRead();
923  // The matrix is not finalized, but multiplication is still possible
924  for (int i = 0; i < height; i++)
925  {
926  RowNode *row = Rows[i];
927  double b = a * xp[i];
928  for ( ; row != NULL; row = row->Prev)
929  {
930  yp[row->Column] += row->Value * b;
931  }
932  }
933  return;
934  }
935 
937  if (At)
938  {
939  At->AddMult(x, y, a);
940  }
941  else
942  {
943  for (int i = 0; i < height; i++)
944  {
945  const double xi = a * x[i];
946  const int end = I[i+1];
947  for (int j = I[i]; j < end; j++)
948  {
949  const int Jj = J[j];
950  y[Jj] += A[j] * xi;
951  }
952  }
953  }
954 }
955 
957 {
958  if (At == NULL)
959  {
960  At = Transpose(*this);
961  }
962 }
963 
965 {
966  delete At;
967  At = NULL;
968 }
969 
971 {
973  {
974  BuildTranspose();
975  }
976 }
977 
979  const Array<int> &rows, const Vector &x, Vector &y) const
980 {
981  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
982 
983  const int n = rows.Size();
984  const int nnz = J.Capacity();
985  auto d_rows = rows.Read();
986  auto d_I = Read(I, height+1);
987  auto d_J = Read(J, nnz);
988  auto d_A = Read(A, nnz);
989  auto d_x = x.Read();
990  auto d_y = y.Write();
991  MFEM_FORALL(i, n,
992  {
993  const int r = d_rows[i];
994  const int end = d_I[r + 1];
995  double a = 0.0;
996  for (int j = d_I[r]; j < end; j++)
997  {
998  a += d_A[j] * d_x[d_J[j]];
999  }
1000  d_y[r] = a;
1001  });
1002 }
1003 
1005  const Array<int> &rows, const Vector &x, Vector &y, const double a) const
1006 {
1007  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
1008 
1009  for (int i = 0; i < rows.Size(); i++)
1010  {
1011  int r = rows[i];
1012  int end = I[r + 1];
1013  double val = 0.0;
1014  for (int j = I[r]; j < end; j++)
1015  {
1016  val += A[j] * x(J[j]);
1017  }
1018  y(r) += a * val;
1019  }
1020 }
1021 
1023 {
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() << ")");
1027 
1029 
1030  const int height = Height();
1031  const int nnz = J.Capacity();
1032  auto d_I = Read(I, height+1);
1033  auto d_J = Read(J, nnz);
1034  auto d_x = Read(x.GetMemory(), x.Size());
1035  auto d_y = Write(y.GetMemory(), y.Size());
1036  MFEM_FORALL(i, height,
1037  {
1038  bool d_yi = false;
1039  const int end = d_I[i+1];
1040  for (int j = d_I[i]; j < end; j++)
1041  {
1042  if (d_x[d_J[j]])
1043  {
1044  d_yi = true;
1045  break;
1046  }
1047  }
1048  d_y[i] = d_yi;
1049  });
1050 }
1051 
1053  Array<int> &y) const
1054 {
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() << ")");
1058 
1059  y.SetSize(Width());
1060  y = 0;
1061 
1062  for (int i = 0; i < Height(); i++)
1063  {
1064  if (x[i])
1065  {
1066  int end = I[i+1];
1067  for (int j = I[i]; j < end; j++)
1068  {
1069  y[J[j]] = x[i];
1070  }
1071  }
1072  }
1073 }
1074 
1075 void SparseMatrix::AbsMult(const Vector &x, Vector &y) const
1076 {
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 << ")");
1081 
1082  if (Finalized()) { y.UseDevice(true); }
1083  y = 0.0;
1084 
1085  if (!Finalized())
1086  {
1087  const double *xp = x.HostRead();
1088  double *yp = y.HostReadWrite();
1089 
1090  // The matrix is not finalized, but multiplication is still possible
1091  for (int i = 0; i < height; i++)
1092  {
1093  RowNode *row = Rows[i];
1094  double b = 0.0;
1095  for ( ; row != NULL; row = row->Prev)
1096  {
1097  b += std::abs(row->Value) * xp[row->Column];
1098  }
1099  *yp += b;
1100  yp++;
1101  }
1102  return;
1103  }
1104 
1105  const int height = this->height;
1106  const int nnz = J.Capacity();
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();
1111  auto d_y = y.ReadWrite();
1112  MFEM_FORALL(i, height,
1113  {
1114  double d = 0.0;
1115  const int end = d_I[i+1];
1116  for (int j = d_I[i]; j < end; j++)
1117  {
1118  d += std::abs(d_A[j]) * d_x[d_J[j]];
1119  }
1120  d_y[i] += d;
1121  });
1122 }
1123 
1125 {
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 << ")");
1130 
1131  y = 0.0;
1132 
1133  if (!Finalized())
1134  {
1135  double *yp = y.GetData();
1136  // The matrix is not finalized, but multiplication is still possible
1137  for (int i = 0; i < height; i++)
1138  {
1139  RowNode *row = Rows[i];
1140  double b = x(i);
1141  for ( ; row != NULL; row = row->Prev)
1142  {
1143  yp[row->Column] += fabs(row->Value) * b;
1144  }
1145  }
1146  return;
1147  }
1148 
1150  if (At)
1151  {
1152  At->AbsMult(x, y);
1153  }
1154  else
1155  {
1156  for (int i = 0; i < height; i++)
1157  {
1158  const double xi = x[i];
1159  const int end = I[i+1];
1160  for (int j = I[i]; j < end; j++)
1161  {
1162  const int Jj = J[j];
1163  y[Jj] += std::abs(A[j]) * xi;
1164  }
1165  }
1166  }
1167 }
1168 
1169 double SparseMatrix::InnerProduct(const Vector &x, const Vector &y) const
1170 {
1171  MFEM_ASSERT(x.Size() == Width(), "x.Size() = " << x.Size()
1172  << " must be equal to Width() = " << Width());
1173  MFEM_ASSERT(y.Size() == Height(), "y.Size() = " << y.Size()
1174  << " must be equal to Height() = " << Height());
1175 
1176  x.HostRead();
1177  y.HostRead();
1178  if (Finalized())
1179  {
1180  const int nnz = J.Capacity();
1181  HostRead(I, height+1);
1182  HostRead(J, nnz);
1183  HostRead(A, nnz);
1184  }
1185 
1186  double prod = 0.0;
1187  for (int i = 0; i < height; i++)
1188  {
1189  double a = 0.0;
1190  if (A)
1191  {
1192  for (int j = I[i], end = I[i+1]; j < end; j++)
1193  {
1194  a += A[j] * x(J[j]);
1195  }
1196  }
1197  else
1198  {
1199  for (RowNode *np = Rows[i]; np != NULL; np = np->Prev)
1200  {
1201  a += np->Value * x(np->Column);
1202  }
1203  }
1204  prod += a * y(i);
1205  }
1206 
1207  return prod;
1208 }
1209 
1211 {
1212  for (int i = 0; i < height; i++)
1213  {
1214  double a = 0.0;
1215  if (A)
1216  {
1217  for (int j = I[i], end = I[i+1]; j < end; j++)
1218  {
1219  a += A[j];
1220  }
1221  }
1222  else
1223  {
1224  for (RowNode *np = Rows[i]; np != NULL; np = np->Prev)
1225  {
1226  a += np->Value;
1227  }
1228  }
1229  x(i) = a;
1230  }
1231 }
1232 
1233 double SparseMatrix::GetRowNorml1(int irow) const
1234 {
1235  MFEM_VERIFY(irow < height,
1236  "row " << irow << " not in matrix with height " << height);
1237 
1238  double a = 0.0;
1239  if (A)
1240  {
1241  for (int j = I[irow], end = I[irow+1]; j < end; j++)
1242  {
1243  a += fabs(A[j]);
1244  }
1245  }
1246  else
1247  {
1248  for (RowNode *np = Rows[irow]; np != NULL; np = np->Prev)
1249  {
1250  a += fabs(np->Value);
1251  }
1252  }
1253 
1254  return a;
1255 }
1256 
1257 void SparseMatrix::Threshold(double tol, bool fix_empty_rows)
1258 {
1259  MFEM_ASSERT(Finalized(), "Matrix must be finalized.");
1260  double atol;
1261  atol = std::abs(tol);
1262 
1263  fix_empty_rows = height == width ? fix_empty_rows : false;
1264 
1265  double *newA;
1266  int *newI, *newJ;
1267  int i, j, nz;
1268 
1269  newI = Memory<int>(height+1);
1270  newI[0] = 0;
1271  for (i = 0, nz = 0; i < height; i++)
1272  {
1273  bool found = false;
1274  for (j = I[i]; j < I[i+1]; j++)
1275  if (std::abs(A[j]) > atol)
1276  {
1277  found = true;
1278  nz++;
1279  }
1280  if (fix_empty_rows && !found) { nz++; }
1281  newI[i+1] = nz;
1282  }
1283 
1284  newJ = Memory<int>(nz);
1285  newA = Memory<double>(nz);
1286  // Assume we're sorted until we find out otherwise
1287  isSorted = true;
1288  for (i = 0, nz = 0; i < height; i++)
1289  {
1290  bool found = false;
1291  int lastCol = -1;
1292  for (j = I[i]; j < I[i+1]; j++)
1293  if (std::abs(A[j]) > atol)
1294  {
1295  found = true;
1296  newJ[nz] = J[j];
1297  newA[nz] = A[j];
1298  if ( lastCol > newJ[nz] )
1299  {
1300  isSorted = false;
1301  }
1302  lastCol = newJ[nz];
1303  nz++;
1304  }
1305  if (fix_empty_rows && !found)
1306  {
1307  newJ[nz] = i;
1308  newA[nz] = 0.0;
1309  nz++;
1310  }
1311  }
1312  Destroy();
1313  I.Wrap(newI, height+1, true);
1314  J.Wrap(newJ, I[height], true);
1315  A.Wrap(newA, I[height], true);
1316 }
1317 
1318 void SparseMatrix::Finalize(int skip_zeros, bool fix_empty_rows)
1319 {
1320  int i, j, nr, nz;
1321  RowNode *aux;
1322 
1323  if (Finalized())
1324  {
1325  return;
1326  }
1327 
1328  delete [] ColPtrNode;
1329  ColPtrNode = NULL;
1330 
1331  I.New(height+1);
1332  I[0] = 0;
1333  for (i = 1; i <= height; i++)
1334  {
1335  nr = 0;
1336  for (aux = Rows[i-1]; aux != NULL; aux = aux->Prev)
1337  {
1338  if (skip_zeros && aux->Value == 0.0)
1339  {
1340  if (skip_zeros == 2) { continue; }
1341  if ((i-1) != aux->Column) { continue; }
1342 
1343  bool found = false;
1344  double found_val;
1345  for (RowNode *other = Rows[aux->Column]; other != NULL; other = other->Prev)
1346  {
1347  if (other->Column == (i-1))
1348  {
1349  found = true;
1350  found_val = other->Value;
1351  break;
1352  }
1353  }
1354  if (found && found_val == 0.0) { continue; }
1355 
1356  }
1357  nr++;
1358  }
1359  if (fix_empty_rows && !nr) { nr = 1; }
1360  I[i] = I[i-1] + nr;
1361  }
1362 
1363  nz = I[height];
1364  J.New(nz);
1365  A.New(nz);
1366  // Assume we're sorted until we find out otherwise
1367  isSorted = true;
1368  for (j = i = 0; i < height; i++)
1369  {
1370  int lastCol = -1;
1371  nr = 0;
1372  for (aux = Rows[i]; aux != NULL; aux = aux->Prev)
1373  {
1374  if (skip_zeros && aux->Value == 0.0)
1375  {
1376  if (skip_zeros == 2) { continue; }
1377  if (i != aux->Column) { continue; }
1378 
1379  bool found = false;
1380  double found_val;
1381  for (RowNode *other = Rows[aux->Column]; other != NULL; other = other->Prev)
1382  {
1383  if (other->Column == i)
1384  {
1385  found = true;
1386  found_val = other->Value;
1387  break;
1388  }
1389  }
1390  if (found && found_val == 0.0) { continue; }
1391  }
1392 
1393  J[j] = aux->Column;
1394  A[j] = aux->Value;
1395 
1396  if ( lastCol > J[j] )
1397  {
1398  isSorted = false;
1399  }
1400  lastCol = J[j];
1401 
1402  j++;
1403  nr++;
1404  }
1405  if (fix_empty_rows && !nr)
1406  {
1407  J[j] = i;
1408  A[j] = 1.0;
1409  j++;
1410  }
1411  }
1412 
1413 #ifdef MFEM_USE_MEMALLOC
1414  delete NodesMem;
1415  NodesMem = NULL;
1416 #else
1417  for (i = 0; i < height; i++)
1418  {
1419  RowNode *node_p = Rows[i];
1420  while (node_p != NULL)
1421  {
1422  aux = node_p;
1423  node_p = node_p->Prev;
1424  delete aux;
1425  }
1426  }
1427 #endif
1428 
1429  delete [] Rows;
1430  Rows = NULL;
1431 }
1432 
1434 {
1435  int br = blocks.NumRows(), bc = blocks.NumCols();
1436  int nr = (height + br - 1)/br, nc = (width + bc - 1)/bc;
1437 
1438  for (int j = 0; j < bc; j++)
1439  {
1440  for (int i = 0; i < br; i++)
1441  {
1442  int *bI = Memory<int>(nr + 1);
1443  for (int k = 0; k <= nr; k++)
1444  {
1445  bI[k] = 0;
1446  }
1447  blocks(i,j) = new SparseMatrix(bI, NULL, NULL, nr, nc);
1448  }
1449  }
1450 
1451  for (int gr = 0; gr < height; gr++)
1452  {
1453  int bi = gr/nr, i = gr%nr + 1;
1454  if (Finalized())
1455  {
1456  for (int j = I[gr]; j < I[gr+1]; j++)
1457  {
1458  if (A[j] != 0.0)
1459  {
1460  blocks(bi, J[j]/nc)->I[i]++;
1461  }
1462  }
1463  }
1464  else
1465  {
1466  for (RowNode *n_p = Rows[gr]; n_p != NULL; n_p = n_p->Prev)
1467  {
1468  if (n_p->Value != 0.0)
1469  {
1470  blocks(bi, n_p->Column/nc)->I[i]++;
1471  }
1472  }
1473  }
1474  }
1475 
1476  for (int j = 0; j < bc; j++)
1477  {
1478  for (int i = 0; i < br; i++)
1479  {
1480  SparseMatrix &b = *blocks(i,j);
1481  int nnz = 0, rs;
1482  for (int k = 1; k <= nr; k++)
1483  {
1484  rs = b.I[k], b.I[k] = nnz, nnz += rs;
1485  }
1486  b.J.New(nnz);
1487  b.A.New(nnz);
1488  }
1489  }
1490 
1491  for (int gr = 0; gr < height; gr++)
1492  {
1493  int bi = gr/nr, i = gr%nr + 1;
1494  if (Finalized())
1495  {
1496  for (int j = I[gr]; j < I[gr+1]; j++)
1497  {
1498  if (A[j] != 0.0)
1499  {
1500  SparseMatrix &b = *blocks(bi, J[j]/nc);
1501  b.J[b.I[i]] = J[j] % nc;
1502  b.A[b.I[i]] = A[j];
1503  b.I[i]++;
1504  }
1505  }
1506  }
1507  else
1508  {
1509  for (RowNode *n_p = Rows[gr]; n_p != NULL; n_p = n_p->Prev)
1510  {
1511  if (n_p->Value != 0.0)
1512  {
1513  SparseMatrix &b = *blocks(bi, n_p->Column/nc);
1514  b.J[b.I[i]] = n_p->Column % nc;
1515  b.A[b.I[i]] = n_p->Value;
1516  b.I[i]++;
1517  }
1518  }
1519  }
1520  }
1521 }
1522 
1524 {
1525  if (height != width)
1526  {
1527  return infinity();
1528  }
1529 
1530  double symm = 0.0;
1531  if (Empty())
1532  {
1533  // return 0.0;
1534  }
1535  else if (Finalized())
1536  {
1537  for (int i = 1; i < height; i++)
1538  {
1539  for (int j = I[i]; j < I[i+1]; j++)
1540  {
1541  if (J[j] < i)
1542  {
1543  symm = std::max(symm, std::abs(A[j]-(*this)(J[j],i)));
1544  }
1545  }
1546  }
1547  }
1548  else
1549  {
1550  for (int i = 0; i < height; i++)
1551  {
1552  for (RowNode *node_p = Rows[i]; node_p != NULL; node_p = node_p->Prev)
1553  {
1554  int col = node_p->Column;
1555  if (col < i)
1556  {
1557  symm = std::max(symm, std::abs(node_p->Value-(*this)(col,i)));
1558  }
1559  }
1560  }
1561  }
1562  return symm;
1563 }
1564 
1566 {
1567  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
1568 
1569  int i, j;
1570  for (i = 1; i < height; i++)
1571  {
1572  for (j = I[i]; j < I[i+1]; j++)
1573  {
1574  if (J[j] < i)
1575  {
1576  A[j] += (*this)(J[j],i);
1577  A[j] *= 0.5;
1578  (*this)(J[j],i) = A[j];
1579  }
1580  }
1581  }
1582 }
1583 
1585 {
1586  if (Finalized())
1587  {
1588  HostReadI();
1589  return I[height];
1590  }
1591  else
1592  {
1593  int nnz = 0;
1594 
1595  for (int i = 0; i < height; i++)
1596  {
1597  for (RowNode *node_p = Rows[i]; node_p != NULL; node_p = node_p->Prev)
1598  {
1599  nnz++;
1600  }
1601  }
1602 
1603  return nnz;
1604  }
1605 }
1606 
1608 {
1609  double m = 0.0;
1610 
1611  if (A)
1612  {
1613  int nnz = I[height];
1614  for (int j = 0; j < nnz; j++)
1615  {
1616  m = std::max(m, std::abs(A[j]));
1617  }
1618  }
1619  else
1620  {
1621  for (int i = 0; i < height; i++)
1622  {
1623  for (RowNode *n_p = Rows[i]; n_p != NULL; n_p = n_p->Prev)
1624  {
1625  m = std::max(m, std::abs(n_p->Value));
1626  }
1627  }
1628  }
1629  return m;
1630 }
1631 
1632 int SparseMatrix::CountSmallElems(double tol) const
1633 {
1634  int counter = 0;
1635 
1636  if (A)
1637  {
1638  const int nz = I[height];
1639  const double *Ap = A;
1640 
1641  for (int i = 0; i < nz; i++)
1642  {
1643  counter += (std::abs(Ap[i]) <= tol);
1644  }
1645  }
1646  else
1647  {
1648  for (int i = 0; i < height; i++)
1649  {
1650  for (RowNode *aux = Rows[i]; aux != NULL; aux = aux->Prev)
1651  {
1652  counter += (std::abs(aux->Value) <= tol);
1653  }
1654  }
1655  }
1656 
1657  return counter;
1658 }
1659 
1661 {
1662  if (Empty())
1663  {
1664  return 0;
1665  }
1666  else if (Finalized())
1667  {
1668  return mfem::CheckFinite(A, I[height]);
1669  }
1670  else
1671  {
1672  int counter = 0;
1673  for (int i = 0; i < height; i++)
1674  {
1675  for (RowNode *aux = Rows[i]; aux != NULL; aux = aux->Prev)
1676  {
1677  counter += !IsFinite(aux->Value);
1678  }
1679  }
1680  return counter;
1681  }
1682 }
1683 
1685 {
1686  return NULL;
1687 }
1688 
1689 void SparseMatrix::EliminateRow(int row, const double sol, Vector &rhs)
1690 {
1691  RowNode *aux;
1692 
1693  MFEM_ASSERT(row < height && row >= 0,
1694  "Row " << row << " not in matrix of height " << height);
1695 
1696  MFEM_VERIFY(!Finalized(), "Matrix must NOT be finalized.");
1697 
1698  for (aux = Rows[row]; aux != NULL; aux = aux->Prev)
1699  {
1700  rhs(aux->Column) -= sol * aux->Value;
1701  aux->Value = 0.0;
1702  }
1703 }
1704 
1706 {
1707  RowNode *aux;
1708 
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");
1712  MFEM_ASSERT(dpolicy != DIAG_ONE || height == width,
1713  "if dpolicy == DIAG_ONE, matrix must be square, not height = "
1714  << height << ", width = " << width);
1715 
1716  if (Rows == NULL)
1717  {
1718  for (int i=I[row]; i < I[row+1]; ++i)
1719  {
1720  A[i]=0.0;
1721  }
1722  }
1723  else
1724  {
1725  for (aux = Rows[row]; aux != NULL; aux = aux->Prev)
1726  {
1727  aux->Value = 0.0;
1728  }
1729  }
1730 
1731  if (dpolicy == DIAG_ONE)
1732  {
1733  SearchRow(row, row) = 1.;
1734  }
1735 }
1736 
1738 {
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");
1742  MFEM_ASSERT(dpolicy != DIAG_ONE || height == width,
1743  "if dpolicy == DIAG_ONE, matrix must be square, not height = "
1744  << height << ", width = " << width);
1745 
1746  if (Rows == NULL)
1747  {
1748  const int nnz = I[height];
1749  for (int jpos = 0; jpos != nnz; ++jpos)
1750  {
1751  if (J[jpos] == col)
1752  {
1753  A[jpos] = 0.0;
1754  }
1755  }
1756  }
1757  else
1758  {
1759  for (int i = 0; i < height; i++)
1760  {
1761  for (RowNode *aux = Rows[i]; aux != NULL; aux = aux->Prev)
1762  {
1763  if (aux->Column == col)
1764  {
1765  aux->Value = 0.0;
1766  }
1767  }
1768  }
1769  }
1770 
1771  if (dpolicy == DIAG_ONE)
1772  {
1773  SearchRow(col, col) = 1.0;
1774  }
1775 }
1776 
1777 void SparseMatrix::EliminateCols(const Array<int> &cols, const Vector *x,
1778  Vector *b)
1779 {
1780  if (Rows == NULL)
1781  {
1782  for (int i = 0; i < height; i++)
1783  {
1784  for (int jpos = I[i]; jpos != I[i+1]; ++jpos)
1785  {
1786  if (cols[ J[jpos]])
1787  {
1788  if (x && b)
1789  {
1790  (*b)(i) -= A[jpos] * (*x)( J[jpos] );
1791  }
1792  A[jpos] = 0.0;
1793  }
1794  }
1795  }
1796  }
1797  else
1798  {
1799  for (int i = 0; i < height; i++)
1800  {
1801  for (RowNode *aux = Rows[i]; aux != NULL; aux = aux->Prev)
1802  {
1803  if (cols[aux -> Column])
1804  {
1805  if (x && b)
1806  {
1807  (*b)(i) -= aux -> Value * (*x)(aux -> Column);
1808  }
1809  aux->Value = 0.0;
1810  }
1811  }
1812  }
1813  }
1814 }
1815 
1817 {
1818  if (Rows)
1819  {
1820  RowNode *nd;
1821  for (int row = 0; row < height; row++)
1822  {
1823  for (nd = Rows[row]; nd != NULL; nd = nd->Prev)
1824  {
1825  if (col_marker[nd->Column])
1826  {
1827  Ae.Add(row, nd->Column, nd->Value);
1828  nd->Value = 0.0;
1829  }
1830  }
1831  }
1832  }
1833  else
1834  {
1835  for (int row = 0; row < height; row++)
1836  {
1837  for (int j = I[row]; j < I[row+1]; j++)
1838  {
1839  if (col_marker[J[j]])
1840  {
1841  Ae.Add(row, J[j], A[j]);
1842  A[j] = 0.0;
1843  }
1844  }
1845  }
1846  }
1847 }
1848 
1849 
1850 void SparseMatrix::EliminateRowCol(int rc, const double sol, Vector &rhs,
1851  DiagonalPolicy dpolicy)
1852 {
1853  MFEM_ASSERT(rc < height && rc >= 0,
1854  "Row " << rc << " not in matrix of height " << height);
1855  HostReadWriteI();
1856  HostReadWriteJ();
1858 
1859  if (Rows == NULL)
1860  {
1861  for (int j = I[rc]; j < I[rc+1]; j++)
1862  {
1863  const int col = J[j];
1864  if (col == rc)
1865  {
1866  switch (dpolicy)
1867  {
1868  case DIAG_KEEP:
1869  rhs(rc) = A[j] * sol;
1870  break;
1871  case DIAG_ONE:
1872  A[j] = 1.0;
1873  rhs(rc) = sol;
1874  break;
1875  case DIAG_ZERO:
1876  A[j] = 0.;
1877  rhs(rc) = 0.;
1878  break;
1879  default:
1880  mfem_error("SparseMatrix::EliminateRowCol () #2");
1881  break;
1882  }
1883  }
1884  else
1885  {
1886  A[j] = 0.0;
1887  for (int k = I[col]; 1; k++)
1888  {
1889  if (k == I[col+1])
1890  {
1891  mfem_error("SparseMatrix::EliminateRowCol () #3");
1892  }
1893  else if (J[k] == rc)
1894  {
1895  rhs(col) -= sol * A[k];
1896  A[k] = 0.0;
1897  break;
1898  }
1899  }
1900  }
1901  }
1902  }
1903  else
1904  {
1905  for (RowNode *aux = Rows[rc]; aux != NULL; aux = aux->Prev)
1906  {
1907  const int col = aux->Column;
1908  if (col == rc)
1909  {
1910  switch (dpolicy)
1911  {
1912  case DIAG_KEEP:
1913  rhs(rc) = aux->Value * sol;
1914  break;
1915  case DIAG_ONE:
1916  aux->Value = 1.0;
1917  rhs(rc) = sol;
1918  break;
1919  case DIAG_ZERO:
1920  aux->Value = 0.;
1921  rhs(rc) = 0.;
1922  break;
1923  default:
1924  mfem_error("SparseMatrix::EliminateRowCol () #4");
1925  break;
1926  }
1927  }
1928  else
1929  {
1930  aux->Value = 0.0;
1931  for (RowNode *node = Rows[col]; 1; node = node->Prev)
1932  {
1933  if (node == NULL)
1934  {
1935  mfem_error("SparseMatrix::EliminateRowCol () #5");
1936  }
1937  else if (node->Column == rc)
1938  {
1939  rhs(col) -= sol * node->Value;
1940  node->Value = 0.0;
1941  break;
1942  }
1943  }
1944  }
1945  }
1946  }
1947 }
1948 
1950  DenseMatrix &rhs,
1951  DiagonalPolicy dpolicy)
1952 {
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() << ")");
1957 
1958  const int num_rhs = rhs.Width();
1959  if (Rows == NULL)
1960  {
1961  for (int j = I[rc]; j < I[rc+1]; j++)
1962  {
1963  const int col = J[j];
1964  if (col == rc)
1965  {
1966  switch (dpolicy)
1967  {
1968  case DIAG_KEEP:
1969  for (int r = 0; r < num_rhs; r++)
1970  {
1971  rhs(rc,r) = A[j] * sol(r);
1972  }
1973  break;
1974  case DIAG_ONE:
1975  A[j] = 1.0;
1976  for (int r = 0; r < num_rhs; r++)
1977  {
1978  rhs(rc,r) = sol(r);
1979  }
1980  break;
1981  case DIAG_ZERO:
1982  A[j] = 0.;
1983  for (int r = 0; r < num_rhs; r++)
1984  {
1985  rhs(rc,r) = 0.;
1986  }
1987  break;
1988  default:
1989  mfem_error("SparseMatrix::EliminateRowColMultipleRHS() #3");
1990  break;
1991  }
1992  }
1993  else
1994  {
1995  A[j] = 0.0;
1996  for (int k = I[col]; 1; k++)
1997  {
1998  if (k == I[col+1])
1999  {
2000  mfem_error("SparseMatrix::EliminateRowColMultipleRHS() #4");
2001  }
2002  else if (J[k] == rc)
2003  {
2004  for (int r = 0; r < num_rhs; r++)
2005  {
2006  rhs(col,r) -= sol(r) * A[k];
2007  }
2008  A[k] = 0.0;
2009  break;
2010  }
2011  }
2012  }
2013  }
2014  }
2015  else
2016  {
2017  for (RowNode *aux = Rows[rc]; aux != NULL; aux = aux->Prev)
2018  {
2019  const int col = aux->Column;
2020  if (col == rc)
2021  {
2022  switch (dpolicy)
2023  {
2024  case DIAG_KEEP:
2025  for (int r = 0; r < num_rhs; r++)
2026  {
2027  rhs(rc,r) = aux->Value * sol(r);
2028  }
2029  break;
2030  case DIAG_ONE:
2031  aux->Value = 1.0;
2032  for (int r = 0; r < num_rhs; r++)
2033  {
2034  rhs(rc,r) = sol(r);
2035  }
2036  break;
2037  case DIAG_ZERO:
2038  aux->Value = 0.;
2039  for (int r = 0; r < num_rhs; r++)
2040  {
2041  rhs(rc,r) = 0.;
2042  }
2043  break;
2044  default:
2045  mfem_error("SparseMatrix::EliminateRowColMultipleRHS() #5");
2046  break;
2047  }
2048  }
2049  else
2050  {
2051  aux->Value = 0.0;
2052  for (RowNode *node = Rows[col]; 1; node = node->Prev)
2053  {
2054  if (node == NULL)
2055  {
2056  mfem_error("SparseMatrix::EliminateRowColMultipleRHS() #6");
2057  }
2058  else if (node->Column == rc)
2059  {
2060  for (int r = 0; r < num_rhs; r++)
2061  {
2062  rhs(col,r) -= sol(r) * node->Value;
2063  }
2064  node->Value = 0.0;
2065  break;
2066  }
2067  }
2068  }
2069  }
2070  }
2071 }
2072 
2074 {
2075  MFEM_ASSERT(rc < height && rc >= 0,
2076  "Row " << rc << " not in matrix of height " << height);
2077 
2078  if (Rows == NULL)
2079  {
2080  const auto &II = this->I; // only use const access for I
2081  const auto &JJ = this->J; // only use const access for J
2082  for (int j = II[rc]; j < II[rc+1]; j++)
2083  {
2084  const int col = JJ[j];
2085  if (col == rc)
2086  {
2087  if (dpolicy == DIAG_ONE)
2088  {
2089  A[j] = 1.0;
2090  }
2091  else if (dpolicy == DIAG_ZERO)
2092  {
2093  A[j] = 0.0;
2094  }
2095  }
2096  else
2097  {
2098  A[j] = 0.0;
2099  for (int k = II[col]; 1; k++)
2100  {
2101  if (k == II[col+1])
2102  {
2103  mfem_error("SparseMatrix::EliminateRowCol() #2");
2104  }
2105  else if (JJ[k] == rc)
2106  {
2107  A[k] = 0.0;
2108  break;
2109  }
2110  }
2111  }
2112  }
2113  }
2114  else
2115  {
2116  RowNode *aux, *node;
2117 
2118  for (aux = Rows[rc]; aux != NULL; aux = aux->Prev)
2119  {
2120  const int col = aux->Column;
2121  if (col == rc)
2122  {
2123  if (dpolicy == DIAG_ONE)
2124  {
2125  aux->Value = 1.0;
2126  }
2127  else if (dpolicy == DIAG_ZERO)
2128  {
2129  aux->Value = 0.;
2130  }
2131  }
2132  else
2133  {
2134  aux->Value = 0.0;
2135  for (node = Rows[col]; 1; node = node->Prev)
2136  {
2137  if (node == NULL)
2138  {
2139  mfem_error("SparseMatrix::EliminateRowCol() #3");
2140  }
2141  else if (node->Column == rc)
2142  {
2143  node->Value = 0.0;
2144  break;
2145  }
2146  }
2147  }
2148  }
2149  }
2150 }
2151 
2152 // This is almost identical to EliminateRowCol(int, int), except for
2153 // the A[j] = value; and aux->Value = value; lines.
2154 void SparseMatrix::EliminateRowColDiag(int rc, double value)
2155 {
2156  MFEM_ASSERT(rc < height && rc >= 0,
2157  "Row " << rc << " not in matrix of height " << height);
2158 
2159  if (Rows == NULL)
2160  {
2161  for (int j = I[rc]; j < I[rc+1]; j++)
2162  {
2163  const int col = J[j];
2164  if (col == rc)
2165  {
2166  A[j] = value;
2167  }
2168  else
2169  {
2170  A[j] = 0.0;
2171  for (int k = I[col]; 1; k++)
2172  {
2173  if (k == I[col+1])
2174  {
2175  mfem_error("SparseMatrix::EliminateRowCol() #2");
2176  }
2177  else if (J[k] == rc)
2178  {
2179  A[k] = 0.0;
2180  break;
2181  }
2182  }
2183  }
2184  }
2185  }
2186  else
2187  {
2188  RowNode *aux, *node;
2189 
2190  for (aux = Rows[rc]; aux != NULL; aux = aux->Prev)
2191  {
2192  const int col = aux->Column;
2193  if (col == rc)
2194  {
2195  aux->Value = value;
2196  }
2197  else
2198  {
2199  aux->Value = 0.0;
2200  for (node = Rows[col]; 1; node = node->Prev)
2201  {
2202  if (node == NULL)
2203  {
2204  mfem_error("SparseMatrix::EliminateRowCol() #3");
2205  }
2206  else if (node->Column == rc)
2207  {
2208  node->Value = 0.0;
2209  break;
2210  }
2211  }
2212  }
2213  }
2214  }
2215 }
2216 
2218  DiagonalPolicy dpolicy)
2219 {
2220  if (Rows)
2221  {
2222  RowNode *nd, *nd2;
2223  for (nd = Rows[rc]; nd != NULL; nd = nd->Prev)
2224  {
2225  const int col = nd->Column;
2226  if (col == rc)
2227  {
2228  switch (dpolicy)
2229  {
2230  case DIAG_ONE:
2231  Ae.Add(rc, rc, nd->Value - 1.0);
2232  nd->Value = 1.0;
2233  break;
2234  case DIAG_ZERO:
2235  Ae.Add(rc, rc, nd->Value);
2236  nd->Value = 0.;
2237  break;
2238  case DIAG_KEEP:
2239  break;
2240  default:
2241  mfem_error("SparseMatrix::EliminateRowCol #1");
2242  break;
2243  }
2244  }
2245  else
2246  {
2247  Ae.Add(rc, col, nd->Value);
2248  nd->Value = 0.0;
2249  for (nd2 = Rows[col]; 1; nd2 = nd2->Prev)
2250  {
2251  if (nd2 == NULL)
2252  {
2253  mfem_error("SparseMatrix::EliminateRowCol #2");
2254  }
2255  else if (nd2->Column == rc)
2256  {
2257  Ae.Add(col, rc, nd2->Value);
2258  nd2->Value = 0.0;
2259  break;
2260  }
2261  }
2262  }
2263  }
2264  }
2265  else
2266  {
2267  for (int j = I[rc]; j < I[rc+1]; j++)
2268  {
2269  const int col = J[j];
2270  if (col == rc)
2271  {
2272  switch (dpolicy)
2273  {
2274  case DIAG_ONE:
2275  Ae.Add(rc, rc, A[j] - 1.0);
2276  A[j] = 1.0;
2277  break;
2278  case DIAG_ZERO:
2279  Ae.Add(rc, rc, A[j]);
2280  A[j] = 0.;
2281  break;
2282  case DIAG_KEEP:
2283  break;
2284  default:
2285  mfem_error("SparseMatrix::EliminateRowCol #3");
2286  break;
2287  }
2288  }
2289  else
2290  {
2291  Ae.Add(rc, col, A[j]);
2292  A[j] = 0.0;
2293  for (int k = I[col]; true; k++)
2294  {
2295  if (k == I[col+1])
2296  {
2297  mfem_error("SparseMatrix::EliminateRowCol #4");
2298  }
2299  else if (J[k] == rc)
2300  {
2301  Ae.Add(col, rc, A[k]);
2302  A[k] = 0.0;
2303  break;
2304  }
2305  }
2306  }
2307  }
2308  }
2309 }
2310 
2312  DiagonalPolicy diag_policy)
2313 {
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();
2318  auto dA = ReadWriteData();
2319 
2320  MFEM_FORALL(i, n_ess_dofs,
2321  {
2322  const int idof = ess_dofs_d[i];
2323  for (int j=dI[idof]; j<dI[idof+1]; ++j)
2324  {
2325  const int jdof = dJ[j];
2326  if (jdof != idof)
2327  {
2328  dA[j] = 0.0;
2329  for (int k=dI[jdof]; k<dI[jdof+1]; ++k)
2330  {
2331  if (dJ[k] == idof)
2332  {
2333  dA[k] = 0.0;
2334  break;
2335  }
2336  }
2337  }
2338  else
2339  {
2340  if (diag_policy == DiagonalPolicy::DIAG_ONE)
2341  {
2342  dA[j] = 1.0;
2343  }
2344  else if (diag_policy == DiagonalPolicy::DIAG_ZERO)
2345  {
2346  dA[j] = 0.0;
2347  }
2348  // else (diag_policy == DiagonalPolicy::DIAG_KEEP)
2349  }
2350  }
2351  });
2352 }
2353 
2355 {
2356  for (int i = 0; i < height; i++)
2357  {
2358  if (I[i+1] == I[i]+1 && fabs(A[I[i]]) < 1e-16)
2359  {
2360  A[I[i]] = 1.0;
2361  }
2362  }
2363 }
2364 
2365 void SparseMatrix::EliminateZeroRows(const double threshold)
2366 {
2367  for (int i = 0; i < height; i++)
2368  {
2369  double zero = 0.0;
2370  for (int j = I[i]; j < I[i+1]; j++)
2371  {
2372  zero += fabs(A[j]);
2373  }
2374  if (zero <= threshold)
2375  {
2376  for (int j = I[i]; j < I[i+1]; j++)
2377  {
2378  A[j] = (J[j] == i) ? 1.0 : 0.0;
2379  }
2380  }
2381  }
2382 }
2383 
2385 {
2386  if (!Finalized())
2387  {
2388  double *yp = y.GetData();
2389  const double *xp = x.GetData();
2390  RowNode *diag_p, *n_p, **R = Rows;
2391 
2392  const int s = height;
2393  for (int i = 0; i < s; i++)
2394  {
2395  double sum = 0.0;
2396  diag_p = NULL;
2397  for (n_p = R[i]; n_p != NULL; n_p = n_p->Prev)
2398  {
2399  const int c = n_p->Column;
2400  if (c == i)
2401  {
2402  diag_p = n_p;
2403  }
2404  else
2405  {
2406  sum += n_p->Value * yp[c];
2407  }
2408  }
2409 
2410  if (diag_p != NULL && diag_p->Value != 0.0)
2411  {
2412  yp[i] = (xp[i] - sum) / diag_p->Value;
2413  }
2414  else if (xp[i] == sum)
2415  {
2416  yp[i] = sum;
2417  }
2418  else
2419  {
2420  mfem_error("SparseMatrix::Gauss_Seidel_forw()");
2421  }
2422  }
2423  }
2424  else
2425  {
2426  const int s = height;
2427  const int nnz = J.Capacity();
2428  const int *Ip = HostRead(I, s+1);
2429  const int *Jp = HostRead(J, nnz);
2430  const double *Ap = HostRead(A, nnz);
2431  double *yp = y.HostReadWrite();
2432  const double *xp = x.HostRead();
2433 
2434  for (int i = 0, j = Ip[0]; i < s; i++)
2435  {
2436  const int end = Ip[i+1];
2437  double sum = 0.0;
2438  int d = -1;
2439  for ( ; j < end; j++)
2440  {
2441  const int c = Jp[j];
2442  if (c == i)
2443  {
2444  d = j;
2445  }
2446  else
2447  {
2448  sum += Ap[j] * yp[c];
2449  }
2450  }
2451 
2452  if (d >= 0 && Ap[d] != 0.0)
2453  {
2454  yp[i] = (xp[i] - sum) / Ap[d];
2455  }
2456  else if (xp[i] == sum)
2457  {
2458  yp[i] = sum;
2459  }
2460  else
2461  {
2462  mfem_error("SparseMatrix::Gauss_Seidel_forw(...) #2");
2463  }
2464  }
2465  }
2466 }
2467 
2469 {
2470  if (!Finalized())
2471  {
2472  double *yp = y.GetData();
2473  const double *xp = x.GetData();
2474  RowNode *diag_p, *n_p, **R = Rows;
2475 
2476  for (int i = height-1; i >= 0; i--)
2477  {
2478  double sum = 0.;
2479  diag_p = NULL;
2480  for (n_p = R[i]; n_p != NULL; n_p = n_p->Prev)
2481  {
2482  const int c = n_p->Column;
2483  if (c == i)
2484  {
2485  diag_p = n_p;
2486  }
2487  else
2488  {
2489  sum += n_p->Value * yp[c];
2490  }
2491  }
2492 
2493  if (diag_p != NULL && diag_p->Value != 0.0)
2494  {
2495  yp[i] = (xp[i] - sum) / diag_p->Value;
2496  }
2497  else if (xp[i] == sum)
2498  {
2499  yp[i] = sum;
2500  }
2501  else
2502  {
2503  mfem_error("SparseMatrix::Gauss_Seidel_back()");
2504  }
2505  }
2506  }
2507  else
2508  {
2509  const int s = height;
2510  const int nnz = J.Capacity();
2511  const int *Ip = HostRead(I, s+1);
2512  const int *Jp = HostRead(J, nnz);
2513  const double *Ap = HostRead(A, nnz);
2514  double *yp = y.HostReadWrite();
2515  const double *xp = x.HostRead();
2516 
2517  for (int i = s-1, j = Ip[s]-1; i >= 0; i--)
2518  {
2519  const int beg = Ip[i];
2520  double sum = 0.;
2521  int d = -1;
2522  for ( ; j >= beg; j--)
2523  {
2524  const int c = Jp[j];
2525  if (c == i)
2526  {
2527  d = j;
2528  }
2529  else
2530  {
2531  sum += Ap[j] * yp[c];
2532  }
2533  }
2534 
2535  if (d >= 0 && Ap[d] != 0.0)
2536  {
2537  yp[i] = (xp[i] - sum) / Ap[d];
2538  }
2539  else if (xp[i] == sum)
2540  {
2541  yp[i] = sum;
2542  }
2543  else
2544  {
2545  mfem_error("SparseMatrix::Gauss_Seidel_back(...) #2");
2546  }
2547  }
2548  }
2549 }
2550 
2552 {
2553  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
2554 
2555  double sc = 1.0;
2556  for (int i = 0; i < height; i++)
2557  {
2558  int d = -1;
2559  double norm = 0.0;
2560  for (int j = I[i]; j < I[i+1]; j++)
2561  {
2562  if (J[j] == i)
2563  {
2564  d = j;
2565  }
2566  norm += fabs(A[j]);
2567  }
2568  if (d >= 0 && A[d] != 0.0)
2569  {
2570  double a = 1.8 * fabs(A[d]) / norm;
2571  if (a < sc)
2572  {
2573  sc = a;
2574  }
2575  }
2576  else
2577  {
2578  mfem_error("SparseMatrix::GetJacobiScaling() #2");
2579  }
2580  }
2581  return sc;
2582 }
2583 
2584 void SparseMatrix::Jacobi(const Vector &b, const Vector &x0, Vector &x1,
2585  double sc, bool use_abs_diag) const
2586 {
2587  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
2588 
2589  for (int i = 0; i < height; i++)
2590  {
2591  int d = -1;
2592  double sum = b(i);
2593  for (int j = I[i]; j < I[i+1]; j++)
2594  {
2595  if (J[j] == i)
2596  {
2597  d = j;
2598  }
2599  else
2600  {
2601  sum -= A[j] * x0(J[j]);
2602  }
2603  }
2604  if (d >= 0 && A[d] != 0.0)
2605  {
2606  const double diag = (use_abs_diag) ? fabs(A[d]) : A[d];
2607  x1(i) = sc * (sum / diag) + (1.0 - sc) * x0(i);
2608  }
2609  else
2610  {
2611  mfem_error("SparseMatrix::Jacobi(...) #2");
2612  }
2613  }
2614 }
2615 
2617  double sc, bool use_abs_diag) const
2618 {
2619  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
2620 
2621  const int H = height;
2622  const int nnz = J.Capacity();
2623  const bool use_dev = b.UseDevice() || x.UseDevice();
2624 
2625  const auto Ap = Read(A, nnz, use_dev);
2626  const auto Ip = Read(I, height+1, use_dev);
2627  const auto Jp = Read(J, nnz, use_dev);
2628 
2629  const auto bp = b.Read(use_dev);
2630  auto xp = x.Write(use_dev);
2631 
2632  MFEM_FORALL_SWITCH(use_dev, i, H,
2633  {
2634  const int end = Ip[i+1];
2635  for (int j = Ip[i]; true; j++)
2636  {
2637  if (j == end)
2638  {
2639  MFEM_ABORT_KERNEL("Diagonal not found in SparseMatrix::DiagScale");
2640  }
2641  if (Jp[j] == i)
2642  {
2643  const double diag = (use_abs_diag) ? fabs(Ap[j]) : Ap[j];
2644  if (diag == 0.0)
2645  {
2646  MFEM_ABORT_KERNEL("Zero diagonal in SparseMatrix::DiagScale");
2647  }
2648  xp[i] = sc * bp[i] / diag;
2649  break;
2650  }
2651  }
2652  });
2653 }
2654 
2655 template <bool useFabs>
2656 static void JacobiDispatch(const Vector &b, const Vector &x0, Vector &x1,
2657  const Memory<int> &I, const Memory<int> &J,
2658  const Memory<double> &A, const int height,
2659  const double sc)
2660 {
2661  const bool useDevice = b.UseDevice() || x0.UseDevice() || x1.UseDevice();
2662 
2663  const auto bp = b.Read(useDevice);
2664  const auto x0p = x0.Read(useDevice);
2665  auto x1p = x1.Write(useDevice);
2666 
2667  const auto Ip = Read(I, height+1, useDevice);
2668  const auto Jp = Read(J, J.Capacity(), useDevice);
2669  const auto Ap = Read(A, J.Capacity(), useDevice);
2670 
2671  MFEM_FORALL_SWITCH(useDevice, i, height,
2672  {
2673  double resi = bp[i], norm = 0.0;
2674  for (int j = Ip[i]; j < Ip[i+1]; j++)
2675  {
2676  resi -= Ap[j] * x0p[Jp[j]];
2677  if (useFabs)
2678  {
2679  norm += fabs(Ap[j]);
2680  }
2681  else
2682  {
2683  norm += Ap[j];
2684  }
2685  }
2686  if (norm > 0.0)
2687  {
2688  x1p[i] = x0p[i] + sc * resi / norm;
2689  }
2690  else
2691  {
2692  if (useFabs)
2693  {
2694  MFEM_ABORT_KERNEL("L1 norm of row is zero.");
2695  }
2696  else
2697  {
2698  MFEM_ABORT_KERNEL("sum of row is zero.");
2699  }
2700  }
2701  });
2702 }
2703 
2704 void SparseMatrix::Jacobi2(const Vector &b, const Vector &x0, Vector &x1,
2705  double sc) const
2706 {
2707  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
2708  JacobiDispatch<true>(b,x0,x1,I,J,A,height,sc);
2709 }
2710 
2711 void SparseMatrix::Jacobi3(const Vector &b, const Vector &x0, Vector &x1,
2712  double sc) const
2713 {
2714  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
2715  JacobiDispatch<false>(b,x0,x1,I,J,A,height,sc);
2716 }
2717 
2718 void SparseMatrix::AddSubMatrix(const Array<int> &rows, const Array<int> &cols,
2719  const DenseMatrix &subm, int skip_zeros)
2720 {
2721  int i, j, gi, gj, s, t;
2722  double a;
2723 
2724  if (Finalized())
2725  {
2726  HostReadI();
2727  HostReadJ();
2729  }
2730 
2731  for (i = 0; i < rows.Size(); i++)
2732  {
2733  if ((gi=rows[i]) < 0) { gi = -1-gi, s = -1; }
2734  else { s = 1; }
2735  MFEM_ASSERT(gi < height,
2736  "Trying to insert a row " << gi << " outside the matrix height "
2737  << height);
2738  SetColPtr(gi);
2739  for (j = 0; j < cols.Size(); j++)
2740  {
2741  if ((gj=cols[j]) < 0) { gj = -1-gj, t = -s; }
2742  else { t = s; }
2743  MFEM_ASSERT(gj < width,
2744  "Trying to insert a column " << gj << " outside the matrix width "
2745  << width);
2746  a = subm(i, j);
2747  if (skip_zeros && a == 0.0)
2748  {
2749  // Skip assembly of zero elements if either:
2750  // (i) user specified to skip zeros regardless of symmetry, or
2751  // (ii) symmetry is not broken.
2752  if (skip_zeros == 2 || &rows != &cols || subm(j, i) == 0.0)
2753  {
2754  continue;
2755  }
2756  }
2757  if (t < 0) { a = -a; }
2758  _Add_(gj, a);
2759  }
2760  ClearColPtr();
2761  }
2762 }
2763 
2764 void SparseMatrix::Set(const int i, const int j, const double val)
2765 {
2766  double a = val;
2767  int gi, gj, s, t;
2768 
2769  if ((gi=i) < 0) { gi = -1-gi, s = -1; }
2770  else { s = 1; }
2771  MFEM_ASSERT(gi < height,
2772  "Trying to set a row " << gi << " outside the matrix height "
2773  << height);
2774  if ((gj=j) < 0) { gj = -1-gj, t = -s; }
2775  else { t = s; }
2776  MFEM_ASSERT(gj < width,
2777  "Trying to set a column " << gj << " outside the matrix width "
2778  << width);
2779  if (t < 0) { a = -a; }
2780  _Set_(gi, gj, a);
2781 }
2782 
2783 void SparseMatrix::Add(const int i, const int j, const double val)
2784 {
2785  int gi, gj, s, t;
2786  double a = val;
2787 
2788  if ((gi=i) < 0) { gi = -1-gi, s = -1; }
2789  else { s = 1; }
2790  MFEM_ASSERT(gi < height,
2791  "Trying to insert a row " << gi << " outside the matrix height "
2792  << height);
2793  if ((gj=j) < 0) { gj = -1-gj, t = -s; }
2794  else { t = s; }
2795  MFEM_ASSERT(gj < width,
2796  "Trying to insert a column " << gj << " outside the matrix width "
2797  << width);
2798  if (t < 0) { a = -a; }
2799  _Add_(gi, gj, a);
2800 }
2801 
2802 void SparseMatrix::SetSubMatrix(const Array<int> &rows, const Array<int> &cols,
2803  const DenseMatrix &subm, int skip_zeros)
2804 {
2805  int i, j, gi, gj, s, t;
2806  double a;
2807 
2808  for (i = 0; i < rows.Size(); i++)
2809  {
2810  if ((gi=rows[i]) < 0) { gi = -1-gi, s = -1; }
2811  else { s = 1; }
2812  MFEM_ASSERT(gi < height,
2813  "Trying to set a row " << gi << " outside the matrix height "
2814  << height);
2815  SetColPtr(gi);
2816  for (j = 0; j < cols.Size(); j++)
2817  {
2818  a = subm(i, j);
2819  if (skip_zeros && a == 0.0)
2820  {
2821  // Skip assembly of zero elements if either:
2822  // (i) user specified to skip zeros regardless of symmetry, or
2823  // (ii) symmetry is not broken.
2824  if (skip_zeros == 2 || &rows != &cols || subm(j, i) == 0.0)
2825  {
2826  continue;
2827  }
2828  }
2829  if ((gj=cols[j]) < 0) { gj = -1-gj, t = -s; }
2830  else { t = s; }
2831  MFEM_ASSERT(gj < width,
2832  "Trying to set a column " << gj << " outside the matrix width "
2833  << width);
2834  if (t < 0) { a = -a; }
2835  _Set_(gj, a);
2836  }
2837  ClearColPtr();
2838  }
2839 }
2840 
2842  const Array<int> &cols,
2843  const DenseMatrix &subm,
2844  int skip_zeros)
2845 {
2846  int i, j, gi, gj, s, t;
2847  double a;
2848 
2849  for (i = 0; i < rows.Size(); i++)
2850  {
2851  if ((gi=rows[i]) < 0) { gi = -1-gi, s = -1; }
2852  else { s = 1; }
2853  MFEM_ASSERT(gi < height,
2854  "Trying to set a row " << gi << " outside the matrix height "
2855  << height);
2856  SetColPtr(gi);
2857  for (j = 0; j < cols.Size(); j++)
2858  {
2859  a = subm(j, i);
2860  if (skip_zeros && a == 0.0)
2861  {
2862  // Skip assembly of zero elements if either:
2863  // (i) user specified to skip zeros regardless of symmetry, or
2864  // (ii) symmetry is not broken.
2865  if (skip_zeros == 2 || &rows != &cols || subm(j, i) == 0.0)
2866  {
2867  continue;
2868  }
2869  }
2870  if ((gj=cols[j]) < 0) { gj = -1-gj, t = -s; }
2871  else { t = s; }
2872  MFEM_ASSERT(gj < width,
2873  "Trying to set a column " << gj << " outside the matrix width "
2874  << width);
2875  if (t < 0) { a = -a; }
2876  _Set_(gj, a);
2877  }
2878  ClearColPtr();
2879  }
2880 }
2881 
2882 void SparseMatrix::GetSubMatrix(const Array<int> &rows, const Array<int> &cols,
2883  DenseMatrix &subm) const
2884 {
2885  int i, j, gi, gj, s, t;
2886  double a;
2887 
2888  for (i = 0; i < rows.Size(); i++)
2889  {
2890  if ((gi=rows[i]) < 0) { gi = -1-gi, s = -1; }
2891  else { s = 1; }
2892  MFEM_ASSERT(gi < height,
2893  "Trying to read a row " << gi << " outside the matrix height "
2894  << height);
2895  SetColPtr(gi);
2896  for (j = 0; j < cols.Size(); j++)
2897  {
2898  if ((gj=cols[j]) < 0) { gj = -1-gj, t = -s; }
2899  else { t = s; }
2900  MFEM_ASSERT(gj < width,
2901  "Trying to read a column " << gj << " outside the matrix width "
2902  << width);
2903  a = _Get_(gj);
2904  subm(i, j) = (t < 0) ? (-a) : (a);
2905  }
2906  ClearColPtr();
2907  }
2908 }
2909 
2910 bool SparseMatrix::RowIsEmpty(const int row) const
2911 {
2912  int gi;
2913 
2914  if ((gi=row) < 0)
2915  {
2916  gi = -1-gi;
2917  }
2918  MFEM_ASSERT(gi < height,
2919  "Trying to query a row " << gi << " outside the matrix height "
2920  << height);
2921  if (Rows)
2922  {
2923  return (Rows[gi] == NULL);
2924  }
2925  else
2926  {
2927  return (I[gi] == I[gi+1]);
2928  }
2929 }
2930 
2931 int SparseMatrix::GetRow(const int row, Array<int> &cols, Vector &srow) const
2932 {
2933  RowNode *n;
2934  int j, gi;
2935 
2936  if ((gi=row) < 0) { gi = -1-gi; }
2937  MFEM_ASSERT(gi < height,
2938  "Trying to read a row " << gi << " outside the matrix height "
2939  << height);
2940  if (Rows)
2941  {
2942  for (n = Rows[gi], j = 0; n; n = n->Prev)
2943  {
2944  j++;
2945  }
2946  cols.SetSize(j);
2947  srow.SetSize(j);
2948  for (n = Rows[gi], j = 0; n; n = n->Prev, j++)
2949  {
2950  cols[j] = n->Column;
2951  srow(j) = n->Value;
2952  }
2953  if (row < 0)
2954  {
2955  srow.Neg();
2956  }
2957 
2958  return 0;
2959  }
2960  else
2961  {
2962  j = I[gi];
2963  cols.MakeRef(const_cast<int*>((const int*)J) + j, I[gi+1]-j);
2964  srow.NewDataAndSize(
2965  const_cast<double*>((const double*)A) + j, cols.Size());
2966  MFEM_ASSERT(row >= 0, "Row not valid: " << row << ", height: " << height);
2967  return 1;
2968  }
2969 }
2970 
2971 void SparseMatrix::SetRow(const int row, const Array<int> &cols,
2972  const Vector &srow)
2973 {
2974  int gi, gj, s, t;
2975  double a;
2976 
2977  if ((gi=row) < 0) { gi = -1-gi, s = -1; }
2978  else { s = 1; }
2979  MFEM_ASSERT(gi < height,
2980  "Trying to set a row " << gi << " outside the matrix height "
2981  << height);
2982 
2983  if (!Finalized())
2984  {
2985  SetColPtr(gi);
2986  for (int j = 0; j < cols.Size(); j++)
2987  {
2988  if ((gj=cols[j]) < 0) { gj = -1-gj, t = -s; }
2989  else { t = s; }
2990  MFEM_ASSERT(gj < width,
2991  "Trying to set a column " << gj << " outside the matrix"
2992  " width " << width);
2993  a = srow(j);
2994  if (t < 0) { a = -a; }
2995  _Set_(gj, a);
2996  }
2997  ClearColPtr();
2998  }
2999  else
3000  {
3001  MFEM_ASSERT(cols.Size() == RowSize(gi), "");
3002  MFEM_ASSERT(cols.Size() == srow.Size(), "");
3003 
3004  for (int i = I[gi], j = 0; j < cols.Size(); j++, i++)
3005  {
3006  if ((gj=cols[j]) < 0) { gj = -1-gj, t = -s; }
3007  else { t = s; }
3008  MFEM_ASSERT(gj < width,
3009  "Trying to set a column " << gj << " outside the matrix"
3010  " width " << width);
3011 
3012  J[i] = gj;
3013  A[i] = srow[j] * t;
3014  }
3015  }
3016 }
3017 
3018 void SparseMatrix::AddRow(const int row, const Array<int> &cols,
3019  const Vector &srow)
3020 {
3021  int j, gi, gj, s, t;
3022  double a;
3023 
3024  MFEM_VERIFY(!Finalized(), "Matrix must NOT be finalized.");
3025 
3026  if ((gi=row) < 0) { gi = -1-gi, s = -1; }
3027  else { s = 1; }
3028  MFEM_ASSERT(gi < height,
3029  "Trying to insert a row " << gi << " outside the matrix height "
3030  << height);
3031  SetColPtr(gi);
3032  for (j = 0; j < cols.Size(); j++)
3033  {
3034  if ((gj=cols[j]) < 0) { gj = -1-gj, t = -s; }
3035  else { t = s; }
3036  MFEM_ASSERT(gj < width,
3037  "Trying to insert a column " << gj << " outside the matrix width "
3038  << width);
3039  a = srow(j);
3040  if (a == 0.0)
3041  {
3042  continue;
3043  }
3044  if (t < 0) { a = -a; }
3045  _Add_(gj, a);
3046  }
3047  ClearColPtr();
3048 }
3049 
3050 void SparseMatrix::ScaleRow(const int row, const double scale)
3051 {
3052  int i;
3053 
3054  if ((i=row) < 0)
3055  {
3056  i = -1-i;
3057  }
3058  if (Rows != NULL)
3059  {
3060  RowNode *aux;
3061 
3062  for (aux = Rows[i]; aux != NULL; aux = aux -> Prev)
3063  {
3064  aux -> Value *= scale;
3065  }
3066  }
3067  else
3068  {
3069  int j, end = I[i+1];
3070 
3071  for (j = I[i]; j < end; j++)
3072  {
3073  A[j] *= scale;
3074  }
3075  }
3076 }
3077 
3079 {
3080  double scale;
3081  if (Rows != NULL)
3082  {
3083  RowNode *aux;
3084  for (int i=0; i < height; ++i)
3085  {
3086  scale = sl(i);
3087  for (aux = Rows[i]; aux != NULL; aux = aux -> Prev)
3088  {
3089  aux -> Value *= scale;
3090  }
3091  }
3092  }
3093  else
3094  {
3095  int j, end;
3096 
3097  for (int i=0; i < height; ++i)
3098  {
3099  end = I[i+1];
3100  scale = sl(i);
3101  for (j = I[i]; j < end; j++)
3102  {
3103  A[j] *= scale;
3104  }
3105  }
3106  }
3107 }
3108 
3110 {
3111  if (Rows != NULL)
3112  {
3113  RowNode *aux;
3114  for (int i=0; i < height; ++i)
3115  {
3116  for (aux = Rows[i]; aux != NULL; aux = aux -> Prev)
3117  {
3118  aux -> Value *= sr(aux->Column);
3119  }
3120  }
3121  }
3122  else
3123  {
3124  int j, end;
3125 
3126  for (int i=0; i < height; ++i)
3127  {
3128  end = I[i+1];
3129  for (j = I[i]; j < end; j++)
3130  {
3131  A[j] *= sr(J[j]);
3132  }
3133  }
3134  }
3135 }
3136 
3138 {
3139  MFEM_ASSERT(height == B.height && width == B.width,
3140  "Mismatch of this matrix size and rhs. This height = "
3141  << height << ", width = " << width << ", B.height = "
3142  << B.height << ", B.width = " << B.width);
3143 
3144  for (int i = 0; i < height; i++)
3145  {
3146  SetColPtr(i);
3147  if (B.Rows)
3148  {
3149  for (RowNode *aux = B.Rows[i]; aux != NULL; aux = aux->Prev)
3150  {
3151  _Add_(aux->Column, aux->Value);
3152  }
3153  }
3154  else
3155  {
3156  for (int j = B.I[i]; j < B.I[i+1]; j++)
3157  {
3158  _Add_(B.J[j], B.A[j]);
3159  }
3160  }
3161  ClearColPtr();
3162  }
3163 
3164  return (*this);
3165 }
3166 
3167 void SparseMatrix::Add(const double a, const SparseMatrix &B)
3168 {
3169  for (int i = 0; i < height; i++)
3170  {
3171  B.SetColPtr(i);
3172  if (Rows)
3173  {
3174  for (RowNode *np = Rows[i]; np != NULL; np = np->Prev)
3175  {
3176  np->Value += a * B._Get_(np->Column);
3177  }
3178  }
3179  else
3180  {
3181  for (int j = I[i]; j < I[i+1]; j++)
3182  {
3183  A[j] += a * B._Get_(J[j]);
3184  }
3185  }
3186  B.ClearColPtr();
3187  }
3188 }
3189 
3191 {
3192  if (Rows == NULL)
3193  {
3194  const int nnz = J.Capacity();
3195  double *h_A = HostWrite(A, nnz);
3196  for (int i = 0; i < nnz; i++)
3197  {
3198  h_A[i] = a;
3199  }
3200  }
3201  else
3202  {
3203  for (int i = 0; i < height; i++)
3204  {
3205  for (RowNode *node_p = Rows[i]; node_p != NULL;
3206  node_p = node_p -> Prev)
3207  {
3208  node_p -> Value = a;
3209  }
3210  }
3211  }
3212 
3213  return (*this);
3214 }
3215 
3217 {
3218  if (Rows == NULL)
3219  {
3220  for (int i = 0, nnz = I[height]; i < nnz; i++)
3221  {
3222  A[i] *= a;
3223  }
3224  }
3225  else
3226  {
3227  for (int i = 0; i < height; i++)
3228  {
3229  for (RowNode *node_p = Rows[i]; node_p != NULL;
3230  node_p = node_p -> Prev)
3231  {
3232  node_p -> Value *= a;
3233  }
3234  }
3235  }
3236 
3237  return (*this);
3238 }
3239 
3240 void SparseMatrix::Print(std::ostream & os, int width_) const
3241 {
3242  int i, j;
3243 
3244  if (A == NULL)
3245  {
3246  RowNode *nd;
3247  for (i = 0; i < height; i++)
3248  {
3249  os << "[row " << i << "]\n";
3250  for (nd = Rows[i], j = 0; nd != NULL; nd = nd->Prev, j++)
3251  {
3252  os << " (" << nd->Column << "," << nd->Value << ")";
3253  if ( !((j+1) % width_) )
3254  {
3255  os << '\n';
3256  }
3257  }
3258  if (j % width_)
3259  {
3260  os << '\n';
3261  }
3262  }
3263  return;
3264  }
3265 
3266  // HostRead forces synchronization
3267  HostReadI();
3268  HostReadJ();
3269  HostReadData();
3270  for (i = 0; i < height; i++)
3271  {
3272  os << "[row " << i << "]\n";
3273  for (j = I[i]; j < I[i+1]; j++)
3274  {
3275  os << " (" << J[j] << "," << A[j] << ")";
3276  if ( !((j+1-I[i]) % width_) )
3277  {
3278  os << '\n';
3279  }
3280  }
3281  if ((j-I[i]) % width_)
3282  {
3283  os << '\n';
3284  }
3285  }
3286 }
3287 
3288 void SparseMatrix::PrintMatlab(std::ostream & os) const
3289 {
3290  os << "% size " << height << " " << width << "\n";
3291  os << "% Non Zeros " << NumNonZeroElems() << "\n";
3292  int i, j;
3293  ios::fmtflags old_fmt = os.flags();
3294  os.setf(ios::scientific);
3295  std::streamsize old_prec = os.precision(14);
3296 
3297  for (i = 0; i < height; i++)
3298  {
3299  for (j = I[i]; j < I[i+1]; j++)
3300  {
3301  os << i+1 << " " << J[j]+1 << " " << A[j] << '\n';
3302  }
3303  }
3304  // Write a zero entry at (m,n) to make sure MATLAB doesn't shrink the matrix
3305  os << height << " " << width << " 0.0\n";
3306  os.precision(old_prec);
3307  os.flags(old_fmt);
3308 }
3309 
3310 void SparseMatrix::PrintMM(std::ostream & os) const
3311 {
3312  int i, j;
3313  ios::fmtflags old_fmt = os.flags();
3314  os.setf(ios::scientific);
3315  std::streamsize old_prec = os.precision(14);
3316 
3317  os << "%%MatrixMarket matrix coordinate real general" << '\n'
3318  << "% Generated by MFEM" << '\n';
3319 
3320  os << height << " " << width << " " << NumNonZeroElems() << '\n';
3321  for (i = 0; i < height; i++)
3322  {
3323  for (j = I[i]; j < I[i+1]; j++)
3324  {
3325  os << i+1 << " " << J[j]+1 << " " << A[j] << '\n';
3326  }
3327  }
3328  os.precision(old_prec);
3329  os.flags(old_fmt);
3330 }
3331 
3332 void SparseMatrix::PrintCSR(std::ostream & os) const
3333 {
3334  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
3335 
3336  int i;
3337 
3338  os << height << '\n'; // number of rows
3339 
3340  for (i = 0; i <= height; i++)
3341  {
3342  os << I[i]+1 << '\n';
3343  }
3344 
3345  for (i = 0; i < I[height]; i++)
3346  {
3347  os << J[i]+1 << '\n';
3348  }
3349 
3350  for (i = 0; i < I[height]; i++)
3351  {
3352  os << A[i] << '\n';
3353  }
3354 }
3355 
3356 void SparseMatrix::PrintCSR2(std::ostream & os) const
3357 {
3358  MFEM_VERIFY(Finalized(), "Matrix must be finalized.");
3359 
3360  int i;
3361 
3362  os << height << '\n'; // number of rows
3363  os << width << '\n'; // number of columns
3364 
3365  for (i = 0; i <= height; i++)
3366  {
3367  os << I[i] << '\n';
3368  }
3369 
3370  for (i = 0; i < I[height]; i++)
3371  {
3372  os << J[i] << '\n';
3373  }
3374 
3375  for (i = 0; i < I[height]; i++)
3376  {
3377  os << A[i] << '\n';
3378  }
3379 }
3380 
3381 void SparseMatrix::PrintInfo(std::ostream &os) const
3382 {
3383  const double MiB = 1024.*1024;
3384  int nnz = NumNonZeroElems();
3385  double pz = 100./nnz;
3386  int nz = CountSmallElems(0.0);
3387  double max_norm = MaxNorm();
3388  double symm = IsSymmetric();
3389  int nnf = CheckFinite();
3390  int ns12 = CountSmallElems(1e-12*max_norm);
3391  int ns15 = CountSmallElems(1e-15*max_norm);
3392  int ns18 = CountSmallElems(1e-18*max_norm);
3393 
3394  os <<
3395  "SparseMatrix statistics:\n"
3396  " Format : " <<
3397  (Empty() ? "(empty)" : (Finalized() ? "CSR" : "LIL")) << "\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";
3409  if (Finalized())
3410  {
3411  os << " Memory used by CSR : " <<
3412  (sizeof(int)*(height+1+nnz)+sizeof(double)*nnz)/MiB << " MiB\n";
3413  }
3414  if (Rows != NULL)
3415  {
3416  size_t used_mem = sizeof(RowNode*)*height;
3417 #ifdef MFEM_USE_MEMALLOC
3418  used_mem += NodesMem->MemoryUsage();
3419 #else
3420  for (int i = 0; i < height; i++)
3421  {
3422  for (RowNode *aux = Rows[i]; aux != NULL; aux = aux->Prev)
3423  {
3424  used_mem += sizeof(RowNode);
3425  }
3426  }
3427 #endif
3428  os << " Memory used by LIL : " << used_mem/MiB << " MiB\n";
3429  }
3430 }
3431 
3433 {
3434  I.Delete();
3435  J.Delete();
3436  A.Delete();
3437 
3438  if (Rows != NULL)
3439  {
3440 #if !defined(MFEM_USE_MEMALLOC)
3441  for (int i = 0; i < height; i++)
3442  {
3443  RowNode *aux, *node_p = Rows[i];
3444  while (node_p != NULL)
3445  {
3446  aux = node_p;
3447  node_p = node_p->Prev;
3448  delete aux;
3449  }
3450  }
3451 #endif
3452  delete [] Rows;
3453  }
3454 
3455  delete [] ColPtrJ;
3456  delete [] ColPtrNode;
3457 #ifdef MFEM_USE_MEMALLOC
3458  delete NodesMem;
3459 #endif
3460  delete At;
3461 
3462  ClearGPUSparse();
3463 }
3464 
3466 {
3467  int awidth = 0;
3468  if (A)
3469  {
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)
3473  {
3474  awidth = std::max(awidth, *jptr + 1);
3475  }
3476  }
3477  else
3478  {
3479  RowNode *aux;
3480  for (int i = 0; i < height; i++)
3481  {
3482  for (aux = Rows[i]; aux != NULL; aux = aux->Prev)
3483  {
3484  awidth = std::max(awidth, aux->Column + 1);
3485  }
3486  }
3487  }
3488  return awidth;
3489 }
3490 
3491 void SparseMatrixFunction (SparseMatrix & S, double (*f)(double))
3492 {
3493  int n = S.NumNonZeroElems();
3494  double * s = S.GetData();
3495 
3496  for (int i = 0; i < n; i++)
3497  {
3498  s[i] = f(s[i]);
3499  }
3500 }
3501 
3503 {
3504  MFEM_VERIFY(
3505  A.Finalized(),
3506  "Finalize must be called before Transpose. Use TransposeRowMatrix instead");
3507 
3508  int i, j, end;
3509  const int *A_i, *A_j;
3510  int m, n, nnz, *At_i, *At_j;
3511  const double *A_data;
3512  double *At_data;
3513 
3514  m = A.Height(); // number of rows of A
3515  n = A.Width(); // number of columns of A
3516  nnz = A.NumNonZeroElems();
3517  A_i = A.HostReadI();
3518  A_j = A.HostReadJ();
3519  A_data = A.HostReadData();
3520 
3521  At_i = Memory<int>(n+1);
3522  At_j = Memory<int>(nnz);
3523  At_data = Memory<double>(nnz);
3524 
3525  for (i = 0; i <= n; i++)
3526  {
3527  At_i[i] = 0;
3528  }
3529  for (i = 0; i < nnz; i++)
3530  {
3531  At_i[A_j[i]+1]++;
3532  }
3533  for (i = 1; i < n; i++)
3534  {
3535  At_i[i+1] += At_i[i];
3536  }
3537 
3538  for (i = j = 0; i < m; i++)
3539  {
3540  end = A_i[i+1];
3541  for ( ; j < end; j++)
3542  {
3543  At_j[At_i[A_j[j]]] = i;
3544  At_data[At_i[A_j[j]]] = A_data[j];
3545  At_i[A_j[j]]++;
3546  }
3547  }
3548 
3549  for (i = n; i > 0; i--)
3550  {
3551  At_i[i] = At_i[i-1];
3552  }
3553  At_i[0] = 0;
3554 
3555  return new SparseMatrix(At_i, At_j, At_data, n, m);
3556 }
3557 
3559  int useActualWidth)
3560 {
3561  int i, j;
3562  int m, n, nnz, *At_i, *At_j;
3563  double *At_data;
3564  Array<int> Acols;
3565  Vector Avals;
3566 
3567  m = A.Height(); // number of rows of A
3568  if (useActualWidth)
3569  {
3570  n = 0;
3571  int tmp;
3572  for (i = 0; i < m; i++)
3573  {
3574  A.GetRow(i, Acols, Avals);
3575  if (Acols.Size())
3576  {
3577  tmp = Acols.Max();
3578  if (tmp > n)
3579  {
3580  n = tmp;
3581  }
3582  }
3583  }
3584  ++n;
3585  }
3586  else
3587  {
3588  n = A.Width(); // number of columns of A
3589  }
3590  nnz = A.NumNonZeroElems();
3591 
3592  At_i = Memory<int>(n+1);
3593  At_j = Memory<int>(nnz);
3594  At_data = Memory<double>(nnz);
3595 
3596  for (i = 0; i <= n; i++)
3597  {
3598  At_i[i] = 0;
3599  }
3600 
3601  for (i = 0; i < m; i++)
3602  {
3603  A.GetRow(i, Acols, Avals);
3604  for (j = 0; j<Acols.Size(); ++j)
3605  {
3606  At_i[Acols[j]+1]++;
3607  }
3608  }
3609  for (i = 1; i < n; i++)
3610  {
3611  At_i[i+1] += At_i[i];
3612  }
3613 
3614  for (i = 0; i < m; i++)
3615  {
3616  A.GetRow(i, Acols, Avals);
3617  for (j = 0; j<Acols.Size(); ++j)
3618  {
3619  At_j[At_i[Acols[j]]] = i;
3620  At_data[At_i[Acols[j]]] = Avals[j];
3621  At_i[Acols[j]]++;
3622  }
3623  }
3624 
3625  for (i = n; i > 0; i--)
3626  {
3627  At_i[i] = At_i[i-1];
3628  }
3629  At_i[0] = 0;
3630 
3631  return new SparseMatrix(At_i, At_j, At_data, n, m);
3632 }
3633 
3634 
3636  SparseMatrix *OAB)
3637 {
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;
3642  double *C_data;
3643  int ia, ib, ic, ja, jb, num_nonzeros;
3644  int row_start, counter;
3645  double a_entry, b_entry;
3646  SparseMatrix *C;
3647 
3648  nrowsA = A.Height();
3649  ncolsA = A.Width();
3650  nrowsB = B.Height();
3651  ncolsB = B.Width();
3652 
3653  MFEM_VERIFY(ncolsA == nrowsB,
3654  "number of columns of A (" << ncolsA
3655  << ") must equal number of rows of B (" << nrowsB << ")");
3656 
3657  A_i = A.HostReadI();
3658  A_j = A.HostReadJ();
3659  A_data = A.HostReadData();
3660  B_i = B.HostReadI();
3661  B_j = B.HostReadJ();
3662  B_data = B.HostReadData();
3663 
3664  B_marker = new int[ncolsB];
3665 
3666  for (ib = 0; ib < ncolsB; ib++)
3667  {
3668  B_marker[ib] = -1;
3669  }
3670 
3671  if (OAB == NULL)
3672  {
3673  C_i = Memory<int>(nrowsA+1);
3674 
3675  C_i[0] = num_nonzeros = 0;
3676  for (ic = 0; ic < nrowsA; ic++)
3677  {
3678  for (ia = A_i[ic]; ia < A_i[ic+1]; ia++)
3679  {
3680  ja = A_j[ia];
3681  for (ib = B_i[ja]; ib < B_i[ja+1]; ib++)
3682  {
3683  jb = B_j[ib];
3684  if (B_marker[jb] != ic)
3685  {
3686  B_marker[jb] = ic;
3687  num_nonzeros++;
3688  }
3689  }
3690  }
3691  C_i[ic+1] = num_nonzeros;
3692  }
3693 
3694  C_j = Memory<int>(num_nonzeros);
3695  C_data = Memory<double>(num_nonzeros);
3696 
3697  C = new SparseMatrix(C_i, C_j, C_data, nrowsA, ncolsB);
3698 
3699  for (ib = 0; ib < ncolsB; ib++)
3700  {
3701  B_marker[ib] = -1;
3702  }
3703  }
3704  else
3705  {
3706  C = OAB;
3707 
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());
3714 
3715  // C_i = C->HostReadI(); // not used
3716  C_j = C->HostWriteJ();
3717  C_data = C->HostWriteData();
3718  }
3719 
3720  counter = 0;
3721  for (ic = 0; ic < nrowsA; ic++)
3722  {
3723  // row_start = C_i[ic];
3724  row_start = counter;
3725  for (ia = A_i[ic]; ia < A_i[ic+1]; ia++)
3726  {
3727  ja = A_j[ia];
3728  a_entry = A_data[ia];
3729  for (ib = B_i[ja]; ib < B_i[ja+1]; ib++)
3730  {
3731  jb = B_j[ib];
3732  b_entry = B_data[ib];
3733  if (B_marker[jb] < row_start)
3734  {
3735  B_marker[jb] = counter;
3736  if (OAB == NULL)
3737  {
3738  C_j[counter] = jb;
3739  }
3740  C_data[counter] = a_entry*b_entry;
3741  counter++;
3742  }
3743  else
3744  {
3745  C_data[B_marker[jb]] += a_entry*b_entry;
3746  }
3747  }
3748  }
3749  }
3750 
3751  MFEM_VERIFY(
3752  OAB == NULL || counter == OAB->NumNonZeroElems(),
3753  "With pre-allocated output matrix, number of non-zeros ("
3754  << OAB->NumNonZeroElems()
3755  << ") did not match number of entries changed from matrix-matrix multiply, "
3756  << counter);
3757 
3758  delete [] B_marker;
3759 
3760  return C;
3761 }
3762 
3764 {
3765  SparseMatrix *At = Transpose(A);
3766  SparseMatrix *AtB = Mult(*At, B);
3767  delete At;
3768  return AtB;
3769 }
3770 
3772  const AbstractSparseMatrix &B)
3773 {
3774  int nrowsA, ncolsA, nrowsB, ncolsB;
3775  int *C_i, *C_j, *B_marker;
3776  double *C_data;
3777  int ia, ib, ic, ja, jb, num_nonzeros;
3778  int row_start, counter;
3779  double a_entry, b_entry;
3780  SparseMatrix *C;
3781 
3782  nrowsA = A.Height();
3783  ncolsA = A.Width();
3784  nrowsB = B.Height();
3785  ncolsB = B.Width();
3786 
3787  MFEM_VERIFY(ncolsA == nrowsB,
3788  "number of columns of A (" << ncolsA
3789  << ") must equal number of rows of B (" << nrowsB << ")");
3790 
3791  B_marker = new int[ncolsB];
3792 
3793  for (ib = 0; ib < ncolsB; ib++)
3794  {
3795  B_marker[ib] = -1;
3796  }
3797 
3798  C_i = Memory<int>(nrowsA+1);
3799 
3800  C_i[0] = num_nonzeros = 0;
3801 
3802  Array<int> colsA, colsB;
3803  Vector dataA, dataB;
3804  for (ic = 0; ic < nrowsA; ic++)
3805  {
3806  A.GetRow(ic, colsA, dataA);
3807  for (ia = 0; ia < colsA.Size(); ia++)
3808  {
3809  ja = colsA[ia];
3810  B.GetRow(ja, colsB, dataB);
3811  for (ib = 0; ib < colsB.Size(); ib++)
3812  {
3813  jb = colsB[ib];
3814  if (B_marker[jb] != ic)
3815  {
3816  B_marker[jb] = ic;
3817  num_nonzeros++;
3818  }
3819  }
3820  }
3821  C_i[ic+1] = num_nonzeros;
3822  }
3823 
3824  C_j = Memory<int>(num_nonzeros);
3825  C_data = Memory<double>(num_nonzeros);
3826 
3827  C = new SparseMatrix(C_i, C_j, C_data, nrowsA, ncolsB);
3828 
3829  for (ib = 0; ib < ncolsB; ib++)
3830  {
3831  B_marker[ib] = -1;
3832  }
3833 
3834  counter = 0;
3835  for (ic = 0; ic < nrowsA; ic++)
3836  {
3837  row_start = counter;
3838  A.GetRow(ic, colsA, dataA);
3839  for (ia = 0; ia < colsA.Size(); ia++)
3840  {
3841  ja = colsA[ia];
3842  a_entry = dataA[ia];
3843  B.GetRow(ja, colsB, dataB);
3844  for (ib = 0; ib < colsB.Size(); ib++)
3845  {
3846  jb = colsB[ib];
3847  b_entry = dataB[ib];
3848  if (B_marker[jb] < row_start)
3849  {
3850  B_marker[jb] = counter;
3851  C_j[counter] = jb;
3852  C_data[counter] = a_entry*b_entry;
3853  counter++;
3854  }
3855  else
3856  {
3857  C_data[B_marker[jb]] += a_entry*b_entry;
3858  }
3859  }
3860  }
3861  }
3862 
3863  delete [] B_marker;
3864 
3865  return C;
3866 }
3867 
3869 {
3870  DenseMatrix *C = new DenseMatrix(A.Height(), B.Width());
3871  Vector columnB, columnC;
3872  for (int j = 0; j < B.Width(); ++j)
3873  {
3874  B.GetColumnReference(j, columnB);
3875  C->GetColumnReference(j, columnC);
3876  A.Mult(columnB, columnC);
3877  }
3878  return C;
3879 }
3880 
3882 {
3883  DenseMatrix R (P, 't'); // R = P^T
3884  DenseMatrix *AP = Mult (A, P);
3885  DenseMatrix *RAP_ = new DenseMatrix(R.Height(), AP->Width());
3886  Mult (R, *AP, *RAP_);
3887  delete AP;
3888  return RAP_;
3889 }
3890 
3892 {
3893  SparseMatrix *R = Transpose(P);
3894  DenseMatrix *RA = Mult(*R, A);
3895  DenseMatrix AtP(*RA, 't');
3896  delete RA;
3897  DenseMatrix *RAtP = Mult(*R, AtP);
3898  delete R;
3899  DenseMatrix * RAP_ = new DenseMatrix(*RAtP, 't');
3900  delete RAtP;
3901  return RAP_;
3902 }
3903 
3905  SparseMatrix *ORAP)
3906 {
3907  SparseMatrix *P = Transpose (R);
3908  SparseMatrix *AP = Mult (A, *P);
3909  delete P;
3910  SparseMatrix *RAP_ = Mult (R, *AP, ORAP);
3911  delete AP;
3912  return RAP_;
3913 }
3914 
3916  const SparseMatrix &P)
3917 {
3918  SparseMatrix * R = Transpose(Rt);
3919  SparseMatrix * RA = Mult(*R,A);
3920  delete R;
3921  SparseMatrix * RAP_ = Mult(*RA, P);
3922  delete RA;
3923  return RAP_;
3924 }
3925 
3927  SparseMatrix *OAtDA)
3928 {
3929  int i, At_nnz, *At_j;
3930  double *At_data;
3931 
3932  SparseMatrix *At = Transpose (A);
3933  At_nnz = At -> NumNonZeroElems();
3934  At_j = At -> GetJ();
3935  At_data = At -> GetData();
3936  for (i = 0; i < At_nnz; i++)
3937  {
3938  At_data[i] *= D(At_j[i]);
3939  }
3940  SparseMatrix *AtDA = Mult (*At, A, OAtDA);
3941  delete At;
3942  return AtDA;
3943 }
3944 
3945 SparseMatrix * Add(double a, const SparseMatrix & A, double b,
3946  const SparseMatrix & B)
3947 {
3948  int nrows = A.Height();
3949  int ncols = A.Width();
3950 
3951  int * C_i = Memory<int>(nrows+1);
3952  int * C_j;
3953  double * C_data;
3954 
3955  const int *A_i = A.GetI();
3956  const int *A_j = A.GetJ();
3957  const double *A_data = A.GetData();
3958 
3959  const int *B_i = B.GetI();
3960  const int *B_j = B.GetJ();
3961  const double *B_data = B.GetData();
3962 
3963  int * marker = new int[ncols];
3964  std::fill(marker, marker+ncols, -1);
3965 
3966  int num_nonzeros = 0, jcol;
3967  C_i[0] = 0;
3968  for (int ic = 0; ic < nrows; ic++)
3969  {
3970  for (int ia = A_i[ic]; ia < A_i[ic+1]; ia++)
3971  {
3972  jcol = A_j[ia];
3973  marker[jcol] = ic;
3974  num_nonzeros++;
3975  }
3976  for (int ib = B_i[ic]; ib < B_i[ic+1]; ib++)
3977  {
3978  jcol = B_j[ib];
3979  if (marker[jcol] != ic)
3980  {
3981  marker[jcol] = ic;
3982  num_nonzeros++;
3983  }
3984  }
3985  C_i[ic+1] = num_nonzeros;
3986  }
3987 
3988  C_j = Memory<int>(num_nonzeros);
3989  C_data = Memory<double>(num_nonzeros);
3990 
3991  for (int ia = 0; ia < ncols; ia++)
3992  {
3993  marker[ia] = -1;
3994  }
3995 
3996  int pos = 0;
3997  for (int ic = 0; ic < nrows; ic++)
3998  {
3999  for (int ia = A_i[ic]; ia < A_i[ic+1]; ia++)
4000  {
4001  jcol = A_j[ia];
4002  C_j[pos] = jcol;
4003  C_data[pos] = a*A_data[ia];
4004  marker[jcol] = pos;
4005  pos++;
4006  }
4007  for (int ib = B_i[ic]; ib < B_i[ic+1]; ib++)
4008  {
4009  jcol = B_j[ib];
4010  if (marker[jcol] < C_i[ic])
4011  {
4012  C_j[pos] = jcol;
4013  C_data[pos] = b*B_data[ib];
4014  marker[jcol] = pos;
4015  pos++;
4016  }
4017  else
4018  {
4019  C_data[marker[jcol]] += b*B_data[ib];
4020  }
4021  }
4022  }
4023 
4024  delete[] marker;
4025  return new SparseMatrix(C_i, C_j, C_data, nrows, ncols);
4026 }
4027 
4029 {
4030  return Add(1.,A,1.,B);
4031 }
4032 
4034 {
4035  MFEM_ASSERT(Ai.Size() > 0, "invalid size Ai.Size() = " << Ai.Size());
4036 
4037  SparseMatrix * accumulate = Ai[0];
4038  SparseMatrix * result = accumulate;
4039 
4040  for (int i=1; i < Ai.Size(); ++i)
4041  {
4042  result = Add(*accumulate, *Ai[i]);
4043  if (i != 1)
4044  {
4045  delete accumulate;
4046  }
4047 
4048  accumulate = result;
4049  }
4050 
4051  return result;
4052 }
4053 
4054 /// B += alpha * A
4055 void Add(const SparseMatrix &A,
4056  double alpha, DenseMatrix &B)
4057 {
4058  for (int r = 0; r < B.Height(); r++)
4059  {
4060  const int * colA = A.GetRowColumns(r);
4061  const double * valA = A.GetRowEntries(r);
4062  for (int i=0; i<A.RowSize(r); i++)
4063  {
4064  B(r, colA[i]) += alpha * valA[i];
4065  }
4066  }
4067 }
4068 
4069 /// Produces a block matrix with blocks A_{ij}*B
4071 {
4072  int mA = A.Height(), nA = A.Width();
4073  int mB = B.Height(), nB = B.Width();
4074 
4075  DenseMatrix *C = new DenseMatrix(mA * mB, nA * nB);
4076  *C = 0.0;
4077  for (int i=0; i<mA; i++)
4078  {
4079  for (int j=0; j<nA; j++)
4080  {
4081  C->AddMatrix(A(i,j), B, i * mB, j * nB);
4082  }
4083  }
4084  return C;
4085 }
4086 
4087 /// Produces a block matrix with blocks A_{ij}*B
4089 {
4090  int mA = A.Height(), nA = A.Width();
4091  int mB = B.Height(), nB = B.Width();
4092 
4093  SparseMatrix *C = new SparseMatrix(mA * mB, nA * nB);
4094 
4095  for (int i=0; i<mA; i++)
4096  {
4097  for (int j=0; j<nA; j++)
4098  {
4099  for (int r=0; r<mB; r++)
4100  {
4101  const int * colB = B.GetRowColumns(r);
4102  const double * valB = B.GetRowEntries(r);
4103 
4104  for (int cj=0; cj<B.RowSize(r); cj++)
4105  {
4106  C->Set(i * mB + r, j * nB + colB[cj], A(i,j) * valB[cj]);
4107  }
4108  }
4109  }
4110  }
4111  C->Finalize();
4112 
4113  return C;
4114 }
4115 
4116 /// Produces a block matrix with blocks A_{ij}*B
4118 {
4119  int mA = A.Height(), nA = A.Width();
4120  int mB = B.Height(), nB = B.Width();
4121 
4122  SparseMatrix *C = new SparseMatrix(mA * mB, nA * nB);
4123 
4124  for (int r=0; r<mA; r++)
4125  {
4126  const int * colA = A.GetRowColumns(r);
4127  const double * valA = A.GetRowEntries(r);
4128 
4129  for (int aj=0; aj<A.RowSize(r); aj++)
4130  {
4131  for (int i=0; i<mB; i++)
4132  {
4133  for (int j=0; j<nB; j++)
4134  {
4135  C->Set(r * mB + i, colA[aj] * nB + j, valA[aj] * B(i, j));
4136  }
4137  }
4138  }
4139  }
4140  C->Finalize();
4141 
4142  return C;
4143 }
4144 
4145 /// Produces a block matrix with blocks A_{ij}*B
4147 {
4148  int mA = A.Height(), nA = A.Width();
4149  int mB = B.Height(), nB = B.Width();
4150 
4151  SparseMatrix *C = new SparseMatrix(mA * mB, nA * nB);
4152 
4153  for (int ar=0; ar<mA; ar++)
4154  {
4155  const int * colA = A.GetRowColumns(ar);
4156  const double * valA = A.GetRowEntries(ar);
4157 
4158  for (int aj=0; aj<A.RowSize(ar); aj++)
4159  {
4160  for (int br=0; br<mB; br++)
4161  {
4162  const int * colB = B.GetRowColumns(br);
4163  const double * valB = B.GetRowEntries(br);
4164 
4165  for (int bj=0; bj<B.RowSize(br); bj++)
4166  {
4167  C->Set(ar * mB + br, colA[aj] * nB + colB[bj],
4168  valA[aj] * valB[bj]);
4169  }
4170  }
4171  }
4172  }
4173  C->Finalize();
4174 
4175  return C;
4176 }
4177 
4179 {
4180  mfem::Swap(width, other.width);
4181  mfem::Swap(height, other.height);
4182  mfem::Swap(I, other.I);
4183  mfem::Swap(J, other.J);
4184  mfem::Swap(A, other.A);
4185  mfem::Swap(Rows, other.Rows);
4187  mfem::Swap(ColPtrJ, other.ColPtrJ);
4189  mfem::Swap(At, other.At);
4190 
4191 #ifdef MFEM_USE_MEMALLOC
4192  mfem::Swap(NodesMem, other.NodesMem);
4193 #endif
4194 
4195  mfem::Swap(isSorted, other.isSorted);
4196 }
4197 
4199 {
4200  Destroy();
4201 #ifdef MFEM_USE_CUDA_OR_HIP
4203  {
4204  if (SparseMatrixCount==1)
4205  {
4206  if (handle)
4207  {
4208  MFEM_cu_or_hip(sparseDestroy)(handle);
4209  handle = nullptr;
4210  }
4211  if (dBuffer)
4212  {
4213  MFEM_Cu_or_Hip(MemFree)(dBuffer);
4214  dBuffer = nullptr;
4215  bufferSize = 0;
4216  }
4217  }
4219  }
4220 #endif // MFEM_USE_CUDA_OR_HIP
4221 }
4222 
4223 }
RowNode ** ColPtrNode
Definition: sparsemat.hpp:77
Memory< int > I
Array with size (height+1) containing the row offsets.
Definition: sparsemat.hpp:62
Elem * Alloc()
Definition: mem_alloc.hpp:166
double GetJacobiScaling() const
Determine appropriate scaling for Jacobi iteration.
Definition: sparsemat.cpp:2551
int CheckFinite() const
Count the number of entries that are NOT finite, i.e. Inf or Nan.
Definition: sparsemat.cpp:1660
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
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.
Definition: sparsemat.cpp:2616
SparseMatrix * TransposeAbstractSparseMatrix(const AbstractSparseMatrix &A, int useActualWidth)
Transpose of a sparse matrix. A does not need to be a CSR matrix.
Definition: sparsemat.cpp:3558
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:344
int CheckFinite(const double *v, const int n)
Definition: vector.hpp:493
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:1584
void * CuMemFree(void *dptr)
Frees device memory and returns destination ptr.
Definition: cuda.cpp:79
void _Add_(const int col, const double a)
Add a value to an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:584
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
Definition: sparsemat.cpp:956
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:163
void Clear()
Clear the contents of the SparseMatrix.
Definition: sparsemat.hpp:209
void EliminateCol(int col, DiagonalPolicy dpolicy=DIAG_ZERO)
Eliminates the column col from the matrix.
Definition: sparsemat.cpp:1737
Memory< T > & GetMemory()
Return a reference to the Memory object used by the Array.
Definition: array.hpp:117
SparseMatrix * At
Transpose of A. Owned. Used to perform MultTranspose() on devices.
Definition: sparsemat.hpp:80
int NumCols() const
Definition: array.hpp:373
static size_t bufferSize
Definition: sparsemat.hpp:101
size_t MemoryUsage() const
Definition: mem_alloc.hpp:212
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:545
void MakeRef(const SparseMatrix &master)
Clear the contents of the SparseMatrix and make it a reference to master.
Definition: sparsemat.cpp:313
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
double & SearchRow(const int col)
Perform a fast search for an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:845
void Delete()
Delete the owned pointers and reset the Memory object.
DenseMatrix * ToDenseMatrix() const
Produces a DenseMatrix from a SparseMatrix.
Definition: sparsemat.cpp:699
bool Empty() const
Check if the SparseMatrix is empty.
Definition: sparsemat.hpp:219
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:391
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:1022
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:73
const double * ReadData(bool on_dev=true) const
Definition: sparsemat.hpp:271
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:227
Biwise-OR of all HIP backends.
Definition: device.hpp:90
Abstract data type for sparse matrices.
Definition: matrix.hpp:73
SparseMatrix & operator+=(const SparseMatrix &B)
Add the sparse matrix &#39;B&#39; to &#39;*this&#39;. This operation will cause an error if &#39;*this&#39; is finalized and ...
Definition: sparsemat.cpp:3137
T * Write(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for write access to mem with the mfem::Device&#39;s DeviceMemoryClass, if on_dev = true...
Definition: device.hpp:336
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
void AddMult(const Vector &x, Vector &y, const double a=1.0) const
y += A * x (default) or y += a * A * x
Definition: sparsemat.cpp:736
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
double MaxNorm() const
Definition: sparsemat.cpp:1607
int * GetI()
Return the array I.
Definition: sparsemat.hpp:222
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
y += At * x (default) or y += a * At * x
Definition: sparsemat.cpp:911
Abstract data type for matrix inverse.
Definition: matrix.hpp:62
static void * dBuffer
Definition: sparsemat.hpp:102
void EliminateRowColMultipleRHS(int rc, const Vector &sol, DenseMatrix &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Similar to EliminateRowCol(int, const double, Vector &amp;, DiagonalPolicy), but multiple values for elim...
Definition: sparsemat.cpp:1949
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...
Definition: sparsemat.cpp:81
void PrintInfo(std::ostream &out) const
Print various sparse matrix statistics.
Definition: sparsemat.cpp:3381
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}.
Definition: sparsemat.cpp:597
void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc, bool use_abs_diag=false) const
Definition: sparsemat.cpp:2584
void AddRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:3018
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
Definition: sparsemat.cpp:405
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:209
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:457
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:118
void PrintMM(std::ostream &out=mfem::out) const
Prints matrix in Matrix Market sparse format.
Definition: sparsemat.cpp:3310
int Capacity() const
Return the size of the allocated memory.
cusparseDnVecDescr_t vecY_descr
Definition: sparsemat.hpp:113
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
Definition: sparsemat.cpp:2882
double * HostReadWriteData()
Definition: sparsemat.hpp:281
void _Set_(const int col, const double a)
Set an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:587
void BooleanMultTranspose(const Array< int > &x, Array< int > &y) const
y = At * x, treating all entries as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:1052
SparseMatrix * Mult_AtDA(const SparseMatrix &A, const Vector &D, SparseMatrix *OAtDA)
Matrix multiplication A^t D A. All matrices must be finalized.
Definition: sparsemat.cpp:3926
void PartAddMult(const Array< int > &rows, const Vector &x, Vector &y, const double a=1.0) const
Definition: sparsemat.cpp:1004
bool RowIsEmpty(const int row) const
Definition: sparsemat.cpp:2910
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
Definition: sparsemat.cpp:1433
const int * ReadJ(bool on_dev=true) const
Definition: sparsemat.hpp:255
void ScaleRow(const int row, const double scale)
Definition: sparsemat.cpp:3050
RowNodeAlloc * NodesMem
Definition: sparsemat.hpp:84
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:232
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
Definition: sparsemat.cpp:1565
bool IsFinite(const double &val)
Definition: vector.hpp:482
void MoveDiagonalFirst()
Move the diagonal entry to the first position in each row, preserving the order of the rest of the co...
Definition: sparsemat.cpp:573
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2282
void Set(const int i, const int j, const double val)
Definition: sparsemat.cpp:2764
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 &quot;current row&quot; set by calling SetColPtr(). This method must be called between any two calls ...
Definition: sparsemat.hpp:826
void ScaleColumns(const Vector &sr)
this = this * diag(sr);
Definition: sparsemat.cpp:3109
virtual ~SparseMatrix()
Destroys sparse matrix.
Definition: sparsemat.cpp:4198
void PrintCSR2(std::ostream &out) const
Prints a sparse matrix to stream out in CSR format.
Definition: sparsemat.cpp:3356
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
virtual int GetRow(const int row, Array< int > &cols, Vector &srow) const
Extract all column indices and values from a given row.
Definition: sparsemat.cpp:2931
bool isSorted
Are the columns sorted already.
Definition: sparsemat.hpp:88
Memory< double > A
Array with size I[height], containing the actual entries of the sparse matrix, as indexed by the I ar...
Definition: sparsemat.hpp:68
SparseMatrix()
Create an empty SparseMatrix.
Definition: sparsemat.hpp:131
Data type sparse matrix.
Definition: sparsemat.hpp:50
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
void EnsureMultTranspose() const
Ensures that the matrix is capable of performing MultTranspose(), AddMultTranspose(), and AbsMultTranspose().
Definition: sparsemat.cpp:970
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:457
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:154
void SetWidth(int width_=-1)
Change the width of a SparseMatrix.
Definition: sparsemat.cpp:419
double b
Definition: lissajous.cpp:42
const double * HostReadData() const
Definition: sparsemat.hpp:277
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition: device.hpp:273
double IsSymmetric() const
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.
Definition: sparsemat.cpp:1523
void * HipMemFree(void *dptr)
Frees device memory.
Definition: hip.cpp:79
const int * HostReadJ() const
Definition: sparsemat.hpp:261
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
static int SparseMatrixCount
Definition: sparsemat.hpp:100
void SetDiagIdentity()
If a row contains only one diag entry of zero, set it to 1.
Definition: sparsemat.cpp:2354
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:304
void ScaleRows(const Vector &sl)
this = diag(sl) * this;
Definition: sparsemat.cpp:3078
static cusparseHandle_t handle
Definition: sparsemat.hpp:107
DenseMatrix * OuterProduct(const DenseMatrix &A, const DenseMatrix &B)
Produces a block matrix with blocks A_{ij}*B.
Definition: sparsemat.cpp:4070
void Sort()
Sorts the array in ascending order. This requires operator&lt; to be defined for T.
Definition: array.hpp:248
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.
Definition: operator.hpp:50
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
Definition: densemat.cpp:3213
virtual void PrintMatlab(std::ostream &out=mfem::out) const
Prints matrix in matlab format.
Definition: sparsemat.cpp:3288
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&#39;s DeviceMemoryClass, if on_dev = true...
Definition: device.hpp:319
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).
Definition: vector.hpp:453
Biwise-OR of all CUDA backends.
Definition: device.hpp:88
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2718
void Add(const int i, const int j, const double val)
Definition: sparsemat.cpp:2783
Dynamic 2D array using row-major layout.
Definition: array.hpp:351
void SetSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2802
bool Finalized() const
Returns whether or not CSR format has been finalized.
Definition: sparsemat.hpp:551
Biwise-OR of all CPU backends.
Definition: device.hpp:86
void * HipMemAlloc(void **dptr, size_t bytes)
Allocates device memory.
Definition: hip.cpp:34
T * HostWrite(Memory< T > &mem, int size)
Shortcut to Write(const Memory&lt;T&gt; &amp;mem, int size, false)
Definition: device.hpp:343
void SparseMatrixFunction(SparseMatrix &S, double(*f)(double))
Applies f() to each element of the matrix (after it is finalized).
Definition: sparsemat.cpp:3491
virtual MatrixInverse * Inverse() const
This virtual method is not supported: it always returns NULL.
Definition: sparsemat.cpp:1684
void Swap(Array< T > &, Array< T > &)
Definition: array.hpp:630
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.
Definition: sparsemat.cpp:2154
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:729
void SetColPtr(const int row) const
Initialize the SparseMatrix for fast access to the entries of the given row which becomes the &quot;curren...
Definition: sparsemat.hpp:791
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
void SetSubMatrixTranspose(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2841
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.
Definition: array.hpp:679
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.
Definition: sparsemat.cpp:1850
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition: table.cpp:413
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:289
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...
Definition: sparsemat.cpp:1257
double * HostWriteData()
Definition: sparsemat.hpp:279
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0&lt;=i&lt;A.Height, 0&lt;=j&lt;A.Width.
Definition: densemat.cpp:1684
void Gauss_Seidel_back(const Vector &x, Vector &y) const
Definition: sparsemat.cpp:2468
const int * ReadI(bool on_dev=true) const
Definition: sparsemat.hpp:239
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
Definition: device.hpp:258
void Gauss_Seidel_forw(const Vector &x, Vector &y) const
Gauss-Seidel forward and backward iterations over a vector x.
Definition: sparsemat.cpp:2384
int * HostReadWriteJ()
Definition: sparsemat.hpp:265
RowNode ** Rows
Array of linked lists, one for every row. This array represents the linked list (LIL) storage format...
Definition: sparsemat.hpp:73
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
void EliminateRow(int row, const double sol, Vector &rhs)
Eliminates a column from the transpose matrix.
Definition: sparsemat.cpp:1689
double a
Definition: lissajous.cpp:41
const T * HostRead(const Memory< T > &mem, int size)
Shortcut to Read(const Memory&lt;T&gt; &amp;mem, int size, false)
Definition: device.hpp:326
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:465
int ActualWidth() const
Returns the actual Width of the matrix.
Definition: sparsemat.cpp:3465
void GetDiag(Vector &d) const
Returns the Diagonal of A.
Definition: sparsemat.cpp:665
MemAlloc< RowNode, 1024 > RowNodeAlloc
Definition: sparsemat.hpp:83
double GetRowNorml1(int irow) const
For i = irow compute .
Definition: sparsemat.cpp:1233
void ResetTranspose() const
Definition: sparsemat.cpp:964
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
Definition: sparsemat.cpp:3763
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
Keep the diagonal value.
Definition: operator.hpp:51
virtual void EliminateZeroRows(const double threshold=1e-12)
If a row contains only zeros, set its diagonal to 1.
Definition: sparsemat.cpp:2365
void AbsMultTranspose(const Vector &x, Vector &y) const
y = |At| * x, using entry-wise absolute values of the transpose of matrix A
Definition: sparsemat.cpp:1124
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition: operator.hpp:47
double & operator()(int i, int j)
Returns reference to A[i][j].
Definition: sparsemat.cpp:607
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:904
void EliminateCols(const Array< int > &cols, const Vector *x=NULL, Vector *b=NULL)
Eliminate all columns i for which cols[i] != 0.
Definition: sparsemat.cpp:1777
void OverrideSize(int height_, int width_)
Sets the height and width of the matrix.
Definition: sparsemat.cpp:297
SparseMatrix & operator=(const SparseMatrix &rhs)
Assignment operator: deep copy.
Definition: sparsemat.cpp:303
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
SparseMatrix & operator*=(double a)
Definition: sparsemat.cpp:3216
RefCoord t[3]
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:2971
cusparseSpMatDescr_t matA_descr
Definition: sparsemat.hpp:111
T * Write(bool on_dev=true)
Shortcut for mfem::Write(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:312
void MakeRef(T *, int)
Make this Array a reference to a pointer.
Definition: array.hpp:864
const double alpha
Definition: ex15.cpp:369
double * ReadWriteData(bool on_dev=true)
Definition: sparsemat.hpp:275
int * ReadWriteJ(bool on_dev=true)
Definition: sparsemat.hpp:259
Vector data type.
Definition: vector.hpp:60
void Print(std::ostream &out=mfem::out, int width_=4) const
Prints matrix to stream out.
Definition: sparsemat.cpp:3240
const int * HostReadI() const
Definition: sparsemat.hpp:245
void Jacobi2(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:2704
cusparseDnVecDescr_t vecX_descr
Definition: sparsemat.hpp:112
int MaxRowSize() const
Returns the maximum number of elements among all rows.
Definition: sparsemat.cpp:367
int * HostReadWriteI()
Definition: sparsemat.hpp:249
RefCoord s[3]
void Jacobi3(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:2711
SparseMatrix * MultAbstractSparseMatrix(const AbstractSparseMatrix &A, const AbstractSparseMatrix &B)
Matrix product of sparse matrices. A and B do not need to be CSR matrices.
Definition: sparsemat.cpp:3771
void AbsMult(const Vector &x, Vector &y) const
y = |A| * x, using entry-wise absolute values of matrix A
Definition: sparsemat.cpp:1075
double _Get_(const int col) const
Read the value of an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:872
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
Memory< int > J
Array with size I[height], containing the column indices for all matrix entries, as indexed by the I ...
Definition: sparsemat.hpp:65
void Swap(SparseMatrix &other)
Definition: sparsemat.cpp:4178
void PrintCSR(std::ostream &out) const
Prints matrix to stream out in hypre_CSRMatrix format.
Definition: sparsemat.cpp:3332
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:1169
void * CuMemAlloc(void **dptr, size_t bytes)
Allocates device memory and returns destination ptr.
Definition: cuda.cpp:34
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:449
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:469
int NumRows() const
Definition: array.hpp:372
Set the diagonal value to zero.
Definition: operator.hpp:49
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
void EliminateBC(const Array< int > &ess_dofs, DiagonalPolicy diag_policy)
Eliminate essential (Dirichlet) boundary conditions.
Definition: sparsemat.cpp:2311
int CountSmallElems(double tol) const
Count the number of entries with |a_ij| &lt;= tol.
Definition: sparsemat.cpp:1632
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
Definition: sparsemat.cpp:978
void GetRowSums(Vector &x) const
For all i compute .
Definition: sparsemat.cpp:1210
void Neg()
(*this) = -(*this)
Definition: vector.cpp:305