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