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