xref: /petsc/src/mat/impls/sell/seq/sell.c (revision 07e43b41f7cb79bba1664e5473a3e2df04f149fa)
1d4002b98SHong Zhang 
2d4002b98SHong Zhang /*
3d4002b98SHong Zhang   Defines the basic matrix operations for the SELL matrix storage format.
4d4002b98SHong Zhang */
5d4002b98SHong Zhang #include <../src/mat/impls/sell/seq/sell.h> /*I   "petscmat.h"  I*/
6d4002b98SHong Zhang #include <petscblaslapack.h>
7d4002b98SHong Zhang #include <petsc/private/kernels/blocktranspose.h>
8ed73aabaSBarry Smith 
9ed73aabaSBarry Smith static PetscBool  cited      = PETSC_FALSE;
109371c9d4SSatish Balay static const char citation[] = "@inproceedings{ZhangELLPACK2018,\n"
11ed73aabaSBarry Smith                                " author = {Hong Zhang and Richard T. Mills and Karl Rupp and Barry F. Smith},\n"
12ed73aabaSBarry Smith                                " title = {Vectorized Parallel Sparse Matrix-Vector Multiplication in {PETSc} Using {AVX-512}},\n"
13ed73aabaSBarry Smith                                " booktitle = {Proceedings of the 47th International Conference on Parallel Processing},\n"
14ed73aabaSBarry Smith                                " year = 2018\n"
15ed73aabaSBarry Smith                                "}\n";
16ed73aabaSBarry Smith 
175f70456aSHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && (defined(__AVX512F__) || (defined(__AVX2__) && defined(__FMA__)) || defined(__AVX__)) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
184243e2ceSHong Zhang 
19d4002b98SHong Zhang   #include <immintrin.h>
20d4002b98SHong Zhang 
21d4002b98SHong Zhang   #if !defined(_MM_SCALE_8)
22d4002b98SHong Zhang     #define _MM_SCALE_8 8
23d4002b98SHong Zhang   #endif
24d4002b98SHong Zhang 
25d4002b98SHong Zhang   #if defined(__AVX512F__)
26d4002b98SHong Zhang     /* these do not work
27d4002b98SHong Zhang    vec_idx  = _mm512_loadunpackhi_epi32(vec_idx,acolidx);
28d4002b98SHong Zhang    vec_vals = _mm512_loadunpackhi_pd(vec_vals,aval);
29d4002b98SHong Zhang   */
30d4002b98SHong Zhang     #define AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y) \
31d4002b98SHong Zhang       /* if the mask bit is set, copy from acolidx, otherwise from vec_idx */ \
32ef588d5cSRichard Tran Mills       vec_idx  = _mm256_loadu_si256((__m256i const *)acolidx); \
33ef588d5cSRichard Tran Mills       vec_vals = _mm512_loadu_pd(aval); \
34d4002b98SHong Zhang       vec_x    = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8); \
35a48a6482SHong Zhang       vec_y    = _mm512_fmadd_pd(vec_x, vec_vals, vec_y)
365f70456aSHong Zhang   #elif defined(__AVX2__) && defined(__FMA__)
37a48a6482SHong Zhang     #define AVX2_Mult_Private(vec_idx, vec_x, vec_vals, vec_y) \
38ef588d5cSRichard Tran Mills       vec_vals = _mm256_loadu_pd(aval); \
39ef588d5cSRichard Tran Mills       vec_idx  = _mm_loadu_si128((__m128i const *)acolidx); /* SSE2 */ \
40a48a6482SHong Zhang       vec_x    = _mm256_i32gather_pd(x, vec_idx, _MM_SCALE_8); \
41a48a6482SHong Zhang       vec_y    = _mm256_fmadd_pd(vec_x, vec_vals, vec_y)
42d4002b98SHong Zhang   #endif
43d4002b98SHong Zhang #endif /* PETSC_HAVE_IMMINTRIN_H */
44d4002b98SHong Zhang 
45d4002b98SHong Zhang /*@C
46d4002b98SHong Zhang  MatSeqSELLSetPreallocation - For good matrix assembly performance
4720f4b53cSBarry Smith  the user should preallocate the matrix storage by setting the parameter `nz`
4820f4b53cSBarry Smith  (or the array `nnz`).
49d4002b98SHong Zhang 
50d083f849SBarry Smith  Collective
51d4002b98SHong Zhang 
52d4002b98SHong Zhang  Input Parameters:
5311a5261eSBarry Smith +  B - The `MATSEQSELL` matrix
5420f4b53cSBarry Smith .  rlenmax - number of nonzeros per row (same for all rows), ignored if `rlen` is provided
5520f4b53cSBarry Smith -  rlen - array containing the number of nonzeros in the various rows (possibly different for each row) or `NULL`
5667be906fSBarry Smith 
5767be906fSBarry Smith  Level: intermediate
58d4002b98SHong Zhang 
59d4002b98SHong Zhang  Notes:
6067be906fSBarry Smith  Specify the preallocated storage with either `rlenmax` or `rlen` (not both).
6167be906fSBarry Smith  Set `rlenmax` = `PETSC_DEFAULT` and `rlen` = `NULL` for PETSc to control dynamic memory
6220f4b53cSBarry Smith  allocation.
63d4002b98SHong Zhang 
6411a5261eSBarry Smith  You can call `MatGetInfo()` to get information on how effective the preallocation was;
65d4002b98SHong Zhang  for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
6667be906fSBarry Smith  You can also run with the option `-info` and look for messages with the string
67d4002b98SHong Zhang  malloc in them to see if additional memory allocation was needed.
68d4002b98SHong Zhang 
6927430b45SBarry Smith  Developer Note:
7067be906fSBarry Smith  Use `rlenmax` of `MAT_SKIP_ALLOCATION` to not allocate any space for the matrix
71d4002b98SHong Zhang  entries or columns indices.
72d4002b98SHong Zhang 
73c7ee91abSRichard Tran Mills  The maximum number of nonzeos in any row should be as accurate as possible.
74c7ee91abSRichard Tran Mills  If it is underestimated, you will get bad performance due to reallocation
7567be906fSBarry Smith  (`MatSeqXSELLReallocateSELL()`).
76d4002b98SHong Zhang 
7767be906fSBarry Smith  .seealso: `Mat`, `MATSEQSELL`, `MATSELL`, `MatCreate()`, `MatCreateSELL()`, `MatSetValues()`, `MatGetInfo()`
78d4002b98SHong Zhang  @*/
79d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLSetPreallocation(Mat B, PetscInt rlenmax, const PetscInt rlen[])
80d71ae5a4SJacob Faibussowitsch {
81d4002b98SHong Zhang   PetscFunctionBegin;
82d4002b98SHong Zhang   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
83d4002b98SHong Zhang   PetscValidType(B, 1);
84cac4c232SBarry Smith   PetscTryMethod(B, "MatSeqSELLSetPreallocation_C", (Mat, PetscInt, const PetscInt[]), (B, rlenmax, rlen));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86d4002b98SHong Zhang }
87d4002b98SHong Zhang 
88d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLSetPreallocation_SeqSELL(Mat B, PetscInt maxallocrow, const PetscInt rlen[])
89d71ae5a4SJacob Faibussowitsch {
90d4002b98SHong Zhang   Mat_SeqSELL *b;
91d4002b98SHong Zhang   PetscInt     i, j, totalslices;
92d4002b98SHong Zhang   PetscBool    skipallocation = PETSC_FALSE, realalloc = PETSC_FALSE;
93d4002b98SHong Zhang 
94d4002b98SHong Zhang   PetscFunctionBegin;
95d4002b98SHong Zhang   if (maxallocrow >= 0 || rlen) realalloc = PETSC_TRUE;
96d4002b98SHong Zhang   if (maxallocrow == MAT_SKIP_ALLOCATION) {
97d4002b98SHong Zhang     skipallocation = PETSC_TRUE;
98d4002b98SHong Zhang     maxallocrow    = 0;
99d4002b98SHong Zhang   }
100d4002b98SHong Zhang 
1019566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
1029566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
103d4002b98SHong Zhang 
104d4002b98SHong Zhang   /* FIXME: if one preallocates more space than needed, the matrix does not shrink automatically, but for best performance it should */
105d4002b98SHong Zhang   if (maxallocrow == PETSC_DEFAULT || maxallocrow == PETSC_DECIDE) maxallocrow = 5;
10608401ef6SPierre Jolivet   PetscCheck(maxallocrow >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "maxallocrow cannot be less than 0: value %" PetscInt_FMT, maxallocrow);
107d4002b98SHong Zhang   if (rlen) {
108d4002b98SHong Zhang     for (i = 0; i < B->rmap->n; i++) {
10908401ef6SPierre Jolivet       PetscCheck(rlen[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "rlen cannot be less than 0: local row %" PetscInt_FMT " value %" PetscInt_FMT, i, rlen[i]);
11008401ef6SPierre Jolivet       PetscCheck(rlen[i] <= B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "rlen cannot be greater than row length: local row %" PetscInt_FMT " value %" PetscInt_FMT " rowlength %" PetscInt_FMT, i, rlen[i], B->cmap->n);
111d4002b98SHong Zhang     }
112d4002b98SHong Zhang   }
113d4002b98SHong Zhang 
114d4002b98SHong Zhang   B->preallocated = PETSC_TRUE;
115d4002b98SHong Zhang 
116d4002b98SHong Zhang   b = (Mat_SeqSELL *)B->data;
117d4002b98SHong Zhang 
118*07e43b41SHong Zhang   if (!b->sliceheight) { /* not set yet */
119*07e43b41SHong Zhang #if defined(PETSC_HAVE_CUDA)
120*07e43b41SHong Zhang     b->sliceheight = 16;
121*07e43b41SHong Zhang #else
122*07e43b41SHong Zhang     b->sliceheight = 8;
123*07e43b41SHong Zhang #endif
124*07e43b41SHong Zhang   }
125*07e43b41SHong Zhang   totalslices    = PetscCeilInt(B->rmap->n, b->sliceheight);
126d4002b98SHong Zhang   b->totalslices = totalslices;
127d4002b98SHong Zhang   if (!skipallocation) {
128*07e43b41SHong Zhang     if (B->rmap->n % b->sliceheight) PetscCall(PetscInfo(B, "Padding rows to the SEQSELL matrix because the number of rows is not the multiple of the slice height (value %" PetscInt_FMT ")\n", B->rmap->n));
129d4002b98SHong Zhang 
130d4002b98SHong Zhang     if (!b->sliidx) { /* sliidx gives the starting index of each slice, the last element is the total space allocated */
1319566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(totalslices + 1, &b->sliidx));
132d4002b98SHong Zhang     }
133d4002b98SHong Zhang     if (!rlen) { /* if rlen is not provided, allocate same space for all the slices */
134d4002b98SHong Zhang       if (maxallocrow == PETSC_DEFAULT || maxallocrow == PETSC_DECIDE) maxallocrow = 10;
135d4002b98SHong Zhang       else if (maxallocrow < 0) maxallocrow = 1;
136*07e43b41SHong Zhang       for (i = 0; i <= totalslices; i++) b->sliidx[i] = b->sliceheight * i * maxallocrow;
137d4002b98SHong Zhang     } else {
138d4002b98SHong Zhang       maxallocrow  = 0;
139d4002b98SHong Zhang       b->sliidx[0] = 0;
140d4002b98SHong Zhang       for (i = 1; i < totalslices; i++) {
141d4002b98SHong Zhang         b->sliidx[i] = 0;
142*07e43b41SHong Zhang         for (j = 0; j < b->sliceheight; j++) { b->sliidx[i] = PetscMax(b->sliidx[i], rlen[b->sliceheight * (i - 1) + j]); }
143d4002b98SHong Zhang         maxallocrow = PetscMax(b->sliidx[i], maxallocrow);
144*07e43b41SHong Zhang         PetscCall(PetscIntSumError(b->sliidx[i - 1], b->sliceheight * b->sliidx[i], &b->sliidx[i]));
145d4002b98SHong Zhang       }
146d4002b98SHong Zhang       /* last slice */
147d4002b98SHong Zhang       b->sliidx[totalslices] = 0;
148*07e43b41SHong Zhang       for (j = b->sliceheight * (totalslices - 1); j < B->rmap->n; j++) b->sliidx[totalslices] = PetscMax(b->sliidx[totalslices], rlen[j]);
149d4002b98SHong Zhang       maxallocrow            = PetscMax(b->sliidx[totalslices], maxallocrow);
150*07e43b41SHong Zhang       b->sliidx[totalslices] = b->sliidx[totalslices - 1] + b->sliceheight * b->sliidx[totalslices];
151d4002b98SHong Zhang     }
152d4002b98SHong Zhang 
153d4002b98SHong Zhang     /* allocate space for val, colidx, rlen */
154d4002b98SHong Zhang     /* FIXME: should B's old memory be unlogged? */
1559566063dSJacob Faibussowitsch     PetscCall(MatSeqXSELLFreeSELL(B, &b->val, &b->colidx));
156d4002b98SHong Zhang     /* FIXME: assuming an element of the bit array takes 8 bits */
1579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(b->sliidx[totalslices], &b->val, b->sliidx[totalslices], &b->colidx));
158d4002b98SHong Zhang     /* b->rlen will count nonzeros in each row so far. We dont copy rlen to b->rlen because the matrix has not been set. */
159*07e43b41SHong Zhang     PetscCall(PetscCalloc1(b->sliceheight * totalslices, &b->rlen));
160d4002b98SHong Zhang 
161d4002b98SHong Zhang     b->singlemalloc = PETSC_TRUE;
162d4002b98SHong Zhang     b->free_val     = PETSC_TRUE;
163d4002b98SHong Zhang     b->free_colidx  = PETSC_TRUE;
164d4002b98SHong Zhang   } else {
165d4002b98SHong Zhang     b->free_val    = PETSC_FALSE;
166d4002b98SHong Zhang     b->free_colidx = PETSC_FALSE;
167d4002b98SHong Zhang   }
168d4002b98SHong Zhang 
169d4002b98SHong Zhang   b->nz               = 0;
170d4002b98SHong Zhang   b->maxallocrow      = maxallocrow;
171d4002b98SHong Zhang   b->rlenmax          = maxallocrow;
172d4002b98SHong Zhang   b->maxallocmat      = b->sliidx[totalslices];
173d4002b98SHong Zhang   B->info.nz_unneeded = (double)b->maxallocmat;
1741baa6e33SBarry Smith   if (realalloc) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
1753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
176d4002b98SHong Zhang }
177d4002b98SHong Zhang 
178d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetRow_SeqSELL(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
179d71ae5a4SJacob Faibussowitsch {
1806108893eSStefano Zampini   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1816108893eSStefano Zampini   PetscInt     shift;
1826108893eSStefano Zampini 
1836108893eSStefano Zampini   PetscFunctionBegin;
184aed4548fSBarry Smith   PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " out of range", row);
1856108893eSStefano Zampini   if (nz) *nz = a->rlen[row];
186*07e43b41SHong Zhang   shift = a->sliidx[row / a->sliceheight] + (row % a->sliceheight);
1872d1451d4SHong Zhang   if (!a->getrowcols) { PetscCall(PetscMalloc2(a->rlenmax, &a->getrowcols, a->rlenmax, &a->getrowvals)); }
1886108893eSStefano Zampini   if (idx) {
1896108893eSStefano Zampini     PetscInt j;
190*07e43b41SHong Zhang     for (j = 0; j < a->rlen[row]; j++) a->getrowcols[j] = a->colidx[shift + a->sliceheight * j];
1916108893eSStefano Zampini     *idx = a->getrowcols;
1926108893eSStefano Zampini   }
1936108893eSStefano Zampini   if (v) {
1946108893eSStefano Zampini     PetscInt j;
195*07e43b41SHong Zhang     for (j = 0; j < a->rlen[row]; j++) a->getrowvals[j] = a->val[shift + a->sliceheight * j];
1966108893eSStefano Zampini     *v = a->getrowvals;
1976108893eSStefano Zampini   }
1983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1996108893eSStefano Zampini }
2006108893eSStefano Zampini 
201d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRestoreRow_SeqSELL(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
202d71ae5a4SJacob Faibussowitsch {
2036108893eSStefano Zampini   PetscFunctionBegin;
2043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2056108893eSStefano Zampini }
2066108893eSStefano Zampini 
207d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_SeqSELL_SeqAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
208d71ae5a4SJacob Faibussowitsch {
209d4002b98SHong Zhang   Mat          B;
210d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
211e3f1f374SStefano Zampini   PetscInt     i;
212d4002b98SHong Zhang 
213d4002b98SHong Zhang   PetscFunctionBegin;
214ad013a7bSRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
215ad013a7bSRichard Tran Mills     B = *newmat;
2169566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(B));
217ad013a7bSRichard Tran Mills   } else {
2189566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2199566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
2209566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQAIJ));
2219566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, a->rlen));
222ad013a7bSRichard Tran Mills   }
223d4002b98SHong Zhang 
224e3f1f374SStefano Zampini   for (i = 0; i < A->rmap->n; i++) {
225e108cb99SStefano Zampini     PetscInt     nz = 0, *cols = NULL;
226e108cb99SStefano Zampini     PetscScalar *vals = NULL;
227e3f1f374SStefano Zampini 
2289566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqSELL(A, i, &nz, &cols, &vals));
2299566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, nz, cols, vals, INSERT_VALUES));
2309566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow_SeqSELL(A, i, &nz, &cols, &vals));
231d4002b98SHong Zhang   }
232e3f1f374SStefano Zampini 
2339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
235d4002b98SHong Zhang   B->rmap->bs = A->rmap->bs;
236d4002b98SHong Zhang 
237d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
2389566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
239d4002b98SHong Zhang   } else {
240d4002b98SHong Zhang     *newmat = B;
241d4002b98SHong Zhang   }
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
243d4002b98SHong Zhang }
244d4002b98SHong Zhang 
245d4002b98SHong Zhang #include <../src/mat/impls/aij/seq/aij.h>
246d4002b98SHong Zhang 
247d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
248d71ae5a4SJacob Faibussowitsch {
249d4002b98SHong Zhang   Mat                B;
250d4002b98SHong Zhang   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *)A->data;
251d4002b98SHong Zhang   PetscInt          *ai = a->i, m = A->rmap->N, n = A->cmap->N, i, *rowlengths, row, ncols;
252d4002b98SHong Zhang   const PetscInt    *cols;
253d4002b98SHong Zhang   const PetscScalar *vals;
254d4002b98SHong Zhang 
255d4002b98SHong Zhang   PetscFunctionBegin;
256ad013a7bSRichard Tran Mills 
257ad013a7bSRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
258ad013a7bSRichard Tran Mills     B = *newmat;
259ad013a7bSRichard Tran Mills   } else {
260d5e5b2e5SBarry Smith     if (PetscDefined(USE_DEBUG) || !a->ilen) {
2619566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m, &rowlengths));
262ad540459SPierre Jolivet       for (i = 0; i < m; i++) rowlengths[i] = ai[i + 1] - ai[i];
263d5e5b2e5SBarry Smith     }
264d5e5b2e5SBarry Smith     if (PetscDefined(USE_DEBUG) && a->ilen) {
265d5e5b2e5SBarry Smith       PetscBool eq;
2669566063dSJacob Faibussowitsch       PetscCall(PetscMemcmp(rowlengths, a->ilen, m * sizeof(PetscInt), &eq));
26728b400f6SJacob Faibussowitsch       PetscCheck(eq, PETSC_COMM_SELF, PETSC_ERR_PLIB, "SeqAIJ ilen array incorrect");
2689566063dSJacob Faibussowitsch       PetscCall(PetscFree(rowlengths));
269d5e5b2e5SBarry Smith       rowlengths = a->ilen;
270d5e5b2e5SBarry Smith     } else if (a->ilen) rowlengths = a->ilen;
2719566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2729566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
2739566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQSELL));
2749566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(B, 0, rowlengths));
2759566063dSJacob Faibussowitsch     if (rowlengths != a->ilen) PetscCall(PetscFree(rowlengths));
276ad013a7bSRichard Tran Mills   }
277d4002b98SHong Zhang 
278d4002b98SHong Zhang   for (row = 0; row < m; row++) {
2799566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqAIJ(A, row, &ncols, (PetscInt **)&cols, (PetscScalar **)&vals));
2809566063dSJacob Faibussowitsch     PetscCall(MatSetValues_SeqSELL(B, 1, &row, ncols, cols, vals, INSERT_VALUES));
2819566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow_SeqAIJ(A, row, &ncols, (PetscInt **)&cols, (PetscScalar **)&vals));
282d4002b98SHong Zhang   }
2839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
285d4002b98SHong Zhang   B->rmap->bs = A->rmap->bs;
286d4002b98SHong Zhang 
287d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
2889566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
289d4002b98SHong Zhang   } else {
290d4002b98SHong Zhang     *newmat = B;
291d4002b98SHong Zhang   }
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
293d4002b98SHong Zhang }
294d4002b98SHong Zhang 
295d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_SeqSELL(Mat A, Vec xx, Vec yy)
296d71ae5a4SJacob Faibussowitsch {
297d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
298d4002b98SHong Zhang   PetscScalar       *y;
299d4002b98SHong Zhang   const PetscScalar *x;
300d4002b98SHong Zhang   const MatScalar   *aval        = a->val;
301d4002b98SHong Zhang   PetscInt           totalslices = a->totalslices;
302d4002b98SHong Zhang   const PetscInt    *acolidx     = a->colidx;
3037285fed1SHong Zhang   PetscInt           i, j;
304d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
305d4002b98SHong Zhang   __m512d  vec_x, vec_y, vec_vals;
306d4002b98SHong Zhang   __m256i  vec_idx;
307d4002b98SHong Zhang   __mmask8 mask;
308d4002b98SHong Zhang   __m512d  vec_x2, vec_y2, vec_vals2, vec_x3, vec_y3, vec_vals3, vec_x4, vec_y4, vec_vals4;
309d4002b98SHong Zhang   __m256i  vec_idx2, vec_idx3, vec_idx4;
3105f70456aSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX2__) && defined(__FMA__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
311a48a6482SHong Zhang   __m128i   vec_idx;
312a48a6482SHong Zhang   __m256d   vec_x, vec_y, vec_y2, vec_vals;
313a48a6482SHong Zhang   MatScalar yval;
314a48a6482SHong Zhang   PetscInt  r, rows_left, row, nnz_in_row;
31521cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
316d4002b98SHong Zhang   __m128d   vec_x_tmp;
317d4002b98SHong Zhang   __m256d   vec_x, vec_y, vec_y2, vec_vals;
318d4002b98SHong Zhang   MatScalar yval;
319d4002b98SHong Zhang   PetscInt  r, rows_left, row, nnz_in_row;
320d4002b98SHong Zhang #else
321*07e43b41SHong Zhang   PetscInt     k, sliceheight = a->sliceheight;
322*07e43b41SHong Zhang   PetscScalar *sum;
323d4002b98SHong Zhang #endif
324d4002b98SHong Zhang 
325d4002b98SHong Zhang #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
326d4002b98SHong Zhang   #pragma disjoint(*x, *y, *aval)
327d4002b98SHong Zhang #endif
328d4002b98SHong Zhang 
329d4002b98SHong Zhang   PetscFunctionBegin;
3309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx, &x));
3319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(yy, &y));
332d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
333*07e43b41SHong Zhang   PetscCheck(a->sliceheight == 8, PETSC_COMM_SELF, PETSC_ERR_SUP, "The kernel requires a slice height of 8, but the input matrix has a slice height of %" PetscInt_FMT, a->sliceheight);
334d4002b98SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
335d4002b98SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
336d4002b98SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
337d4002b98SHong Zhang 
338d4002b98SHong Zhang     vec_y  = _mm512_setzero_pd();
339d4002b98SHong Zhang     vec_y2 = _mm512_setzero_pd();
340d4002b98SHong Zhang     vec_y3 = _mm512_setzero_pd();
341d4002b98SHong Zhang     vec_y4 = _mm512_setzero_pd();
342d4002b98SHong Zhang 
343da81f932SPierre Jolivet     j = a->sliidx[i] >> 3; /* 8 bytes are read at each time, corresponding to a slice column */
344d4002b98SHong Zhang     switch ((a->sliidx[i + 1] - a->sliidx[i]) / 8 & 3) {
345d4002b98SHong Zhang     case 3:
346d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3479371c9d4SSatish Balay       acolidx += 8;
3489371c9d4SSatish Balay       aval += 8;
349d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
3509371c9d4SSatish Balay       acolidx += 8;
3519371c9d4SSatish Balay       aval += 8;
352d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
3539371c9d4SSatish Balay       acolidx += 8;
3549371c9d4SSatish Balay       aval += 8;
355d4002b98SHong Zhang       j += 3;
356d4002b98SHong Zhang       break;
357d4002b98SHong Zhang     case 2:
358d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3599371c9d4SSatish Balay       acolidx += 8;
3609371c9d4SSatish Balay       aval += 8;
361d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
3629371c9d4SSatish Balay       acolidx += 8;
3639371c9d4SSatish Balay       aval += 8;
364d4002b98SHong Zhang       j += 2;
365d4002b98SHong Zhang       break;
366d4002b98SHong Zhang     case 1:
367d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3689371c9d4SSatish Balay       acolidx += 8;
3699371c9d4SSatish Balay       aval += 8;
370d4002b98SHong Zhang       j += 1;
371d4002b98SHong Zhang       break;
372d4002b98SHong Zhang     }
373d4002b98SHong Zhang   #pragma novector
374d4002b98SHong Zhang     for (; j < (a->sliidx[i + 1] >> 3); j += 4) {
375d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3769371c9d4SSatish Balay       acolidx += 8;
3779371c9d4SSatish Balay       aval += 8;
378d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
3799371c9d4SSatish Balay       acolidx += 8;
3809371c9d4SSatish Balay       aval += 8;
381d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
3829371c9d4SSatish Balay       acolidx += 8;
3839371c9d4SSatish Balay       aval += 8;
384d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx4, vec_x4, vec_vals4, vec_y4);
3859371c9d4SSatish Balay       acolidx += 8;
3869371c9d4SSatish Balay       aval += 8;
387d4002b98SHong Zhang     }
388d4002b98SHong Zhang 
389d4002b98SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y2);
390d4002b98SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y3);
391d4002b98SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y4);
392d4002b98SHong Zhang     if (i == totalslices - 1 && A->rmap->n & 0x07) { /* if last slice has padding rows */
393d4002b98SHong Zhang       mask = (__mmask8)(0xff >> (8 - (A->rmap->n & 0x07)));
394ef588d5cSRichard Tran Mills       _mm512_mask_storeu_pd(&y[8 * i], mask, vec_y);
395d4002b98SHong Zhang     } else {
396ef588d5cSRichard Tran Mills       _mm512_storeu_pd(&y[8 * i], vec_y);
397d4002b98SHong Zhang     }
398d4002b98SHong Zhang   }
3995f70456aSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX2__) && defined(__FMA__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
400*07e43b41SHong Zhang   PetscCheck(a->sliceheight == 8, PETSC_COMM_SELF, PETSC_ERR_SUP, "The kernel requires a slice height of 8, but the input matrix has a slice height of %" PetscInt_FMT, a->sliceheight);
401a48a6482SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over full slices */
402a48a6482SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
403a48a6482SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
404a48a6482SHong Zhang 
405a48a6482SHong Zhang     /* last slice may have padding rows. Don't use vectorization. */
406a48a6482SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
407a48a6482SHong Zhang       rows_left = A->rmap->n - 8 * i;
408a48a6482SHong Zhang       for (r = 0; r < rows_left; ++r) {
409a48a6482SHong Zhang         yval       = (MatScalar)0;
410a48a6482SHong Zhang         row        = 8 * i + r;
411a48a6482SHong Zhang         nnz_in_row = a->rlen[row];
412a48a6482SHong Zhang         for (j = 0; j < nnz_in_row; ++j) yval += aval[8 * j + r] * x[acolidx[8 * j + r]];
413a48a6482SHong Zhang         y[row] = yval;
414a48a6482SHong Zhang       }
415a48a6482SHong Zhang       break;
416a48a6482SHong Zhang     }
417a48a6482SHong Zhang 
418a48a6482SHong Zhang     vec_y  = _mm256_setzero_pd();
419a48a6482SHong Zhang     vec_y2 = _mm256_setzero_pd();
420a48a6482SHong Zhang 
421a48a6482SHong Zhang   /* Process slice of height 8 (512 bits) via two subslices of height 4 (256 bits) via AVX */
422a48a6482SHong Zhang   #pragma novector
423a48a6482SHong Zhang   #pragma unroll(2)
424a48a6482SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
425a48a6482SHong Zhang       AVX2_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
4269371c9d4SSatish Balay       aval += 4;
4279371c9d4SSatish Balay       acolidx += 4;
428a48a6482SHong Zhang       AVX2_Mult_Private(vec_idx, vec_x, vec_vals, vec_y2);
4299371c9d4SSatish Balay       aval += 4;
4309371c9d4SSatish Balay       acolidx += 4;
431a48a6482SHong Zhang     }
432a48a6482SHong Zhang 
433ef588d5cSRichard Tran Mills     _mm256_storeu_pd(y + i * 8, vec_y);
434ef588d5cSRichard Tran Mills     _mm256_storeu_pd(y + i * 8 + 4, vec_y2);
435a48a6482SHong Zhang   }
43621cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
437*07e43b41SHong Zhang   PetscCheck(a->sliceheight == 8, PETSC_COMM_SELF, PETSC_ERR_SUP, "The kernel requires a slice height of 8, but the input matrix has a slice height of %" PetscInt_FMT, a->sliceheight);
438d4002b98SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over full slices */
439d4002b98SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
440d4002b98SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
441d4002b98SHong Zhang 
442d4002b98SHong Zhang     vec_y  = _mm256_setzero_pd();
443d4002b98SHong Zhang     vec_y2 = _mm256_setzero_pd();
444d4002b98SHong Zhang 
445d4002b98SHong Zhang     /* last slice may have padding rows. Don't use vectorization. */
446d4002b98SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
447d4002b98SHong Zhang       rows_left = A->rmap->n - 8 * i;
448d4002b98SHong Zhang       for (r = 0; r < rows_left; ++r) {
449d4002b98SHong Zhang         yval       = (MatScalar)0;
450d4002b98SHong Zhang         row        = 8 * i + r;
451d4002b98SHong Zhang         nnz_in_row = a->rlen[row];
452d4002b98SHong Zhang         for (j = 0; j < nnz_in_row; ++j) yval += aval[8 * j + r] * x[acolidx[8 * j + r]];
453d4002b98SHong Zhang         y[row] = yval;
454d4002b98SHong Zhang       }
455d4002b98SHong Zhang       break;
456d4002b98SHong Zhang     }
457d4002b98SHong Zhang 
458d4002b98SHong Zhang   /* Process slice of height 8 (512 bits) via two subslices of height 4 (256 bits) via AVX */
459a48a6482SHong Zhang   #pragma novector
460a48a6482SHong Zhang   #pragma unroll(2)
4617285fed1SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
462d4002b98SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
463165f9cc3SJed Brown       vec_x_tmp = _mm_setzero_pd();
464d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
465d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
466d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
467d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
468d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
469d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
470d4002b98SHong Zhang       vec_y     = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y);
471d4002b98SHong Zhang       aval += 4;
472d4002b98SHong Zhang 
473d4002b98SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
474d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
475d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
476d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
477d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
478d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
479d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
480d4002b98SHong Zhang       vec_y2    = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y2);
481d4002b98SHong Zhang       aval += 4;
482d4002b98SHong Zhang     }
483d4002b98SHong Zhang 
484d4002b98SHong Zhang     _mm256_storeu_pd(y + i * 8, vec_y);
485d4002b98SHong Zhang     _mm256_storeu_pd(y + i * 8 + 4, vec_y2);
486d4002b98SHong Zhang   }
487d4002b98SHong Zhang #else
488*07e43b41SHong Zhang   PetscCall(PetscMalloc1(sliceheight, &sum));
489d4002b98SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
490*07e43b41SHong Zhang     for (j = 0; j < sliceheight; j++) {
4912d1451d4SHong Zhang       sum[j] = 0.0;
492*07e43b41SHong Zhang       for (k = a->sliidx[i] + j; k < a->sliidx[i + 1]; k += sliceheight) sum[j] += aval[k] * x[acolidx[k]];
493d4002b98SHong Zhang     }
494*07e43b41SHong Zhang     if (i == totalslices - 1 && (A->rmap->n % sliceheight)) { /* if last slice has padding rows */
495*07e43b41SHong Zhang       for (j = 0; j < (A->rmap->n % sliceheight); j++) y[sliceheight * i + j] = sum[j];
496d4002b98SHong Zhang     } else {
497*07e43b41SHong Zhang       for (j = 0; j < sliceheight; j++) y[sliceheight * i + j] = sum[j];
498d4002b98SHong Zhang     }
499d4002b98SHong Zhang   }
500*07e43b41SHong Zhang   PetscCall(PetscFree(sum));
501d4002b98SHong Zhang #endif
502d4002b98SHong Zhang 
5039566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->nz - a->nonzerorowcnt)); /* theoretical minimal FLOPs */
5049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx, &x));
5059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(yy, &y));
5063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
507d4002b98SHong Zhang }
508d4002b98SHong Zhang 
509d4002b98SHong Zhang #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h>
510d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAdd_SeqSELL(Mat A, Vec xx, Vec yy, Vec zz)
511d71ae5a4SJacob Faibussowitsch {
512d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
513d4002b98SHong Zhang   PetscScalar       *y, *z;
514d4002b98SHong Zhang   const PetscScalar *x;
515d4002b98SHong Zhang   const MatScalar   *aval        = a->val;
516d4002b98SHong Zhang   PetscInt           totalslices = a->totalslices;
517d4002b98SHong Zhang   const PetscInt    *acolidx     = a->colidx;
518d4002b98SHong Zhang   PetscInt           i, j;
519d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
5207285fed1SHong Zhang   __m512d  vec_x, vec_y, vec_vals;
521d4002b98SHong Zhang   __m256i  vec_idx;
522d4002b98SHong Zhang   __mmask8 mask;
5237285fed1SHong Zhang   __m512d  vec_x2, vec_y2, vec_vals2, vec_x3, vec_y3, vec_vals3, vec_x4, vec_y4, vec_vals4;
5247285fed1SHong Zhang   __m256i  vec_idx2, vec_idx3, vec_idx4;
52521cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
5267285fed1SHong Zhang   __m128d   vec_x_tmp;
5277285fed1SHong Zhang   __m256d   vec_x, vec_y, vec_y2, vec_vals;
5287285fed1SHong Zhang   MatScalar yval;
5297285fed1SHong Zhang   PetscInt  r, row, nnz_in_row;
530d4002b98SHong Zhang #else
531*07e43b41SHong Zhang   PetscInt     k, sliceheight = a->sliceheight;
532*07e43b41SHong Zhang   PetscScalar *sum;
533d4002b98SHong Zhang #endif
534d4002b98SHong Zhang 
535d4002b98SHong Zhang #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
536d4002b98SHong Zhang   #pragma disjoint(*x, *y, *aval)
537d4002b98SHong Zhang #endif
538d4002b98SHong Zhang 
539d4002b98SHong Zhang   PetscFunctionBegin;
5402d1451d4SHong Zhang   if (!a->nz) {
5412d1451d4SHong Zhang     PetscCall(VecCopy(yy, zz));
5422d1451d4SHong Zhang     PetscFunctionReturn(PETSC_SUCCESS);
5432d1451d4SHong Zhang   }
5449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx, &x));
5459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayPair(yy, zz, &y, &z));
546d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
547*07e43b41SHong Zhang   PetscCheck(a->sliceheight == 8, PETSC_COMM_SELF, PETSC_ERR_SUP, "The kernel requires a slice height of 8, but the input matrix has a slice height of %" PetscInt_FMT, a->sliceheight);
5487285fed1SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
5497285fed1SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
5507285fed1SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
5517285fed1SHong Zhang 
552d4002b98SHong Zhang     if (i == totalslices - 1 && A->rmap->n & 0x07) { /* if last slice has padding rows */
553d4002b98SHong Zhang       mask  = (__mmask8)(0xff >> (8 - (A->rmap->n & 0x07)));
554ef588d5cSRichard Tran Mills       vec_y = _mm512_mask_loadu_pd(vec_y, mask, &y[8 * i]);
5557285fed1SHong Zhang     } else {
556ef588d5cSRichard Tran Mills       vec_y = _mm512_loadu_pd(&y[8 * i]);
5577285fed1SHong Zhang     }
5587285fed1SHong Zhang     vec_y2 = _mm512_setzero_pd();
5597285fed1SHong Zhang     vec_y3 = _mm512_setzero_pd();
5607285fed1SHong Zhang     vec_y4 = _mm512_setzero_pd();
5617285fed1SHong Zhang 
562da81f932SPierre Jolivet     j = a->sliidx[i] >> 3; /* 8 bytes are read at each time, corresponding to a slice column */
5637285fed1SHong Zhang     switch ((a->sliidx[i + 1] - a->sliidx[i]) / 8 & 3) {
5647285fed1SHong Zhang     case 3:
5657285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
5669371c9d4SSatish Balay       acolidx += 8;
5679371c9d4SSatish Balay       aval += 8;
5687285fed1SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
5699371c9d4SSatish Balay       acolidx += 8;
5709371c9d4SSatish Balay       aval += 8;
5717285fed1SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
5729371c9d4SSatish Balay       acolidx += 8;
5739371c9d4SSatish Balay       aval += 8;
5747285fed1SHong Zhang       j += 3;
5757285fed1SHong Zhang       break;
5767285fed1SHong Zhang     case 2:
5777285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
5789371c9d4SSatish Balay       acolidx += 8;
5799371c9d4SSatish Balay       aval += 8;
5807285fed1SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
5819371c9d4SSatish Balay       acolidx += 8;
5829371c9d4SSatish Balay       aval += 8;
5837285fed1SHong Zhang       j += 2;
5847285fed1SHong Zhang       break;
5857285fed1SHong Zhang     case 1:
5867285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
5879371c9d4SSatish Balay       acolidx += 8;
5889371c9d4SSatish Balay       aval += 8;
5897285fed1SHong Zhang       j += 1;
5907285fed1SHong Zhang       break;
5917285fed1SHong Zhang     }
5927285fed1SHong Zhang   #pragma novector
5937285fed1SHong Zhang     for (; j < (a->sliidx[i + 1] >> 3); j += 4) {
5947285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
5959371c9d4SSatish Balay       acolidx += 8;
5969371c9d4SSatish Balay       aval += 8;
5977285fed1SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
5989371c9d4SSatish Balay       acolidx += 8;
5999371c9d4SSatish Balay       aval += 8;
6007285fed1SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
6019371c9d4SSatish Balay       acolidx += 8;
6029371c9d4SSatish Balay       aval += 8;
6037285fed1SHong Zhang       AVX512_Mult_Private(vec_idx4, vec_x4, vec_vals4, vec_y4);
6049371c9d4SSatish Balay       acolidx += 8;
6059371c9d4SSatish Balay       aval += 8;
6067285fed1SHong Zhang     }
6077285fed1SHong Zhang 
6087285fed1SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y2);
6097285fed1SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y3);
6107285fed1SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y4);
6117285fed1SHong Zhang     if (i == totalslices - 1 && A->rmap->n & 0x07) { /* if last slice has padding rows */
612ef588d5cSRichard Tran Mills       _mm512_mask_storeu_pd(&z[8 * i], mask, vec_y);
613d4002b98SHong Zhang     } else {
614ef588d5cSRichard Tran Mills       _mm512_storeu_pd(&z[8 * i], vec_y);
615d4002b98SHong Zhang     }
6167285fed1SHong Zhang   }
61721cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
618*07e43b41SHong Zhang   PetscCheck(a->sliceheight == 8, PETSC_COMM_SELF, PETSC_ERR_SUP, "The kernel requires a slice height of 8, but the input matrix has a slice height of %" PetscInt_FMT, a->sliceheight);
6197285fed1SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over full slices */
6207285fed1SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
6217285fed1SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
6227285fed1SHong Zhang 
6237285fed1SHong Zhang     /* last slice may have padding rows. Don't use vectorization. */
6247285fed1SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
6257285fed1SHong Zhang       for (r = 0; r < (A->rmap->n & 0x07); ++r) {
6267285fed1SHong Zhang         row        = 8 * i + r;
6277285fed1SHong Zhang         yval       = (MatScalar)0.0;
6287285fed1SHong Zhang         nnz_in_row = a->rlen[row];
6297285fed1SHong Zhang         for (j = 0; j < nnz_in_row; ++j) yval += aval[8 * j + r] * x[acolidx[8 * j + r]];
6307285fed1SHong Zhang         z[row] = y[row] + yval;
6317285fed1SHong Zhang       }
6327285fed1SHong Zhang       break;
6337285fed1SHong Zhang     }
6347285fed1SHong Zhang 
6357285fed1SHong Zhang     vec_y  = _mm256_loadu_pd(y + 8 * i);
6367285fed1SHong Zhang     vec_y2 = _mm256_loadu_pd(y + 8 * i + 4);
6377285fed1SHong Zhang 
6387285fed1SHong Zhang     /* Process slice of height 8 (512 bits) via two subslices of height 4 (256 bits) via AVX */
6397285fed1SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
6407285fed1SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
641165f9cc3SJed Brown       vec_x_tmp = _mm_setzero_pd();
6427285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6437285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
644165f9cc3SJed Brown       vec_x     = _mm256_setzero_pd();
6457285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
6467285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6477285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
6487285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
6497285fed1SHong Zhang       vec_y     = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y);
6507285fed1SHong Zhang       aval += 4;
6517285fed1SHong Zhang 
6527285fed1SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
6537285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6547285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
6557285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
6567285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6577285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
6587285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
6597285fed1SHong Zhang       vec_y2    = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y2);
6607285fed1SHong Zhang       aval += 4;
6617285fed1SHong Zhang     }
6627285fed1SHong Zhang 
6637285fed1SHong Zhang     _mm256_storeu_pd(z + i * 8, vec_y);
6647285fed1SHong Zhang     _mm256_storeu_pd(z + i * 8 + 4, vec_y2);
6657285fed1SHong Zhang   }
666d4002b98SHong Zhang #else
667*07e43b41SHong Zhang   PetscCall(PetscMalloc1(sliceheight, &sum));
6687285fed1SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
669*07e43b41SHong Zhang     for (j = 0; j < sliceheight; j++) {
6702d1451d4SHong Zhang       sum[j] = 0.0;
671*07e43b41SHong Zhang       for (k = a->sliidx[i] + j; k < a->sliidx[i + 1]; k += sliceheight) sum[j] += aval[k] * x[acolidx[k]];
672d4002b98SHong Zhang     }
673*07e43b41SHong Zhang     if (i == totalslices - 1 && (A->rmap->n % sliceheight)) {
674*07e43b41SHong Zhang       for (j = 0; j < (A->rmap->n % sliceheight); j++) z[sliceheight * i + j] = y[sliceheight * i + j] + sum[j];
675d4002b98SHong Zhang     } else {
676*07e43b41SHong Zhang       for (j = 0; j < sliceheight; j++) z[sliceheight * i + j] = y[sliceheight * i + j] + sum[j];
6777285fed1SHong Zhang     }
678d4002b98SHong Zhang   }
679*07e43b41SHong Zhang   PetscCall(PetscFree(sum));
680d4002b98SHong Zhang #endif
681d4002b98SHong Zhang 
6829566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->nz));
6839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx, &x));
6849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayPair(yy, zz, &y, &z));
6853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
686d4002b98SHong Zhang }
687d4002b98SHong Zhang 
688d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAdd_SeqSELL(Mat A, Vec xx, Vec zz, Vec yy)
689d71ae5a4SJacob Faibussowitsch {
690d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
691d4002b98SHong Zhang   PetscScalar       *y;
692d4002b98SHong Zhang   const PetscScalar *x;
693d4002b98SHong Zhang   const MatScalar   *aval    = a->val;
694d4002b98SHong Zhang   const PetscInt    *acolidx = a->colidx;
695*07e43b41SHong Zhang   PetscInt           i, j, r, row, nnz_in_row, totalslices = a->totalslices, sliceheight = a->sliceheight;
696d4002b98SHong Zhang 
697d4002b98SHong Zhang #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
698d4002b98SHong Zhang   #pragma disjoint(*x, *y, *aval)
699d4002b98SHong Zhang #endif
700d4002b98SHong Zhang 
701d4002b98SHong Zhang   PetscFunctionBegin;
702b94d7dedSBarry Smith   if (A->symmetric == PETSC_BOOL3_TRUE) {
7039566063dSJacob Faibussowitsch     PetscCall(MatMultAdd_SeqSELL(A, xx, zz, yy));
7043ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
7059fc32365SStefano Zampini   }
7069566063dSJacob Faibussowitsch   if (zz != yy) PetscCall(VecCopy(zz, yy));
7072d1451d4SHong Zhang 
7082d1451d4SHong Zhang   if (a->nz) {
7099566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(xx, &x));
7109566063dSJacob Faibussowitsch     PetscCall(VecGetArray(yy, &y));
711d4002b98SHong Zhang     for (i = 0; i < a->totalslices; i++) { /* loop over slices */
712*07e43b41SHong Zhang       if (i == totalslices - 1 && (A->rmap->n % sliceheight)) {
713*07e43b41SHong Zhang         for (r = 0; r < (A->rmap->n % sliceheight); ++r) {
714*07e43b41SHong Zhang           row        = sliceheight * i + r;
7157285fed1SHong Zhang           nnz_in_row = a->rlen[row];
716*07e43b41SHong Zhang           for (j = 0; j < nnz_in_row; ++j) y[acolidx[sliceheight * j + r]] += aval[sliceheight * j + r] * x[row];
7177285fed1SHong Zhang         }
7187285fed1SHong Zhang         break;
7197285fed1SHong Zhang       }
720*07e43b41SHong Zhang       for (r = 0; r < sliceheight; ++r)
721*07e43b41SHong Zhang         for (j = a->sliidx[i] + r; j < a->sliidx[i + 1]; j += sliceheight) y[acolidx[j]] += aval[j] * x[sliceheight * i + r];
722d4002b98SHong Zhang     }
7232d1451d4SHong Zhang     PetscCall(PetscLogFlops(2.0 * a->nz));
7249566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(xx, &x));
7259566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(yy, &y));
7262d1451d4SHong Zhang   }
7273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
728d4002b98SHong Zhang }
729d4002b98SHong Zhang 
730d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTranspose_SeqSELL(Mat A, Vec xx, Vec yy)
731d71ae5a4SJacob Faibussowitsch {
732d4002b98SHong Zhang   PetscFunctionBegin;
733b94d7dedSBarry Smith   if (A->symmetric == PETSC_BOOL3_TRUE) {
7349566063dSJacob Faibussowitsch     PetscCall(MatMult_SeqSELL(A, xx, yy));
7359fc32365SStefano Zampini   } else {
7369566063dSJacob Faibussowitsch     PetscCall(VecSet(yy, 0.0));
7379566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd_SeqSELL(A, xx, yy, yy));
7389fc32365SStefano Zampini   }
7393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
740d4002b98SHong Zhang }
741d4002b98SHong Zhang 
742d4002b98SHong Zhang /*
743d4002b98SHong Zhang      Checks for missing diagonals
744d4002b98SHong Zhang */
745d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMissingDiagonal_SeqSELL(Mat A, PetscBool *missing, PetscInt *d)
746d71ae5a4SJacob Faibussowitsch {
747d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
748d4002b98SHong Zhang   PetscInt    *diag, i;
749d4002b98SHong Zhang 
750d4002b98SHong Zhang   PetscFunctionBegin;
751d4002b98SHong Zhang   *missing = PETSC_FALSE;
752d4002b98SHong Zhang   if (A->rmap->n > 0 && !(a->colidx)) {
753d4002b98SHong Zhang     *missing = PETSC_TRUE;
754d4002b98SHong Zhang     if (d) *d = 0;
7559566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Matrix has no entries therefore is missing diagonal\n"));
756d4002b98SHong Zhang   } else {
757d4002b98SHong Zhang     diag = a->diag;
758d4002b98SHong Zhang     for (i = 0; i < A->rmap->n; i++) {
759d4002b98SHong Zhang       if (diag[i] == -1) {
760d4002b98SHong Zhang         *missing = PETSC_TRUE;
761d4002b98SHong Zhang         if (d) *d = i;
7629566063dSJacob Faibussowitsch         PetscCall(PetscInfo(A, "Matrix is missing diagonal number %" PetscInt_FMT "\n", i));
763d4002b98SHong Zhang         break;
764d4002b98SHong Zhang       }
765d4002b98SHong Zhang     }
766d4002b98SHong Zhang   }
7673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
768d4002b98SHong Zhang }
769d4002b98SHong Zhang 
770d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMarkDiagonal_SeqSELL(Mat A)
771d71ae5a4SJacob Faibussowitsch {
772d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
773d4002b98SHong Zhang   PetscInt     i, j, m = A->rmap->n, shift;
774d4002b98SHong Zhang 
775d4002b98SHong Zhang   PetscFunctionBegin;
776d4002b98SHong Zhang   if (!a->diag) {
7779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &a->diag));
778d4002b98SHong Zhang     a->free_diag = PETSC_TRUE;
779d4002b98SHong Zhang   }
780d4002b98SHong Zhang   for (i = 0; i < m; i++) {                                          /* loop over rows */
781*07e43b41SHong Zhang     shift      = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
782d4002b98SHong Zhang     a->diag[i] = -1;
783d4002b98SHong Zhang     for (j = 0; j < a->rlen[i]; j++) {
784*07e43b41SHong Zhang       if (a->colidx[shift + a->sliceheight * j] == i) {
785*07e43b41SHong Zhang         a->diag[i] = shift + a->sliceheight * j;
786d4002b98SHong Zhang         break;
787d4002b98SHong Zhang       }
788d4002b98SHong Zhang     }
789d4002b98SHong Zhang   }
7903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
791d4002b98SHong Zhang }
792d4002b98SHong Zhang 
793d4002b98SHong Zhang /*
794d4002b98SHong Zhang   Negative shift indicates do not generate an error if there is a zero diagonal, just invert it anyways
795d4002b98SHong Zhang */
796d71ae5a4SJacob Faibussowitsch PetscErrorCode MatInvertDiagonal_SeqSELL(Mat A, PetscScalar omega, PetscScalar fshift)
797d71ae5a4SJacob Faibussowitsch {
798d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
799d4002b98SHong Zhang   PetscInt     i, *diag, m = A->rmap->n;
800d4002b98SHong Zhang   MatScalar   *val = a->val;
801d4002b98SHong Zhang   PetscScalar *idiag, *mdiag;
802d4002b98SHong Zhang 
803d4002b98SHong Zhang   PetscFunctionBegin;
8043ba16761SJacob Faibussowitsch   if (a->idiagvalid) PetscFunctionReturn(PETSC_SUCCESS);
8059566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_SeqSELL(A));
806d4002b98SHong Zhang   diag = a->diag;
807d4002b98SHong Zhang   if (!a->idiag) {
8089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(m, &a->idiag, m, &a->mdiag, m, &a->ssor_work));
809d4002b98SHong Zhang     val = a->val;
810d4002b98SHong Zhang   }
811d4002b98SHong Zhang   mdiag = a->mdiag;
812d4002b98SHong Zhang   idiag = a->idiag;
813d4002b98SHong Zhang 
814d4002b98SHong Zhang   if (omega == 1.0 && PetscRealPart(fshift) <= 0.0) {
815d4002b98SHong Zhang     for (i = 0; i < m; i++) {
816d4002b98SHong Zhang       mdiag[i] = val[diag[i]];
817d4002b98SHong Zhang       if (!PetscAbsScalar(mdiag[i])) { /* zero diagonal */
8180fdf79fbSJacob Faibussowitsch         PetscCheck(PetscRealPart(fshift), PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Zero diagonal on row %" PetscInt_FMT, i);
8199566063dSJacob Faibussowitsch         PetscCall(PetscInfo(A, "Zero diagonal on row %" PetscInt_FMT "\n", i));
820d4002b98SHong Zhang         A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
821d4002b98SHong Zhang         A->factorerror_zeropivot_value = 0.0;
822d4002b98SHong Zhang         A->factorerror_zeropivot_row   = i;
823d4002b98SHong Zhang       }
824d4002b98SHong Zhang       idiag[i] = 1.0 / val[diag[i]];
825d4002b98SHong Zhang     }
8269566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(m));
827d4002b98SHong Zhang   } else {
828d4002b98SHong Zhang     for (i = 0; i < m; i++) {
829d4002b98SHong Zhang       mdiag[i] = val[diag[i]];
830d4002b98SHong Zhang       idiag[i] = omega / (fshift + val[diag[i]]);
831d4002b98SHong Zhang     }
8329566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(2.0 * m));
833d4002b98SHong Zhang   }
834d4002b98SHong Zhang   a->idiagvalid = PETSC_TRUE;
8353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
836d4002b98SHong Zhang }
837d4002b98SHong Zhang 
838d71ae5a4SJacob Faibussowitsch PetscErrorCode MatZeroEntries_SeqSELL(Mat A)
839d71ae5a4SJacob Faibussowitsch {
840d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
841d4002b98SHong Zhang 
842d4002b98SHong Zhang   PetscFunctionBegin;
8439566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(a->val, a->sliidx[a->totalslices]));
8449566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
8453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
846d4002b98SHong Zhang }
847d4002b98SHong Zhang 
848d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_SeqSELL(Mat A)
849d71ae5a4SJacob Faibussowitsch {
850d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
851d4002b98SHong Zhang 
852d4002b98SHong Zhang   PetscFunctionBegin;
853d4002b98SHong Zhang #if defined(PETSC_USE_LOG)
8543ba16761SJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)A, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, A->rmap->n, A->cmap->n, a->nz));
855d4002b98SHong Zhang #endif
8569566063dSJacob Faibussowitsch   PetscCall(MatSeqXSELLFreeSELL(A, &a->val, &a->colidx));
8579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->row));
8589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->col));
8599566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->diag));
8609566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->rlen));
8619566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->sliidx));
8629566063dSJacob Faibussowitsch   PetscCall(PetscFree3(a->idiag, a->mdiag, a->ssor_work));
8639566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->solve_work));
8649566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->icol));
8659566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->saved_values));
8669566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a->getrowcols, a->getrowvals));
8679566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
868d4002b98SHong Zhang 
8699566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
8709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatStoreValues_C", NULL));
8719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatRetrieveValues_C", NULL));
8729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLSetPreallocation_C", NULL));
8732e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetArray_C", NULL));
8742e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLRestoreArray_C", NULL));
8752d1451d4SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_SeqSELL_SeqAIJ_C", NULL));
8762d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
8772d1451d4SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_SeqSELL_SeqSELLCUDA_C", NULL));
8782d1451d4SHong Zhang #endif
879*07e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetFillRatio_C", NULL));
880*07e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetMaxSliceWidth_C", NULL));
881*07e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetAvgSliceWidth_C", NULL));
882*07e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLSetSliceHeight_C", NULL));
8833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
884d4002b98SHong Zhang }
885d4002b98SHong Zhang 
886d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetOption_SeqSELL(Mat A, MatOption op, PetscBool flg)
887d71ae5a4SJacob Faibussowitsch {
888d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
889d4002b98SHong Zhang 
890d4002b98SHong Zhang   PetscFunctionBegin;
891d4002b98SHong Zhang   switch (op) {
892d71ae5a4SJacob Faibussowitsch   case MAT_ROW_ORIENTED:
893d71ae5a4SJacob Faibussowitsch     a->roworiented = flg;
894d71ae5a4SJacob Faibussowitsch     break;
895d71ae5a4SJacob Faibussowitsch   case MAT_KEEP_NONZERO_PATTERN:
896d71ae5a4SJacob Faibussowitsch     a->keepnonzeropattern = flg;
897d71ae5a4SJacob Faibussowitsch     break;
898d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_LOCATIONS:
899d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? 0 : 1);
900d71ae5a4SJacob Faibussowitsch     break;
901d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_LOCATION_ERR:
902d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? -1 : 0);
903d71ae5a4SJacob Faibussowitsch     break;
904d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_ALLOCATION_ERR:
905d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? -2 : 0);
906d71ae5a4SJacob Faibussowitsch     break;
907d71ae5a4SJacob Faibussowitsch   case MAT_UNUSED_NONZERO_LOCATION_ERR:
908d71ae5a4SJacob Faibussowitsch     a->nounused = (flg ? -1 : 0);
909d71ae5a4SJacob Faibussowitsch     break;
9108c78258cSHong Zhang   case MAT_FORCE_DIAGONAL_ENTRIES:
911d4002b98SHong Zhang   case MAT_IGNORE_OFF_PROC_ENTRIES:
912d4002b98SHong Zhang   case MAT_USE_HASH_TABLE:
913d71ae5a4SJacob Faibussowitsch   case MAT_SORTED_FULL:
914d71ae5a4SJacob Faibussowitsch     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
915d71ae5a4SJacob Faibussowitsch     break;
916d4002b98SHong Zhang   case MAT_SPD:
917d4002b98SHong Zhang   case MAT_SYMMETRIC:
918d4002b98SHong Zhang   case MAT_STRUCTURALLY_SYMMETRIC:
919d4002b98SHong Zhang   case MAT_HERMITIAN:
920d4002b98SHong Zhang   case MAT_SYMMETRY_ETERNAL:
921b94d7dedSBarry Smith   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
922b94d7dedSBarry Smith   case MAT_SPD_ETERNAL:
923d4002b98SHong Zhang     /* These options are handled directly by MatSetOption() */
924d4002b98SHong Zhang     break;
925d71ae5a4SJacob Faibussowitsch   default:
926d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
927d4002b98SHong Zhang   }
9283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
929d4002b98SHong Zhang }
930d4002b98SHong Zhang 
931d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_SeqSELL(Mat A, Vec v)
932d71ae5a4SJacob Faibussowitsch {
933d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
934d4002b98SHong Zhang   PetscInt     i, j, n, shift;
935d4002b98SHong Zhang   PetscScalar *x, zero = 0.0;
936d4002b98SHong Zhang 
937d4002b98SHong Zhang   PetscFunctionBegin;
9389566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
93908401ef6SPierre Jolivet   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
940d4002b98SHong Zhang 
941d4002b98SHong Zhang   if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) {
942d4002b98SHong Zhang     PetscInt *diag = a->diag;
9439566063dSJacob Faibussowitsch     PetscCall(VecGetArray(v, &x));
944d4002b98SHong Zhang     for (i = 0; i < n; i++) x[i] = 1.0 / a->val[diag[i]];
9459566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(v, &x));
9463ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
947d4002b98SHong Zhang   }
948d4002b98SHong Zhang 
9499566063dSJacob Faibussowitsch   PetscCall(VecSet(v, zero));
9509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &x));
951d4002b98SHong Zhang   for (i = 0; i < n; i++) {                                     /* loop over rows */
952*07e43b41SHong Zhang     shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
953d4002b98SHong Zhang     x[i]  = 0;
954d4002b98SHong Zhang     for (j = 0; j < a->rlen[i]; j++) {
955*07e43b41SHong Zhang       if (a->colidx[shift + a->sliceheight * j] == i) {
956*07e43b41SHong Zhang         x[i] = a->val[shift + a->sliceheight * j];
957d4002b98SHong Zhang         break;
958d4002b98SHong Zhang       }
959d4002b98SHong Zhang     }
960d4002b98SHong Zhang   }
9619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &x));
9623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
963d4002b98SHong Zhang }
964d4002b98SHong Zhang 
965d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalScale_SeqSELL(Mat A, Vec ll, Vec rr)
966d71ae5a4SJacob Faibussowitsch {
967d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
968d4002b98SHong Zhang   const PetscScalar *l, *r;
969d4002b98SHong Zhang   PetscInt           i, j, m, n, row;
970d4002b98SHong Zhang 
971d4002b98SHong Zhang   PetscFunctionBegin;
972d4002b98SHong Zhang   if (ll) {
973d4002b98SHong Zhang     /* The local size is used so that VecMPI can be passed to this routine
974d4002b98SHong Zhang        by MatDiagonalScale_MPISELL */
9759566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(ll, &m));
97608401ef6SPierre Jolivet     PetscCheck(m == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Left scaling vector wrong length");
9779566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(ll, &l));
978d4002b98SHong Zhang     for (i = 0; i < a->totalslices; i++) {                            /* loop over slices */
979*07e43b41SHong Zhang       if (i == a->totalslices - 1 && (A->rmap->n % a->sliceheight)) { /* if last slice has padding rows */
980*07e43b41SHong Zhang         for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = (row + 1) % a->sliceheight) {
981*07e43b41SHong Zhang           if (row < (A->rmap->n % a->sliceheight)) a->val[j] *= l[a->sliceheight * i + row];
982dab86139SHong Zhang         }
983dab86139SHong Zhang       } else {
984*07e43b41SHong Zhang         for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = (row + 1) % a->sliceheight) { a->val[j] *= l[a->sliceheight * i + row]; }
985d4002b98SHong Zhang       }
986dab86139SHong Zhang     }
9879566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(ll, &l));
9889566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(a->nz));
989d4002b98SHong Zhang   }
990d4002b98SHong Zhang   if (rr) {
9919566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rr, &n));
99208401ef6SPierre Jolivet     PetscCheck(n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Right scaling vector wrong length");
9939566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(rr, &r));
994d4002b98SHong Zhang     for (i = 0; i < a->totalslices; i++) {                            /* loop over slices */
995*07e43b41SHong Zhang       if (i == a->totalslices - 1 && (A->rmap->n % a->sliceheight)) { /* if last slice has padding rows */
996*07e43b41SHong Zhang         for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = ((row + 1) % a->sliceheight)) {
997*07e43b41SHong Zhang           if (row < (A->rmap->n % a->sliceheight)) a->val[j] *= r[a->colidx[j]];
998dab86139SHong Zhang         }
999dab86139SHong Zhang       } else {
1000ad540459SPierre Jolivet         for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j++) a->val[j] *= r[a->colidx[j]];
1001d4002b98SHong Zhang       }
1002dab86139SHong Zhang     }
10039566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(rr, &r));
10049566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(a->nz));
1005d4002b98SHong Zhang   }
10069566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
10072d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
10082d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
10092d1451d4SHong Zhang #endif
10103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1011d4002b98SHong Zhang }
1012d4002b98SHong Zhang 
1013d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetValues_SeqSELL(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], PetscScalar v[])
1014d71ae5a4SJacob Faibussowitsch {
1015d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1016d4002b98SHong Zhang   PetscInt    *cp, i, k, low, high, t, row, col, l;
1017d4002b98SHong Zhang   PetscInt     shift;
1018d4002b98SHong Zhang   MatScalar   *vp;
1019d4002b98SHong Zhang 
1020d4002b98SHong Zhang   PetscFunctionBegin;
102168aafef3SStefano Zampini   for (k = 0; k < m; k++) { /* loop over requested rows */
1022d4002b98SHong Zhang     row = im[k];
1023d4002b98SHong Zhang     if (row < 0) continue;
10246bdcaf15SBarry Smith     PetscCheck(row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, A->rmap->n - 1);
1025*07e43b41SHong Zhang     shift = a->sliidx[row / a->sliceheight] + (row % a->sliceheight); /* starting index of the row */
1026d4002b98SHong Zhang     cp    = a->colidx + shift;                                        /* pointer to the row */
1027d4002b98SHong Zhang     vp    = a->val + shift;                                           /* pointer to the row */
102868aafef3SStefano Zampini     for (l = 0; l < n; l++) {                                         /* loop over requested columns */
1029d4002b98SHong Zhang       col = in[l];
1030d4002b98SHong Zhang       if (col < 0) continue;
10316bdcaf15SBarry Smith       PetscCheck(col < A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: row %" PetscInt_FMT " max %" PetscInt_FMT, col, A->cmap->n - 1);
10329371c9d4SSatish Balay       high = a->rlen[row];
10339371c9d4SSatish Balay       low  = 0; /* assume unsorted */
1034d4002b98SHong Zhang       while (high - low > 5) {
1035d4002b98SHong Zhang         t = (low + high) / 2;
1036*07e43b41SHong Zhang         if (*(cp + a->sliceheight * t) > col) high = t;
1037d4002b98SHong Zhang         else low = t;
1038d4002b98SHong Zhang       }
1039d4002b98SHong Zhang       for (i = low; i < high; i++) {
1040*07e43b41SHong Zhang         if (*(cp + a->sliceheight * i) > col) break;
1041*07e43b41SHong Zhang         if (*(cp + a->sliceheight * i) == col) {
1042*07e43b41SHong Zhang           *v++ = *(vp + a->sliceheight * i);
1043d4002b98SHong Zhang           goto finished;
1044d4002b98SHong Zhang         }
1045d4002b98SHong Zhang       }
1046d4002b98SHong Zhang       *v++ = 0.0;
1047d4002b98SHong Zhang     finished:;
1048d4002b98SHong Zhang     }
1049d4002b98SHong Zhang   }
10503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1051d4002b98SHong Zhang }
1052d4002b98SHong Zhang 
1053d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL_ASCII(Mat A, PetscViewer viewer)
1054d71ae5a4SJacob Faibussowitsch {
1055d4002b98SHong Zhang   Mat_SeqSELL      *a = (Mat_SeqSELL *)A->data;
1056d4002b98SHong Zhang   PetscInt          i, j, m = A->rmap->n, shift;
1057d4002b98SHong Zhang   const char       *name;
1058d4002b98SHong Zhang   PetscViewerFormat format;
1059d4002b98SHong Zhang 
1060d4002b98SHong Zhang   PetscFunctionBegin;
10619566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
1062d4002b98SHong Zhang   if (format == PETSC_VIEWER_ASCII_MATLAB) {
1063d4002b98SHong Zhang     PetscInt nofinalvalue = 0;
1064d4002b98SHong Zhang     /*
1065d4002b98SHong Zhang     if (m && ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap->n-1))) {
1066d4002b98SHong Zhang       nofinalvalue = 1;
1067d4002b98SHong Zhang     }
1068d4002b98SHong Zhang     */
10699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
10709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Size = %" PetscInt_FMT " %" PetscInt_FMT " \n", m, A->cmap->n));
10719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Nonzeros = %" PetscInt_FMT " \n", a->nz));
1072d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
10739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",4);\n", a->nz + nofinalvalue));
1074d4002b98SHong Zhang #else
10759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",3);\n", a->nz + nofinalvalue));
1076d4002b98SHong Zhang #endif
10779566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = [\n"));
1078d4002b98SHong Zhang 
1079d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1080*07e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1081d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1082d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1083*07e43b41SHong Zhang         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e %18.16e\n", i + 1, a->colidx[shift + a->sliceheight * j] + 1, (double)PetscRealPart(a->val[shift + a->sliceheight * j]), (double)PetscImaginaryPart(a->val[shift + a->sliceheight * j])));
1084d4002b98SHong Zhang #else
1085*07e43b41SHong Zhang         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n", i + 1, a->colidx[shift + a->sliceheight * j] + 1, (double)a->val[shift + a->sliceheight * j]));
1086d4002b98SHong Zhang #endif
1087d4002b98SHong Zhang       }
1088d4002b98SHong Zhang     }
1089d4002b98SHong Zhang     /*
1090d4002b98SHong Zhang     if (nofinalvalue) {
1091d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
10929566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e %18.16e\n",m,A->cmap->n,0.,0.));
1093d4002b98SHong Zhang #else
10949566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n",m,A->cmap->n,0.0));
1095d4002b98SHong Zhang #endif
1096d4002b98SHong Zhang     }
1097d4002b98SHong Zhang     */
10989566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)A, &name));
10999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "];\n %s = spconvert(zzz);\n", name));
11009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1101d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
11023ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1103d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
11049566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1105d4002b98SHong Zhang     for (i = 0; i < m; i++) {
11069566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
1107*07e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1108d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1109d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1110*07e43b41SHong Zhang         if (PetscImaginaryPart(a->val[shift + a->sliceheight * j]) > 0.0 && PetscRealPart(a->val[shift + a->sliceheight * j]) != 0.0) {
1111*07e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[shift + a->sliceheight * j], (double)PetscRealPart(a->val[shift + a->sliceheight * j]), (double)PetscImaginaryPart(a->val[shift + a->sliceheight * j])));
1112*07e43b41SHong Zhang         } else if (PetscImaginaryPart(a->val[shift + a->sliceheight * j]) < 0.0 && PetscRealPart(a->val[shift + a->sliceheight * j]) != 0.0) {
1113*07e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[shift + a->sliceheight * j], (double)PetscRealPart(a->val[shift + a->sliceheight * j]), (double)-PetscImaginaryPart(a->val[shift + a->sliceheight * j])));
1114*07e43b41SHong Zhang         } else if (PetscRealPart(a->val[shift + a->sliceheight * j]) != 0.0) {
1115*07e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + a->sliceheight * j], (double)PetscRealPart(a->val[shift + a->sliceheight * j])));
1116d4002b98SHong Zhang         }
1117d4002b98SHong Zhang #else
1118*07e43b41SHong Zhang         if (a->val[shift + a->sliceheight * j] != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + a->sliceheight * j], (double)a->val[shift + a->sliceheight * j]));
1119d4002b98SHong Zhang #endif
1120d4002b98SHong Zhang       }
11219566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1122d4002b98SHong Zhang     }
11239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1124d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_DENSE) {
1125d4002b98SHong Zhang     PetscInt    cnt = 0, jcnt;
1126d4002b98SHong Zhang     PetscScalar value;
1127d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1128d4002b98SHong Zhang     PetscBool realonly = PETSC_TRUE;
1129d4002b98SHong Zhang     for (i = 0; i < a->sliidx[a->totalslices]; i++) {
1130d4002b98SHong Zhang       if (PetscImaginaryPart(a->val[i]) != 0.0) {
1131d4002b98SHong Zhang         realonly = PETSC_FALSE;
1132d4002b98SHong Zhang         break;
1133d4002b98SHong Zhang       }
1134d4002b98SHong Zhang     }
1135d4002b98SHong Zhang #endif
1136d4002b98SHong Zhang 
11379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1138d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1139d4002b98SHong Zhang       jcnt  = 0;
1140*07e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1141d4002b98SHong Zhang       for (j = 0; j < A->cmap->n; j++) {
1142*07e43b41SHong Zhang         if (jcnt < a->rlen[i] && j == a->colidx[shift + a->sliceheight * j]) {
1143d4002b98SHong Zhang           value = a->val[cnt++];
1144d4002b98SHong Zhang           jcnt++;
1145d4002b98SHong Zhang         } else {
1146d4002b98SHong Zhang           value = 0.0;
1147d4002b98SHong Zhang         }
1148d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1149d4002b98SHong Zhang         if (realonly) {
11509566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)PetscRealPart(value)));
1151d4002b98SHong Zhang         } else {
11529566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e+%7.5e i ", (double)PetscRealPart(value), (double)PetscImaginaryPart(value)));
1153d4002b98SHong Zhang         }
1154d4002b98SHong Zhang #else
11559566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)value));
1156d4002b98SHong Zhang #endif
1157d4002b98SHong Zhang       }
11589566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1159d4002b98SHong Zhang     }
11609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1161d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_MATRIXMARKET) {
1162d4002b98SHong Zhang     PetscInt fshift = 1;
11639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1164d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
11659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate complex general\n"));
1166d4002b98SHong Zhang #else
11679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate real general\n"));
1168d4002b98SHong Zhang #endif
11699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", m, A->cmap->n, a->nz));
1170d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1171*07e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1172d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1173d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1174*07e43b41SHong Zhang         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %g %g\n", i + fshift, a->colidx[shift + a->sliceheight * j] + fshift, (double)PetscRealPart(a->val[shift + a->sliceheight * j]), (double)PetscImaginaryPart(a->val[shift + a->sliceheight * j])));
1175d4002b98SHong Zhang #else
1176*07e43b41SHong Zhang         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %g\n", i + fshift, a->colidx[shift + a->sliceheight * j] + fshift, (double)a->val[shift + a->sliceheight * j]));
1177d4002b98SHong Zhang #endif
1178d4002b98SHong Zhang       }
1179d4002b98SHong Zhang     }
11809566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
118168aafef3SStefano Zampini   } else if (format == PETSC_VIEWER_NATIVE) {
118268aafef3SStefano Zampini     for (i = 0; i < a->totalslices; i++) { /* loop over slices */
118368aafef3SStefano Zampini       PetscInt row;
11849566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "slice %" PetscInt_FMT ": %" PetscInt_FMT " %" PetscInt_FMT "\n", i, a->sliidx[i], a->sliidx[i + 1]));
1185*07e43b41SHong Zhang       for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = (row + 1) % a->sliceheight) {
118668aafef3SStefano Zampini #if defined(PETSC_USE_COMPLEX)
118768aafef3SStefano Zampini         if (PetscImaginaryPart(a->val[j]) > 0.0) {
1188*07e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g + %g i\n", a->sliceheight * i + row, a->colidx[j], (double)PetscRealPart(a->val[j]), (double)PetscImaginaryPart(a->val[j])));
118968aafef3SStefano Zampini         } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
1190*07e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g - %g i\n", a->sliceheight * i + row, a->colidx[j], (double)PetscRealPart(a->val[j]), -(double)PetscImaginaryPart(a->val[j])));
119168aafef3SStefano Zampini         } else {
1192*07e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g\n", a->sliceheight * i + row, a->colidx[j], (double)PetscRealPart(a->val[j])));
119368aafef3SStefano Zampini         }
119468aafef3SStefano Zampini #else
1195*07e43b41SHong Zhang         PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g\n", a->sliceheight * i + row, a->colidx[j], (double)a->val[j]));
119668aafef3SStefano Zampini #endif
119768aafef3SStefano Zampini       }
119868aafef3SStefano Zampini     }
1199d4002b98SHong Zhang   } else {
12009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1201d4002b98SHong Zhang     if (A->factortype) {
1202d4002b98SHong Zhang       for (i = 0; i < m; i++) {
1203*07e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
12049566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
1205d4002b98SHong Zhang         /* L part */
1206*07e43b41SHong Zhang         for (j = shift; j < a->diag[i]; j += a->sliceheight) {
1207d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1208*07e43b41SHong Zhang           if (PetscImaginaryPart(a->val[shift + a->sliceheight * j]) > 0.0) {
12099566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)PetscImaginaryPart(a->val[j])));
1210*07e43b41SHong Zhang           } else if (PetscImaginaryPart(a->val[shift + a->sliceheight * j]) < 0.0) {
12119566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)(-PetscImaginaryPart(a->val[j]))));
1212d4002b98SHong Zhang           } else {
12139566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)PetscRealPart(a->val[j])));
1214d4002b98SHong Zhang           }
1215d4002b98SHong Zhang #else
12169566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)a->val[j]));
1217d4002b98SHong Zhang #endif
1218d4002b98SHong Zhang         }
1219d4002b98SHong Zhang         /* diagonal */
1220d4002b98SHong Zhang         j = a->diag[i];
1221d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1222d4002b98SHong Zhang         if (PetscImaginaryPart(a->val[j]) > 0.0) {
12239566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[j], (double)PetscRealPart(1.0 / a->val[j]), (double)PetscImaginaryPart(1.0 / a->val[j])));
1224d4002b98SHong Zhang         } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
12259566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[j], (double)PetscRealPart(1.0 / a->val[j]), (double)(-PetscImaginaryPart(1.0 / a->val[j]))));
1226d4002b98SHong Zhang         } else {
12279566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)PetscRealPart(1.0 / a->val[j])));
1228d4002b98SHong Zhang         }
1229d4002b98SHong Zhang #else
12309566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)(1.0 / a->val[j])));
1231d4002b98SHong Zhang #endif
1232d4002b98SHong Zhang 
1233d4002b98SHong Zhang         /* U part */
1234*07e43b41SHong Zhang         for (j = a->diag[i] + 1; j < shift + a->sliceheight * a->rlen[i]; j += a->sliceheight) {
1235d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1236d4002b98SHong Zhang           if (PetscImaginaryPart(a->val[j]) > 0.0) {
12379566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)PetscImaginaryPart(a->val[j])));
1238d4002b98SHong Zhang           } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
12399566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)(-PetscImaginaryPart(a->val[j]))));
1240d4002b98SHong Zhang           } else {
12419566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)PetscRealPart(a->val[j])));
1242d4002b98SHong Zhang           }
1243d4002b98SHong Zhang #else
12449566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)a->val[j]));
1245d4002b98SHong Zhang #endif
1246d4002b98SHong Zhang         }
12479566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1248d4002b98SHong Zhang       }
1249d4002b98SHong Zhang     } else {
1250d4002b98SHong Zhang       for (i = 0; i < m; i++) {
1251*07e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
12529566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
1253d4002b98SHong Zhang         for (j = 0; j < a->rlen[i]; j++) {
1254d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1255d4002b98SHong Zhang           if (PetscImaginaryPart(a->val[j]) > 0.0) {
1256*07e43b41SHong Zhang             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[shift + a->sliceheight * j], (double)PetscRealPart(a->val[shift + a->sliceheight * j]), (double)PetscImaginaryPart(a->val[shift + a->sliceheight * j])));
1257d4002b98SHong Zhang           } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
1258*07e43b41SHong Zhang             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[shift + a->sliceheight * j], (double)PetscRealPart(a->val[shift + a->sliceheight * j]), (double)-PetscImaginaryPart(a->val[shift + a->sliceheight * j])));
1259d4002b98SHong Zhang           } else {
1260*07e43b41SHong Zhang             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + a->sliceheight * j], (double)PetscRealPart(a->val[shift + a->sliceheight * j])));
1261d4002b98SHong Zhang           }
1262d4002b98SHong Zhang #else
1263*07e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + a->sliceheight * j], (double)a->val[shift + a->sliceheight * j]));
1264d4002b98SHong Zhang #endif
1265d4002b98SHong Zhang         }
12669566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1267d4002b98SHong Zhang       }
1268d4002b98SHong Zhang     }
12699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1270d4002b98SHong Zhang   }
12719566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
12723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1273d4002b98SHong Zhang }
1274d4002b98SHong Zhang 
1275d4002b98SHong Zhang #include <petscdraw.h>
1276d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL_Draw_Zoom(PetscDraw draw, void *Aa)
1277d71ae5a4SJacob Faibussowitsch {
1278d4002b98SHong Zhang   Mat               A = (Mat)Aa;
1279d4002b98SHong Zhang   Mat_SeqSELL      *a = (Mat_SeqSELL *)A->data;
1280d4002b98SHong Zhang   PetscInt          i, j, m = A->rmap->n, shift;
1281d4002b98SHong Zhang   int               color;
1282d4002b98SHong Zhang   PetscReal         xl, yl, xr, yr, x_l, x_r, y_l, y_r;
1283d4002b98SHong Zhang   PetscViewer       viewer;
1284d4002b98SHong Zhang   PetscViewerFormat format;
1285d4002b98SHong Zhang 
1286d4002b98SHong Zhang   PetscFunctionBegin;
12879566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "Zoomviewer", (PetscObject *)&viewer));
12889566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
12899566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
1290d4002b98SHong Zhang 
1291d4002b98SHong Zhang   /* loop over matrix elements drawing boxes */
1292d4002b98SHong Zhang 
1293d4002b98SHong Zhang   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1294d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
1295d4002b98SHong Zhang     /* Blue for negative, Cyan for zero and  Red for positive */
1296d4002b98SHong Zhang     color = PETSC_DRAW_BLUE;
1297d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1298*07e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
12999371c9d4SSatish Balay       y_l   = m - i - 1.0;
13009371c9d4SSatish Balay       y_r   = y_l + 1.0;
1301d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1302*07e43b41SHong Zhang         x_l = a->colidx[shift + a->sliceheight * j];
13039371c9d4SSatish Balay         x_r = x_l + 1.0;
1304*07e43b41SHong Zhang         if (PetscRealPart(a->val[shift + a->sliceheight * j]) >= 0.) continue;
13059566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1306d4002b98SHong Zhang       }
1307d4002b98SHong Zhang     }
1308d4002b98SHong Zhang     color = PETSC_DRAW_CYAN;
1309d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1310*07e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
13119371c9d4SSatish Balay       y_l   = m - i - 1.0;
13129371c9d4SSatish Balay       y_r   = y_l + 1.0;
1313d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1314*07e43b41SHong Zhang         x_l = a->colidx[shift + a->sliceheight * j];
13159371c9d4SSatish Balay         x_r = x_l + 1.0;
1316*07e43b41SHong Zhang         if (a->val[shift + a->sliceheight * j] != 0.) continue;
13179566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1318d4002b98SHong Zhang       }
1319d4002b98SHong Zhang     }
1320d4002b98SHong Zhang     color = PETSC_DRAW_RED;
1321d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1322*07e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
13239371c9d4SSatish Balay       y_l   = m - i - 1.0;
13249371c9d4SSatish Balay       y_r   = y_l + 1.0;
1325d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1326*07e43b41SHong Zhang         x_l = a->colidx[shift + a->sliceheight * j];
13279371c9d4SSatish Balay         x_r = x_l + 1.0;
1328*07e43b41SHong Zhang         if (PetscRealPart(a->val[shift + a->sliceheight * j]) <= 0.) continue;
13299566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1330d4002b98SHong Zhang       }
1331d4002b98SHong Zhang     }
1332d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1333d4002b98SHong Zhang   } else {
1334d4002b98SHong Zhang     /* use contour shading to indicate magnitude of values */
1335d4002b98SHong Zhang     /* first determine max of all nonzero values */
1336d4002b98SHong Zhang     PetscReal minv = 0.0, maxv = 0.0;
1337d4002b98SHong Zhang     PetscInt  count = 0;
1338d4002b98SHong Zhang     PetscDraw popup;
1339d4002b98SHong Zhang     for (i = 0; i < a->sliidx[a->totalslices]; i++) {
1340d4002b98SHong Zhang       if (PetscAbsScalar(a->val[i]) > maxv) maxv = PetscAbsScalar(a->val[i]);
1341d4002b98SHong Zhang     }
1342d4002b98SHong Zhang     if (minv >= maxv) maxv = minv + PETSC_SMALL;
13439566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetPopup(draw, &popup));
13449566063dSJacob Faibussowitsch     PetscCall(PetscDrawScalePopup(popup, minv, maxv));
1345d4002b98SHong Zhang 
1346d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
1347d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1348*07e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1349d4002b98SHong Zhang       y_l   = m - i - 1.0;
1350d4002b98SHong Zhang       y_r   = y_l + 1.0;
1351d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1352*07e43b41SHong Zhang         x_l   = a->colidx[shift + a->sliceheight * j];
1353d4002b98SHong Zhang         x_r   = x_l + 1.0;
1354d4002b98SHong Zhang         color = PetscDrawRealToColor(PetscAbsScalar(a->val[count]), minv, maxv);
13559566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1356d4002b98SHong Zhang         count++;
1357d4002b98SHong Zhang       }
1358d4002b98SHong Zhang     }
1359d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1360d4002b98SHong Zhang   }
13613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1362d4002b98SHong Zhang }
1363d4002b98SHong Zhang 
1364d4002b98SHong Zhang #include <petscdraw.h>
1365d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL_Draw(Mat A, PetscViewer viewer)
1366d71ae5a4SJacob Faibussowitsch {
1367d4002b98SHong Zhang   PetscDraw draw;
1368d4002b98SHong Zhang   PetscReal xr, yr, xl, yl, h, w;
1369d4002b98SHong Zhang   PetscBool isnull;
1370d4002b98SHong Zhang 
1371d4002b98SHong Zhang   PetscFunctionBegin;
13729566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
13739566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
13743ba16761SJacob Faibussowitsch   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
1375d4002b98SHong Zhang 
13769371c9d4SSatish Balay   xr = A->cmap->n;
13779371c9d4SSatish Balay   yr = A->rmap->n;
13789371c9d4SSatish Balay   h  = yr / 10.0;
13799371c9d4SSatish Balay   w  = xr / 10.0;
13809371c9d4SSatish Balay   xr += w;
13819371c9d4SSatish Balay   yr += h;
13829371c9d4SSatish Balay   xl = -w;
13839371c9d4SSatish Balay   yl = -h;
13849566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
13859566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", (PetscObject)viewer));
13869566063dSJacob Faibussowitsch   PetscCall(PetscDrawZoom(draw, MatView_SeqSELL_Draw_Zoom, A));
13879566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", NULL));
13889566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
13893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1390d4002b98SHong Zhang }
1391d4002b98SHong Zhang 
1392d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL(Mat A, PetscViewer viewer)
1393d71ae5a4SJacob Faibussowitsch {
1394d4002b98SHong Zhang   PetscBool iascii, isbinary, isdraw;
1395d4002b98SHong Zhang 
1396d4002b98SHong Zhang   PetscFunctionBegin;
13979566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
13989566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
13999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1400d4002b98SHong Zhang   if (iascii) {
14019566063dSJacob Faibussowitsch     PetscCall(MatView_SeqSELL_ASCII(A, viewer));
1402d4002b98SHong Zhang   } else if (isbinary) {
14039566063dSJacob Faibussowitsch     /* PetscCall(MatView_SeqSELL_Binary(A,viewer)); */
14041baa6e33SBarry Smith   } else if (isdraw) PetscCall(MatView_SeqSELL_Draw(A, viewer));
14053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1406d4002b98SHong Zhang }
1407d4002b98SHong Zhang 
1408d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_SeqSELL(Mat A, MatAssemblyType mode)
1409d71ae5a4SJacob Faibussowitsch {
1410d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1411d4002b98SHong Zhang   PetscInt     i, shift, row_in_slice, row, nrow, *cp, lastcol, j, k;
1412d4002b98SHong Zhang   MatScalar   *vp;
1413d4002b98SHong Zhang 
1414d4002b98SHong Zhang   PetscFunctionBegin;
14153ba16761SJacob Faibussowitsch   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
1416d4002b98SHong Zhang   /* To do: compress out the unused elements */
14179566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_SeqSELL(A));
14189566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Matrix size: %" PetscInt_FMT " X %" PetscInt_FMT "; storage space: %" PetscInt_FMT " allocated %" PetscInt_FMT " used (%" PetscInt_FMT " nonzeros+%" PetscInt_FMT " paddedzeros)\n", A->rmap->n, A->cmap->n, a->maxallocmat, a->sliidx[a->totalslices], a->nz, a->sliidx[a->totalslices] - a->nz));
14199566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs));
14209566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", a->rlenmax));
14212d1451d4SHong Zhang   a->nonzerorowcnt = 0;
1422d4002b98SHong Zhang   /* Set unused slots for column indices to last valid column index. Set unused slots for values to zero. This allows for a use of unmasked intrinsics -> higher performance */
1423d4002b98SHong Zhang   for (i = 0; i < a->totalslices; ++i) {
1424d4002b98SHong Zhang     shift = a->sliidx[i];                                                   /* starting index of the slice */
1425d4002b98SHong Zhang     cp    = a->colidx + shift;                                              /* pointer to the column indices of the slice */
1426d4002b98SHong Zhang     vp    = a->val + shift;                                                 /* pointer to the nonzero values of the slice */
1427*07e43b41SHong Zhang     for (row_in_slice = 0; row_in_slice < a->sliceheight; ++row_in_slice) { /* loop over rows in the slice */
1428*07e43b41SHong Zhang       row  = a->sliceheight * i + row_in_slice;
1429d4002b98SHong Zhang       nrow = a->rlen[row]; /* number of nonzeros in row */
1430d4002b98SHong Zhang       /*
1431d4002b98SHong Zhang         Search for the nearest nonzero. Normally setting the index to zero may cause extra communication.
1432d4002b98SHong Zhang         But if the entire slice are empty, it is fine to use 0 since the index will not be loaded.
1433d4002b98SHong Zhang       */
1434d4002b98SHong Zhang       lastcol = 0;
1435d4002b98SHong Zhang       if (nrow > 0) { /* nonempty row */
14362d1451d4SHong Zhang         a->nonzerorowcnt++;
1437*07e43b41SHong Zhang         lastcol = cp[a->sliceheight * (nrow - 1) + row_in_slice]; /* use the index from the last nonzero at current row */
1438aaa8cc7dSPierre Jolivet       } else if (!row_in_slice) {                                 /* first row of the correct slice is empty */
1439*07e43b41SHong Zhang         for (j = 1; j < a->sliceheight; j++) {
1440*07e43b41SHong Zhang           if (a->rlen[a->sliceheight * i + j]) {
1441d4002b98SHong Zhang             lastcol = cp[j];
1442d4002b98SHong Zhang             break;
1443d4002b98SHong Zhang           }
1444d4002b98SHong Zhang         }
1445d4002b98SHong Zhang       } else {
1446d4002b98SHong Zhang         if (a->sliidx[i + 1] != shift) lastcol = cp[row_in_slice - 1]; /* use the index from the previous row */
1447d4002b98SHong Zhang       }
1448d4002b98SHong Zhang 
1449*07e43b41SHong Zhang       for (k = nrow; k < (a->sliidx[i + 1] - shift) / a->sliceheight; ++k) {
1450*07e43b41SHong Zhang         cp[a->sliceheight * k + row_in_slice] = lastcol;
1451*07e43b41SHong Zhang         vp[a->sliceheight * k + row_in_slice] = (MatScalar)0;
1452d4002b98SHong Zhang       }
1453d4002b98SHong Zhang     }
1454d4002b98SHong Zhang   }
1455d4002b98SHong Zhang 
1456d4002b98SHong Zhang   A->info.mallocs += a->reallocs;
1457d4002b98SHong Zhang   a->reallocs = 0;
1458d4002b98SHong Zhang 
14599566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
14603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1461d4002b98SHong Zhang }
1462d4002b98SHong Zhang 
1463d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetInfo_SeqSELL(Mat A, MatInfoType flag, MatInfo *info)
1464d71ae5a4SJacob Faibussowitsch {
1465d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1466d4002b98SHong Zhang 
1467d4002b98SHong Zhang   PetscFunctionBegin;
1468d4002b98SHong Zhang   info->block_size   = 1.0;
14693966268fSBarry Smith   info->nz_allocated = a->maxallocmat;
14703966268fSBarry Smith   info->nz_used      = a->sliidx[a->totalslices]; /* include padding zeros */
14713966268fSBarry Smith   info->nz_unneeded  = (a->maxallocmat - a->sliidx[a->totalslices]);
14723966268fSBarry Smith   info->assemblies   = A->num_ass;
14733966268fSBarry Smith   info->mallocs      = A->info.mallocs;
14744dfa11a4SJacob Faibussowitsch   info->memory       = 0; /* REVIEW ME */
1475d4002b98SHong Zhang   if (A->factortype) {
1476d4002b98SHong Zhang     info->fill_ratio_given  = A->info.fill_ratio_given;
1477d4002b98SHong Zhang     info->fill_ratio_needed = A->info.fill_ratio_needed;
1478d4002b98SHong Zhang     info->factor_mallocs    = A->info.factor_mallocs;
1479d4002b98SHong Zhang   } else {
1480d4002b98SHong Zhang     info->fill_ratio_given  = 0;
1481d4002b98SHong Zhang     info->fill_ratio_needed = 0;
1482d4002b98SHong Zhang     info->factor_mallocs    = 0;
1483d4002b98SHong Zhang   }
14843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1485d4002b98SHong Zhang }
1486d4002b98SHong Zhang 
1487d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetValues_SeqSELL(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
1488d71ae5a4SJacob Faibussowitsch {
1489d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1490d4002b98SHong Zhang   PetscInt     shift, i, k, l, low, high, t, ii, row, col, nrow;
1491d4002b98SHong Zhang   PetscInt    *cp, nonew = a->nonew, lastcol = -1;
1492d4002b98SHong Zhang   MatScalar   *vp, value;
14932d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
14942d1451d4SHong Zhang   PetscBool inserted = PETSC_FALSE;
14952d1451d4SHong Zhang #endif
1496d4002b98SHong Zhang 
1497d4002b98SHong Zhang   PetscFunctionBegin;
1498d4002b98SHong Zhang   for (k = 0; k < m; k++) { /* loop over added rows */
1499d4002b98SHong Zhang     row = im[k];
1500d4002b98SHong Zhang     if (row < 0) continue;
15016bdcaf15SBarry Smith     PetscCheck(row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, A->rmap->n - 1);
1502*07e43b41SHong Zhang     shift = a->sliidx[row / a->sliceheight] + row % a->sliceheight; /* starting index of the row */
1503d4002b98SHong Zhang     cp    = a->colidx + shift;                                      /* pointer to the row */
1504d4002b98SHong Zhang     vp    = a->val + shift;                                         /* pointer to the row */
1505d4002b98SHong Zhang     nrow  = a->rlen[row];
1506d4002b98SHong Zhang     low   = 0;
1507d4002b98SHong Zhang     high  = nrow;
1508d4002b98SHong Zhang 
1509d4002b98SHong Zhang     for (l = 0; l < n; l++) { /* loop over added columns */
1510d4002b98SHong Zhang       col = in[l];
1511d4002b98SHong Zhang       if (col < 0) continue;
15126bdcaf15SBarry Smith       PetscCheck(col < A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Col too large: row %" PetscInt_FMT " max %" PetscInt_FMT, col, A->cmap->n - 1);
1513d4002b98SHong Zhang       if (a->roworiented) {
1514d4002b98SHong Zhang         value = v[l + k * n];
1515d4002b98SHong Zhang       } else {
1516d4002b98SHong Zhang         value = v[k + l * m];
1517d4002b98SHong Zhang       }
1518d4002b98SHong Zhang       if ((value == 0.0 && a->ignorezeroentries) && (is == ADD_VALUES)) continue;
1519d4002b98SHong Zhang 
1520ed73aabaSBarry Smith       /* search in this row for the specified column, i indicates the column to be set */
1521d4002b98SHong Zhang       if (col <= lastcol) low = 0;
1522d4002b98SHong Zhang       else high = nrow;
1523d4002b98SHong Zhang       lastcol = col;
1524d4002b98SHong Zhang       while (high - low > 5) {
1525d4002b98SHong Zhang         t = (low + high) / 2;
1526*07e43b41SHong Zhang         if (*(cp + a->sliceheight * t) > col) high = t;
1527d4002b98SHong Zhang         else low = t;
1528d4002b98SHong Zhang       }
1529d4002b98SHong Zhang       for (i = low; i < high; i++) {
1530*07e43b41SHong Zhang         if (*(cp + a->sliceheight * i) > col) break;
1531*07e43b41SHong Zhang         if (*(cp + a->sliceheight * i) == col) {
1532*07e43b41SHong Zhang           if (is == ADD_VALUES) *(vp + a->sliceheight * i) += value;
1533*07e43b41SHong Zhang           else *(vp + a->sliceheight * i) = value;
15342d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
15352d1451d4SHong Zhang           inserted = PETSC_TRUE;
15362d1451d4SHong Zhang #endif
1537d4002b98SHong Zhang           low = i + 1;
1538d4002b98SHong Zhang           goto noinsert;
1539d4002b98SHong Zhang         }
1540d4002b98SHong Zhang       }
1541d4002b98SHong Zhang       if (value == 0.0 && a->ignorezeroentries) goto noinsert;
1542d4002b98SHong Zhang       if (nonew == 1) goto noinsert;
154308401ef6SPierre Jolivet       PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
1544d4002b98SHong Zhang       /* If the current row length exceeds the slice width (e.g. nrow==slice_width), allocate a new space, otherwise do nothing */
1545*07e43b41SHong Zhang       MatSeqXSELLReallocateSELL(A, A->rmap->n, 1, nrow, a->sliidx, a->sliceheight, row / a->sliceheight, row, col, a->colidx, a->val, cp, vp, nonew, MatScalar);
1546d4002b98SHong Zhang       /* add the new nonzero to the high position, shift the remaining elements in current row to the right by one slot */
1547d4002b98SHong Zhang       for (ii = nrow - 1; ii >= i; ii--) {
1548*07e43b41SHong Zhang         *(cp + a->sliceheight * (ii + 1)) = *(cp + a->sliceheight * ii);
1549*07e43b41SHong Zhang         *(vp + a->sliceheight * (ii + 1)) = *(vp + a->sliceheight * ii);
1550d4002b98SHong Zhang       }
1551d4002b98SHong Zhang       a->rlen[row]++;
1552*07e43b41SHong Zhang       *(cp + a->sliceheight * i) = col;
1553*07e43b41SHong Zhang       *(vp + a->sliceheight * i) = value;
1554d4002b98SHong Zhang       a->nz++;
1555d4002b98SHong Zhang       A->nonzerostate++;
15562d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
15572d1451d4SHong Zhang       inserted = PETSC_TRUE;
15582d1451d4SHong Zhang #endif
15599371c9d4SSatish Balay       low = i + 1;
15609371c9d4SSatish Balay       high++;
15619371c9d4SSatish Balay       nrow++;
1562d4002b98SHong Zhang     noinsert:;
1563d4002b98SHong Zhang     }
1564d4002b98SHong Zhang     a->rlen[row] = nrow;
1565d4002b98SHong Zhang   }
15662d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
15672d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED && inserted) A->offloadmask = PETSC_OFFLOAD_CPU;
15682d1451d4SHong Zhang #endif
15693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1570d4002b98SHong Zhang }
1571d4002b98SHong Zhang 
1572d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCopy_SeqSELL(Mat A, Mat B, MatStructure str)
1573d71ae5a4SJacob Faibussowitsch {
1574d4002b98SHong Zhang   PetscFunctionBegin;
1575d4002b98SHong Zhang   /* If the two matrices have the same copy implementation, use fast copy. */
1576d4002b98SHong Zhang   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1577d4002b98SHong Zhang     Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1578d4002b98SHong Zhang     Mat_SeqSELL *b = (Mat_SeqSELL *)B->data;
1579d4002b98SHong Zhang 
158008401ef6SPierre Jolivet     PetscCheck(a->sliidx[a->totalslices] == b->sliidx[b->totalslices], PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of nonzeros in two matrices are different");
15819566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(b->val, a->val, a->sliidx[a->totalslices]));
1582d4002b98SHong Zhang   } else {
15839566063dSJacob Faibussowitsch     PetscCall(MatCopy_Basic(A, B, str));
1584d4002b98SHong Zhang   }
15853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1586d4002b98SHong Zhang }
1587d4002b98SHong Zhang 
1588d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUp_SeqSELL(Mat A)
1589d71ae5a4SJacob Faibussowitsch {
1590d4002b98SHong Zhang   PetscFunctionBegin;
15919566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(A, PETSC_DEFAULT, NULL));
15923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1593d4002b98SHong Zhang }
1594d4002b98SHong Zhang 
1595d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLGetArray_SeqSELL(Mat A, PetscScalar *array[])
1596d71ae5a4SJacob Faibussowitsch {
1597d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1598d4002b98SHong Zhang 
1599d4002b98SHong Zhang   PetscFunctionBegin;
1600d4002b98SHong Zhang   *array = a->val;
16013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1602d4002b98SHong Zhang }
1603d4002b98SHong Zhang 
1604d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLRestoreArray_SeqSELL(Mat A, PetscScalar *array[])
1605d71ae5a4SJacob Faibussowitsch {
1606d4002b98SHong Zhang   PetscFunctionBegin;
16073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1608d4002b98SHong Zhang }
1609d4002b98SHong Zhang 
1610d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRealPart_SeqSELL(Mat A)
1611d71ae5a4SJacob Faibussowitsch {
1612d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1613d4002b98SHong Zhang   PetscInt     i;
1614d4002b98SHong Zhang   MatScalar   *aval = a->val;
1615d4002b98SHong Zhang 
1616d4002b98SHong Zhang   PetscFunctionBegin;
1617d4002b98SHong Zhang   for (i = 0; i < a->sliidx[a->totalslices]; i++) aval[i] = PetscRealPart(aval[i]);
16182d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
16192d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
16202d1451d4SHong Zhang #endif
16213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1622d4002b98SHong Zhang }
1623d4002b98SHong Zhang 
1624d71ae5a4SJacob Faibussowitsch PetscErrorCode MatImaginaryPart_SeqSELL(Mat A)
1625d71ae5a4SJacob Faibussowitsch {
1626d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1627d4002b98SHong Zhang   PetscInt     i;
1628d4002b98SHong Zhang   MatScalar   *aval = a->val;
1629d4002b98SHong Zhang 
1630d4002b98SHong Zhang   PetscFunctionBegin;
1631d4002b98SHong Zhang   for (i = 0; i < a->sliidx[a->totalslices]; i++) aval[i] = PetscImaginaryPart(aval[i]);
16329566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
16332d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
16342d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
16352d1451d4SHong Zhang #endif
16363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1637d4002b98SHong Zhang }
1638d4002b98SHong Zhang 
1639d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_SeqSELL(Mat inA, PetscScalar alpha)
1640d71ae5a4SJacob Faibussowitsch {
1641d4002b98SHong Zhang   Mat_SeqSELL *a      = (Mat_SeqSELL *)inA->data;
1642d4002b98SHong Zhang   MatScalar   *aval   = a->val;
1643d4002b98SHong Zhang   PetscScalar  oalpha = alpha;
1644d4002b98SHong Zhang   PetscBLASInt one    = 1, size;
1645d4002b98SHong Zhang 
1646d4002b98SHong Zhang   PetscFunctionBegin;
16479566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(a->sliidx[a->totalslices], &size));
1648792fecdfSBarry Smith   PetscCallBLAS("BLASscal", BLASscal_(&size, &oalpha, aval, &one));
16499566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(a->nz));
16509566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(inA));
16512d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
16522d1451d4SHong Zhang   if (inA->offloadmask != PETSC_OFFLOAD_UNALLOCATED) inA->offloadmask = PETSC_OFFLOAD_CPU;
16532d1451d4SHong Zhang #endif
16543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1655d4002b98SHong Zhang }
1656d4002b98SHong Zhang 
1657d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift_SeqSELL(Mat Y, PetscScalar a)
1658d71ae5a4SJacob Faibussowitsch {
1659d4002b98SHong Zhang   Mat_SeqSELL *y = (Mat_SeqSELL *)Y->data;
1660d4002b98SHong Zhang 
1661d4002b98SHong Zhang   PetscFunctionBegin;
166248a46eb9SPierre Jolivet   if (!Y->preallocated || !y->nz) PetscCall(MatSeqSELLSetPreallocation(Y, 1, NULL));
16639566063dSJacob Faibussowitsch   PetscCall(MatShift_Basic(Y, a));
16643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1665d4002b98SHong Zhang }
1666d4002b98SHong Zhang 
1667d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSOR_SeqSELL(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
1668d71ae5a4SJacob Faibussowitsch {
1669d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
1670d4002b98SHong Zhang   PetscScalar       *x, sum, *t;
1671f4259b30SLisandro Dalcin   const MatScalar   *idiag = NULL, *mdiag;
1672d4002b98SHong Zhang   const PetscScalar *b, *xb;
1673d4002b98SHong Zhang   PetscInt           n, m = A->rmap->n, i, j, shift;
1674d4002b98SHong Zhang   const PetscInt    *diag;
1675d4002b98SHong Zhang 
1676d4002b98SHong Zhang   PetscFunctionBegin;
1677d4002b98SHong Zhang   its = its * lits;
1678d4002b98SHong Zhang 
1679d4002b98SHong Zhang   if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
16809566063dSJacob Faibussowitsch   if (!a->idiagvalid) PetscCall(MatInvertDiagonal_SeqSELL(A, omega, fshift));
1681d4002b98SHong Zhang   a->fshift = fshift;
1682d4002b98SHong Zhang   a->omega  = omega;
1683d4002b98SHong Zhang 
1684d4002b98SHong Zhang   diag  = a->diag;
1685d4002b98SHong Zhang   t     = a->ssor_work;
1686d4002b98SHong Zhang   idiag = a->idiag;
1687d4002b98SHong Zhang   mdiag = a->mdiag;
1688d4002b98SHong Zhang 
16899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
16909566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1691d4002b98SHong Zhang   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
169208401ef6SPierre Jolivet   PetscCheck(flag != SOR_APPLY_UPPER, PETSC_COMM_SELF, PETSC_ERR_SUP, "SOR_APPLY_UPPER is not implemented");
169308401ef6SPierre Jolivet   PetscCheck(flag != SOR_APPLY_LOWER, PETSC_COMM_SELF, PETSC_ERR_SUP, "SOR_APPLY_LOWER is not implemented");
1694aed4548fSBarry Smith   PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat");
1695d4002b98SHong Zhang 
1696d4002b98SHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
1697d4002b98SHong Zhang     if ((flag & SOR_FORWARD_SWEEP) || (flag & SOR_LOCAL_FORWARD_SWEEP)) {
1698d4002b98SHong Zhang       for (i = 0; i < m; i++) {
1699*07e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
1700d4002b98SHong Zhang         sum   = b[i];
1701*07e43b41SHong Zhang         n     = (diag[i] - shift) / a->sliceheight;
1702*07e43b41SHong Zhang         for (j = 0; j < n; j++) sum -= a->val[shift + a->sliceheight * j] * x[a->colidx[shift + a->sliceheight * j]];
1703d4002b98SHong Zhang         t[i] = sum;
1704d4002b98SHong Zhang         x[i] = sum * idiag[i];
1705d4002b98SHong Zhang       }
1706d4002b98SHong Zhang       xb = t;
17079566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(a->nz));
1708d4002b98SHong Zhang     } else xb = b;
1709d4002b98SHong Zhang     if ((flag & SOR_BACKWARD_SWEEP) || (flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1710d4002b98SHong Zhang       for (i = m - 1; i >= 0; i--) {
1711*07e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
1712d4002b98SHong Zhang         sum   = xb[i];
1713*07e43b41SHong Zhang         n     = a->rlen[i] - (diag[i] - shift) / a->sliceheight - 1;
1714*07e43b41SHong Zhang         for (j = 1; j <= n; j++) sum -= a->val[diag[i] + a->sliceheight * j] * x[a->colidx[diag[i] + a->sliceheight * j]];
1715d4002b98SHong Zhang         if (xb == b) {
1716d4002b98SHong Zhang           x[i] = sum * idiag[i];
1717d4002b98SHong Zhang         } else {
1718d4002b98SHong Zhang           x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
1719d4002b98SHong Zhang         }
1720d4002b98SHong Zhang       }
17219566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
1722d4002b98SHong Zhang     }
1723d4002b98SHong Zhang     its--;
1724d4002b98SHong Zhang   }
1725d4002b98SHong Zhang   while (its--) {
1726d4002b98SHong Zhang     if ((flag & SOR_FORWARD_SWEEP) || (flag & SOR_LOCAL_FORWARD_SWEEP)) {
1727d4002b98SHong Zhang       for (i = 0; i < m; i++) {
1728d4002b98SHong Zhang         /* lower */
1729*07e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
1730d4002b98SHong Zhang         sum   = b[i];
1731*07e43b41SHong Zhang         n     = (diag[i] - shift) / a->sliceheight;
1732*07e43b41SHong Zhang         for (j = 0; j < n; j++) sum -= a->val[shift + a->sliceheight * j] * x[a->colidx[shift + a->sliceheight * j]];
1733d4002b98SHong Zhang         t[i] = sum; /* save application of the lower-triangular part */
1734d4002b98SHong Zhang         /* upper */
1735*07e43b41SHong Zhang         n = a->rlen[i] - (diag[i] - shift) / a->sliceheight - 1;
1736*07e43b41SHong Zhang         for (j = 1; j <= n; j++) sum -= a->val[diag[i] + a->sliceheight * j] * x[a->colidx[diag[i] + a->sliceheight * j]];
1737d4002b98SHong Zhang         x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
1738d4002b98SHong Zhang       }
1739d4002b98SHong Zhang       xb = t;
17409566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(2.0 * a->nz));
1741d4002b98SHong Zhang     } else xb = b;
1742d4002b98SHong Zhang     if ((flag & SOR_BACKWARD_SWEEP) || (flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1743d4002b98SHong Zhang       for (i = m - 1; i >= 0; i--) {
1744*07e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
1745d4002b98SHong Zhang         sum   = xb[i];
1746d4002b98SHong Zhang         if (xb == b) {
1747d4002b98SHong Zhang           /* whole matrix (no checkpointing available) */
1748d4002b98SHong Zhang           n = a->rlen[i];
1749*07e43b41SHong Zhang           for (j = 0; j < n; j++) sum -= a->val[shift + a->sliceheight * j] * x[a->colidx[shift + a->sliceheight * j]];
1750d4002b98SHong Zhang           x[i] = (1. - omega) * x[i] + (sum + mdiag[i] * x[i]) * idiag[i];
1751d4002b98SHong Zhang         } else { /* lower-triangular part has been saved, so only apply upper-triangular */
1752*07e43b41SHong Zhang           n = a->rlen[i] - (diag[i] - shift) / a->sliceheight - 1;
1753*07e43b41SHong Zhang           for (j = 1; j <= n; j++) sum -= a->val[diag[i] + a->sliceheight * j] * x[a->colidx[diag[i] + a->sliceheight * j]];
1754d4002b98SHong Zhang           x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
1755d4002b98SHong Zhang         }
1756d4002b98SHong Zhang       }
1757d4002b98SHong Zhang       if (xb == b) {
17589566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(2.0 * a->nz));
1759d4002b98SHong Zhang       } else {
17609566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
1761d4002b98SHong Zhang       }
1762d4002b98SHong Zhang     }
1763d4002b98SHong Zhang   }
17649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
17659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
17663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1767d4002b98SHong Zhang }
1768d4002b98SHong Zhang 
1769d4002b98SHong Zhang static struct _MatOps MatOps_Values = {MatSetValues_SeqSELL,
17706108893eSStefano Zampini                                        MatGetRow_SeqSELL,
17716108893eSStefano Zampini                                        MatRestoreRow_SeqSELL,
1772d4002b98SHong Zhang                                        MatMult_SeqSELL,
1773d4002b98SHong Zhang                                        /* 4*/ MatMultAdd_SeqSELL,
1774d4002b98SHong Zhang                                        MatMultTranspose_SeqSELL,
1775d4002b98SHong Zhang                                        MatMultTransposeAdd_SeqSELL,
1776f4259b30SLisandro Dalcin                                        NULL,
1777f4259b30SLisandro Dalcin                                        NULL,
1778f4259b30SLisandro Dalcin                                        NULL,
1779f4259b30SLisandro Dalcin                                        /* 10*/ NULL,
1780f4259b30SLisandro Dalcin                                        NULL,
1781f4259b30SLisandro Dalcin                                        NULL,
1782d4002b98SHong Zhang                                        MatSOR_SeqSELL,
1783f4259b30SLisandro Dalcin                                        NULL,
1784d4002b98SHong Zhang                                        /* 15*/ MatGetInfo_SeqSELL,
1785d4002b98SHong Zhang                                        MatEqual_SeqSELL,
1786d4002b98SHong Zhang                                        MatGetDiagonal_SeqSELL,
1787d4002b98SHong Zhang                                        MatDiagonalScale_SeqSELL,
1788f4259b30SLisandro Dalcin                                        NULL,
1789f4259b30SLisandro Dalcin                                        /* 20*/ NULL,
1790d4002b98SHong Zhang                                        MatAssemblyEnd_SeqSELL,
1791d4002b98SHong Zhang                                        MatSetOption_SeqSELL,
1792d4002b98SHong Zhang                                        MatZeroEntries_SeqSELL,
1793f4259b30SLisandro Dalcin                                        /* 24*/ NULL,
1794f4259b30SLisandro Dalcin                                        NULL,
1795f4259b30SLisandro Dalcin                                        NULL,
1796f4259b30SLisandro Dalcin                                        NULL,
1797f4259b30SLisandro Dalcin                                        NULL,
1798d4002b98SHong Zhang                                        /* 29*/ MatSetUp_SeqSELL,
1799f4259b30SLisandro Dalcin                                        NULL,
1800f4259b30SLisandro Dalcin                                        NULL,
1801f4259b30SLisandro Dalcin                                        NULL,
1802f4259b30SLisandro Dalcin                                        NULL,
1803d4002b98SHong Zhang                                        /* 34*/ MatDuplicate_SeqSELL,
1804f4259b30SLisandro Dalcin                                        NULL,
1805f4259b30SLisandro Dalcin                                        NULL,
1806f4259b30SLisandro Dalcin                                        NULL,
1807f4259b30SLisandro Dalcin                                        NULL,
1808f4259b30SLisandro Dalcin                                        /* 39*/ NULL,
1809f4259b30SLisandro Dalcin                                        NULL,
1810f4259b30SLisandro Dalcin                                        NULL,
1811d4002b98SHong Zhang                                        MatGetValues_SeqSELL,
1812d4002b98SHong Zhang                                        MatCopy_SeqSELL,
1813f4259b30SLisandro Dalcin                                        /* 44*/ NULL,
1814d4002b98SHong Zhang                                        MatScale_SeqSELL,
1815d4002b98SHong Zhang                                        MatShift_SeqSELL,
1816f4259b30SLisandro Dalcin                                        NULL,
1817f4259b30SLisandro Dalcin                                        NULL,
1818f4259b30SLisandro Dalcin                                        /* 49*/ NULL,
1819f4259b30SLisandro Dalcin                                        NULL,
1820f4259b30SLisandro Dalcin                                        NULL,
1821f4259b30SLisandro Dalcin                                        NULL,
1822f4259b30SLisandro Dalcin                                        NULL,
1823d4002b98SHong Zhang                                        /* 54*/ MatFDColoringCreate_SeqXAIJ,
1824f4259b30SLisandro Dalcin                                        NULL,
1825f4259b30SLisandro Dalcin                                        NULL,
1826f4259b30SLisandro Dalcin                                        NULL,
1827f4259b30SLisandro Dalcin                                        NULL,
1828f4259b30SLisandro Dalcin                                        /* 59*/ NULL,
1829d4002b98SHong Zhang                                        MatDestroy_SeqSELL,
1830d4002b98SHong Zhang                                        MatView_SeqSELL,
1831f4259b30SLisandro Dalcin                                        NULL,
1832f4259b30SLisandro Dalcin                                        NULL,
1833f4259b30SLisandro Dalcin                                        /* 64*/ NULL,
1834f4259b30SLisandro Dalcin                                        NULL,
1835f4259b30SLisandro Dalcin                                        NULL,
1836f4259b30SLisandro Dalcin                                        NULL,
1837f4259b30SLisandro Dalcin                                        NULL,
1838f4259b30SLisandro Dalcin                                        /* 69*/ NULL,
1839f4259b30SLisandro Dalcin                                        NULL,
1840f4259b30SLisandro Dalcin                                        NULL,
1841f4259b30SLisandro Dalcin                                        NULL,
1842f4259b30SLisandro Dalcin                                        NULL,
1843f4259b30SLisandro Dalcin                                        /* 74*/ NULL,
1844d4002b98SHong Zhang                                        MatFDColoringApply_AIJ, /* reuse the FDColoring function for AIJ */
1845f4259b30SLisandro Dalcin                                        NULL,
1846f4259b30SLisandro Dalcin                                        NULL,
1847f4259b30SLisandro Dalcin                                        NULL,
1848f4259b30SLisandro Dalcin                                        /* 79*/ NULL,
1849f4259b30SLisandro Dalcin                                        NULL,
1850f4259b30SLisandro Dalcin                                        NULL,
1851f4259b30SLisandro Dalcin                                        NULL,
1852f4259b30SLisandro Dalcin                                        NULL,
1853f4259b30SLisandro Dalcin                                        /* 84*/ NULL,
1854f4259b30SLisandro Dalcin                                        NULL,
1855f4259b30SLisandro Dalcin                                        NULL,
1856f4259b30SLisandro Dalcin                                        NULL,
1857f4259b30SLisandro Dalcin                                        NULL,
1858f4259b30SLisandro Dalcin                                        /* 89*/ NULL,
1859f4259b30SLisandro Dalcin                                        NULL,
1860f4259b30SLisandro Dalcin                                        NULL,
1861f4259b30SLisandro Dalcin                                        NULL,
1862f4259b30SLisandro Dalcin                                        NULL,
1863f4259b30SLisandro Dalcin                                        /* 94*/ NULL,
1864f4259b30SLisandro Dalcin                                        NULL,
1865f4259b30SLisandro Dalcin                                        NULL,
1866f4259b30SLisandro Dalcin                                        NULL,
1867f4259b30SLisandro Dalcin                                        NULL,
1868f4259b30SLisandro Dalcin                                        /* 99*/ NULL,
1869f4259b30SLisandro Dalcin                                        NULL,
1870f4259b30SLisandro Dalcin                                        NULL,
1871d4002b98SHong Zhang                                        MatConjugate_SeqSELL,
1872f4259b30SLisandro Dalcin                                        NULL,
1873f4259b30SLisandro Dalcin                                        /*104*/ NULL,
1874f4259b30SLisandro Dalcin                                        NULL,
1875f4259b30SLisandro Dalcin                                        NULL,
1876f4259b30SLisandro Dalcin                                        NULL,
1877f4259b30SLisandro Dalcin                                        NULL,
1878f4259b30SLisandro Dalcin                                        /*109*/ NULL,
1879f4259b30SLisandro Dalcin                                        NULL,
1880f4259b30SLisandro Dalcin                                        NULL,
1881f4259b30SLisandro Dalcin                                        NULL,
1882d4002b98SHong Zhang                                        MatMissingDiagonal_SeqSELL,
1883f4259b30SLisandro Dalcin                                        /*114*/ NULL,
1884f4259b30SLisandro Dalcin                                        NULL,
1885f4259b30SLisandro Dalcin                                        NULL,
1886f4259b30SLisandro Dalcin                                        NULL,
1887f4259b30SLisandro Dalcin                                        NULL,
1888f4259b30SLisandro Dalcin                                        /*119*/ NULL,
1889f4259b30SLisandro Dalcin                                        NULL,
1890f4259b30SLisandro Dalcin                                        NULL,
1891f4259b30SLisandro Dalcin                                        NULL,
1892f4259b30SLisandro Dalcin                                        NULL,
1893f4259b30SLisandro Dalcin                                        /*124*/ NULL,
1894f4259b30SLisandro Dalcin                                        NULL,
1895f4259b30SLisandro Dalcin                                        NULL,
1896f4259b30SLisandro Dalcin                                        NULL,
1897f4259b30SLisandro Dalcin                                        NULL,
1898f4259b30SLisandro Dalcin                                        /*129*/ NULL,
1899f4259b30SLisandro Dalcin                                        NULL,
1900f4259b30SLisandro Dalcin                                        NULL,
1901f4259b30SLisandro Dalcin                                        NULL,
1902f4259b30SLisandro Dalcin                                        NULL,
1903f4259b30SLisandro Dalcin                                        /*134*/ NULL,
1904f4259b30SLisandro Dalcin                                        NULL,
1905f4259b30SLisandro Dalcin                                        NULL,
1906f4259b30SLisandro Dalcin                                        NULL,
1907f4259b30SLisandro Dalcin                                        NULL,
1908f4259b30SLisandro Dalcin                                        /*139*/ NULL,
1909f4259b30SLisandro Dalcin                                        NULL,
1910f4259b30SLisandro Dalcin                                        NULL,
1911d4002b98SHong Zhang                                        MatFDColoringSetUp_SeqXAIJ,
1912f4259b30SLisandro Dalcin                                        NULL,
1913d70f29a3SPierre Jolivet                                        /*144*/ NULL,
1914d70f29a3SPierre Jolivet                                        NULL,
1915d70f29a3SPierre Jolivet                                        NULL,
191699a7f59eSMark Adams                                        NULL,
191799a7f59eSMark Adams                                        NULL,
19187fb60732SBarry Smith                                        NULL,
1919dec0b466SHong Zhang                                        /*150*/ NULL,
1920dec0b466SHong Zhang                                        NULL};
1921d4002b98SHong Zhang 
1922d71ae5a4SJacob Faibussowitsch PetscErrorCode MatStoreValues_SeqSELL(Mat mat)
1923d71ae5a4SJacob Faibussowitsch {
1924d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
1925d4002b98SHong Zhang 
1926d4002b98SHong Zhang   PetscFunctionBegin;
192728b400f6SJacob Faibussowitsch   PetscCheck(a->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1928d4002b98SHong Zhang 
1929d4002b98SHong Zhang   /* allocate space for values if not already there */
1930aa624791SPierre Jolivet   if (!a->saved_values) PetscCall(PetscMalloc1(a->sliidx[a->totalslices] + 1, &a->saved_values));
1931d4002b98SHong Zhang 
1932d4002b98SHong Zhang   /* copy values over */
19339566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->saved_values, a->val, a->sliidx[a->totalslices]));
19343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1935d4002b98SHong Zhang }
1936d4002b98SHong Zhang 
1937d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRetrieveValues_SeqSELL(Mat mat)
1938d71ae5a4SJacob Faibussowitsch {
1939d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
1940d4002b98SHong Zhang 
1941d4002b98SHong Zhang   PetscFunctionBegin;
194228b400f6SJacob Faibussowitsch   PetscCheck(a->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
194328b400f6SJacob Faibussowitsch   PetscCheck(a->saved_values, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatStoreValues(A);first");
19449566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->val, a->saved_values, a->sliidx[a->totalslices]));
19453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1946d4002b98SHong Zhang }
1947d4002b98SHong Zhang 
1948*07e43b41SHong Zhang PetscErrorCode MatSeqSELLGetFillRatio_SeqSELL(Mat mat, PetscReal *ratio)
1949*07e43b41SHong Zhang {
1950*07e43b41SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
1951*07e43b41SHong Zhang 
1952*07e43b41SHong Zhang   PetscFunctionBegin;
1953*07e43b41SHong Zhang   if (a->totalslices && a->sliidx[a->totalslices]) {
1954*07e43b41SHong Zhang     *ratio = (PetscReal)(a->sliidx[a->totalslices] - a->nz) / a->sliidx[a->totalslices];
1955*07e43b41SHong Zhang   } else {
1956*07e43b41SHong Zhang     *ratio = 0.0;
1957*07e43b41SHong Zhang   }
1958*07e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
1959*07e43b41SHong Zhang }
1960*07e43b41SHong Zhang 
1961*07e43b41SHong Zhang PetscErrorCode MatSeqSELLGetMaxSliceWidth_SeqSELL(Mat mat, PetscInt *slicewidth)
1962*07e43b41SHong Zhang {
1963*07e43b41SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
1964*07e43b41SHong Zhang   PetscInt     i, current_slicewidth;
1965*07e43b41SHong Zhang 
1966*07e43b41SHong Zhang   PetscFunctionBegin;
1967*07e43b41SHong Zhang   *slicewidth = 0;
1968*07e43b41SHong Zhang   for (i = 0; i < a->totalslices; i++) {
1969*07e43b41SHong Zhang     current_slicewidth = (a->sliidx[i + 1] - a->sliidx[i]) / a->sliceheight;
1970*07e43b41SHong Zhang     if (current_slicewidth > *slicewidth) *slicewidth = current_slicewidth;
1971*07e43b41SHong Zhang   }
1972*07e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
1973*07e43b41SHong Zhang }
1974*07e43b41SHong Zhang 
1975*07e43b41SHong Zhang PetscErrorCode MatSeqSELLGetAvgSliceWidth_SeqSELL(Mat mat, PetscReal *slicewidth)
1976*07e43b41SHong Zhang {
1977*07e43b41SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
1978*07e43b41SHong Zhang 
1979*07e43b41SHong Zhang   PetscFunctionBegin;
1980*07e43b41SHong Zhang   *slicewidth = 0;
1981*07e43b41SHong Zhang   if (a->totalslices) { *slicewidth = (PetscReal)a->sliidx[a->totalslices] / a->sliceheight / a->totalslices; }
1982*07e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
1983*07e43b41SHong Zhang }
1984*07e43b41SHong Zhang 
1985*07e43b41SHong Zhang PetscErrorCode MatSeqSELLSetSliceHeight_SeqSELL(Mat A, PetscInt sliceheight)
1986*07e43b41SHong Zhang {
1987*07e43b41SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1988*07e43b41SHong Zhang 
1989*07e43b41SHong Zhang   PetscFunctionBegin;
1990*07e43b41SHong Zhang   if (A->preallocated) PetscFunctionReturn(PETSC_SUCCESS);
1991*07e43b41SHong Zhang   PetscCheck(a->sliceheight <= 0 || a->sliceheight == sliceheight, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot change slice height %" PetscInt_FMT " to %" PetscInt_FMT, a->sliceheight, sliceheight);
1992*07e43b41SHong Zhang   a->sliceheight = sliceheight;
1993*07e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
1994*07e43b41SHong Zhang }
1995*07e43b41SHong Zhang 
1996d4002b98SHong Zhang /*@C
199711a5261eSBarry Smith   MatSeqSELLRestoreArray - returns access to the array where the data for a `MATSEQSELL` matrix is stored obtained by `MatSeqSELLGetArray()`
1998d4002b98SHong Zhang 
1999d4002b98SHong Zhang   Not Collective
2000d4002b98SHong Zhang 
2001d4002b98SHong Zhang   Input Parameters:
2002*07e43b41SHong Zhang +  A - a `MATSEQSELL` matrix
200320f4b53cSBarry Smith -  array - pointer to the data
2004d4002b98SHong Zhang 
2005d4002b98SHong Zhang   Level: intermediate
2006d4002b98SHong Zhang 
200767be906fSBarry Smith .seealso: `Mat`, `MATSEQSELL`, `MatSeqSELLGetArray()`, `MatSeqSELLRestoreArrayF90()`
2008d4002b98SHong Zhang @*/
2009d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLRestoreArray(Mat A, PetscScalar **array)
2010d71ae5a4SJacob Faibussowitsch {
2011d4002b98SHong Zhang   PetscFunctionBegin;
2012cac4c232SBarry Smith   PetscUseMethod(A, "MatSeqSELLRestoreArray_C", (Mat, PetscScalar **), (A, array));
20133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2014d4002b98SHong Zhang }
2015d4002b98SHong Zhang 
2016*07e43b41SHong Zhang /*@C
2017*07e43b41SHong Zhang   MatSeqSELLGetFillRatio - returns a ratio that indicates the irregularity of the matrix.
2018*07e43b41SHong Zhang 
2019*07e43b41SHong Zhang   Not Collective
2020*07e43b41SHong Zhang 
2021*07e43b41SHong Zhang   Input Parameter:
2022*07e43b41SHong Zhang .  A - a MATSEQSELL matrix
2023*07e43b41SHong Zhang 
2024*07e43b41SHong Zhang   Output Parameter:
2025*07e43b41SHong Zhang .  ratio - ratio of number of padded zeros to number of allocated elements
2026*07e43b41SHong Zhang 
2027*07e43b41SHong Zhang   Level: intermediate
2028*07e43b41SHong Zhang @*/
2029*07e43b41SHong Zhang PetscErrorCode MatSeqSELLGetFillRatio(Mat A, PetscReal *ratio)
2030*07e43b41SHong Zhang {
2031*07e43b41SHong Zhang   PetscFunctionBegin;
2032*07e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLGetFillRatio_C", (Mat, PetscReal *), (A, ratio));
2033*07e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
2034*07e43b41SHong Zhang }
2035*07e43b41SHong Zhang 
2036*07e43b41SHong Zhang /*@C
2037*07e43b41SHong Zhang   MatSeqSELLGetMaxSliceWidth - returns the maximum slice width.
2038*07e43b41SHong Zhang 
2039*07e43b41SHong Zhang   Not Collective
2040*07e43b41SHong Zhang 
2041*07e43b41SHong Zhang   Input Parameter:
2042*07e43b41SHong Zhang .  A - a MATSEQSELL matrix
2043*07e43b41SHong Zhang 
2044*07e43b41SHong Zhang   Output Parameter:
2045*07e43b41SHong Zhang .  slicewidth - maximum slice width
2046*07e43b41SHong Zhang 
2047*07e43b41SHong Zhang  Level: intermediate
2048*07e43b41SHong Zhang @*/
2049*07e43b41SHong Zhang PetscErrorCode MatSeqSELLGetMaxSliceWidth(Mat A, PetscInt *slicewidth)
2050*07e43b41SHong Zhang {
2051*07e43b41SHong Zhang   PetscFunctionBegin;
2052*07e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLGetMaxSliceWidth_C", (Mat, PetscInt *), (A, slicewidth));
2053*07e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
2054*07e43b41SHong Zhang }
2055*07e43b41SHong Zhang 
2056*07e43b41SHong Zhang /*@C
2057*07e43b41SHong Zhang   MatSeqSELLGetAvgSliceWidth - returns the average slice width.
2058*07e43b41SHong Zhang 
2059*07e43b41SHong Zhang   Not Collective
2060*07e43b41SHong Zhang 
2061*07e43b41SHong Zhang   Input Parameter:
2062*07e43b41SHong Zhang .  A - a MATSEQSELL matrix
2063*07e43b41SHong Zhang 
2064*07e43b41SHong Zhang   Output Parameter:
2065*07e43b41SHong Zhang .  slicewidth - average slice width
2066*07e43b41SHong Zhang 
2067*07e43b41SHong Zhang   Level: intermediate
2068*07e43b41SHong Zhang @*/
2069*07e43b41SHong Zhang PetscErrorCode MatSeqSELLGetAvgSliceWidth(Mat A, PetscReal *slicewidth)
2070*07e43b41SHong Zhang {
2071*07e43b41SHong Zhang   PetscFunctionBegin;
2072*07e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLGetAvgSliceWidth_C", (Mat, PetscReal *), (A, slicewidth));
2073*07e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
2074*07e43b41SHong Zhang }
2075*07e43b41SHong Zhang 
2076*07e43b41SHong Zhang /*@C
2077*07e43b41SHong Zhang   MatSeqSELLSetSliceHeight - sets the slice height.
2078*07e43b41SHong Zhang 
2079*07e43b41SHong Zhang   Not Collective
2080*07e43b41SHong Zhang 
2081*07e43b41SHong Zhang   Input Parameters:
2082*07e43b41SHong Zhang +  A - a MATSEQSELL matrix
2083*07e43b41SHong Zhang -  sliceheight - slice height
2084*07e43b41SHong Zhang 
2085*07e43b41SHong Zhang   Notes:
2086*07e43b41SHong Zhang   You cannot change the slice height once it have been set.
2087*07e43b41SHong Zhang 
2088*07e43b41SHong Zhang   The slice height must be set before MatSetUp() or MatXXXSetPreallocation() is called.
2089*07e43b41SHong Zhang 
2090*07e43b41SHong Zhang   Level: intermediate
2091*07e43b41SHong Zhang @*/
2092*07e43b41SHong Zhang PetscErrorCode MatSeqSELLSetSliceHeight(Mat A, PetscInt sliceheight)
2093*07e43b41SHong Zhang {
2094*07e43b41SHong Zhang   PetscFunctionBegin;
2095*07e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLSetSliceHeight_C", (Mat, PetscInt), (A, sliceheight));
2096*07e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
2097*07e43b41SHong Zhang }
2098*07e43b41SHong Zhang 
2099*07e43b41SHong Zhang /*@C
2100*07e43b41SHong Zhang   MatSeqSELLGetVarSliceSize - returns the variance of the slice size.
2101*07e43b41SHong Zhang 
2102*07e43b41SHong Zhang   Not Collective
2103*07e43b41SHong Zhang 
2104*07e43b41SHong Zhang   Input Parameter:
2105*07e43b41SHong Zhang .  A - a MATSEQSELL matrix
2106*07e43b41SHong Zhang 
2107*07e43b41SHong Zhang   Output Parameter:
2108*07e43b41SHong Zhang .  variance - variance of the slice size
2109*07e43b41SHong Zhang 
2110*07e43b41SHong Zhang   Level: intermediate
2111*07e43b41SHong Zhang @*/
2112*07e43b41SHong Zhang PetscErrorCode MatSeqSELLGetVarSliceSize(Mat A, PetscReal *variance)
2113*07e43b41SHong Zhang {
2114*07e43b41SHong Zhang   PetscFunctionBegin;
2115*07e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLGetVarSliceSize_C", (Mat, PetscReal *), (A, variance));
2116*07e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
2117*07e43b41SHong Zhang }
2118*07e43b41SHong Zhang 
21192d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
21202d1451d4SHong Zhang PETSC_EXTERN PetscErrorCode MatConvert_SeqSELL_SeqSELLCUDA(Mat);
21212d1451d4SHong Zhang #endif
21222d1451d4SHong Zhang 
2123d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_SeqSELL(Mat B)
2124d71ae5a4SJacob Faibussowitsch {
2125d4002b98SHong Zhang   Mat_SeqSELL *b;
2126d4002b98SHong Zhang   PetscMPIInt  size;
2127d4002b98SHong Zhang 
2128d4002b98SHong Zhang   PetscFunctionBegin;
21299566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
21309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
213108401ef6SPierre Jolivet   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Comm must be of size 1");
2132d4002b98SHong Zhang 
21334dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
2134d4002b98SHong Zhang 
2135d4002b98SHong Zhang   B->data = (void *)b;
2136d4002b98SHong Zhang 
21379566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(B->ops, &MatOps_Values, sizeof(struct _MatOps)));
2138d4002b98SHong Zhang 
2139f4259b30SLisandro Dalcin   b->row                = NULL;
2140f4259b30SLisandro Dalcin   b->col                = NULL;
2141f4259b30SLisandro Dalcin   b->icol               = NULL;
2142d4002b98SHong Zhang   b->reallocs           = 0;
2143d4002b98SHong Zhang   b->ignorezeroentries  = PETSC_FALSE;
2144d4002b98SHong Zhang   b->roworiented        = PETSC_TRUE;
2145d4002b98SHong Zhang   b->nonew              = 0;
2146f4259b30SLisandro Dalcin   b->diag               = NULL;
2147f4259b30SLisandro Dalcin   b->solve_work         = NULL;
2148f4259b30SLisandro Dalcin   B->spptr              = NULL;
2149f4259b30SLisandro Dalcin   b->saved_values       = NULL;
2150f4259b30SLisandro Dalcin   b->idiag              = NULL;
2151f4259b30SLisandro Dalcin   b->mdiag              = NULL;
2152f4259b30SLisandro Dalcin   b->ssor_work          = NULL;
2153d4002b98SHong Zhang   b->omega              = 1.0;
2154d4002b98SHong Zhang   b->fshift             = 0.0;
2155d4002b98SHong Zhang   b->idiagvalid         = PETSC_FALSE;
2156d4002b98SHong Zhang   b->keepnonzeropattern = PETSC_FALSE;
2157*07e43b41SHong Zhang   b->sliceheight        = 0;
2158d4002b98SHong Zhang 
21599566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQSELL));
21609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetArray_C", MatSeqSELLGetArray_SeqSELL));
21619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLRestoreArray_C", MatSeqSELLRestoreArray_SeqSELL));
21629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_SeqSELL));
21639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_SeqSELL));
21649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLSetPreallocation_C", MatSeqSELLSetPreallocation_SeqSELL));
21652d1451d4SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_SeqSELL_SeqAIJ_C", MatConvert_SeqSELL_SeqAIJ));
21662d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
21672d1451d4SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_SeqSELL_SeqSELLCUDA_C", MatConvert_SeqSELL_SeqSELLCUDA));
21682d1451d4SHong Zhang #endif
2169*07e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetFillRatio_C", MatSeqSELLGetFillRatio_SeqSELL));
2170*07e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetMaxSliceWidth_C", MatSeqSELLGetMaxSliceWidth_SeqSELL));
2171*07e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetAvgSliceWidth_C", MatSeqSELLGetAvgSliceWidth_SeqSELL));
2172*07e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLSetSliceHeight_C", MatSeqSELLSetSliceHeight_SeqSELL));
2173*07e43b41SHong Zhang 
2174*07e43b41SHong Zhang   PetscObjectOptionsBegin((PetscObject)B);
2175*07e43b41SHong Zhang   {
2176*07e43b41SHong Zhang     PetscInt  newsh = -1;
2177*07e43b41SHong Zhang     PetscBool flg;
2178*07e43b41SHong Zhang 
2179*07e43b41SHong Zhang     PetscCall(PetscOptionsInt("-mat_sell_slice_height", "Set the slice height used to store SELL matrix", "MatSELLSetSliceHeight", newsh, &newsh, &flg));
2180*07e43b41SHong Zhang     if (flg) { PetscCall(MatSeqSELLSetSliceHeight(B, newsh)); }
2181*07e43b41SHong Zhang   }
2182*07e43b41SHong Zhang   PetscOptionsEnd();
21833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2184d4002b98SHong Zhang }
2185d4002b98SHong Zhang 
2186d4002b98SHong Zhang /*
2187d4002b98SHong Zhang  Given a matrix generated with MatGetFactor() duplicates all the information in A into B
2188d4002b98SHong Zhang  */
2189d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicateNoCreate_SeqSELL(Mat C, Mat A, MatDuplicateOption cpvalues, PetscBool mallocmatspace)
2190d71ae5a4SJacob Faibussowitsch {
2191ed73aabaSBarry Smith   Mat_SeqSELL *c = (Mat_SeqSELL *)C->data, *a = (Mat_SeqSELL *)A->data;
2192d4002b98SHong Zhang   PetscInt     i, m                           = A->rmap->n;
2193d4002b98SHong Zhang   PetscInt     totalslices = a->totalslices;
2194d4002b98SHong Zhang 
2195d4002b98SHong Zhang   PetscFunctionBegin;
2196d4002b98SHong Zhang   C->factortype = A->factortype;
2197f4259b30SLisandro Dalcin   c->row        = NULL;
2198f4259b30SLisandro Dalcin   c->col        = NULL;
2199f4259b30SLisandro Dalcin   c->icol       = NULL;
2200d4002b98SHong Zhang   c->reallocs   = 0;
2201d4002b98SHong Zhang   C->assembled  = PETSC_TRUE;
2202d4002b98SHong Zhang 
22039566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->rmap, &C->rmap));
22049566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &C->cmap));
2205d4002b98SHong Zhang 
2206*07e43b41SHong Zhang   PetscCall(PetscMalloc1(a->sliceheight * totalslices, &c->rlen));
22079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(totalslices + 1, &c->sliidx));
2208d4002b98SHong Zhang 
2209d4002b98SHong Zhang   for (i = 0; i < m; i++) c->rlen[i] = a->rlen[i];
2210d4002b98SHong Zhang   for (i = 0; i < totalslices + 1; i++) c->sliidx[i] = a->sliidx[i];
2211d4002b98SHong Zhang 
2212d4002b98SHong Zhang   /* allocate the matrix space */
2213d4002b98SHong Zhang   if (mallocmatspace) {
22149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(a->maxallocmat, &c->val, a->maxallocmat, &c->colidx));
2215d4002b98SHong Zhang 
2216d4002b98SHong Zhang     c->singlemalloc = PETSC_TRUE;
2217d4002b98SHong Zhang 
2218d4002b98SHong Zhang     if (m > 0) {
22199566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->colidx, a->colidx, a->maxallocmat));
2220d4002b98SHong Zhang       if (cpvalues == MAT_COPY_VALUES) {
22219566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(c->val, a->val, a->maxallocmat));
2222d4002b98SHong Zhang       } else {
22239566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(c->val, a->maxallocmat));
2224d4002b98SHong Zhang       }
2225d4002b98SHong Zhang     }
2226d4002b98SHong Zhang   }
2227d4002b98SHong Zhang 
2228d4002b98SHong Zhang   c->ignorezeroentries = a->ignorezeroentries;
2229d4002b98SHong Zhang   c->roworiented       = a->roworiented;
2230d4002b98SHong Zhang   c->nonew             = a->nonew;
2231d4002b98SHong Zhang   if (a->diag) {
22329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &c->diag));
2233ad540459SPierre Jolivet     for (i = 0; i < m; i++) c->diag[i] = a->diag[i];
2234f4259b30SLisandro Dalcin   } else c->diag = NULL;
2235d4002b98SHong Zhang 
2236f4259b30SLisandro Dalcin   c->solve_work         = NULL;
2237f4259b30SLisandro Dalcin   c->saved_values       = NULL;
2238f4259b30SLisandro Dalcin   c->idiag              = NULL;
2239f4259b30SLisandro Dalcin   c->ssor_work          = NULL;
2240d4002b98SHong Zhang   c->keepnonzeropattern = a->keepnonzeropattern;
2241d4002b98SHong Zhang   c->free_val           = PETSC_TRUE;
2242d4002b98SHong Zhang   c->free_colidx        = PETSC_TRUE;
2243d4002b98SHong Zhang 
2244d4002b98SHong Zhang   c->maxallocmat  = a->maxallocmat;
2245d4002b98SHong Zhang   c->maxallocrow  = a->maxallocrow;
2246d4002b98SHong Zhang   c->rlenmax      = a->rlenmax;
2247d4002b98SHong Zhang   c->nz           = a->nz;
2248d4002b98SHong Zhang   C->preallocated = PETSC_TRUE;
2249d4002b98SHong Zhang 
2250d4002b98SHong Zhang   c->nonzerorowcnt = a->nonzerorowcnt;
2251d4002b98SHong Zhang   C->nonzerostate  = A->nonzerostate;
2252d4002b98SHong Zhang 
22539566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist, &((PetscObject)C)->qlist));
22543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2255d4002b98SHong Zhang }
2256d4002b98SHong Zhang 
2257d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_SeqSELL(Mat A, MatDuplicateOption cpvalues, Mat *B)
2258d71ae5a4SJacob Faibussowitsch {
2259d4002b98SHong Zhang   PetscFunctionBegin;
22609566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
22619566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n));
226248a46eb9SPierre Jolivet   if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) PetscCall(MatSetBlockSizesFromMats(*B, A, A));
22639566063dSJacob Faibussowitsch   PetscCall(MatSetType(*B, ((PetscObject)A)->type_name));
22649566063dSJacob Faibussowitsch   PetscCall(MatDuplicateNoCreate_SeqSELL(*B, A, cpvalues, PETSC_TRUE));
22653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2266d4002b98SHong Zhang }
2267d4002b98SHong Zhang 
2268ed73aabaSBarry Smith /*MC
2269ed73aabaSBarry Smith    MATSEQSELL - MATSEQSELL = "seqsell" - A matrix type to be used for sequential sparse matrices,
2270ed73aabaSBarry Smith    based on the sliced Ellpack format
2271ed73aabaSBarry Smith 
227220f4b53cSBarry Smith    Options Database Key:
227311a5261eSBarry Smith . -mat_type seqsell - sets the matrix type to "`MATSEQELL` during a call to `MatSetFromOptions()`
2274ed73aabaSBarry Smith 
2275ed73aabaSBarry Smith    Level: beginner
2276ed73aabaSBarry Smith 
227767be906fSBarry Smith .seealso: `Mat`, `MatCreateSeqSell()`, `MATSELL`, `MATMPISELL`, `MATSEQAIJ`, `MATAIJ`, `MATMPIAIJ`
2278ed73aabaSBarry Smith M*/
2279ed73aabaSBarry Smith 
2280ed73aabaSBarry Smith /*MC
2281ed73aabaSBarry Smith    MATSELL - MATSELL = "sell" - A matrix type to be used for sparse matrices.
2282ed73aabaSBarry Smith 
228311a5261eSBarry Smith    This matrix type is identical to `MATSEQSELL` when constructed with a single process communicator,
228411a5261eSBarry Smith    and `MATMPISELL` otherwise.  As a result, for single process communicators,
228511a5261eSBarry Smith   `MatSeqSELLSetPreallocation()` is supported, and similarly `MatMPISELLSetPreallocation()` is supported
2286ed73aabaSBarry Smith   for communicators controlling multiple processes.  It is recommended that you call both of
2287ed73aabaSBarry Smith   the above preallocation routines for simplicity.
2288ed73aabaSBarry Smith 
228920f4b53cSBarry Smith    Options Database Key:
2290ed73aabaSBarry Smith . -mat_type sell - sets the matrix type to "sell" during a call to MatSetFromOptions()
2291ed73aabaSBarry Smith 
2292ed73aabaSBarry Smith   Level: beginner
2293ed73aabaSBarry Smith 
2294ed73aabaSBarry Smith   Notes:
2295ed73aabaSBarry Smith    This format is only supported for real scalars, double precision, and 32 bit indices (the defaults).
2296ed73aabaSBarry Smith 
2297ed73aabaSBarry Smith    It can provide better performance on Intel and AMD processes with AVX2 or AVX512 support for matrices that have a similar number of
2298ed73aabaSBarry Smith    non-zeros in contiguous groups of rows. However if the computation is memory bandwidth limited it may not provide much improvement.
2299ed73aabaSBarry Smith 
2300ed73aabaSBarry Smith   Developer Notes:
2301ed73aabaSBarry Smith    On Intel (and AMD) systems some of the matrix operations use SIMD (AVX) instructions to achieve higher performance.
2302ed73aabaSBarry Smith 
2303ed73aabaSBarry Smith    The sparse matrix format is as follows. For simplicity we assume a slice size of 2, it is actually 8
2304ed73aabaSBarry Smith .vb
2305ed73aabaSBarry Smith                             (2 0  3 4)
2306ed73aabaSBarry Smith    Consider the matrix A =  (5 0  6 0)
2307ed73aabaSBarry Smith                             (0 0  7 8)
2308ed73aabaSBarry Smith                             (0 0  9 9)
2309ed73aabaSBarry Smith 
2310ed73aabaSBarry Smith    symbolically the Ellpack format can be written as
2311ed73aabaSBarry Smith 
2312ed73aabaSBarry Smith         (2 3 4 |)           (0 2 3 |)
2313ed73aabaSBarry Smith    v =  (5 6 0 |)  colidx = (0 2 2 |)
2314ed73aabaSBarry Smith         --------            ---------
2315ed73aabaSBarry Smith         (7 8 |)             (2 3 |)
2316ed73aabaSBarry Smith         (9 9 |)             (2 3 |)
2317ed73aabaSBarry Smith 
2318ed73aabaSBarry Smith     The data for 2 contiguous rows of the matrix are stored together (in column-major format) (with any left-over rows handled as a special case).
2319ed73aabaSBarry Smith     Any of the rows in a slice fewer columns than the rest of the slice (row 1 above) are padded with a previous valid column in their "extra" colidx[] locations and
2320ed73aabaSBarry Smith     zeros in their "extra" v locations so that the matrix operations do not need special code to handle different length rows within the 2 rows in a slice.
2321ed73aabaSBarry Smith 
2322ed73aabaSBarry Smith     The one-dimensional representation of v used in the code is (2 5 3 6 4 0 7 9 8 9)  and for colidx is (0 0 2 2 3 2 2 2 3 3)
2323ed73aabaSBarry Smith 
2324ed73aabaSBarry Smith .ve
2325ed73aabaSBarry Smith 
2326ed73aabaSBarry Smith       See MatMult_SeqSELL() for how this format is used with the SIMD operations to achieve high performance.
2327ed73aabaSBarry Smith 
2328ed73aabaSBarry Smith  References:
2329606c0280SSatish Balay . * - Hong Zhang, Richard T. Mills, Karl Rupp, and Barry F. Smith, Vectorized Parallel Sparse Matrix-Vector Multiplication in {PETSc} Using {AVX-512},
2330ed73aabaSBarry Smith    Proceedings of the 47th International Conference on Parallel Processing, 2018.
2331ed73aabaSBarry Smith 
233267be906fSBarry Smith .seealso: `Mat`, `MatCreateSeqSELL()`, `MatCreateSeqAIJ()`, `MatCreateSell()`, `MATSEQSELL`, `MATMPISELL`, `MATSEQAIJ`, `MATMPIAIJ`, `MATAIJ`
2333ed73aabaSBarry Smith M*/
2334ed73aabaSBarry Smith 
2335d4002b98SHong Zhang /*@C
233611a5261eSBarry Smith        MatCreateSeqSELL - Creates a sparse matrix in `MATSEQSELL` format.
2337d4002b98SHong Zhang 
23382ef1f0ffSBarry Smith  Collective
2339d4002b98SHong Zhang 
2340d4002b98SHong Zhang  Input Parameters:
234111a5261eSBarry Smith +  comm - MPI communicator, set to `PETSC_COMM_SELF`
2342d4002b98SHong Zhang .  m - number of rows
2343d4002b98SHong Zhang .  n - number of columns
234420f4b53cSBarry Smith .  rlenmax - maximum number of nonzeros in a row, ignored if `rlen` is provided
234520f4b53cSBarry Smith -  rlen - array containing the number of nonzeros in the various rows (possibly different for each row) or NULL
2346d4002b98SHong Zhang 
2347d4002b98SHong Zhang  Output Parameter:
2348d4002b98SHong Zhang .  A - the matrix
2349d4002b98SHong Zhang 
235020f4b53cSBarry Smith  Level: intermediate
235120f4b53cSBarry Smith 
235220f4b53cSBarry Smith  Notes:
235311a5261eSBarry Smith  It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`,
2354f6f02116SRichard Tran Mills  MatXXXXSetPreallocation() paradigm instead of this routine directly.
235511a5261eSBarry Smith  [MatXXXXSetPreallocation() is, for example, `MatSeqSELLSetPreallocation()`]
2356d4002b98SHong Zhang 
235720f4b53cSBarry Smith  Specify the preallocated storage with either `rlenmax` or `rlen` (not both).
235820f4b53cSBarry Smith  Set `rlenmax` = `PETSC_DEFAULT` and `rlen` = `NULL` for PETSc to control dynamic memory
235920f4b53cSBarry Smith  allocation.
2360d4002b98SHong Zhang 
236167be906fSBarry Smith  .seealso: `Mat`, `MATSEQSELL`, `MatCreate()`, `MatCreateSELL()`, `MatSetValues()`, `MatSeqSELLSetPreallocation()`, `MATSELL`, `MATSEQSELL`, `MATMPISELL`
2362d4002b98SHong Zhang  @*/
236320f4b53cSBarry Smith PetscErrorCode MatCreateSeqSELL(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt rlenmax, const PetscInt rlen[], Mat *A)
2364d71ae5a4SJacob Faibussowitsch {
2365d4002b98SHong Zhang   PetscFunctionBegin;
23669566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
23679566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, n, m, n));
23689566063dSJacob Faibussowitsch   PetscCall(MatSetType(*A, MATSEQSELL));
236920f4b53cSBarry Smith   PetscCall(MatSeqSELLSetPreallocation_SeqSELL(*A, rlenmax, rlen));
23703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2371d4002b98SHong Zhang }
2372d4002b98SHong Zhang 
2373d71ae5a4SJacob Faibussowitsch PetscErrorCode MatEqual_SeqSELL(Mat A, Mat B, PetscBool *flg)
2374d71ae5a4SJacob Faibussowitsch {
2375d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data, *b = (Mat_SeqSELL *)B->data;
2376d4002b98SHong Zhang   PetscInt     totalslices = a->totalslices;
2377d4002b98SHong Zhang 
2378d4002b98SHong Zhang   PetscFunctionBegin;
2379d4002b98SHong Zhang   /* If the  matrix dimensions are not equal,or no of nonzeros */
2380d4002b98SHong Zhang   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) || (a->nz != b->nz) || (a->rlenmax != b->rlenmax)) {
2381d4002b98SHong Zhang     *flg = PETSC_FALSE;
23823ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2383d4002b98SHong Zhang   }
2384d4002b98SHong Zhang   /* if the a->colidx are the same */
23859566063dSJacob Faibussowitsch   PetscCall(PetscArraycmp(a->colidx, b->colidx, a->sliidx[totalslices], flg));
23863ba16761SJacob Faibussowitsch   if (!*flg) PetscFunctionReturn(PETSC_SUCCESS);
2387d4002b98SHong Zhang   /* if a->val are the same */
23889566063dSJacob Faibussowitsch   PetscCall(PetscArraycmp(a->val, b->val, a->sliidx[totalslices], flg));
23893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2390d4002b98SHong Zhang }
2391d4002b98SHong Zhang 
2392d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLInvalidateDiagonal(Mat A)
2393d71ae5a4SJacob Faibussowitsch {
2394d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
2395d4002b98SHong Zhang 
2396d4002b98SHong Zhang   PetscFunctionBegin;
2397d4002b98SHong Zhang   a->idiagvalid = PETSC_FALSE;
23983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2399d4002b98SHong Zhang }
2400d4002b98SHong Zhang 
2401d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConjugate_SeqSELL(Mat A)
2402d71ae5a4SJacob Faibussowitsch {
2403d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
2404d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
2405d4002b98SHong Zhang   PetscInt     i;
2406d4002b98SHong Zhang   PetscScalar *val = a->val;
2407d4002b98SHong Zhang 
2408d4002b98SHong Zhang   PetscFunctionBegin;
24092d1451d4SHong Zhang   for (i = 0; i < a->sliidx[a->totalslices]; i++) { val[i] = PetscConj(val[i]); }
24102d1451d4SHong Zhang   #if defined(PETSC_HAVE_CUDA)
24112d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
24122d1451d4SHong Zhang   #endif
2413d4002b98SHong Zhang #else
2414d4002b98SHong Zhang   PetscFunctionBegin;
2415d4002b98SHong Zhang #endif
24163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2417d4002b98SHong Zhang }
2418