xref: /petsc/src/mat/impls/sell/seq/sell.c (revision b921024ee4b1b8e0b257f8784c7a4ecadbc34b4e)
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;
924e58db63SHong Zhang #if defined(PETSC_HAVE_CUDA)
934e58db63SHong Zhang   PetscInt rlenmax = 0;
944e58db63SHong Zhang #endif
95d4002b98SHong Zhang   PetscBool skipallocation = PETSC_FALSE, realalloc = PETSC_FALSE;
96d4002b98SHong Zhang 
97d4002b98SHong Zhang   PetscFunctionBegin;
98d4002b98SHong Zhang   if (maxallocrow >= 0 || rlen) realalloc = PETSC_TRUE;
99d4002b98SHong Zhang   if (maxallocrow == MAT_SKIP_ALLOCATION) {
100d4002b98SHong Zhang     skipallocation = PETSC_TRUE;
101d4002b98SHong Zhang     maxallocrow    = 0;
102d4002b98SHong Zhang   }
103d4002b98SHong Zhang 
1049566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
1059566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
106d4002b98SHong Zhang 
107d4002b98SHong Zhang   /* FIXME: if one preallocates more space than needed, the matrix does not shrink automatically, but for best performance it should */
108d4002b98SHong Zhang   if (maxallocrow == PETSC_DEFAULT || maxallocrow == PETSC_DECIDE) maxallocrow = 5;
10908401ef6SPierre Jolivet   PetscCheck(maxallocrow >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "maxallocrow cannot be less than 0: value %" PetscInt_FMT, maxallocrow);
110d4002b98SHong Zhang   if (rlen) {
111d4002b98SHong Zhang     for (i = 0; i < B->rmap->n; i++) {
11208401ef6SPierre 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]);
11308401ef6SPierre 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);
114d4002b98SHong Zhang     }
115d4002b98SHong Zhang   }
116d4002b98SHong Zhang 
117d4002b98SHong Zhang   B->preallocated = PETSC_TRUE;
118d4002b98SHong Zhang 
119d4002b98SHong Zhang   b = (Mat_SeqSELL *)B->data;
120d4002b98SHong Zhang 
12107e43b41SHong Zhang   if (!b->sliceheight) { /* not set yet */
12207e43b41SHong Zhang #if defined(PETSC_HAVE_CUDA)
12307e43b41SHong Zhang     b->sliceheight = 16;
12407e43b41SHong Zhang #else
12507e43b41SHong Zhang     b->sliceheight = 8;
12607e43b41SHong Zhang #endif
12707e43b41SHong Zhang   }
12807e43b41SHong Zhang   totalslices    = PetscCeilInt(B->rmap->n, b->sliceheight);
129d4002b98SHong Zhang   b->totalslices = totalslices;
130d4002b98SHong Zhang   if (!skipallocation) {
13107e43b41SHong 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));
132d4002b98SHong Zhang 
133d4002b98SHong Zhang     if (!b->sliidx) { /* sliidx gives the starting index of each slice, the last element is the total space allocated */
1349566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(totalslices + 1, &b->sliidx));
135d4002b98SHong Zhang     }
136d4002b98SHong Zhang     if (!rlen) { /* if rlen is not provided, allocate same space for all the slices */
137d4002b98SHong Zhang       if (maxallocrow == PETSC_DEFAULT || maxallocrow == PETSC_DECIDE) maxallocrow = 10;
138d4002b98SHong Zhang       else if (maxallocrow < 0) maxallocrow = 1;
1394e58db63SHong Zhang #if defined(PETSC_HAVE_CUDA)
1404e58db63SHong Zhang       rlenmax = maxallocrow;
1414e58db63SHong Zhang       /* Pad the slice to DEVICE_MEM_ALIGN */
1424e58db63SHong Zhang       while (b->sliceheight * maxallocrow % DEVICE_MEM_ALIGN) maxallocrow++;
1434e58db63SHong Zhang #endif
14407e43b41SHong Zhang       for (i = 0; i <= totalslices; i++) b->sliidx[i] = b->sliceheight * i * maxallocrow;
145d4002b98SHong Zhang     } else {
1464e58db63SHong Zhang #if defined(PETSC_HAVE_CUDA)
1474e58db63SHong Zhang       PetscInt mul = DEVICE_MEM_ALIGN / b->sliceheight;
1484e58db63SHong Zhang #endif
149d4002b98SHong Zhang       maxallocrow  = 0;
150d4002b98SHong Zhang       b->sliidx[0] = 0;
151d4002b98SHong Zhang       for (i = 1; i < totalslices; i++) {
152d4002b98SHong Zhang         b->sliidx[i] = 0;
15307e43b41SHong Zhang         for (j = 0; j < b->sliceheight; j++) { b->sliidx[i] = PetscMax(b->sliidx[i], rlen[b->sliceheight * (i - 1) + j]); }
1544e58db63SHong Zhang #if defined(PETSC_HAVE_CUDA)
1554e58db63SHong Zhang         rlenmax = PetscMax(b->sliidx[i], rlenmax);
1564e58db63SHong Zhang         /* Pad the slice to DEVICE_MEM_ALIGN */
1574e58db63SHong Zhang         b->sliidx[i] = ((b->sliidx[i] - 1) / mul + 1) * mul;
1584e58db63SHong Zhang #endif
159d4002b98SHong Zhang         maxallocrow = PetscMax(b->sliidx[i], maxallocrow);
16007e43b41SHong Zhang         PetscCall(PetscIntSumError(b->sliidx[i - 1], b->sliceheight * b->sliidx[i], &b->sliidx[i]));
161d4002b98SHong Zhang       }
162d4002b98SHong Zhang       /* last slice */
163d4002b98SHong Zhang       b->sliidx[totalslices] = 0;
16407e43b41SHong Zhang       for (j = b->sliceheight * (totalslices - 1); j < B->rmap->n; j++) b->sliidx[totalslices] = PetscMax(b->sliidx[totalslices], rlen[j]);
1654e58db63SHong Zhang #if defined(PETSC_HAVE_CUDA)
1664e58db63SHong Zhang       rlenmax                = PetscMax(b->sliidx[i], rlenmax);
1674e58db63SHong Zhang       b->sliidx[totalslices] = ((b->sliidx[totalslices] - 1) / mul + 1) * mul;
1684e58db63SHong Zhang #endif
169d4002b98SHong Zhang       maxallocrow            = PetscMax(b->sliidx[totalslices], maxallocrow);
17007e43b41SHong Zhang       b->sliidx[totalslices] = b->sliidx[totalslices - 1] + b->sliceheight * b->sliidx[totalslices];
171d4002b98SHong Zhang     }
172d4002b98SHong Zhang 
173d4002b98SHong Zhang     /* allocate space for val, colidx, rlen */
174d4002b98SHong Zhang     /* FIXME: should B's old memory be unlogged? */
1759566063dSJacob Faibussowitsch     PetscCall(MatSeqXSELLFreeSELL(B, &b->val, &b->colidx));
176d4002b98SHong Zhang     /* FIXME: assuming an element of the bit array takes 8 bits */
1779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(b->sliidx[totalslices], &b->val, b->sliidx[totalslices], &b->colidx));
178d4002b98SHong 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. */
17907e43b41SHong Zhang     PetscCall(PetscCalloc1(b->sliceheight * totalslices, &b->rlen));
180d4002b98SHong Zhang 
181d4002b98SHong Zhang     b->singlemalloc = PETSC_TRUE;
182d4002b98SHong Zhang     b->free_val     = PETSC_TRUE;
183d4002b98SHong Zhang     b->free_colidx  = PETSC_TRUE;
184d4002b98SHong Zhang   } else {
185d4002b98SHong Zhang     b->free_val    = PETSC_FALSE;
186d4002b98SHong Zhang     b->free_colidx = PETSC_FALSE;
187d4002b98SHong Zhang   }
188d4002b98SHong Zhang 
189d4002b98SHong Zhang   b->nz          = 0;
190d4002b98SHong Zhang   b->maxallocrow = maxallocrow;
1914e58db63SHong Zhang #if defined(PETSC_HAVE_CUDA)
1924e58db63SHong Zhang   b->rlenmax = rlenmax;
1934e58db63SHong Zhang #else
194d4002b98SHong Zhang   b->rlenmax = maxallocrow;
1954e58db63SHong Zhang #endif
196d4002b98SHong Zhang   b->maxallocmat      = b->sliidx[totalslices];
197d4002b98SHong Zhang   B->info.nz_unneeded = (double)b->maxallocmat;
1981baa6e33SBarry Smith   if (realalloc) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
200d4002b98SHong Zhang }
201d4002b98SHong Zhang 
202d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetRow_SeqSELL(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
203d71ae5a4SJacob Faibussowitsch {
2046108893eSStefano Zampini   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
2056108893eSStefano Zampini   PetscInt     shift;
2066108893eSStefano Zampini 
2076108893eSStefano Zampini   PetscFunctionBegin;
208aed4548fSBarry Smith   PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " out of range", row);
2096108893eSStefano Zampini   if (nz) *nz = a->rlen[row];
21007e43b41SHong Zhang   shift = a->sliidx[row / a->sliceheight] + (row % a->sliceheight);
2112d1451d4SHong Zhang   if (!a->getrowcols) { PetscCall(PetscMalloc2(a->rlenmax, &a->getrowcols, a->rlenmax, &a->getrowvals)); }
2126108893eSStefano Zampini   if (idx) {
2136108893eSStefano Zampini     PetscInt j;
21407e43b41SHong Zhang     for (j = 0; j < a->rlen[row]; j++) a->getrowcols[j] = a->colidx[shift + a->sliceheight * j];
2156108893eSStefano Zampini     *idx = a->getrowcols;
2166108893eSStefano Zampini   }
2176108893eSStefano Zampini   if (v) {
2186108893eSStefano Zampini     PetscInt j;
21907e43b41SHong Zhang     for (j = 0; j < a->rlen[row]; j++) a->getrowvals[j] = a->val[shift + a->sliceheight * j];
2206108893eSStefano Zampini     *v = a->getrowvals;
2216108893eSStefano Zampini   }
2223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2236108893eSStefano Zampini }
2246108893eSStefano Zampini 
225d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRestoreRow_SeqSELL(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
226d71ae5a4SJacob Faibussowitsch {
2276108893eSStefano Zampini   PetscFunctionBegin;
2283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2296108893eSStefano Zampini }
2306108893eSStefano Zampini 
231d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_SeqSELL_SeqAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
232d71ae5a4SJacob Faibussowitsch {
233d4002b98SHong Zhang   Mat          B;
234d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
235e3f1f374SStefano Zampini   PetscInt     i;
236d4002b98SHong Zhang 
237d4002b98SHong Zhang   PetscFunctionBegin;
238ad013a7bSRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
239ad013a7bSRichard Tran Mills     B = *newmat;
2409566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(B));
241ad013a7bSRichard Tran Mills   } else {
2429566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2439566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
2449566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQAIJ));
2459566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, a->rlen));
246ad013a7bSRichard Tran Mills   }
247d4002b98SHong Zhang 
248e3f1f374SStefano Zampini   for (i = 0; i < A->rmap->n; i++) {
249e108cb99SStefano Zampini     PetscInt     nz = 0, *cols = NULL;
250e108cb99SStefano Zampini     PetscScalar *vals = NULL;
251e3f1f374SStefano Zampini 
2529566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqSELL(A, i, &nz, &cols, &vals));
2539566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, nz, cols, vals, INSERT_VALUES));
2549566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow_SeqSELL(A, i, &nz, &cols, &vals));
255d4002b98SHong Zhang   }
256e3f1f374SStefano Zampini 
2579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
259d4002b98SHong Zhang   B->rmap->bs = A->rmap->bs;
260d4002b98SHong Zhang 
261d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
2629566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
263d4002b98SHong Zhang   } else {
264d4002b98SHong Zhang     *newmat = B;
265d4002b98SHong Zhang   }
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
267d4002b98SHong Zhang }
268d4002b98SHong Zhang 
269d4002b98SHong Zhang #include <../src/mat/impls/aij/seq/aij.h>
270d4002b98SHong Zhang 
271d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
272d71ae5a4SJacob Faibussowitsch {
273d4002b98SHong Zhang   Mat                B;
274d4002b98SHong Zhang   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *)A->data;
275d4002b98SHong Zhang   PetscInt          *ai = a->i, m = A->rmap->N, n = A->cmap->N, i, *rowlengths, row, ncols;
276d4002b98SHong Zhang   const PetscInt    *cols;
277d4002b98SHong Zhang   const PetscScalar *vals;
278d4002b98SHong Zhang 
279d4002b98SHong Zhang   PetscFunctionBegin;
280ad013a7bSRichard Tran Mills 
281ad013a7bSRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
282ad013a7bSRichard Tran Mills     B = *newmat;
283ad013a7bSRichard Tran Mills   } else {
284d5e5b2e5SBarry Smith     if (PetscDefined(USE_DEBUG) || !a->ilen) {
2859566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m, &rowlengths));
286ad540459SPierre Jolivet       for (i = 0; i < m; i++) rowlengths[i] = ai[i + 1] - ai[i];
287d5e5b2e5SBarry Smith     }
288d5e5b2e5SBarry Smith     if (PetscDefined(USE_DEBUG) && a->ilen) {
289d5e5b2e5SBarry Smith       PetscBool eq;
2909566063dSJacob Faibussowitsch       PetscCall(PetscMemcmp(rowlengths, a->ilen, m * sizeof(PetscInt), &eq));
29128b400f6SJacob Faibussowitsch       PetscCheck(eq, PETSC_COMM_SELF, PETSC_ERR_PLIB, "SeqAIJ ilen array incorrect");
2929566063dSJacob Faibussowitsch       PetscCall(PetscFree(rowlengths));
293d5e5b2e5SBarry Smith       rowlengths = a->ilen;
294d5e5b2e5SBarry Smith     } else if (a->ilen) rowlengths = a->ilen;
2959566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2969566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
2979566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQSELL));
2989566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(B, 0, rowlengths));
2999566063dSJacob Faibussowitsch     if (rowlengths != a->ilen) PetscCall(PetscFree(rowlengths));
300ad013a7bSRichard Tran Mills   }
301d4002b98SHong Zhang 
302d4002b98SHong Zhang   for (row = 0; row < m; row++) {
3039566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqAIJ(A, row, &ncols, (PetscInt **)&cols, (PetscScalar **)&vals));
3049566063dSJacob Faibussowitsch     PetscCall(MatSetValues_SeqSELL(B, 1, &row, ncols, cols, vals, INSERT_VALUES));
3059566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow_SeqAIJ(A, row, &ncols, (PetscInt **)&cols, (PetscScalar **)&vals));
306d4002b98SHong Zhang   }
3079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
309d4002b98SHong Zhang   B->rmap->bs = A->rmap->bs;
310d4002b98SHong Zhang 
311d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
3129566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
313d4002b98SHong Zhang   } else {
314d4002b98SHong Zhang     *newmat = B;
315d4002b98SHong Zhang   }
3163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
317d4002b98SHong Zhang }
318d4002b98SHong Zhang 
319d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_SeqSELL(Mat A, Vec xx, Vec yy)
320d71ae5a4SJacob Faibussowitsch {
321d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
322d4002b98SHong Zhang   PetscScalar       *y;
323d4002b98SHong Zhang   const PetscScalar *x;
324d4002b98SHong Zhang   const MatScalar   *aval        = a->val;
325d4002b98SHong Zhang   PetscInt           totalslices = a->totalslices;
326d4002b98SHong Zhang   const PetscInt    *acolidx     = a->colidx;
3277285fed1SHong Zhang   PetscInt           i, j;
328d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
329d4002b98SHong Zhang   __m512d  vec_x, vec_y, vec_vals;
330d4002b98SHong Zhang   __m256i  vec_idx;
331d4002b98SHong Zhang   __mmask8 mask;
332d4002b98SHong Zhang   __m512d  vec_x2, vec_y2, vec_vals2, vec_x3, vec_y3, vec_vals3, vec_x4, vec_y4, vec_vals4;
333d4002b98SHong Zhang   __m256i  vec_idx2, vec_idx3, vec_idx4;
3345f70456aSHong 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)
335a48a6482SHong Zhang   __m128i   vec_idx;
336a48a6482SHong Zhang   __m256d   vec_x, vec_y, vec_y2, vec_vals;
337a48a6482SHong Zhang   MatScalar yval;
338a48a6482SHong Zhang   PetscInt  r, rows_left, row, nnz_in_row;
33921cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
340d4002b98SHong Zhang   __m128d   vec_x_tmp;
341d4002b98SHong Zhang   __m256d   vec_x, vec_y, vec_y2, vec_vals;
342d4002b98SHong Zhang   MatScalar yval;
343d4002b98SHong Zhang   PetscInt  r, rows_left, row, nnz_in_row;
344d4002b98SHong Zhang #else
34507e43b41SHong Zhang   PetscInt     k, sliceheight = a->sliceheight;
34607e43b41SHong Zhang   PetscScalar *sum;
347d4002b98SHong Zhang #endif
348d4002b98SHong Zhang 
349d4002b98SHong Zhang #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
350d4002b98SHong Zhang   #pragma disjoint(*x, *y, *aval)
351d4002b98SHong Zhang #endif
352d4002b98SHong Zhang 
353d4002b98SHong Zhang   PetscFunctionBegin;
3549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx, &x));
3559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(yy, &y));
356d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
35707e43b41SHong 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);
358d4002b98SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
359d4002b98SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
360d4002b98SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
361d4002b98SHong Zhang 
362d4002b98SHong Zhang     vec_y  = _mm512_setzero_pd();
363d4002b98SHong Zhang     vec_y2 = _mm512_setzero_pd();
364d4002b98SHong Zhang     vec_y3 = _mm512_setzero_pd();
365d4002b98SHong Zhang     vec_y4 = _mm512_setzero_pd();
366d4002b98SHong Zhang 
367da81f932SPierre Jolivet     j = a->sliidx[i] >> 3; /* 8 bytes are read at each time, corresponding to a slice column */
368d4002b98SHong Zhang     switch ((a->sliidx[i + 1] - a->sliidx[i]) / 8 & 3) {
369d4002b98SHong Zhang     case 3:
370d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3719371c9d4SSatish Balay       acolidx += 8;
3729371c9d4SSatish Balay       aval += 8;
373d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
3749371c9d4SSatish Balay       acolidx += 8;
3759371c9d4SSatish Balay       aval += 8;
376d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
3779371c9d4SSatish Balay       acolidx += 8;
3789371c9d4SSatish Balay       aval += 8;
379d4002b98SHong Zhang       j += 3;
380d4002b98SHong Zhang       break;
381d4002b98SHong Zhang     case 2:
382d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3839371c9d4SSatish Balay       acolidx += 8;
3849371c9d4SSatish Balay       aval += 8;
385d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
3869371c9d4SSatish Balay       acolidx += 8;
3879371c9d4SSatish Balay       aval += 8;
388d4002b98SHong Zhang       j += 2;
389d4002b98SHong Zhang       break;
390d4002b98SHong Zhang     case 1:
391d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3929371c9d4SSatish Balay       acolidx += 8;
3939371c9d4SSatish Balay       aval += 8;
394d4002b98SHong Zhang       j += 1;
395d4002b98SHong Zhang       break;
396d4002b98SHong Zhang     }
397d4002b98SHong Zhang   #pragma novector
398d4002b98SHong Zhang     for (; j < (a->sliidx[i + 1] >> 3); j += 4) {
399d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
4009371c9d4SSatish Balay       acolidx += 8;
4019371c9d4SSatish Balay       aval += 8;
402d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
4039371c9d4SSatish Balay       acolidx += 8;
4049371c9d4SSatish Balay       aval += 8;
405d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
4069371c9d4SSatish Balay       acolidx += 8;
4079371c9d4SSatish Balay       aval += 8;
408d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx4, vec_x4, vec_vals4, vec_y4);
4099371c9d4SSatish Balay       acolidx += 8;
4109371c9d4SSatish Balay       aval += 8;
411d4002b98SHong Zhang     }
412d4002b98SHong Zhang 
413d4002b98SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y2);
414d4002b98SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y3);
415d4002b98SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y4);
416d4002b98SHong Zhang     if (i == totalslices - 1 && A->rmap->n & 0x07) { /* if last slice has padding rows */
417d4002b98SHong Zhang       mask = (__mmask8)(0xff >> (8 - (A->rmap->n & 0x07)));
418ef588d5cSRichard Tran Mills       _mm512_mask_storeu_pd(&y[8 * i], mask, vec_y);
419d4002b98SHong Zhang     } else {
420ef588d5cSRichard Tran Mills       _mm512_storeu_pd(&y[8 * i], vec_y);
421d4002b98SHong Zhang     }
422d4002b98SHong Zhang   }
4235f70456aSHong 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)
42407e43b41SHong 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);
425a48a6482SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over full slices */
426a48a6482SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
427a48a6482SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
428a48a6482SHong Zhang 
429a48a6482SHong Zhang     /* last slice may have padding rows. Don't use vectorization. */
430a48a6482SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
431a48a6482SHong Zhang       rows_left = A->rmap->n - 8 * i;
432a48a6482SHong Zhang       for (r = 0; r < rows_left; ++r) {
433a48a6482SHong Zhang         yval       = (MatScalar)0;
434a48a6482SHong Zhang         row        = 8 * i + r;
435a48a6482SHong Zhang         nnz_in_row = a->rlen[row];
436a48a6482SHong Zhang         for (j = 0; j < nnz_in_row; ++j) yval += aval[8 * j + r] * x[acolidx[8 * j + r]];
437a48a6482SHong Zhang         y[row] = yval;
438a48a6482SHong Zhang       }
439a48a6482SHong Zhang       break;
440a48a6482SHong Zhang     }
441a48a6482SHong Zhang 
442a48a6482SHong Zhang     vec_y  = _mm256_setzero_pd();
443a48a6482SHong Zhang     vec_y2 = _mm256_setzero_pd();
444a48a6482SHong Zhang 
445a48a6482SHong Zhang   /* Process slice of height 8 (512 bits) via two subslices of height 4 (256 bits) via AVX */
446a48a6482SHong Zhang   #pragma novector
447a48a6482SHong Zhang   #pragma unroll(2)
448a48a6482SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
449a48a6482SHong Zhang       AVX2_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
4509371c9d4SSatish Balay       aval += 4;
4519371c9d4SSatish Balay       acolidx += 4;
452a48a6482SHong Zhang       AVX2_Mult_Private(vec_idx, vec_x, vec_vals, vec_y2);
4539371c9d4SSatish Balay       aval += 4;
4549371c9d4SSatish Balay       acolidx += 4;
455a48a6482SHong Zhang     }
456a48a6482SHong Zhang 
457ef588d5cSRichard Tran Mills     _mm256_storeu_pd(y + i * 8, vec_y);
458ef588d5cSRichard Tran Mills     _mm256_storeu_pd(y + i * 8 + 4, vec_y2);
459a48a6482SHong Zhang   }
46021cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
46107e43b41SHong 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);
462d4002b98SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over full slices */
463d4002b98SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
464d4002b98SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
465d4002b98SHong Zhang 
466d4002b98SHong Zhang     vec_y  = _mm256_setzero_pd();
467d4002b98SHong Zhang     vec_y2 = _mm256_setzero_pd();
468d4002b98SHong Zhang 
469d4002b98SHong Zhang     /* last slice may have padding rows. Don't use vectorization. */
470d4002b98SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
471d4002b98SHong Zhang       rows_left = A->rmap->n - 8 * i;
472d4002b98SHong Zhang       for (r = 0; r < rows_left; ++r) {
473d4002b98SHong Zhang         yval       = (MatScalar)0;
474d4002b98SHong Zhang         row        = 8 * i + r;
475d4002b98SHong Zhang         nnz_in_row = a->rlen[row];
476d4002b98SHong Zhang         for (j = 0; j < nnz_in_row; ++j) yval += aval[8 * j + r] * x[acolidx[8 * j + r]];
477d4002b98SHong Zhang         y[row] = yval;
478d4002b98SHong Zhang       }
479d4002b98SHong Zhang       break;
480d4002b98SHong Zhang     }
481d4002b98SHong Zhang 
482d4002b98SHong Zhang   /* Process slice of height 8 (512 bits) via two subslices of height 4 (256 bits) via AVX */
483a48a6482SHong Zhang   #pragma novector
484a48a6482SHong Zhang   #pragma unroll(2)
4857285fed1SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
486d4002b98SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
487165f9cc3SJed Brown       vec_x_tmp = _mm_setzero_pd();
488d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
489d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
490d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
491d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
492d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
493d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
494d4002b98SHong Zhang       vec_y     = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y);
495d4002b98SHong Zhang       aval += 4;
496d4002b98SHong Zhang 
497d4002b98SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
498d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
499d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
500d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
501d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
502d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
503d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
504d4002b98SHong Zhang       vec_y2    = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y2);
505d4002b98SHong Zhang       aval += 4;
506d4002b98SHong Zhang     }
507d4002b98SHong Zhang 
508d4002b98SHong Zhang     _mm256_storeu_pd(y + i * 8, vec_y);
509d4002b98SHong Zhang     _mm256_storeu_pd(y + i * 8 + 4, vec_y2);
510d4002b98SHong Zhang   }
511d4002b98SHong Zhang #else
51207e43b41SHong Zhang   PetscCall(PetscMalloc1(sliceheight, &sum));
513d4002b98SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
51407e43b41SHong Zhang     for (j = 0; j < sliceheight; j++) {
5152d1451d4SHong Zhang       sum[j] = 0.0;
51607e43b41SHong Zhang       for (k = a->sliidx[i] + j; k < a->sliidx[i + 1]; k += sliceheight) sum[j] += aval[k] * x[acolidx[k]];
517d4002b98SHong Zhang     }
51807e43b41SHong Zhang     if (i == totalslices - 1 && (A->rmap->n % sliceheight)) { /* if last slice has padding rows */
51907e43b41SHong Zhang       for (j = 0; j < (A->rmap->n % sliceheight); j++) y[sliceheight * i + j] = sum[j];
520d4002b98SHong Zhang     } else {
52107e43b41SHong Zhang       for (j = 0; j < sliceheight; j++) y[sliceheight * i + j] = sum[j];
522d4002b98SHong Zhang     }
523d4002b98SHong Zhang   }
52407e43b41SHong Zhang   PetscCall(PetscFree(sum));
525d4002b98SHong Zhang #endif
526d4002b98SHong Zhang 
5279566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->nz - a->nonzerorowcnt)); /* theoretical minimal FLOPs */
5289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx, &x));
5299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(yy, &y));
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
531d4002b98SHong Zhang }
532d4002b98SHong Zhang 
533d4002b98SHong Zhang #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h>
534d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAdd_SeqSELL(Mat A, Vec xx, Vec yy, Vec zz)
535d71ae5a4SJacob Faibussowitsch {
536d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
537d4002b98SHong Zhang   PetscScalar       *y, *z;
538d4002b98SHong Zhang   const PetscScalar *x;
539d4002b98SHong Zhang   const MatScalar   *aval        = a->val;
540d4002b98SHong Zhang   PetscInt           totalslices = a->totalslices;
541d4002b98SHong Zhang   const PetscInt    *acolidx     = a->colidx;
542d4002b98SHong Zhang   PetscInt           i, j;
543d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
5447285fed1SHong Zhang   __m512d  vec_x, vec_y, vec_vals;
545d4002b98SHong Zhang   __m256i  vec_idx;
546d4002b98SHong Zhang   __mmask8 mask;
5477285fed1SHong Zhang   __m512d  vec_x2, vec_y2, vec_vals2, vec_x3, vec_y3, vec_vals3, vec_x4, vec_y4, vec_vals4;
5487285fed1SHong Zhang   __m256i  vec_idx2, vec_idx3, vec_idx4;
54921cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
5507285fed1SHong Zhang   __m128d   vec_x_tmp;
5517285fed1SHong Zhang   __m256d   vec_x, vec_y, vec_y2, vec_vals;
5527285fed1SHong Zhang   MatScalar yval;
5537285fed1SHong Zhang   PetscInt  r, row, nnz_in_row;
554d4002b98SHong Zhang #else
55507e43b41SHong Zhang   PetscInt     k, sliceheight = a->sliceheight;
55607e43b41SHong Zhang   PetscScalar *sum;
557d4002b98SHong Zhang #endif
558d4002b98SHong Zhang 
559d4002b98SHong Zhang #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
560d4002b98SHong Zhang   #pragma disjoint(*x, *y, *aval)
561d4002b98SHong Zhang #endif
562d4002b98SHong Zhang 
563d4002b98SHong Zhang   PetscFunctionBegin;
5642d1451d4SHong Zhang   if (!a->nz) {
5652d1451d4SHong Zhang     PetscCall(VecCopy(yy, zz));
5662d1451d4SHong Zhang     PetscFunctionReturn(PETSC_SUCCESS);
5672d1451d4SHong Zhang   }
5689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx, &x));
5699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayPair(yy, zz, &y, &z));
570d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
57107e43b41SHong 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);
5727285fed1SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
5737285fed1SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
5747285fed1SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
5757285fed1SHong Zhang 
576d4002b98SHong Zhang     if (i == totalslices - 1 && A->rmap->n & 0x07) { /* if last slice has padding rows */
577d4002b98SHong Zhang       mask  = (__mmask8)(0xff >> (8 - (A->rmap->n & 0x07)));
578ef588d5cSRichard Tran Mills       vec_y = _mm512_mask_loadu_pd(vec_y, mask, &y[8 * i]);
5797285fed1SHong Zhang     } else {
580ef588d5cSRichard Tran Mills       vec_y = _mm512_loadu_pd(&y[8 * i]);
5817285fed1SHong Zhang     }
5827285fed1SHong Zhang     vec_y2 = _mm512_setzero_pd();
5837285fed1SHong Zhang     vec_y3 = _mm512_setzero_pd();
5847285fed1SHong Zhang     vec_y4 = _mm512_setzero_pd();
5857285fed1SHong Zhang 
586da81f932SPierre Jolivet     j = a->sliidx[i] >> 3; /* 8 bytes are read at each time, corresponding to a slice column */
5877285fed1SHong Zhang     switch ((a->sliidx[i + 1] - a->sliidx[i]) / 8 & 3) {
5887285fed1SHong Zhang     case 3:
5897285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
5909371c9d4SSatish Balay       acolidx += 8;
5919371c9d4SSatish Balay       aval += 8;
5927285fed1SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
5939371c9d4SSatish Balay       acolidx += 8;
5949371c9d4SSatish Balay       aval += 8;
5957285fed1SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
5969371c9d4SSatish Balay       acolidx += 8;
5979371c9d4SSatish Balay       aval += 8;
5987285fed1SHong Zhang       j += 3;
5997285fed1SHong Zhang       break;
6007285fed1SHong Zhang     case 2:
6017285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
6029371c9d4SSatish Balay       acolidx += 8;
6039371c9d4SSatish Balay       aval += 8;
6047285fed1SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
6059371c9d4SSatish Balay       acolidx += 8;
6069371c9d4SSatish Balay       aval += 8;
6077285fed1SHong Zhang       j += 2;
6087285fed1SHong Zhang       break;
6097285fed1SHong Zhang     case 1:
6107285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
6119371c9d4SSatish Balay       acolidx += 8;
6129371c9d4SSatish Balay       aval += 8;
6137285fed1SHong Zhang       j += 1;
6147285fed1SHong Zhang       break;
6157285fed1SHong Zhang     }
6167285fed1SHong Zhang   #pragma novector
6177285fed1SHong Zhang     for (; j < (a->sliidx[i + 1] >> 3); j += 4) {
6187285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
6199371c9d4SSatish Balay       acolidx += 8;
6209371c9d4SSatish Balay       aval += 8;
6217285fed1SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
6229371c9d4SSatish Balay       acolidx += 8;
6239371c9d4SSatish Balay       aval += 8;
6247285fed1SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
6259371c9d4SSatish Balay       acolidx += 8;
6269371c9d4SSatish Balay       aval += 8;
6277285fed1SHong Zhang       AVX512_Mult_Private(vec_idx4, vec_x4, vec_vals4, vec_y4);
6289371c9d4SSatish Balay       acolidx += 8;
6299371c9d4SSatish Balay       aval += 8;
6307285fed1SHong Zhang     }
6317285fed1SHong Zhang 
6327285fed1SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y2);
6337285fed1SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y3);
6347285fed1SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y4);
6357285fed1SHong Zhang     if (i == totalslices - 1 && A->rmap->n & 0x07) { /* if last slice has padding rows */
636ef588d5cSRichard Tran Mills       _mm512_mask_storeu_pd(&z[8 * i], mask, vec_y);
637d4002b98SHong Zhang     } else {
638ef588d5cSRichard Tran Mills       _mm512_storeu_pd(&z[8 * i], vec_y);
639d4002b98SHong Zhang     }
6407285fed1SHong Zhang   }
64121cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
64207e43b41SHong 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);
6437285fed1SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over full slices */
6447285fed1SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
6457285fed1SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
6467285fed1SHong Zhang 
6477285fed1SHong Zhang     /* last slice may have padding rows. Don't use vectorization. */
6487285fed1SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
6497285fed1SHong Zhang       for (r = 0; r < (A->rmap->n & 0x07); ++r) {
6507285fed1SHong Zhang         row        = 8 * i + r;
6517285fed1SHong Zhang         yval       = (MatScalar)0.0;
6527285fed1SHong Zhang         nnz_in_row = a->rlen[row];
6537285fed1SHong Zhang         for (j = 0; j < nnz_in_row; ++j) yval += aval[8 * j + r] * x[acolidx[8 * j + r]];
6547285fed1SHong Zhang         z[row] = y[row] + yval;
6557285fed1SHong Zhang       }
6567285fed1SHong Zhang       break;
6577285fed1SHong Zhang     }
6587285fed1SHong Zhang 
6597285fed1SHong Zhang     vec_y  = _mm256_loadu_pd(y + 8 * i);
6607285fed1SHong Zhang     vec_y2 = _mm256_loadu_pd(y + 8 * i + 4);
6617285fed1SHong Zhang 
6627285fed1SHong Zhang     /* Process slice of height 8 (512 bits) via two subslices of height 4 (256 bits) via AVX */
6637285fed1SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
6647285fed1SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
665165f9cc3SJed Brown       vec_x_tmp = _mm_setzero_pd();
6667285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6677285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
668165f9cc3SJed Brown       vec_x     = _mm256_setzero_pd();
6697285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
6707285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6717285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
6727285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
6737285fed1SHong Zhang       vec_y     = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y);
6747285fed1SHong Zhang       aval += 4;
6757285fed1SHong Zhang 
6767285fed1SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
6777285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6787285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
6797285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
6807285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6817285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
6827285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
6837285fed1SHong Zhang       vec_y2    = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y2);
6847285fed1SHong Zhang       aval += 4;
6857285fed1SHong Zhang     }
6867285fed1SHong Zhang 
6877285fed1SHong Zhang     _mm256_storeu_pd(z + i * 8, vec_y);
6887285fed1SHong Zhang     _mm256_storeu_pd(z + i * 8 + 4, vec_y2);
6897285fed1SHong Zhang   }
690d4002b98SHong Zhang #else
69107e43b41SHong Zhang   PetscCall(PetscMalloc1(sliceheight, &sum));
6927285fed1SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
69307e43b41SHong Zhang     for (j = 0; j < sliceheight; j++) {
6942d1451d4SHong Zhang       sum[j] = 0.0;
69507e43b41SHong Zhang       for (k = a->sliidx[i] + j; k < a->sliidx[i + 1]; k += sliceheight) sum[j] += aval[k] * x[acolidx[k]];
696d4002b98SHong Zhang     }
69707e43b41SHong Zhang     if (i == totalslices - 1 && (A->rmap->n % sliceheight)) {
69807e43b41SHong Zhang       for (j = 0; j < (A->rmap->n % sliceheight); j++) z[sliceheight * i + j] = y[sliceheight * i + j] + sum[j];
699d4002b98SHong Zhang     } else {
70007e43b41SHong Zhang       for (j = 0; j < sliceheight; j++) z[sliceheight * i + j] = y[sliceheight * i + j] + sum[j];
7017285fed1SHong Zhang     }
702d4002b98SHong Zhang   }
70307e43b41SHong Zhang   PetscCall(PetscFree(sum));
704d4002b98SHong Zhang #endif
705d4002b98SHong Zhang 
7069566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->nz));
7079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx, &x));
7089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayPair(yy, zz, &y, &z));
7093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
710d4002b98SHong Zhang }
711d4002b98SHong Zhang 
712d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAdd_SeqSELL(Mat A, Vec xx, Vec zz, Vec yy)
713d71ae5a4SJacob Faibussowitsch {
714d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
715d4002b98SHong Zhang   PetscScalar       *y;
716d4002b98SHong Zhang   const PetscScalar *x;
717d4002b98SHong Zhang   const MatScalar   *aval    = a->val;
718d4002b98SHong Zhang   const PetscInt    *acolidx = a->colidx;
71907e43b41SHong Zhang   PetscInt           i, j, r, row, nnz_in_row, totalslices = a->totalslices, sliceheight = a->sliceheight;
720d4002b98SHong Zhang 
721d4002b98SHong Zhang #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
722d4002b98SHong Zhang   #pragma disjoint(*x, *y, *aval)
723d4002b98SHong Zhang #endif
724d4002b98SHong Zhang 
725d4002b98SHong Zhang   PetscFunctionBegin;
726b94d7dedSBarry Smith   if (A->symmetric == PETSC_BOOL3_TRUE) {
7279566063dSJacob Faibussowitsch     PetscCall(MatMultAdd_SeqSELL(A, xx, zz, yy));
7283ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
7299fc32365SStefano Zampini   }
7309566063dSJacob Faibussowitsch   if (zz != yy) PetscCall(VecCopy(zz, yy));
7312d1451d4SHong Zhang 
7322d1451d4SHong Zhang   if (a->nz) {
7339566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(xx, &x));
7349566063dSJacob Faibussowitsch     PetscCall(VecGetArray(yy, &y));
735d4002b98SHong Zhang     for (i = 0; i < a->totalslices; i++) { /* loop over slices */
73607e43b41SHong Zhang       if (i == totalslices - 1 && (A->rmap->n % sliceheight)) {
73707e43b41SHong Zhang         for (r = 0; r < (A->rmap->n % sliceheight); ++r) {
73807e43b41SHong Zhang           row        = sliceheight * i + r;
7397285fed1SHong Zhang           nnz_in_row = a->rlen[row];
74007e43b41SHong Zhang           for (j = 0; j < nnz_in_row; ++j) y[acolidx[sliceheight * j + r]] += aval[sliceheight * j + r] * x[row];
7417285fed1SHong Zhang         }
7427285fed1SHong Zhang         break;
7437285fed1SHong Zhang       }
74407e43b41SHong Zhang       for (r = 0; r < sliceheight; ++r)
74507e43b41SHong Zhang         for (j = a->sliidx[i] + r; j < a->sliidx[i + 1]; j += sliceheight) y[acolidx[j]] += aval[j] * x[sliceheight * i + r];
746d4002b98SHong Zhang     }
7472d1451d4SHong Zhang     PetscCall(PetscLogFlops(2.0 * a->nz));
7489566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(xx, &x));
7499566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(yy, &y));
7502d1451d4SHong Zhang   }
7513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
752d4002b98SHong Zhang }
753d4002b98SHong Zhang 
754d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTranspose_SeqSELL(Mat A, Vec xx, Vec yy)
755d71ae5a4SJacob Faibussowitsch {
756d4002b98SHong Zhang   PetscFunctionBegin;
757b94d7dedSBarry Smith   if (A->symmetric == PETSC_BOOL3_TRUE) {
7589566063dSJacob Faibussowitsch     PetscCall(MatMult_SeqSELL(A, xx, yy));
7599fc32365SStefano Zampini   } else {
7609566063dSJacob Faibussowitsch     PetscCall(VecSet(yy, 0.0));
7619566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd_SeqSELL(A, xx, yy, yy));
7629fc32365SStefano Zampini   }
7633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
764d4002b98SHong Zhang }
765d4002b98SHong Zhang 
766d4002b98SHong Zhang /*
767d4002b98SHong Zhang      Checks for missing diagonals
768d4002b98SHong Zhang */
769d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMissingDiagonal_SeqSELL(Mat A, PetscBool *missing, PetscInt *d)
770d71ae5a4SJacob Faibussowitsch {
771d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
772d4002b98SHong Zhang   PetscInt    *diag, i;
773d4002b98SHong Zhang 
774d4002b98SHong Zhang   PetscFunctionBegin;
775d4002b98SHong Zhang   *missing = PETSC_FALSE;
776d4002b98SHong Zhang   if (A->rmap->n > 0 && !(a->colidx)) {
777d4002b98SHong Zhang     *missing = PETSC_TRUE;
778d4002b98SHong Zhang     if (d) *d = 0;
7799566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Matrix has no entries therefore is missing diagonal\n"));
780d4002b98SHong Zhang   } else {
781d4002b98SHong Zhang     diag = a->diag;
782d4002b98SHong Zhang     for (i = 0; i < A->rmap->n; i++) {
783d4002b98SHong Zhang       if (diag[i] == -1) {
784d4002b98SHong Zhang         *missing = PETSC_TRUE;
785d4002b98SHong Zhang         if (d) *d = i;
7869566063dSJacob Faibussowitsch         PetscCall(PetscInfo(A, "Matrix is missing diagonal number %" PetscInt_FMT "\n", i));
787d4002b98SHong Zhang         break;
788d4002b98SHong Zhang       }
789d4002b98SHong Zhang     }
790d4002b98SHong Zhang   }
7913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
792d4002b98SHong Zhang }
793d4002b98SHong Zhang 
794d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMarkDiagonal_SeqSELL(Mat A)
795d71ae5a4SJacob Faibussowitsch {
796d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
797d4002b98SHong Zhang   PetscInt     i, j, m = A->rmap->n, shift;
798d4002b98SHong Zhang 
799d4002b98SHong Zhang   PetscFunctionBegin;
800d4002b98SHong Zhang   if (!a->diag) {
8019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &a->diag));
802d4002b98SHong Zhang     a->free_diag = PETSC_TRUE;
803d4002b98SHong Zhang   }
804d4002b98SHong Zhang   for (i = 0; i < m; i++) {                                          /* loop over rows */
80507e43b41SHong Zhang     shift      = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
806d4002b98SHong Zhang     a->diag[i] = -1;
807d4002b98SHong Zhang     for (j = 0; j < a->rlen[i]; j++) {
80807e43b41SHong Zhang       if (a->colidx[shift + a->sliceheight * j] == i) {
80907e43b41SHong Zhang         a->diag[i] = shift + a->sliceheight * j;
810d4002b98SHong Zhang         break;
811d4002b98SHong Zhang       }
812d4002b98SHong Zhang     }
813d4002b98SHong Zhang   }
8143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
815d4002b98SHong Zhang }
816d4002b98SHong Zhang 
817d4002b98SHong Zhang /*
818d4002b98SHong Zhang   Negative shift indicates do not generate an error if there is a zero diagonal, just invert it anyways
819d4002b98SHong Zhang */
820d71ae5a4SJacob Faibussowitsch PetscErrorCode MatInvertDiagonal_SeqSELL(Mat A, PetscScalar omega, PetscScalar fshift)
821d71ae5a4SJacob Faibussowitsch {
822d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
823d4002b98SHong Zhang   PetscInt     i, *diag, m = A->rmap->n;
824d4002b98SHong Zhang   MatScalar   *val = a->val;
825d4002b98SHong Zhang   PetscScalar *idiag, *mdiag;
826d4002b98SHong Zhang 
827d4002b98SHong Zhang   PetscFunctionBegin;
8283ba16761SJacob Faibussowitsch   if (a->idiagvalid) PetscFunctionReturn(PETSC_SUCCESS);
8299566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_SeqSELL(A));
830d4002b98SHong Zhang   diag = a->diag;
831d4002b98SHong Zhang   if (!a->idiag) {
8329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(m, &a->idiag, m, &a->mdiag, m, &a->ssor_work));
833d4002b98SHong Zhang     val = a->val;
834d4002b98SHong Zhang   }
835d4002b98SHong Zhang   mdiag = a->mdiag;
836d4002b98SHong Zhang   idiag = a->idiag;
837d4002b98SHong Zhang 
838d4002b98SHong Zhang   if (omega == 1.0 && PetscRealPart(fshift) <= 0.0) {
839d4002b98SHong Zhang     for (i = 0; i < m; i++) {
840d4002b98SHong Zhang       mdiag[i] = val[diag[i]];
841d4002b98SHong Zhang       if (!PetscAbsScalar(mdiag[i])) { /* zero diagonal */
8420fdf79fbSJacob Faibussowitsch         PetscCheck(PetscRealPart(fshift), PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Zero diagonal on row %" PetscInt_FMT, i);
8439566063dSJacob Faibussowitsch         PetscCall(PetscInfo(A, "Zero diagonal on row %" PetscInt_FMT "\n", i));
844d4002b98SHong Zhang         A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
845d4002b98SHong Zhang         A->factorerror_zeropivot_value = 0.0;
846d4002b98SHong Zhang         A->factorerror_zeropivot_row   = i;
847d4002b98SHong Zhang       }
848d4002b98SHong Zhang       idiag[i] = 1.0 / val[diag[i]];
849d4002b98SHong Zhang     }
8509566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(m));
851d4002b98SHong Zhang   } else {
852d4002b98SHong Zhang     for (i = 0; i < m; i++) {
853d4002b98SHong Zhang       mdiag[i] = val[diag[i]];
854d4002b98SHong Zhang       idiag[i] = omega / (fshift + val[diag[i]]);
855d4002b98SHong Zhang     }
8569566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(2.0 * m));
857d4002b98SHong Zhang   }
858d4002b98SHong Zhang   a->idiagvalid = PETSC_TRUE;
8593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
860d4002b98SHong Zhang }
861d4002b98SHong Zhang 
862d71ae5a4SJacob Faibussowitsch PetscErrorCode MatZeroEntries_SeqSELL(Mat A)
863d71ae5a4SJacob Faibussowitsch {
864d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
865d4002b98SHong Zhang 
866d4002b98SHong Zhang   PetscFunctionBegin;
8679566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(a->val, a->sliidx[a->totalslices]));
8689566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
8693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
870d4002b98SHong Zhang }
871d4002b98SHong Zhang 
872d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_SeqSELL(Mat A)
873d71ae5a4SJacob Faibussowitsch {
874d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
875d4002b98SHong Zhang 
876d4002b98SHong Zhang   PetscFunctionBegin;
877d4002b98SHong Zhang #if defined(PETSC_USE_LOG)
8783ba16761SJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)A, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, A->rmap->n, A->cmap->n, a->nz));
879d4002b98SHong Zhang #endif
8809566063dSJacob Faibussowitsch   PetscCall(MatSeqXSELLFreeSELL(A, &a->val, &a->colidx));
8819566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->row));
8829566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->col));
8839566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->diag));
8849566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->rlen));
8859566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->sliidx));
8869566063dSJacob Faibussowitsch   PetscCall(PetscFree3(a->idiag, a->mdiag, a->ssor_work));
8879566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->solve_work));
8889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->icol));
8899566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->saved_values));
8909566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a->getrowcols, a->getrowvals));
8919566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
89290d2215bSHong Zhang #if defined(PETSC_HAVE_CUDA)
89390d2215bSHong Zhang   PetscCall(PetscFree(a->chunk_slice_map));
89490d2215bSHong Zhang #endif
895d4002b98SHong Zhang 
8969566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
8979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatStoreValues_C", NULL));
8989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatRetrieveValues_C", NULL));
8999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLSetPreallocation_C", NULL));
9002e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetArray_C", NULL));
9012e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLRestoreArray_C", NULL));
9022d1451d4SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_SeqSELL_SeqAIJ_C", NULL));
9032d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
9042d1451d4SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_SeqSELL_SeqSELLCUDA_C", NULL));
9052d1451d4SHong Zhang #endif
90607e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetFillRatio_C", NULL));
90707e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetMaxSliceWidth_C", NULL));
90807e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetAvgSliceWidth_C", NULL));
909*b921024eSHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetVarSliceSize_C", NULL));
91007e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLSetSliceHeight_C", NULL));
9113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
912d4002b98SHong Zhang }
913d4002b98SHong Zhang 
914d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetOption_SeqSELL(Mat A, MatOption op, PetscBool flg)
915d71ae5a4SJacob Faibussowitsch {
916d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
917d4002b98SHong Zhang 
918d4002b98SHong Zhang   PetscFunctionBegin;
919d4002b98SHong Zhang   switch (op) {
920d71ae5a4SJacob Faibussowitsch   case MAT_ROW_ORIENTED:
921d71ae5a4SJacob Faibussowitsch     a->roworiented = flg;
922d71ae5a4SJacob Faibussowitsch     break;
923d71ae5a4SJacob Faibussowitsch   case MAT_KEEP_NONZERO_PATTERN:
924d71ae5a4SJacob Faibussowitsch     a->keepnonzeropattern = flg;
925d71ae5a4SJacob Faibussowitsch     break;
926d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_LOCATIONS:
927d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? 0 : 1);
928d71ae5a4SJacob Faibussowitsch     break;
929d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_LOCATION_ERR:
930d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? -1 : 0);
931d71ae5a4SJacob Faibussowitsch     break;
932d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_ALLOCATION_ERR:
933d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? -2 : 0);
934d71ae5a4SJacob Faibussowitsch     break;
935d71ae5a4SJacob Faibussowitsch   case MAT_UNUSED_NONZERO_LOCATION_ERR:
936d71ae5a4SJacob Faibussowitsch     a->nounused = (flg ? -1 : 0);
937d71ae5a4SJacob Faibussowitsch     break;
9388c78258cSHong Zhang   case MAT_FORCE_DIAGONAL_ENTRIES:
939d4002b98SHong Zhang   case MAT_IGNORE_OFF_PROC_ENTRIES:
940d4002b98SHong Zhang   case MAT_USE_HASH_TABLE:
941d71ae5a4SJacob Faibussowitsch   case MAT_SORTED_FULL:
942d71ae5a4SJacob Faibussowitsch     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
943d71ae5a4SJacob Faibussowitsch     break;
944d4002b98SHong Zhang   case MAT_SPD:
945d4002b98SHong Zhang   case MAT_SYMMETRIC:
946d4002b98SHong Zhang   case MAT_STRUCTURALLY_SYMMETRIC:
947d4002b98SHong Zhang   case MAT_HERMITIAN:
948d4002b98SHong Zhang   case MAT_SYMMETRY_ETERNAL:
949b94d7dedSBarry Smith   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
950b94d7dedSBarry Smith   case MAT_SPD_ETERNAL:
951d4002b98SHong Zhang     /* These options are handled directly by MatSetOption() */
952d4002b98SHong Zhang     break;
953d71ae5a4SJacob Faibussowitsch   default:
954d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
955d4002b98SHong Zhang   }
9563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
957d4002b98SHong Zhang }
958d4002b98SHong Zhang 
959d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_SeqSELL(Mat A, Vec v)
960d71ae5a4SJacob Faibussowitsch {
961d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
962d4002b98SHong Zhang   PetscInt     i, j, n, shift;
963d4002b98SHong Zhang   PetscScalar *x, zero = 0.0;
964d4002b98SHong Zhang 
965d4002b98SHong Zhang   PetscFunctionBegin;
9669566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
96708401ef6SPierre Jolivet   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
968d4002b98SHong Zhang 
969d4002b98SHong Zhang   if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) {
970d4002b98SHong Zhang     PetscInt *diag = a->diag;
9719566063dSJacob Faibussowitsch     PetscCall(VecGetArray(v, &x));
972d4002b98SHong Zhang     for (i = 0; i < n; i++) x[i] = 1.0 / a->val[diag[i]];
9739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(v, &x));
9743ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
975d4002b98SHong Zhang   }
976d4002b98SHong Zhang 
9779566063dSJacob Faibussowitsch   PetscCall(VecSet(v, zero));
9789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &x));
979d4002b98SHong Zhang   for (i = 0; i < n; i++) {                                     /* loop over rows */
98007e43b41SHong Zhang     shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
981d4002b98SHong Zhang     x[i]  = 0;
982d4002b98SHong Zhang     for (j = 0; j < a->rlen[i]; j++) {
98307e43b41SHong Zhang       if (a->colidx[shift + a->sliceheight * j] == i) {
98407e43b41SHong Zhang         x[i] = a->val[shift + a->sliceheight * j];
985d4002b98SHong Zhang         break;
986d4002b98SHong Zhang       }
987d4002b98SHong Zhang     }
988d4002b98SHong Zhang   }
9899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &x));
9903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
991d4002b98SHong Zhang }
992d4002b98SHong Zhang 
993d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalScale_SeqSELL(Mat A, Vec ll, Vec rr)
994d71ae5a4SJacob Faibussowitsch {
995d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
996d4002b98SHong Zhang   const PetscScalar *l, *r;
997d4002b98SHong Zhang   PetscInt           i, j, m, n, row;
998d4002b98SHong Zhang 
999d4002b98SHong Zhang   PetscFunctionBegin;
1000d4002b98SHong Zhang   if (ll) {
1001d4002b98SHong Zhang     /* The local size is used so that VecMPI can be passed to this routine
1002d4002b98SHong Zhang        by MatDiagonalScale_MPISELL */
10039566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(ll, &m));
100408401ef6SPierre Jolivet     PetscCheck(m == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Left scaling vector wrong length");
10059566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(ll, &l));
1006d4002b98SHong Zhang     for (i = 0; i < a->totalslices; i++) {                            /* loop over slices */
100707e43b41SHong Zhang       if (i == a->totalslices - 1 && (A->rmap->n % a->sliceheight)) { /* if last slice has padding rows */
100807e43b41SHong Zhang         for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = (row + 1) % a->sliceheight) {
100907e43b41SHong Zhang           if (row < (A->rmap->n % a->sliceheight)) a->val[j] *= l[a->sliceheight * i + row];
1010dab86139SHong Zhang         }
1011dab86139SHong Zhang       } else {
101207e43b41SHong 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]; }
1013d4002b98SHong Zhang       }
1014dab86139SHong Zhang     }
10159566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(ll, &l));
10169566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(a->nz));
1017d4002b98SHong Zhang   }
1018d4002b98SHong Zhang   if (rr) {
10199566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rr, &n));
102008401ef6SPierre Jolivet     PetscCheck(n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Right scaling vector wrong length");
10219566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(rr, &r));
1022d4002b98SHong Zhang     for (i = 0; i < a->totalslices; i++) {                            /* loop over slices */
102307e43b41SHong Zhang       if (i == a->totalslices - 1 && (A->rmap->n % a->sliceheight)) { /* if last slice has padding rows */
102407e43b41SHong Zhang         for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = ((row + 1) % a->sliceheight)) {
102507e43b41SHong Zhang           if (row < (A->rmap->n % a->sliceheight)) a->val[j] *= r[a->colidx[j]];
1026dab86139SHong Zhang         }
1027dab86139SHong Zhang       } else {
1028ad540459SPierre Jolivet         for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j++) a->val[j] *= r[a->colidx[j]];
1029d4002b98SHong Zhang       }
1030dab86139SHong Zhang     }
10319566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(rr, &r));
10329566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(a->nz));
1033d4002b98SHong Zhang   }
10349566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
10352d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
10362d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
10372d1451d4SHong Zhang #endif
10383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1039d4002b98SHong Zhang }
1040d4002b98SHong Zhang 
1041d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetValues_SeqSELL(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], PetscScalar v[])
1042d71ae5a4SJacob Faibussowitsch {
1043d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1044d4002b98SHong Zhang   PetscInt    *cp, i, k, low, high, t, row, col, l;
1045d4002b98SHong Zhang   PetscInt     shift;
1046d4002b98SHong Zhang   MatScalar   *vp;
1047d4002b98SHong Zhang 
1048d4002b98SHong Zhang   PetscFunctionBegin;
104968aafef3SStefano Zampini   for (k = 0; k < m; k++) { /* loop over requested rows */
1050d4002b98SHong Zhang     row = im[k];
1051d4002b98SHong Zhang     if (row < 0) continue;
10526bdcaf15SBarry 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);
105307e43b41SHong Zhang     shift = a->sliidx[row / a->sliceheight] + (row % a->sliceheight); /* starting index of the row */
1054d4002b98SHong Zhang     cp    = a->colidx + shift;                                        /* pointer to the row */
1055d4002b98SHong Zhang     vp    = a->val + shift;                                           /* pointer to the row */
105668aafef3SStefano Zampini     for (l = 0; l < n; l++) {                                         /* loop over requested columns */
1057d4002b98SHong Zhang       col = in[l];
1058d4002b98SHong Zhang       if (col < 0) continue;
10596bdcaf15SBarry 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);
10609371c9d4SSatish Balay       high = a->rlen[row];
10619371c9d4SSatish Balay       low  = 0; /* assume unsorted */
1062d4002b98SHong Zhang       while (high - low > 5) {
1063d4002b98SHong Zhang         t = (low + high) / 2;
106407e43b41SHong Zhang         if (*(cp + a->sliceheight * t) > col) high = t;
1065d4002b98SHong Zhang         else low = t;
1066d4002b98SHong Zhang       }
1067d4002b98SHong Zhang       for (i = low; i < high; i++) {
106807e43b41SHong Zhang         if (*(cp + a->sliceheight * i) > col) break;
106907e43b41SHong Zhang         if (*(cp + a->sliceheight * i) == col) {
107007e43b41SHong Zhang           *v++ = *(vp + a->sliceheight * i);
1071d4002b98SHong Zhang           goto finished;
1072d4002b98SHong Zhang         }
1073d4002b98SHong Zhang       }
1074d4002b98SHong Zhang       *v++ = 0.0;
1075d4002b98SHong Zhang     finished:;
1076d4002b98SHong Zhang     }
1077d4002b98SHong Zhang   }
10783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1079d4002b98SHong Zhang }
1080d4002b98SHong Zhang 
1081d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL_ASCII(Mat A, PetscViewer viewer)
1082d71ae5a4SJacob Faibussowitsch {
1083d4002b98SHong Zhang   Mat_SeqSELL      *a = (Mat_SeqSELL *)A->data;
1084d4002b98SHong Zhang   PetscInt          i, j, m = A->rmap->n, shift;
1085d4002b98SHong Zhang   const char       *name;
1086d4002b98SHong Zhang   PetscViewerFormat format;
1087d4002b98SHong Zhang 
1088d4002b98SHong Zhang   PetscFunctionBegin;
10899566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
1090d4002b98SHong Zhang   if (format == PETSC_VIEWER_ASCII_MATLAB) {
1091d4002b98SHong Zhang     PetscInt nofinalvalue = 0;
1092d4002b98SHong Zhang     /*
1093d4002b98SHong Zhang     if (m && ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap->n-1))) {
1094d4002b98SHong Zhang       nofinalvalue = 1;
1095d4002b98SHong Zhang     }
1096d4002b98SHong Zhang     */
10979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
10989566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Size = %" PetscInt_FMT " %" PetscInt_FMT " \n", m, A->cmap->n));
10999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Nonzeros = %" PetscInt_FMT " \n", a->nz));
1100d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
11019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",4);\n", a->nz + nofinalvalue));
1102d4002b98SHong Zhang #else
11039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",3);\n", a->nz + nofinalvalue));
1104d4002b98SHong Zhang #endif
11059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = [\n"));
1106d4002b98SHong Zhang 
1107d4002b98SHong Zhang     for (i = 0; i < m; i++) {
110807e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1109d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1110d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
111107e43b41SHong 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])));
1112d4002b98SHong Zhang #else
111307e43b41SHong 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]));
1114d4002b98SHong Zhang #endif
1115d4002b98SHong Zhang       }
1116d4002b98SHong Zhang     }
1117d4002b98SHong Zhang     /*
1118d4002b98SHong Zhang     if (nofinalvalue) {
1119d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
11209566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e %18.16e\n",m,A->cmap->n,0.,0.));
1121d4002b98SHong Zhang #else
11229566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n",m,A->cmap->n,0.0));
1123d4002b98SHong Zhang #endif
1124d4002b98SHong Zhang     }
1125d4002b98SHong Zhang     */
11269566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)A, &name));
11279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "];\n %s = spconvert(zzz);\n", name));
11289566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1129d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
11303ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1131d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
11329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1133d4002b98SHong Zhang     for (i = 0; i < m; i++) {
11349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
113507e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1136d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1137d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
113807e43b41SHong Zhang         if (PetscImaginaryPart(a->val[shift + a->sliceheight * j]) > 0.0 && PetscRealPart(a->val[shift + a->sliceheight * j]) != 0.0) {
113907e43b41SHong 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])));
114007e43b41SHong Zhang         } else if (PetscImaginaryPart(a->val[shift + a->sliceheight * j]) < 0.0 && PetscRealPart(a->val[shift + a->sliceheight * j]) != 0.0) {
114107e43b41SHong 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])));
114207e43b41SHong Zhang         } else if (PetscRealPart(a->val[shift + a->sliceheight * j]) != 0.0) {
114307e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + a->sliceheight * j], (double)PetscRealPart(a->val[shift + a->sliceheight * j])));
1144d4002b98SHong Zhang         }
1145d4002b98SHong Zhang #else
114607e43b41SHong 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]));
1147d4002b98SHong Zhang #endif
1148d4002b98SHong Zhang       }
11499566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1150d4002b98SHong Zhang     }
11519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1152d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_DENSE) {
1153d4002b98SHong Zhang     PetscInt    cnt = 0, jcnt;
1154d4002b98SHong Zhang     PetscScalar value;
1155d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1156d4002b98SHong Zhang     PetscBool realonly = PETSC_TRUE;
1157d4002b98SHong Zhang     for (i = 0; i < a->sliidx[a->totalslices]; i++) {
1158d4002b98SHong Zhang       if (PetscImaginaryPart(a->val[i]) != 0.0) {
1159d4002b98SHong Zhang         realonly = PETSC_FALSE;
1160d4002b98SHong Zhang         break;
1161d4002b98SHong Zhang       }
1162d4002b98SHong Zhang     }
1163d4002b98SHong Zhang #endif
1164d4002b98SHong Zhang 
11659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1166d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1167d4002b98SHong Zhang       jcnt  = 0;
116807e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1169d4002b98SHong Zhang       for (j = 0; j < A->cmap->n; j++) {
117007e43b41SHong Zhang         if (jcnt < a->rlen[i] && j == a->colidx[shift + a->sliceheight * j]) {
1171d4002b98SHong Zhang           value = a->val[cnt++];
1172d4002b98SHong Zhang           jcnt++;
1173d4002b98SHong Zhang         } else {
1174d4002b98SHong Zhang           value = 0.0;
1175d4002b98SHong Zhang         }
1176d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1177d4002b98SHong Zhang         if (realonly) {
11789566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)PetscRealPart(value)));
1179d4002b98SHong Zhang         } else {
11809566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e+%7.5e i ", (double)PetscRealPart(value), (double)PetscImaginaryPart(value)));
1181d4002b98SHong Zhang         }
1182d4002b98SHong Zhang #else
11839566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)value));
1184d4002b98SHong Zhang #endif
1185d4002b98SHong Zhang       }
11869566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1187d4002b98SHong Zhang     }
11889566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1189d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_MATRIXMARKET) {
1190d4002b98SHong Zhang     PetscInt fshift = 1;
11919566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1192d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
11939566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate complex general\n"));
1194d4002b98SHong Zhang #else
11959566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate real general\n"));
1196d4002b98SHong Zhang #endif
11979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", m, A->cmap->n, a->nz));
1198d4002b98SHong Zhang     for (i = 0; i < m; i++) {
119907e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1200d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1201d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
120207e43b41SHong 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])));
1203d4002b98SHong Zhang #else
120407e43b41SHong 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]));
1205d4002b98SHong Zhang #endif
1206d4002b98SHong Zhang       }
1207d4002b98SHong Zhang     }
12089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
120968aafef3SStefano Zampini   } else if (format == PETSC_VIEWER_NATIVE) {
121068aafef3SStefano Zampini     for (i = 0; i < a->totalslices; i++) { /* loop over slices */
121168aafef3SStefano Zampini       PetscInt row;
12129566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "slice %" PetscInt_FMT ": %" PetscInt_FMT " %" PetscInt_FMT "\n", i, a->sliidx[i], a->sliidx[i + 1]));
121307e43b41SHong Zhang       for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = (row + 1) % a->sliceheight) {
121468aafef3SStefano Zampini #if defined(PETSC_USE_COMPLEX)
121568aafef3SStefano Zampini         if (PetscImaginaryPart(a->val[j]) > 0.0) {
121607e43b41SHong 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])));
121768aafef3SStefano Zampini         } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
121807e43b41SHong 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])));
121968aafef3SStefano Zampini         } else {
122007e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g\n", a->sliceheight * i + row, a->colidx[j], (double)PetscRealPart(a->val[j])));
122168aafef3SStefano Zampini         }
122268aafef3SStefano Zampini #else
122307e43b41SHong Zhang         PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g\n", a->sliceheight * i + row, a->colidx[j], (double)a->val[j]));
122468aafef3SStefano Zampini #endif
122568aafef3SStefano Zampini       }
122668aafef3SStefano Zampini     }
1227d4002b98SHong Zhang   } else {
12289566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1229d4002b98SHong Zhang     if (A->factortype) {
1230d4002b98SHong Zhang       for (i = 0; i < m; i++) {
123107e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
12329566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
1233d4002b98SHong Zhang         /* L part */
123407e43b41SHong Zhang         for (j = shift; j < a->diag[i]; j += a->sliceheight) {
1235d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
123607e43b41SHong Zhang           if (PetscImaginaryPart(a->val[shift + a->sliceheight * 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])));
123807e43b41SHong Zhang           } else if (PetscImaginaryPart(a->val[shift + a->sliceheight * 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         }
1247d4002b98SHong Zhang         /* diagonal */
1248d4002b98SHong Zhang         j = a->diag[i];
1249d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1250d4002b98SHong Zhang         if (PetscImaginaryPart(a->val[j]) > 0.0) {
12519566063dSJacob 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])));
1252d4002b98SHong Zhang         } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
12539566063dSJacob 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]))));
1254d4002b98SHong Zhang         } else {
12559566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)PetscRealPart(1.0 / a->val[j])));
1256d4002b98SHong Zhang         }
1257d4002b98SHong Zhang #else
12589566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)(1.0 / a->val[j])));
1259d4002b98SHong Zhang #endif
1260d4002b98SHong Zhang 
1261d4002b98SHong Zhang         /* U part */
126207e43b41SHong Zhang         for (j = a->diag[i] + 1; j < shift + a->sliceheight * a->rlen[i]; j += a->sliceheight) {
1263d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1264d4002b98SHong Zhang           if (PetscImaginaryPart(a->val[j]) > 0.0) {
12659566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)PetscImaginaryPart(a->val[j])));
1266d4002b98SHong Zhang           } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
12679566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)(-PetscImaginaryPart(a->val[j]))));
1268d4002b98SHong Zhang           } else {
12699566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)PetscRealPart(a->val[j])));
1270d4002b98SHong Zhang           }
1271d4002b98SHong Zhang #else
12729566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)a->val[j]));
1273d4002b98SHong Zhang #endif
1274d4002b98SHong Zhang         }
12759566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1276d4002b98SHong Zhang       }
1277d4002b98SHong Zhang     } else {
1278d4002b98SHong Zhang       for (i = 0; i < m; i++) {
127907e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
12809566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
1281d4002b98SHong Zhang         for (j = 0; j < a->rlen[i]; j++) {
1282d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1283d4002b98SHong Zhang           if (PetscImaginaryPart(a->val[j]) > 0.0) {
128407e43b41SHong 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])));
1285d4002b98SHong Zhang           } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
128607e43b41SHong 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])));
1287d4002b98SHong Zhang           } else {
128807e43b41SHong Zhang             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + a->sliceheight * j], (double)PetscRealPart(a->val[shift + a->sliceheight * j])));
1289d4002b98SHong Zhang           }
1290d4002b98SHong Zhang #else
129107e43b41SHong Zhang           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + a->sliceheight * j], (double)a->val[shift + a->sliceheight * j]));
1292d4002b98SHong Zhang #endif
1293d4002b98SHong Zhang         }
12949566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1295d4002b98SHong Zhang       }
1296d4002b98SHong Zhang     }
12979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1298d4002b98SHong Zhang   }
12999566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
13003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1301d4002b98SHong Zhang }
1302d4002b98SHong Zhang 
1303d4002b98SHong Zhang #include <petscdraw.h>
1304d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL_Draw_Zoom(PetscDraw draw, void *Aa)
1305d71ae5a4SJacob Faibussowitsch {
1306d4002b98SHong Zhang   Mat               A = (Mat)Aa;
1307d4002b98SHong Zhang   Mat_SeqSELL      *a = (Mat_SeqSELL *)A->data;
1308d4002b98SHong Zhang   PetscInt          i, j, m = A->rmap->n, shift;
1309d4002b98SHong Zhang   int               color;
1310d4002b98SHong Zhang   PetscReal         xl, yl, xr, yr, x_l, x_r, y_l, y_r;
1311d4002b98SHong Zhang   PetscViewer       viewer;
1312d4002b98SHong Zhang   PetscViewerFormat format;
1313d4002b98SHong Zhang 
1314d4002b98SHong Zhang   PetscFunctionBegin;
13159566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "Zoomviewer", (PetscObject *)&viewer));
13169566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
13179566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
1318d4002b98SHong Zhang 
1319d4002b98SHong Zhang   /* loop over matrix elements drawing boxes */
1320d4002b98SHong Zhang 
1321d4002b98SHong Zhang   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1322d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
1323d4002b98SHong Zhang     /* Blue for negative, Cyan for zero and  Red for positive */
1324d4002b98SHong Zhang     color = PETSC_DRAW_BLUE;
1325d4002b98SHong Zhang     for (i = 0; i < m; i++) {
132607e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
13279371c9d4SSatish Balay       y_l   = m - i - 1.0;
13289371c9d4SSatish Balay       y_r   = y_l + 1.0;
1329d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
133007e43b41SHong Zhang         x_l = a->colidx[shift + a->sliceheight * j];
13319371c9d4SSatish Balay         x_r = x_l + 1.0;
133207e43b41SHong Zhang         if (PetscRealPart(a->val[shift + a->sliceheight * j]) >= 0.) continue;
13339566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1334d4002b98SHong Zhang       }
1335d4002b98SHong Zhang     }
1336d4002b98SHong Zhang     color = PETSC_DRAW_CYAN;
1337d4002b98SHong Zhang     for (i = 0; i < m; i++) {
133807e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
13399371c9d4SSatish Balay       y_l   = m - i - 1.0;
13409371c9d4SSatish Balay       y_r   = y_l + 1.0;
1341d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
134207e43b41SHong Zhang         x_l = a->colidx[shift + a->sliceheight * j];
13439371c9d4SSatish Balay         x_r = x_l + 1.0;
134407e43b41SHong Zhang         if (a->val[shift + a->sliceheight * j] != 0.) continue;
13459566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1346d4002b98SHong Zhang       }
1347d4002b98SHong Zhang     }
1348d4002b98SHong Zhang     color = PETSC_DRAW_RED;
1349d4002b98SHong Zhang     for (i = 0; i < m; i++) {
135007e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
13519371c9d4SSatish Balay       y_l   = m - i - 1.0;
13529371c9d4SSatish Balay       y_r   = y_l + 1.0;
1353d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
135407e43b41SHong Zhang         x_l = a->colidx[shift + a->sliceheight * j];
13559371c9d4SSatish Balay         x_r = x_l + 1.0;
135607e43b41SHong Zhang         if (PetscRealPart(a->val[shift + a->sliceheight * j]) <= 0.) continue;
13579566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1358d4002b98SHong Zhang       }
1359d4002b98SHong Zhang     }
1360d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1361d4002b98SHong Zhang   } else {
1362d4002b98SHong Zhang     /* use contour shading to indicate magnitude of values */
1363d4002b98SHong Zhang     /* first determine max of all nonzero values */
1364d4002b98SHong Zhang     PetscReal minv = 0.0, maxv = 0.0;
1365d4002b98SHong Zhang     PetscInt  count = 0;
1366d4002b98SHong Zhang     PetscDraw popup;
1367d4002b98SHong Zhang     for (i = 0; i < a->sliidx[a->totalslices]; i++) {
1368d4002b98SHong Zhang       if (PetscAbsScalar(a->val[i]) > maxv) maxv = PetscAbsScalar(a->val[i]);
1369d4002b98SHong Zhang     }
1370d4002b98SHong Zhang     if (minv >= maxv) maxv = minv + PETSC_SMALL;
13719566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetPopup(draw, &popup));
13729566063dSJacob Faibussowitsch     PetscCall(PetscDrawScalePopup(popup, minv, maxv));
1373d4002b98SHong Zhang 
1374d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
1375d4002b98SHong Zhang     for (i = 0; i < m; i++) {
137607e43b41SHong Zhang       shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight;
1377d4002b98SHong Zhang       y_l   = m - i - 1.0;
1378d4002b98SHong Zhang       y_r   = y_l + 1.0;
1379d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
138007e43b41SHong Zhang         x_l   = a->colidx[shift + a->sliceheight * j];
1381d4002b98SHong Zhang         x_r   = x_l + 1.0;
1382d4002b98SHong Zhang         color = PetscDrawRealToColor(PetscAbsScalar(a->val[count]), minv, maxv);
13839566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1384d4002b98SHong Zhang         count++;
1385d4002b98SHong Zhang       }
1386d4002b98SHong Zhang     }
1387d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1388d4002b98SHong Zhang   }
13893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1390d4002b98SHong Zhang }
1391d4002b98SHong Zhang 
1392d4002b98SHong Zhang #include <petscdraw.h>
1393d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL_Draw(Mat A, PetscViewer viewer)
1394d71ae5a4SJacob Faibussowitsch {
1395d4002b98SHong Zhang   PetscDraw draw;
1396d4002b98SHong Zhang   PetscReal xr, yr, xl, yl, h, w;
1397d4002b98SHong Zhang   PetscBool isnull;
1398d4002b98SHong Zhang 
1399d4002b98SHong Zhang   PetscFunctionBegin;
14009566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
14019566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
14023ba16761SJacob Faibussowitsch   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
1403d4002b98SHong Zhang 
14049371c9d4SSatish Balay   xr = A->cmap->n;
14059371c9d4SSatish Balay   yr = A->rmap->n;
14069371c9d4SSatish Balay   h  = yr / 10.0;
14079371c9d4SSatish Balay   w  = xr / 10.0;
14089371c9d4SSatish Balay   xr += w;
14099371c9d4SSatish Balay   yr += h;
14109371c9d4SSatish Balay   xl = -w;
14119371c9d4SSatish Balay   yl = -h;
14129566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
14139566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", (PetscObject)viewer));
14149566063dSJacob Faibussowitsch   PetscCall(PetscDrawZoom(draw, MatView_SeqSELL_Draw_Zoom, A));
14159566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", NULL));
14169566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
14173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1418d4002b98SHong Zhang }
1419d4002b98SHong Zhang 
1420d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL(Mat A, PetscViewer viewer)
1421d71ae5a4SJacob Faibussowitsch {
1422d4002b98SHong Zhang   PetscBool iascii, isbinary, isdraw;
1423d4002b98SHong Zhang 
1424d4002b98SHong Zhang   PetscFunctionBegin;
14259566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
14269566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
14279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1428d4002b98SHong Zhang   if (iascii) {
14299566063dSJacob Faibussowitsch     PetscCall(MatView_SeqSELL_ASCII(A, viewer));
1430d4002b98SHong Zhang   } else if (isbinary) {
14319566063dSJacob Faibussowitsch     /* PetscCall(MatView_SeqSELL_Binary(A,viewer)); */
14321baa6e33SBarry Smith   } else if (isdraw) PetscCall(MatView_SeqSELL_Draw(A, viewer));
14333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1434d4002b98SHong Zhang }
1435d4002b98SHong Zhang 
1436d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_SeqSELL(Mat A, MatAssemblyType mode)
1437d71ae5a4SJacob Faibussowitsch {
1438d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1439d4002b98SHong Zhang   PetscInt     i, shift, row_in_slice, row, nrow, *cp, lastcol, j, k;
1440d4002b98SHong Zhang   MatScalar   *vp;
144190d2215bSHong Zhang #if defined(PETSC_HAVE_CUDA)
144290d2215bSHong Zhang   PetscInt totalchunks = 0;
144390d2215bSHong Zhang #endif
1444d4002b98SHong Zhang 
1445d4002b98SHong Zhang   PetscFunctionBegin;
14463ba16761SJacob Faibussowitsch   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
1447d4002b98SHong Zhang   /* To do: compress out the unused elements */
14489566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_SeqSELL(A));
14499566063dSJacob 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));
14509566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs));
14519566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", a->rlenmax));
14522d1451d4SHong Zhang   a->nonzerorowcnt = 0;
1453d4002b98SHong 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 */
1454d4002b98SHong Zhang   for (i = 0; i < a->totalslices; ++i) {
1455d4002b98SHong Zhang     shift = a->sliidx[i];                                                   /* starting index of the slice */
1456d4002b98SHong Zhang     cp    = a->colidx + shift;                                              /* pointer to the column indices of the slice */
1457d4002b98SHong Zhang     vp    = a->val + shift;                                                 /* pointer to the nonzero values of the slice */
145807e43b41SHong Zhang     for (row_in_slice = 0; row_in_slice < a->sliceheight; ++row_in_slice) { /* loop over rows in the slice */
145907e43b41SHong Zhang       row  = a->sliceheight * i + row_in_slice;
1460d4002b98SHong Zhang       nrow = a->rlen[row]; /* number of nonzeros in row */
1461d4002b98SHong Zhang       /*
1462d4002b98SHong Zhang         Search for the nearest nonzero. Normally setting the index to zero may cause extra communication.
1463d4002b98SHong Zhang         But if the entire slice are empty, it is fine to use 0 since the index will not be loaded.
1464d4002b98SHong Zhang       */
1465d4002b98SHong Zhang       lastcol = 0;
1466d4002b98SHong Zhang       if (nrow > 0) { /* nonempty row */
14672d1451d4SHong Zhang         a->nonzerorowcnt++;
146807e43b41SHong Zhang         lastcol = cp[a->sliceheight * (nrow - 1) + row_in_slice]; /* use the index from the last nonzero at current row */
1469aaa8cc7dSPierre Jolivet       } else if (!row_in_slice) {                                 /* first row of the correct slice is empty */
147007e43b41SHong Zhang         for (j = 1; j < a->sliceheight; j++) {
147107e43b41SHong Zhang           if (a->rlen[a->sliceheight * i + j]) {
1472d4002b98SHong Zhang             lastcol = cp[j];
1473d4002b98SHong Zhang             break;
1474d4002b98SHong Zhang           }
1475d4002b98SHong Zhang         }
1476d4002b98SHong Zhang       } else {
1477d4002b98SHong Zhang         if (a->sliidx[i + 1] != shift) lastcol = cp[row_in_slice - 1]; /* use the index from the previous row */
1478d4002b98SHong Zhang       }
1479d4002b98SHong Zhang 
148007e43b41SHong Zhang       for (k = nrow; k < (a->sliidx[i + 1] - shift) / a->sliceheight; ++k) {
148107e43b41SHong Zhang         cp[a->sliceheight * k + row_in_slice] = lastcol;
148207e43b41SHong Zhang         vp[a->sliceheight * k + row_in_slice] = (MatScalar)0;
1483d4002b98SHong Zhang       }
1484d4002b98SHong Zhang     }
1485d4002b98SHong Zhang   }
1486d4002b98SHong Zhang 
1487d4002b98SHong Zhang   A->info.mallocs += a->reallocs;
1488d4002b98SHong Zhang   a->reallocs = 0;
1489d4002b98SHong Zhang 
14909566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
149190d2215bSHong Zhang #if defined(PETSC_HAVE_CUDA)
149290d2215bSHong Zhang   if (!a->chunksize && a->totalslices) {
149390d2215bSHong Zhang     a->chunksize = 64;
149490d2215bSHong Zhang     while (a->chunksize < 1024 && 2 * a->chunksize <= a->sliidx[a->totalslices] / a->totalslices) a->chunksize *= 2;
149590d2215bSHong Zhang     totalchunks = 1 + (a->sliidx[a->totalslices] - 1) / a->chunksize;
149690d2215bSHong Zhang   }
149790d2215bSHong Zhang   if (totalchunks != a->totalchunks) {
149890d2215bSHong Zhang     PetscCall(PetscFree(a->chunk_slice_map));
149990d2215bSHong Zhang     PetscCall(PetscMalloc1(totalchunks, &a->chunk_slice_map));
150090d2215bSHong Zhang     a->totalchunks = totalchunks;
150190d2215bSHong Zhang   }
150290d2215bSHong Zhang   j = 0;
150390d2215bSHong Zhang   for (i = 0; i < totalchunks; i++) {
150490d2215bSHong Zhang     while (a->sliidx[j + 1] <= i * a->chunksize && j < a->totalslices) j++;
150590d2215bSHong Zhang     a->chunk_slice_map[i] = j;
150690d2215bSHong Zhang   }
150790d2215bSHong Zhang #endif
15083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1509d4002b98SHong Zhang }
1510d4002b98SHong Zhang 
1511d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetInfo_SeqSELL(Mat A, MatInfoType flag, MatInfo *info)
1512d71ae5a4SJacob Faibussowitsch {
1513d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1514d4002b98SHong Zhang 
1515d4002b98SHong Zhang   PetscFunctionBegin;
1516d4002b98SHong Zhang   info->block_size   = 1.0;
15173966268fSBarry Smith   info->nz_allocated = a->maxallocmat;
15183966268fSBarry Smith   info->nz_used      = a->sliidx[a->totalslices]; /* include padding zeros */
15193966268fSBarry Smith   info->nz_unneeded  = (a->maxallocmat - a->sliidx[a->totalslices]);
15203966268fSBarry Smith   info->assemblies   = A->num_ass;
15213966268fSBarry Smith   info->mallocs      = A->info.mallocs;
15224dfa11a4SJacob Faibussowitsch   info->memory       = 0; /* REVIEW ME */
1523d4002b98SHong Zhang   if (A->factortype) {
1524d4002b98SHong Zhang     info->fill_ratio_given  = A->info.fill_ratio_given;
1525d4002b98SHong Zhang     info->fill_ratio_needed = A->info.fill_ratio_needed;
1526d4002b98SHong Zhang     info->factor_mallocs    = A->info.factor_mallocs;
1527d4002b98SHong Zhang   } else {
1528d4002b98SHong Zhang     info->fill_ratio_given  = 0;
1529d4002b98SHong Zhang     info->fill_ratio_needed = 0;
1530d4002b98SHong Zhang     info->factor_mallocs    = 0;
1531d4002b98SHong Zhang   }
15323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1533d4002b98SHong Zhang }
1534d4002b98SHong Zhang 
1535d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetValues_SeqSELL(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
1536d71ae5a4SJacob Faibussowitsch {
1537d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1538d4002b98SHong Zhang   PetscInt     shift, i, k, l, low, high, t, ii, row, col, nrow;
1539d4002b98SHong Zhang   PetscInt    *cp, nonew = a->nonew, lastcol = -1;
1540d4002b98SHong Zhang   MatScalar   *vp, value;
15412d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
15422d1451d4SHong Zhang   PetscBool inserted = PETSC_FALSE;
15434e58db63SHong Zhang   PetscInt  mul      = DEVICE_MEM_ALIGN / a->sliceheight;
15442d1451d4SHong Zhang #endif
1545d4002b98SHong Zhang 
1546d4002b98SHong Zhang   PetscFunctionBegin;
1547d4002b98SHong Zhang   for (k = 0; k < m; k++) { /* loop over added rows */
1548d4002b98SHong Zhang     row = im[k];
1549d4002b98SHong Zhang     if (row < 0) continue;
15506bdcaf15SBarry 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);
155107e43b41SHong Zhang     shift = a->sliidx[row / a->sliceheight] + row % a->sliceheight; /* starting index of the row */
1552d4002b98SHong Zhang     cp    = a->colidx + shift;                                      /* pointer to the row */
1553d4002b98SHong Zhang     vp    = a->val + shift;                                         /* pointer to the row */
1554d4002b98SHong Zhang     nrow  = a->rlen[row];
1555d4002b98SHong Zhang     low   = 0;
1556d4002b98SHong Zhang     high  = nrow;
1557d4002b98SHong Zhang 
1558d4002b98SHong Zhang     for (l = 0; l < n; l++) { /* loop over added columns */
1559d4002b98SHong Zhang       col = in[l];
1560d4002b98SHong Zhang       if (col < 0) continue;
15616bdcaf15SBarry 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);
1562d4002b98SHong Zhang       if (a->roworiented) {
1563d4002b98SHong Zhang         value = v[l + k * n];
1564d4002b98SHong Zhang       } else {
1565d4002b98SHong Zhang         value = v[k + l * m];
1566d4002b98SHong Zhang       }
1567d4002b98SHong Zhang       if ((value == 0.0 && a->ignorezeroentries) && (is == ADD_VALUES)) continue;
1568d4002b98SHong Zhang 
1569ed73aabaSBarry Smith       /* search in this row for the specified column, i indicates the column to be set */
1570d4002b98SHong Zhang       if (col <= lastcol) low = 0;
1571d4002b98SHong Zhang       else high = nrow;
1572d4002b98SHong Zhang       lastcol = col;
1573d4002b98SHong Zhang       while (high - low > 5) {
1574d4002b98SHong Zhang         t = (low + high) / 2;
157507e43b41SHong Zhang         if (*(cp + a->sliceheight * t) > col) high = t;
1576d4002b98SHong Zhang         else low = t;
1577d4002b98SHong Zhang       }
1578d4002b98SHong Zhang       for (i = low; i < high; i++) {
157907e43b41SHong Zhang         if (*(cp + a->sliceheight * i) > col) break;
158007e43b41SHong Zhang         if (*(cp + a->sliceheight * i) == col) {
158107e43b41SHong Zhang           if (is == ADD_VALUES) *(vp + a->sliceheight * i) += value;
158207e43b41SHong Zhang           else *(vp + a->sliceheight * i) = value;
15832d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
15842d1451d4SHong Zhang           inserted = PETSC_TRUE;
15852d1451d4SHong Zhang #endif
1586d4002b98SHong Zhang           low = i + 1;
1587d4002b98SHong Zhang           goto noinsert;
1588d4002b98SHong Zhang         }
1589d4002b98SHong Zhang       }
1590d4002b98SHong Zhang       if (value == 0.0 && a->ignorezeroentries) goto noinsert;
1591d4002b98SHong Zhang       if (nonew == 1) goto noinsert;
159208401ef6SPierre Jolivet       PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
15934e58db63SHong Zhang #if defined(PETSC_HAVE_CUDA)
15944e58db63SHong 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, mul);
15954e58db63SHong Zhang #else
1596d4002b98SHong Zhang       /* If the current row length exceeds the slice width (e.g. nrow==slice_width), allocate a new space, otherwise do nothing */
15974e58db63SHong 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, 1);
15984e58db63SHong Zhang #endif
1599d4002b98SHong Zhang       /* add the new nonzero to the high position, shift the remaining elements in current row to the right by one slot */
1600d4002b98SHong Zhang       for (ii = nrow - 1; ii >= i; ii--) {
160107e43b41SHong Zhang         *(cp + a->sliceheight * (ii + 1)) = *(cp + a->sliceheight * ii);
160207e43b41SHong Zhang         *(vp + a->sliceheight * (ii + 1)) = *(vp + a->sliceheight * ii);
1603d4002b98SHong Zhang       }
1604d4002b98SHong Zhang       a->rlen[row]++;
160507e43b41SHong Zhang       *(cp + a->sliceheight * i) = col;
160607e43b41SHong Zhang       *(vp + a->sliceheight * i) = value;
1607d4002b98SHong Zhang       a->nz++;
1608d4002b98SHong Zhang       A->nonzerostate++;
16092d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
16102d1451d4SHong Zhang       inserted = PETSC_TRUE;
16112d1451d4SHong Zhang #endif
16129371c9d4SSatish Balay       low = i + 1;
16139371c9d4SSatish Balay       high++;
16149371c9d4SSatish Balay       nrow++;
1615d4002b98SHong Zhang     noinsert:;
1616d4002b98SHong Zhang     }
1617d4002b98SHong Zhang     a->rlen[row] = nrow;
1618d4002b98SHong Zhang   }
16192d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
16202d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED && inserted) A->offloadmask = PETSC_OFFLOAD_CPU;
16212d1451d4SHong Zhang #endif
16223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1623d4002b98SHong Zhang }
1624d4002b98SHong Zhang 
1625d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCopy_SeqSELL(Mat A, Mat B, MatStructure str)
1626d71ae5a4SJacob Faibussowitsch {
1627d4002b98SHong Zhang   PetscFunctionBegin;
1628d4002b98SHong Zhang   /* If the two matrices have the same copy implementation, use fast copy. */
1629d4002b98SHong Zhang   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1630d4002b98SHong Zhang     Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1631d4002b98SHong Zhang     Mat_SeqSELL *b = (Mat_SeqSELL *)B->data;
1632d4002b98SHong Zhang 
163308401ef6SPierre 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");
16349566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(b->val, a->val, a->sliidx[a->totalslices]));
1635d4002b98SHong Zhang   } else {
16369566063dSJacob Faibussowitsch     PetscCall(MatCopy_Basic(A, B, str));
1637d4002b98SHong Zhang   }
16383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1639d4002b98SHong Zhang }
1640d4002b98SHong Zhang 
1641d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUp_SeqSELL(Mat A)
1642d71ae5a4SJacob Faibussowitsch {
1643d4002b98SHong Zhang   PetscFunctionBegin;
16449566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(A, PETSC_DEFAULT, NULL));
16453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1646d4002b98SHong Zhang }
1647d4002b98SHong Zhang 
1648d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLGetArray_SeqSELL(Mat A, PetscScalar *array[])
1649d71ae5a4SJacob Faibussowitsch {
1650d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1651d4002b98SHong Zhang 
1652d4002b98SHong Zhang   PetscFunctionBegin;
1653d4002b98SHong Zhang   *array = a->val;
16543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1655d4002b98SHong Zhang }
1656d4002b98SHong Zhang 
1657d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLRestoreArray_SeqSELL(Mat A, PetscScalar *array[])
1658d71ae5a4SJacob Faibussowitsch {
1659d4002b98SHong Zhang   PetscFunctionBegin;
16603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1661d4002b98SHong Zhang }
1662d4002b98SHong Zhang 
1663d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRealPart_SeqSELL(Mat A)
1664d71ae5a4SJacob Faibussowitsch {
1665d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1666d4002b98SHong Zhang   PetscInt     i;
1667d4002b98SHong Zhang   MatScalar   *aval = a->val;
1668d4002b98SHong Zhang 
1669d4002b98SHong Zhang   PetscFunctionBegin;
1670d4002b98SHong Zhang   for (i = 0; i < a->sliidx[a->totalslices]; i++) aval[i] = PetscRealPart(aval[i]);
16712d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
16722d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
16732d1451d4SHong Zhang #endif
16743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1675d4002b98SHong Zhang }
1676d4002b98SHong Zhang 
1677d71ae5a4SJacob Faibussowitsch PetscErrorCode MatImaginaryPart_SeqSELL(Mat A)
1678d71ae5a4SJacob Faibussowitsch {
1679d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1680d4002b98SHong Zhang   PetscInt     i;
1681d4002b98SHong Zhang   MatScalar   *aval = a->val;
1682d4002b98SHong Zhang 
1683d4002b98SHong Zhang   PetscFunctionBegin;
1684d4002b98SHong Zhang   for (i = 0; i < a->sliidx[a->totalslices]; i++) aval[i] = PetscImaginaryPart(aval[i]);
16859566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
16862d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
16872d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
16882d1451d4SHong Zhang #endif
16893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1690d4002b98SHong Zhang }
1691d4002b98SHong Zhang 
1692d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_SeqSELL(Mat inA, PetscScalar alpha)
1693d71ae5a4SJacob Faibussowitsch {
1694d4002b98SHong Zhang   Mat_SeqSELL *a      = (Mat_SeqSELL *)inA->data;
1695d4002b98SHong Zhang   MatScalar   *aval   = a->val;
1696d4002b98SHong Zhang   PetscScalar  oalpha = alpha;
1697d4002b98SHong Zhang   PetscBLASInt one    = 1, size;
1698d4002b98SHong Zhang 
1699d4002b98SHong Zhang   PetscFunctionBegin;
17009566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(a->sliidx[a->totalslices], &size));
1701792fecdfSBarry Smith   PetscCallBLAS("BLASscal", BLASscal_(&size, &oalpha, aval, &one));
17029566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(a->nz));
17039566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(inA));
17042d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
17052d1451d4SHong Zhang   if (inA->offloadmask != PETSC_OFFLOAD_UNALLOCATED) inA->offloadmask = PETSC_OFFLOAD_CPU;
17062d1451d4SHong Zhang #endif
17073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1708d4002b98SHong Zhang }
1709d4002b98SHong Zhang 
1710d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift_SeqSELL(Mat Y, PetscScalar a)
1711d71ae5a4SJacob Faibussowitsch {
1712d4002b98SHong Zhang   Mat_SeqSELL *y = (Mat_SeqSELL *)Y->data;
1713d4002b98SHong Zhang 
1714d4002b98SHong Zhang   PetscFunctionBegin;
171548a46eb9SPierre Jolivet   if (!Y->preallocated || !y->nz) PetscCall(MatSeqSELLSetPreallocation(Y, 1, NULL));
17169566063dSJacob Faibussowitsch   PetscCall(MatShift_Basic(Y, a));
17173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1718d4002b98SHong Zhang }
1719d4002b98SHong Zhang 
1720d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSOR_SeqSELL(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
1721d71ae5a4SJacob Faibussowitsch {
1722d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
1723d4002b98SHong Zhang   PetscScalar       *x, sum, *t;
1724f4259b30SLisandro Dalcin   const MatScalar   *idiag = NULL, *mdiag;
1725d4002b98SHong Zhang   const PetscScalar *b, *xb;
1726d4002b98SHong Zhang   PetscInt           n, m = A->rmap->n, i, j, shift;
1727d4002b98SHong Zhang   const PetscInt    *diag;
1728d4002b98SHong Zhang 
1729d4002b98SHong Zhang   PetscFunctionBegin;
1730d4002b98SHong Zhang   its = its * lits;
1731d4002b98SHong Zhang 
1732d4002b98SHong Zhang   if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
17339566063dSJacob Faibussowitsch   if (!a->idiagvalid) PetscCall(MatInvertDiagonal_SeqSELL(A, omega, fshift));
1734d4002b98SHong Zhang   a->fshift = fshift;
1735d4002b98SHong Zhang   a->omega  = omega;
1736d4002b98SHong Zhang 
1737d4002b98SHong Zhang   diag  = a->diag;
1738d4002b98SHong Zhang   t     = a->ssor_work;
1739d4002b98SHong Zhang   idiag = a->idiag;
1740d4002b98SHong Zhang   mdiag = a->mdiag;
1741d4002b98SHong Zhang 
17429566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
17439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1744d4002b98SHong Zhang   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
174508401ef6SPierre Jolivet   PetscCheck(flag != SOR_APPLY_UPPER, PETSC_COMM_SELF, PETSC_ERR_SUP, "SOR_APPLY_UPPER is not implemented");
174608401ef6SPierre Jolivet   PetscCheck(flag != SOR_APPLY_LOWER, PETSC_COMM_SELF, PETSC_ERR_SUP, "SOR_APPLY_LOWER is not implemented");
1747aed4548fSBarry Smith   PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat");
1748d4002b98SHong Zhang 
1749d4002b98SHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
1750d4002b98SHong Zhang     if ((flag & SOR_FORWARD_SWEEP) || (flag & SOR_LOCAL_FORWARD_SWEEP)) {
1751d4002b98SHong Zhang       for (i = 0; i < m; i++) {
175207e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
1753d4002b98SHong Zhang         sum   = b[i];
175407e43b41SHong Zhang         n     = (diag[i] - shift) / a->sliceheight;
175507e43b41SHong Zhang         for (j = 0; j < n; j++) sum -= a->val[shift + a->sliceheight * j] * x[a->colidx[shift + a->sliceheight * j]];
1756d4002b98SHong Zhang         t[i] = sum;
1757d4002b98SHong Zhang         x[i] = sum * idiag[i];
1758d4002b98SHong Zhang       }
1759d4002b98SHong Zhang       xb = t;
17609566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(a->nz));
1761d4002b98SHong Zhang     } else xb = b;
1762d4002b98SHong Zhang     if ((flag & SOR_BACKWARD_SWEEP) || (flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1763d4002b98SHong Zhang       for (i = m - 1; i >= 0; i--) {
176407e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
1765d4002b98SHong Zhang         sum   = xb[i];
176607e43b41SHong Zhang         n     = a->rlen[i] - (diag[i] - shift) / a->sliceheight - 1;
176707e43b41SHong Zhang         for (j = 1; j <= n; j++) sum -= a->val[diag[i] + a->sliceheight * j] * x[a->colidx[diag[i] + a->sliceheight * j]];
1768d4002b98SHong Zhang         if (xb == b) {
1769d4002b98SHong Zhang           x[i] = sum * idiag[i];
1770d4002b98SHong Zhang         } else {
1771d4002b98SHong Zhang           x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
1772d4002b98SHong Zhang         }
1773d4002b98SHong Zhang       }
17749566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
1775d4002b98SHong Zhang     }
1776d4002b98SHong Zhang     its--;
1777d4002b98SHong Zhang   }
1778d4002b98SHong Zhang   while (its--) {
1779d4002b98SHong Zhang     if ((flag & SOR_FORWARD_SWEEP) || (flag & SOR_LOCAL_FORWARD_SWEEP)) {
1780d4002b98SHong Zhang       for (i = 0; i < m; i++) {
1781d4002b98SHong Zhang         /* lower */
178207e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
1783d4002b98SHong Zhang         sum   = b[i];
178407e43b41SHong Zhang         n     = (diag[i] - shift) / a->sliceheight;
178507e43b41SHong Zhang         for (j = 0; j < n; j++) sum -= a->val[shift + a->sliceheight * j] * x[a->colidx[shift + a->sliceheight * j]];
1786d4002b98SHong Zhang         t[i] = sum; /* save application of the lower-triangular part */
1787d4002b98SHong Zhang         /* upper */
178807e43b41SHong Zhang         n = a->rlen[i] - (diag[i] - shift) / a->sliceheight - 1;
178907e43b41SHong Zhang         for (j = 1; j <= n; j++) sum -= a->val[diag[i] + a->sliceheight * j] * x[a->colidx[diag[i] + a->sliceheight * j]];
1790d4002b98SHong Zhang         x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
1791d4002b98SHong Zhang       }
1792d4002b98SHong Zhang       xb = t;
17939566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(2.0 * a->nz));
1794d4002b98SHong Zhang     } else xb = b;
1795d4002b98SHong Zhang     if ((flag & SOR_BACKWARD_SWEEP) || (flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1796d4002b98SHong Zhang       for (i = m - 1; i >= 0; i--) {
179707e43b41SHong Zhang         shift = a->sliidx[i / a->sliceheight] + i % a->sliceheight; /* starting index of the row i */
1798d4002b98SHong Zhang         sum   = xb[i];
1799d4002b98SHong Zhang         if (xb == b) {
1800d4002b98SHong Zhang           /* whole matrix (no checkpointing available) */
1801d4002b98SHong Zhang           n = a->rlen[i];
180207e43b41SHong Zhang           for (j = 0; j < n; j++) sum -= a->val[shift + a->sliceheight * j] * x[a->colidx[shift + a->sliceheight * j]];
1803d4002b98SHong Zhang           x[i] = (1. - omega) * x[i] + (sum + mdiag[i] * x[i]) * idiag[i];
1804d4002b98SHong Zhang         } else { /* lower-triangular part has been saved, so only apply upper-triangular */
180507e43b41SHong Zhang           n = a->rlen[i] - (diag[i] - shift) / a->sliceheight - 1;
180607e43b41SHong Zhang           for (j = 1; j <= n; j++) sum -= a->val[diag[i] + a->sliceheight * j] * x[a->colidx[diag[i] + a->sliceheight * j]];
1807d4002b98SHong Zhang           x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
1808d4002b98SHong Zhang         }
1809d4002b98SHong Zhang       }
1810d4002b98SHong Zhang       if (xb == b) {
18119566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(2.0 * a->nz));
1812d4002b98SHong Zhang       } else {
18139566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
1814d4002b98SHong Zhang       }
1815d4002b98SHong Zhang     }
1816d4002b98SHong Zhang   }
18179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
18189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
18193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1820d4002b98SHong Zhang }
1821d4002b98SHong Zhang 
1822d4002b98SHong Zhang static struct _MatOps MatOps_Values = {MatSetValues_SeqSELL,
18236108893eSStefano Zampini                                        MatGetRow_SeqSELL,
18246108893eSStefano Zampini                                        MatRestoreRow_SeqSELL,
1825d4002b98SHong Zhang                                        MatMult_SeqSELL,
1826d4002b98SHong Zhang                                        /* 4*/ MatMultAdd_SeqSELL,
1827d4002b98SHong Zhang                                        MatMultTranspose_SeqSELL,
1828d4002b98SHong Zhang                                        MatMultTransposeAdd_SeqSELL,
1829f4259b30SLisandro Dalcin                                        NULL,
1830f4259b30SLisandro Dalcin                                        NULL,
1831f4259b30SLisandro Dalcin                                        NULL,
1832f4259b30SLisandro Dalcin                                        /* 10*/ NULL,
1833f4259b30SLisandro Dalcin                                        NULL,
1834f4259b30SLisandro Dalcin                                        NULL,
1835d4002b98SHong Zhang                                        MatSOR_SeqSELL,
1836f4259b30SLisandro Dalcin                                        NULL,
1837d4002b98SHong Zhang                                        /* 15*/ MatGetInfo_SeqSELL,
1838d4002b98SHong Zhang                                        MatEqual_SeqSELL,
1839d4002b98SHong Zhang                                        MatGetDiagonal_SeqSELL,
1840d4002b98SHong Zhang                                        MatDiagonalScale_SeqSELL,
1841f4259b30SLisandro Dalcin                                        NULL,
1842f4259b30SLisandro Dalcin                                        /* 20*/ NULL,
1843d4002b98SHong Zhang                                        MatAssemblyEnd_SeqSELL,
1844d4002b98SHong Zhang                                        MatSetOption_SeqSELL,
1845d4002b98SHong Zhang                                        MatZeroEntries_SeqSELL,
1846f4259b30SLisandro Dalcin                                        /* 24*/ NULL,
1847f4259b30SLisandro Dalcin                                        NULL,
1848f4259b30SLisandro Dalcin                                        NULL,
1849f4259b30SLisandro Dalcin                                        NULL,
1850f4259b30SLisandro Dalcin                                        NULL,
1851d4002b98SHong Zhang                                        /* 29*/ MatSetUp_SeqSELL,
1852f4259b30SLisandro Dalcin                                        NULL,
1853f4259b30SLisandro Dalcin                                        NULL,
1854f4259b30SLisandro Dalcin                                        NULL,
1855f4259b30SLisandro Dalcin                                        NULL,
1856d4002b98SHong Zhang                                        /* 34*/ MatDuplicate_SeqSELL,
1857f4259b30SLisandro Dalcin                                        NULL,
1858f4259b30SLisandro Dalcin                                        NULL,
1859f4259b30SLisandro Dalcin                                        NULL,
1860f4259b30SLisandro Dalcin                                        NULL,
1861f4259b30SLisandro Dalcin                                        /* 39*/ NULL,
1862f4259b30SLisandro Dalcin                                        NULL,
1863f4259b30SLisandro Dalcin                                        NULL,
1864d4002b98SHong Zhang                                        MatGetValues_SeqSELL,
1865d4002b98SHong Zhang                                        MatCopy_SeqSELL,
1866f4259b30SLisandro Dalcin                                        /* 44*/ NULL,
1867d4002b98SHong Zhang                                        MatScale_SeqSELL,
1868d4002b98SHong Zhang                                        MatShift_SeqSELL,
1869f4259b30SLisandro Dalcin                                        NULL,
1870f4259b30SLisandro Dalcin                                        NULL,
1871f4259b30SLisandro Dalcin                                        /* 49*/ NULL,
1872f4259b30SLisandro Dalcin                                        NULL,
1873f4259b30SLisandro Dalcin                                        NULL,
1874f4259b30SLisandro Dalcin                                        NULL,
1875f4259b30SLisandro Dalcin                                        NULL,
1876d4002b98SHong Zhang                                        /* 54*/ MatFDColoringCreate_SeqXAIJ,
1877f4259b30SLisandro Dalcin                                        NULL,
1878f4259b30SLisandro Dalcin                                        NULL,
1879f4259b30SLisandro Dalcin                                        NULL,
1880f4259b30SLisandro Dalcin                                        NULL,
1881f4259b30SLisandro Dalcin                                        /* 59*/ NULL,
1882d4002b98SHong Zhang                                        MatDestroy_SeqSELL,
1883d4002b98SHong Zhang                                        MatView_SeqSELL,
1884f4259b30SLisandro Dalcin                                        NULL,
1885f4259b30SLisandro Dalcin                                        NULL,
1886f4259b30SLisandro Dalcin                                        /* 64*/ NULL,
1887f4259b30SLisandro Dalcin                                        NULL,
1888f4259b30SLisandro Dalcin                                        NULL,
1889f4259b30SLisandro Dalcin                                        NULL,
1890f4259b30SLisandro Dalcin                                        NULL,
1891f4259b30SLisandro Dalcin                                        /* 69*/ NULL,
1892f4259b30SLisandro Dalcin                                        NULL,
1893f4259b30SLisandro Dalcin                                        NULL,
1894f4259b30SLisandro Dalcin                                        NULL,
1895f4259b30SLisandro Dalcin                                        NULL,
1896f4259b30SLisandro Dalcin                                        /* 74*/ NULL,
1897d4002b98SHong Zhang                                        MatFDColoringApply_AIJ, /* reuse the FDColoring function for AIJ */
1898f4259b30SLisandro Dalcin                                        NULL,
1899f4259b30SLisandro Dalcin                                        NULL,
1900f4259b30SLisandro Dalcin                                        NULL,
1901f4259b30SLisandro Dalcin                                        /* 79*/ NULL,
1902f4259b30SLisandro Dalcin                                        NULL,
1903f4259b30SLisandro Dalcin                                        NULL,
1904f4259b30SLisandro Dalcin                                        NULL,
1905f4259b30SLisandro Dalcin                                        NULL,
1906f4259b30SLisandro Dalcin                                        /* 84*/ NULL,
1907f4259b30SLisandro Dalcin                                        NULL,
1908f4259b30SLisandro Dalcin                                        NULL,
1909f4259b30SLisandro Dalcin                                        NULL,
1910f4259b30SLisandro Dalcin                                        NULL,
1911f4259b30SLisandro Dalcin                                        /* 89*/ NULL,
1912f4259b30SLisandro Dalcin                                        NULL,
1913f4259b30SLisandro Dalcin                                        NULL,
1914f4259b30SLisandro Dalcin                                        NULL,
1915f4259b30SLisandro Dalcin                                        NULL,
1916f4259b30SLisandro Dalcin                                        /* 94*/ NULL,
1917f4259b30SLisandro Dalcin                                        NULL,
1918f4259b30SLisandro Dalcin                                        NULL,
1919f4259b30SLisandro Dalcin                                        NULL,
1920f4259b30SLisandro Dalcin                                        NULL,
1921f4259b30SLisandro Dalcin                                        /* 99*/ NULL,
1922f4259b30SLisandro Dalcin                                        NULL,
1923f4259b30SLisandro Dalcin                                        NULL,
1924d4002b98SHong Zhang                                        MatConjugate_SeqSELL,
1925f4259b30SLisandro Dalcin                                        NULL,
1926f4259b30SLisandro Dalcin                                        /*104*/ NULL,
1927f4259b30SLisandro Dalcin                                        NULL,
1928f4259b30SLisandro Dalcin                                        NULL,
1929f4259b30SLisandro Dalcin                                        NULL,
1930f4259b30SLisandro Dalcin                                        NULL,
1931f4259b30SLisandro Dalcin                                        /*109*/ NULL,
1932f4259b30SLisandro Dalcin                                        NULL,
1933f4259b30SLisandro Dalcin                                        NULL,
1934f4259b30SLisandro Dalcin                                        NULL,
1935d4002b98SHong Zhang                                        MatMissingDiagonal_SeqSELL,
1936f4259b30SLisandro Dalcin                                        /*114*/ NULL,
1937f4259b30SLisandro Dalcin                                        NULL,
1938f4259b30SLisandro Dalcin                                        NULL,
1939f4259b30SLisandro Dalcin                                        NULL,
1940f4259b30SLisandro Dalcin                                        NULL,
1941f4259b30SLisandro Dalcin                                        /*119*/ NULL,
1942f4259b30SLisandro Dalcin                                        NULL,
1943f4259b30SLisandro Dalcin                                        NULL,
1944f4259b30SLisandro Dalcin                                        NULL,
1945f4259b30SLisandro Dalcin                                        NULL,
1946f4259b30SLisandro Dalcin                                        /*124*/ NULL,
1947f4259b30SLisandro Dalcin                                        NULL,
1948f4259b30SLisandro Dalcin                                        NULL,
1949f4259b30SLisandro Dalcin                                        NULL,
1950f4259b30SLisandro Dalcin                                        NULL,
1951f4259b30SLisandro Dalcin                                        /*129*/ NULL,
1952f4259b30SLisandro Dalcin                                        NULL,
1953f4259b30SLisandro Dalcin                                        NULL,
1954f4259b30SLisandro Dalcin                                        NULL,
1955f4259b30SLisandro Dalcin                                        NULL,
1956f4259b30SLisandro Dalcin                                        /*134*/ NULL,
1957f4259b30SLisandro Dalcin                                        NULL,
1958f4259b30SLisandro Dalcin                                        NULL,
1959f4259b30SLisandro Dalcin                                        NULL,
1960f4259b30SLisandro Dalcin                                        NULL,
1961f4259b30SLisandro Dalcin                                        /*139*/ NULL,
1962f4259b30SLisandro Dalcin                                        NULL,
1963f4259b30SLisandro Dalcin                                        NULL,
1964d4002b98SHong Zhang                                        MatFDColoringSetUp_SeqXAIJ,
1965f4259b30SLisandro Dalcin                                        NULL,
1966d70f29a3SPierre Jolivet                                        /*144*/ NULL,
1967d70f29a3SPierre Jolivet                                        NULL,
1968d70f29a3SPierre Jolivet                                        NULL,
196999a7f59eSMark Adams                                        NULL,
197099a7f59eSMark Adams                                        NULL,
19717fb60732SBarry Smith                                        NULL,
1972dec0b466SHong Zhang                                        /*150*/ NULL,
1973dec0b466SHong Zhang                                        NULL};
1974d4002b98SHong Zhang 
1975d71ae5a4SJacob Faibussowitsch PetscErrorCode MatStoreValues_SeqSELL(Mat mat)
1976d71ae5a4SJacob Faibussowitsch {
1977d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
1978d4002b98SHong Zhang 
1979d4002b98SHong Zhang   PetscFunctionBegin;
198028b400f6SJacob Faibussowitsch   PetscCheck(a->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1981d4002b98SHong Zhang 
1982d4002b98SHong Zhang   /* allocate space for values if not already there */
1983aa624791SPierre Jolivet   if (!a->saved_values) PetscCall(PetscMalloc1(a->sliidx[a->totalslices] + 1, &a->saved_values));
1984d4002b98SHong Zhang 
1985d4002b98SHong Zhang   /* copy values over */
19869566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->saved_values, a->val, a->sliidx[a->totalslices]));
19873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1988d4002b98SHong Zhang }
1989d4002b98SHong Zhang 
1990d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRetrieveValues_SeqSELL(Mat mat)
1991d71ae5a4SJacob Faibussowitsch {
1992d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
1993d4002b98SHong Zhang 
1994d4002b98SHong Zhang   PetscFunctionBegin;
199528b400f6SJacob Faibussowitsch   PetscCheck(a->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
199628b400f6SJacob Faibussowitsch   PetscCheck(a->saved_values, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatStoreValues(A);first");
19979566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->val, a->saved_values, a->sliidx[a->totalslices]));
19983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1999d4002b98SHong Zhang }
2000d4002b98SHong Zhang 
200107e43b41SHong Zhang PetscErrorCode MatSeqSELLGetFillRatio_SeqSELL(Mat mat, PetscReal *ratio)
200207e43b41SHong Zhang {
200307e43b41SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
200407e43b41SHong Zhang 
200507e43b41SHong Zhang   PetscFunctionBegin;
200607e43b41SHong Zhang   if (a->totalslices && a->sliidx[a->totalslices]) {
200707e43b41SHong Zhang     *ratio = (PetscReal)(a->sliidx[a->totalslices] - a->nz) / a->sliidx[a->totalslices];
200807e43b41SHong Zhang   } else {
200907e43b41SHong Zhang     *ratio = 0.0;
201007e43b41SHong Zhang   }
201107e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
201207e43b41SHong Zhang }
201307e43b41SHong Zhang 
201407e43b41SHong Zhang PetscErrorCode MatSeqSELLGetMaxSliceWidth_SeqSELL(Mat mat, PetscInt *slicewidth)
201507e43b41SHong Zhang {
201607e43b41SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
201707e43b41SHong Zhang   PetscInt     i, current_slicewidth;
201807e43b41SHong Zhang 
201907e43b41SHong Zhang   PetscFunctionBegin;
202007e43b41SHong Zhang   *slicewidth = 0;
202107e43b41SHong Zhang   for (i = 0; i < a->totalslices; i++) {
202207e43b41SHong Zhang     current_slicewidth = (a->sliidx[i + 1] - a->sliidx[i]) / a->sliceheight;
202307e43b41SHong Zhang     if (current_slicewidth > *slicewidth) *slicewidth = current_slicewidth;
202407e43b41SHong Zhang   }
202507e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
202607e43b41SHong Zhang }
202707e43b41SHong Zhang 
202807e43b41SHong Zhang PetscErrorCode MatSeqSELLGetAvgSliceWidth_SeqSELL(Mat mat, PetscReal *slicewidth)
202907e43b41SHong Zhang {
203007e43b41SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
203107e43b41SHong Zhang 
203207e43b41SHong Zhang   PetscFunctionBegin;
203307e43b41SHong Zhang   *slicewidth = 0;
203407e43b41SHong Zhang   if (a->totalslices) { *slicewidth = (PetscReal)a->sliidx[a->totalslices] / a->sliceheight / a->totalslices; }
203507e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
203607e43b41SHong Zhang }
203707e43b41SHong Zhang 
2038*b921024eSHong Zhang PetscErrorCode MatSeqSELLGetVarSliceSize_SeqSELL(Mat mat, PetscReal *variance)
2039*b921024eSHong Zhang {
2040*b921024eSHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
2041*b921024eSHong Zhang   PetscReal    mean;
2042*b921024eSHong Zhang   PetscInt     i, totalslices = a->totalslices, *sliidx = a->sliidx;
2043*b921024eSHong Zhang 
2044*b921024eSHong Zhang   PetscFunctionBegin;
2045*b921024eSHong Zhang   *variance = 0;
2046*b921024eSHong Zhang   if (totalslices) {
2047*b921024eSHong Zhang     mean = (PetscReal)sliidx[totalslices] / totalslices;
2048*b921024eSHong Zhang     for (i = 1; i <= totalslices; i++) { *variance += ((PetscReal)(sliidx[i] - sliidx[i - 1]) - mean) * ((PetscReal)(sliidx[i] - sliidx[i - 1]) - mean) / totalslices; }
2049*b921024eSHong Zhang   }
2050*b921024eSHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
2051*b921024eSHong Zhang }
2052*b921024eSHong Zhang 
205307e43b41SHong Zhang PetscErrorCode MatSeqSELLSetSliceHeight_SeqSELL(Mat A, PetscInt sliceheight)
205407e43b41SHong Zhang {
205507e43b41SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
205607e43b41SHong Zhang 
205707e43b41SHong Zhang   PetscFunctionBegin;
205807e43b41SHong Zhang   if (A->preallocated) PetscFunctionReturn(PETSC_SUCCESS);
205907e43b41SHong 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);
206007e43b41SHong Zhang   a->sliceheight = sliceheight;
20614e58db63SHong Zhang #if defined(PETSC_HAVE_CUDA)
20624e58db63SHong Zhang   PetscCheck(DEVICE_MEM_ALIGN % sliceheight == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "DEVICE_MEM_ALIGN is not divisible by the slice height %" PetscInt_FMT, sliceheight);
20634e58db63SHong Zhang #endif
206407e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
206507e43b41SHong Zhang }
206607e43b41SHong Zhang 
2067d4002b98SHong Zhang /*@C
206811a5261eSBarry Smith   MatSeqSELLRestoreArray - returns access to the array where the data for a `MATSEQSELL` matrix is stored obtained by `MatSeqSELLGetArray()`
2069d4002b98SHong Zhang 
2070d4002b98SHong Zhang   Not Collective
2071d4002b98SHong Zhang 
2072d4002b98SHong Zhang   Input Parameters:
207307e43b41SHong Zhang +  A - a `MATSEQSELL` matrix
207420f4b53cSBarry Smith -  array - pointer to the data
2075d4002b98SHong Zhang 
2076d4002b98SHong Zhang   Level: intermediate
2077d4002b98SHong Zhang 
207867be906fSBarry Smith .seealso: `Mat`, `MATSEQSELL`, `MatSeqSELLGetArray()`, `MatSeqSELLRestoreArrayF90()`
2079d4002b98SHong Zhang @*/
2080d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLRestoreArray(Mat A, PetscScalar **array)
2081d71ae5a4SJacob Faibussowitsch {
2082d4002b98SHong Zhang   PetscFunctionBegin;
2083cac4c232SBarry Smith   PetscUseMethod(A, "MatSeqSELLRestoreArray_C", (Mat, PetscScalar **), (A, array));
20843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2085d4002b98SHong Zhang }
2086d4002b98SHong Zhang 
208707e43b41SHong Zhang /*@C
208807e43b41SHong Zhang   MatSeqSELLGetFillRatio - returns a ratio that indicates the irregularity of the matrix.
208907e43b41SHong Zhang 
209007e43b41SHong Zhang   Not Collective
209107e43b41SHong Zhang 
209207e43b41SHong Zhang   Input Parameter:
209307e43b41SHong Zhang .  A - a MATSEQSELL matrix
209407e43b41SHong Zhang 
209507e43b41SHong Zhang   Output Parameter:
209607e43b41SHong Zhang .  ratio - ratio of number of padded zeros to number of allocated elements
209707e43b41SHong Zhang 
209807e43b41SHong Zhang   Level: intermediate
209907e43b41SHong Zhang @*/
210007e43b41SHong Zhang PetscErrorCode MatSeqSELLGetFillRatio(Mat A, PetscReal *ratio)
210107e43b41SHong Zhang {
210207e43b41SHong Zhang   PetscFunctionBegin;
210307e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLGetFillRatio_C", (Mat, PetscReal *), (A, ratio));
210407e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
210507e43b41SHong Zhang }
210607e43b41SHong Zhang 
210707e43b41SHong Zhang /*@C
210807e43b41SHong Zhang   MatSeqSELLGetMaxSliceWidth - returns the maximum slice width.
210907e43b41SHong Zhang 
211007e43b41SHong Zhang   Not Collective
211107e43b41SHong Zhang 
211207e43b41SHong Zhang   Input Parameter:
211307e43b41SHong Zhang .  A - a MATSEQSELL matrix
211407e43b41SHong Zhang 
211507e43b41SHong Zhang   Output Parameter:
211607e43b41SHong Zhang .  slicewidth - maximum slice width
211707e43b41SHong Zhang 
211807e43b41SHong Zhang  Level: intermediate
211907e43b41SHong Zhang @*/
212007e43b41SHong Zhang PetscErrorCode MatSeqSELLGetMaxSliceWidth(Mat A, PetscInt *slicewidth)
212107e43b41SHong Zhang {
212207e43b41SHong Zhang   PetscFunctionBegin;
212307e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLGetMaxSliceWidth_C", (Mat, PetscInt *), (A, slicewidth));
212407e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
212507e43b41SHong Zhang }
212607e43b41SHong Zhang 
212707e43b41SHong Zhang /*@C
212807e43b41SHong Zhang   MatSeqSELLGetAvgSliceWidth - returns the average slice width.
212907e43b41SHong Zhang 
213007e43b41SHong Zhang   Not Collective
213107e43b41SHong Zhang 
213207e43b41SHong Zhang   Input Parameter:
213307e43b41SHong Zhang .  A - a MATSEQSELL matrix
213407e43b41SHong Zhang 
213507e43b41SHong Zhang   Output Parameter:
213607e43b41SHong Zhang .  slicewidth - average slice width
213707e43b41SHong Zhang 
213807e43b41SHong Zhang   Level: intermediate
213907e43b41SHong Zhang @*/
214007e43b41SHong Zhang PetscErrorCode MatSeqSELLGetAvgSliceWidth(Mat A, PetscReal *slicewidth)
214107e43b41SHong Zhang {
214207e43b41SHong Zhang   PetscFunctionBegin;
214307e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLGetAvgSliceWidth_C", (Mat, PetscReal *), (A, slicewidth));
214407e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
214507e43b41SHong Zhang }
214607e43b41SHong Zhang 
214707e43b41SHong Zhang /*@C
214807e43b41SHong Zhang   MatSeqSELLSetSliceHeight - sets the slice height.
214907e43b41SHong Zhang 
215007e43b41SHong Zhang   Not Collective
215107e43b41SHong Zhang 
215207e43b41SHong Zhang   Input Parameters:
215307e43b41SHong Zhang +  A - a MATSEQSELL matrix
215407e43b41SHong Zhang -  sliceheight - slice height
215507e43b41SHong Zhang 
215607e43b41SHong Zhang   Notes:
215707e43b41SHong Zhang   You cannot change the slice height once it have been set.
215807e43b41SHong Zhang 
215907e43b41SHong Zhang   The slice height must be set before MatSetUp() or MatXXXSetPreallocation() is called.
216007e43b41SHong Zhang 
216107e43b41SHong Zhang   Level: intermediate
216207e43b41SHong Zhang @*/
216307e43b41SHong Zhang PetscErrorCode MatSeqSELLSetSliceHeight(Mat A, PetscInt sliceheight)
216407e43b41SHong Zhang {
216507e43b41SHong Zhang   PetscFunctionBegin;
216607e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLSetSliceHeight_C", (Mat, PetscInt), (A, sliceheight));
216707e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
216807e43b41SHong Zhang }
216907e43b41SHong Zhang 
217007e43b41SHong Zhang /*@C
217107e43b41SHong Zhang   MatSeqSELLGetVarSliceSize - returns the variance of the slice size.
217207e43b41SHong Zhang 
217307e43b41SHong Zhang   Not Collective
217407e43b41SHong Zhang 
217507e43b41SHong Zhang   Input Parameter:
217607e43b41SHong Zhang .  A - a MATSEQSELL matrix
217707e43b41SHong Zhang 
217807e43b41SHong Zhang   Output Parameter:
217907e43b41SHong Zhang .  variance - variance of the slice size
218007e43b41SHong Zhang 
218107e43b41SHong Zhang   Level: intermediate
218207e43b41SHong Zhang @*/
218307e43b41SHong Zhang PetscErrorCode MatSeqSELLGetVarSliceSize(Mat A, PetscReal *variance)
218407e43b41SHong Zhang {
218507e43b41SHong Zhang   PetscFunctionBegin;
218607e43b41SHong Zhang   PetscUseMethod(A, "MatSeqSELLGetVarSliceSize_C", (Mat, PetscReal *), (A, variance));
218707e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
218807e43b41SHong Zhang }
218907e43b41SHong Zhang 
21902d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
21912d1451d4SHong Zhang PETSC_EXTERN PetscErrorCode MatConvert_SeqSELL_SeqSELLCUDA(Mat);
21922d1451d4SHong Zhang #endif
21932d1451d4SHong Zhang 
2194d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_SeqSELL(Mat B)
2195d71ae5a4SJacob Faibussowitsch {
2196d4002b98SHong Zhang   Mat_SeqSELL *b;
2197d4002b98SHong Zhang   PetscMPIInt  size;
2198d4002b98SHong Zhang 
2199d4002b98SHong Zhang   PetscFunctionBegin;
22009566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
22019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
220208401ef6SPierre Jolivet   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Comm must be of size 1");
2203d4002b98SHong Zhang 
22044dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
2205d4002b98SHong Zhang 
2206d4002b98SHong Zhang   B->data = (void *)b;
2207d4002b98SHong Zhang 
22089566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(B->ops, &MatOps_Values, sizeof(struct _MatOps)));
2209d4002b98SHong Zhang 
2210f4259b30SLisandro Dalcin   b->row                = NULL;
2211f4259b30SLisandro Dalcin   b->col                = NULL;
2212f4259b30SLisandro Dalcin   b->icol               = NULL;
2213d4002b98SHong Zhang   b->reallocs           = 0;
2214d4002b98SHong Zhang   b->ignorezeroentries  = PETSC_FALSE;
2215d4002b98SHong Zhang   b->roworiented        = PETSC_TRUE;
2216d4002b98SHong Zhang   b->nonew              = 0;
2217f4259b30SLisandro Dalcin   b->diag               = NULL;
2218f4259b30SLisandro Dalcin   b->solve_work         = NULL;
2219f4259b30SLisandro Dalcin   B->spptr              = NULL;
2220f4259b30SLisandro Dalcin   b->saved_values       = NULL;
2221f4259b30SLisandro Dalcin   b->idiag              = NULL;
2222f4259b30SLisandro Dalcin   b->mdiag              = NULL;
2223f4259b30SLisandro Dalcin   b->ssor_work          = NULL;
2224d4002b98SHong Zhang   b->omega              = 1.0;
2225d4002b98SHong Zhang   b->fshift             = 0.0;
2226d4002b98SHong Zhang   b->idiagvalid         = PETSC_FALSE;
2227d4002b98SHong Zhang   b->keepnonzeropattern = PETSC_FALSE;
222807e43b41SHong Zhang   b->sliceheight        = 0;
2229d4002b98SHong Zhang 
22309566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQSELL));
22319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetArray_C", MatSeqSELLGetArray_SeqSELL));
22329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLRestoreArray_C", MatSeqSELLRestoreArray_SeqSELL));
22339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_SeqSELL));
22349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_SeqSELL));
22359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLSetPreallocation_C", MatSeqSELLSetPreallocation_SeqSELL));
22362d1451d4SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_SeqSELL_SeqAIJ_C", MatConvert_SeqSELL_SeqAIJ));
22372d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
22382d1451d4SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_SeqSELL_SeqSELLCUDA_C", MatConvert_SeqSELL_SeqSELLCUDA));
22392d1451d4SHong Zhang #endif
224007e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetFillRatio_C", MatSeqSELLGetFillRatio_SeqSELL));
224107e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetMaxSliceWidth_C", MatSeqSELLGetMaxSliceWidth_SeqSELL));
224207e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetAvgSliceWidth_C", MatSeqSELLGetAvgSliceWidth_SeqSELL));
2243*b921024eSHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetVarSliceSize_C", MatSeqSELLGetVarSliceSize_SeqSELL));
224407e43b41SHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLSetSliceHeight_C", MatSeqSELLSetSliceHeight_SeqSELL));
224507e43b41SHong Zhang 
224607e43b41SHong Zhang   PetscObjectOptionsBegin((PetscObject)B);
224707e43b41SHong Zhang   {
224807e43b41SHong Zhang     PetscInt  newsh = -1;
224907e43b41SHong Zhang     PetscBool flg;
225090d2215bSHong Zhang #if defined(PETSC_HAVE_CUDA)
225190d2215bSHong Zhang     PetscInt chunksize = 0;
225290d2215bSHong Zhang #endif
225307e43b41SHong Zhang 
225407e43b41SHong Zhang     PetscCall(PetscOptionsInt("-mat_sell_slice_height", "Set the slice height used to store SELL matrix", "MatSELLSetSliceHeight", newsh, &newsh, &flg));
225507e43b41SHong Zhang     if (flg) { PetscCall(MatSeqSELLSetSliceHeight(B, newsh)); }
225690d2215bSHong Zhang #if defined(PETSC_HAVE_CUDA)
225790d2215bSHong Zhang     PetscCall(PetscOptionsInt("-mat_sell_chunk_size", "Set the chunksize for load-balanced CUDA kernels. Choices include 64,128,256,512,1024", NULL, chunksize, &chunksize, &flg));
225890d2215bSHong Zhang     if (flg) {
225990d2215bSHong Zhang       PetscCheck(chunksize >= 64 && chunksize <= 1024 && chunksize % 64 == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "chunksize must be a number in {64,128,256,512,1024}: value %" PetscInt_FMT, chunksize);
226090d2215bSHong Zhang       b->chunksize = chunksize;
226190d2215bSHong Zhang     }
226290d2215bSHong Zhang #endif
226307e43b41SHong Zhang   }
226407e43b41SHong Zhang   PetscOptionsEnd();
22653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2266d4002b98SHong Zhang }
2267d4002b98SHong Zhang 
2268d4002b98SHong Zhang /*
2269d4002b98SHong Zhang  Given a matrix generated with MatGetFactor() duplicates all the information in A into B
2270d4002b98SHong Zhang  */
2271d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicateNoCreate_SeqSELL(Mat C, Mat A, MatDuplicateOption cpvalues, PetscBool mallocmatspace)
2272d71ae5a4SJacob Faibussowitsch {
2273ed73aabaSBarry Smith   Mat_SeqSELL *c = (Mat_SeqSELL *)C->data, *a = (Mat_SeqSELL *)A->data;
2274d4002b98SHong Zhang   PetscInt     i, m                           = A->rmap->n;
2275d4002b98SHong Zhang   PetscInt     totalslices = a->totalslices;
2276d4002b98SHong Zhang 
2277d4002b98SHong Zhang   PetscFunctionBegin;
2278d4002b98SHong Zhang   C->factortype = A->factortype;
2279f4259b30SLisandro Dalcin   c->row        = NULL;
2280f4259b30SLisandro Dalcin   c->col        = NULL;
2281f4259b30SLisandro Dalcin   c->icol       = NULL;
2282d4002b98SHong Zhang   c->reallocs   = 0;
2283d4002b98SHong Zhang   C->assembled  = PETSC_TRUE;
2284d4002b98SHong Zhang 
22859566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->rmap, &C->rmap));
22869566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &C->cmap));
2287d4002b98SHong Zhang 
228807e43b41SHong Zhang   PetscCall(PetscMalloc1(a->sliceheight * totalslices, &c->rlen));
22899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(totalslices + 1, &c->sliidx));
2290d4002b98SHong Zhang 
2291d4002b98SHong Zhang   for (i = 0; i < m; i++) c->rlen[i] = a->rlen[i];
2292d4002b98SHong Zhang   for (i = 0; i < totalslices + 1; i++) c->sliidx[i] = a->sliidx[i];
2293d4002b98SHong Zhang 
2294d4002b98SHong Zhang   /* allocate the matrix space */
2295d4002b98SHong Zhang   if (mallocmatspace) {
22969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(a->maxallocmat, &c->val, a->maxallocmat, &c->colidx));
2297d4002b98SHong Zhang 
2298d4002b98SHong Zhang     c->singlemalloc = PETSC_TRUE;
2299d4002b98SHong Zhang 
2300d4002b98SHong Zhang     if (m > 0) {
23019566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->colidx, a->colidx, a->maxallocmat));
2302d4002b98SHong Zhang       if (cpvalues == MAT_COPY_VALUES) {
23039566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(c->val, a->val, a->maxallocmat));
2304d4002b98SHong Zhang       } else {
23059566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(c->val, a->maxallocmat));
2306d4002b98SHong Zhang       }
2307d4002b98SHong Zhang     }
2308d4002b98SHong Zhang   }
2309d4002b98SHong Zhang 
2310d4002b98SHong Zhang   c->ignorezeroentries = a->ignorezeroentries;
2311d4002b98SHong Zhang   c->roworiented       = a->roworiented;
2312d4002b98SHong Zhang   c->nonew             = a->nonew;
2313d4002b98SHong Zhang   if (a->diag) {
23149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &c->diag));
2315ad540459SPierre Jolivet     for (i = 0; i < m; i++) c->diag[i] = a->diag[i];
2316f4259b30SLisandro Dalcin   } else c->diag = NULL;
2317d4002b98SHong Zhang 
2318f4259b30SLisandro Dalcin   c->solve_work         = NULL;
2319f4259b30SLisandro Dalcin   c->saved_values       = NULL;
2320f4259b30SLisandro Dalcin   c->idiag              = NULL;
2321f4259b30SLisandro Dalcin   c->ssor_work          = NULL;
2322d4002b98SHong Zhang   c->keepnonzeropattern = a->keepnonzeropattern;
2323d4002b98SHong Zhang   c->free_val           = PETSC_TRUE;
2324d4002b98SHong Zhang   c->free_colidx        = PETSC_TRUE;
2325d4002b98SHong Zhang 
2326d4002b98SHong Zhang   c->maxallocmat  = a->maxallocmat;
2327d4002b98SHong Zhang   c->maxallocrow  = a->maxallocrow;
2328d4002b98SHong Zhang   c->rlenmax      = a->rlenmax;
2329d4002b98SHong Zhang   c->nz           = a->nz;
2330d4002b98SHong Zhang   C->preallocated = PETSC_TRUE;
2331d4002b98SHong Zhang 
2332d4002b98SHong Zhang   c->nonzerorowcnt = a->nonzerorowcnt;
2333d4002b98SHong Zhang   C->nonzerostate  = A->nonzerostate;
2334d4002b98SHong Zhang 
23359566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist, &((PetscObject)C)->qlist));
23363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2337d4002b98SHong Zhang }
2338d4002b98SHong Zhang 
2339d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_SeqSELL(Mat A, MatDuplicateOption cpvalues, Mat *B)
2340d71ae5a4SJacob Faibussowitsch {
2341d4002b98SHong Zhang   PetscFunctionBegin;
23429566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
23439566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n));
234448a46eb9SPierre Jolivet   if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) PetscCall(MatSetBlockSizesFromMats(*B, A, A));
23459566063dSJacob Faibussowitsch   PetscCall(MatSetType(*B, ((PetscObject)A)->type_name));
23469566063dSJacob Faibussowitsch   PetscCall(MatDuplicateNoCreate_SeqSELL(*B, A, cpvalues, PETSC_TRUE));
23473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2348d4002b98SHong Zhang }
2349d4002b98SHong Zhang 
2350ed73aabaSBarry Smith /*MC
2351ed73aabaSBarry Smith    MATSEQSELL - MATSEQSELL = "seqsell" - A matrix type to be used for sequential sparse matrices,
2352ed73aabaSBarry Smith    based on the sliced Ellpack format
2353ed73aabaSBarry Smith 
235420f4b53cSBarry Smith    Options Database Key:
235511a5261eSBarry Smith . -mat_type seqsell - sets the matrix type to "`MATSEQELL` during a call to `MatSetFromOptions()`
2356ed73aabaSBarry Smith 
2357ed73aabaSBarry Smith    Level: beginner
2358ed73aabaSBarry Smith 
235967be906fSBarry Smith .seealso: `Mat`, `MatCreateSeqSell()`, `MATSELL`, `MATMPISELL`, `MATSEQAIJ`, `MATAIJ`, `MATMPIAIJ`
2360ed73aabaSBarry Smith M*/
2361ed73aabaSBarry Smith 
2362ed73aabaSBarry Smith /*MC
2363ed73aabaSBarry Smith    MATSELL - MATSELL = "sell" - A matrix type to be used for sparse matrices.
2364ed73aabaSBarry Smith 
236511a5261eSBarry Smith    This matrix type is identical to `MATSEQSELL` when constructed with a single process communicator,
236611a5261eSBarry Smith    and `MATMPISELL` otherwise.  As a result, for single process communicators,
236711a5261eSBarry Smith   `MatSeqSELLSetPreallocation()` is supported, and similarly `MatMPISELLSetPreallocation()` is supported
2368ed73aabaSBarry Smith   for communicators controlling multiple processes.  It is recommended that you call both of
2369ed73aabaSBarry Smith   the above preallocation routines for simplicity.
2370ed73aabaSBarry Smith 
237120f4b53cSBarry Smith    Options Database Key:
2372ed73aabaSBarry Smith . -mat_type sell - sets the matrix type to "sell" during a call to MatSetFromOptions()
2373ed73aabaSBarry Smith 
2374ed73aabaSBarry Smith   Level: beginner
2375ed73aabaSBarry Smith 
2376ed73aabaSBarry Smith   Notes:
2377ed73aabaSBarry Smith    This format is only supported for real scalars, double precision, and 32 bit indices (the defaults).
2378ed73aabaSBarry Smith 
2379ed73aabaSBarry Smith    It can provide better performance on Intel and AMD processes with AVX2 or AVX512 support for matrices that have a similar number of
2380ed73aabaSBarry Smith    non-zeros in contiguous groups of rows. However if the computation is memory bandwidth limited it may not provide much improvement.
2381ed73aabaSBarry Smith 
2382ed73aabaSBarry Smith   Developer Notes:
2383ed73aabaSBarry Smith    On Intel (and AMD) systems some of the matrix operations use SIMD (AVX) instructions to achieve higher performance.
2384ed73aabaSBarry Smith 
2385ed73aabaSBarry Smith    The sparse matrix format is as follows. For simplicity we assume a slice size of 2, it is actually 8
2386ed73aabaSBarry Smith .vb
2387ed73aabaSBarry Smith                             (2 0  3 4)
2388ed73aabaSBarry Smith    Consider the matrix A =  (5 0  6 0)
2389ed73aabaSBarry Smith                             (0 0  7 8)
2390ed73aabaSBarry Smith                             (0 0  9 9)
2391ed73aabaSBarry Smith 
2392ed73aabaSBarry Smith    symbolically the Ellpack format can be written as
2393ed73aabaSBarry Smith 
2394ed73aabaSBarry Smith         (2 3 4 |)           (0 2 3 |)
2395ed73aabaSBarry Smith    v =  (5 6 0 |)  colidx = (0 2 2 |)
2396ed73aabaSBarry Smith         --------            ---------
2397ed73aabaSBarry Smith         (7 8 |)             (2 3 |)
2398ed73aabaSBarry Smith         (9 9 |)             (2 3 |)
2399ed73aabaSBarry Smith 
2400ed73aabaSBarry 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).
2401ed73aabaSBarry 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
2402ed73aabaSBarry 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.
2403ed73aabaSBarry Smith 
2404ed73aabaSBarry 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)
2405ed73aabaSBarry Smith 
2406ed73aabaSBarry Smith .ve
2407ed73aabaSBarry Smith 
2408ed73aabaSBarry Smith       See MatMult_SeqSELL() for how this format is used with the SIMD operations to achieve high performance.
2409ed73aabaSBarry Smith 
2410ed73aabaSBarry Smith  References:
2411606c0280SSatish Balay . * - Hong Zhang, Richard T. Mills, Karl Rupp, and Barry F. Smith, Vectorized Parallel Sparse Matrix-Vector Multiplication in {PETSc} Using {AVX-512},
2412ed73aabaSBarry Smith    Proceedings of the 47th International Conference on Parallel Processing, 2018.
2413ed73aabaSBarry Smith 
241467be906fSBarry Smith .seealso: `Mat`, `MatCreateSeqSELL()`, `MatCreateSeqAIJ()`, `MatCreateSell()`, `MATSEQSELL`, `MATMPISELL`, `MATSEQAIJ`, `MATMPIAIJ`, `MATAIJ`
2415ed73aabaSBarry Smith M*/
2416ed73aabaSBarry Smith 
2417d4002b98SHong Zhang /*@C
241811a5261eSBarry Smith        MatCreateSeqSELL - Creates a sparse matrix in `MATSEQSELL` format.
2419d4002b98SHong Zhang 
24202ef1f0ffSBarry Smith  Collective
2421d4002b98SHong Zhang 
2422d4002b98SHong Zhang  Input Parameters:
242311a5261eSBarry Smith +  comm - MPI communicator, set to `PETSC_COMM_SELF`
2424d4002b98SHong Zhang .  m - number of rows
2425d4002b98SHong Zhang .  n - number of columns
242620f4b53cSBarry Smith .  rlenmax - maximum number of nonzeros in a row, ignored if `rlen` is provided
242720f4b53cSBarry Smith -  rlen - array containing the number of nonzeros in the various rows (possibly different for each row) or NULL
2428d4002b98SHong Zhang 
2429d4002b98SHong Zhang  Output Parameter:
2430d4002b98SHong Zhang .  A - the matrix
2431d4002b98SHong Zhang 
243220f4b53cSBarry Smith  Level: intermediate
243320f4b53cSBarry Smith 
243420f4b53cSBarry Smith  Notes:
243511a5261eSBarry Smith  It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`,
2436f6f02116SRichard Tran Mills  MatXXXXSetPreallocation() paradigm instead of this routine directly.
243711a5261eSBarry Smith  [MatXXXXSetPreallocation() is, for example, `MatSeqSELLSetPreallocation()`]
2438d4002b98SHong Zhang 
243920f4b53cSBarry Smith  Specify the preallocated storage with either `rlenmax` or `rlen` (not both).
244020f4b53cSBarry Smith  Set `rlenmax` = `PETSC_DEFAULT` and `rlen` = `NULL` for PETSc to control dynamic memory
244120f4b53cSBarry Smith  allocation.
2442d4002b98SHong Zhang 
244367be906fSBarry Smith  .seealso: `Mat`, `MATSEQSELL`, `MatCreate()`, `MatCreateSELL()`, `MatSetValues()`, `MatSeqSELLSetPreallocation()`, `MATSELL`, `MATSEQSELL`, `MATMPISELL`
2444d4002b98SHong Zhang  @*/
244520f4b53cSBarry Smith PetscErrorCode MatCreateSeqSELL(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt rlenmax, const PetscInt rlen[], Mat *A)
2446d71ae5a4SJacob Faibussowitsch {
2447d4002b98SHong Zhang   PetscFunctionBegin;
24489566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
24499566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, n, m, n));
24509566063dSJacob Faibussowitsch   PetscCall(MatSetType(*A, MATSEQSELL));
245120f4b53cSBarry Smith   PetscCall(MatSeqSELLSetPreallocation_SeqSELL(*A, rlenmax, rlen));
24523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2453d4002b98SHong Zhang }
2454d4002b98SHong Zhang 
2455d71ae5a4SJacob Faibussowitsch PetscErrorCode MatEqual_SeqSELL(Mat A, Mat B, PetscBool *flg)
2456d71ae5a4SJacob Faibussowitsch {
2457d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data, *b = (Mat_SeqSELL *)B->data;
2458d4002b98SHong Zhang   PetscInt     totalslices = a->totalslices;
2459d4002b98SHong Zhang 
2460d4002b98SHong Zhang   PetscFunctionBegin;
2461d4002b98SHong Zhang   /* If the  matrix dimensions are not equal,or no of nonzeros */
2462d4002b98SHong Zhang   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) || (a->nz != b->nz) || (a->rlenmax != b->rlenmax)) {
2463d4002b98SHong Zhang     *flg = PETSC_FALSE;
24643ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2465d4002b98SHong Zhang   }
2466d4002b98SHong Zhang   /* if the a->colidx are the same */
24679566063dSJacob Faibussowitsch   PetscCall(PetscArraycmp(a->colidx, b->colidx, a->sliidx[totalslices], flg));
24683ba16761SJacob Faibussowitsch   if (!*flg) PetscFunctionReturn(PETSC_SUCCESS);
2469d4002b98SHong Zhang   /* if a->val are the same */
24709566063dSJacob Faibussowitsch   PetscCall(PetscArraycmp(a->val, b->val, a->sliidx[totalslices], flg));
24713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2472d4002b98SHong Zhang }
2473d4002b98SHong Zhang 
2474d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLInvalidateDiagonal(Mat A)
2475d71ae5a4SJacob Faibussowitsch {
2476d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
2477d4002b98SHong Zhang 
2478d4002b98SHong Zhang   PetscFunctionBegin;
2479d4002b98SHong Zhang   a->idiagvalid = PETSC_FALSE;
24803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2481d4002b98SHong Zhang }
2482d4002b98SHong Zhang 
2483d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConjugate_SeqSELL(Mat A)
2484d71ae5a4SJacob Faibussowitsch {
2485d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
2486d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
2487d4002b98SHong Zhang   PetscInt     i;
2488d4002b98SHong Zhang   PetscScalar *val = a->val;
2489d4002b98SHong Zhang 
2490d4002b98SHong Zhang   PetscFunctionBegin;
24912d1451d4SHong Zhang   for (i = 0; i < a->sliidx[a->totalslices]; i++) { val[i] = PetscConj(val[i]); }
24922d1451d4SHong Zhang   #if defined(PETSC_HAVE_CUDA)
24932d1451d4SHong Zhang   if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) A->offloadmask = PETSC_OFFLOAD_CPU;
24942d1451d4SHong Zhang   #endif
2495d4002b98SHong Zhang #else
2496d4002b98SHong Zhang   PetscFunctionBegin;
2497d4002b98SHong Zhang #endif
24983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2499d4002b98SHong Zhang }
2500