xref: /petsc/src/mat/impls/sell/seq/sell.c (revision aaa8cc7d2a5c3913edcbb923e20f154fe9c4aa65)
1d4002b98SHong Zhang 
2d4002b98SHong Zhang /*
3d4002b98SHong Zhang   Defines the basic matrix operations for the SELL matrix storage format.
4d4002b98SHong Zhang */
5d4002b98SHong Zhang #include <../src/mat/impls/sell/seq/sell.h> /*I   "petscmat.h"  I*/
6d4002b98SHong Zhang #include <petscblaslapack.h>
7d4002b98SHong Zhang #include <petsc/private/kernels/blocktranspose.h>
8ed73aabaSBarry Smith 
9ed73aabaSBarry Smith static PetscBool  cited      = PETSC_FALSE;
109371c9d4SSatish Balay static const char citation[] = "@inproceedings{ZhangELLPACK2018,\n"
11ed73aabaSBarry Smith                                " author = {Hong Zhang and Richard T. Mills and Karl Rupp and Barry F. Smith},\n"
12ed73aabaSBarry Smith                                " title = {Vectorized Parallel Sparse Matrix-Vector Multiplication in {PETSc} Using {AVX-512}},\n"
13ed73aabaSBarry Smith                                " booktitle = {Proceedings of the 47th International Conference on Parallel Processing},\n"
14ed73aabaSBarry Smith                                " year = 2018\n"
15ed73aabaSBarry Smith                                "}\n";
16ed73aabaSBarry Smith 
175f70456aSHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && (defined(__AVX512F__) || (defined(__AVX2__) && defined(__FMA__)) || defined(__AVX__)) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
184243e2ceSHong Zhang 
19d4002b98SHong Zhang   #include <immintrin.h>
20d4002b98SHong Zhang 
21d4002b98SHong Zhang   #if !defined(_MM_SCALE_8)
22d4002b98SHong Zhang     #define _MM_SCALE_8 8
23d4002b98SHong Zhang   #endif
24d4002b98SHong Zhang 
25d4002b98SHong Zhang   #if defined(__AVX512F__)
26d4002b98SHong Zhang     /* these do not work
27d4002b98SHong Zhang    vec_idx  = _mm512_loadunpackhi_epi32(vec_idx,acolidx);
28d4002b98SHong Zhang    vec_vals = _mm512_loadunpackhi_pd(vec_vals,aval);
29d4002b98SHong Zhang   */
30d4002b98SHong Zhang     #define AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y) \
31d4002b98SHong Zhang       /* if the mask bit is set, copy from acolidx, otherwise from vec_idx */ \
32ef588d5cSRichard Tran Mills       vec_idx  = _mm256_loadu_si256((__m256i const *)acolidx); \
33ef588d5cSRichard Tran Mills       vec_vals = _mm512_loadu_pd(aval); \
34d4002b98SHong Zhang       vec_x    = _mm512_i32gather_pd(vec_idx, x, _MM_SCALE_8); \
35a48a6482SHong Zhang       vec_y    = _mm512_fmadd_pd(vec_x, vec_vals, vec_y)
365f70456aSHong Zhang   #elif defined(__AVX2__) && defined(__FMA__)
37a48a6482SHong Zhang     #define AVX2_Mult_Private(vec_idx, vec_x, vec_vals, vec_y) \
38ef588d5cSRichard Tran Mills       vec_vals = _mm256_loadu_pd(aval); \
39ef588d5cSRichard Tran Mills       vec_idx  = _mm_loadu_si128((__m128i const *)acolidx); /* SSE2 */ \
40a48a6482SHong Zhang       vec_x    = _mm256_i32gather_pd(x, vec_idx, _MM_SCALE_8); \
41a48a6482SHong Zhang       vec_y    = _mm256_fmadd_pd(vec_x, vec_vals, vec_y)
42d4002b98SHong Zhang   #endif
43d4002b98SHong Zhang #endif /* PETSC_HAVE_IMMINTRIN_H */
44d4002b98SHong Zhang 
45d4002b98SHong Zhang /*@C
46d4002b98SHong Zhang  MatSeqSELLSetPreallocation - For good matrix assembly performance
4720f4b53cSBarry Smith  the user should preallocate the matrix storage by setting the parameter `nz`
4820f4b53cSBarry Smith  (or the array `nnz`).
49d4002b98SHong Zhang 
50d083f849SBarry Smith  Collective
51d4002b98SHong Zhang 
52d4002b98SHong Zhang  Input Parameters:
5311a5261eSBarry Smith +  B - The `MATSEQSELL` matrix
5420f4b53cSBarry Smith .  rlenmax - number of nonzeros per row (same for all rows), ignored if `rlen` is provided
5520f4b53cSBarry Smith -  rlen - array containing the number of nonzeros in the various rows (possibly different for each row) or `NULL`
5667be906fSBarry Smith 
5767be906fSBarry Smith  Level: intermediate
58d4002b98SHong Zhang 
59d4002b98SHong Zhang  Notes:
6067be906fSBarry Smith  Specify the preallocated storage with either `rlenmax` or `rlen` (not both).
6167be906fSBarry Smith  Set `rlenmax` = `PETSC_DEFAULT` and `rlen` = `NULL` for PETSc to control dynamic memory
6220f4b53cSBarry Smith  allocation.
63d4002b98SHong Zhang 
6411a5261eSBarry Smith  You can call `MatGetInfo()` to get information on how effective the preallocation was;
65d4002b98SHong Zhang  for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
6667be906fSBarry Smith  You can also run with the option `-info` and look for messages with the string
67d4002b98SHong Zhang  malloc in them to see if additional memory allocation was needed.
68d4002b98SHong Zhang 
6927430b45SBarry Smith  Developer Note:
7067be906fSBarry Smith  Use `rlenmax` of `MAT_SKIP_ALLOCATION` to not allocate any space for the matrix
71d4002b98SHong Zhang  entries or columns indices.
72d4002b98SHong Zhang 
73c7ee91abSRichard Tran Mills  The maximum number of nonzeos in any row should be as accurate as possible.
74c7ee91abSRichard Tran Mills  If it is underestimated, you will get bad performance due to reallocation
7567be906fSBarry Smith  (`MatSeqXSELLReallocateSELL()`).
76d4002b98SHong Zhang 
7767be906fSBarry Smith  .seealso: `Mat`, `MATSEQSELL`, `MATSELL`, `MatCreate()`, `MatCreateSELL()`, `MatSetValues()`, `MatGetInfo()`
78d4002b98SHong Zhang  @*/
79d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLSetPreallocation(Mat B, PetscInt rlenmax, const PetscInt rlen[])
80d71ae5a4SJacob Faibussowitsch {
81d4002b98SHong Zhang   PetscFunctionBegin;
82d4002b98SHong Zhang   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
83d4002b98SHong Zhang   PetscValidType(B, 1);
84cac4c232SBarry Smith   PetscTryMethod(B, "MatSeqSELLSetPreallocation_C", (Mat, PetscInt, const PetscInt[]), (B, rlenmax, rlen));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86d4002b98SHong Zhang }
87d4002b98SHong Zhang 
88d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLSetPreallocation_SeqSELL(Mat B, PetscInt maxallocrow, const PetscInt rlen[])
89d71ae5a4SJacob Faibussowitsch {
90d4002b98SHong Zhang   Mat_SeqSELL *b;
91d4002b98SHong Zhang   PetscInt     i, j, totalslices;
92d4002b98SHong Zhang   PetscBool    skipallocation = PETSC_FALSE, realalloc = PETSC_FALSE;
93d4002b98SHong Zhang 
94d4002b98SHong Zhang   PetscFunctionBegin;
95d4002b98SHong Zhang   if (maxallocrow >= 0 || rlen) realalloc = PETSC_TRUE;
96d4002b98SHong Zhang   if (maxallocrow == MAT_SKIP_ALLOCATION) {
97d4002b98SHong Zhang     skipallocation = PETSC_TRUE;
98d4002b98SHong Zhang     maxallocrow    = 0;
99d4002b98SHong Zhang   }
100d4002b98SHong Zhang 
1019566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
1029566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
103d4002b98SHong Zhang 
104d4002b98SHong Zhang   /* FIXME: if one preallocates more space than needed, the matrix does not shrink automatically, but for best performance it should */
105d4002b98SHong Zhang   if (maxallocrow == PETSC_DEFAULT || maxallocrow == PETSC_DECIDE) maxallocrow = 5;
10608401ef6SPierre Jolivet   PetscCheck(maxallocrow >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "maxallocrow cannot be less than 0: value %" PetscInt_FMT, maxallocrow);
107d4002b98SHong Zhang   if (rlen) {
108d4002b98SHong Zhang     for (i = 0; i < B->rmap->n; i++) {
10908401ef6SPierre Jolivet       PetscCheck(rlen[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "rlen cannot be less than 0: local row %" PetscInt_FMT " value %" PetscInt_FMT, i, rlen[i]);
11008401ef6SPierre Jolivet       PetscCheck(rlen[i] <= B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "rlen cannot be greater than row length: local row %" PetscInt_FMT " value %" PetscInt_FMT " rowlength %" PetscInt_FMT, i, rlen[i], B->cmap->n);
111d4002b98SHong Zhang     }
112d4002b98SHong Zhang   }
113d4002b98SHong Zhang 
114d4002b98SHong Zhang   B->preallocated = PETSC_TRUE;
115d4002b98SHong Zhang 
116d4002b98SHong Zhang   b = (Mat_SeqSELL *)B->data;
117d4002b98SHong Zhang 
118faa75363SBarry Smith   totalslices    = PetscCeilInt(B->rmap->n, 8);
119d4002b98SHong Zhang   b->totalslices = totalslices;
120d4002b98SHong Zhang   if (!skipallocation) {
1219566063dSJacob Faibussowitsch     if (B->rmap->n & 0x07) PetscCall(PetscInfo(B, "Padding rows to the SEQSELL matrix because the number of rows is not the multiple of 8 (value %" PetscInt_FMT ")\n", B->rmap->n));
122d4002b98SHong Zhang 
123d4002b98SHong Zhang     if (!b->sliidx) { /* sliidx gives the starting index of each slice, the last element is the total space allocated */
1249566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(totalslices + 1, &b->sliidx));
125d4002b98SHong Zhang     }
126d4002b98SHong Zhang     if (!rlen) { /* if rlen is not provided, allocate same space for all the slices */
127d4002b98SHong Zhang       if (maxallocrow == PETSC_DEFAULT || maxallocrow == PETSC_DECIDE) maxallocrow = 10;
128d4002b98SHong Zhang       else if (maxallocrow < 0) maxallocrow = 1;
129d4002b98SHong Zhang       for (i = 0; i <= totalslices; i++) b->sliidx[i] = i * 8 * maxallocrow;
130d4002b98SHong Zhang     } else {
131d4002b98SHong Zhang       maxallocrow  = 0;
132d4002b98SHong Zhang       b->sliidx[0] = 0;
133d4002b98SHong Zhang       for (i = 1; i < totalslices; i++) {
134d4002b98SHong Zhang         b->sliidx[i] = 0;
135ad540459SPierre Jolivet         for (j = 0; j < 8; j++) b->sliidx[i] = PetscMax(b->sliidx[i], rlen[8 * (i - 1) + j]);
136d4002b98SHong Zhang         maxallocrow = PetscMax(b->sliidx[i], maxallocrow);
1379566063dSJacob Faibussowitsch         PetscCall(PetscIntSumError(b->sliidx[i - 1], 8 * b->sliidx[i], &b->sliidx[i]));
138d4002b98SHong Zhang       }
139d4002b98SHong Zhang       /* last slice */
140d4002b98SHong Zhang       b->sliidx[totalslices] = 0;
141d4002b98SHong Zhang       for (j = (totalslices - 1) * 8; j < B->rmap->n; j++) b->sliidx[totalslices] = PetscMax(b->sliidx[totalslices], rlen[j]);
142d4002b98SHong Zhang       maxallocrow            = PetscMax(b->sliidx[totalslices], maxallocrow);
143d4002b98SHong Zhang       b->sliidx[totalslices] = b->sliidx[totalslices - 1] + 8 * b->sliidx[totalslices];
144d4002b98SHong Zhang     }
145d4002b98SHong Zhang 
146d4002b98SHong Zhang     /* allocate space for val, colidx, rlen */
147d4002b98SHong Zhang     /* FIXME: should B's old memory be unlogged? */
1489566063dSJacob Faibussowitsch     PetscCall(MatSeqXSELLFreeSELL(B, &b->val, &b->colidx));
149d4002b98SHong Zhang     /* FIXME: assuming an element of the bit array takes 8 bits */
1509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(b->sliidx[totalslices], &b->val, b->sliidx[totalslices], &b->colidx));
151d4002b98SHong 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. */
1529566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(8 * totalslices, &b->rlen));
153d4002b98SHong Zhang 
154d4002b98SHong Zhang     b->singlemalloc = PETSC_TRUE;
155d4002b98SHong Zhang     b->free_val     = PETSC_TRUE;
156d4002b98SHong Zhang     b->free_colidx  = PETSC_TRUE;
157d4002b98SHong Zhang   } else {
158d4002b98SHong Zhang     b->free_val    = PETSC_FALSE;
159d4002b98SHong Zhang     b->free_colidx = PETSC_FALSE;
160d4002b98SHong Zhang   }
161d4002b98SHong Zhang 
162d4002b98SHong Zhang   b->nz               = 0;
163d4002b98SHong Zhang   b->maxallocrow      = maxallocrow;
164d4002b98SHong Zhang   b->rlenmax          = maxallocrow;
165d4002b98SHong Zhang   b->maxallocmat      = b->sliidx[totalslices];
166d4002b98SHong Zhang   B->info.nz_unneeded = (double)b->maxallocmat;
1671baa6e33SBarry Smith   if (realalloc) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
169d4002b98SHong Zhang }
170d4002b98SHong Zhang 
171d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetRow_SeqSELL(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
172d71ae5a4SJacob Faibussowitsch {
1736108893eSStefano Zampini   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1746108893eSStefano Zampini   PetscInt     shift;
1756108893eSStefano Zampini 
1766108893eSStefano Zampini   PetscFunctionBegin;
177aed4548fSBarry Smith   PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row %" PetscInt_FMT " out of range", row);
1786108893eSStefano Zampini   if (nz) *nz = a->rlen[row];
1796108893eSStefano Zampini   shift = a->sliidx[row >> 3] + (row & 0x07);
18048a46eb9SPierre Jolivet   if (!a->getrowcols) PetscCall(PetscMalloc2(a->rlenmax, &a->getrowcols, a->rlenmax, &a->getrowvals));
1816108893eSStefano Zampini   if (idx) {
1826108893eSStefano Zampini     PetscInt j;
1836108893eSStefano Zampini     for (j = 0; j < a->rlen[row]; j++) a->getrowcols[j] = a->colidx[shift + 8 * j];
1846108893eSStefano Zampini     *idx = a->getrowcols;
1856108893eSStefano Zampini   }
1866108893eSStefano Zampini   if (v) {
1876108893eSStefano Zampini     PetscInt j;
1886108893eSStefano Zampini     for (j = 0; j < a->rlen[row]; j++) a->getrowvals[j] = a->val[shift + 8 * j];
1896108893eSStefano Zampini     *v = a->getrowvals;
1906108893eSStefano Zampini   }
1913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1926108893eSStefano Zampini }
1936108893eSStefano Zampini 
194d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRestoreRow_SeqSELL(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
195d71ae5a4SJacob Faibussowitsch {
1966108893eSStefano Zampini   PetscFunctionBegin;
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1986108893eSStefano Zampini }
1996108893eSStefano Zampini 
200d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_SeqSELL_SeqAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
201d71ae5a4SJacob Faibussowitsch {
202d4002b98SHong Zhang   Mat          B;
203d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
204e3f1f374SStefano Zampini   PetscInt     i;
205d4002b98SHong Zhang 
206d4002b98SHong Zhang   PetscFunctionBegin;
207ad013a7bSRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
208ad013a7bSRichard Tran Mills     B = *newmat;
2099566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(B));
210ad013a7bSRichard Tran Mills   } else {
2119566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2129566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
2139566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQAIJ));
2149566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, a->rlen));
215ad013a7bSRichard Tran Mills   }
216d4002b98SHong Zhang 
217e3f1f374SStefano Zampini   for (i = 0; i < A->rmap->n; i++) {
218e108cb99SStefano Zampini     PetscInt     nz = 0, *cols = NULL;
219e108cb99SStefano Zampini     PetscScalar *vals = NULL;
220e3f1f374SStefano Zampini 
2219566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqSELL(A, i, &nz, &cols, &vals));
2229566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &i, nz, cols, vals, INSERT_VALUES));
2239566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow_SeqSELL(A, i, &nz, &cols, &vals));
224d4002b98SHong Zhang   }
225e3f1f374SStefano Zampini 
2269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
228d4002b98SHong Zhang   B->rmap->bs = A->rmap->bs;
229d4002b98SHong Zhang 
230d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
2319566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
232d4002b98SHong Zhang   } else {
233d4002b98SHong Zhang     *newmat = B;
234d4002b98SHong Zhang   }
2353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
236d4002b98SHong Zhang }
237d4002b98SHong Zhang 
238d4002b98SHong Zhang #include <../src/mat/impls/aij/seq/aij.h>
239d4002b98SHong Zhang 
240d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_SeqAIJ_SeqSELL(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
241d71ae5a4SJacob Faibussowitsch {
242d4002b98SHong Zhang   Mat                B;
243d4002b98SHong Zhang   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *)A->data;
244d4002b98SHong Zhang   PetscInt          *ai = a->i, m = A->rmap->N, n = A->cmap->N, i, *rowlengths, row, ncols;
245d4002b98SHong Zhang   const PetscInt    *cols;
246d4002b98SHong Zhang   const PetscScalar *vals;
247d4002b98SHong Zhang 
248d4002b98SHong Zhang   PetscFunctionBegin;
249ad013a7bSRichard Tran Mills 
250ad013a7bSRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
251ad013a7bSRichard Tran Mills     B = *newmat;
252ad013a7bSRichard Tran Mills   } else {
253d5e5b2e5SBarry Smith     if (PetscDefined(USE_DEBUG) || !a->ilen) {
2549566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m, &rowlengths));
255ad540459SPierre Jolivet       for (i = 0; i < m; i++) rowlengths[i] = ai[i + 1] - ai[i];
256d5e5b2e5SBarry Smith     }
257d5e5b2e5SBarry Smith     if (PetscDefined(USE_DEBUG) && a->ilen) {
258d5e5b2e5SBarry Smith       PetscBool eq;
2599566063dSJacob Faibussowitsch       PetscCall(PetscMemcmp(rowlengths, a->ilen, m * sizeof(PetscInt), &eq));
26028b400f6SJacob Faibussowitsch       PetscCheck(eq, PETSC_COMM_SELF, PETSC_ERR_PLIB, "SeqAIJ ilen array incorrect");
2619566063dSJacob Faibussowitsch       PetscCall(PetscFree(rowlengths));
262d5e5b2e5SBarry Smith       rowlengths = a->ilen;
263d5e5b2e5SBarry Smith     } else if (a->ilen) rowlengths = a->ilen;
2649566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
2659566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, m, n, m, n));
2669566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATSEQSELL));
2679566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(B, 0, rowlengths));
2689566063dSJacob Faibussowitsch     if (rowlengths != a->ilen) PetscCall(PetscFree(rowlengths));
269ad013a7bSRichard Tran Mills   }
270d4002b98SHong Zhang 
271d4002b98SHong Zhang   for (row = 0; row < m; row++) {
2729566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqAIJ(A, row, &ncols, (PetscInt **)&cols, (PetscScalar **)&vals));
2739566063dSJacob Faibussowitsch     PetscCall(MatSetValues_SeqSELL(B, 1, &row, ncols, cols, vals, INSERT_VALUES));
2749566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow_SeqAIJ(A, row, &ncols, (PetscInt **)&cols, (PetscScalar **)&vals));
275d4002b98SHong Zhang   }
2769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
278d4002b98SHong Zhang   B->rmap->bs = A->rmap->bs;
279d4002b98SHong Zhang 
280d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
2819566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
282d4002b98SHong Zhang   } else {
283d4002b98SHong Zhang     *newmat = B;
284d4002b98SHong Zhang   }
2853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
286d4002b98SHong Zhang }
287d4002b98SHong Zhang 
288d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_SeqSELL(Mat A, Vec xx, Vec yy)
289d71ae5a4SJacob Faibussowitsch {
290d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
291d4002b98SHong Zhang   PetscScalar       *y;
292d4002b98SHong Zhang   const PetscScalar *x;
293d4002b98SHong Zhang   const MatScalar   *aval        = a->val;
294d4002b98SHong Zhang   PetscInt           totalslices = a->totalslices;
295d4002b98SHong Zhang   const PetscInt    *acolidx     = a->colidx;
2967285fed1SHong Zhang   PetscInt           i, j;
297d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
298d4002b98SHong Zhang   __m512d  vec_x, vec_y, vec_vals;
299d4002b98SHong Zhang   __m256i  vec_idx;
300d4002b98SHong Zhang   __mmask8 mask;
301d4002b98SHong Zhang   __m512d  vec_x2, vec_y2, vec_vals2, vec_x3, vec_y3, vec_vals3, vec_x4, vec_y4, vec_vals4;
302d4002b98SHong Zhang   __m256i  vec_idx2, vec_idx3, vec_idx4;
3035f70456aSHong 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)
304a48a6482SHong Zhang   __m128i   vec_idx;
305a48a6482SHong Zhang   __m256d   vec_x, vec_y, vec_y2, vec_vals;
306a48a6482SHong Zhang   MatScalar yval;
307a48a6482SHong Zhang   PetscInt  r, rows_left, row, nnz_in_row;
30821cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
309d4002b98SHong Zhang   __m128d   vec_x_tmp;
310d4002b98SHong Zhang   __m256d   vec_x, vec_y, vec_y2, vec_vals;
311d4002b98SHong Zhang   MatScalar yval;
312d4002b98SHong Zhang   PetscInt  r, rows_left, row, nnz_in_row;
313d4002b98SHong Zhang #else
314d4002b98SHong Zhang   PetscScalar sum[8];
315d4002b98SHong Zhang #endif
316d4002b98SHong Zhang 
317d4002b98SHong Zhang #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
318d4002b98SHong Zhang   #pragma disjoint(*x, *y, *aval)
319d4002b98SHong Zhang #endif
320d4002b98SHong Zhang 
321d4002b98SHong Zhang   PetscFunctionBegin;
3229566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx, &x));
3239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(yy, &y));
324d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
325d4002b98SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
326d4002b98SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
327d4002b98SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
328d4002b98SHong Zhang 
329d4002b98SHong Zhang     vec_y  = _mm512_setzero_pd();
330d4002b98SHong Zhang     vec_y2 = _mm512_setzero_pd();
331d4002b98SHong Zhang     vec_y3 = _mm512_setzero_pd();
332d4002b98SHong Zhang     vec_y4 = _mm512_setzero_pd();
333d4002b98SHong Zhang 
334da81f932SPierre Jolivet     j = a->sliidx[i] >> 3; /* 8 bytes are read at each time, corresponding to a slice column */
335d4002b98SHong Zhang     switch ((a->sliidx[i + 1] - a->sliidx[i]) / 8 & 3) {
336d4002b98SHong Zhang     case 3:
337d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3389371c9d4SSatish Balay       acolidx += 8;
3399371c9d4SSatish Balay       aval += 8;
340d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
3419371c9d4SSatish Balay       acolidx += 8;
3429371c9d4SSatish Balay       aval += 8;
343d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
3449371c9d4SSatish Balay       acolidx += 8;
3459371c9d4SSatish Balay       aval += 8;
346d4002b98SHong Zhang       j += 3;
347d4002b98SHong Zhang       break;
348d4002b98SHong Zhang     case 2:
349d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3509371c9d4SSatish Balay       acolidx += 8;
3519371c9d4SSatish Balay       aval += 8;
352d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
3539371c9d4SSatish Balay       acolidx += 8;
3549371c9d4SSatish Balay       aval += 8;
355d4002b98SHong Zhang       j += 2;
356d4002b98SHong Zhang       break;
357d4002b98SHong Zhang     case 1:
358d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3599371c9d4SSatish Balay       acolidx += 8;
3609371c9d4SSatish Balay       aval += 8;
361d4002b98SHong Zhang       j += 1;
362d4002b98SHong Zhang       break;
363d4002b98SHong Zhang     }
364d4002b98SHong Zhang   #pragma novector
365d4002b98SHong Zhang     for (; j < (a->sliidx[i + 1] >> 3); j += 4) {
366d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
3679371c9d4SSatish Balay       acolidx += 8;
3689371c9d4SSatish Balay       aval += 8;
369d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
3709371c9d4SSatish Balay       acolidx += 8;
3719371c9d4SSatish Balay       aval += 8;
372d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
3739371c9d4SSatish Balay       acolidx += 8;
3749371c9d4SSatish Balay       aval += 8;
375d4002b98SHong Zhang       AVX512_Mult_Private(vec_idx4, vec_x4, vec_vals4, vec_y4);
3769371c9d4SSatish Balay       acolidx += 8;
3779371c9d4SSatish Balay       aval += 8;
378d4002b98SHong Zhang     }
379d4002b98SHong Zhang 
380d4002b98SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y2);
381d4002b98SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y3);
382d4002b98SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y4);
383d4002b98SHong Zhang     if (i == totalslices - 1 && A->rmap->n & 0x07) { /* if last slice has padding rows */
384d4002b98SHong Zhang       mask = (__mmask8)(0xff >> (8 - (A->rmap->n & 0x07)));
385ef588d5cSRichard Tran Mills       _mm512_mask_storeu_pd(&y[8 * i], mask, vec_y);
386d4002b98SHong Zhang     } else {
387ef588d5cSRichard Tran Mills       _mm512_storeu_pd(&y[8 * i], vec_y);
388d4002b98SHong Zhang     }
389d4002b98SHong Zhang   }
3905f70456aSHong 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)
391a48a6482SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over full slices */
392a48a6482SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
393a48a6482SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
394a48a6482SHong Zhang 
395a48a6482SHong Zhang     /* last slice may have padding rows. Don't use vectorization. */
396a48a6482SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
397a48a6482SHong Zhang       rows_left = A->rmap->n - 8 * i;
398a48a6482SHong Zhang       for (r = 0; r < rows_left; ++r) {
399a48a6482SHong Zhang         yval       = (MatScalar)0;
400a48a6482SHong Zhang         row        = 8 * i + r;
401a48a6482SHong Zhang         nnz_in_row = a->rlen[row];
402a48a6482SHong Zhang         for (j = 0; j < nnz_in_row; ++j) yval += aval[8 * j + r] * x[acolidx[8 * j + r]];
403a48a6482SHong Zhang         y[row] = yval;
404a48a6482SHong Zhang       }
405a48a6482SHong Zhang       break;
406a48a6482SHong Zhang     }
407a48a6482SHong Zhang 
408a48a6482SHong Zhang     vec_y  = _mm256_setzero_pd();
409a48a6482SHong Zhang     vec_y2 = _mm256_setzero_pd();
410a48a6482SHong Zhang 
411a48a6482SHong Zhang   /* Process slice of height 8 (512 bits) via two subslices of height 4 (256 bits) via AVX */
412a48a6482SHong Zhang   #pragma novector
413a48a6482SHong Zhang   #pragma unroll(2)
414a48a6482SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
415a48a6482SHong Zhang       AVX2_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
4169371c9d4SSatish Balay       aval += 4;
4179371c9d4SSatish Balay       acolidx += 4;
418a48a6482SHong Zhang       AVX2_Mult_Private(vec_idx, vec_x, vec_vals, vec_y2);
4199371c9d4SSatish Balay       aval += 4;
4209371c9d4SSatish Balay       acolidx += 4;
421a48a6482SHong Zhang     }
422a48a6482SHong Zhang 
423ef588d5cSRichard Tran Mills     _mm256_storeu_pd(y + i * 8, vec_y);
424ef588d5cSRichard Tran Mills     _mm256_storeu_pd(y + i * 8 + 4, vec_y2);
425a48a6482SHong Zhang   }
42621cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
427d4002b98SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over full slices */
428d4002b98SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
429d4002b98SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
430d4002b98SHong Zhang 
431d4002b98SHong Zhang     vec_y  = _mm256_setzero_pd();
432d4002b98SHong Zhang     vec_y2 = _mm256_setzero_pd();
433d4002b98SHong Zhang 
434d4002b98SHong Zhang     /* last slice may have padding rows. Don't use vectorization. */
435d4002b98SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
436d4002b98SHong Zhang       rows_left = A->rmap->n - 8 * i;
437d4002b98SHong Zhang       for (r = 0; r < rows_left; ++r) {
438d4002b98SHong Zhang         yval       = (MatScalar)0;
439d4002b98SHong Zhang         row        = 8 * i + r;
440d4002b98SHong Zhang         nnz_in_row = a->rlen[row];
441d4002b98SHong Zhang         for (j = 0; j < nnz_in_row; ++j) yval += aval[8 * j + r] * x[acolidx[8 * j + r]];
442d4002b98SHong Zhang         y[row] = yval;
443d4002b98SHong Zhang       }
444d4002b98SHong Zhang       break;
445d4002b98SHong Zhang     }
446d4002b98SHong Zhang 
447d4002b98SHong Zhang   /* Process slice of height 8 (512 bits) via two subslices of height 4 (256 bits) via AVX */
448a48a6482SHong Zhang   #pragma novector
449a48a6482SHong Zhang   #pragma unroll(2)
4507285fed1SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
451d4002b98SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
452165f9cc3SJed Brown       vec_x_tmp = _mm_setzero_pd();
453d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
454d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
455d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
456d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
457d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
458d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
459d4002b98SHong Zhang       vec_y     = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y);
460d4002b98SHong Zhang       aval += 4;
461d4002b98SHong Zhang 
462d4002b98SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
463d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
464d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
465d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
466d4002b98SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
467d4002b98SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
468d4002b98SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
469d4002b98SHong Zhang       vec_y2    = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y2);
470d4002b98SHong Zhang       aval += 4;
471d4002b98SHong Zhang     }
472d4002b98SHong Zhang 
473d4002b98SHong Zhang     _mm256_storeu_pd(y + i * 8, vec_y);
474d4002b98SHong Zhang     _mm256_storeu_pd(y + i * 8 + 4, vec_y2);
475d4002b98SHong Zhang   }
476d4002b98SHong Zhang #else
477d4002b98SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
478d4002b98SHong Zhang     for (j = 0; j < 8; j++) sum[j] = 0.0;
479d4002b98SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
480d4002b98SHong Zhang       sum[0] += aval[j] * x[acolidx[j]];
481d4002b98SHong Zhang       sum[1] += aval[j + 1] * x[acolidx[j + 1]];
482d4002b98SHong Zhang       sum[2] += aval[j + 2] * x[acolidx[j + 2]];
483d4002b98SHong Zhang       sum[3] += aval[j + 3] * x[acolidx[j + 3]];
484d4002b98SHong Zhang       sum[4] += aval[j + 4] * x[acolidx[j + 4]];
485d4002b98SHong Zhang       sum[5] += aval[j + 5] * x[acolidx[j + 5]];
486d4002b98SHong Zhang       sum[6] += aval[j + 6] * x[acolidx[j + 6]];
487d4002b98SHong Zhang       sum[7] += aval[j + 7] * x[acolidx[j + 7]];
488d4002b98SHong Zhang     }
489d4002b98SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) { /* if last slice has padding rows */
490d4002b98SHong Zhang       for (j = 0; j < (A->rmap->n & 0x07); j++) y[8 * i + j] = sum[j];
491d4002b98SHong Zhang     } else {
4927285fed1SHong Zhang       for (j = 0; j < 8; j++) y[8 * i + j] = sum[j];
493d4002b98SHong Zhang     }
494d4002b98SHong Zhang   }
495d4002b98SHong Zhang #endif
496d4002b98SHong Zhang 
4979566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->nz - a->nonzerorowcnt)); /* theoretical minimal FLOPs */
4989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx, &x));
4999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(yy, &y));
5003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
501d4002b98SHong Zhang }
502d4002b98SHong Zhang 
503d4002b98SHong Zhang #include <../src/mat/impls/aij/seq/ftn-kernels/fmultadd.h>
504d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAdd_SeqSELL(Mat A, Vec xx, Vec yy, Vec zz)
505d71ae5a4SJacob Faibussowitsch {
506d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
507d4002b98SHong Zhang   PetscScalar       *y, *z;
508d4002b98SHong Zhang   const PetscScalar *x;
509d4002b98SHong Zhang   const MatScalar   *aval        = a->val;
510d4002b98SHong Zhang   PetscInt           totalslices = a->totalslices;
511d4002b98SHong Zhang   const PetscInt    *acolidx     = a->colidx;
512d4002b98SHong Zhang   PetscInt           i, j;
513d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
5147285fed1SHong Zhang   __m512d  vec_x, vec_y, vec_vals;
515d4002b98SHong Zhang   __m256i  vec_idx;
516d4002b98SHong Zhang   __mmask8 mask;
5177285fed1SHong Zhang   __m512d  vec_x2, vec_y2, vec_vals2, vec_x3, vec_y3, vec_vals3, vec_x4, vec_y4, vec_vals4;
5187285fed1SHong Zhang   __m256i  vec_idx2, vec_idx3, vec_idx4;
51921cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
5207285fed1SHong Zhang   __m128d   vec_x_tmp;
5217285fed1SHong Zhang   __m256d   vec_x, vec_y, vec_y2, vec_vals;
5227285fed1SHong Zhang   MatScalar yval;
5237285fed1SHong Zhang   PetscInt  r, row, nnz_in_row;
524d4002b98SHong Zhang #else
525d4002b98SHong Zhang   PetscScalar sum[8];
526d4002b98SHong Zhang #endif
527d4002b98SHong Zhang 
528d4002b98SHong Zhang #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
529d4002b98SHong Zhang   #pragma disjoint(*x, *y, *aval)
530d4002b98SHong Zhang #endif
531d4002b98SHong Zhang 
532d4002b98SHong Zhang   PetscFunctionBegin;
5339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx, &x));
5349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayPair(yy, zz, &y, &z));
535d4002b98SHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX512F__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
5367285fed1SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
5377285fed1SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
5387285fed1SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
5397285fed1SHong Zhang 
540d4002b98SHong Zhang     if (i == totalslices - 1 && A->rmap->n & 0x07) { /* if last slice has padding rows */
541d4002b98SHong Zhang       mask  = (__mmask8)(0xff >> (8 - (A->rmap->n & 0x07)));
542ef588d5cSRichard Tran Mills       vec_y = _mm512_mask_loadu_pd(vec_y, mask, &y[8 * i]);
5437285fed1SHong Zhang     } else {
544ef588d5cSRichard Tran Mills       vec_y = _mm512_loadu_pd(&y[8 * i]);
5457285fed1SHong Zhang     }
5467285fed1SHong Zhang     vec_y2 = _mm512_setzero_pd();
5477285fed1SHong Zhang     vec_y3 = _mm512_setzero_pd();
5487285fed1SHong Zhang     vec_y4 = _mm512_setzero_pd();
5497285fed1SHong Zhang 
550da81f932SPierre Jolivet     j = a->sliidx[i] >> 3; /* 8 bytes are read at each time, corresponding to a slice column */
5517285fed1SHong Zhang     switch ((a->sliidx[i + 1] - a->sliidx[i]) / 8 & 3) {
5527285fed1SHong Zhang     case 3:
5537285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
5549371c9d4SSatish Balay       acolidx += 8;
5559371c9d4SSatish Balay       aval += 8;
5567285fed1SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
5579371c9d4SSatish Balay       acolidx += 8;
5589371c9d4SSatish Balay       aval += 8;
5597285fed1SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
5609371c9d4SSatish Balay       acolidx += 8;
5619371c9d4SSatish Balay       aval += 8;
5627285fed1SHong Zhang       j += 3;
5637285fed1SHong Zhang       break;
5647285fed1SHong Zhang     case 2:
5657285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
5669371c9d4SSatish Balay       acolidx += 8;
5679371c9d4SSatish Balay       aval += 8;
5687285fed1SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
5699371c9d4SSatish Balay       acolidx += 8;
5709371c9d4SSatish Balay       aval += 8;
5717285fed1SHong Zhang       j += 2;
5727285fed1SHong Zhang       break;
5737285fed1SHong Zhang     case 1:
5747285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
5759371c9d4SSatish Balay       acolidx += 8;
5769371c9d4SSatish Balay       aval += 8;
5777285fed1SHong Zhang       j += 1;
5787285fed1SHong Zhang       break;
5797285fed1SHong Zhang     }
5807285fed1SHong Zhang   #pragma novector
5817285fed1SHong Zhang     for (; j < (a->sliidx[i + 1] >> 3); j += 4) {
5827285fed1SHong Zhang       AVX512_Mult_Private(vec_idx, vec_x, vec_vals, vec_y);
5839371c9d4SSatish Balay       acolidx += 8;
5849371c9d4SSatish Balay       aval += 8;
5857285fed1SHong Zhang       AVX512_Mult_Private(vec_idx2, vec_x2, vec_vals2, vec_y2);
5869371c9d4SSatish Balay       acolidx += 8;
5879371c9d4SSatish Balay       aval += 8;
5887285fed1SHong Zhang       AVX512_Mult_Private(vec_idx3, vec_x3, vec_vals3, vec_y3);
5899371c9d4SSatish Balay       acolidx += 8;
5909371c9d4SSatish Balay       aval += 8;
5917285fed1SHong Zhang       AVX512_Mult_Private(vec_idx4, vec_x4, vec_vals4, vec_y4);
5929371c9d4SSatish Balay       acolidx += 8;
5939371c9d4SSatish Balay       aval += 8;
5947285fed1SHong Zhang     }
5957285fed1SHong Zhang 
5967285fed1SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y2);
5977285fed1SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y3);
5987285fed1SHong Zhang     vec_y = _mm512_add_pd(vec_y, vec_y4);
5997285fed1SHong Zhang     if (i == totalslices - 1 && A->rmap->n & 0x07) { /* if last slice has padding rows */
600ef588d5cSRichard Tran Mills       _mm512_mask_storeu_pd(&z[8 * i], mask, vec_y);
601d4002b98SHong Zhang     } else {
602ef588d5cSRichard Tran Mills       _mm512_storeu_pd(&z[8 * i], vec_y);
603d4002b98SHong Zhang     }
6047285fed1SHong Zhang   }
60521cec45eSHong Zhang #elif defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
6067285fed1SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over full slices */
6077285fed1SHong Zhang     PetscPrefetchBlock(acolidx, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
6087285fed1SHong Zhang     PetscPrefetchBlock(aval, a->sliidx[i + 1] - a->sliidx[i], 0, PETSC_PREFETCH_HINT_T0);
6097285fed1SHong Zhang 
6107285fed1SHong Zhang     /* last slice may have padding rows. Don't use vectorization. */
6117285fed1SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
6127285fed1SHong Zhang       for (r = 0; r < (A->rmap->n & 0x07); ++r) {
6137285fed1SHong Zhang         row        = 8 * i + r;
6147285fed1SHong Zhang         yval       = (MatScalar)0.0;
6157285fed1SHong Zhang         nnz_in_row = a->rlen[row];
6167285fed1SHong Zhang         for (j = 0; j < nnz_in_row; ++j) yval += aval[8 * j + r] * x[acolidx[8 * j + r]];
6177285fed1SHong Zhang         z[row] = y[row] + yval;
6187285fed1SHong Zhang       }
6197285fed1SHong Zhang       break;
6207285fed1SHong Zhang     }
6217285fed1SHong Zhang 
6227285fed1SHong Zhang     vec_y  = _mm256_loadu_pd(y + 8 * i);
6237285fed1SHong Zhang     vec_y2 = _mm256_loadu_pd(y + 8 * i + 4);
6247285fed1SHong Zhang 
6257285fed1SHong Zhang     /* Process slice of height 8 (512 bits) via two subslices of height 4 (256 bits) via AVX */
6267285fed1SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
6277285fed1SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
628165f9cc3SJed Brown       vec_x_tmp = _mm_setzero_pd();
6297285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6307285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
631165f9cc3SJed Brown       vec_x     = _mm256_setzero_pd();
6327285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
6337285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6347285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
6357285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
6367285fed1SHong Zhang       vec_y     = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y);
6377285fed1SHong Zhang       aval += 4;
6387285fed1SHong Zhang 
6397285fed1SHong Zhang       vec_vals  = _mm256_loadu_pd(aval);
6407285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6417285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
6427285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 0);
6437285fed1SHong Zhang       vec_x_tmp = _mm_loadl_pd(vec_x_tmp, x + *acolidx++);
6447285fed1SHong Zhang       vec_x_tmp = _mm_loadh_pd(vec_x_tmp, x + *acolidx++);
6457285fed1SHong Zhang       vec_x     = _mm256_insertf128_pd(vec_x, vec_x_tmp, 1);
6467285fed1SHong Zhang       vec_y2    = _mm256_add_pd(_mm256_mul_pd(vec_x, vec_vals), vec_y2);
6477285fed1SHong Zhang       aval += 4;
6487285fed1SHong Zhang     }
6497285fed1SHong Zhang 
6507285fed1SHong Zhang     _mm256_storeu_pd(z + i * 8, vec_y);
6517285fed1SHong Zhang     _mm256_storeu_pd(z + i * 8 + 4, vec_y2);
6527285fed1SHong Zhang   }
653d4002b98SHong Zhang #else
6547285fed1SHong Zhang   for (i = 0; i < totalslices; i++) { /* loop over slices */
6557285fed1SHong Zhang     for (j = 0; j < 8; j++) sum[j] = 0.0;
656d4002b98SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
657d4002b98SHong Zhang       sum[0] += aval[j] * x[acolidx[j]];
658d4002b98SHong Zhang       sum[1] += aval[j + 1] * x[acolidx[j + 1]];
659d4002b98SHong Zhang       sum[2] += aval[j + 2] * x[acolidx[j + 2]];
660d4002b98SHong Zhang       sum[3] += aval[j + 3] * x[acolidx[j + 3]];
661d4002b98SHong Zhang       sum[4] += aval[j + 4] * x[acolidx[j + 4]];
662d4002b98SHong Zhang       sum[5] += aval[j + 5] * x[acolidx[j + 5]];
663d4002b98SHong Zhang       sum[6] += aval[j + 6] * x[acolidx[j + 6]];
664d4002b98SHong Zhang       sum[7] += aval[j + 7] * x[acolidx[j + 7]];
665d4002b98SHong Zhang     }
6667285fed1SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
6677285fed1SHong Zhang       for (j = 0; j < (A->rmap->n & 0x07); j++) z[8 * i + j] = y[8 * i + j] + sum[j];
668d4002b98SHong Zhang     } else {
6697285fed1SHong Zhang       for (j = 0; j < 8; j++) z[8 * i + j] = y[8 * i + j] + sum[j];
6707285fed1SHong Zhang     }
671d4002b98SHong Zhang   }
672d4002b98SHong Zhang #endif
673d4002b98SHong Zhang 
6749566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->nz));
6759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx, &x));
6769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayPair(yy, zz, &y, &z));
6773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
678d4002b98SHong Zhang }
679d4002b98SHong Zhang 
680d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAdd_SeqSELL(Mat A, Vec xx, Vec zz, Vec yy)
681d71ae5a4SJacob Faibussowitsch {
682d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
683d4002b98SHong Zhang   PetscScalar       *y;
684d4002b98SHong Zhang   const PetscScalar *x;
685d4002b98SHong Zhang   const MatScalar   *aval    = a->val;
686d4002b98SHong Zhang   const PetscInt    *acolidx = a->colidx;
6877285fed1SHong Zhang   PetscInt           i, j, r, row, nnz_in_row, totalslices = a->totalslices;
688d4002b98SHong Zhang 
689d4002b98SHong Zhang #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
690d4002b98SHong Zhang   #pragma disjoint(*x, *y, *aval)
691d4002b98SHong Zhang #endif
692d4002b98SHong Zhang 
693d4002b98SHong Zhang   PetscFunctionBegin;
694b94d7dedSBarry Smith   if (A->symmetric == PETSC_BOOL3_TRUE) {
6959566063dSJacob Faibussowitsch     PetscCall(MatMultAdd_SeqSELL(A, xx, zz, yy));
6963ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6979fc32365SStefano Zampini   }
6989566063dSJacob Faibussowitsch   if (zz != yy) PetscCall(VecCopy(zz, yy));
6999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xx, &x));
7009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(yy, &y));
701d4002b98SHong Zhang   for (i = 0; i < a->totalslices; i++) { /* loop over slices */
7027285fed1SHong Zhang     if (i == totalslices - 1 && (A->rmap->n & 0x07)) {
7037285fed1SHong Zhang       for (r = 0; r < (A->rmap->n & 0x07); ++r) {
7047285fed1SHong Zhang         row        = 8 * i + r;
7057285fed1SHong Zhang         nnz_in_row = a->rlen[row];
7067285fed1SHong Zhang         for (j = 0; j < nnz_in_row; ++j) y[acolidx[8 * j + r]] += aval[8 * j + r] * x[row];
7077285fed1SHong Zhang       }
7087285fed1SHong Zhang       break;
7097285fed1SHong Zhang     }
7107285fed1SHong Zhang     for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j += 8) {
7117285fed1SHong Zhang       y[acolidx[j]] += aval[j] * x[8 * i];
7127285fed1SHong Zhang       y[acolidx[j + 1]] += aval[j + 1] * x[8 * i + 1];
7137285fed1SHong Zhang       y[acolidx[j + 2]] += aval[j + 2] * x[8 * i + 2];
7147285fed1SHong Zhang       y[acolidx[j + 3]] += aval[j + 3] * x[8 * i + 3];
7157285fed1SHong Zhang       y[acolidx[j + 4]] += aval[j + 4] * x[8 * i + 4];
7167285fed1SHong Zhang       y[acolidx[j + 5]] += aval[j + 5] * x[8 * i + 5];
7177285fed1SHong Zhang       y[acolidx[j + 6]] += aval[j + 6] * x[8 * i + 6];
7187285fed1SHong Zhang       y[acolidx[j + 7]] += aval[j + 7] * x[8 * i + 7];
719d4002b98SHong Zhang     }
720d4002b98SHong Zhang   }
7219566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.0 * a->sliidx[a->totalslices]));
7229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xx, &x));
7239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(yy, &y));
7243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
725d4002b98SHong Zhang }
726d4002b98SHong Zhang 
727d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTranspose_SeqSELL(Mat A, Vec xx, Vec yy)
728d71ae5a4SJacob Faibussowitsch {
729d4002b98SHong Zhang   PetscFunctionBegin;
730b94d7dedSBarry Smith   if (A->symmetric == PETSC_BOOL3_TRUE) {
7319566063dSJacob Faibussowitsch     PetscCall(MatMult_SeqSELL(A, xx, yy));
7329fc32365SStefano Zampini   } else {
7339566063dSJacob Faibussowitsch     PetscCall(VecSet(yy, 0.0));
7349566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd_SeqSELL(A, xx, yy, yy));
7359fc32365SStefano Zampini   }
7363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
737d4002b98SHong Zhang }
738d4002b98SHong Zhang 
739d4002b98SHong Zhang /*
740d4002b98SHong Zhang      Checks for missing diagonals
741d4002b98SHong Zhang */
742d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMissingDiagonal_SeqSELL(Mat A, PetscBool *missing, PetscInt *d)
743d71ae5a4SJacob Faibussowitsch {
744d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
745d4002b98SHong Zhang   PetscInt    *diag, i;
746d4002b98SHong Zhang 
747d4002b98SHong Zhang   PetscFunctionBegin;
748d4002b98SHong Zhang   *missing = PETSC_FALSE;
749d4002b98SHong Zhang   if (A->rmap->n > 0 && !(a->colidx)) {
750d4002b98SHong Zhang     *missing = PETSC_TRUE;
751d4002b98SHong Zhang     if (d) *d = 0;
7529566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Matrix has no entries therefore is missing diagonal\n"));
753d4002b98SHong Zhang   } else {
754d4002b98SHong Zhang     diag = a->diag;
755d4002b98SHong Zhang     for (i = 0; i < A->rmap->n; i++) {
756d4002b98SHong Zhang       if (diag[i] == -1) {
757d4002b98SHong Zhang         *missing = PETSC_TRUE;
758d4002b98SHong Zhang         if (d) *d = i;
7599566063dSJacob Faibussowitsch         PetscCall(PetscInfo(A, "Matrix is missing diagonal number %" PetscInt_FMT "\n", i));
760d4002b98SHong Zhang         break;
761d4002b98SHong Zhang       }
762d4002b98SHong Zhang     }
763d4002b98SHong Zhang   }
7643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
765d4002b98SHong Zhang }
766d4002b98SHong Zhang 
767d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMarkDiagonal_SeqSELL(Mat A)
768d71ae5a4SJacob Faibussowitsch {
769d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
770d4002b98SHong Zhang   PetscInt     i, j, m = A->rmap->n, shift;
771d4002b98SHong Zhang 
772d4002b98SHong Zhang   PetscFunctionBegin;
773d4002b98SHong Zhang   if (!a->diag) {
7749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &a->diag));
775d4002b98SHong Zhang     a->free_diag = PETSC_TRUE;
776d4002b98SHong Zhang   }
777d4002b98SHong Zhang   for (i = 0; i < m; i++) {                      /* loop over rows */
778d4002b98SHong Zhang     shift      = a->sliidx[i >> 3] + (i & 0x07); /* starting index of the row i */
779d4002b98SHong Zhang     a->diag[i] = -1;
780d4002b98SHong Zhang     for (j = 0; j < a->rlen[i]; j++) {
781d4002b98SHong Zhang       if (a->colidx[shift + j * 8] == i) {
782d4002b98SHong Zhang         a->diag[i] = shift + j * 8;
783d4002b98SHong Zhang         break;
784d4002b98SHong Zhang       }
785d4002b98SHong Zhang     }
786d4002b98SHong Zhang   }
7873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
788d4002b98SHong Zhang }
789d4002b98SHong Zhang 
790d4002b98SHong Zhang /*
791d4002b98SHong Zhang   Negative shift indicates do not generate an error if there is a zero diagonal, just invert it anyways
792d4002b98SHong Zhang */
793d71ae5a4SJacob Faibussowitsch PetscErrorCode MatInvertDiagonal_SeqSELL(Mat A, PetscScalar omega, PetscScalar fshift)
794d71ae5a4SJacob Faibussowitsch {
795d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
796d4002b98SHong Zhang   PetscInt     i, *diag, m = A->rmap->n;
797d4002b98SHong Zhang   MatScalar   *val = a->val;
798d4002b98SHong Zhang   PetscScalar *idiag, *mdiag;
799d4002b98SHong Zhang 
800d4002b98SHong Zhang   PetscFunctionBegin;
8013ba16761SJacob Faibussowitsch   if (a->idiagvalid) PetscFunctionReturn(PETSC_SUCCESS);
8029566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_SeqSELL(A));
803d4002b98SHong Zhang   diag = a->diag;
804d4002b98SHong Zhang   if (!a->idiag) {
8059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(m, &a->idiag, m, &a->mdiag, m, &a->ssor_work));
806d4002b98SHong Zhang     val = a->val;
807d4002b98SHong Zhang   }
808d4002b98SHong Zhang   mdiag = a->mdiag;
809d4002b98SHong Zhang   idiag = a->idiag;
810d4002b98SHong Zhang 
811d4002b98SHong Zhang   if (omega == 1.0 && PetscRealPart(fshift) <= 0.0) {
812d4002b98SHong Zhang     for (i = 0; i < m; i++) {
813d4002b98SHong Zhang       mdiag[i] = val[diag[i]];
814d4002b98SHong Zhang       if (!PetscAbsScalar(mdiag[i])) { /* zero diagonal */
8150fdf79fbSJacob Faibussowitsch         PetscCheck(PetscRealPart(fshift), PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Zero diagonal on row %" PetscInt_FMT, i);
8169566063dSJacob Faibussowitsch         PetscCall(PetscInfo(A, "Zero diagonal on row %" PetscInt_FMT "\n", i));
817d4002b98SHong Zhang         A->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
818d4002b98SHong Zhang         A->factorerror_zeropivot_value = 0.0;
819d4002b98SHong Zhang         A->factorerror_zeropivot_row   = i;
820d4002b98SHong Zhang       }
821d4002b98SHong Zhang       idiag[i] = 1.0 / val[diag[i]];
822d4002b98SHong Zhang     }
8239566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(m));
824d4002b98SHong Zhang   } else {
825d4002b98SHong Zhang     for (i = 0; i < m; i++) {
826d4002b98SHong Zhang       mdiag[i] = val[diag[i]];
827d4002b98SHong Zhang       idiag[i] = omega / (fshift + val[diag[i]]);
828d4002b98SHong Zhang     }
8299566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(2.0 * m));
830d4002b98SHong Zhang   }
831d4002b98SHong Zhang   a->idiagvalid = PETSC_TRUE;
8323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
833d4002b98SHong Zhang }
834d4002b98SHong Zhang 
835d71ae5a4SJacob Faibussowitsch PetscErrorCode MatZeroEntries_SeqSELL(Mat A)
836d71ae5a4SJacob Faibussowitsch {
837d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
838d4002b98SHong Zhang 
839d4002b98SHong Zhang   PetscFunctionBegin;
8409566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(a->val, a->sliidx[a->totalslices]));
8419566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
8423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
843d4002b98SHong Zhang }
844d4002b98SHong Zhang 
845d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_SeqSELL(Mat A)
846d71ae5a4SJacob Faibussowitsch {
847d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
848d4002b98SHong Zhang 
849d4002b98SHong Zhang   PetscFunctionBegin;
850d4002b98SHong Zhang #if defined(PETSC_USE_LOG)
8513ba16761SJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)A, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, A->rmap->n, A->cmap->n, a->nz));
852d4002b98SHong Zhang #endif
8539566063dSJacob Faibussowitsch   PetscCall(MatSeqXSELLFreeSELL(A, &a->val, &a->colidx));
8549566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->row));
8559566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->col));
8569566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->diag));
8579566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->rlen));
8589566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->sliidx));
8599566063dSJacob Faibussowitsch   PetscCall(PetscFree3(a->idiag, a->mdiag, a->ssor_work));
8609566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->solve_work));
8619566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&a->icol));
8629566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->saved_values));
8639566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a->getrowcols, a->getrowvals));
864d4002b98SHong Zhang 
8659566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
866d4002b98SHong Zhang 
8679566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
8689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatStoreValues_C", NULL));
8699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatRetrieveValues_C", NULL));
8709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLSetPreallocation_C", NULL));
8712e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLGetArray_C", NULL));
8722e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatSeqSELLRestoreArray_C", NULL));
8732e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_seqsell_seqaij_C", NULL));
8743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
875d4002b98SHong Zhang }
876d4002b98SHong Zhang 
877d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetOption_SeqSELL(Mat A, MatOption op, PetscBool flg)
878d71ae5a4SJacob Faibussowitsch {
879d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
880d4002b98SHong Zhang 
881d4002b98SHong Zhang   PetscFunctionBegin;
882d4002b98SHong Zhang   switch (op) {
883d71ae5a4SJacob Faibussowitsch   case MAT_ROW_ORIENTED:
884d71ae5a4SJacob Faibussowitsch     a->roworiented = flg;
885d71ae5a4SJacob Faibussowitsch     break;
886d71ae5a4SJacob Faibussowitsch   case MAT_KEEP_NONZERO_PATTERN:
887d71ae5a4SJacob Faibussowitsch     a->keepnonzeropattern = flg;
888d71ae5a4SJacob Faibussowitsch     break;
889d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_LOCATIONS:
890d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? 0 : 1);
891d71ae5a4SJacob Faibussowitsch     break;
892d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_LOCATION_ERR:
893d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? -1 : 0);
894d71ae5a4SJacob Faibussowitsch     break;
895d71ae5a4SJacob Faibussowitsch   case MAT_NEW_NONZERO_ALLOCATION_ERR:
896d71ae5a4SJacob Faibussowitsch     a->nonew = (flg ? -2 : 0);
897d71ae5a4SJacob Faibussowitsch     break;
898d71ae5a4SJacob Faibussowitsch   case MAT_UNUSED_NONZERO_LOCATION_ERR:
899d71ae5a4SJacob Faibussowitsch     a->nounused = (flg ? -1 : 0);
900d71ae5a4SJacob Faibussowitsch     break;
9018c78258cSHong Zhang   case MAT_FORCE_DIAGONAL_ENTRIES:
902d4002b98SHong Zhang   case MAT_IGNORE_OFF_PROC_ENTRIES:
903d4002b98SHong Zhang   case MAT_USE_HASH_TABLE:
904d71ae5a4SJacob Faibussowitsch   case MAT_SORTED_FULL:
905d71ae5a4SJacob Faibussowitsch     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
906d71ae5a4SJacob Faibussowitsch     break;
907d4002b98SHong Zhang   case MAT_SPD:
908d4002b98SHong Zhang   case MAT_SYMMETRIC:
909d4002b98SHong Zhang   case MAT_STRUCTURALLY_SYMMETRIC:
910d4002b98SHong Zhang   case MAT_HERMITIAN:
911d4002b98SHong Zhang   case MAT_SYMMETRY_ETERNAL:
912b94d7dedSBarry Smith   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
913b94d7dedSBarry Smith   case MAT_SPD_ETERNAL:
914d4002b98SHong Zhang     /* These options are handled directly by MatSetOption() */
915d4002b98SHong Zhang     break;
916d71ae5a4SJacob Faibussowitsch   default:
917d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
918d4002b98SHong Zhang   }
9193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
920d4002b98SHong Zhang }
921d4002b98SHong Zhang 
922d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_SeqSELL(Mat A, Vec v)
923d71ae5a4SJacob Faibussowitsch {
924d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
925d4002b98SHong Zhang   PetscInt     i, j, n, shift;
926d4002b98SHong Zhang   PetscScalar *x, zero = 0.0;
927d4002b98SHong Zhang 
928d4002b98SHong Zhang   PetscFunctionBegin;
9299566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
93008401ef6SPierre Jolivet   PetscCheck(n == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming matrix and vector");
931d4002b98SHong Zhang 
932d4002b98SHong Zhang   if (A->factortype == MAT_FACTOR_ILU || A->factortype == MAT_FACTOR_LU) {
933d4002b98SHong Zhang     PetscInt *diag = a->diag;
9349566063dSJacob Faibussowitsch     PetscCall(VecGetArray(v, &x));
935d4002b98SHong Zhang     for (i = 0; i < n; i++) x[i] = 1.0 / a->val[diag[i]];
9369566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(v, &x));
9373ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
938d4002b98SHong Zhang   }
939d4002b98SHong Zhang 
9409566063dSJacob Faibussowitsch   PetscCall(VecSet(v, zero));
9419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &x));
942d4002b98SHong Zhang   for (i = 0; i < n; i++) {                 /* loop over rows */
943d4002b98SHong Zhang     shift = a->sliidx[i >> 3] + (i & 0x07); /* starting index of the row i */
944d4002b98SHong Zhang     x[i]  = 0;
945d4002b98SHong Zhang     for (j = 0; j < a->rlen[i]; j++) {
946d4002b98SHong Zhang       if (a->colidx[shift + j * 8] == i) {
947d4002b98SHong Zhang         x[i] = a->val[shift + j * 8];
948d4002b98SHong Zhang         break;
949d4002b98SHong Zhang       }
950d4002b98SHong Zhang     }
951d4002b98SHong Zhang   }
9529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &x));
9533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
954d4002b98SHong Zhang }
955d4002b98SHong Zhang 
956d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalScale_SeqSELL(Mat A, Vec ll, Vec rr)
957d71ae5a4SJacob Faibussowitsch {
958d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
959d4002b98SHong Zhang   const PetscScalar *l, *r;
960d4002b98SHong Zhang   PetscInt           i, j, m, n, row;
961d4002b98SHong Zhang 
962d4002b98SHong Zhang   PetscFunctionBegin;
963d4002b98SHong Zhang   if (ll) {
964d4002b98SHong Zhang     /* The local size is used so that VecMPI can be passed to this routine
965d4002b98SHong Zhang        by MatDiagonalScale_MPISELL */
9669566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(ll, &m));
96708401ef6SPierre Jolivet     PetscCheck(m == A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Left scaling vector wrong length");
9689566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(ll, &l));
969d4002b98SHong Zhang     for (i = 0; i < a->totalslices; i++) {                  /* loop over slices */
970dab86139SHong Zhang       if (i == a->totalslices - 1 && (A->rmap->n & 0x07)) { /* if last slice has padding rows */
971dab86139SHong Zhang         for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = ((row + 1) & 0x07)) {
972dab86139SHong Zhang           if (row < (A->rmap->n & 0x07)) a->val[j] *= l[8 * i + row];
973dab86139SHong Zhang         }
974dab86139SHong Zhang       } else {
975ad540459SPierre Jolivet         for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = ((row + 1) & 0x07)) a->val[j] *= l[8 * i + row];
976d4002b98SHong Zhang       }
977dab86139SHong Zhang     }
9789566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(ll, &l));
9799566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(a->nz));
980d4002b98SHong Zhang   }
981d4002b98SHong Zhang   if (rr) {
9829566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rr, &n));
98308401ef6SPierre Jolivet     PetscCheck(n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Right scaling vector wrong length");
9849566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(rr, &r));
985d4002b98SHong Zhang     for (i = 0; i < a->totalslices; i++) {                  /* loop over slices */
986dab86139SHong Zhang       if (i == a->totalslices - 1 && (A->rmap->n & 0x07)) { /* if last slice has padding rows */
987dab86139SHong Zhang         for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = ((row + 1) & 0x07)) {
988dab86139SHong Zhang           if (row < (A->rmap->n & 0x07)) a->val[j] *= r[a->colidx[j]];
989dab86139SHong Zhang         }
990dab86139SHong Zhang       } else {
991ad540459SPierre Jolivet         for (j = a->sliidx[i]; j < a->sliidx[i + 1]; j++) a->val[j] *= r[a->colidx[j]];
992d4002b98SHong Zhang       }
993dab86139SHong Zhang     }
9949566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(rr, &r));
9959566063dSJacob Faibussowitsch     PetscCall(PetscLogFlops(a->nz));
996d4002b98SHong Zhang   }
9979566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
9983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
999d4002b98SHong Zhang }
1000d4002b98SHong Zhang 
1001d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetValues_SeqSELL(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], PetscScalar v[])
1002d71ae5a4SJacob Faibussowitsch {
1003d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1004d4002b98SHong Zhang   PetscInt    *cp, i, k, low, high, t, row, col, l;
1005d4002b98SHong Zhang   PetscInt     shift;
1006d4002b98SHong Zhang   MatScalar   *vp;
1007d4002b98SHong Zhang 
1008d4002b98SHong Zhang   PetscFunctionBegin;
100968aafef3SStefano Zampini   for (k = 0; k < m; k++) { /* loop over requested rows */
1010d4002b98SHong Zhang     row = im[k];
1011d4002b98SHong Zhang     if (row < 0) continue;
10126bdcaf15SBarry 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);
1013d4002b98SHong Zhang     shift = a->sliidx[row >> 3] + (row & 0x07); /* starting index of the row */
1014d4002b98SHong Zhang     cp    = a->colidx + shift;                  /* pointer to the row */
1015d4002b98SHong Zhang     vp    = a->val + shift;                     /* pointer to the row */
101668aafef3SStefano Zampini     for (l = 0; l < n; l++) {                   /* loop over requested columns */
1017d4002b98SHong Zhang       col = in[l];
1018d4002b98SHong Zhang       if (col < 0) continue;
10196bdcaf15SBarry 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);
10209371c9d4SSatish Balay       high = a->rlen[row];
10219371c9d4SSatish Balay       low  = 0; /* assume unsorted */
1022d4002b98SHong Zhang       while (high - low > 5) {
1023d4002b98SHong Zhang         t = (low + high) / 2;
1024d4002b98SHong Zhang         if (*(cp + t * 8) > col) high = t;
1025d4002b98SHong Zhang         else low = t;
1026d4002b98SHong Zhang       }
1027d4002b98SHong Zhang       for (i = low; i < high; i++) {
1028d4002b98SHong Zhang         if (*(cp + 8 * i) > col) break;
1029d4002b98SHong Zhang         if (*(cp + 8 * i) == col) {
1030d4002b98SHong Zhang           *v++ = *(vp + 8 * i);
1031d4002b98SHong Zhang           goto finished;
1032d4002b98SHong Zhang         }
1033d4002b98SHong Zhang       }
1034d4002b98SHong Zhang       *v++ = 0.0;
1035d4002b98SHong Zhang     finished:;
1036d4002b98SHong Zhang     }
1037d4002b98SHong Zhang   }
10383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1039d4002b98SHong Zhang }
1040d4002b98SHong Zhang 
1041d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL_ASCII(Mat A, PetscViewer viewer)
1042d71ae5a4SJacob Faibussowitsch {
1043d4002b98SHong Zhang   Mat_SeqSELL      *a = (Mat_SeqSELL *)A->data;
1044d4002b98SHong Zhang   PetscInt          i, j, m = A->rmap->n, shift;
1045d4002b98SHong Zhang   const char       *name;
1046d4002b98SHong Zhang   PetscViewerFormat format;
1047d4002b98SHong Zhang 
1048d4002b98SHong Zhang   PetscFunctionBegin;
10499566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
1050d4002b98SHong Zhang   if (format == PETSC_VIEWER_ASCII_MATLAB) {
1051d4002b98SHong Zhang     PetscInt nofinalvalue = 0;
1052d4002b98SHong Zhang     /*
1053d4002b98SHong Zhang     if (m && ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap->n-1))) {
1054d4002b98SHong Zhang       nofinalvalue = 1;
1055d4002b98SHong Zhang     }
1056d4002b98SHong Zhang     */
10579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
10589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Size = %" PetscInt_FMT " %" PetscInt_FMT " \n", m, A->cmap->n));
10599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%% Nonzeros = %" PetscInt_FMT " \n", a->nz));
1060d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
10619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",4);\n", a->nz + nofinalvalue));
1062d4002b98SHong Zhang #else
10639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = zeros(%" PetscInt_FMT ",3);\n", a->nz + nofinalvalue));
1064d4002b98SHong Zhang #endif
10659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "zzz = [\n"));
1066d4002b98SHong Zhang 
1067d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1068d4002b98SHong Zhang       shift = a->sliidx[i >> 3] + (i & 0x07);
1069d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1070d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
10719566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e %18.16e\n", i + 1, a->colidx[shift + 8 * j] + 1, (double)PetscRealPart(a->val[shift + 8 * j]), (double)PetscImaginaryPart(a->val[shift + 8 * j])));
1072d4002b98SHong Zhang #else
10739566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n", i + 1, a->colidx[shift + 8 * j] + 1, (double)a->val[shift + 8 * j]));
1074d4002b98SHong Zhang #endif
1075d4002b98SHong Zhang       }
1076d4002b98SHong Zhang     }
1077d4002b98SHong Zhang     /*
1078d4002b98SHong Zhang     if (nofinalvalue) {
1079d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
10809566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e %18.16e\n",m,A->cmap->n,0.,0.));
1081d4002b98SHong Zhang #else
10829566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT "  %18.16e\n",m,A->cmap->n,0.0));
1083d4002b98SHong Zhang #endif
1084d4002b98SHong Zhang     }
1085d4002b98SHong Zhang     */
10869566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)A, &name));
10879566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "];\n %s = spconvert(zzz);\n", name));
10889566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1089d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) {
10903ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1091d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
10929566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1093d4002b98SHong Zhang     for (i = 0; i < m; i++) {
10949566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
1095d4002b98SHong Zhang       shift = a->sliidx[i >> 3] + (i & 0x07);
1096d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1097d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1098d4002b98SHong Zhang         if (PetscImaginaryPart(a->val[shift + 8 * j]) > 0.0 && PetscRealPart(a->val[shift + 8 * j]) != 0.0) {
10999566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[shift + 8 * j], (double)PetscRealPart(a->val[shift + 8 * j]), (double)PetscImaginaryPart(a->val[shift + 8 * j])));
1100d4002b98SHong Zhang         } else if (PetscImaginaryPart(a->val[shift + 8 * j]) < 0.0 && PetscRealPart(a->val[shift + 8 * j]) != 0.0) {
11019566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[shift + 8 * j], (double)PetscRealPart(a->val[shift + 8 * j]), (double)-PetscImaginaryPart(a->val[shift + 8 * j])));
1102d4002b98SHong Zhang         } else if (PetscRealPart(a->val[shift + 8 * j]) != 0.0) {
11039566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + 8 * j], (double)PetscRealPart(a->val[shift + 8 * j])));
1104d4002b98SHong Zhang         }
1105d4002b98SHong Zhang #else
11069566063dSJacob Faibussowitsch         if (a->val[shift + 8 * j] != 0.0) PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + 8 * j], (double)a->val[shift + 8 * j]));
1107d4002b98SHong Zhang #endif
1108d4002b98SHong Zhang       }
11099566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1110d4002b98SHong Zhang     }
11119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1112d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_DENSE) {
1113d4002b98SHong Zhang     PetscInt    cnt = 0, jcnt;
1114d4002b98SHong Zhang     PetscScalar value;
1115d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1116d4002b98SHong Zhang     PetscBool realonly = PETSC_TRUE;
1117d4002b98SHong Zhang     for (i = 0; i < a->sliidx[a->totalslices]; i++) {
1118d4002b98SHong Zhang       if (PetscImaginaryPart(a->val[i]) != 0.0) {
1119d4002b98SHong Zhang         realonly = PETSC_FALSE;
1120d4002b98SHong Zhang         break;
1121d4002b98SHong Zhang       }
1122d4002b98SHong Zhang     }
1123d4002b98SHong Zhang #endif
1124d4002b98SHong Zhang 
11259566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1126d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1127d4002b98SHong Zhang       jcnt  = 0;
1128d4002b98SHong Zhang       shift = a->sliidx[i >> 3] + (i & 0x07);
1129d4002b98SHong Zhang       for (j = 0; j < A->cmap->n; j++) {
1130d4002b98SHong Zhang         if (jcnt < a->rlen[i] && j == a->colidx[shift + 8 * j]) {
1131d4002b98SHong Zhang           value = a->val[cnt++];
1132d4002b98SHong Zhang           jcnt++;
1133d4002b98SHong Zhang         } else {
1134d4002b98SHong Zhang           value = 0.0;
1135d4002b98SHong Zhang         }
1136d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1137d4002b98SHong Zhang         if (realonly) {
11389566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)PetscRealPart(value)));
1139d4002b98SHong Zhang         } else {
11409566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e+%7.5e i ", (double)PetscRealPart(value), (double)PetscImaginaryPart(value)));
1141d4002b98SHong Zhang         }
1142d4002b98SHong Zhang #else
11439566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " %7.5e ", (double)value));
1144d4002b98SHong Zhang #endif
1145d4002b98SHong Zhang       }
11469566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1147d4002b98SHong Zhang     }
11489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1149d4002b98SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_MATRIXMARKET) {
1150d4002b98SHong Zhang     PetscInt fshift = 1;
11519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1152d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
11539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate complex general\n"));
1154d4002b98SHong Zhang #else
11559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%%%%MatrixMarket matrix coordinate real general\n"));
1156d4002b98SHong Zhang #endif
11579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", m, A->cmap->n, a->nz));
1158d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1159d4002b98SHong Zhang       shift = a->sliidx[i >> 3] + (i & 0x07);
1160d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1161d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
11629566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %g %g\n", i + fshift, a->colidx[shift + 8 * j] + fshift, (double)PetscRealPart(a->val[shift + 8 * j]), (double)PetscImaginaryPart(a->val[shift + 8 * j])));
1163d4002b98SHong Zhang #else
11649566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %g\n", i + fshift, a->colidx[shift + 8 * j] + fshift, (double)a->val[shift + 8 * j]));
1165d4002b98SHong Zhang #endif
1166d4002b98SHong Zhang       }
1167d4002b98SHong Zhang     }
11689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
116968aafef3SStefano Zampini   } else if (format == PETSC_VIEWER_NATIVE) {
117068aafef3SStefano Zampini     for (i = 0; i < a->totalslices; i++) { /* loop over slices */
117168aafef3SStefano Zampini       PetscInt row;
11729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "slice %" PetscInt_FMT ": %" PetscInt_FMT " %" PetscInt_FMT "\n", i, a->sliidx[i], a->sliidx[i + 1]));
117368aafef3SStefano Zampini       for (j = a->sliidx[i], row = 0; j < a->sliidx[i + 1]; j++, row = ((row + 1) & 0x07)) {
117468aafef3SStefano Zampini #if defined(PETSC_USE_COMPLEX)
117568aafef3SStefano Zampini         if (PetscImaginaryPart(a->val[j]) > 0.0) {
11769566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g + %g i\n", 8 * i + row, a->colidx[j], (double)PetscRealPart(a->val[j]), (double)PetscImaginaryPart(a->val[j])));
117768aafef3SStefano Zampini         } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
11789566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g - %g i\n", 8 * i + row, a->colidx[j], (double)PetscRealPart(a->val[j]), -(double)PetscImaginaryPart(a->val[j])));
117968aafef3SStefano Zampini         } else {
11809566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g\n", 8 * i + row, a->colidx[j], (double)PetscRealPart(a->val[j])));
118168aafef3SStefano Zampini         }
118268aafef3SStefano Zampini #else
11839566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " %" PetscInt_FMT " %g\n", 8 * i + row, a->colidx[j], (double)a->val[j]));
118468aafef3SStefano Zampini #endif
118568aafef3SStefano Zampini       }
118668aafef3SStefano Zampini     }
1187d4002b98SHong Zhang   } else {
11889566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1189d4002b98SHong Zhang     if (A->factortype) {
1190d4002b98SHong Zhang       for (i = 0; i < m; i++) {
1191d4002b98SHong Zhang         shift = a->sliidx[i >> 3] + (i & 0x07);
11929566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
1193d4002b98SHong Zhang         /* L part */
1194d4002b98SHong Zhang         for (j = shift; j < a->diag[i]; j += 8) {
1195d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1196d4002b98SHong Zhang           if (PetscImaginaryPart(a->val[shift + 8 * j]) > 0.0) {
11979566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)PetscImaginaryPart(a->val[j])));
1198d4002b98SHong Zhang           } else if (PetscImaginaryPart(a->val[shift + 8 * j]) < 0.0) {
11999566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)(-PetscImaginaryPart(a->val[j]))));
1200d4002b98SHong Zhang           } else {
12019566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)PetscRealPart(a->val[j])));
1202d4002b98SHong Zhang           }
1203d4002b98SHong Zhang #else
12049566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)a->val[j]));
1205d4002b98SHong Zhang #endif
1206d4002b98SHong Zhang         }
1207d4002b98SHong Zhang         /* diagonal */
1208d4002b98SHong Zhang         j = a->diag[i];
1209d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1210d4002b98SHong Zhang         if (PetscImaginaryPart(a->val[j]) > 0.0) {
12119566063dSJacob 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])));
1212d4002b98SHong Zhang         } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
12139566063dSJacob 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]))));
1214d4002b98SHong Zhang         } else {
12159566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)PetscRealPart(1.0 / a->val[j])));
1216d4002b98SHong Zhang         }
1217d4002b98SHong Zhang #else
12189566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)(1.0 / a->val[j])));
1219d4002b98SHong Zhang #endif
1220d4002b98SHong Zhang 
1221d4002b98SHong Zhang         /* U part */
1222d4002b98SHong Zhang         for (j = a->diag[i] + 1; j < shift + 8 * a->rlen[i]; j += 8) {
1223d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1224d4002b98SHong Zhang           if (PetscImaginaryPart(a->val[j]) > 0.0) {
12259566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)PetscImaginaryPart(a->val[j])));
1226d4002b98SHong Zhang           } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
12279566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[j], (double)PetscRealPart(a->val[j]), (double)(-PetscImaginaryPart(a->val[j]))));
1228d4002b98SHong Zhang           } else {
12299566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)PetscRealPart(a->val[j])));
1230d4002b98SHong Zhang           }
1231d4002b98SHong Zhang #else
12329566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[j], (double)a->val[j]));
1233d4002b98SHong Zhang #endif
1234d4002b98SHong Zhang         }
12359566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1236d4002b98SHong Zhang       }
1237d4002b98SHong Zhang     } else {
1238d4002b98SHong Zhang       for (i = 0; i < m; i++) {
1239d4002b98SHong Zhang         shift = a->sliidx[i >> 3] + (i & 0x07);
12409566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "row %" PetscInt_FMT ":", i));
1241d4002b98SHong Zhang         for (j = 0; j < a->rlen[i]; j++) {
1242d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
1243d4002b98SHong Zhang           if (PetscImaginaryPart(a->val[j]) > 0.0) {
12449566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g + %g i)", a->colidx[shift + 8 * j], (double)PetscRealPart(a->val[shift + 8 * j]), (double)PetscImaginaryPart(a->val[shift + 8 * j])));
1245d4002b98SHong Zhang           } else if (PetscImaginaryPart(a->val[j]) < 0.0) {
12469566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g - %g i)", a->colidx[shift + 8 * j], (double)PetscRealPart(a->val[shift + 8 * j]), (double)-PetscImaginaryPart(a->val[shift + 8 * j])));
1247d4002b98SHong Zhang           } else {
12489566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + 8 * j], (double)PetscRealPart(a->val[shift + 8 * j])));
1249d4002b98SHong Zhang           }
1250d4002b98SHong Zhang #else
12519566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, " (%" PetscInt_FMT ", %g) ", a->colidx[shift + 8 * j], (double)a->val[shift + 8 * j]));
1252d4002b98SHong Zhang #endif
1253d4002b98SHong Zhang         }
12549566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1255d4002b98SHong Zhang       }
1256d4002b98SHong Zhang     }
12579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1258d4002b98SHong Zhang   }
12599566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
12603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1261d4002b98SHong Zhang }
1262d4002b98SHong Zhang 
1263d4002b98SHong Zhang #include <petscdraw.h>
1264d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL_Draw_Zoom(PetscDraw draw, void *Aa)
1265d71ae5a4SJacob Faibussowitsch {
1266d4002b98SHong Zhang   Mat               A = (Mat)Aa;
1267d4002b98SHong Zhang   Mat_SeqSELL      *a = (Mat_SeqSELL *)A->data;
1268d4002b98SHong Zhang   PetscInt          i, j, m = A->rmap->n, shift;
1269d4002b98SHong Zhang   int               color;
1270d4002b98SHong Zhang   PetscReal         xl, yl, xr, yr, x_l, x_r, y_l, y_r;
1271d4002b98SHong Zhang   PetscViewer       viewer;
1272d4002b98SHong Zhang   PetscViewerFormat format;
1273d4002b98SHong Zhang 
1274d4002b98SHong Zhang   PetscFunctionBegin;
12759566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, "Zoomviewer", (PetscObject *)&viewer));
12769566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
12779566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
1278d4002b98SHong Zhang 
1279d4002b98SHong Zhang   /* loop over matrix elements drawing boxes */
1280d4002b98SHong Zhang 
1281d4002b98SHong Zhang   if (format != PETSC_VIEWER_DRAW_CONTOUR) {
1282d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
1283d4002b98SHong Zhang     /* Blue for negative, Cyan for zero and  Red for positive */
1284d4002b98SHong Zhang     color = PETSC_DRAW_BLUE;
1285d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1286d4002b98SHong Zhang       shift = a->sliidx[i >> 3] + (i & 0x07); /* starting index of the row i */
12879371c9d4SSatish Balay       y_l   = m - i - 1.0;
12889371c9d4SSatish Balay       y_r   = y_l + 1.0;
1289d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
12909371c9d4SSatish Balay         x_l = a->colidx[shift + j * 8];
12919371c9d4SSatish Balay         x_r = x_l + 1.0;
1292d4002b98SHong Zhang         if (PetscRealPart(a->val[shift + 8 * j]) >= 0.) continue;
12939566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1294d4002b98SHong Zhang       }
1295d4002b98SHong Zhang     }
1296d4002b98SHong Zhang     color = PETSC_DRAW_CYAN;
1297d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1298d4002b98SHong Zhang       shift = a->sliidx[i >> 3] + (i & 0x07);
12999371c9d4SSatish Balay       y_l   = m - i - 1.0;
13009371c9d4SSatish Balay       y_r   = y_l + 1.0;
1301d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
13029371c9d4SSatish Balay         x_l = a->colidx[shift + j * 8];
13039371c9d4SSatish Balay         x_r = x_l + 1.0;
1304d4002b98SHong Zhang         if (a->val[shift + 8 * j] != 0.) continue;
13059566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1306d4002b98SHong Zhang       }
1307d4002b98SHong Zhang     }
1308d4002b98SHong Zhang     color = PETSC_DRAW_RED;
1309d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1310d4002b98SHong Zhang       shift = a->sliidx[i >> 3] + (i & 0x07);
13119371c9d4SSatish Balay       y_l   = m - i - 1.0;
13129371c9d4SSatish Balay       y_r   = y_l + 1.0;
1313d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
13149371c9d4SSatish Balay         x_l = a->colidx[shift + j * 8];
13159371c9d4SSatish Balay         x_r = x_l + 1.0;
1316d4002b98SHong Zhang         if (PetscRealPart(a->val[shift + 8 * j]) <= 0.) continue;
13179566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1318d4002b98SHong Zhang       }
1319d4002b98SHong Zhang     }
1320d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1321d4002b98SHong Zhang   } else {
1322d4002b98SHong Zhang     /* use contour shading to indicate magnitude of values */
1323d4002b98SHong Zhang     /* first determine max of all nonzero values */
1324d4002b98SHong Zhang     PetscReal minv = 0.0, maxv = 0.0;
1325d4002b98SHong Zhang     PetscInt  count = 0;
1326d4002b98SHong Zhang     PetscDraw popup;
1327d4002b98SHong Zhang     for (i = 0; i < a->sliidx[a->totalslices]; i++) {
1328d4002b98SHong Zhang       if (PetscAbsScalar(a->val[i]) > maxv) maxv = PetscAbsScalar(a->val[i]);
1329d4002b98SHong Zhang     }
1330d4002b98SHong Zhang     if (minv >= maxv) maxv = minv + PETSC_SMALL;
13319566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetPopup(draw, &popup));
13329566063dSJacob Faibussowitsch     PetscCall(PetscDrawScalePopup(popup, minv, maxv));
1333d4002b98SHong Zhang 
1334d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
1335d4002b98SHong Zhang     for (i = 0; i < m; i++) {
1336d4002b98SHong Zhang       shift = a->sliidx[i >> 3] + (i & 0x07);
1337d4002b98SHong Zhang       y_l   = m - i - 1.0;
1338d4002b98SHong Zhang       y_r   = y_l + 1.0;
1339d4002b98SHong Zhang       for (j = 0; j < a->rlen[i]; j++) {
1340d4002b98SHong Zhang         x_l   = a->colidx[shift + j * 8];
1341d4002b98SHong Zhang         x_r   = x_l + 1.0;
1342d4002b98SHong Zhang         color = PetscDrawRealToColor(PetscAbsScalar(a->val[count]), minv, maxv);
13439566063dSJacob Faibussowitsch         PetscCall(PetscDrawRectangle(draw, x_l, y_l, x_r, y_r, color, color, color, color));
1344d4002b98SHong Zhang         count++;
1345d4002b98SHong Zhang       }
1346d4002b98SHong Zhang     }
1347d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
1348d4002b98SHong Zhang   }
13493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1350d4002b98SHong Zhang }
1351d4002b98SHong Zhang 
1352d4002b98SHong Zhang #include <petscdraw.h>
1353d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL_Draw(Mat A, PetscViewer viewer)
1354d71ae5a4SJacob Faibussowitsch {
1355d4002b98SHong Zhang   PetscDraw draw;
1356d4002b98SHong Zhang   PetscReal xr, yr, xl, yl, h, w;
1357d4002b98SHong Zhang   PetscBool isnull;
1358d4002b98SHong Zhang 
1359d4002b98SHong Zhang   PetscFunctionBegin;
13609566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
13619566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
13623ba16761SJacob Faibussowitsch   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
1363d4002b98SHong Zhang 
13649371c9d4SSatish Balay   xr = A->cmap->n;
13659371c9d4SSatish Balay   yr = A->rmap->n;
13669371c9d4SSatish Balay   h  = yr / 10.0;
13679371c9d4SSatish Balay   w  = xr / 10.0;
13689371c9d4SSatish Balay   xr += w;
13699371c9d4SSatish Balay   yr += h;
13709371c9d4SSatish Balay   xl = -w;
13719371c9d4SSatish Balay   yl = -h;
13729566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetCoordinates(draw, xl, yl, xr, yr));
13739566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", (PetscObject)viewer));
13749566063dSJacob Faibussowitsch   PetscCall(PetscDrawZoom(draw, MatView_SeqSELL_Draw_Zoom, A));
13759566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, "Zoomviewer", NULL));
13769566063dSJacob Faibussowitsch   PetscCall(PetscDrawSave(draw));
13773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1378d4002b98SHong Zhang }
1379d4002b98SHong Zhang 
1380d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_SeqSELL(Mat A, PetscViewer viewer)
1381d71ae5a4SJacob Faibussowitsch {
1382d4002b98SHong Zhang   PetscBool iascii, isbinary, isdraw;
1383d4002b98SHong Zhang 
1384d4002b98SHong Zhang   PetscFunctionBegin;
13859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
13869566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
13879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1388d4002b98SHong Zhang   if (iascii) {
13899566063dSJacob Faibussowitsch     PetscCall(MatView_SeqSELL_ASCII(A, viewer));
1390d4002b98SHong Zhang   } else if (isbinary) {
13919566063dSJacob Faibussowitsch     /* PetscCall(MatView_SeqSELL_Binary(A,viewer)); */
13921baa6e33SBarry Smith   } else if (isdraw) PetscCall(MatView_SeqSELL_Draw(A, viewer));
13933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1394d4002b98SHong Zhang }
1395d4002b98SHong Zhang 
1396d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_SeqSELL(Mat A, MatAssemblyType mode)
1397d71ae5a4SJacob Faibussowitsch {
1398d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1399d4002b98SHong Zhang   PetscInt     i, shift, row_in_slice, row, nrow, *cp, lastcol, j, k;
1400d4002b98SHong Zhang   MatScalar   *vp;
1401d4002b98SHong Zhang 
1402d4002b98SHong Zhang   PetscFunctionBegin;
14033ba16761SJacob Faibussowitsch   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
1404d4002b98SHong Zhang   /* To do: compress out the unused elements */
14059566063dSJacob Faibussowitsch   PetscCall(MatMarkDiagonal_SeqSELL(A));
14069566063dSJacob 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));
14079566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Number of mallocs during MatSetValues() is %" PetscInt_FMT "\n", a->reallocs));
14089566063dSJacob Faibussowitsch   PetscCall(PetscInfo(A, "Maximum nonzeros in any row is %" PetscInt_FMT "\n", a->rlenmax));
1409d4002b98SHong 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 */
1410d4002b98SHong Zhang   for (i = 0; i < a->totalslices; ++i) {
1411d4002b98SHong Zhang     shift = a->sliidx[i];                                      /* starting index of the slice */
1412d4002b98SHong Zhang     cp    = a->colidx + shift;                                 /* pointer to the column indices of the slice */
1413d4002b98SHong Zhang     vp    = a->val + shift;                                    /* pointer to the nonzero values of the slice */
1414d4002b98SHong Zhang     for (row_in_slice = 0; row_in_slice < 8; ++row_in_slice) { /* loop over rows in the slice */
1415d4002b98SHong Zhang       row  = 8 * i + row_in_slice;
1416d4002b98SHong Zhang       nrow = a->rlen[row]; /* number of nonzeros in row */
1417d4002b98SHong Zhang       /*
1418d4002b98SHong Zhang         Search for the nearest nonzero. Normally setting the index to zero may cause extra communication.
1419d4002b98SHong Zhang         But if the entire slice are empty, it is fine to use 0 since the index will not be loaded.
1420d4002b98SHong Zhang       */
1421d4002b98SHong Zhang       lastcol = 0;
1422d4002b98SHong Zhang       if (nrow > 0) {                                /* nonempty row */
1423d4002b98SHong Zhang         lastcol = cp[8 * (nrow - 1) + row_in_slice]; /* use the index from the last nonzero at current row */
1424*aaa8cc7dSPierre Jolivet       } else if (!row_in_slice) {                    /* first row of the correct slice is empty */
1425d4002b98SHong Zhang         for (j = 1; j < 8; j++) {
1426d4002b98SHong Zhang           if (a->rlen[8 * i + j]) {
1427d4002b98SHong Zhang             lastcol = cp[j];
1428d4002b98SHong Zhang             break;
1429d4002b98SHong Zhang           }
1430d4002b98SHong Zhang         }
1431d4002b98SHong Zhang       } else {
1432d4002b98SHong Zhang         if (a->sliidx[i + 1] != shift) lastcol = cp[row_in_slice - 1]; /* use the index from the previous row */
1433d4002b98SHong Zhang       }
1434d4002b98SHong Zhang 
1435d4002b98SHong Zhang       for (k = nrow; k < (a->sliidx[i + 1] - shift) / 8; ++k) {
1436d4002b98SHong Zhang         cp[8 * k + row_in_slice] = lastcol;
1437d4002b98SHong Zhang         vp[8 * k + row_in_slice] = (MatScalar)0;
1438d4002b98SHong Zhang       }
1439d4002b98SHong Zhang     }
1440d4002b98SHong Zhang   }
1441d4002b98SHong Zhang 
1442d4002b98SHong Zhang   A->info.mallocs += a->reallocs;
1443d4002b98SHong Zhang   a->reallocs = 0;
1444d4002b98SHong Zhang 
14459566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
14463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1447d4002b98SHong Zhang }
1448d4002b98SHong Zhang 
1449d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetInfo_SeqSELL(Mat A, MatInfoType flag, MatInfo *info)
1450d71ae5a4SJacob Faibussowitsch {
1451d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1452d4002b98SHong Zhang 
1453d4002b98SHong Zhang   PetscFunctionBegin;
1454d4002b98SHong Zhang   info->block_size   = 1.0;
14553966268fSBarry Smith   info->nz_allocated = a->maxallocmat;
14563966268fSBarry Smith   info->nz_used      = a->sliidx[a->totalslices]; /* include padding zeros */
14573966268fSBarry Smith   info->nz_unneeded  = (a->maxallocmat - a->sliidx[a->totalslices]);
14583966268fSBarry Smith   info->assemblies   = A->num_ass;
14593966268fSBarry Smith   info->mallocs      = A->info.mallocs;
14604dfa11a4SJacob Faibussowitsch   info->memory       = 0; /* REVIEW ME */
1461d4002b98SHong Zhang   if (A->factortype) {
1462d4002b98SHong Zhang     info->fill_ratio_given  = A->info.fill_ratio_given;
1463d4002b98SHong Zhang     info->fill_ratio_needed = A->info.fill_ratio_needed;
1464d4002b98SHong Zhang     info->factor_mallocs    = A->info.factor_mallocs;
1465d4002b98SHong Zhang   } else {
1466d4002b98SHong Zhang     info->fill_ratio_given  = 0;
1467d4002b98SHong Zhang     info->fill_ratio_needed = 0;
1468d4002b98SHong Zhang     info->factor_mallocs    = 0;
1469d4002b98SHong Zhang   }
14703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1471d4002b98SHong Zhang }
1472d4002b98SHong Zhang 
1473d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetValues_SeqSELL(Mat A, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode is)
1474d71ae5a4SJacob Faibussowitsch {
1475d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1476d4002b98SHong Zhang   PetscInt     shift, i, k, l, low, high, t, ii, row, col, nrow;
1477d4002b98SHong Zhang   PetscInt    *cp, nonew = a->nonew, lastcol = -1;
1478d4002b98SHong Zhang   MatScalar   *vp, value;
1479d4002b98SHong Zhang 
1480d4002b98SHong Zhang   PetscFunctionBegin;
1481d4002b98SHong Zhang   for (k = 0; k < m; k++) { /* loop over added rows */
1482d4002b98SHong Zhang     row = im[k];
1483d4002b98SHong Zhang     if (row < 0) continue;
14846bdcaf15SBarry 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);
1485d4002b98SHong Zhang     shift = a->sliidx[row >> 3] + (row & 0x07); /* starting index of the row */
1486d4002b98SHong Zhang     cp    = a->colidx + shift;                  /* pointer to the row */
1487d4002b98SHong Zhang     vp    = a->val + shift;                     /* pointer to the row */
1488d4002b98SHong Zhang     nrow  = a->rlen[row];
1489d4002b98SHong Zhang     low   = 0;
1490d4002b98SHong Zhang     high  = nrow;
1491d4002b98SHong Zhang 
1492d4002b98SHong Zhang     for (l = 0; l < n; l++) { /* loop over added columns */
1493d4002b98SHong Zhang       col = in[l];
1494d4002b98SHong Zhang       if (col < 0) continue;
14956bdcaf15SBarry 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);
1496d4002b98SHong Zhang       if (a->roworiented) {
1497d4002b98SHong Zhang         value = v[l + k * n];
1498d4002b98SHong Zhang       } else {
1499d4002b98SHong Zhang         value = v[k + l * m];
1500d4002b98SHong Zhang       }
1501d4002b98SHong Zhang       if ((value == 0.0 && a->ignorezeroentries) && (is == ADD_VALUES)) continue;
1502d4002b98SHong Zhang 
1503ed73aabaSBarry Smith       /* search in this row for the specified column, i indicates the column to be set */
1504d4002b98SHong Zhang       if (col <= lastcol) low = 0;
1505d4002b98SHong Zhang       else high = nrow;
1506d4002b98SHong Zhang       lastcol = col;
1507d4002b98SHong Zhang       while (high - low > 5) {
1508d4002b98SHong Zhang         t = (low + high) / 2;
1509d4002b98SHong Zhang         if (*(cp + t * 8) > col) high = t;
1510d4002b98SHong Zhang         else low = t;
1511d4002b98SHong Zhang       }
1512d4002b98SHong Zhang       for (i = low; i < high; i++) {
1513d4002b98SHong Zhang         if (*(cp + i * 8) > col) break;
1514d4002b98SHong Zhang         if (*(cp + i * 8) == col) {
1515d4002b98SHong Zhang           if (is == ADD_VALUES) *(vp + i * 8) += value;
1516d4002b98SHong Zhang           else *(vp + i * 8) = value;
1517d4002b98SHong Zhang           low = i + 1;
1518d4002b98SHong Zhang           goto noinsert;
1519d4002b98SHong Zhang         }
1520d4002b98SHong Zhang       }
1521d4002b98SHong Zhang       if (value == 0.0 && a->ignorezeroentries) goto noinsert;
1522d4002b98SHong Zhang       if (nonew == 1) goto noinsert;
152308401ef6SPierre Jolivet       PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero (%" PetscInt_FMT ", %" PetscInt_FMT ") in the matrix", row, col);
1524d4002b98SHong Zhang       /* If the current row length exceeds the slice width (e.g. nrow==slice_width), allocate a new space, otherwise do nothing */
1525d4002b98SHong Zhang       MatSeqXSELLReallocateSELL(A, A->rmap->n, 1, nrow, a->sliidx, row / 8, row, col, a->colidx, a->val, cp, vp, nonew, MatScalar);
1526d4002b98SHong Zhang       /* add the new nonzero to the high position, shift the remaining elements in current row to the right by one slot */
1527d4002b98SHong Zhang       for (ii = nrow - 1; ii >= i; ii--) {
1528d4002b98SHong Zhang         *(cp + (ii + 1) * 8) = *(cp + ii * 8);
1529d4002b98SHong Zhang         *(vp + (ii + 1) * 8) = *(vp + ii * 8);
1530d4002b98SHong Zhang       }
1531d4002b98SHong Zhang       a->rlen[row]++;
1532d4002b98SHong Zhang       *(cp + i * 8) = col;
1533d4002b98SHong Zhang       *(vp + i * 8) = value;
1534d4002b98SHong Zhang       a->nz++;
1535d4002b98SHong Zhang       A->nonzerostate++;
15369371c9d4SSatish Balay       low = i + 1;
15379371c9d4SSatish Balay       high++;
15389371c9d4SSatish Balay       nrow++;
1539d4002b98SHong Zhang     noinsert:;
1540d4002b98SHong Zhang     }
1541d4002b98SHong Zhang     a->rlen[row] = nrow;
1542d4002b98SHong Zhang   }
15433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1544d4002b98SHong Zhang }
1545d4002b98SHong Zhang 
1546d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCopy_SeqSELL(Mat A, Mat B, MatStructure str)
1547d71ae5a4SJacob Faibussowitsch {
1548d4002b98SHong Zhang   PetscFunctionBegin;
1549d4002b98SHong Zhang   /* If the two matrices have the same copy implementation, use fast copy. */
1550d4002b98SHong Zhang   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1551d4002b98SHong Zhang     Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1552d4002b98SHong Zhang     Mat_SeqSELL *b = (Mat_SeqSELL *)B->data;
1553d4002b98SHong Zhang 
155408401ef6SPierre 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");
15559566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(b->val, a->val, a->sliidx[a->totalslices]));
1556d4002b98SHong Zhang   } else {
15579566063dSJacob Faibussowitsch     PetscCall(MatCopy_Basic(A, B, str));
1558d4002b98SHong Zhang   }
15593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1560d4002b98SHong Zhang }
1561d4002b98SHong Zhang 
1562d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUp_SeqSELL(Mat A)
1563d71ae5a4SJacob Faibussowitsch {
1564d4002b98SHong Zhang   PetscFunctionBegin;
15659566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(A, PETSC_DEFAULT, NULL));
15663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1567d4002b98SHong Zhang }
1568d4002b98SHong Zhang 
1569d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLGetArray_SeqSELL(Mat A, PetscScalar *array[])
1570d71ae5a4SJacob Faibussowitsch {
1571d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1572d4002b98SHong Zhang 
1573d4002b98SHong Zhang   PetscFunctionBegin;
1574d4002b98SHong Zhang   *array = a->val;
15753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1576d4002b98SHong Zhang }
1577d4002b98SHong Zhang 
1578d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLRestoreArray_SeqSELL(Mat A, PetscScalar *array[])
1579d71ae5a4SJacob Faibussowitsch {
1580d4002b98SHong Zhang   PetscFunctionBegin;
15813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1582d4002b98SHong Zhang }
1583d4002b98SHong Zhang 
1584d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRealPart_SeqSELL(Mat A)
1585d71ae5a4SJacob Faibussowitsch {
1586d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1587d4002b98SHong Zhang   PetscInt     i;
1588d4002b98SHong Zhang   MatScalar   *aval = a->val;
1589d4002b98SHong Zhang 
1590d4002b98SHong Zhang   PetscFunctionBegin;
1591d4002b98SHong Zhang   for (i = 0; i < a->sliidx[a->totalslices]; i++) aval[i] = PetscRealPart(aval[i]);
15923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1593d4002b98SHong Zhang }
1594d4002b98SHong Zhang 
1595d71ae5a4SJacob Faibussowitsch PetscErrorCode MatImaginaryPart_SeqSELL(Mat A)
1596d71ae5a4SJacob Faibussowitsch {
1597d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
1598d4002b98SHong Zhang   PetscInt     i;
1599d4002b98SHong Zhang   MatScalar   *aval = a->val;
1600d4002b98SHong Zhang 
1601d4002b98SHong Zhang   PetscFunctionBegin;
1602d4002b98SHong Zhang   for (i = 0; i < a->sliidx[a->totalslices]; i++) aval[i] = PetscImaginaryPart(aval[i]);
16039566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(A));
16043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1605d4002b98SHong Zhang }
1606d4002b98SHong Zhang 
1607d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_SeqSELL(Mat inA, PetscScalar alpha)
1608d71ae5a4SJacob Faibussowitsch {
1609d4002b98SHong Zhang   Mat_SeqSELL *a      = (Mat_SeqSELL *)inA->data;
1610d4002b98SHong Zhang   MatScalar   *aval   = a->val;
1611d4002b98SHong Zhang   PetscScalar  oalpha = alpha;
1612d4002b98SHong Zhang   PetscBLASInt one    = 1, size;
1613d4002b98SHong Zhang 
1614d4002b98SHong Zhang   PetscFunctionBegin;
16159566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(a->sliidx[a->totalslices], &size));
1616792fecdfSBarry Smith   PetscCallBLAS("BLASscal", BLASscal_(&size, &oalpha, aval, &one));
16179566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(a->nz));
16189566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLInvalidateDiagonal(inA));
16193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1620d4002b98SHong Zhang }
1621d4002b98SHong Zhang 
1622d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift_SeqSELL(Mat Y, PetscScalar a)
1623d71ae5a4SJacob Faibussowitsch {
1624d4002b98SHong Zhang   Mat_SeqSELL *y = (Mat_SeqSELL *)Y->data;
1625d4002b98SHong Zhang 
1626d4002b98SHong Zhang   PetscFunctionBegin;
162748a46eb9SPierre Jolivet   if (!Y->preallocated || !y->nz) PetscCall(MatSeqSELLSetPreallocation(Y, 1, NULL));
16289566063dSJacob Faibussowitsch   PetscCall(MatShift_Basic(Y, a));
16293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1630d4002b98SHong Zhang }
1631d4002b98SHong Zhang 
1632d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSOR_SeqSELL(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
1633d71ae5a4SJacob Faibussowitsch {
1634d4002b98SHong Zhang   Mat_SeqSELL       *a = (Mat_SeqSELL *)A->data;
1635d4002b98SHong Zhang   PetscScalar       *x, sum, *t;
1636f4259b30SLisandro Dalcin   const MatScalar   *idiag = NULL, *mdiag;
1637d4002b98SHong Zhang   const PetscScalar *b, *xb;
1638d4002b98SHong Zhang   PetscInt           n, m = A->rmap->n, i, j, shift;
1639d4002b98SHong Zhang   const PetscInt    *diag;
1640d4002b98SHong Zhang 
1641d4002b98SHong Zhang   PetscFunctionBegin;
1642d4002b98SHong Zhang   its = its * lits;
1643d4002b98SHong Zhang 
1644d4002b98SHong Zhang   if (fshift != a->fshift || omega != a->omega) a->idiagvalid = PETSC_FALSE; /* must recompute idiag[] */
16459566063dSJacob Faibussowitsch   if (!a->idiagvalid) PetscCall(MatInvertDiagonal_SeqSELL(A, omega, fshift));
1646d4002b98SHong Zhang   a->fshift = fshift;
1647d4002b98SHong Zhang   a->omega  = omega;
1648d4002b98SHong Zhang 
1649d4002b98SHong Zhang   diag  = a->diag;
1650d4002b98SHong Zhang   t     = a->ssor_work;
1651d4002b98SHong Zhang   idiag = a->idiag;
1652d4002b98SHong Zhang   mdiag = a->mdiag;
1653d4002b98SHong Zhang 
16549566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
16559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(bb, &b));
1656d4002b98SHong Zhang   /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */
165708401ef6SPierre Jolivet   PetscCheck(flag != SOR_APPLY_UPPER, PETSC_COMM_SELF, PETSC_ERR_SUP, "SOR_APPLY_UPPER is not implemented");
165808401ef6SPierre Jolivet   PetscCheck(flag != SOR_APPLY_LOWER, PETSC_COMM_SELF, PETSC_ERR_SUP, "SOR_APPLY_LOWER is not implemented");
1659aed4548fSBarry Smith   PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat");
1660d4002b98SHong Zhang 
1661d4002b98SHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
1662d4002b98SHong Zhang     if ((flag & SOR_FORWARD_SWEEP) || (flag & SOR_LOCAL_FORWARD_SWEEP)) {
1663d4002b98SHong Zhang       for (i = 0; i < m; i++) {
1664d4002b98SHong Zhang         shift = a->sliidx[i >> 3] + (i & 0x07); /* starting index of the row i */
1665d4002b98SHong Zhang         sum   = b[i];
1666d4002b98SHong Zhang         n     = (diag[i] - shift) / 8;
1667d4002b98SHong Zhang         for (j = 0; j < n; j++) sum -= a->val[shift + j * 8] * x[a->colidx[shift + j * 8]];
1668d4002b98SHong Zhang         t[i] = sum;
1669d4002b98SHong Zhang         x[i] = sum * idiag[i];
1670d4002b98SHong Zhang       }
1671d4002b98SHong Zhang       xb = t;
16729566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(a->nz));
1673d4002b98SHong Zhang     } else xb = b;
1674d4002b98SHong Zhang     if ((flag & SOR_BACKWARD_SWEEP) || (flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1675d4002b98SHong Zhang       for (i = m - 1; i >= 0; i--) {
1676d4002b98SHong Zhang         shift = a->sliidx[i >> 3] + (i & 0x07); /* starting index of the row i */
1677d4002b98SHong Zhang         sum   = xb[i];
1678d4002b98SHong Zhang         n     = a->rlen[i] - (diag[i] - shift) / 8 - 1;
1679d4002b98SHong Zhang         for (j = 1; j <= n; j++) sum -= a->val[diag[i] + j * 8] * x[a->colidx[diag[i] + j * 8]];
1680d4002b98SHong Zhang         if (xb == b) {
1681d4002b98SHong Zhang           x[i] = sum * idiag[i];
1682d4002b98SHong Zhang         } else {
1683d4002b98SHong Zhang           x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
1684d4002b98SHong Zhang         }
1685d4002b98SHong Zhang       }
16869566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
1687d4002b98SHong Zhang     }
1688d4002b98SHong Zhang     its--;
1689d4002b98SHong Zhang   }
1690d4002b98SHong Zhang   while (its--) {
1691d4002b98SHong Zhang     if ((flag & SOR_FORWARD_SWEEP) || (flag & SOR_LOCAL_FORWARD_SWEEP)) {
1692d4002b98SHong Zhang       for (i = 0; i < m; i++) {
1693d4002b98SHong Zhang         /* lower */
1694d4002b98SHong Zhang         shift = a->sliidx[i >> 3] + (i & 0x07); /* starting index of the row i */
1695d4002b98SHong Zhang         sum   = b[i];
1696d4002b98SHong Zhang         n     = (diag[i] - shift) / 8;
1697d4002b98SHong Zhang         for (j = 0; j < n; j++) sum -= a->val[shift + j * 8] * x[a->colidx[shift + j * 8]];
1698d4002b98SHong Zhang         t[i] = sum; /* save application of the lower-triangular part */
1699d4002b98SHong Zhang         /* upper */
1700d4002b98SHong Zhang         n = a->rlen[i] - (diag[i] - shift) / 8 - 1;
1701d4002b98SHong Zhang         for (j = 1; j <= n; j++) sum -= a->val[diag[i] + j * 8] * x[a->colidx[diag[i] + j * 8]];
1702d4002b98SHong Zhang         x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
1703d4002b98SHong Zhang       }
1704d4002b98SHong Zhang       xb = t;
17059566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(2.0 * a->nz));
1706d4002b98SHong Zhang     } else xb = b;
1707d4002b98SHong Zhang     if ((flag & SOR_BACKWARD_SWEEP) || (flag & SOR_LOCAL_BACKWARD_SWEEP)) {
1708d4002b98SHong Zhang       for (i = m - 1; i >= 0; i--) {
1709d4002b98SHong Zhang         shift = a->sliidx[i >> 3] + (i & 0x07); /* starting index of the row i */
1710d4002b98SHong Zhang         sum   = xb[i];
1711d4002b98SHong Zhang         if (xb == b) {
1712d4002b98SHong Zhang           /* whole matrix (no checkpointing available) */
1713d4002b98SHong Zhang           n = a->rlen[i];
1714d4002b98SHong Zhang           for (j = 0; j < n; j++) sum -= a->val[shift + j * 8] * x[a->colidx[shift + j * 8]];
1715d4002b98SHong Zhang           x[i] = (1. - omega) * x[i] + (sum + mdiag[i] * x[i]) * idiag[i];
1716d4002b98SHong Zhang         } else { /* lower-triangular part has been saved, so only apply upper-triangular */
1717d4002b98SHong Zhang           n = a->rlen[i] - (diag[i] - shift) / 8 - 1;
1718d4002b98SHong Zhang           for (j = 1; j <= n; j++) sum -= a->val[diag[i] + j * 8] * x[a->colidx[diag[i] + j * 8]];
1719d4002b98SHong Zhang           x[i] = (1. - omega) * x[i] + sum * idiag[i]; /* omega in idiag */
1720d4002b98SHong Zhang         }
1721d4002b98SHong Zhang       }
1722d4002b98SHong Zhang       if (xb == b) {
17239566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(2.0 * a->nz));
1724d4002b98SHong Zhang       } else {
17259566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(a->nz)); /* assumes 1/2 in upper */
1726d4002b98SHong Zhang       }
1727d4002b98SHong Zhang     }
1728d4002b98SHong Zhang   }
17299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
17309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(bb, &b));
17313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1732d4002b98SHong Zhang }
1733d4002b98SHong Zhang 
1734d4002b98SHong Zhang static struct _MatOps MatOps_Values = {MatSetValues_SeqSELL,
17356108893eSStefano Zampini                                        MatGetRow_SeqSELL,
17366108893eSStefano Zampini                                        MatRestoreRow_SeqSELL,
1737d4002b98SHong Zhang                                        MatMult_SeqSELL,
1738d4002b98SHong Zhang                                        /* 4*/ MatMultAdd_SeqSELL,
1739d4002b98SHong Zhang                                        MatMultTranspose_SeqSELL,
1740d4002b98SHong Zhang                                        MatMultTransposeAdd_SeqSELL,
1741f4259b30SLisandro Dalcin                                        NULL,
1742f4259b30SLisandro Dalcin                                        NULL,
1743f4259b30SLisandro Dalcin                                        NULL,
1744f4259b30SLisandro Dalcin                                        /* 10*/ NULL,
1745f4259b30SLisandro Dalcin                                        NULL,
1746f4259b30SLisandro Dalcin                                        NULL,
1747d4002b98SHong Zhang                                        MatSOR_SeqSELL,
1748f4259b30SLisandro Dalcin                                        NULL,
1749d4002b98SHong Zhang                                        /* 15*/ MatGetInfo_SeqSELL,
1750d4002b98SHong Zhang                                        MatEqual_SeqSELL,
1751d4002b98SHong Zhang                                        MatGetDiagonal_SeqSELL,
1752d4002b98SHong Zhang                                        MatDiagonalScale_SeqSELL,
1753f4259b30SLisandro Dalcin                                        NULL,
1754f4259b30SLisandro Dalcin                                        /* 20*/ NULL,
1755d4002b98SHong Zhang                                        MatAssemblyEnd_SeqSELL,
1756d4002b98SHong Zhang                                        MatSetOption_SeqSELL,
1757d4002b98SHong Zhang                                        MatZeroEntries_SeqSELL,
1758f4259b30SLisandro Dalcin                                        /* 24*/ NULL,
1759f4259b30SLisandro Dalcin                                        NULL,
1760f4259b30SLisandro Dalcin                                        NULL,
1761f4259b30SLisandro Dalcin                                        NULL,
1762f4259b30SLisandro Dalcin                                        NULL,
1763d4002b98SHong Zhang                                        /* 29*/ MatSetUp_SeqSELL,
1764f4259b30SLisandro Dalcin                                        NULL,
1765f4259b30SLisandro Dalcin                                        NULL,
1766f4259b30SLisandro Dalcin                                        NULL,
1767f4259b30SLisandro Dalcin                                        NULL,
1768d4002b98SHong Zhang                                        /* 34*/ MatDuplicate_SeqSELL,
1769f4259b30SLisandro Dalcin                                        NULL,
1770f4259b30SLisandro Dalcin                                        NULL,
1771f4259b30SLisandro Dalcin                                        NULL,
1772f4259b30SLisandro Dalcin                                        NULL,
1773f4259b30SLisandro Dalcin                                        /* 39*/ NULL,
1774f4259b30SLisandro Dalcin                                        NULL,
1775f4259b30SLisandro Dalcin                                        NULL,
1776d4002b98SHong Zhang                                        MatGetValues_SeqSELL,
1777d4002b98SHong Zhang                                        MatCopy_SeqSELL,
1778f4259b30SLisandro Dalcin                                        /* 44*/ NULL,
1779d4002b98SHong Zhang                                        MatScale_SeqSELL,
1780d4002b98SHong Zhang                                        MatShift_SeqSELL,
1781f4259b30SLisandro Dalcin                                        NULL,
1782f4259b30SLisandro Dalcin                                        NULL,
1783f4259b30SLisandro Dalcin                                        /* 49*/ NULL,
1784f4259b30SLisandro Dalcin                                        NULL,
1785f4259b30SLisandro Dalcin                                        NULL,
1786f4259b30SLisandro Dalcin                                        NULL,
1787f4259b30SLisandro Dalcin                                        NULL,
1788d4002b98SHong Zhang                                        /* 54*/ MatFDColoringCreate_SeqXAIJ,
1789f4259b30SLisandro Dalcin                                        NULL,
1790f4259b30SLisandro Dalcin                                        NULL,
1791f4259b30SLisandro Dalcin                                        NULL,
1792f4259b30SLisandro Dalcin                                        NULL,
1793f4259b30SLisandro Dalcin                                        /* 59*/ NULL,
1794d4002b98SHong Zhang                                        MatDestroy_SeqSELL,
1795d4002b98SHong Zhang                                        MatView_SeqSELL,
1796f4259b30SLisandro Dalcin                                        NULL,
1797f4259b30SLisandro Dalcin                                        NULL,
1798f4259b30SLisandro Dalcin                                        /* 64*/ NULL,
1799f4259b30SLisandro Dalcin                                        NULL,
1800f4259b30SLisandro Dalcin                                        NULL,
1801f4259b30SLisandro Dalcin                                        NULL,
1802f4259b30SLisandro Dalcin                                        NULL,
1803f4259b30SLisandro Dalcin                                        /* 69*/ NULL,
1804f4259b30SLisandro Dalcin                                        NULL,
1805f4259b30SLisandro Dalcin                                        NULL,
1806f4259b30SLisandro Dalcin                                        NULL,
1807f4259b30SLisandro Dalcin                                        NULL,
1808f4259b30SLisandro Dalcin                                        /* 74*/ NULL,
1809d4002b98SHong Zhang                                        MatFDColoringApply_AIJ, /* reuse the FDColoring function for AIJ */
1810f4259b30SLisandro Dalcin                                        NULL,
1811f4259b30SLisandro Dalcin                                        NULL,
1812f4259b30SLisandro Dalcin                                        NULL,
1813f4259b30SLisandro Dalcin                                        /* 79*/ NULL,
1814f4259b30SLisandro Dalcin                                        NULL,
1815f4259b30SLisandro Dalcin                                        NULL,
1816f4259b30SLisandro Dalcin                                        NULL,
1817f4259b30SLisandro Dalcin                                        NULL,
1818f4259b30SLisandro Dalcin                                        /* 84*/ NULL,
1819f4259b30SLisandro Dalcin                                        NULL,
1820f4259b30SLisandro Dalcin                                        NULL,
1821f4259b30SLisandro Dalcin                                        NULL,
1822f4259b30SLisandro Dalcin                                        NULL,
1823f4259b30SLisandro Dalcin                                        /* 89*/ NULL,
1824f4259b30SLisandro Dalcin                                        NULL,
1825f4259b30SLisandro Dalcin                                        NULL,
1826f4259b30SLisandro Dalcin                                        NULL,
1827f4259b30SLisandro Dalcin                                        NULL,
1828f4259b30SLisandro Dalcin                                        /* 94*/ NULL,
1829f4259b30SLisandro Dalcin                                        NULL,
1830f4259b30SLisandro Dalcin                                        NULL,
1831f4259b30SLisandro Dalcin                                        NULL,
1832f4259b30SLisandro Dalcin                                        NULL,
1833f4259b30SLisandro Dalcin                                        /* 99*/ NULL,
1834f4259b30SLisandro Dalcin                                        NULL,
1835f4259b30SLisandro Dalcin                                        NULL,
1836d4002b98SHong Zhang                                        MatConjugate_SeqSELL,
1837f4259b30SLisandro Dalcin                                        NULL,
1838f4259b30SLisandro Dalcin                                        /*104*/ NULL,
1839f4259b30SLisandro Dalcin                                        NULL,
1840f4259b30SLisandro Dalcin                                        NULL,
1841f4259b30SLisandro Dalcin                                        NULL,
1842f4259b30SLisandro Dalcin                                        NULL,
1843f4259b30SLisandro Dalcin                                        /*109*/ NULL,
1844f4259b30SLisandro Dalcin                                        NULL,
1845f4259b30SLisandro Dalcin                                        NULL,
1846f4259b30SLisandro Dalcin                                        NULL,
1847d4002b98SHong Zhang                                        MatMissingDiagonal_SeqSELL,
1848f4259b30SLisandro Dalcin                                        /*114*/ NULL,
1849f4259b30SLisandro Dalcin                                        NULL,
1850f4259b30SLisandro Dalcin                                        NULL,
1851f4259b30SLisandro Dalcin                                        NULL,
1852f4259b30SLisandro Dalcin                                        NULL,
1853f4259b30SLisandro Dalcin                                        /*119*/ NULL,
1854f4259b30SLisandro Dalcin                                        NULL,
1855f4259b30SLisandro Dalcin                                        NULL,
1856f4259b30SLisandro Dalcin                                        NULL,
1857f4259b30SLisandro Dalcin                                        NULL,
1858f4259b30SLisandro Dalcin                                        /*124*/ NULL,
1859f4259b30SLisandro Dalcin                                        NULL,
1860f4259b30SLisandro Dalcin                                        NULL,
1861f4259b30SLisandro Dalcin                                        NULL,
1862f4259b30SLisandro Dalcin                                        NULL,
1863f4259b30SLisandro Dalcin                                        /*129*/ NULL,
1864f4259b30SLisandro Dalcin                                        NULL,
1865f4259b30SLisandro Dalcin                                        NULL,
1866f4259b30SLisandro Dalcin                                        NULL,
1867f4259b30SLisandro Dalcin                                        NULL,
1868f4259b30SLisandro Dalcin                                        /*134*/ NULL,
1869f4259b30SLisandro Dalcin                                        NULL,
1870f4259b30SLisandro Dalcin                                        NULL,
1871f4259b30SLisandro Dalcin                                        NULL,
1872f4259b30SLisandro Dalcin                                        NULL,
1873f4259b30SLisandro Dalcin                                        /*139*/ NULL,
1874f4259b30SLisandro Dalcin                                        NULL,
1875f4259b30SLisandro Dalcin                                        NULL,
1876d4002b98SHong Zhang                                        MatFDColoringSetUp_SeqXAIJ,
1877f4259b30SLisandro Dalcin                                        NULL,
1878d70f29a3SPierre Jolivet                                        /*144*/ NULL,
1879d70f29a3SPierre Jolivet                                        NULL,
1880d70f29a3SPierre Jolivet                                        NULL,
188199a7f59eSMark Adams                                        NULL,
188299a7f59eSMark Adams                                        NULL,
18837fb60732SBarry Smith                                        NULL,
1884dec0b466SHong Zhang                                        /*150*/ NULL,
1885dec0b466SHong Zhang                                        NULL};
1886d4002b98SHong Zhang 
1887d71ae5a4SJacob Faibussowitsch PetscErrorCode MatStoreValues_SeqSELL(Mat mat)
1888d71ae5a4SJacob Faibussowitsch {
1889d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
1890d4002b98SHong Zhang 
1891d4002b98SHong Zhang   PetscFunctionBegin;
189228b400f6SJacob Faibussowitsch   PetscCheck(a->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1893d4002b98SHong Zhang 
1894d4002b98SHong Zhang   /* allocate space for values if not already there */
1895aa624791SPierre Jolivet   if (!a->saved_values) PetscCall(PetscMalloc1(a->sliidx[a->totalslices] + 1, &a->saved_values));
1896d4002b98SHong Zhang 
1897d4002b98SHong Zhang   /* copy values over */
18989566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->saved_values, a->val, a->sliidx[a->totalslices]));
18993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1900d4002b98SHong Zhang }
1901d4002b98SHong Zhang 
1902d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRetrieveValues_SeqSELL(Mat mat)
1903d71ae5a4SJacob Faibussowitsch {
1904d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)mat->data;
1905d4002b98SHong Zhang 
1906d4002b98SHong Zhang   PetscFunctionBegin;
190728b400f6SJacob Faibussowitsch   PetscCheck(a->nonew, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
190828b400f6SJacob Faibussowitsch   PetscCheck(a->saved_values, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call MatStoreValues(A);first");
19099566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(a->val, a->saved_values, a->sliidx[a->totalslices]));
19103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1911d4002b98SHong Zhang }
1912d4002b98SHong Zhang 
1913d4002b98SHong Zhang /*@C
191411a5261eSBarry Smith  MatSeqSELLRestoreArray - returns access to the array where the data for a `MATSEQSELL` matrix is stored obtained by `MatSeqSELLGetArray()`
1915d4002b98SHong Zhang 
1916d4002b98SHong Zhang  Not Collective
1917d4002b98SHong Zhang 
1918d4002b98SHong Zhang  Input Parameters:
191920f4b53cSBarry Smith +  mat - a `MATSEQSELL` matrix
192020f4b53cSBarry Smith -  array - pointer to the data
1921d4002b98SHong Zhang 
1922d4002b98SHong Zhang  Level: intermediate
1923d4002b98SHong Zhang 
192467be906fSBarry Smith  .seealso: `Mat`, `MATSEQSELL`, `MatSeqSELLGetArray()`, `MatSeqSELLRestoreArrayF90()`
1925d4002b98SHong Zhang  @*/
1926d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLRestoreArray(Mat A, PetscScalar **array)
1927d71ae5a4SJacob Faibussowitsch {
1928d4002b98SHong Zhang   PetscFunctionBegin;
1929cac4c232SBarry Smith   PetscUseMethod(A, "MatSeqSELLRestoreArray_C", (Mat, PetscScalar **), (A, array));
19303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1931d4002b98SHong Zhang }
1932d4002b98SHong Zhang 
1933d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_SeqSELL(Mat B)
1934d71ae5a4SJacob Faibussowitsch {
1935d4002b98SHong Zhang   Mat_SeqSELL *b;
1936d4002b98SHong Zhang   PetscMPIInt  size;
1937d4002b98SHong Zhang 
1938d4002b98SHong Zhang   PetscFunctionBegin;
19399566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
19409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
194108401ef6SPierre Jolivet   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Comm must be of size 1");
1942d4002b98SHong Zhang 
19434dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
1944d4002b98SHong Zhang 
1945d4002b98SHong Zhang   B->data = (void *)b;
1946d4002b98SHong Zhang 
19479566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(B->ops, &MatOps_Values, sizeof(struct _MatOps)));
1948d4002b98SHong Zhang 
1949f4259b30SLisandro Dalcin   b->row                = NULL;
1950f4259b30SLisandro Dalcin   b->col                = NULL;
1951f4259b30SLisandro Dalcin   b->icol               = NULL;
1952d4002b98SHong Zhang   b->reallocs           = 0;
1953d4002b98SHong Zhang   b->ignorezeroentries  = PETSC_FALSE;
1954d4002b98SHong Zhang   b->roworiented        = PETSC_TRUE;
1955d4002b98SHong Zhang   b->nonew              = 0;
1956f4259b30SLisandro Dalcin   b->diag               = NULL;
1957f4259b30SLisandro Dalcin   b->solve_work         = NULL;
1958f4259b30SLisandro Dalcin   B->spptr              = NULL;
1959f4259b30SLisandro Dalcin   b->saved_values       = NULL;
1960f4259b30SLisandro Dalcin   b->idiag              = NULL;
1961f4259b30SLisandro Dalcin   b->mdiag              = NULL;
1962f4259b30SLisandro Dalcin   b->ssor_work          = NULL;
1963d4002b98SHong Zhang   b->omega              = 1.0;
1964d4002b98SHong Zhang   b->fshift             = 0.0;
1965d4002b98SHong Zhang   b->idiagvalid         = PETSC_FALSE;
1966d4002b98SHong Zhang   b->keepnonzeropattern = PETSC_FALSE;
1967d4002b98SHong Zhang 
19689566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQSELL));
19699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLGetArray_C", MatSeqSELLGetArray_SeqSELL));
19709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLRestoreArray_C", MatSeqSELLRestoreArray_SeqSELL));
19719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_SeqSELL));
19729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_SeqSELL));
19739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatSeqSELLSetPreallocation_C", MatSeqSELLSetPreallocation_SeqSELL));
19749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqsell_seqaij_C", MatConvert_SeqSELL_SeqAIJ));
19753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1976d4002b98SHong Zhang }
1977d4002b98SHong Zhang 
1978d4002b98SHong Zhang /*
1979d4002b98SHong Zhang  Given a matrix generated with MatGetFactor() duplicates all the information in A into B
1980d4002b98SHong Zhang  */
1981d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicateNoCreate_SeqSELL(Mat C, Mat A, MatDuplicateOption cpvalues, PetscBool mallocmatspace)
1982d71ae5a4SJacob Faibussowitsch {
1983ed73aabaSBarry Smith   Mat_SeqSELL *c = (Mat_SeqSELL *)C->data, *a = (Mat_SeqSELL *)A->data;
1984d4002b98SHong Zhang   PetscInt     i, m                           = A->rmap->n;
1985d4002b98SHong Zhang   PetscInt     totalslices = a->totalslices;
1986d4002b98SHong Zhang 
1987d4002b98SHong Zhang   PetscFunctionBegin;
1988d4002b98SHong Zhang   C->factortype = A->factortype;
1989f4259b30SLisandro Dalcin   c->row        = NULL;
1990f4259b30SLisandro Dalcin   c->col        = NULL;
1991f4259b30SLisandro Dalcin   c->icol       = NULL;
1992d4002b98SHong Zhang   c->reallocs   = 0;
1993d4002b98SHong Zhang   C->assembled  = PETSC_TRUE;
1994d4002b98SHong Zhang 
19959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->rmap, &C->rmap));
19969566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &C->cmap));
1997d4002b98SHong Zhang 
19989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(8 * totalslices, &c->rlen));
19999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(totalslices + 1, &c->sliidx));
2000d4002b98SHong Zhang 
2001d4002b98SHong Zhang   for (i = 0; i < m; i++) c->rlen[i] = a->rlen[i];
2002d4002b98SHong Zhang   for (i = 0; i < totalslices + 1; i++) c->sliidx[i] = a->sliidx[i];
2003d4002b98SHong Zhang 
2004d4002b98SHong Zhang   /* allocate the matrix space */
2005d4002b98SHong Zhang   if (mallocmatspace) {
20069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(a->maxallocmat, &c->val, a->maxallocmat, &c->colidx));
2007d4002b98SHong Zhang 
2008d4002b98SHong Zhang     c->singlemalloc = PETSC_TRUE;
2009d4002b98SHong Zhang 
2010d4002b98SHong Zhang     if (m > 0) {
20119566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(c->colidx, a->colidx, a->maxallocmat));
2012d4002b98SHong Zhang       if (cpvalues == MAT_COPY_VALUES) {
20139566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(c->val, a->val, a->maxallocmat));
2014d4002b98SHong Zhang       } else {
20159566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(c->val, a->maxallocmat));
2016d4002b98SHong Zhang       }
2017d4002b98SHong Zhang     }
2018d4002b98SHong Zhang   }
2019d4002b98SHong Zhang 
2020d4002b98SHong Zhang   c->ignorezeroentries = a->ignorezeroentries;
2021d4002b98SHong Zhang   c->roworiented       = a->roworiented;
2022d4002b98SHong Zhang   c->nonew             = a->nonew;
2023d4002b98SHong Zhang   if (a->diag) {
20249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &c->diag));
2025ad540459SPierre Jolivet     for (i = 0; i < m; i++) c->diag[i] = a->diag[i];
2026f4259b30SLisandro Dalcin   } else c->diag = NULL;
2027d4002b98SHong Zhang 
2028f4259b30SLisandro Dalcin   c->solve_work         = NULL;
2029f4259b30SLisandro Dalcin   c->saved_values       = NULL;
2030f4259b30SLisandro Dalcin   c->idiag              = NULL;
2031f4259b30SLisandro Dalcin   c->ssor_work          = NULL;
2032d4002b98SHong Zhang   c->keepnonzeropattern = a->keepnonzeropattern;
2033d4002b98SHong Zhang   c->free_val           = PETSC_TRUE;
2034d4002b98SHong Zhang   c->free_colidx        = PETSC_TRUE;
2035d4002b98SHong Zhang 
2036d4002b98SHong Zhang   c->maxallocmat  = a->maxallocmat;
2037d4002b98SHong Zhang   c->maxallocrow  = a->maxallocrow;
2038d4002b98SHong Zhang   c->rlenmax      = a->rlenmax;
2039d4002b98SHong Zhang   c->nz           = a->nz;
2040d4002b98SHong Zhang   C->preallocated = PETSC_TRUE;
2041d4002b98SHong Zhang 
2042d4002b98SHong Zhang   c->nonzerorowcnt = a->nonzerorowcnt;
2043d4002b98SHong Zhang   C->nonzerostate  = A->nonzerostate;
2044d4002b98SHong Zhang 
20459566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDuplicate(((PetscObject)A)->qlist, &((PetscObject)C)->qlist));
20463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2047d4002b98SHong Zhang }
2048d4002b98SHong Zhang 
2049d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_SeqSELL(Mat A, MatDuplicateOption cpvalues, Mat *B)
2050d71ae5a4SJacob Faibussowitsch {
2051d4002b98SHong Zhang   PetscFunctionBegin;
20529566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
20539566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->n, A->cmap->n));
205448a46eb9SPierre Jolivet   if (!(A->rmap->n % A->rmap->bs) && !(A->cmap->n % A->cmap->bs)) PetscCall(MatSetBlockSizesFromMats(*B, A, A));
20559566063dSJacob Faibussowitsch   PetscCall(MatSetType(*B, ((PetscObject)A)->type_name));
20569566063dSJacob Faibussowitsch   PetscCall(MatDuplicateNoCreate_SeqSELL(*B, A, cpvalues, PETSC_TRUE));
20573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2058d4002b98SHong Zhang }
2059d4002b98SHong Zhang 
2060ed73aabaSBarry Smith /*MC
2061ed73aabaSBarry Smith    MATSEQSELL - MATSEQSELL = "seqsell" - A matrix type to be used for sequential sparse matrices,
2062ed73aabaSBarry Smith    based on the sliced Ellpack format
2063ed73aabaSBarry Smith 
206420f4b53cSBarry Smith    Options Database Key:
206511a5261eSBarry Smith . -mat_type seqsell - sets the matrix type to "`MATSEQELL` during a call to `MatSetFromOptions()`
2066ed73aabaSBarry Smith 
2067ed73aabaSBarry Smith    Level: beginner
2068ed73aabaSBarry Smith 
206967be906fSBarry Smith .seealso: `Mat`, `MatCreateSeqSell()`, `MATSELL`, `MATMPISELL`, `MATSEQAIJ`, `MATAIJ`, `MATMPIAIJ`
2070ed73aabaSBarry Smith M*/
2071ed73aabaSBarry Smith 
2072ed73aabaSBarry Smith /*MC
2073ed73aabaSBarry Smith    MATSELL - MATSELL = "sell" - A matrix type to be used for sparse matrices.
2074ed73aabaSBarry Smith 
207511a5261eSBarry Smith    This matrix type is identical to `MATSEQSELL` when constructed with a single process communicator,
207611a5261eSBarry Smith    and `MATMPISELL` otherwise.  As a result, for single process communicators,
207711a5261eSBarry Smith   `MatSeqSELLSetPreallocation()` is supported, and similarly `MatMPISELLSetPreallocation()` is supported
2078ed73aabaSBarry Smith   for communicators controlling multiple processes.  It is recommended that you call both of
2079ed73aabaSBarry Smith   the above preallocation routines for simplicity.
2080ed73aabaSBarry Smith 
208120f4b53cSBarry Smith    Options Database Key:
2082ed73aabaSBarry Smith . -mat_type sell - sets the matrix type to "sell" during a call to MatSetFromOptions()
2083ed73aabaSBarry Smith 
2084ed73aabaSBarry Smith   Level: beginner
2085ed73aabaSBarry Smith 
2086ed73aabaSBarry Smith   Notes:
2087ed73aabaSBarry Smith    This format is only supported for real scalars, double precision, and 32 bit indices (the defaults).
2088ed73aabaSBarry Smith 
2089ed73aabaSBarry Smith    It can provide better performance on Intel and AMD processes with AVX2 or AVX512 support for matrices that have a similar number of
2090ed73aabaSBarry Smith    non-zeros in contiguous groups of rows. However if the computation is memory bandwidth limited it may not provide much improvement.
2091ed73aabaSBarry Smith 
2092ed73aabaSBarry Smith   Developer Notes:
2093ed73aabaSBarry Smith    On Intel (and AMD) systems some of the matrix operations use SIMD (AVX) instructions to achieve higher performance.
2094ed73aabaSBarry Smith 
2095ed73aabaSBarry Smith    The sparse matrix format is as follows. For simplicity we assume a slice size of 2, it is actually 8
2096ed73aabaSBarry Smith .vb
2097ed73aabaSBarry Smith                             (2 0  3 4)
2098ed73aabaSBarry Smith    Consider the matrix A =  (5 0  6 0)
2099ed73aabaSBarry Smith                             (0 0  7 8)
2100ed73aabaSBarry Smith                             (0 0  9 9)
2101ed73aabaSBarry Smith 
2102ed73aabaSBarry Smith    symbolically the Ellpack format can be written as
2103ed73aabaSBarry Smith 
2104ed73aabaSBarry Smith         (2 3 4 |)           (0 2 3 |)
2105ed73aabaSBarry Smith    v =  (5 6 0 |)  colidx = (0 2 2 |)
2106ed73aabaSBarry Smith         --------            ---------
2107ed73aabaSBarry Smith         (7 8 |)             (2 3 |)
2108ed73aabaSBarry Smith         (9 9 |)             (2 3 |)
2109ed73aabaSBarry Smith 
2110ed73aabaSBarry 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).
2111ed73aabaSBarry 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
2112ed73aabaSBarry 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.
2113ed73aabaSBarry Smith 
2114ed73aabaSBarry 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)
2115ed73aabaSBarry Smith 
2116ed73aabaSBarry Smith .ve
2117ed73aabaSBarry Smith 
2118ed73aabaSBarry Smith       See MatMult_SeqSELL() for how this format is used with the SIMD operations to achieve high performance.
2119ed73aabaSBarry Smith 
2120ed73aabaSBarry Smith  References:
2121606c0280SSatish Balay . * - Hong Zhang, Richard T. Mills, Karl Rupp, and Barry F. Smith, Vectorized Parallel Sparse Matrix-Vector Multiplication in {PETSc} Using {AVX-512},
2122ed73aabaSBarry Smith    Proceedings of the 47th International Conference on Parallel Processing, 2018.
2123ed73aabaSBarry Smith 
212467be906fSBarry Smith .seealso: `Mat`, `MatCreateSeqSELL()`, `MatCreateSeqAIJ()`, `MatCreateSell()`, `MATSEQSELL`, `MATMPISELL`, `MATSEQAIJ`, `MATMPIAIJ`, `MATAIJ`
2125ed73aabaSBarry Smith M*/
2126ed73aabaSBarry Smith 
2127d4002b98SHong Zhang /*@C
212811a5261eSBarry Smith        MatCreateSeqSELL - Creates a sparse matrix in `MATSEQSELL` format.
2129d4002b98SHong Zhang 
21302ef1f0ffSBarry Smith  Collective
2131d4002b98SHong Zhang 
2132d4002b98SHong Zhang  Input Parameters:
213311a5261eSBarry Smith +  comm - MPI communicator, set to `PETSC_COMM_SELF`
2134d4002b98SHong Zhang .  m - number of rows
2135d4002b98SHong Zhang .  n - number of columns
213620f4b53cSBarry Smith .  rlenmax - maximum number of nonzeros in a row, ignored if `rlen` is provided
213720f4b53cSBarry Smith -  rlen - array containing the number of nonzeros in the various rows (possibly different for each row) or NULL
2138d4002b98SHong Zhang 
2139d4002b98SHong Zhang  Output Parameter:
2140d4002b98SHong Zhang .  A - the matrix
2141d4002b98SHong Zhang 
214220f4b53cSBarry Smith  Level: intermediate
214320f4b53cSBarry Smith 
214420f4b53cSBarry Smith  Notes:
214511a5261eSBarry Smith  It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`,
2146f6f02116SRichard Tran Mills  MatXXXXSetPreallocation() paradigm instead of this routine directly.
214711a5261eSBarry Smith  [MatXXXXSetPreallocation() is, for example, `MatSeqSELLSetPreallocation()`]
2148d4002b98SHong Zhang 
214920f4b53cSBarry Smith  Specify the preallocated storage with either `rlenmax` or `rlen` (not both).
215020f4b53cSBarry Smith  Set `rlenmax` = `PETSC_DEFAULT` and `rlen` = `NULL` for PETSc to control dynamic memory
215120f4b53cSBarry Smith  allocation.
2152d4002b98SHong Zhang 
215367be906fSBarry Smith  .seealso: `Mat`, `MATSEQSELL`, `MatCreate()`, `MatCreateSELL()`, `MatSetValues()`, `MatSeqSELLSetPreallocation()`, `MATSELL`, `MATSEQSELL`, `MATMPISELL`
2154d4002b98SHong Zhang  @*/
215520f4b53cSBarry Smith PetscErrorCode MatCreateSeqSELL(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt rlenmax, const PetscInt rlen[], Mat *A)
2156d71ae5a4SJacob Faibussowitsch {
2157d4002b98SHong Zhang   PetscFunctionBegin;
21589566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
21599566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, n, m, n));
21609566063dSJacob Faibussowitsch   PetscCall(MatSetType(*A, MATSEQSELL));
216120f4b53cSBarry Smith   PetscCall(MatSeqSELLSetPreallocation_SeqSELL(*A, rlenmax, rlen));
21623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2163d4002b98SHong Zhang }
2164d4002b98SHong Zhang 
2165d71ae5a4SJacob Faibussowitsch PetscErrorCode MatEqual_SeqSELL(Mat A, Mat B, PetscBool *flg)
2166d71ae5a4SJacob Faibussowitsch {
2167d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data, *b = (Mat_SeqSELL *)B->data;
2168d4002b98SHong Zhang   PetscInt     totalslices = a->totalslices;
2169d4002b98SHong Zhang 
2170d4002b98SHong Zhang   PetscFunctionBegin;
2171d4002b98SHong Zhang   /* If the  matrix dimensions are not equal,or no of nonzeros */
2172d4002b98SHong Zhang   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n) || (a->nz != b->nz) || (a->rlenmax != b->rlenmax)) {
2173d4002b98SHong Zhang     *flg = PETSC_FALSE;
21743ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2175d4002b98SHong Zhang   }
2176d4002b98SHong Zhang   /* if the a->colidx are the same */
21779566063dSJacob Faibussowitsch   PetscCall(PetscArraycmp(a->colidx, b->colidx, a->sliidx[totalslices], flg));
21783ba16761SJacob Faibussowitsch   if (!*flg) PetscFunctionReturn(PETSC_SUCCESS);
2179d4002b98SHong Zhang   /* if a->val are the same */
21809566063dSJacob Faibussowitsch   PetscCall(PetscArraycmp(a->val, b->val, a->sliidx[totalslices], flg));
21813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2182d4002b98SHong Zhang }
2183d4002b98SHong Zhang 
2184d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSeqSELLInvalidateDiagonal(Mat A)
2185d71ae5a4SJacob Faibussowitsch {
2186d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
2187d4002b98SHong Zhang 
2188d4002b98SHong Zhang   PetscFunctionBegin;
2189d4002b98SHong Zhang   a->idiagvalid = PETSC_FALSE;
21903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2191d4002b98SHong Zhang }
2192d4002b98SHong Zhang 
2193d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConjugate_SeqSELL(Mat A)
2194d71ae5a4SJacob Faibussowitsch {
2195d4002b98SHong Zhang #if defined(PETSC_USE_COMPLEX)
2196d4002b98SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
2197d4002b98SHong Zhang   PetscInt     i;
2198d4002b98SHong Zhang   PetscScalar *val = a->val;
2199d4002b98SHong Zhang 
2200d4002b98SHong Zhang   PetscFunctionBegin;
2201ad540459SPierre Jolivet   for (i = 0; i < a->sliidx[a->totalslices]; i++) val[i] = PetscConj(val[i]);
2202d4002b98SHong Zhang #else
2203d4002b98SHong Zhang   PetscFunctionBegin;
2204d4002b98SHong Zhang #endif
22053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2206d4002b98SHong Zhang }
2207