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