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