xref: /petsc/src/mat/impls/sell/mpi/mpisell.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1d4002b98SHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>   /*I "petscmat.h" I*/
2d4002b98SHong Zhang #include <../src/mat/impls/sell/mpi/mpisell.h> /*I "petscmat.h" I*/
3d4002b98SHong Zhang #include <petsc/private/vecimpl.h>
4d4002b98SHong Zhang #include <petsc/private/isimpl.h>
5d4002b98SHong Zhang #include <petscblaslapack.h>
6d4002b98SHong Zhang #include <petscsf.h>
7d4002b98SHong Zhang 
8d4002b98SHong Zhang /*MC
9d4002b98SHong Zhang    MATSELL - MATSELL = "sell" - A matrix type to be used for sparse matrices.
10d4002b98SHong Zhang 
11d4002b98SHong Zhang    This matrix type is identical to MATSEQSELL when constructed with a single process communicator,
12d4002b98SHong Zhang    and MATMPISELL otherwise.  As a result, for single process communicators,
13d4002b98SHong Zhang   MatSeqSELLSetPreallocation is supported, and similarly MatMPISELLSetPreallocation is supported
14d4002b98SHong Zhang   for communicators controlling multiple processes.  It is recommended that you call both of
15d4002b98SHong Zhang   the above preallocation routines for simplicity.
16d4002b98SHong Zhang 
17d4002b98SHong Zhang    Options Database Keys:
18d4002b98SHong Zhang . -mat_type sell - sets the matrix type to "sell" during a call to MatSetFromOptions()
19d4002b98SHong Zhang 
20d4002b98SHong Zhang   Level: beginner
21d4002b98SHong Zhang 
22db781477SPatrick Sanan .seealso: `MatCreateSELL()`, `MatCreateSeqSELL()`, `MATSEQSELL`, `MATMPISELL`
23d4002b98SHong Zhang M*/
24d4002b98SHong Zhang 
259371c9d4SSatish Balay PetscErrorCode MatDiagonalSet_MPISELL(Mat Y, Vec D, InsertMode is) {
26d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)Y->data;
27d4002b98SHong Zhang 
28d4002b98SHong Zhang   PetscFunctionBegin;
29d4002b98SHong Zhang   if (Y->assembled && Y->rmap->rstart == Y->cmap->rstart && Y->rmap->rend == Y->cmap->rend) {
309566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(sell->A, D, is));
31d4002b98SHong Zhang   } else {
329566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet_Default(Y, D, is));
33d4002b98SHong Zhang   }
34d4002b98SHong Zhang   PetscFunctionReturn(0);
35d4002b98SHong Zhang }
36d4002b98SHong Zhang 
37d4002b98SHong Zhang /*
38d4002b98SHong Zhang   Local utility routine that creates a mapping from the global column
39d4002b98SHong Zhang number to the local number in the off-diagonal part of the local
40d4002b98SHong Zhang storage of the matrix.  When PETSC_USE_CTABLE is used this is scalable at
41d4002b98SHong Zhang a slightly higher hash table cost; without it it is not scalable (each processor
42d4002b98SHong Zhang has an order N integer array but is fast to acess.
43d4002b98SHong Zhang */
449371c9d4SSatish Balay PetscErrorCode MatCreateColmap_MPISELL_Private(Mat mat) {
45d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
46d4002b98SHong Zhang   PetscInt     n    = sell->B->cmap->n, i;
47d4002b98SHong Zhang 
48d4002b98SHong Zhang   PetscFunctionBegin;
4928b400f6SJacob Faibussowitsch   PetscCheck(sell->garray, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPISELL Matrix was assembled but is missing garray");
50d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
519566063dSJacob Faibussowitsch   PetscCall(PetscTableCreate(n, mat->cmap->N + 1, &sell->colmap));
52*48a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall(PetscTableAdd(sell->colmap, sell->garray[i] + 1, i + 1, INSERT_VALUES));
53d4002b98SHong Zhang #else
549566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(mat->cmap->N + 1, &sell->colmap));
559566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)mat, (mat->cmap->N + 1) * sizeof(PetscInt)));
56d4002b98SHong Zhang   for (i = 0; i < n; i++) sell->colmap[sell->garray[i]] = i + 1;
57d4002b98SHong Zhang #endif
58d4002b98SHong Zhang   PetscFunctionReturn(0);
59d4002b98SHong Zhang }
60d4002b98SHong Zhang 
61d4002b98SHong Zhang #define MatSetValues_SeqSELL_A_Private(row, col, value, addv, orow, ocol) \
62d4002b98SHong Zhang   { \
63d4002b98SHong Zhang     if (col <= lastcol1) low1 = 0; \
64d4002b98SHong Zhang     else high1 = nrow1; \
65d4002b98SHong Zhang     lastcol1 = col; \
66d4002b98SHong Zhang     while (high1 - low1 > 5) { \
67d4002b98SHong Zhang       t = (low1 + high1) / 2; \
68d4002b98SHong Zhang       if (*(cp1 + 8 * t) > col) high1 = t; \
69d4002b98SHong Zhang       else low1 = t; \
70d4002b98SHong Zhang     } \
71d4002b98SHong Zhang     for (_i = low1; _i < high1; _i++) { \
72d4002b98SHong Zhang       if (*(cp1 + 8 * _i) > col) break; \
73d4002b98SHong Zhang       if (*(cp1 + 8 * _i) == col) { \
74d4002b98SHong Zhang         if (addv == ADD_VALUES) *(vp1 + 8 * _i) += value; \
75d4002b98SHong Zhang         else *(vp1 + 8 * _i) = value; \
76d4002b98SHong Zhang         goto a_noinsert; \
77d4002b98SHong Zhang       } \
78d4002b98SHong Zhang     } \
799371c9d4SSatish Balay     if (value == 0.0 && ignorezeroentries) { \
809371c9d4SSatish Balay       low1  = 0; \
819371c9d4SSatish Balay       high1 = nrow1; \
829371c9d4SSatish Balay       goto a_noinsert; \
839371c9d4SSatish Balay     } \
849371c9d4SSatish Balay     if (nonew == 1) { \
859371c9d4SSatish Balay       low1  = 0; \
869371c9d4SSatish Balay       high1 = nrow1; \
879371c9d4SSatish Balay       goto a_noinsert; \
889371c9d4SSatish Balay     } \
8908401ef6SPierre Jolivet     PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", orow, ocol); \
90d4002b98SHong Zhang     MatSeqXSELLReallocateSELL(A, am, 1, nrow1, a->sliidx, row / 8, row, col, a->colidx, a->val, cp1, vp1, nonew, MatScalar); \
91d4002b98SHong Zhang     /* shift up all the later entries in this row */ \
92d4002b98SHong Zhang     for (ii = nrow1 - 1; ii >= _i; ii--) { \
93d4002b98SHong Zhang       *(cp1 + 8 * (ii + 1)) = *(cp1 + 8 * ii); \
94d4002b98SHong Zhang       *(vp1 + 8 * (ii + 1)) = *(vp1 + 8 * ii); \
95d4002b98SHong Zhang     } \
96d4002b98SHong Zhang     *(cp1 + 8 * _i) = col; \
97d4002b98SHong Zhang     *(vp1 + 8 * _i) = value; \
989371c9d4SSatish Balay     a->nz++; \
999371c9d4SSatish Balay     nrow1++; \
1009371c9d4SSatish Balay     A->nonzerostate++; \
101d4002b98SHong Zhang   a_noinsert:; \
102d4002b98SHong Zhang     a->rlen[row] = nrow1; \
103d4002b98SHong Zhang   }
104d4002b98SHong Zhang 
105d4002b98SHong Zhang #define MatSetValues_SeqSELL_B_Private(row, col, value, addv, orow, ocol) \
106d4002b98SHong Zhang   { \
107d4002b98SHong Zhang     if (col <= lastcol2) low2 = 0; \
108d4002b98SHong Zhang     else high2 = nrow2; \
109d4002b98SHong Zhang     lastcol2 = col; \
110d4002b98SHong Zhang     while (high2 - low2 > 5) { \
111d4002b98SHong Zhang       t = (low2 + high2) / 2; \
112d4002b98SHong Zhang       if (*(cp2 + 8 * t) > col) high2 = t; \
113d4002b98SHong Zhang       else low2 = t; \
114d4002b98SHong Zhang     } \
115d4002b98SHong Zhang     for (_i = low2; _i < high2; _i++) { \
116d4002b98SHong Zhang       if (*(cp2 + 8 * _i) > col) break; \
117d4002b98SHong Zhang       if (*(cp2 + 8 * _i) == col) { \
118d4002b98SHong Zhang         if (addv == ADD_VALUES) *(vp2 + 8 * _i) += value; \
119d4002b98SHong Zhang         else *(vp2 + 8 * _i) = value; \
120d4002b98SHong Zhang         goto b_noinsert; \
121d4002b98SHong Zhang       } \
122d4002b98SHong Zhang     } \
1239371c9d4SSatish Balay     if (value == 0.0 && ignorezeroentries) { \
1249371c9d4SSatish Balay       low2  = 0; \
1259371c9d4SSatish Balay       high2 = nrow2; \
1269371c9d4SSatish Balay       goto b_noinsert; \
1279371c9d4SSatish Balay     } \
1289371c9d4SSatish Balay     if (nonew == 1) { \
1299371c9d4SSatish Balay       low2  = 0; \
1309371c9d4SSatish Balay       high2 = nrow2; \
1319371c9d4SSatish Balay       goto b_noinsert; \
1329371c9d4SSatish Balay     } \
13308401ef6SPierre Jolivet     PetscCheck(nonew != -1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", orow, ocol); \
134d4002b98SHong Zhang     MatSeqXSELLReallocateSELL(B, bm, 1, nrow2, b->sliidx, row / 8, row, col, b->colidx, b->val, cp2, vp2, nonew, MatScalar); \
135d4002b98SHong Zhang     /* shift up all the later entries in this row */ \
136d4002b98SHong Zhang     for (ii = nrow2 - 1; ii >= _i; ii--) { \
137d4002b98SHong Zhang       *(cp2 + 8 * (ii + 1)) = *(cp2 + 8 * ii); \
138d4002b98SHong Zhang       *(vp2 + 8 * (ii + 1)) = *(vp2 + 8 * ii); \
139d4002b98SHong Zhang     } \
140d4002b98SHong Zhang     *(cp2 + 8 * _i) = col; \
141d4002b98SHong Zhang     *(vp2 + 8 * _i) = value; \
1429371c9d4SSatish Balay     b->nz++; \
1439371c9d4SSatish Balay     nrow2++; \
1449371c9d4SSatish Balay     B->nonzerostate++; \
145d4002b98SHong Zhang   b_noinsert:; \
146d4002b98SHong Zhang     b->rlen[row] = nrow2; \
147d4002b98SHong Zhang   }
148d4002b98SHong Zhang 
1499371c9d4SSatish Balay PetscErrorCode MatSetValues_MPISELL(Mat mat, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode addv) {
150d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
151d4002b98SHong Zhang   PetscScalar  value;
152d4002b98SHong Zhang   PetscInt     i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend, shift1, shift2;
153d4002b98SHong Zhang   PetscInt     cstart = mat->cmap->rstart, cend = mat->cmap->rend, row, col;
154d4002b98SHong Zhang   PetscBool    roworiented = sell->roworiented;
155d4002b98SHong Zhang 
156d4002b98SHong Zhang   /* Some Variables required in the macro */
157d4002b98SHong Zhang   Mat          A                 = sell->A;
158d4002b98SHong Zhang   Mat_SeqSELL *a                 = (Mat_SeqSELL *)A->data;
159d4002b98SHong Zhang   PetscBool    ignorezeroentries = a->ignorezeroentries, found;
160d4002b98SHong Zhang   Mat          B                 = sell->B;
161d4002b98SHong Zhang   Mat_SeqSELL *b                 = (Mat_SeqSELL *)B->data;
162d4002b98SHong Zhang   PetscInt    *cp1, *cp2, ii, _i, nrow1, nrow2, low1, high1, low2, high2, t, lastcol1, lastcol2;
163d4002b98SHong Zhang   MatScalar   *vp1, *vp2;
164d4002b98SHong Zhang 
165d4002b98SHong Zhang   PetscFunctionBegin;
166d4002b98SHong Zhang   for (i = 0; i < m; i++) {
167d4002b98SHong Zhang     if (im[i] < 0) continue;
1686bdcaf15SBarry Smith     PetscCheck(im[i] < mat->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, im[i], mat->rmap->N - 1);
169d4002b98SHong Zhang     if (im[i] >= rstart && im[i] < rend) {
170d4002b98SHong Zhang       row      = im[i] - rstart;
171d4002b98SHong Zhang       lastcol1 = -1;
172d4002b98SHong Zhang       shift1   = a->sliidx[row >> 3] + (row & 0x07); /* starting index of the row */
173d4002b98SHong Zhang       cp1      = a->colidx + shift1;
174d4002b98SHong Zhang       vp1      = a->val + shift1;
175d4002b98SHong Zhang       nrow1    = a->rlen[row];
176d4002b98SHong Zhang       low1     = 0;
177d4002b98SHong Zhang       high1    = nrow1;
178d4002b98SHong Zhang       lastcol2 = -1;
179d4002b98SHong Zhang       shift2   = b->sliidx[row >> 3] + (row & 0x07); /* starting index of the row */
180d4002b98SHong Zhang       cp2      = b->colidx + shift2;
181d4002b98SHong Zhang       vp2      = b->val + shift2;
182d4002b98SHong Zhang       nrow2    = b->rlen[row];
183d4002b98SHong Zhang       low2     = 0;
184d4002b98SHong Zhang       high2    = nrow2;
185d4002b98SHong Zhang 
186d4002b98SHong Zhang       for (j = 0; j < n; j++) {
187d4002b98SHong Zhang         if (roworiented) value = v[i * n + j];
188d4002b98SHong Zhang         else value = v[i + j * m];
189d4002b98SHong Zhang         if (ignorezeroentries && value == 0.0 && (addv == ADD_VALUES)) continue;
190d4002b98SHong Zhang         if (in[j] >= cstart && in[j] < cend) {
191d4002b98SHong Zhang           col = in[j] - cstart;
192d4002b98SHong Zhang           MatSetValue_SeqSELL_Private(A, row, col, value, addv, im[i], in[j], cp1, vp1, lastcol1, low1, high1); /* set one value */
193f7d195e4SLawrence Mitchell         } else if (in[j] < 0) {
194f7d195e4SLawrence Mitchell           continue;
195f7d195e4SLawrence Mitchell         } else {
196f7d195e4SLawrence Mitchell           PetscCheck(in[j] < mat->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, in[j], mat->cmap->N - 1);
197d4002b98SHong Zhang           if (mat->was_assembled) {
198*48a46eb9SPierre Jolivet             if (!sell->colmap) PetscCall(MatCreateColmap_MPISELL_Private(mat));
199d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
2009566063dSJacob Faibussowitsch             PetscCall(PetscTableFind(sell->colmap, in[j] + 1, &col));
201d4002b98SHong Zhang             col--;
202d4002b98SHong Zhang #else
203d4002b98SHong Zhang             col = sell->colmap[in[j]] - 1;
204d4002b98SHong Zhang #endif
205d4002b98SHong Zhang             if (col < 0 && !((Mat_SeqSELL *)(sell->B->data))->nonew) {
2069566063dSJacob Faibussowitsch               PetscCall(MatDisAssemble_MPISELL(mat));
207d4002b98SHong Zhang               col    = in[j];
208d4002b98SHong Zhang               /* Reinitialize the variables required by MatSetValues_SeqSELL_B_Private() */
209d4002b98SHong Zhang               B      = sell->B;
210d4002b98SHong Zhang               b      = (Mat_SeqSELL *)B->data;
211d4002b98SHong Zhang               shift2 = b->sliidx[row >> 3] + (row & 0x07); /* starting index of the row */
212d4002b98SHong Zhang               cp2    = b->colidx + shift2;
213d4002b98SHong Zhang               vp2    = b->val + shift2;
214d4002b98SHong Zhang               nrow2  = b->rlen[row];
215d4002b98SHong Zhang               low2   = 0;
216d4002b98SHong Zhang               high2  = nrow2;
217f7d195e4SLawrence Mitchell             } else {
218f7d195e4SLawrence Mitchell               PetscCheck(col >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Inserting a new nonzero at global row/column (%" PetscInt_FMT ", %" PetscInt_FMT ") into matrix", im[i], in[j]);
219f7d195e4SLawrence Mitchell             }
220d4002b98SHong Zhang           } else col = in[j];
221d4002b98SHong Zhang           MatSetValue_SeqSELL_Private(B, row, col, value, addv, im[i], in[j], cp2, vp2, lastcol2, low2, high2); /* set one value */
222d4002b98SHong Zhang         }
223d4002b98SHong Zhang       }
224d4002b98SHong Zhang     } else {
22528b400f6SJacob Faibussowitsch       PetscCheck(!mat->nooffprocentries, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Setting off process row %" PetscInt_FMT " even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set", im[i]);
226d4002b98SHong Zhang       if (!sell->donotstash) {
227d4002b98SHong Zhang         mat->assembled = PETSC_FALSE;
228d4002b98SHong Zhang         if (roworiented) {
2299566063dSJacob Faibussowitsch           PetscCall(MatStashValuesRow_Private(&mat->stash, im[i], n, in, v + i * n, (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
230d4002b98SHong Zhang         } else {
2319566063dSJacob Faibussowitsch           PetscCall(MatStashValuesCol_Private(&mat->stash, im[i], n, in, v + i, m, (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
232d4002b98SHong Zhang         }
233d4002b98SHong Zhang       }
234d4002b98SHong Zhang     }
235d4002b98SHong Zhang   }
236d4002b98SHong Zhang   PetscFunctionReturn(0);
237d4002b98SHong Zhang }
238d4002b98SHong Zhang 
2399371c9d4SSatish Balay PetscErrorCode MatGetValues_MPISELL(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[]) {
240d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
241d4002b98SHong Zhang   PetscInt     i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend;
242d4002b98SHong Zhang   PetscInt     cstart = mat->cmap->rstart, cend = mat->cmap->rend, row, col;
243d4002b98SHong Zhang 
244d4002b98SHong Zhang   PetscFunctionBegin;
245d4002b98SHong Zhang   for (i = 0; i < m; i++) {
24654c59aa7SJacob Faibussowitsch     if (idxm[i] < 0) continue; /* negative row */
24754c59aa7SJacob Faibussowitsch     PetscCheck(idxm[i] < mat->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, idxm[i], mat->rmap->N - 1);
248d4002b98SHong Zhang     if (idxm[i] >= rstart && idxm[i] < rend) {
249d4002b98SHong Zhang       row = idxm[i] - rstart;
250d4002b98SHong Zhang       for (j = 0; j < n; j++) {
25154c59aa7SJacob Faibussowitsch         if (idxn[j] < 0) continue; /* negative column */
25254c59aa7SJacob Faibussowitsch         PetscCheck(idxn[j] < mat->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT, idxn[j], mat->cmap->N - 1);
253d4002b98SHong Zhang         if (idxn[j] >= cstart && idxn[j] < cend) {
254d4002b98SHong Zhang           col = idxn[j] - cstart;
2559566063dSJacob Faibussowitsch           PetscCall(MatGetValues(sell->A, 1, &row, 1, &col, v + i * n + j));
256d4002b98SHong Zhang         } else {
257*48a46eb9SPierre Jolivet           if (!sell->colmap) PetscCall(MatCreateColmap_MPISELL_Private(mat));
258d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
2599566063dSJacob Faibussowitsch           PetscCall(PetscTableFind(sell->colmap, idxn[j] + 1, &col));
260d4002b98SHong Zhang           col--;
261d4002b98SHong Zhang #else
262d4002b98SHong Zhang           col = sell->colmap[idxn[j]] - 1;
263d4002b98SHong Zhang #endif
264d4002b98SHong Zhang           if ((col < 0) || (sell->garray[col] != idxn[j])) *(v + i * n + j) = 0.0;
265*48a46eb9SPierre Jolivet           else PetscCall(MatGetValues(sell->B, 1, &row, 1, &col, v + i * n + j));
266d4002b98SHong Zhang         }
267d4002b98SHong Zhang       }
268d4002b98SHong Zhang     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only local values currently supported");
269d4002b98SHong Zhang   }
270d4002b98SHong Zhang   PetscFunctionReturn(0);
271d4002b98SHong Zhang }
272d4002b98SHong Zhang 
273d4002b98SHong Zhang extern PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat, Vec, Vec);
274d4002b98SHong Zhang 
2759371c9d4SSatish Balay PetscErrorCode MatAssemblyBegin_MPISELL(Mat mat, MatAssemblyType mode) {
276d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
277d4002b98SHong Zhang   PetscInt     nstash, reallocs;
278d4002b98SHong Zhang 
279d4002b98SHong Zhang   PetscFunctionBegin;
280d4002b98SHong Zhang   if (sell->donotstash || mat->nooffprocentries) PetscFunctionReturn(0);
281d4002b98SHong Zhang 
2829566063dSJacob Faibussowitsch   PetscCall(MatStashScatterBegin_Private(mat, &mat->stash, mat->rmap->range));
2839566063dSJacob Faibussowitsch   PetscCall(MatStashGetInfo_Private(&mat->stash, &nstash, &reallocs));
2849566063dSJacob Faibussowitsch   PetscCall(PetscInfo(sell->A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
285d4002b98SHong Zhang   PetscFunctionReturn(0);
286d4002b98SHong Zhang }
287d4002b98SHong Zhang 
2889371c9d4SSatish Balay PetscErrorCode MatAssemblyEnd_MPISELL(Mat mat, MatAssemblyType mode) {
289d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
290d4002b98SHong Zhang   PetscMPIInt  n;
291d4002b98SHong Zhang   PetscInt     i, flg;
292d4002b98SHong Zhang   PetscInt    *row, *col;
293d4002b98SHong Zhang   PetscScalar *val;
294d4002b98SHong Zhang   PetscBool    other_disassembled;
295d4002b98SHong Zhang   /* do not use 'b = (Mat_SeqSELL*)sell->B->data' as B can be reset in disassembly */
296d4002b98SHong Zhang   PetscFunctionBegin;
297d4002b98SHong Zhang   if (!sell->donotstash && !mat->nooffprocentries) {
298d4002b98SHong Zhang     while (1) {
2999566063dSJacob Faibussowitsch       PetscCall(MatStashScatterGetMesg_Private(&mat->stash, &n, &row, &col, &val, &flg));
300d4002b98SHong Zhang       if (!flg) break;
301d4002b98SHong Zhang 
302d4002b98SHong Zhang       for (i = 0; i < n; i++) { /* assemble one by one */
3039566063dSJacob Faibussowitsch         PetscCall(MatSetValues_MPISELL(mat, 1, row + i, 1, col + i, val + i, mat->insertmode));
304d4002b98SHong Zhang       }
305d4002b98SHong Zhang     }
3069566063dSJacob Faibussowitsch     PetscCall(MatStashScatterEnd_Private(&mat->stash));
307d4002b98SHong Zhang   }
3089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sell->A, mode));
3099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sell->A, mode));
310d4002b98SHong Zhang 
311d4002b98SHong Zhang   /*
312d4002b98SHong Zhang      determine if any processor has disassembled, if so we must
3136aad120cSJose E. Roman      also disassemble ourselves, in order that we may reassemble.
314d4002b98SHong Zhang   */
315d4002b98SHong Zhang   /*
316d4002b98SHong Zhang      if nonzero structure of submatrix B cannot change then we know that
317d4002b98SHong Zhang      no processor disassembled thus we can skip this stuff
318d4002b98SHong Zhang   */
319d4002b98SHong Zhang   if (!((Mat_SeqSELL *)sell->B->data)->nonew) {
3205f9db2b2SJunchao Zhang     PetscCall(MPIU_Allreduce(&mat->was_assembled, &other_disassembled, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)mat)));
32108401ef6SPierre Jolivet     PetscCheck(!mat->was_assembled || other_disassembled, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatDisAssemble not implemented yet");
322d4002b98SHong Zhang   }
323*48a46eb9SPierre Jolivet   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) PetscCall(MatSetUpMultiply_MPISELL(mat));
324d4002b98SHong Zhang   /*
3259566063dSJacob Faibussowitsch   PetscCall(MatSetOption(sell->B,MAT_USE_INODES,PETSC_FALSE));
326d4002b98SHong Zhang   */
3279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sell->B, mode));
3289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sell->B, mode));
3299566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sell->rowvalues, sell->rowindices));
330f4259b30SLisandro Dalcin   sell->rowvalues = NULL;
3319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->diag));
332d4002b98SHong Zhang 
333d4002b98SHong Zhang   /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
334d4002b98SHong Zhang   if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqSELL *)(sell->A->data))->nonew) {
335d4002b98SHong Zhang     PetscObjectState state = sell->A->nonzerostate + sell->B->nonzerostate;
3361c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&state, &mat->nonzerostate, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)mat)));
337d4002b98SHong Zhang   }
338d4002b98SHong Zhang   PetscFunctionReturn(0);
339d4002b98SHong Zhang }
340d4002b98SHong Zhang 
3419371c9d4SSatish Balay PetscErrorCode MatZeroEntries_MPISELL(Mat A) {
342d4002b98SHong Zhang   Mat_MPISELL *l = (Mat_MPISELL *)A->data;
343d4002b98SHong Zhang 
344d4002b98SHong Zhang   PetscFunctionBegin;
3459566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->A));
3469566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->B));
347d4002b98SHong Zhang   PetscFunctionReturn(0);
348d4002b98SHong Zhang }
349d4002b98SHong Zhang 
3509371c9d4SSatish Balay PetscErrorCode MatMult_MPISELL(Mat A, Vec xx, Vec yy) {
351d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
352d4002b98SHong Zhang   PetscInt     nt;
353d4002b98SHong Zhang 
354d4002b98SHong Zhang   PetscFunctionBegin;
3559566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xx, &nt));
35608401ef6SPierre Jolivet   PetscCheck(nt == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible partition of A (%" PetscInt_FMT ") and xx (%" PetscInt_FMT ")", A->cmap->n, nt);
3579566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3589566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->mult)(a->A, xx, yy));
3599566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3609566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B, a->lvec, yy, yy));
361d4002b98SHong Zhang   PetscFunctionReturn(0);
362d4002b98SHong Zhang }
363d4002b98SHong Zhang 
3649371c9d4SSatish Balay PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat A, Vec bb, Vec xx) {
365d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
366d4002b98SHong Zhang 
367d4002b98SHong Zhang   PetscFunctionBegin;
3689566063dSJacob Faibussowitsch   PetscCall(MatMultDiagonalBlock(a->A, bb, xx));
369d4002b98SHong Zhang   PetscFunctionReturn(0);
370d4002b98SHong Zhang }
371d4002b98SHong Zhang 
3729371c9d4SSatish Balay PetscErrorCode MatMultAdd_MPISELL(Mat A, Vec xx, Vec yy, Vec zz) {
373d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
374d4002b98SHong Zhang 
375d4002b98SHong Zhang   PetscFunctionBegin;
3769566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3779566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multadd)(a->A, xx, yy, zz));
3789566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3799566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B, a->lvec, zz, zz));
380d4002b98SHong Zhang   PetscFunctionReturn(0);
381d4002b98SHong Zhang }
382d4002b98SHong Zhang 
3839371c9d4SSatish Balay PetscErrorCode MatMultTranspose_MPISELL(Mat A, Vec xx, Vec yy) {
384d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
385d4002b98SHong Zhang 
386d4002b98SHong Zhang   PetscFunctionBegin;
387d4002b98SHong Zhang   /* do nondiagonal part */
3889566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B, xx, a->lvec));
389d4002b98SHong Zhang   /* do local part */
3909566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multtranspose)(a->A, xx, yy));
391a29b4eb7SJunchao Zhang   /* add partial results together */
3929566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, a->lvec, yy, ADD_VALUES, SCATTER_REVERSE));
3939566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, a->lvec, yy, ADD_VALUES, SCATTER_REVERSE));
394d4002b98SHong Zhang   PetscFunctionReturn(0);
395d4002b98SHong Zhang }
396d4002b98SHong Zhang 
3979371c9d4SSatish Balay PetscErrorCode MatIsTranspose_MPISELL(Mat Amat, Mat Bmat, PetscReal tol, PetscBool *f) {
398d4002b98SHong Zhang   MPI_Comm     comm;
399d4002b98SHong Zhang   Mat_MPISELL *Asell = (Mat_MPISELL *)Amat->data, *Bsell;
400d4002b98SHong Zhang   Mat          Adia  = Asell->A, Bdia, Aoff, Boff, *Aoffs, *Boffs;
401d4002b98SHong Zhang   IS           Me, Notme;
402d4002b98SHong Zhang   PetscInt     M, N, first, last, *notme, i;
403d4002b98SHong Zhang   PetscMPIInt  size;
404d4002b98SHong Zhang 
405d4002b98SHong Zhang   PetscFunctionBegin;
406d4002b98SHong Zhang   /* Easy test: symmetric diagonal block */
4079371c9d4SSatish Balay   Bsell = (Mat_MPISELL *)Bmat->data;
4089371c9d4SSatish Balay   Bdia  = Bsell->A;
4099566063dSJacob Faibussowitsch   PetscCall(MatIsTranspose(Adia, Bdia, tol, f));
410d4002b98SHong Zhang   if (!*f) PetscFunctionReturn(0);
4119566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
4129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
413d4002b98SHong Zhang   if (size == 1) PetscFunctionReturn(0);
414d4002b98SHong Zhang 
415d4002b98SHong Zhang   /* Hard test: off-diagonal block. This takes a MatCreateSubMatrix. */
4169566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Amat, &M, &N));
4179566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &first, &last));
4189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N - last + first, &notme));
419d4002b98SHong Zhang   for (i = 0; i < first; i++) notme[i] = i;
420d4002b98SHong Zhang   for (i = last; i < M; i++) notme[i - last + first] = i;
4219566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(MPI_COMM_SELF, N - last + first, notme, PETSC_COPY_VALUES, &Notme));
4229566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(MPI_COMM_SELF, last - first, first, 1, &Me));
4239566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Amat, 1, &Me, &Notme, MAT_INITIAL_MATRIX, &Aoffs));
424d4002b98SHong Zhang   Aoff = Aoffs[0];
4259566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Bmat, 1, &Notme, &Me, MAT_INITIAL_MATRIX, &Boffs));
426d4002b98SHong Zhang   Boff = Boffs[0];
4279566063dSJacob Faibussowitsch   PetscCall(MatIsTranspose(Aoff, Boff, tol, f));
4289566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1, &Aoffs));
4299566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1, &Boffs));
4309566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Me));
4319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Notme));
4329566063dSJacob Faibussowitsch   PetscCall(PetscFree(notme));
433d4002b98SHong Zhang   PetscFunctionReturn(0);
434d4002b98SHong Zhang }
435d4002b98SHong Zhang 
4369371c9d4SSatish Balay PetscErrorCode MatMultTransposeAdd_MPISELL(Mat A, Vec xx, Vec yy, Vec zz) {
437d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
438d4002b98SHong Zhang 
439d4002b98SHong Zhang   PetscFunctionBegin;
440d4002b98SHong Zhang   /* do nondiagonal part */
4419566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B, xx, a->lvec));
442d4002b98SHong Zhang   /* do local part */
4439566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multtransposeadd)(a->A, xx, yy, zz));
444e4a140f6SJunchao Zhang   /* add partial results together */
4459566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, a->lvec, zz, ADD_VALUES, SCATTER_REVERSE));
4469566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, a->lvec, zz, ADD_VALUES, SCATTER_REVERSE));
447d4002b98SHong Zhang   PetscFunctionReturn(0);
448d4002b98SHong Zhang }
449d4002b98SHong Zhang 
450d4002b98SHong Zhang /*
451d4002b98SHong Zhang   This only works correctly for square matrices where the subblock A->A is the
452d4002b98SHong Zhang    diagonal block
453d4002b98SHong Zhang */
4549371c9d4SSatish Balay PetscErrorCode MatGetDiagonal_MPISELL(Mat A, Vec v) {
455d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
456d4002b98SHong Zhang 
457d4002b98SHong Zhang   PetscFunctionBegin;
45808401ef6SPierre Jolivet   PetscCheck(A->rmap->N == A->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Supports only square matrix where A->A is diag block");
459aed4548fSBarry Smith   PetscCheck(A->rmap->rstart == A->cmap->rstart && A->rmap->rend == A->cmap->rend, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "row partition must equal col partition");
4609566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(a->A, v));
461d4002b98SHong Zhang   PetscFunctionReturn(0);
462d4002b98SHong Zhang }
463d4002b98SHong Zhang 
4649371c9d4SSatish Balay PetscErrorCode MatScale_MPISELL(Mat A, PetscScalar aa) {
465d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
466d4002b98SHong Zhang 
467d4002b98SHong Zhang   PetscFunctionBegin;
4689566063dSJacob Faibussowitsch   PetscCall(MatScale(a->A, aa));
4699566063dSJacob Faibussowitsch   PetscCall(MatScale(a->B, aa));
470d4002b98SHong Zhang   PetscFunctionReturn(0);
471d4002b98SHong Zhang }
472d4002b98SHong Zhang 
4739371c9d4SSatish Balay PetscErrorCode MatDestroy_MPISELL(Mat mat) {
474d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
475d4002b98SHong Zhang 
476d4002b98SHong Zhang   PetscFunctionBegin;
477d4002b98SHong Zhang #if defined(PETSC_USE_LOG)
478c0aa6a63SJacob Faibussowitsch   PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT, mat->rmap->N, mat->cmap->N);
479d4002b98SHong Zhang #endif
4809566063dSJacob Faibussowitsch   PetscCall(MatStashDestroy_Private(&mat->stash));
4819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->diag));
4829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sell->A));
4839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sell->B));
484d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
4859566063dSJacob Faibussowitsch   PetscCall(PetscTableDestroy(&sell->colmap));
486d4002b98SHong Zhang #else
4879566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->colmap));
488d4002b98SHong Zhang #endif
4899566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->garray));
4909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->lvec));
4919566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sell->Mvctx));
4929566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sell->rowvalues, sell->rowindices));
4939566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->ld));
4949566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
495d4002b98SHong Zhang 
4969566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL));
4979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatStoreValues_C", NULL));
4989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatRetrieveValues_C", NULL));
4999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatIsTranspose_C", NULL));
5009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPISELLSetPreallocation_C", NULL));
5019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisell_mpiaij_C", NULL));
5029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalScaleLocal_C", NULL));
503d4002b98SHong Zhang   PetscFunctionReturn(0);
504d4002b98SHong Zhang }
505d4002b98SHong Zhang 
506d4002b98SHong Zhang #include <petscdraw.h>
5079371c9d4SSatish Balay PetscErrorCode MatView_MPISELL_ASCIIorDraworSocket(Mat mat, PetscViewer viewer) {
508d4002b98SHong Zhang   Mat_MPISELL      *sell = (Mat_MPISELL *)mat->data;
509d4002b98SHong Zhang   PetscMPIInt       rank = sell->rank, size = sell->size;
510d4002b98SHong Zhang   PetscBool         isdraw, iascii, isbinary;
511d4002b98SHong Zhang   PetscViewer       sviewer;
512d4002b98SHong Zhang   PetscViewerFormat format;
513d4002b98SHong Zhang 
514d4002b98SHong Zhang   PetscFunctionBegin;
5159566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
5169566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
5179566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
518d4002b98SHong Zhang   if (iascii) {
5199566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
520d4002b98SHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
521d4002b98SHong Zhang       MatInfo   info;
5226335e310SSatish Balay       PetscInt *inodes;
523d4002b98SHong Zhang 
5249566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat), &rank));
5259566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(mat, MAT_LOCAL, &info));
5269566063dSJacob Faibussowitsch       PetscCall(MatInodeGetInodeSizes(sell->A, NULL, &inodes, NULL));
5279566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
528d4002b98SHong Zhang       if (!inodes) {
5299371c9d4SSatish Balay         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Local rows %" PetscInt_FMT " nz %" PetscInt_FMT " nz alloced %" PetscInt_FMT " mem %" PetscInt_FMT ", not using I-node routines\n", rank, mat->rmap->n, (PetscInt)info.nz_used,
5309371c9d4SSatish Balay                                                      (PetscInt)info.nz_allocated, (PetscInt)info.memory));
531d4002b98SHong Zhang       } else {
5329371c9d4SSatish Balay         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Local rows %" PetscInt_FMT " nz %" PetscInt_FMT " nz alloced %" PetscInt_FMT " mem %" PetscInt_FMT ", using I-node routines\n", rank, mat->rmap->n, (PetscInt)info.nz_used,
5339371c9d4SSatish Balay                                                      (PetscInt)info.nz_allocated, (PetscInt)info.memory));
534d4002b98SHong Zhang       }
5359566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(sell->A, MAT_LOCAL, &info));
5369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] on-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
5379566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(sell->B, MAT_LOCAL, &info));
5389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] off-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
5399566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
5409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
5419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Information on VecScatter used in matrix-vector product: \n"));
5429566063dSJacob Faibussowitsch       PetscCall(VecScatterView(sell->Mvctx, viewer));
543d4002b98SHong Zhang       PetscFunctionReturn(0);
544d4002b98SHong Zhang     } else if (format == PETSC_VIEWER_ASCII_INFO) {
545d4002b98SHong Zhang       PetscInt inodecount, inodelimit, *inodes;
5469566063dSJacob Faibussowitsch       PetscCall(MatInodeGetInodeSizes(sell->A, &inodecount, &inodes, &inodelimit));
547d4002b98SHong Zhang       if (inodes) {
5489566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "using I-node (on process 0) routines: found %" PetscInt_FMT " nodes, limit used is %" PetscInt_FMT "\n", inodecount, inodelimit));
549d4002b98SHong Zhang       } else {
5509566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "not using I-node (on process 0) routines\n"));
551d4002b98SHong Zhang       }
552d4002b98SHong Zhang       PetscFunctionReturn(0);
553d4002b98SHong Zhang     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
554d4002b98SHong Zhang       PetscFunctionReturn(0);
555d4002b98SHong Zhang     }
556d4002b98SHong Zhang   } else if (isbinary) {
557d4002b98SHong Zhang     if (size == 1) {
5589566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)sell->A, ((PetscObject)mat)->name));
5599566063dSJacob Faibussowitsch       PetscCall(MatView(sell->A, viewer));
560d4002b98SHong Zhang     } else {
5619566063dSJacob Faibussowitsch       /* PetscCall(MatView_MPISELL_Binary(mat,viewer)); */
562d4002b98SHong Zhang     }
563d4002b98SHong Zhang     PetscFunctionReturn(0);
564d4002b98SHong Zhang   } else if (isdraw) {
565d4002b98SHong Zhang     PetscDraw draw;
566d4002b98SHong Zhang     PetscBool isnull;
5679566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
5689566063dSJacob Faibussowitsch     PetscCall(PetscDrawIsNull(draw, &isnull));
569d4002b98SHong Zhang     if (isnull) PetscFunctionReturn(0);
570d4002b98SHong Zhang   }
571d4002b98SHong Zhang 
572d4002b98SHong Zhang   {
573d4002b98SHong Zhang     /* assemble the entire matrix onto first processor. */
574d4002b98SHong Zhang     Mat          A;
575d4002b98SHong Zhang     Mat_SeqSELL *Aloc;
576d4002b98SHong Zhang     PetscInt     M = mat->rmap->N, N = mat->cmap->N, *acolidx, row, col, i, j;
577d4002b98SHong Zhang     MatScalar   *aval;
578d4002b98SHong Zhang     PetscBool    isnonzero;
579d4002b98SHong Zhang 
5809566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &A));
581dd400576SPatrick Sanan     if (rank == 0) {
5829566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A, M, N, M, N));
583d4002b98SHong Zhang     } else {
5849566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A, 0, 0, M, N));
585d4002b98SHong Zhang     }
586d4002b98SHong Zhang     /* This is just a temporary matrix, so explicitly using MATMPISELL is probably best */
5879566063dSJacob Faibussowitsch     PetscCall(MatSetType(A, MATMPISELL));
5889566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(A, 0, NULL, 0, NULL));
5899566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
5909566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)mat, (PetscObject)A));
591d4002b98SHong Zhang 
592d4002b98SHong Zhang     /* copy over the A part */
593d4002b98SHong Zhang     Aloc    = (Mat_SeqSELL *)sell->A->data;
5949371c9d4SSatish Balay     acolidx = Aloc->colidx;
5959371c9d4SSatish Balay     aval    = Aloc->val;
596d4002b98SHong Zhang     for (i = 0; i < Aloc->totalslices; i++) { /* loop over slices */
597d4002b98SHong Zhang       for (j = Aloc->sliidx[i]; j < Aloc->sliidx[i + 1]; j++) {
598d4002b98SHong Zhang         isnonzero = (PetscBool)((j - Aloc->sliidx[i]) / 8 < Aloc->rlen[(i << 3) + (j & 0x07)]);
599d4002b98SHong Zhang         if (isnonzero) {                                   /* check the mask bit */
600d4002b98SHong Zhang           row = (i << 3) + (j & 0x07) + mat->rmap->rstart; /* i<<3 is the starting row of this slice */
601d4002b98SHong Zhang           col = *acolidx + mat->rmap->rstart;
6029566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A, 1, &row, 1, &col, aval, INSERT_VALUES));
603d4002b98SHong Zhang         }
6049371c9d4SSatish Balay         aval++;
6059371c9d4SSatish Balay         acolidx++;
606d4002b98SHong Zhang       }
607d4002b98SHong Zhang     }
608d4002b98SHong Zhang 
609d4002b98SHong Zhang     /* copy over the B part */
610d4002b98SHong Zhang     Aloc    = (Mat_SeqSELL *)sell->B->data;
6119371c9d4SSatish Balay     acolidx = Aloc->colidx;
6129371c9d4SSatish Balay     aval    = Aloc->val;
613d4002b98SHong Zhang     for (i = 0; i < Aloc->totalslices; i++) {
614d4002b98SHong Zhang       for (j = Aloc->sliidx[i]; j < Aloc->sliidx[i + 1]; j++) {
615d4002b98SHong Zhang         isnonzero = (PetscBool)((j - Aloc->sliidx[i]) / 8 < Aloc->rlen[(i << 3) + (j & 0x07)]);
616d4002b98SHong Zhang         if (isnonzero) {
617d4002b98SHong Zhang           row = (i << 3) + (j & 0x07) + mat->rmap->rstart;
618d4002b98SHong Zhang           col = sell->garray[*acolidx];
6199566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A, 1, &row, 1, &col, aval, INSERT_VALUES));
620d4002b98SHong Zhang         }
6219371c9d4SSatish Balay         aval++;
6229371c9d4SSatish Balay         acolidx++;
623d4002b98SHong Zhang       }
624d4002b98SHong Zhang     }
625d4002b98SHong Zhang 
6269566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
6279566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
628d4002b98SHong Zhang     /*
629d4002b98SHong Zhang        Everyone has to call to draw the matrix since the graphics waits are
630d4002b98SHong Zhang        synchronized across all processors that share the PetscDraw object
631d4002b98SHong Zhang     */
6329566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
633dd400576SPatrick Sanan     if (rank == 0) {
6349566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)((Mat_MPISELL *)(A->data))->A, ((PetscObject)mat)->name));
6359566063dSJacob Faibussowitsch       PetscCall(MatView_SeqSELL(((Mat_MPISELL *)(A->data))->A, sviewer));
636d4002b98SHong Zhang     }
6379566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
6389566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
6399566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
640d4002b98SHong Zhang   }
641d4002b98SHong Zhang   PetscFunctionReturn(0);
642d4002b98SHong Zhang }
643d4002b98SHong Zhang 
6449371c9d4SSatish Balay PetscErrorCode MatView_MPISELL(Mat mat, PetscViewer viewer) {
645d4002b98SHong Zhang   PetscBool iascii, isdraw, issocket, isbinary;
646d4002b98SHong Zhang 
647d4002b98SHong Zhang   PetscFunctionBegin;
6489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
6499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
6509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
6519566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSOCKET, &issocket));
652*48a46eb9SPierre Jolivet   if (iascii || isdraw || isbinary || issocket) PetscCall(MatView_MPISELL_ASCIIorDraworSocket(mat, viewer));
653d4002b98SHong Zhang   PetscFunctionReturn(0);
654d4002b98SHong Zhang }
655d4002b98SHong Zhang 
6569371c9d4SSatish Balay PetscErrorCode MatGetGhosts_MPISELL(Mat mat, PetscInt *nghosts, const PetscInt *ghosts[]) {
657d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
658d4002b98SHong Zhang 
659d4002b98SHong Zhang   PetscFunctionBegin;
6609566063dSJacob Faibussowitsch   PetscCall(MatGetSize(sell->B, NULL, nghosts));
661d4002b98SHong Zhang   if (ghosts) *ghosts = sell->garray;
662d4002b98SHong Zhang   PetscFunctionReturn(0);
663d4002b98SHong Zhang }
664d4002b98SHong Zhang 
6659371c9d4SSatish Balay PetscErrorCode MatGetInfo_MPISELL(Mat matin, MatInfoType flag, MatInfo *info) {
666d4002b98SHong Zhang   Mat_MPISELL   *mat = (Mat_MPISELL *)matin->data;
667d4002b98SHong Zhang   Mat            A = mat->A, B = mat->B;
6683966268fSBarry Smith   PetscLogDouble isend[5], irecv[5];
669d4002b98SHong Zhang 
670d4002b98SHong Zhang   PetscFunctionBegin;
671d4002b98SHong Zhang   info->block_size = 1.0;
6729566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_LOCAL, info));
673d4002b98SHong Zhang 
6749371c9d4SSatish Balay   isend[0] = info->nz_used;
6759371c9d4SSatish Balay   isend[1] = info->nz_allocated;
6769371c9d4SSatish Balay   isend[2] = info->nz_unneeded;
6779371c9d4SSatish Balay   isend[3] = info->memory;
6789371c9d4SSatish Balay   isend[4] = info->mallocs;
679d4002b98SHong Zhang 
6809566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(B, MAT_LOCAL, info));
681d4002b98SHong Zhang 
6829371c9d4SSatish Balay   isend[0] += info->nz_used;
6839371c9d4SSatish Balay   isend[1] += info->nz_allocated;
6849371c9d4SSatish Balay   isend[2] += info->nz_unneeded;
6859371c9d4SSatish Balay   isend[3] += info->memory;
6869371c9d4SSatish Balay   isend[4] += info->mallocs;
687d4002b98SHong Zhang   if (flag == MAT_LOCAL) {
688d4002b98SHong Zhang     info->nz_used      = isend[0];
689d4002b98SHong Zhang     info->nz_allocated = isend[1];
690d4002b98SHong Zhang     info->nz_unneeded  = isend[2];
691d4002b98SHong Zhang     info->memory       = isend[3];
692d4002b98SHong Zhang     info->mallocs      = isend[4];
693d4002b98SHong Zhang   } else if (flag == MAT_GLOBAL_MAX) {
6941c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_MAX, PetscObjectComm((PetscObject)matin)));
695d4002b98SHong Zhang 
696d4002b98SHong Zhang     info->nz_used      = irecv[0];
697d4002b98SHong Zhang     info->nz_allocated = irecv[1];
698d4002b98SHong Zhang     info->nz_unneeded  = irecv[2];
699d4002b98SHong Zhang     info->memory       = irecv[3];
700d4002b98SHong Zhang     info->mallocs      = irecv[4];
701d4002b98SHong Zhang   } else if (flag == MAT_GLOBAL_SUM) {
7021c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_SUM, PetscObjectComm((PetscObject)matin)));
703d4002b98SHong Zhang 
704d4002b98SHong Zhang     info->nz_used      = irecv[0];
705d4002b98SHong Zhang     info->nz_allocated = irecv[1];
706d4002b98SHong Zhang     info->nz_unneeded  = irecv[2];
707d4002b98SHong Zhang     info->memory       = irecv[3];
708d4002b98SHong Zhang     info->mallocs      = irecv[4];
709d4002b98SHong Zhang   }
710d4002b98SHong Zhang   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
711d4002b98SHong Zhang   info->fill_ratio_needed = 0;
712d4002b98SHong Zhang   info->factor_mallocs    = 0;
713d4002b98SHong Zhang   PetscFunctionReturn(0);
714d4002b98SHong Zhang }
715d4002b98SHong Zhang 
7169371c9d4SSatish Balay PetscErrorCode MatSetOption_MPISELL(Mat A, MatOption op, PetscBool flg) {
717d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
718d4002b98SHong Zhang 
719d4002b98SHong Zhang   PetscFunctionBegin;
720d4002b98SHong Zhang   switch (op) {
721d4002b98SHong Zhang   case MAT_NEW_NONZERO_LOCATIONS:
722d4002b98SHong Zhang   case MAT_NEW_NONZERO_ALLOCATION_ERR:
723d4002b98SHong Zhang   case MAT_UNUSED_NONZERO_LOCATION_ERR:
724d4002b98SHong Zhang   case MAT_KEEP_NONZERO_PATTERN:
725d4002b98SHong Zhang   case MAT_NEW_NONZERO_LOCATION_ERR:
726d4002b98SHong Zhang   case MAT_USE_INODES:
727d4002b98SHong Zhang   case MAT_IGNORE_ZERO_ENTRIES:
728d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7299566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
7309566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B, op, flg));
731d4002b98SHong Zhang     break;
732d4002b98SHong Zhang   case MAT_ROW_ORIENTED:
733d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
734d4002b98SHong Zhang     a->roworiented = flg;
735d4002b98SHong Zhang 
7369566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
7379566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B, op, flg));
738d4002b98SHong Zhang     break;
7398c78258cSHong Zhang   case MAT_FORCE_DIAGONAL_ENTRIES:
7409371c9d4SSatish Balay   case MAT_SORTED_FULL: PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op])); break;
7419371c9d4SSatish Balay   case MAT_IGNORE_OFF_PROC_ENTRIES: a->donotstash = flg; break;
742d4002b98SHong Zhang   case MAT_SPD:
7439371c9d4SSatish Balay   case MAT_SPD_ETERNAL: break;
744d4002b98SHong Zhang   case MAT_SYMMETRIC:
745d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7469566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
747d4002b98SHong Zhang     break;
748d4002b98SHong Zhang   case MAT_STRUCTURALLY_SYMMETRIC:
749d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7509566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
751d4002b98SHong Zhang     break;
752d4002b98SHong Zhang   case MAT_HERMITIAN:
753d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7549566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
755d4002b98SHong Zhang     break;
756d4002b98SHong Zhang   case MAT_SYMMETRY_ETERNAL:
757d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7589566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
759d4002b98SHong Zhang     break;
760b94d7dedSBarry Smith   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
761b94d7dedSBarry Smith     MatCheckPreallocated(A, 1);
762b94d7dedSBarry Smith     PetscCall(MatSetOption(a->A, op, flg));
763b94d7dedSBarry Smith     break;
7649371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
765d4002b98SHong Zhang   }
766d4002b98SHong Zhang   PetscFunctionReturn(0);
767d4002b98SHong Zhang }
768d4002b98SHong Zhang 
7699371c9d4SSatish Balay PetscErrorCode MatDiagonalScale_MPISELL(Mat mat, Vec ll, Vec rr) {
770d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
771d4002b98SHong Zhang   Mat          a = sell->A, b = sell->B;
772d4002b98SHong Zhang   PetscInt     s1, s2, s3;
773d4002b98SHong Zhang 
774d4002b98SHong Zhang   PetscFunctionBegin;
7759566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mat, &s2, &s3));
776d4002b98SHong Zhang   if (rr) {
7779566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rr, &s1));
77808401ef6SPierre Jolivet     PetscCheck(s1 == s3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "right vector non-conforming local size");
779d4002b98SHong Zhang     /* Overlap communication with computation. */
7809566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(sell->Mvctx, rr, sell->lvec, INSERT_VALUES, SCATTER_FORWARD));
781d4002b98SHong Zhang   }
782d4002b98SHong Zhang   if (ll) {
7839566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(ll, &s1));
78408401ef6SPierre Jolivet     PetscCheck(s1 == s2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "left vector non-conforming local size");
785dbbe0bcdSBarry Smith     PetscUseTypeMethod(b, diagonalscale, ll, NULL);
786d4002b98SHong Zhang   }
787d4002b98SHong Zhang   /* scale  the diagonal block */
788dbbe0bcdSBarry Smith   PetscUseTypeMethod(a, diagonalscale, ll, rr);
789d4002b98SHong Zhang 
790d4002b98SHong Zhang   if (rr) {
791d4002b98SHong Zhang     /* Do a scatter end and then right scale the off-diagonal block */
7929566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(sell->Mvctx, rr, sell->lvec, INSERT_VALUES, SCATTER_FORWARD));
793dbbe0bcdSBarry Smith     PetscUseTypeMethod(b, diagonalscale, NULL, sell->lvec);
794d4002b98SHong Zhang   }
795d4002b98SHong Zhang   PetscFunctionReturn(0);
796d4002b98SHong Zhang }
797d4002b98SHong Zhang 
7989371c9d4SSatish Balay PetscErrorCode MatSetUnfactored_MPISELL(Mat A) {
799d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
800d4002b98SHong Zhang 
801d4002b98SHong Zhang   PetscFunctionBegin;
8029566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(a->A));
803d4002b98SHong Zhang   PetscFunctionReturn(0);
804d4002b98SHong Zhang }
805d4002b98SHong Zhang 
8069371c9d4SSatish Balay PetscErrorCode MatEqual_MPISELL(Mat A, Mat B, PetscBool *flag) {
807d4002b98SHong Zhang   Mat_MPISELL *matB = (Mat_MPISELL *)B->data, *matA = (Mat_MPISELL *)A->data;
808d4002b98SHong Zhang   Mat          a, b, c, d;
809d4002b98SHong Zhang   PetscBool    flg;
810d4002b98SHong Zhang 
811d4002b98SHong Zhang   PetscFunctionBegin;
8129371c9d4SSatish Balay   a = matA->A;
8139371c9d4SSatish Balay   b = matA->B;
8149371c9d4SSatish Balay   c = matB->A;
8159371c9d4SSatish Balay   d = matB->B;
816d4002b98SHong Zhang 
8179566063dSJacob Faibussowitsch   PetscCall(MatEqual(a, c, &flg));
818*48a46eb9SPierre Jolivet   if (flg) PetscCall(MatEqual(b, d, &flg));
8191c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&flg, flag, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
820d4002b98SHong Zhang   PetscFunctionReturn(0);
821d4002b98SHong Zhang }
822d4002b98SHong Zhang 
8239371c9d4SSatish Balay PetscErrorCode MatCopy_MPISELL(Mat A, Mat B, MatStructure str) {
824d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
825d4002b98SHong Zhang   Mat_MPISELL *b = (Mat_MPISELL *)B->data;
826d4002b98SHong Zhang 
827d4002b98SHong Zhang   PetscFunctionBegin;
828d4002b98SHong Zhang   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
829d4002b98SHong Zhang   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
830d4002b98SHong Zhang     /* because of the column compression in the off-processor part of the matrix a->B,
831d4002b98SHong Zhang        the number of columns in a->B and b->B may be different, hence we cannot call
832d4002b98SHong Zhang        the MatCopy() directly on the two parts. If need be, we can provide a more
833d4002b98SHong Zhang        efficient copy than the MatCopy_Basic() by first uncompressing the a->B matrices
834d4002b98SHong Zhang        then copying the submatrices */
8359566063dSJacob Faibussowitsch     PetscCall(MatCopy_Basic(A, B, str));
836d4002b98SHong Zhang   } else {
8379566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->A, b->A, str));
8389566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->B, b->B, str));
839d4002b98SHong Zhang   }
840d4002b98SHong Zhang   PetscFunctionReturn(0);
841d4002b98SHong Zhang }
842d4002b98SHong Zhang 
8439371c9d4SSatish Balay PetscErrorCode MatSetUp_MPISELL(Mat A) {
844d4002b98SHong Zhang   PetscFunctionBegin;
8459566063dSJacob Faibussowitsch   PetscCall(MatMPISELLSetPreallocation(A, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
846d4002b98SHong Zhang   PetscFunctionReturn(0);
847d4002b98SHong Zhang }
848d4002b98SHong Zhang 
849d4002b98SHong Zhang extern PetscErrorCode MatConjugate_SeqSELL(Mat);
850d4002b98SHong Zhang 
8519371c9d4SSatish Balay PetscErrorCode MatConjugate_MPISELL(Mat mat) {
8525f80ce2aSJacob Faibussowitsch   PetscFunctionBegin;
8535f80ce2aSJacob Faibussowitsch   if (PetscDefined(USE_COMPLEX)) {
854d4002b98SHong Zhang     Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
855d4002b98SHong Zhang 
8569566063dSJacob Faibussowitsch     PetscCall(MatConjugate_SeqSELL(sell->A));
8579566063dSJacob Faibussowitsch     PetscCall(MatConjugate_SeqSELL(sell->B));
8585f80ce2aSJacob Faibussowitsch   }
859d4002b98SHong Zhang   PetscFunctionReturn(0);
860d4002b98SHong Zhang }
861d4002b98SHong Zhang 
8629371c9d4SSatish Balay PetscErrorCode MatRealPart_MPISELL(Mat A) {
863d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
864d4002b98SHong Zhang 
865d4002b98SHong Zhang   PetscFunctionBegin;
8669566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->A));
8679566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->B));
868d4002b98SHong Zhang   PetscFunctionReturn(0);
869d4002b98SHong Zhang }
870d4002b98SHong Zhang 
8719371c9d4SSatish Balay PetscErrorCode MatImaginaryPart_MPISELL(Mat A) {
872d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
873d4002b98SHong Zhang 
874d4002b98SHong Zhang   PetscFunctionBegin;
8759566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->A));
8769566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->B));
877d4002b98SHong Zhang   PetscFunctionReturn(0);
878d4002b98SHong Zhang }
879d4002b98SHong Zhang 
8809371c9d4SSatish Balay PetscErrorCode MatInvertBlockDiagonal_MPISELL(Mat A, const PetscScalar **values) {
881d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
882d4002b98SHong Zhang 
883d4002b98SHong Zhang   PetscFunctionBegin;
8849566063dSJacob Faibussowitsch   PetscCall(MatInvertBlockDiagonal(a->A, values));
885d4002b98SHong Zhang   A->factorerrortype = a->A->factorerrortype;
886d4002b98SHong Zhang   PetscFunctionReturn(0);
887d4002b98SHong Zhang }
888d4002b98SHong Zhang 
8899371c9d4SSatish Balay static PetscErrorCode MatSetRandom_MPISELL(Mat x, PetscRandom rctx) {
890d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)x->data;
891d4002b98SHong Zhang 
892d4002b98SHong Zhang   PetscFunctionBegin;
8939566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(sell->A, rctx));
8949566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(sell->B, rctx));
8959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
8969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
897d4002b98SHong Zhang   PetscFunctionReturn(0);
898d4002b98SHong Zhang }
899d4002b98SHong Zhang 
9009371c9d4SSatish Balay PetscErrorCode MatSetFromOptions_MPISELL(Mat A, PetscOptionItems *PetscOptionsObject) {
901d4002b98SHong Zhang   PetscFunctionBegin;
902d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "MPISELL options");
903d0609cedSBarry Smith   PetscOptionsHeadEnd();
904d4002b98SHong Zhang   PetscFunctionReturn(0);
905d4002b98SHong Zhang }
906d4002b98SHong Zhang 
9079371c9d4SSatish Balay PetscErrorCode MatShift_MPISELL(Mat Y, PetscScalar a) {
908d4002b98SHong Zhang   Mat_MPISELL *msell = (Mat_MPISELL *)Y->data;
909d4002b98SHong Zhang   Mat_SeqSELL *sell  = (Mat_SeqSELL *)msell->A->data;
910d4002b98SHong Zhang 
911d4002b98SHong Zhang   PetscFunctionBegin;
912d4002b98SHong Zhang   if (!Y->preallocated) {
9139566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(Y, 1, NULL, 0, NULL));
914d4002b98SHong Zhang   } else if (!sell->nz) {
915d4002b98SHong Zhang     PetscInt nonew = sell->nonew;
9169566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(msell->A, 1, NULL));
917d4002b98SHong Zhang     sell->nonew = nonew;
918d4002b98SHong Zhang   }
9199566063dSJacob Faibussowitsch   PetscCall(MatShift_Basic(Y, a));
920d4002b98SHong Zhang   PetscFunctionReturn(0);
921d4002b98SHong Zhang }
922d4002b98SHong Zhang 
9239371c9d4SSatish Balay PetscErrorCode MatMissingDiagonal_MPISELL(Mat A, PetscBool *missing, PetscInt *d) {
924d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
925d4002b98SHong Zhang 
926d4002b98SHong Zhang   PetscFunctionBegin;
92708401ef6SPierre Jolivet   PetscCheck(A->rmap->n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only works for square matrices");
9289566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(a->A, missing, d));
929d4002b98SHong Zhang   if (d) {
930d4002b98SHong Zhang     PetscInt rstart;
9319566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
932d4002b98SHong Zhang     *d += rstart;
933d4002b98SHong Zhang   }
934d4002b98SHong Zhang   PetscFunctionReturn(0);
935d4002b98SHong Zhang }
936d4002b98SHong Zhang 
9379371c9d4SSatish Balay PetscErrorCode MatGetDiagonalBlock_MPISELL(Mat A, Mat *a) {
938d4002b98SHong Zhang   PetscFunctionBegin;
939d4002b98SHong Zhang   *a = ((Mat_MPISELL *)A->data)->A;
940d4002b98SHong Zhang   PetscFunctionReturn(0);
941d4002b98SHong Zhang }
942d4002b98SHong Zhang 
943d4002b98SHong Zhang /* -------------------------------------------------------------------*/
944d4002b98SHong Zhang static struct _MatOps MatOps_Values = {MatSetValues_MPISELL,
945f4259b30SLisandro Dalcin                                        NULL,
946f4259b30SLisandro Dalcin                                        NULL,
947d4002b98SHong Zhang                                        MatMult_MPISELL,
948d4002b98SHong Zhang                                        /* 4*/ MatMultAdd_MPISELL,
949d4002b98SHong Zhang                                        MatMultTranspose_MPISELL,
950d4002b98SHong Zhang                                        MatMultTransposeAdd_MPISELL,
951f4259b30SLisandro Dalcin                                        NULL,
952f4259b30SLisandro Dalcin                                        NULL,
953f4259b30SLisandro Dalcin                                        NULL,
954f4259b30SLisandro Dalcin                                        /*10*/ NULL,
955f4259b30SLisandro Dalcin                                        NULL,
956f4259b30SLisandro Dalcin                                        NULL,
957d4002b98SHong Zhang                                        MatSOR_MPISELL,
958f4259b30SLisandro Dalcin                                        NULL,
959d4002b98SHong Zhang                                        /*15*/ MatGetInfo_MPISELL,
960d4002b98SHong Zhang                                        MatEqual_MPISELL,
961d4002b98SHong Zhang                                        MatGetDiagonal_MPISELL,
962d4002b98SHong Zhang                                        MatDiagonalScale_MPISELL,
963f4259b30SLisandro Dalcin                                        NULL,
964d4002b98SHong Zhang                                        /*20*/ MatAssemblyBegin_MPISELL,
965d4002b98SHong Zhang                                        MatAssemblyEnd_MPISELL,
966d4002b98SHong Zhang                                        MatSetOption_MPISELL,
967d4002b98SHong Zhang                                        MatZeroEntries_MPISELL,
968f4259b30SLisandro Dalcin                                        /*24*/ NULL,
969f4259b30SLisandro Dalcin                                        NULL,
970f4259b30SLisandro Dalcin                                        NULL,
971f4259b30SLisandro Dalcin                                        NULL,
972f4259b30SLisandro Dalcin                                        NULL,
973d4002b98SHong Zhang                                        /*29*/ MatSetUp_MPISELL,
974f4259b30SLisandro Dalcin                                        NULL,
975f4259b30SLisandro Dalcin                                        NULL,
976d4002b98SHong Zhang                                        MatGetDiagonalBlock_MPISELL,
977f4259b30SLisandro Dalcin                                        NULL,
978d4002b98SHong Zhang                                        /*34*/ MatDuplicate_MPISELL,
979f4259b30SLisandro Dalcin                                        NULL,
980f4259b30SLisandro Dalcin                                        NULL,
981f4259b30SLisandro Dalcin                                        NULL,
982f4259b30SLisandro Dalcin                                        NULL,
983f4259b30SLisandro Dalcin                                        /*39*/ NULL,
984f4259b30SLisandro Dalcin                                        NULL,
985f4259b30SLisandro Dalcin                                        NULL,
986d4002b98SHong Zhang                                        MatGetValues_MPISELL,
987d4002b98SHong Zhang                                        MatCopy_MPISELL,
988f4259b30SLisandro Dalcin                                        /*44*/ NULL,
989d4002b98SHong Zhang                                        MatScale_MPISELL,
990d4002b98SHong Zhang                                        MatShift_MPISELL,
991d4002b98SHong Zhang                                        MatDiagonalSet_MPISELL,
992f4259b30SLisandro Dalcin                                        NULL,
993d4002b98SHong Zhang                                        /*49*/ MatSetRandom_MPISELL,
994f4259b30SLisandro Dalcin                                        NULL,
995f4259b30SLisandro Dalcin                                        NULL,
996f4259b30SLisandro Dalcin                                        NULL,
997f4259b30SLisandro Dalcin                                        NULL,
998d4002b98SHong Zhang                                        /*54*/ MatFDColoringCreate_MPIXAIJ,
999f4259b30SLisandro Dalcin                                        NULL,
1000d4002b98SHong Zhang                                        MatSetUnfactored_MPISELL,
1001f4259b30SLisandro Dalcin                                        NULL,
1002f4259b30SLisandro Dalcin                                        NULL,
1003f4259b30SLisandro Dalcin                                        /*59*/ NULL,
1004d4002b98SHong Zhang                                        MatDestroy_MPISELL,
1005d4002b98SHong Zhang                                        MatView_MPISELL,
1006f4259b30SLisandro Dalcin                                        NULL,
1007f4259b30SLisandro Dalcin                                        NULL,
1008f4259b30SLisandro Dalcin                                        /*64*/ NULL,
1009f4259b30SLisandro Dalcin                                        NULL,
1010f4259b30SLisandro Dalcin                                        NULL,
1011f4259b30SLisandro Dalcin                                        NULL,
1012f4259b30SLisandro Dalcin                                        NULL,
1013f4259b30SLisandro Dalcin                                        /*69*/ NULL,
1014f4259b30SLisandro Dalcin                                        NULL,
1015f4259b30SLisandro Dalcin                                        NULL,
1016f4259b30SLisandro Dalcin                                        NULL,
1017f4259b30SLisandro Dalcin                                        NULL,
1018f4259b30SLisandro Dalcin                                        NULL,
1019d4002b98SHong Zhang                                        /*75*/ MatFDColoringApply_AIJ, /* reuse AIJ function */
1020d4002b98SHong Zhang                                        MatSetFromOptions_MPISELL,
1021f4259b30SLisandro Dalcin                                        NULL,
1022f4259b30SLisandro Dalcin                                        NULL,
1023f4259b30SLisandro Dalcin                                        NULL,
1024f4259b30SLisandro Dalcin                                        /*80*/ NULL,
1025f4259b30SLisandro Dalcin                                        NULL,
1026f4259b30SLisandro Dalcin                                        NULL,
1027f4259b30SLisandro Dalcin                                        /*83*/ NULL,
1028f4259b30SLisandro Dalcin                                        NULL,
1029f4259b30SLisandro Dalcin                                        NULL,
1030f4259b30SLisandro Dalcin                                        NULL,
1031f4259b30SLisandro Dalcin                                        NULL,
1032f4259b30SLisandro Dalcin                                        NULL,
1033f4259b30SLisandro Dalcin                                        /*89*/ NULL,
1034f4259b30SLisandro Dalcin                                        NULL,
1035f4259b30SLisandro Dalcin                                        NULL,
1036f4259b30SLisandro Dalcin                                        NULL,
1037f4259b30SLisandro Dalcin                                        NULL,
1038f4259b30SLisandro Dalcin                                        /*94*/ NULL,
1039f4259b30SLisandro Dalcin                                        NULL,
1040f4259b30SLisandro Dalcin                                        NULL,
1041f4259b30SLisandro Dalcin                                        NULL,
1042f4259b30SLisandro Dalcin                                        NULL,
1043f4259b30SLisandro Dalcin                                        /*99*/ NULL,
1044f4259b30SLisandro Dalcin                                        NULL,
1045f4259b30SLisandro Dalcin                                        NULL,
1046d4002b98SHong Zhang                                        MatConjugate_MPISELL,
1047f4259b30SLisandro Dalcin                                        NULL,
1048f4259b30SLisandro Dalcin                                        /*104*/ NULL,
1049d4002b98SHong Zhang                                        MatRealPart_MPISELL,
1050d4002b98SHong Zhang                                        MatImaginaryPart_MPISELL,
1051f4259b30SLisandro Dalcin                                        NULL,
1052f4259b30SLisandro Dalcin                                        NULL,
1053f4259b30SLisandro Dalcin                                        /*109*/ NULL,
1054f4259b30SLisandro Dalcin                                        NULL,
1055f4259b30SLisandro Dalcin                                        NULL,
1056f4259b30SLisandro Dalcin                                        NULL,
1057d4002b98SHong Zhang                                        MatMissingDiagonal_MPISELL,
1058f4259b30SLisandro Dalcin                                        /*114*/ NULL,
1059f4259b30SLisandro Dalcin                                        NULL,
1060d4002b98SHong Zhang                                        MatGetGhosts_MPISELL,
1061f4259b30SLisandro Dalcin                                        NULL,
1062f4259b30SLisandro Dalcin                                        NULL,
1063f4259b30SLisandro Dalcin                                        /*119*/ NULL,
1064f4259b30SLisandro Dalcin                                        NULL,
1065f4259b30SLisandro Dalcin                                        NULL,
1066f4259b30SLisandro Dalcin                                        NULL,
1067f4259b30SLisandro Dalcin                                        NULL,
1068f4259b30SLisandro Dalcin                                        /*124*/ NULL,
1069f4259b30SLisandro Dalcin                                        NULL,
1070d4002b98SHong Zhang                                        MatInvertBlockDiagonal_MPISELL,
1071f4259b30SLisandro Dalcin                                        NULL,
1072f4259b30SLisandro Dalcin                                        NULL,
1073f4259b30SLisandro Dalcin                                        /*129*/ NULL,
1074f4259b30SLisandro Dalcin                                        NULL,
1075f4259b30SLisandro Dalcin                                        NULL,
1076f4259b30SLisandro Dalcin                                        NULL,
1077f4259b30SLisandro Dalcin                                        NULL,
1078f4259b30SLisandro Dalcin                                        /*134*/ NULL,
1079f4259b30SLisandro Dalcin                                        NULL,
1080f4259b30SLisandro Dalcin                                        NULL,
1081f4259b30SLisandro Dalcin                                        NULL,
1082f4259b30SLisandro Dalcin                                        NULL,
1083f4259b30SLisandro Dalcin                                        /*139*/ NULL,
1084f4259b30SLisandro Dalcin                                        NULL,
1085f4259b30SLisandro Dalcin                                        NULL,
1086d4002b98SHong Zhang                                        MatFDColoringSetUp_MPIXAIJ,
1087f4259b30SLisandro Dalcin                                        NULL,
1088d70f29a3SPierre Jolivet                                        /*144*/ NULL,
1089d70f29a3SPierre Jolivet                                        NULL,
1090d70f29a3SPierre Jolivet                                        NULL,
109199a7f59eSMark Adams                                        NULL,
109299a7f59eSMark Adams                                        NULL,
10937fb60732SBarry Smith                                        NULL,
10949371c9d4SSatish Balay                                        /*150*/ NULL};
1095d4002b98SHong Zhang 
1096d4002b98SHong Zhang /* ----------------------------------------------------------------------------------------*/
1097d4002b98SHong Zhang 
10989371c9d4SSatish Balay PetscErrorCode MatStoreValues_MPISELL(Mat mat) {
1099d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
1100d4002b98SHong Zhang 
1101d4002b98SHong Zhang   PetscFunctionBegin;
11029566063dSJacob Faibussowitsch   PetscCall(MatStoreValues(sell->A));
11039566063dSJacob Faibussowitsch   PetscCall(MatStoreValues(sell->B));
1104d4002b98SHong Zhang   PetscFunctionReturn(0);
1105d4002b98SHong Zhang }
1106d4002b98SHong Zhang 
11079371c9d4SSatish Balay PetscErrorCode MatRetrieveValues_MPISELL(Mat mat) {
1108d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
1109d4002b98SHong Zhang 
1110d4002b98SHong Zhang   PetscFunctionBegin;
11119566063dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(sell->A));
11129566063dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(sell->B));
1113d4002b98SHong Zhang   PetscFunctionReturn(0);
1114d4002b98SHong Zhang }
1115d4002b98SHong Zhang 
11169371c9d4SSatish Balay PetscErrorCode MatMPISELLSetPreallocation_MPISELL(Mat B, PetscInt d_rlenmax, const PetscInt d_rlen[], PetscInt o_rlenmax, const PetscInt o_rlen[]) {
1117d4002b98SHong Zhang   Mat_MPISELL *b;
1118d4002b98SHong Zhang 
1119d4002b98SHong Zhang   PetscFunctionBegin;
11209566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
11219566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
1122d4002b98SHong Zhang   b = (Mat_MPISELL *)B->data;
1123d4002b98SHong Zhang 
1124d4002b98SHong Zhang   if (!B->preallocated) {
1125d4002b98SHong Zhang     /* Explicitly create 2 MATSEQSELL matrices. */
11269566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &b->A));
11279566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(b->A, B->rmap->n, B->cmap->n, B->rmap->n, B->cmap->n));
11289566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(b->A, B, B));
11299566063dSJacob Faibussowitsch     PetscCall(MatSetType(b->A, MATSEQSELL));
11309566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)B, (PetscObject)b->A));
11319566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &b->B));
11329566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(b->B, B->rmap->n, B->cmap->N, B->rmap->n, B->cmap->N));
11339566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(b->B, B, B));
11349566063dSJacob Faibussowitsch     PetscCall(MatSetType(b->B, MATSEQSELL));
11359566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)B, (PetscObject)b->B));
1136d4002b98SHong Zhang   }
1137d4002b98SHong Zhang 
11389566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(b->A, d_rlenmax, d_rlen));
11399566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(b->B, o_rlenmax, o_rlen));
1140d4002b98SHong Zhang   B->preallocated  = PETSC_TRUE;
1141d4002b98SHong Zhang   B->was_assembled = PETSC_FALSE;
1142d4002b98SHong Zhang   /*
1143d4002b98SHong Zhang     critical for MatAssemblyEnd to work.
1144d4002b98SHong Zhang     MatAssemblyBegin checks it to set up was_assembled
1145d4002b98SHong Zhang     and MatAssemblyEnd checks was_assembled to determine whether to build garray
1146d4002b98SHong Zhang   */
1147d4002b98SHong Zhang   B->assembled     = PETSC_FALSE;
1148d4002b98SHong Zhang   PetscFunctionReturn(0);
1149d4002b98SHong Zhang }
1150d4002b98SHong Zhang 
11519371c9d4SSatish Balay PetscErrorCode MatDuplicate_MPISELL(Mat matin, MatDuplicateOption cpvalues, Mat *newmat) {
1152d4002b98SHong Zhang   Mat          mat;
1153d4002b98SHong Zhang   Mat_MPISELL *a, *oldmat = (Mat_MPISELL *)matin->data;
1154d4002b98SHong Zhang 
1155d4002b98SHong Zhang   PetscFunctionBegin;
1156f4259b30SLisandro Dalcin   *newmat = NULL;
11579566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)matin), &mat));
11589566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, matin->rmap->n, matin->cmap->n, matin->rmap->N, matin->cmap->N));
11599566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizesFromMats(mat, matin, matin));
11609566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, ((PetscObject)matin)->type_name));
1161d4002b98SHong Zhang   a = (Mat_MPISELL *)mat->data;
1162d4002b98SHong Zhang 
1163d4002b98SHong Zhang   mat->factortype   = matin->factortype;
1164d4002b98SHong Zhang   mat->assembled    = PETSC_TRUE;
1165d4002b98SHong Zhang   mat->insertmode   = NOT_SET_VALUES;
1166d4002b98SHong Zhang   mat->preallocated = PETSC_TRUE;
1167d4002b98SHong Zhang 
1168d4002b98SHong Zhang   a->size         = oldmat->size;
1169d4002b98SHong Zhang   a->rank         = oldmat->rank;
1170d4002b98SHong Zhang   a->donotstash   = oldmat->donotstash;
1171d4002b98SHong Zhang   a->roworiented  = oldmat->roworiented;
1172f4259b30SLisandro Dalcin   a->rowindices   = NULL;
1173f4259b30SLisandro Dalcin   a->rowvalues    = NULL;
1174d4002b98SHong Zhang   a->getrowactive = PETSC_FALSE;
1175d4002b98SHong Zhang 
11769566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->rmap, &mat->rmap));
11779566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->cmap, &mat->cmap));
1178d4002b98SHong Zhang 
1179d4002b98SHong Zhang   if (oldmat->colmap) {
1180d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
11819566063dSJacob Faibussowitsch     PetscCall(PetscTableCreateCopy(oldmat->colmap, &a->colmap));
1182d4002b98SHong Zhang #else
11839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mat->cmap->N, &a->colmap));
11849566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)mat, (mat->cmap->N) * sizeof(PetscInt)));
11859566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->colmap, oldmat->colmap, mat->cmap->N));
1186d4002b98SHong Zhang #endif
1187f4259b30SLisandro Dalcin   } else a->colmap = NULL;
1188d4002b98SHong Zhang   if (oldmat->garray) {
1189d4002b98SHong Zhang     PetscInt len;
1190d4002b98SHong Zhang     len = oldmat->B->cmap->n;
11919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len + 1, &a->garray));
11929566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)mat, len * sizeof(PetscInt)));
11939566063dSJacob Faibussowitsch     if (len) PetscCall(PetscArraycpy(a->garray, oldmat->garray, len));
1194f4259b30SLisandro Dalcin   } else a->garray = NULL;
1195d4002b98SHong Zhang 
11969566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(oldmat->lvec, &a->lvec));
11979566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat, (PetscObject)a->lvec));
11989566063dSJacob Faibussowitsch   PetscCall(VecScatterCopy(oldmat->Mvctx, &a->Mvctx));
11999566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat, (PetscObject)a->Mvctx));
12009566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->A, cpvalues, &a->A));
12019566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat, (PetscObject)a->A));
12029566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->B, cpvalues, &a->B));
12039566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat, (PetscObject)a->B));
12049566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDuplicate(((PetscObject)matin)->qlist, &((PetscObject)mat)->qlist));
1205d4002b98SHong Zhang   *newmat = mat;
1206d4002b98SHong Zhang   PetscFunctionReturn(0);
1207d4002b98SHong Zhang }
1208d4002b98SHong Zhang 
1209d4002b98SHong Zhang /*@C
1210d4002b98SHong Zhang    MatMPISELLSetPreallocation - Preallocates memory for a sparse parallel matrix in sell format.
1211d4002b98SHong Zhang    For good matrix assembly performance the user should preallocate the matrix storage by
1212d4002b98SHong Zhang    setting the parameters d_nz (or d_nnz) and o_nz (or o_nnz).
1213d4002b98SHong Zhang 
1214d083f849SBarry Smith    Collective
1215d4002b98SHong Zhang 
1216d4002b98SHong Zhang    Input Parameters:
1217d4002b98SHong Zhang +  B - the matrix
1218d4002b98SHong Zhang .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1219d4002b98SHong Zhang            (same value is used for all local rows)
1220d4002b98SHong Zhang .  d_nnz - array containing the number of nonzeros in the various rows of the
1221d4002b98SHong Zhang            DIAGONAL portion of the local submatrix (possibly different for each row)
1222d4002b98SHong Zhang            or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
1223d4002b98SHong Zhang            The size of this array is equal to the number of local rows, i.e 'm'.
1224d4002b98SHong Zhang            For matrices that will be factored, you must leave room for (and set)
1225d4002b98SHong Zhang            the diagonal entry even if it is zero.
1226d4002b98SHong Zhang .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1227d4002b98SHong Zhang            submatrix (same value is used for all local rows).
1228d4002b98SHong Zhang -  o_nnz - array containing the number of nonzeros in the various rows of the
1229d4002b98SHong Zhang            OFF-DIAGONAL portion of the local submatrix (possibly different for
1230d4002b98SHong Zhang            each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
1231d4002b98SHong Zhang            structure. The size of this array is equal to the number
1232d4002b98SHong Zhang            of local rows, i.e 'm'.
1233d4002b98SHong Zhang 
1234d4002b98SHong Zhang    If the *_nnz parameter is given then the *_nz parameter is ignored
1235d4002b98SHong Zhang 
1236d4002b98SHong Zhang    The stored row and column indices begin with zero.
1237d4002b98SHong Zhang 
1238d4002b98SHong Zhang    The parallel matrix is partitioned such that the first m0 rows belong to
1239d4002b98SHong Zhang    process 0, the next m1 rows belong to process 1, the next m2 rows belong
1240d4002b98SHong Zhang    to process 2 etc.. where m0,m1,m2... are the input parameter 'm'.
1241d4002b98SHong Zhang 
1242d4002b98SHong Zhang    The DIAGONAL portion of the local submatrix of a processor can be defined
1243d4002b98SHong Zhang    as the submatrix which is obtained by extraction the part corresponding to
1244d4002b98SHong Zhang    the rows r1-r2 and columns c1-c2 of the global matrix, where r1 is the
1245d4002b98SHong Zhang    first row that belongs to the processor, r2 is the last row belonging to
1246d4002b98SHong Zhang    the this processor, and c1-c2 is range of indices of the local part of a
1247d4002b98SHong Zhang    vector suitable for applying the matrix to.  This is an mxn matrix.  In the
1248d4002b98SHong Zhang    common case of a square matrix, the row and column ranges are the same and
1249d4002b98SHong Zhang    the DIAGONAL part is also square. The remaining portion of the local
1250d4002b98SHong Zhang    submatrix (mxN) constitute the OFF-DIAGONAL portion.
1251d4002b98SHong Zhang 
1252d4002b98SHong Zhang    If o_nnz, d_nnz are specified, then o_nz, and d_nz are ignored.
1253d4002b98SHong Zhang 
1254d4002b98SHong Zhang    You can call MatGetInfo() to get information on how effective the preallocation was;
1255d4002b98SHong Zhang    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
1256d4002b98SHong Zhang    You can also run with the option -info and look for messages with the string
1257d4002b98SHong Zhang    malloc in them to see if additional memory allocation was needed.
1258d4002b98SHong Zhang 
1259d4002b98SHong Zhang    Example usage:
1260d4002b98SHong Zhang 
1261d4002b98SHong Zhang    Consider the following 8x8 matrix with 34 non-zero values, that is
1262d4002b98SHong Zhang    assembled across 3 processors. Lets assume that proc0 owns 3 rows,
1263d4002b98SHong Zhang    proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
1264d4002b98SHong Zhang    as follows:
1265d4002b98SHong Zhang 
1266d4002b98SHong Zhang .vb
1267d4002b98SHong Zhang             1  2  0  |  0  3  0  |  0  4
1268d4002b98SHong Zhang     Proc0   0  5  6  |  7  0  0  |  8  0
1269d4002b98SHong Zhang             9  0 10  | 11  0  0  | 12  0
1270d4002b98SHong Zhang     -------------------------------------
1271d4002b98SHong Zhang            13  0 14  | 15 16 17  |  0  0
1272d4002b98SHong Zhang     Proc1   0 18  0  | 19 20 21  |  0  0
1273d4002b98SHong Zhang             0  0  0  | 22 23  0  | 24  0
1274d4002b98SHong Zhang     -------------------------------------
1275d4002b98SHong Zhang     Proc2  25 26 27  |  0  0 28  | 29  0
1276d4002b98SHong Zhang            30  0  0  | 31 32 33  |  0 34
1277d4002b98SHong Zhang .ve
1278d4002b98SHong Zhang 
1279d4002b98SHong Zhang    This can be represented as a collection of submatrices as:
1280d4002b98SHong Zhang 
1281d4002b98SHong Zhang .vb
1282d4002b98SHong Zhang       A B C
1283d4002b98SHong Zhang       D E F
1284d4002b98SHong Zhang       G H I
1285d4002b98SHong Zhang .ve
1286d4002b98SHong Zhang 
1287d4002b98SHong Zhang    Where the submatrices A,B,C are owned by proc0, D,E,F are
1288d4002b98SHong Zhang    owned by proc1, G,H,I are owned by proc2.
1289d4002b98SHong Zhang 
1290d4002b98SHong Zhang    The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1291d4002b98SHong Zhang    The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1292d4002b98SHong Zhang    The 'M','N' parameters are 8,8, and have the same values on all procs.
1293d4002b98SHong Zhang 
1294d4002b98SHong Zhang    The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
1295d4002b98SHong Zhang    submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
1296d4002b98SHong Zhang    corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
1297d4002b98SHong Zhang    Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
1298d4002b98SHong Zhang    part as SeqSELL matrices. for eg: proc1 will store [E] as a SeqSELL
1299d4002b98SHong Zhang    matrix, ans [DF] as another SeqSELL matrix.
1300d4002b98SHong Zhang 
1301d4002b98SHong Zhang    When d_nz, o_nz parameters are specified, d_nz storage elements are
1302d4002b98SHong Zhang    allocated for every row of the local diagonal submatrix, and o_nz
1303d4002b98SHong Zhang    storage locations are allocated for every row of the OFF-DIAGONAL submat.
1304d4002b98SHong Zhang    One way to choose d_nz and o_nz is to use the max nonzerors per local
1305d4002b98SHong Zhang    rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
1306d4002b98SHong Zhang    In this case, the values of d_nz,o_nz are:
1307d4002b98SHong Zhang .vb
1308d4002b98SHong Zhang      proc0 : dnz = 2, o_nz = 2
1309d4002b98SHong Zhang      proc1 : dnz = 3, o_nz = 2
1310d4002b98SHong Zhang      proc2 : dnz = 1, o_nz = 4
1311d4002b98SHong Zhang .ve
1312d4002b98SHong Zhang    We are allocating m*(d_nz+o_nz) storage locations for every proc. This
1313d4002b98SHong Zhang    translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
1314d4002b98SHong Zhang    for proc3. i.e we are using 12+15+10=37 storage locations to store
1315d4002b98SHong Zhang    34 values.
1316d4002b98SHong Zhang 
1317d4002b98SHong Zhang    When d_nnz, o_nnz parameters are specified, the storage is specified
1318a5b23f4aSJose E. Roman    for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
1319d4002b98SHong Zhang    In the above case the values for d_nnz,o_nnz are:
1320d4002b98SHong Zhang .vb
1321d4002b98SHong Zhang      proc0: d_nnz = [2,2,2] and o_nnz = [2,2,2]
1322d4002b98SHong Zhang      proc1: d_nnz = [3,3,2] and o_nnz = [2,1,1]
1323d4002b98SHong Zhang      proc2: d_nnz = [1,1]   and o_nnz = [4,4]
1324d4002b98SHong Zhang .ve
1325d4002b98SHong Zhang    Here the space allocated is according to nz (or maximum values in the nnz
1326d4002b98SHong Zhang    if nnz is provided) for DIAGONAL and OFF-DIAGONAL submatrices, i.e (2+2+3+2)*3+(1+4)*2=37
1327d4002b98SHong Zhang 
1328d4002b98SHong Zhang    Level: intermediate
1329d4002b98SHong Zhang 
1330db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateSeqSELL()`, `MatSetValues()`, `MatCreatesell()`,
1331db781477SPatrick Sanan           `MATMPISELL`, `MatGetInfo()`, `PetscSplitOwnership()`
1332d4002b98SHong Zhang @*/
13339371c9d4SSatish Balay PetscErrorCode MatMPISELLSetPreallocation(Mat B, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[]) {
1334d4002b98SHong Zhang   PetscFunctionBegin;
1335d4002b98SHong Zhang   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
1336d4002b98SHong Zhang   PetscValidType(B, 1);
1337cac4c232SBarry Smith   PetscTryMethod(B, "MatMPISELLSetPreallocation_C", (Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]), (B, d_nz, d_nnz, o_nz, o_nnz));
1338d4002b98SHong Zhang   PetscFunctionReturn(0);
1339d4002b98SHong Zhang }
1340d4002b98SHong Zhang 
1341ed73aabaSBarry Smith /*MC
1342ed73aabaSBarry Smith    MATMPISELL - MATMPISELL = "mpisell" - A matrix type to be used for MPI sparse matrices,
1343ed73aabaSBarry Smith    based on the sliced Ellpack format
1344ed73aabaSBarry Smith 
1345ed73aabaSBarry Smith    Options Database Keys:
1346ed73aabaSBarry Smith . -mat_type sell - sets the matrix type to "seqsell" during a call to MatSetFromOptions()
1347ed73aabaSBarry Smith 
1348ed73aabaSBarry Smith    Level: beginner
1349ed73aabaSBarry Smith 
1350db781477SPatrick Sanan .seealso: `MatCreateSell()`, `MATSEQSELL`, `MATSELL`, `MATSEQAIJ`, `MATAIJ`, `MATMPIAIJ`
1351ed73aabaSBarry Smith M*/
1352ed73aabaSBarry Smith 
1353d4002b98SHong Zhang /*@C
1354d4002b98SHong Zhang    MatCreateSELL - Creates a sparse parallel matrix in SELL format.
1355d4002b98SHong Zhang 
1356d083f849SBarry Smith    Collective
1357d4002b98SHong Zhang 
1358d4002b98SHong Zhang    Input Parameters:
1359d4002b98SHong Zhang +  comm - MPI communicator
1360d4002b98SHong Zhang .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1361d4002b98SHong Zhang            This value should be the same as the local size used in creating the
1362d4002b98SHong Zhang            y vector for the matrix-vector product y = Ax.
1363d4002b98SHong Zhang .  n - This value should be the same as the local size used in creating the
1364d4002b98SHong Zhang        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
1365d4002b98SHong Zhang        calculated if N is given) For square matrices n is almost always m.
1366d4002b98SHong Zhang .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
1367d4002b98SHong Zhang .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
1368d4002b98SHong Zhang .  d_rlenmax - max number of nonzeros per row in DIAGONAL portion of local submatrix
1369d4002b98SHong Zhang                (same value is used for all local rows)
1370d4002b98SHong Zhang .  d_rlen - array containing the number of nonzeros in the various rows of the
1371d4002b98SHong Zhang             DIAGONAL portion of the local submatrix (possibly different for each row)
1372d4002b98SHong Zhang             or NULL, if d_rlenmax is used to specify the nonzero structure.
1373d4002b98SHong Zhang             The size of this array is equal to the number of local rows, i.e 'm'.
1374d4002b98SHong Zhang .  o_rlenmax - max number of nonzeros per row in the OFF-DIAGONAL portion of local
1375d4002b98SHong Zhang                submatrix (same value is used for all local rows).
1376d4002b98SHong Zhang -  o_rlen - array containing the number of nonzeros in the various rows of the
1377d4002b98SHong Zhang             OFF-DIAGONAL portion of the local submatrix (possibly different for
1378d4002b98SHong Zhang             each row) or NULL, if o_rlenmax is used to specify the nonzero
1379d4002b98SHong Zhang             structure. The size of this array is equal to the number
1380d4002b98SHong Zhang             of local rows, i.e 'm'.
1381d4002b98SHong Zhang 
1382d4002b98SHong Zhang    Output Parameter:
1383d4002b98SHong Zhang .  A - the matrix
1384d4002b98SHong Zhang 
1385d4002b98SHong Zhang    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
1386f6f02116SRichard Tran Mills    MatXXXXSetPreallocation() paradigm instead of this routine directly.
1387d4002b98SHong Zhang    [MatXXXXSetPreallocation() is, for example, MatSeqSELLSetPreallocation]
1388d4002b98SHong Zhang 
1389d4002b98SHong Zhang    Notes:
1390d4002b98SHong Zhang    If the *_rlen parameter is given then the *_rlenmax parameter is ignored
1391d4002b98SHong Zhang 
1392d4002b98SHong Zhang    m,n,M,N parameters specify the size of the matrix, and its partitioning across
1393d4002b98SHong Zhang    processors, while d_rlenmax,d_rlen,o_rlenmax,o_rlen parameters specify the approximate
1394d4002b98SHong Zhang    storage requirements for this matrix.
1395d4002b98SHong Zhang 
1396d4002b98SHong Zhang    If PETSC_DECIDE or  PETSC_DETERMINE is used for a particular argument on one
1397d4002b98SHong Zhang    processor than it must be used on all processors that share the object for
1398d4002b98SHong Zhang    that argument.
1399d4002b98SHong Zhang 
1400d4002b98SHong Zhang    The user MUST specify either the local or global matrix dimensions
1401d4002b98SHong Zhang    (possibly both).
1402d4002b98SHong Zhang 
1403d4002b98SHong Zhang    The parallel matrix is partitioned across processors such that the
1404d4002b98SHong Zhang    first m0 rows belong to process 0, the next m1 rows belong to
1405d4002b98SHong Zhang    process 1, the next m2 rows belong to process 2 etc.. where
1406d4002b98SHong Zhang    m0,m1,m2,.. are the input parameter 'm'. i.e each processor stores
1407d4002b98SHong Zhang    values corresponding to [m x N] submatrix.
1408d4002b98SHong Zhang 
1409d4002b98SHong Zhang    The columns are logically partitioned with the n0 columns belonging
1410d4002b98SHong Zhang    to 0th partition, the next n1 columns belonging to the next
1411d4002b98SHong Zhang    partition etc.. where n0,n1,n2... are the input parameter 'n'.
1412d4002b98SHong Zhang 
1413d4002b98SHong Zhang    The DIAGONAL portion of the local submatrix on any given processor
1414d4002b98SHong Zhang    is the submatrix corresponding to the rows and columns m,n
1415d4002b98SHong Zhang    corresponding to the given processor. i.e diagonal matrix on
1416d4002b98SHong Zhang    process 0 is [m0 x n0], diagonal matrix on process 1 is [m1 x n1]
1417d4002b98SHong Zhang    etc. The remaining portion of the local submatrix [m x (N-n)]
1418d4002b98SHong Zhang    constitute the OFF-DIAGONAL portion. The example below better
1419d4002b98SHong Zhang    illustrates this concept.
1420d4002b98SHong Zhang 
1421d4002b98SHong Zhang    For a square global matrix we define each processor's diagonal portion
1422d4002b98SHong Zhang    to be its local rows and the corresponding columns (a square submatrix);
1423d4002b98SHong Zhang    each processor's off-diagonal portion encompasses the remainder of the
1424d4002b98SHong Zhang    local matrix (a rectangular submatrix).
1425d4002b98SHong Zhang 
1426d4002b98SHong Zhang    If o_rlen, d_rlen are specified, then o_rlenmax, and d_rlenmax are ignored.
1427d4002b98SHong Zhang 
1428d4002b98SHong Zhang    When calling this routine with a single process communicator, a matrix of
1429d4002b98SHong Zhang    type SEQSELL is returned.  If a matrix of type MATMPISELL is desired for this
1430d4002b98SHong Zhang    type of communicator, use the construction mechanism:
1431d4002b98SHong Zhang      MatCreate(...,&A); MatSetType(A,MATMPISELL); MatSetSizes(A, m,n,M,N); MatMPISELLSetPreallocation(A,...);
1432d4002b98SHong Zhang 
1433d4002b98SHong Zhang    Options Database Keys:
1434d4002b98SHong Zhang -  -mat_sell_oneindex - Internally use indexing starting at 1
1435d4002b98SHong Zhang         rather than 0.  Note that when calling MatSetValues(),
1436d4002b98SHong Zhang         the user still MUST index entries starting at 0!
1437d4002b98SHong Zhang 
1438d4002b98SHong Zhang    Example usage:
1439d4002b98SHong Zhang 
1440d4002b98SHong Zhang    Consider the following 8x8 matrix with 34 non-zero values, that is
1441d4002b98SHong Zhang    assembled across 3 processors. Lets assume that proc0 owns 3 rows,
1442d4002b98SHong Zhang    proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
1443d4002b98SHong Zhang    as follows:
1444d4002b98SHong Zhang 
1445d4002b98SHong Zhang .vb
1446d4002b98SHong Zhang             1  2  0  |  0  3  0  |  0  4
1447d4002b98SHong Zhang     Proc0   0  5  6  |  7  0  0  |  8  0
1448d4002b98SHong Zhang             9  0 10  | 11  0  0  | 12  0
1449d4002b98SHong Zhang     -------------------------------------
1450d4002b98SHong Zhang            13  0 14  | 15 16 17  |  0  0
1451d4002b98SHong Zhang     Proc1   0 18  0  | 19 20 21  |  0  0
1452d4002b98SHong Zhang             0  0  0  | 22 23  0  | 24  0
1453d4002b98SHong Zhang     -------------------------------------
1454d4002b98SHong Zhang     Proc2  25 26 27  |  0  0 28  | 29  0
1455d4002b98SHong Zhang            30  0  0  | 31 32 33  |  0 34
1456d4002b98SHong Zhang .ve
1457d4002b98SHong Zhang 
1458d4002b98SHong Zhang    This can be represented as a collection of submatrices as:
1459d4002b98SHong Zhang 
1460d4002b98SHong Zhang .vb
1461d4002b98SHong Zhang       A B C
1462d4002b98SHong Zhang       D E F
1463d4002b98SHong Zhang       G H I
1464d4002b98SHong Zhang .ve
1465d4002b98SHong Zhang 
1466d4002b98SHong Zhang    Where the submatrices A,B,C are owned by proc0, D,E,F are
1467d4002b98SHong Zhang    owned by proc1, G,H,I are owned by proc2.
1468d4002b98SHong Zhang 
1469d4002b98SHong Zhang    The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1470d4002b98SHong Zhang    The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1471d4002b98SHong Zhang    The 'M','N' parameters are 8,8, and have the same values on all procs.
1472d4002b98SHong Zhang 
1473d4002b98SHong Zhang    The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
1474d4002b98SHong Zhang    submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
1475d4002b98SHong Zhang    corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
1476d4002b98SHong Zhang    Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
1477d4002b98SHong Zhang    part as SeqSELL matrices. for eg: proc1 will store [E] as a SeqSELL
1478d4002b98SHong Zhang    matrix, ans [DF] as another SeqSELL matrix.
1479d4002b98SHong Zhang 
1480d4002b98SHong Zhang    When d_rlenmax, o_rlenmax parameters are specified, d_rlenmax storage elements are
1481d4002b98SHong Zhang    allocated for every row of the local diagonal submatrix, and o_rlenmax
1482d4002b98SHong Zhang    storage locations are allocated for every row of the OFF-DIAGONAL submat.
1483d4002b98SHong Zhang    One way to choose d_rlenmax and o_rlenmax is to use the max nonzerors per local
1484d4002b98SHong Zhang    rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
1485d4002b98SHong Zhang    In this case, the values of d_rlenmax,o_rlenmax are:
1486d4002b98SHong Zhang .vb
1487d4002b98SHong Zhang      proc0 : d_rlenmax = 2, o_rlenmax = 2
1488d4002b98SHong Zhang      proc1 : d_rlenmax = 3, o_rlenmax = 2
1489d4002b98SHong Zhang      proc2 : d_rlenmax = 1, o_rlenmax = 4
1490d4002b98SHong Zhang .ve
1491d4002b98SHong Zhang    We are allocating m*(d_rlenmax+o_rlenmax) storage locations for every proc. This
1492d4002b98SHong Zhang    translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
1493d4002b98SHong Zhang    for proc3. i.e we are using 12+15+10=37 storage locations to store
1494d4002b98SHong Zhang    34 values.
1495d4002b98SHong Zhang 
1496d4002b98SHong Zhang    When d_rlen, o_rlen parameters are specified, the storage is specified
1497a5b23f4aSJose E. Roman    for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
1498d4002b98SHong Zhang    In the above case the values for d_nnz,o_nnz are:
1499d4002b98SHong Zhang .vb
1500d4002b98SHong Zhang      proc0: d_nnz = [2,2,2] and o_nnz = [2,2,2]
1501d4002b98SHong Zhang      proc1: d_nnz = [3,3,2] and o_nnz = [2,1,1]
1502d4002b98SHong Zhang      proc2: d_nnz = [1,1]   and o_nnz = [4,4]
1503d4002b98SHong Zhang .ve
1504d4002b98SHong Zhang    Here the space allocated is still 37 though there are 34 nonzeros because
1505d4002b98SHong Zhang    the allocation is always done according to rlenmax.
1506d4002b98SHong Zhang 
1507d4002b98SHong Zhang    Level: intermediate
1508d4002b98SHong Zhang 
1509db781477SPatrick Sanan .seealso: `MatCreate()`, `MatCreateSeqSELL()`, `MatSetValues()`, `MatMPISELLSetPreallocation()`, `MatMPISELLSetPreallocationSELL()`,
1510db781477SPatrick Sanan           `MATMPISELL`, `MatCreateMPISELLWithArrays()`
1511d4002b98SHong Zhang @*/
15129371c9d4SSatish Balay PetscErrorCode MatCreateSELL(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt d_rlenmax, const PetscInt d_rlen[], PetscInt o_rlenmax, const PetscInt o_rlen[], Mat *A) {
1513d4002b98SHong Zhang   PetscMPIInt size;
1514d4002b98SHong Zhang 
1515d4002b98SHong Zhang   PetscFunctionBegin;
15169566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
15179566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, n, M, N));
15189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1519d4002b98SHong Zhang   if (size > 1) {
15209566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A, MATMPISELL));
15219566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(*A, d_rlenmax, d_rlen, o_rlenmax, o_rlen));
1522d4002b98SHong Zhang   } else {
15239566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A, MATSEQSELL));
15249566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(*A, d_rlenmax, d_rlen));
1525d4002b98SHong Zhang   }
1526d4002b98SHong Zhang   PetscFunctionReturn(0);
1527d4002b98SHong Zhang }
1528d4002b98SHong Zhang 
15299371c9d4SSatish Balay PetscErrorCode MatMPISELLGetSeqSELL(Mat A, Mat *Ad, Mat *Ao, const PetscInt *colmap[]) {
1530d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1531d4002b98SHong Zhang   PetscBool    flg;
1532d4002b98SHong Zhang 
1533d4002b98SHong Zhang   PetscFunctionBegin;
15349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISELL, &flg));
153528b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "This function requires a MATMPISELL matrix as input");
1536d4002b98SHong Zhang   if (Ad) *Ad = a->A;
1537d4002b98SHong Zhang   if (Ao) *Ao = a->B;
1538d4002b98SHong Zhang   if (colmap) *colmap = a->garray;
1539d4002b98SHong Zhang   PetscFunctionReturn(0);
1540d4002b98SHong Zhang }
1541d4002b98SHong Zhang 
1542d4002b98SHong Zhang /*@C
1543d4002b98SHong Zhang      MatMPISELLGetLocalMatCondensed - Creates a SeqSELL matrix from an MATMPISELL matrix by taking all its local rows and NON-ZERO columns
1544d4002b98SHong Zhang 
1545d4002b98SHong Zhang     Not Collective
1546d4002b98SHong Zhang 
1547d4002b98SHong Zhang    Input Parameters:
1548d4002b98SHong Zhang +    A - the matrix
1549d4002b98SHong Zhang .    scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
1550d4002b98SHong Zhang -    row, col - index sets of rows and columns to extract (or NULL)
1551d4002b98SHong Zhang 
1552d4002b98SHong Zhang    Output Parameter:
1553d4002b98SHong Zhang .    A_loc - the local sequential matrix generated
1554d4002b98SHong Zhang 
1555d4002b98SHong Zhang     Level: developer
1556d4002b98SHong Zhang 
1557db781477SPatrick Sanan .seealso: `MatGetOwnershipRange()`, `MatMPISELLGetLocalMat()`
1558d4002b98SHong Zhang 
1559d4002b98SHong Zhang @*/
15609371c9d4SSatish Balay PetscErrorCode MatMPISELLGetLocalMatCondensed(Mat A, MatReuse scall, IS *row, IS *col, Mat *A_loc) {
1561d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1562d4002b98SHong Zhang   PetscInt     i, start, end, ncols, nzA, nzB, *cmap, imark, *idx;
1563d4002b98SHong Zhang   IS           isrowa, iscola;
1564d4002b98SHong Zhang   Mat         *aloc;
1565d4002b98SHong Zhang   PetscBool    match;
1566d4002b98SHong Zhang 
1567d4002b98SHong Zhang   PetscFunctionBegin;
15689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISELL, &match));
156928b400f6SJacob Faibussowitsch   PetscCheck(match, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Requires MATMPISELL matrix as input");
15709566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_Getlocalmatcondensed, A, 0, 0, 0));
1571d4002b98SHong Zhang   if (!row) {
15729371c9d4SSatish Balay     start = A->rmap->rstart;
15739371c9d4SSatish Balay     end   = A->rmap->rend;
15749566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, end - start, start, 1, &isrowa));
1575d4002b98SHong Zhang   } else {
1576d4002b98SHong Zhang     isrowa = *row;
1577d4002b98SHong Zhang   }
1578d4002b98SHong Zhang   if (!col) {
1579d4002b98SHong Zhang     start = A->cmap->rstart;
1580d4002b98SHong Zhang     cmap  = a->garray;
1581d4002b98SHong Zhang     nzA   = a->A->cmap->n;
1582d4002b98SHong Zhang     nzB   = a->B->cmap->n;
15839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nzA + nzB, &idx));
1584d4002b98SHong Zhang     ncols = 0;
1585d4002b98SHong Zhang     for (i = 0; i < nzB; i++) {
1586d4002b98SHong Zhang       if (cmap[i] < start) idx[ncols++] = cmap[i];
1587d4002b98SHong Zhang       else break;
1588d4002b98SHong Zhang     }
1589d4002b98SHong Zhang     imark = i;
1590d4002b98SHong Zhang     for (i = 0; i < nzA; i++) idx[ncols++] = start + i;
1591d4002b98SHong Zhang     for (i = imark; i < nzB; i++) idx[ncols++] = cmap[i];
15929566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncols, idx, PETSC_OWN_POINTER, &iscola));
1593d4002b98SHong Zhang   } else {
1594d4002b98SHong Zhang     iscola = *col;
1595d4002b98SHong Zhang   }
1596d4002b98SHong Zhang   if (scall != MAT_INITIAL_MATRIX) {
15979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &aloc));
1598d4002b98SHong Zhang     aloc[0] = *A_loc;
1599d4002b98SHong Zhang   }
16009566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, 1, &isrowa, &iscola, scall, &aloc));
1601d4002b98SHong Zhang   *A_loc = aloc[0];
16029566063dSJacob Faibussowitsch   PetscCall(PetscFree(aloc));
1603*48a46eb9SPierre Jolivet   if (!row) PetscCall(ISDestroy(&isrowa));
1604*48a46eb9SPierre Jolivet   if (!col) PetscCall(ISDestroy(&iscola));
16059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_Getlocalmatcondensed, A, 0, 0, 0));
1606d4002b98SHong Zhang   PetscFunctionReturn(0);
1607d4002b98SHong Zhang }
1608d4002b98SHong Zhang 
1609d4002b98SHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>
1610d4002b98SHong Zhang 
16119371c9d4SSatish Balay PetscErrorCode MatConvert_MPISELL_MPIAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) {
1612d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1613d4002b98SHong Zhang   Mat          B;
1614d4002b98SHong Zhang   Mat_MPIAIJ  *b;
1615d4002b98SHong Zhang 
1616d4002b98SHong Zhang   PetscFunctionBegin;
161728b400f6SJacob Faibussowitsch   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Matrix must be assembled");
1618d4002b98SHong Zhang 
161994a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
162094a8b381SRichard Tran Mills     B = *newmat;
162194a8b381SRichard Tran Mills   } else {
16229566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
16239566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATMPIAIJ));
16249566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
16259566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
16269566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL));
16279566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL));
162894a8b381SRichard Tran Mills   }
1629d4002b98SHong Zhang   b = (Mat_MPIAIJ *)B->data;
163094a8b381SRichard Tran Mills 
163194a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
16329566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_REUSE_MATRIX, &b->A));
16339566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_REUSE_MATRIX, &b->B));
163494a8b381SRichard Tran Mills   } else {
16359566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->A));
16369566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->B));
16379566063dSJacob Faibussowitsch     PetscCall(MatDisAssemble_MPISELL(A));
16389566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->A));
16399566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->B));
16409566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
16419566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
16429566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
16439566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
164494a8b381SRichard Tran Mills   }
1645d4002b98SHong Zhang 
1646d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
16479566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
1648d4002b98SHong Zhang   } else {
1649d4002b98SHong Zhang     *newmat = B;
1650d4002b98SHong Zhang   }
1651d4002b98SHong Zhang   PetscFunctionReturn(0);
1652d4002b98SHong Zhang }
1653d4002b98SHong Zhang 
16549371c9d4SSatish Balay PetscErrorCode MatConvert_MPIAIJ_MPISELL(Mat A, MatType newtype, MatReuse reuse, Mat *newmat) {
1655d4002b98SHong Zhang   Mat_MPIAIJ  *a = (Mat_MPIAIJ *)A->data;
1656d4002b98SHong Zhang   Mat          B;
1657d4002b98SHong Zhang   Mat_MPISELL *b;
1658d4002b98SHong Zhang 
1659d4002b98SHong Zhang   PetscFunctionBegin;
166028b400f6SJacob Faibussowitsch   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Matrix must be assembled");
1661d4002b98SHong Zhang 
166294a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
166394a8b381SRichard Tran Mills     B = *newmat;
166494a8b381SRichard Tran Mills   } else {
16659566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
16669566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATMPISELL));
16679566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
16689566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
16699566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL));
16709566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL));
167194a8b381SRichard Tran Mills   }
1672d4002b98SHong Zhang   b = (Mat_MPISELL *)B->data;
167394a8b381SRichard Tran Mills 
167494a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
16759566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_REUSE_MATRIX, &b->A));
16769566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_REUSE_MATRIX, &b->B));
167794a8b381SRichard Tran Mills   } else {
16789566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->A));
16799566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->B));
16809566063dSJacob Faibussowitsch     PetscCall(MatDisAssemble_MPIAIJ(A));
16819566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_INITIAL_MATRIX, &b->A));
16829566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_INITIAL_MATRIX, &b->B));
16839566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
16849566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
16859566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
16869566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
168794a8b381SRichard Tran Mills   }
1688d4002b98SHong Zhang 
1689d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
16909566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
1691d4002b98SHong Zhang   } else {
1692d4002b98SHong Zhang     *newmat = B;
1693d4002b98SHong Zhang   }
1694d4002b98SHong Zhang   PetscFunctionReturn(0);
1695d4002b98SHong Zhang }
1696d4002b98SHong Zhang 
16979371c9d4SSatish Balay PetscErrorCode MatSOR_MPISELL(Mat matin, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) {
1698d4002b98SHong Zhang   Mat_MPISELL *mat = (Mat_MPISELL *)matin->data;
1699f4259b30SLisandro Dalcin   Vec          bb1 = NULL;
1700d4002b98SHong Zhang 
1701d4002b98SHong Zhang   PetscFunctionBegin;
1702d4002b98SHong Zhang   if (flag == SOR_APPLY_UPPER) {
17039566063dSJacob Faibussowitsch     PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1704d4002b98SHong Zhang     PetscFunctionReturn(0);
1705d4002b98SHong Zhang   }
1706d4002b98SHong Zhang 
1707*48a46eb9SPierre Jolivet   if (its > 1 || ~flag & SOR_ZERO_INITIAL_GUESS || flag & SOR_EISENSTAT) PetscCall(VecDuplicate(bb, &bb1));
1708d4002b98SHong Zhang 
1709d4002b98SHong Zhang   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
1710d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17119566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1712d4002b98SHong Zhang       its--;
1713d4002b98SHong Zhang     }
1714d4002b98SHong Zhang 
1715d4002b98SHong Zhang     while (its--) {
17169566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
17179566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1718d4002b98SHong Zhang 
1719d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
17209566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec, -1.0));
17219566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));
1722d4002b98SHong Zhang 
1723d4002b98SHong Zhang       /* local sweep */
17249566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_SYMMETRIC_SWEEP, fshift, lits, 1, xx));
1725d4002b98SHong Zhang     }
1726d4002b98SHong Zhang   } else if (flag & SOR_LOCAL_FORWARD_SWEEP) {
1727d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17289566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1729d4002b98SHong Zhang       its--;
1730d4002b98SHong Zhang     }
1731d4002b98SHong Zhang     while (its--) {
17329566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
17339566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1734d4002b98SHong Zhang 
1735d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
17369566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec, -1.0));
17379566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));
1738d4002b98SHong Zhang 
1739d4002b98SHong Zhang       /* local sweep */
17409566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_FORWARD_SWEEP, fshift, lits, 1, xx));
1741d4002b98SHong Zhang     }
1742d4002b98SHong Zhang   } else if (flag & SOR_LOCAL_BACKWARD_SWEEP) {
1743d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17449566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1745d4002b98SHong Zhang       its--;
1746d4002b98SHong Zhang     }
1747d4002b98SHong Zhang     while (its--) {
17489566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
17499566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1750d4002b98SHong Zhang 
1751d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
17529566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec, -1.0));
17539566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));
1754d4002b98SHong Zhang 
1755d4002b98SHong Zhang       /* local sweep */
17569566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_BACKWARD_SWEEP, fshift, lits, 1, xx));
1757d4002b98SHong Zhang     }
1758d4002b98SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)matin), PETSC_ERR_SUP, "Parallel SOR not supported");
1759d4002b98SHong Zhang 
17609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bb1));
1761d4002b98SHong Zhang 
1762d4002b98SHong Zhang   matin->factorerrortype = mat->A->factorerrortype;
1763d4002b98SHong Zhang   PetscFunctionReturn(0);
1764d4002b98SHong Zhang }
1765d4002b98SHong Zhang 
1766d4002b98SHong Zhang /*MC
1767d4002b98SHong Zhang    MATMPISELL - MATMPISELL = "MPISELL" - A matrix type to be used for parallel sparse matrices.
1768d4002b98SHong Zhang 
1769d4002b98SHong Zhang    Options Database Keys:
1770d4002b98SHong Zhang . -mat_type MPISELL - sets the matrix type to "MPISELL" during a call to MatSetFromOptions()
1771d4002b98SHong Zhang 
1772d4002b98SHong Zhang   Level: beginner
1773d4002b98SHong Zhang 
1774db781477SPatrick Sanan .seealso: `MatCreateSELL()`
1775d4002b98SHong Zhang M*/
17769371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat B) {
1777d4002b98SHong Zhang   Mat_MPISELL *b;
1778d4002b98SHong Zhang   PetscMPIInt  size;
1779d4002b98SHong Zhang 
1780d4002b98SHong Zhang   PetscFunctionBegin;
17819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
17829566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B, &b));
1783d4002b98SHong Zhang   B->data = (void *)b;
17849566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(B->ops, &MatOps_Values, sizeof(struct _MatOps)));
1785d4002b98SHong Zhang   B->assembled  = PETSC_FALSE;
1786d4002b98SHong Zhang   B->insertmode = NOT_SET_VALUES;
1787d4002b98SHong Zhang   b->size       = size;
17889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)B), &b->rank));
1789d4002b98SHong Zhang   /* build cache for off array entries formed */
17909566063dSJacob Faibussowitsch   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)B), 1, &B->stash));
1791d4002b98SHong Zhang 
1792d4002b98SHong Zhang   b->donotstash  = PETSC_FALSE;
1793f4259b30SLisandro Dalcin   b->colmap      = NULL;
1794f4259b30SLisandro Dalcin   b->garray      = NULL;
1795d4002b98SHong Zhang   b->roworiented = PETSC_TRUE;
1796d4002b98SHong Zhang 
1797d4002b98SHong Zhang   /* stuff used for matrix vector multiply */
1798d4002b98SHong Zhang   b->lvec  = NULL;
1799d4002b98SHong Zhang   b->Mvctx = NULL;
1800d4002b98SHong Zhang 
1801d4002b98SHong Zhang   /* stuff for MatGetRow() */
1802f4259b30SLisandro Dalcin   b->rowindices   = NULL;
1803f4259b30SLisandro Dalcin   b->rowvalues    = NULL;
1804d4002b98SHong Zhang   b->getrowactive = PETSC_FALSE;
1805d4002b98SHong Zhang 
18069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_MPISELL));
18079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_MPISELL));
18089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsTranspose_C", MatIsTranspose_MPISELL));
18099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPISELLSetPreallocation_C", MatMPISELLSetPreallocation_MPISELL));
18109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisell_mpiaij_C", MatConvert_MPISELL_MPIAIJ));
18119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatDiagonalScaleLocal_C", MatDiagonalScaleLocal_MPISELL));
18129566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPISELL));
1813d4002b98SHong Zhang   PetscFunctionReturn(0);
1814d4002b98SHong Zhang }
1815