xref: /petsc/src/mat/impls/sell/mpi/mpisell.c (revision b5917f1bc12a7822acea9597c8db5406cdee658a)
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 
1111a5261eSBarry Smith    This matrix type is identical to `MATSEQSELL` when constructed with a single process communicator,
1211a5261eSBarry Smith    and `MATMPISELL` otherwise.  As a result, for single process communicators,
1311a5261eSBarry Smith   `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:
1820f4b53cSBarry Smith . -mat_type sell - sets the matrix type to `MATSELL` during a call to `MatSetFromOptions()`
19d4002b98SHong Zhang 
20d4002b98SHong Zhang   Level: beginner
21d4002b98SHong Zhang 
2267be906fSBarry Smith .seealso: `Mat`, `MATAIJ`, `MATBAIJ`, `MATSBAIJ`, `MatCreateSELL()`, `MatCreateSeqSELL()`, `MATSEQSELL`, `MATMPISELL`
23d4002b98SHong Zhang M*/
24d4002b98SHong Zhang 
25d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalSet_MPISELL(Mat Y, Vec D, InsertMode is)
26d71ae5a4SJacob Faibussowitsch {
27d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)Y->data;
28d4002b98SHong Zhang 
29d4002b98SHong Zhang   PetscFunctionBegin;
30d4002b98SHong Zhang   if (Y->assembled && Y->rmap->rstart == Y->cmap->rstart && Y->rmap->rend == Y->cmap->rend) {
319566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet(sell->A, D, is));
32d4002b98SHong Zhang   } else {
339566063dSJacob Faibussowitsch     PetscCall(MatDiagonalSet_Default(Y, D, is));
34d4002b98SHong Zhang   }
353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36d4002b98SHong Zhang }
37d4002b98SHong Zhang 
38d4002b98SHong Zhang /*
39d4002b98SHong Zhang   Local utility routine that creates a mapping from the global column
40d4002b98SHong Zhang number to the local number in the off-diagonal part of the local
41d4002b98SHong Zhang storage of the matrix.  When PETSC_USE_CTABLE is used this is scalable at
42d4002b98SHong Zhang a slightly higher hash table cost; without it it is not scalable (each processor
43da81f932SPierre Jolivet has an order N integer array but is fast to access.
44d4002b98SHong Zhang */
45d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateColmap_MPISELL_Private(Mat mat)
46d71ae5a4SJacob Faibussowitsch {
47d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
48d4002b98SHong Zhang   PetscInt     n    = sell->B->cmap->n, i;
49d4002b98SHong Zhang 
50d4002b98SHong Zhang   PetscFunctionBegin;
5128b400f6SJacob Faibussowitsch   PetscCheck(sell->garray, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPISELL Matrix was assembled but is missing garray");
52d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
53eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapICreateWithSize(n, &sell->colmap));
54c76ffc5fSJacob Faibussowitsch   for (i = 0; i < n; i++) PetscCall(PetscHMapISet(sell->colmap, sell->garray[i] + 1, i + 1));
55d4002b98SHong Zhang #else
569566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(mat->cmap->N + 1, &sell->colmap));
57d4002b98SHong Zhang   for (i = 0; i < n; i++) sell->colmap[sell->garray[i]] = i + 1;
58d4002b98SHong Zhang #endif
593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60d4002b98SHong Zhang }
61d4002b98SHong Zhang 
62d4002b98SHong Zhang #define MatSetValues_SeqSELL_A_Private(row, col, value, addv, orow, ocol) \
63d4002b98SHong Zhang   { \
64d4002b98SHong Zhang     if (col <= lastcol1) low1 = 0; \
65d4002b98SHong Zhang     else high1 = nrow1; \
66d4002b98SHong Zhang     lastcol1 = col; \
67d4002b98SHong Zhang     while (high1 - low1 > 5) { \
68d4002b98SHong Zhang       t = (low1 + high1) / 2; \
694f8ff0b3SHong Zhang       if (cp1[sliceheight * t] > col) high1 = t; \
70d4002b98SHong Zhang       else low1 = t; \
71d4002b98SHong Zhang     } \
72d4002b98SHong Zhang     for (_i = low1; _i < high1; _i++) { \
734f8ff0b3SHong Zhang       if (cp1[sliceheight * _i] > col) break; \
744f8ff0b3SHong Zhang       if (cp1[sliceheight * _i] == col) { \
754f8ff0b3SHong Zhang         if (addv == ADD_VALUES) vp1[sliceheight * _i] += value; \
764f8ff0b3SHong Zhang         else vp1[sliceheight * _i] = value; \
772d1451d4SHong Zhang         inserted = PETSC_TRUE; \
78d4002b98SHong Zhang         goto a_noinsert; \
79d4002b98SHong Zhang       } \
80d4002b98SHong Zhang     } \
819371c9d4SSatish Balay     if (value == 0.0 && ignorezeroentries) { \
829371c9d4SSatish Balay       low1  = 0; \
839371c9d4SSatish Balay       high1 = nrow1; \
849371c9d4SSatish Balay       goto a_noinsert; \
859371c9d4SSatish Balay     } \
869371c9d4SSatish Balay     if (nonew == 1) { \
879371c9d4SSatish Balay       low1  = 0; \
889371c9d4SSatish Balay       high1 = nrow1; \
899371c9d4SSatish Balay       goto a_noinsert; \
909371c9d4SSatish Balay     } \
9108401ef6SPierre 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); \
9207e43b41SHong Zhang     MatSeqXSELLReallocateSELL(A, am, 1, nrow1, a->sliidx, a->sliceheight, row / sliceheight, row, col, a->colidx, a->val, cp1, vp1, nonew, MatScalar); \
93d4002b98SHong Zhang     /* shift up all the later entries in this row */ \
94d4002b98SHong Zhang     for (ii = nrow1 - 1; ii >= _i; ii--) { \
954f8ff0b3SHong Zhang       cp1[sliceheight * (ii + 1)] = cp1[sliceheight * ii]; \
964f8ff0b3SHong Zhang       vp1[sliceheight * (ii + 1)] = vp1[sliceheight * ii]; \
97d4002b98SHong Zhang     } \
984f8ff0b3SHong Zhang     cp1[sliceheight * _i] = col; \
994f8ff0b3SHong Zhang     vp1[sliceheight * _i] = value; \
1009371c9d4SSatish Balay     a->nz++; \
1019371c9d4SSatish Balay     nrow1++; \
1029371c9d4SSatish Balay     A->nonzerostate++; \
103d4002b98SHong Zhang   a_noinsert:; \
104d4002b98SHong Zhang     a->rlen[row] = nrow1; \
105d4002b98SHong Zhang   }
106d4002b98SHong Zhang 
107d4002b98SHong Zhang #define MatSetValues_SeqSELL_B_Private(row, col, value, addv, orow, ocol) \
108d4002b98SHong Zhang   { \
109d4002b98SHong Zhang     if (col <= lastcol2) low2 = 0; \
110d4002b98SHong Zhang     else high2 = nrow2; \
111d4002b98SHong Zhang     lastcol2 = col; \
112d4002b98SHong Zhang     while (high2 - low2 > 5) { \
113d4002b98SHong Zhang       t = (low2 + high2) / 2; \
1144f8ff0b3SHong Zhang       if (cp2[sliceheight * t] > col) high2 = t; \
115d4002b98SHong Zhang       else low2 = t; \
116d4002b98SHong Zhang     } \
117d4002b98SHong Zhang     for (_i = low2; _i < high2; _i++) { \
1184f8ff0b3SHong Zhang       if (cp2[sliceheight * _i] > col) break; \
1194f8ff0b3SHong Zhang       if (cp2[sliceheight * _i] == col) { \
1204f8ff0b3SHong Zhang         if (addv == ADD_VALUES) vp2[sliceheight * _i] += value; \
1214f8ff0b3SHong Zhang         else vp2[sliceheight * _i] = value; \
1222d1451d4SHong Zhang         inserted = PETSC_TRUE; \
123d4002b98SHong Zhang         goto b_noinsert; \
124d4002b98SHong Zhang       } \
125d4002b98SHong Zhang     } \
1269371c9d4SSatish Balay     if (value == 0.0 && ignorezeroentries) { \
1279371c9d4SSatish Balay       low2  = 0; \
1289371c9d4SSatish Balay       high2 = nrow2; \
1299371c9d4SSatish Balay       goto b_noinsert; \
1309371c9d4SSatish Balay     } \
1319371c9d4SSatish Balay     if (nonew == 1) { \
1329371c9d4SSatish Balay       low2  = 0; \
1339371c9d4SSatish Balay       high2 = nrow2; \
1349371c9d4SSatish Balay       goto b_noinsert; \
1359371c9d4SSatish Balay     } \
13608401ef6SPierre 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); \
13707e43b41SHong Zhang     MatSeqXSELLReallocateSELL(B, bm, 1, nrow2, b->sliidx, b->sliceheight, row / sliceheight, row, col, b->colidx, b->val, cp2, vp2, nonew, MatScalar); \
138d4002b98SHong Zhang     /* shift up all the later entries in this row */ \
139d4002b98SHong Zhang     for (ii = nrow2 - 1; ii >= _i; ii--) { \
1404f8ff0b3SHong Zhang       cp2[sliceheight * (ii + 1)] = cp2[sliceheight * ii]; \
1414f8ff0b3SHong Zhang       vp2[sliceheight * (ii + 1)] = vp2[sliceheight * ii]; \
142d4002b98SHong Zhang     } \
1434f8ff0b3SHong Zhang     cp2[sliceheight * _i] = col; \
1444f8ff0b3SHong Zhang     vp2[sliceheight * _i] = value; \
1459371c9d4SSatish Balay     b->nz++; \
1469371c9d4SSatish Balay     nrow2++; \
1479371c9d4SSatish Balay     B->nonzerostate++; \
148d4002b98SHong Zhang   b_noinsert:; \
149d4002b98SHong Zhang     b->rlen[row] = nrow2; \
150d4002b98SHong Zhang   }
151d4002b98SHong Zhang 
152d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetValues_MPISELL(Mat mat, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode addv)
153d71ae5a4SJacob Faibussowitsch {
154d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
155d4002b98SHong Zhang   PetscScalar  value;
156d4002b98SHong Zhang   PetscInt     i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend, shift1, shift2;
157d4002b98SHong Zhang   PetscInt     cstart = mat->cmap->rstart, cend = mat->cmap->rend, row, col;
158d4002b98SHong Zhang   PetscBool    roworiented = sell->roworiented;
159d4002b98SHong Zhang 
160d4002b98SHong Zhang   /* Some Variables required in the macro */
161d4002b98SHong Zhang   Mat          A                 = sell->A;
162d4002b98SHong Zhang   Mat_SeqSELL *a                 = (Mat_SeqSELL *)A->data;
163d4002b98SHong Zhang   PetscBool    ignorezeroentries = a->ignorezeroentries, found;
164d4002b98SHong Zhang   Mat          B                 = sell->B;
165d4002b98SHong Zhang   Mat_SeqSELL *b                 = (Mat_SeqSELL *)B->data;
16607e43b41SHong Zhang   PetscInt    *cp1, *cp2, ii, _i, nrow1, nrow2, low1, high1, low2, high2, t, lastcol1, lastcol2, sliceheight = a->sliceheight;
167d4002b98SHong Zhang   MatScalar   *vp1, *vp2;
168d4002b98SHong Zhang 
169d4002b98SHong Zhang   PetscFunctionBegin;
170d4002b98SHong Zhang   for (i = 0; i < m; i++) {
171d4002b98SHong Zhang     if (im[i] < 0) continue;
1726bdcaf15SBarry 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);
173d4002b98SHong Zhang     if (im[i] >= rstart && im[i] < rend) {
174d4002b98SHong Zhang       row      = im[i] - rstart;
175d4002b98SHong Zhang       lastcol1 = -1;
17607e43b41SHong Zhang       shift1   = a->sliidx[row / sliceheight] + (row % sliceheight); /* starting index of the row */
177d4002b98SHong Zhang       cp1      = a->colidx + shift1;
178d4002b98SHong Zhang       vp1      = a->val + shift1;
179d4002b98SHong Zhang       nrow1    = a->rlen[row];
180d4002b98SHong Zhang       low1     = 0;
181d4002b98SHong Zhang       high1    = nrow1;
182d4002b98SHong Zhang       lastcol2 = -1;
18307e43b41SHong Zhang       shift2   = b->sliidx[row / sliceheight] + (row % sliceheight); /* starting index of the row */
184d4002b98SHong Zhang       cp2      = b->colidx + shift2;
185d4002b98SHong Zhang       vp2      = b->val + shift2;
186d4002b98SHong Zhang       nrow2    = b->rlen[row];
187d4002b98SHong Zhang       low2     = 0;
188d4002b98SHong Zhang       high2    = nrow2;
189d4002b98SHong Zhang 
190d4002b98SHong Zhang       for (j = 0; j < n; j++) {
191d4002b98SHong Zhang         if (roworiented) value = v[i * n + j];
192d4002b98SHong Zhang         else value = v[i + j * m];
193d4002b98SHong Zhang         if (ignorezeroentries && value == 0.0 && (addv == ADD_VALUES)) continue;
194d4002b98SHong Zhang         if (in[j] >= cstart && in[j] < cend) {
195d4002b98SHong Zhang           col = in[j] - cstart;
196d4002b98SHong Zhang           MatSetValue_SeqSELL_Private(A, row, col, value, addv, im[i], in[j], cp1, vp1, lastcol1, low1, high1); /* set one value */
1972d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
1982d1451d4SHong Zhang           if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED && found) A->offloadmask = PETSC_OFFLOAD_CPU;
1992d1451d4SHong Zhang #endif
200f7d195e4SLawrence Mitchell         } else if (in[j] < 0) {
201f7d195e4SLawrence Mitchell           continue;
202f7d195e4SLawrence Mitchell         } else {
203f7d195e4SLawrence 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);
204d4002b98SHong Zhang           if (mat->was_assembled) {
20548a46eb9SPierre Jolivet             if (!sell->colmap) PetscCall(MatCreateColmap_MPISELL_Private(mat));
206d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
207eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapIGetWithDefault(sell->colmap, in[j] + 1, 0, &col));
208d4002b98SHong Zhang             col--;
209d4002b98SHong Zhang #else
210d4002b98SHong Zhang             col = sell->colmap[in[j]] - 1;
211d4002b98SHong Zhang #endif
212d4002b98SHong Zhang             if (col < 0 && !((Mat_SeqSELL *)(sell->B->data))->nonew) {
2139566063dSJacob Faibussowitsch               PetscCall(MatDisAssemble_MPISELL(mat));
214d4002b98SHong Zhang               col = in[j];
215d4002b98SHong Zhang               /* Reinitialize the variables required by MatSetValues_SeqSELL_B_Private() */
216d4002b98SHong Zhang               B      = sell->B;
217d4002b98SHong Zhang               b      = (Mat_SeqSELL *)B->data;
21807e43b41SHong Zhang               shift2 = b->sliidx[row / sliceheight] + (row % sliceheight); /* starting index of the row */
219d4002b98SHong Zhang               cp2    = b->colidx + shift2;
220d4002b98SHong Zhang               vp2    = b->val + shift2;
221d4002b98SHong Zhang               nrow2  = b->rlen[row];
222d4002b98SHong Zhang               low2   = 0;
223d4002b98SHong Zhang               high2  = nrow2;
2242d1451d4SHong Zhang               found  = PETSC_FALSE;
225f7d195e4SLawrence Mitchell             } else {
226f7d195e4SLawrence 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]);
227f7d195e4SLawrence Mitchell             }
228d4002b98SHong Zhang           } else col = in[j];
229d4002b98SHong Zhang           MatSetValue_SeqSELL_Private(B, row, col, value, addv, im[i], in[j], cp2, vp2, lastcol2, low2, high2); /* set one value */
2302d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
2312d1451d4SHong Zhang           if (B->offloadmask != PETSC_OFFLOAD_UNALLOCATED && found) B->offloadmask = PETSC_OFFLOAD_CPU;
2322d1451d4SHong Zhang #endif
233d4002b98SHong Zhang         }
234d4002b98SHong Zhang       }
235d4002b98SHong Zhang     } else {
23628b400f6SJacob 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]);
237d4002b98SHong Zhang       if (!sell->donotstash) {
238d4002b98SHong Zhang         mat->assembled = PETSC_FALSE;
239d4002b98SHong Zhang         if (roworiented) {
2409566063dSJacob Faibussowitsch           PetscCall(MatStashValuesRow_Private(&mat->stash, im[i], n, in, v + i * n, (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
241d4002b98SHong Zhang         } else {
2429566063dSJacob Faibussowitsch           PetscCall(MatStashValuesCol_Private(&mat->stash, im[i], n, in, v + i, m, (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
243d4002b98SHong Zhang         }
244d4002b98SHong Zhang       }
245d4002b98SHong Zhang     }
246d4002b98SHong Zhang   }
2473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
248d4002b98SHong Zhang }
249d4002b98SHong Zhang 
250d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetValues_MPISELL(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[])
251d71ae5a4SJacob Faibussowitsch {
252d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
253d4002b98SHong Zhang   PetscInt     i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend;
254d4002b98SHong Zhang   PetscInt     cstart = mat->cmap->rstart, cend = mat->cmap->rend, row, col;
255d4002b98SHong Zhang 
256d4002b98SHong Zhang   PetscFunctionBegin;
257d4002b98SHong Zhang   for (i = 0; i < m; i++) {
25854c59aa7SJacob Faibussowitsch     if (idxm[i] < 0) continue; /* negative row */
25954c59aa7SJacob 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);
260d4002b98SHong Zhang     if (idxm[i] >= rstart && idxm[i] < rend) {
261d4002b98SHong Zhang       row = idxm[i] - rstart;
262d4002b98SHong Zhang       for (j = 0; j < n; j++) {
26354c59aa7SJacob Faibussowitsch         if (idxn[j] < 0) continue; /* negative column */
26454c59aa7SJacob 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);
265d4002b98SHong Zhang         if (idxn[j] >= cstart && idxn[j] < cend) {
266d4002b98SHong Zhang           col = idxn[j] - cstart;
2679566063dSJacob Faibussowitsch           PetscCall(MatGetValues(sell->A, 1, &row, 1, &col, v + i * n + j));
268d4002b98SHong Zhang         } else {
26948a46eb9SPierre Jolivet           if (!sell->colmap) PetscCall(MatCreateColmap_MPISELL_Private(mat));
270d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
271eec179cfSJacob Faibussowitsch           PetscCall(PetscHMapIGetWithDefault(sell->colmap, idxn[j] + 1, 0, &col));
272d4002b98SHong Zhang           col--;
273d4002b98SHong Zhang #else
274d4002b98SHong Zhang           col = sell->colmap[idxn[j]] - 1;
275d4002b98SHong Zhang #endif
276d4002b98SHong Zhang           if ((col < 0) || (sell->garray[col] != idxn[j])) *(v + i * n + j) = 0.0;
27748a46eb9SPierre Jolivet           else PetscCall(MatGetValues(sell->B, 1, &row, 1, &col, v + i * n + j));
278d4002b98SHong Zhang         }
279d4002b98SHong Zhang       }
280d4002b98SHong Zhang     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only local values currently supported");
281d4002b98SHong Zhang   }
2823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
283d4002b98SHong Zhang }
284d4002b98SHong Zhang 
285d4002b98SHong Zhang extern PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat, Vec, Vec);
286d4002b98SHong Zhang 
287d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyBegin_MPISELL(Mat mat, MatAssemblyType mode)
288d71ae5a4SJacob Faibussowitsch {
289d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
290d4002b98SHong Zhang   PetscInt     nstash, reallocs;
291d4002b98SHong Zhang 
292d4002b98SHong Zhang   PetscFunctionBegin;
2933ba16761SJacob Faibussowitsch   if (sell->donotstash || mat->nooffprocentries) PetscFunctionReturn(PETSC_SUCCESS);
294d4002b98SHong Zhang 
2959566063dSJacob Faibussowitsch   PetscCall(MatStashScatterBegin_Private(mat, &mat->stash, mat->rmap->range));
2969566063dSJacob Faibussowitsch   PetscCall(MatStashGetInfo_Private(&mat->stash, &nstash, &reallocs));
2979566063dSJacob Faibussowitsch   PetscCall(PetscInfo(sell->A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
299d4002b98SHong Zhang }
300d4002b98SHong Zhang 
301d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_MPISELL(Mat mat, MatAssemblyType mode)
302d71ae5a4SJacob Faibussowitsch {
303d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
304d4002b98SHong Zhang   PetscMPIInt  n;
305d4002b98SHong Zhang   PetscInt     i, flg;
306d4002b98SHong Zhang   PetscInt    *row, *col;
307d4002b98SHong Zhang   PetscScalar *val;
308d4002b98SHong Zhang   PetscBool    other_disassembled;
309d4002b98SHong Zhang   /* do not use 'b = (Mat_SeqSELL*)sell->B->data' as B can be reset in disassembly */
310d4002b98SHong Zhang   PetscFunctionBegin;
311d4002b98SHong Zhang   if (!sell->donotstash && !mat->nooffprocentries) {
312d4002b98SHong Zhang     while (1) {
3139566063dSJacob Faibussowitsch       PetscCall(MatStashScatterGetMesg_Private(&mat->stash, &n, &row, &col, &val, &flg));
314d4002b98SHong Zhang       if (!flg) break;
315d4002b98SHong Zhang 
316d4002b98SHong Zhang       for (i = 0; i < n; i++) { /* assemble one by one */
3179566063dSJacob Faibussowitsch         PetscCall(MatSetValues_MPISELL(mat, 1, row + i, 1, col + i, val + i, mat->insertmode));
318d4002b98SHong Zhang       }
319d4002b98SHong Zhang     }
3209566063dSJacob Faibussowitsch     PetscCall(MatStashScatterEnd_Private(&mat->stash));
321d4002b98SHong Zhang   }
3222d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
3232d1451d4SHong Zhang   if (mat->offloadmask == PETSC_OFFLOAD_CPU) sell->A->offloadmask = PETSC_OFFLOAD_CPU;
3242d1451d4SHong Zhang #endif
3259566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sell->A, mode));
3269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sell->A, mode));
327d4002b98SHong Zhang 
328d4002b98SHong Zhang   /*
329d4002b98SHong Zhang      determine if any processor has disassembled, if so we must
3306aad120cSJose E. Roman      also disassemble ourselves, in order that we may reassemble.
331d4002b98SHong Zhang   */
332d4002b98SHong Zhang   /*
333d4002b98SHong Zhang      if nonzero structure of submatrix B cannot change then we know that
334d4002b98SHong Zhang      no processor disassembled thus we can skip this stuff
335d4002b98SHong Zhang   */
336d4002b98SHong Zhang   if (!((Mat_SeqSELL *)sell->B->data)->nonew) {
3375f9db2b2SJunchao Zhang     PetscCall(MPIU_Allreduce(&mat->was_assembled, &other_disassembled, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)mat)));
3382d1451d4SHong Zhang     if (mat->was_assembled && !other_disassembled) PetscCall(MatDisAssemble_MPISELL(mat));
339d4002b98SHong Zhang   }
34048a46eb9SPierre Jolivet   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) PetscCall(MatSetUpMultiply_MPISELL(mat));
3412d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
3422d1451d4SHong Zhang   if (mat->offloadmask == PETSC_OFFLOAD_CPU && sell->B->offloadmask != PETSC_OFFLOAD_UNALLOCATED) sell->B->offloadmask = PETSC_OFFLOAD_CPU;
3432d1451d4SHong Zhang #endif
3449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sell->B, mode));
3459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sell->B, mode));
3469566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sell->rowvalues, sell->rowindices));
347f4259b30SLisandro Dalcin   sell->rowvalues = NULL;
3489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->diag));
349d4002b98SHong Zhang 
350d4002b98SHong Zhang   /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
351d4002b98SHong Zhang   if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqSELL *)(sell->A->data))->nonew) {
352d4002b98SHong Zhang     PetscObjectState state = sell->A->nonzerostate + sell->B->nonzerostate;
3531c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&state, &mat->nonzerostate, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)mat)));
354d4002b98SHong Zhang   }
3552d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
3562d1451d4SHong Zhang   mat->offloadmask = PETSC_OFFLOAD_BOTH;
3572d1451d4SHong Zhang #endif
3583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
359d4002b98SHong Zhang }
360d4002b98SHong Zhang 
361d71ae5a4SJacob Faibussowitsch PetscErrorCode MatZeroEntries_MPISELL(Mat A)
362d71ae5a4SJacob Faibussowitsch {
363d4002b98SHong Zhang   Mat_MPISELL *l = (Mat_MPISELL *)A->data;
364d4002b98SHong Zhang 
365d4002b98SHong Zhang   PetscFunctionBegin;
3669566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->A));
3679566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->B));
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
369d4002b98SHong Zhang }
370d4002b98SHong Zhang 
371d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_MPISELL(Mat A, Vec xx, Vec yy)
372d71ae5a4SJacob Faibussowitsch {
373d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
374d4002b98SHong Zhang   PetscInt     nt;
375d4002b98SHong Zhang 
376d4002b98SHong Zhang   PetscFunctionBegin;
3779566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xx, &nt));
37808401ef6SPierre 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);
3799566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3809566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->mult)(a->A, xx, yy));
3819566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3829566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B, a->lvec, yy, yy));
3833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
384d4002b98SHong Zhang }
385d4002b98SHong Zhang 
386d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat A, Vec bb, Vec xx)
387d71ae5a4SJacob Faibussowitsch {
388d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
389d4002b98SHong Zhang 
390d4002b98SHong Zhang   PetscFunctionBegin;
3919566063dSJacob Faibussowitsch   PetscCall(MatMultDiagonalBlock(a->A, bb, xx));
3923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
393d4002b98SHong Zhang }
394d4002b98SHong Zhang 
395d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAdd_MPISELL(Mat A, Vec xx, Vec yy, Vec zz)
396d71ae5a4SJacob Faibussowitsch {
397d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
398d4002b98SHong Zhang 
399d4002b98SHong Zhang   PetscFunctionBegin;
4009566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
4019566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multadd)(a->A, xx, yy, zz));
4029566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
4039566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B, a->lvec, zz, zz));
4043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
405d4002b98SHong Zhang }
406d4002b98SHong Zhang 
407d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTranspose_MPISELL(Mat A, Vec xx, Vec yy)
408d71ae5a4SJacob Faibussowitsch {
409d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
410d4002b98SHong Zhang 
411d4002b98SHong Zhang   PetscFunctionBegin;
412d4002b98SHong Zhang   /* do nondiagonal part */
4139566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B, xx, a->lvec));
414d4002b98SHong Zhang   /* do local part */
4159566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multtranspose)(a->A, xx, yy));
416a29b4eb7SJunchao Zhang   /* add partial results together */
4179566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, a->lvec, yy, ADD_VALUES, SCATTER_REVERSE));
4189566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, a->lvec, yy, ADD_VALUES, SCATTER_REVERSE));
4193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
420d4002b98SHong Zhang }
421d4002b98SHong Zhang 
422d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIsTranspose_MPISELL(Mat Amat, Mat Bmat, PetscReal tol, PetscBool *f)
423d71ae5a4SJacob Faibussowitsch {
424d4002b98SHong Zhang   MPI_Comm     comm;
425d4002b98SHong Zhang   Mat_MPISELL *Asell = (Mat_MPISELL *)Amat->data, *Bsell;
426d4002b98SHong Zhang   Mat          Adia  = Asell->A, Bdia, Aoff, Boff, *Aoffs, *Boffs;
427d4002b98SHong Zhang   IS           Me, Notme;
428d4002b98SHong Zhang   PetscInt     M, N, first, last, *notme, i;
429d4002b98SHong Zhang   PetscMPIInt  size;
430d4002b98SHong Zhang 
431d4002b98SHong Zhang   PetscFunctionBegin;
432d4002b98SHong Zhang   /* Easy test: symmetric diagonal block */
4339371c9d4SSatish Balay   Bsell = (Mat_MPISELL *)Bmat->data;
4349371c9d4SSatish Balay   Bdia  = Bsell->A;
4359566063dSJacob Faibussowitsch   PetscCall(MatIsTranspose(Adia, Bdia, tol, f));
4363ba16761SJacob Faibussowitsch   if (!*f) PetscFunctionReturn(PETSC_SUCCESS);
4379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
4389566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
4393ba16761SJacob Faibussowitsch   if (size == 1) PetscFunctionReturn(PETSC_SUCCESS);
440d4002b98SHong Zhang 
441d4002b98SHong Zhang   /* Hard test: off-diagonal block. This takes a MatCreateSubMatrix. */
4429566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Amat, &M, &N));
4439566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &first, &last));
4449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N - last + first, &notme));
445d4002b98SHong Zhang   for (i = 0; i < first; i++) notme[i] = i;
446d4002b98SHong Zhang   for (i = last; i < M; i++) notme[i - last + first] = i;
4479566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(MPI_COMM_SELF, N - last + first, notme, PETSC_COPY_VALUES, &Notme));
4489566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(MPI_COMM_SELF, last - first, first, 1, &Me));
4499566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Amat, 1, &Me, &Notme, MAT_INITIAL_MATRIX, &Aoffs));
450d4002b98SHong Zhang   Aoff = Aoffs[0];
4519566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Bmat, 1, &Notme, &Me, MAT_INITIAL_MATRIX, &Boffs));
452d4002b98SHong Zhang   Boff = Boffs[0];
4539566063dSJacob Faibussowitsch   PetscCall(MatIsTranspose(Aoff, Boff, tol, f));
4549566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1, &Aoffs));
4559566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1, &Boffs));
4569566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Me));
4579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Notme));
4589566063dSJacob Faibussowitsch   PetscCall(PetscFree(notme));
4593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
460d4002b98SHong Zhang }
461d4002b98SHong Zhang 
462d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultTransposeAdd_MPISELL(Mat A, Vec xx, Vec yy, Vec zz)
463d71ae5a4SJacob Faibussowitsch {
464d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
465d4002b98SHong Zhang 
466d4002b98SHong Zhang   PetscFunctionBegin;
467d4002b98SHong Zhang   /* do nondiagonal part */
4689566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B, xx, a->lvec));
469d4002b98SHong Zhang   /* do local part */
4709566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multtransposeadd)(a->A, xx, yy, zz));
471e4a140f6SJunchao Zhang   /* add partial results together */
4729566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, a->lvec, zz, ADD_VALUES, SCATTER_REVERSE));
4739566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, a->lvec, zz, ADD_VALUES, SCATTER_REVERSE));
4743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
475d4002b98SHong Zhang }
476d4002b98SHong Zhang 
477d4002b98SHong Zhang /*
478d4002b98SHong Zhang   This only works correctly for square matrices where the subblock A->A is the
479d4002b98SHong Zhang    diagonal block
480d4002b98SHong Zhang */
481d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_MPISELL(Mat A, Vec v)
482d71ae5a4SJacob Faibussowitsch {
483d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
484d4002b98SHong Zhang 
485d4002b98SHong Zhang   PetscFunctionBegin;
48608401ef6SPierre Jolivet   PetscCheck(A->rmap->N == A->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Supports only square matrix where A->A is diag block");
487aed4548fSBarry 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");
4889566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(a->A, v));
4893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
490d4002b98SHong Zhang }
491d4002b98SHong Zhang 
492d71ae5a4SJacob Faibussowitsch PetscErrorCode MatScale_MPISELL(Mat A, PetscScalar aa)
493d71ae5a4SJacob Faibussowitsch {
494d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
495d4002b98SHong Zhang 
496d4002b98SHong Zhang   PetscFunctionBegin;
4979566063dSJacob Faibussowitsch   PetscCall(MatScale(a->A, aa));
4989566063dSJacob Faibussowitsch   PetscCall(MatScale(a->B, aa));
4993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
500d4002b98SHong Zhang }
501d4002b98SHong Zhang 
502d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_MPISELL(Mat mat)
503d71ae5a4SJacob Faibussowitsch {
504d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
505d4002b98SHong Zhang 
506d4002b98SHong Zhang   PetscFunctionBegin;
507d4002b98SHong Zhang #if defined(PETSC_USE_LOG)
5083ba16761SJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT, mat->rmap->N, mat->cmap->N));
509d4002b98SHong Zhang #endif
5109566063dSJacob Faibussowitsch   PetscCall(MatStashDestroy_Private(&mat->stash));
5119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->diag));
5129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sell->A));
5139566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sell->B));
514d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
515eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&sell->colmap));
516d4002b98SHong Zhang #else
5179566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->colmap));
518d4002b98SHong Zhang #endif
5199566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->garray));
5209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->lvec));
5219566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sell->Mvctx));
5229566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sell->rowvalues, sell->rowindices));
5239566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->ld));
5249566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
525d4002b98SHong Zhang 
5269566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL));
5279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatStoreValues_C", NULL));
5289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatRetrieveValues_C", NULL));
5299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatIsTranspose_C", NULL));
5309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPISELLSetPreallocation_C", NULL));
5319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisell_mpiaij_C", NULL));
532*b5917f1bSHong Zhang #if defined(PETSC_HAVE_CUDA)
533*b5917f1bSHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisell_mpisellcuda_C", NULL));
534*b5917f1bSHong Zhang #endif
5359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalScaleLocal_C", NULL));
5363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
537d4002b98SHong Zhang }
538d4002b98SHong Zhang 
539d4002b98SHong Zhang #include <petscdraw.h>
540d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_MPISELL_ASCIIorDraworSocket(Mat mat, PetscViewer viewer)
541d71ae5a4SJacob Faibussowitsch {
542d4002b98SHong Zhang   Mat_MPISELL      *sell = (Mat_MPISELL *)mat->data;
543d4002b98SHong Zhang   PetscMPIInt       rank = sell->rank, size = sell->size;
544d4002b98SHong Zhang   PetscBool         isdraw, iascii, isbinary;
545d4002b98SHong Zhang   PetscViewer       sviewer;
546d4002b98SHong Zhang   PetscViewerFormat format;
547d4002b98SHong Zhang 
548d4002b98SHong Zhang   PetscFunctionBegin;
5499566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
5509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
5519566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
552d4002b98SHong Zhang   if (iascii) {
5539566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
554d4002b98SHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
555d4002b98SHong Zhang       MatInfo   info;
5566335e310SSatish Balay       PetscInt *inodes;
557d4002b98SHong Zhang 
5589566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat), &rank));
5599566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(mat, MAT_LOCAL, &info));
5609566063dSJacob Faibussowitsch       PetscCall(MatInodeGetInodeSizes(sell->A, NULL, &inodes, NULL));
5619566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
562d4002b98SHong Zhang       if (!inodes) {
5639371c9d4SSatish 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,
5649371c9d4SSatish Balay                                                      (PetscInt)info.nz_allocated, (PetscInt)info.memory));
565d4002b98SHong Zhang       } else {
5669371c9d4SSatish 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,
5679371c9d4SSatish Balay                                                      (PetscInt)info.nz_allocated, (PetscInt)info.memory));
568d4002b98SHong Zhang       }
5699566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(sell->A, MAT_LOCAL, &info));
5709566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] on-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
5719566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(sell->B, MAT_LOCAL, &info));
5729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] off-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
5739566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
5749566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
5759566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Information on VecScatter used in matrix-vector product: \n"));
5769566063dSJacob Faibussowitsch       PetscCall(VecScatterView(sell->Mvctx, viewer));
5773ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
578d4002b98SHong Zhang     } else if (format == PETSC_VIEWER_ASCII_INFO) {
579d4002b98SHong Zhang       PetscInt inodecount, inodelimit, *inodes;
5809566063dSJacob Faibussowitsch       PetscCall(MatInodeGetInodeSizes(sell->A, &inodecount, &inodes, &inodelimit));
581d4002b98SHong Zhang       if (inodes) {
5829566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "using I-node (on process 0) routines: found %" PetscInt_FMT " nodes, limit used is %" PetscInt_FMT "\n", inodecount, inodelimit));
583d4002b98SHong Zhang       } else {
5849566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "not using I-node (on process 0) routines\n"));
585d4002b98SHong Zhang       }
5863ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
587d4002b98SHong Zhang     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
5883ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
589d4002b98SHong Zhang     }
590d4002b98SHong Zhang   } else if (isbinary) {
591d4002b98SHong Zhang     if (size == 1) {
5929566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)sell->A, ((PetscObject)mat)->name));
5939566063dSJacob Faibussowitsch       PetscCall(MatView(sell->A, viewer));
594d4002b98SHong Zhang     } else {
5959566063dSJacob Faibussowitsch       /* PetscCall(MatView_MPISELL_Binary(mat,viewer)); */
596d4002b98SHong Zhang     }
5973ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
598d4002b98SHong Zhang   } else if (isdraw) {
599d4002b98SHong Zhang     PetscDraw draw;
600d4002b98SHong Zhang     PetscBool isnull;
6019566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
6029566063dSJacob Faibussowitsch     PetscCall(PetscDrawIsNull(draw, &isnull));
6033ba16761SJacob Faibussowitsch     if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
604d4002b98SHong Zhang   }
605d4002b98SHong Zhang 
606d4002b98SHong Zhang   {
607d4002b98SHong Zhang     /* assemble the entire matrix onto first processor. */
608d4002b98SHong Zhang     Mat          A;
609d4002b98SHong Zhang     Mat_SeqSELL *Aloc;
610d4002b98SHong Zhang     PetscInt     M = mat->rmap->N, N = mat->cmap->N, *acolidx, row, col, i, j;
611d4002b98SHong Zhang     MatScalar   *aval;
612d4002b98SHong Zhang     PetscBool    isnonzero;
613d4002b98SHong Zhang 
6149566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &A));
615dd400576SPatrick Sanan     if (rank == 0) {
6169566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A, M, N, M, N));
617d4002b98SHong Zhang     } else {
6189566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A, 0, 0, M, N));
619d4002b98SHong Zhang     }
620d4002b98SHong Zhang     /* This is just a temporary matrix, so explicitly using MATMPISELL is probably best */
6219566063dSJacob Faibussowitsch     PetscCall(MatSetType(A, MATMPISELL));
6229566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(A, 0, NULL, 0, NULL));
6239566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
624d4002b98SHong Zhang 
625d4002b98SHong Zhang     /* copy over the A part */
626d4002b98SHong Zhang     Aloc    = (Mat_SeqSELL *)sell->A->data;
6279371c9d4SSatish Balay     acolidx = Aloc->colidx;
6289371c9d4SSatish Balay     aval    = Aloc->val;
629d4002b98SHong Zhang     for (i = 0; i < Aloc->totalslices; i++) { /* loop over slices */
630d4002b98SHong Zhang       for (j = Aloc->sliidx[i]; j < Aloc->sliidx[i + 1]; j++) {
63107e43b41SHong Zhang         isnonzero = (PetscBool)((j - Aloc->sliidx[i]) / Aloc->sliceheight < Aloc->rlen[i * Aloc->sliceheight + j % Aloc->sliceheight]);
632d4002b98SHong Zhang         if (isnonzero) { /* check the mask bit */
63307e43b41SHong Zhang           row = i * Aloc->sliceheight + j % Aloc->sliceheight + mat->rmap->rstart;
634d4002b98SHong Zhang           col = *acolidx + mat->rmap->rstart;
6359566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A, 1, &row, 1, &col, aval, INSERT_VALUES));
636d4002b98SHong Zhang         }
6379371c9d4SSatish Balay         aval++;
6389371c9d4SSatish Balay         acolidx++;
639d4002b98SHong Zhang       }
640d4002b98SHong Zhang     }
641d4002b98SHong Zhang 
642d4002b98SHong Zhang     /* copy over the B part */
643d4002b98SHong Zhang     Aloc    = (Mat_SeqSELL *)sell->B->data;
6449371c9d4SSatish Balay     acolidx = Aloc->colidx;
6459371c9d4SSatish Balay     aval    = Aloc->val;
646d4002b98SHong Zhang     for (i = 0; i < Aloc->totalslices; i++) {
647d4002b98SHong Zhang       for (j = Aloc->sliidx[i]; j < Aloc->sliidx[i + 1]; j++) {
64807e43b41SHong Zhang         isnonzero = (PetscBool)((j - Aloc->sliidx[i]) / Aloc->sliceheight < Aloc->rlen[i * Aloc->sliceheight + j % Aloc->sliceheight]);
649d4002b98SHong Zhang         if (isnonzero) {
65007e43b41SHong Zhang           row = i * Aloc->sliceheight + j % Aloc->sliceheight + mat->rmap->rstart;
651d4002b98SHong Zhang           col = sell->garray[*acolidx];
6529566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A, 1, &row, 1, &col, aval, INSERT_VALUES));
653d4002b98SHong Zhang         }
6549371c9d4SSatish Balay         aval++;
6559371c9d4SSatish Balay         acolidx++;
656d4002b98SHong Zhang       }
657d4002b98SHong Zhang     }
658d4002b98SHong Zhang 
6599566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
6609566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
661d4002b98SHong Zhang     /*
662d4002b98SHong Zhang        Everyone has to call to draw the matrix since the graphics waits are
663d4002b98SHong Zhang        synchronized across all processors that share the PetscDraw object
664d4002b98SHong Zhang     */
6659566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
666dd400576SPatrick Sanan     if (rank == 0) {
6679566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)((Mat_MPISELL *)(A->data))->A, ((PetscObject)mat)->name));
6689566063dSJacob Faibussowitsch       PetscCall(MatView_SeqSELL(((Mat_MPISELL *)(A->data))->A, sviewer));
669d4002b98SHong Zhang     }
6709566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
6719566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
6729566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
673d4002b98SHong Zhang   }
6743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
675d4002b98SHong Zhang }
676d4002b98SHong Zhang 
677d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_MPISELL(Mat mat, PetscViewer viewer)
678d71ae5a4SJacob Faibussowitsch {
679d4002b98SHong Zhang   PetscBool iascii, isdraw, issocket, isbinary;
680d4002b98SHong Zhang 
681d4002b98SHong Zhang   PetscFunctionBegin;
6829566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
6839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
6849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
6859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSOCKET, &issocket));
68648a46eb9SPierre Jolivet   if (iascii || isdraw || isbinary || issocket) PetscCall(MatView_MPISELL_ASCIIorDraworSocket(mat, viewer));
6873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
688d4002b98SHong Zhang }
689d4002b98SHong Zhang 
690d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetGhosts_MPISELL(Mat mat, PetscInt *nghosts, const PetscInt *ghosts[])
691d71ae5a4SJacob Faibussowitsch {
692d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
693d4002b98SHong Zhang 
694d4002b98SHong Zhang   PetscFunctionBegin;
6959566063dSJacob Faibussowitsch   PetscCall(MatGetSize(sell->B, NULL, nghosts));
696d4002b98SHong Zhang   if (ghosts) *ghosts = sell->garray;
6973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
698d4002b98SHong Zhang }
699d4002b98SHong Zhang 
700d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetInfo_MPISELL(Mat matin, MatInfoType flag, MatInfo *info)
701d71ae5a4SJacob Faibussowitsch {
702d4002b98SHong Zhang   Mat_MPISELL   *mat = (Mat_MPISELL *)matin->data;
703d4002b98SHong Zhang   Mat            A = mat->A, B = mat->B;
7043966268fSBarry Smith   PetscLogDouble isend[5], irecv[5];
705d4002b98SHong Zhang 
706d4002b98SHong Zhang   PetscFunctionBegin;
707d4002b98SHong Zhang   info->block_size = 1.0;
7089566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_LOCAL, info));
709d4002b98SHong Zhang 
7109371c9d4SSatish Balay   isend[0] = info->nz_used;
7119371c9d4SSatish Balay   isend[1] = info->nz_allocated;
7129371c9d4SSatish Balay   isend[2] = info->nz_unneeded;
7139371c9d4SSatish Balay   isend[3] = info->memory;
7149371c9d4SSatish Balay   isend[4] = info->mallocs;
715d4002b98SHong Zhang 
7169566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(B, MAT_LOCAL, info));
717d4002b98SHong Zhang 
7189371c9d4SSatish Balay   isend[0] += info->nz_used;
7199371c9d4SSatish Balay   isend[1] += info->nz_allocated;
7209371c9d4SSatish Balay   isend[2] += info->nz_unneeded;
7219371c9d4SSatish Balay   isend[3] += info->memory;
7229371c9d4SSatish Balay   isend[4] += info->mallocs;
723d4002b98SHong Zhang   if (flag == MAT_LOCAL) {
724d4002b98SHong Zhang     info->nz_used      = isend[0];
725d4002b98SHong Zhang     info->nz_allocated = isend[1];
726d4002b98SHong Zhang     info->nz_unneeded  = isend[2];
727d4002b98SHong Zhang     info->memory       = isend[3];
728d4002b98SHong Zhang     info->mallocs      = isend[4];
729d4002b98SHong Zhang   } else if (flag == MAT_GLOBAL_MAX) {
7301c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_MAX, PetscObjectComm((PetscObject)matin)));
731d4002b98SHong Zhang 
732d4002b98SHong Zhang     info->nz_used      = irecv[0];
733d4002b98SHong Zhang     info->nz_allocated = irecv[1];
734d4002b98SHong Zhang     info->nz_unneeded  = irecv[2];
735d4002b98SHong Zhang     info->memory       = irecv[3];
736d4002b98SHong Zhang     info->mallocs      = irecv[4];
737d4002b98SHong Zhang   } else if (flag == MAT_GLOBAL_SUM) {
7381c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_SUM, PetscObjectComm((PetscObject)matin)));
739d4002b98SHong Zhang 
740d4002b98SHong Zhang     info->nz_used      = irecv[0];
741d4002b98SHong Zhang     info->nz_allocated = irecv[1];
742d4002b98SHong Zhang     info->nz_unneeded  = irecv[2];
743d4002b98SHong Zhang     info->memory       = irecv[3];
744d4002b98SHong Zhang     info->mallocs      = irecv[4];
745d4002b98SHong Zhang   }
746d4002b98SHong Zhang   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
747d4002b98SHong Zhang   info->fill_ratio_needed = 0;
748d4002b98SHong Zhang   info->factor_mallocs    = 0;
7493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
750d4002b98SHong Zhang }
751d4002b98SHong Zhang 
752d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetOption_MPISELL(Mat A, MatOption op, PetscBool flg)
753d71ae5a4SJacob Faibussowitsch {
754d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
755d4002b98SHong Zhang 
756d4002b98SHong Zhang   PetscFunctionBegin;
757d4002b98SHong Zhang   switch (op) {
758d4002b98SHong Zhang   case MAT_NEW_NONZERO_LOCATIONS:
759d4002b98SHong Zhang   case MAT_NEW_NONZERO_ALLOCATION_ERR:
760d4002b98SHong Zhang   case MAT_UNUSED_NONZERO_LOCATION_ERR:
761d4002b98SHong Zhang   case MAT_KEEP_NONZERO_PATTERN:
762d4002b98SHong Zhang   case MAT_NEW_NONZERO_LOCATION_ERR:
763d4002b98SHong Zhang   case MAT_USE_INODES:
764d4002b98SHong Zhang   case MAT_IGNORE_ZERO_ENTRIES:
765d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7669566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
7679566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B, op, flg));
768d4002b98SHong Zhang     break;
769d4002b98SHong Zhang   case MAT_ROW_ORIENTED:
770d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
771d4002b98SHong Zhang     a->roworiented = flg;
772d4002b98SHong Zhang 
7739566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
7749566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B, op, flg));
775d4002b98SHong Zhang     break;
7768c78258cSHong Zhang   case MAT_FORCE_DIAGONAL_ENTRIES:
777d71ae5a4SJacob Faibussowitsch   case MAT_SORTED_FULL:
778d71ae5a4SJacob Faibussowitsch     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
779d71ae5a4SJacob Faibussowitsch     break;
780d71ae5a4SJacob Faibussowitsch   case MAT_IGNORE_OFF_PROC_ENTRIES:
781d71ae5a4SJacob Faibussowitsch     a->donotstash = flg;
782d71ae5a4SJacob Faibussowitsch     break;
783d4002b98SHong Zhang   case MAT_SPD:
784d71ae5a4SJacob Faibussowitsch   case MAT_SPD_ETERNAL:
785d71ae5a4SJacob Faibussowitsch     break;
786d4002b98SHong Zhang   case MAT_SYMMETRIC:
787d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7889566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
789d4002b98SHong Zhang     break;
790d4002b98SHong Zhang   case MAT_STRUCTURALLY_SYMMETRIC:
791d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7929566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
793d4002b98SHong Zhang     break;
794d4002b98SHong Zhang   case MAT_HERMITIAN:
795d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7969566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
797d4002b98SHong Zhang     break;
798d4002b98SHong Zhang   case MAT_SYMMETRY_ETERNAL:
799d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
8009566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
801d4002b98SHong Zhang     break;
802b94d7dedSBarry Smith   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
803b94d7dedSBarry Smith     MatCheckPreallocated(A, 1);
804b94d7dedSBarry Smith     PetscCall(MatSetOption(a->A, op, flg));
805b94d7dedSBarry Smith     break;
806d71ae5a4SJacob Faibussowitsch   default:
807d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unknown option %d", op);
808d4002b98SHong Zhang   }
8093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
810d4002b98SHong Zhang }
811d4002b98SHong Zhang 
812d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDiagonalScale_MPISELL(Mat mat, Vec ll, Vec rr)
813d71ae5a4SJacob Faibussowitsch {
814d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
815d4002b98SHong Zhang   Mat          a = sell->A, b = sell->B;
816d4002b98SHong Zhang   PetscInt     s1, s2, s3;
817d4002b98SHong Zhang 
818d4002b98SHong Zhang   PetscFunctionBegin;
8199566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mat, &s2, &s3));
820d4002b98SHong Zhang   if (rr) {
8219566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rr, &s1));
82208401ef6SPierre Jolivet     PetscCheck(s1 == s3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "right vector non-conforming local size");
823d4002b98SHong Zhang     /* Overlap communication with computation. */
8249566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(sell->Mvctx, rr, sell->lvec, INSERT_VALUES, SCATTER_FORWARD));
825d4002b98SHong Zhang   }
826d4002b98SHong Zhang   if (ll) {
8279566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(ll, &s1));
82808401ef6SPierre Jolivet     PetscCheck(s1 == s2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "left vector non-conforming local size");
829dbbe0bcdSBarry Smith     PetscUseTypeMethod(b, diagonalscale, ll, NULL);
830d4002b98SHong Zhang   }
831d4002b98SHong Zhang   /* scale  the diagonal block */
832dbbe0bcdSBarry Smith   PetscUseTypeMethod(a, diagonalscale, ll, rr);
833d4002b98SHong Zhang 
834d4002b98SHong Zhang   if (rr) {
835d4002b98SHong Zhang     /* Do a scatter end and then right scale the off-diagonal block */
8369566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(sell->Mvctx, rr, sell->lvec, INSERT_VALUES, SCATTER_FORWARD));
837dbbe0bcdSBarry Smith     PetscUseTypeMethod(b, diagonalscale, NULL, sell->lvec);
838d4002b98SHong Zhang   }
8393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
840d4002b98SHong Zhang }
841d4002b98SHong Zhang 
842d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUnfactored_MPISELL(Mat A)
843d71ae5a4SJacob Faibussowitsch {
844d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
845d4002b98SHong Zhang 
846d4002b98SHong Zhang   PetscFunctionBegin;
8479566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(a->A));
8483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
849d4002b98SHong Zhang }
850d4002b98SHong Zhang 
851d71ae5a4SJacob Faibussowitsch PetscErrorCode MatEqual_MPISELL(Mat A, Mat B, PetscBool *flag)
852d71ae5a4SJacob Faibussowitsch {
853d4002b98SHong Zhang   Mat_MPISELL *matB = (Mat_MPISELL *)B->data, *matA = (Mat_MPISELL *)A->data;
854d4002b98SHong Zhang   Mat          a, b, c, d;
855d4002b98SHong Zhang   PetscBool    flg;
856d4002b98SHong Zhang 
857d4002b98SHong Zhang   PetscFunctionBegin;
8589371c9d4SSatish Balay   a = matA->A;
8599371c9d4SSatish Balay   b = matA->B;
8609371c9d4SSatish Balay   c = matB->A;
8619371c9d4SSatish Balay   d = matB->B;
862d4002b98SHong Zhang 
8639566063dSJacob Faibussowitsch   PetscCall(MatEqual(a, c, &flg));
86448a46eb9SPierre Jolivet   if (flg) PetscCall(MatEqual(b, d, &flg));
8651c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&flg, flag, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
8663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
867d4002b98SHong Zhang }
868d4002b98SHong Zhang 
869d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCopy_MPISELL(Mat A, Mat B, MatStructure str)
870d71ae5a4SJacob Faibussowitsch {
871d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
872d4002b98SHong Zhang   Mat_MPISELL *b = (Mat_MPISELL *)B->data;
873d4002b98SHong Zhang 
874d4002b98SHong Zhang   PetscFunctionBegin;
875d4002b98SHong Zhang   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
876d4002b98SHong Zhang   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
877d4002b98SHong Zhang     /* because of the column compression in the off-processor part of the matrix a->B,
878d4002b98SHong Zhang        the number of columns in a->B and b->B may be different, hence we cannot call
879d4002b98SHong Zhang        the MatCopy() directly on the two parts. If need be, we can provide a more
880d4002b98SHong Zhang        efficient copy than the MatCopy_Basic() by first uncompressing the a->B matrices
881d4002b98SHong Zhang        then copying the submatrices */
8829566063dSJacob Faibussowitsch     PetscCall(MatCopy_Basic(A, B, str));
883d4002b98SHong Zhang   } else {
8849566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->A, b->A, str));
8859566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->B, b->B, str));
886d4002b98SHong Zhang   }
8873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
888d4002b98SHong Zhang }
889d4002b98SHong Zhang 
890d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUp_MPISELL(Mat A)
891d71ae5a4SJacob Faibussowitsch {
892d4002b98SHong Zhang   PetscFunctionBegin;
8939566063dSJacob Faibussowitsch   PetscCall(MatMPISELLSetPreallocation(A, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
8943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
895d4002b98SHong Zhang }
896d4002b98SHong Zhang 
897d4002b98SHong Zhang extern PetscErrorCode MatConjugate_SeqSELL(Mat);
898d4002b98SHong Zhang 
899d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConjugate_MPISELL(Mat mat)
900d71ae5a4SJacob Faibussowitsch {
9015f80ce2aSJacob Faibussowitsch   PetscFunctionBegin;
9025f80ce2aSJacob Faibussowitsch   if (PetscDefined(USE_COMPLEX)) {
903d4002b98SHong Zhang     Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
904d4002b98SHong Zhang 
9059566063dSJacob Faibussowitsch     PetscCall(MatConjugate_SeqSELL(sell->A));
9069566063dSJacob Faibussowitsch     PetscCall(MatConjugate_SeqSELL(sell->B));
9075f80ce2aSJacob Faibussowitsch   }
9083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
909d4002b98SHong Zhang }
910d4002b98SHong Zhang 
911d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRealPart_MPISELL(Mat A)
912d71ae5a4SJacob Faibussowitsch {
913d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
914d4002b98SHong Zhang 
915d4002b98SHong Zhang   PetscFunctionBegin;
9169566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->A));
9179566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->B));
9183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
919d4002b98SHong Zhang }
920d4002b98SHong Zhang 
921d71ae5a4SJacob Faibussowitsch PetscErrorCode MatImaginaryPart_MPISELL(Mat A)
922d71ae5a4SJacob Faibussowitsch {
923d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
924d4002b98SHong Zhang 
925d4002b98SHong Zhang   PetscFunctionBegin;
9269566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->A));
9279566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->B));
9283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
929d4002b98SHong Zhang }
930d4002b98SHong Zhang 
931d71ae5a4SJacob Faibussowitsch PetscErrorCode MatInvertBlockDiagonal_MPISELL(Mat A, const PetscScalar **values)
932d71ae5a4SJacob Faibussowitsch {
933d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
934d4002b98SHong Zhang 
935d4002b98SHong Zhang   PetscFunctionBegin;
9369566063dSJacob Faibussowitsch   PetscCall(MatInvertBlockDiagonal(a->A, values));
937d4002b98SHong Zhang   A->factorerrortype = a->A->factorerrortype;
9383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
939d4002b98SHong Zhang }
940d4002b98SHong Zhang 
941d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetRandom_MPISELL(Mat x, PetscRandom rctx)
942d71ae5a4SJacob Faibussowitsch {
943d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)x->data;
944d4002b98SHong Zhang 
945d4002b98SHong Zhang   PetscFunctionBegin;
9469566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(sell->A, rctx));
9479566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(sell->B, rctx));
9489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
9499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
9503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
951d4002b98SHong Zhang }
952d4002b98SHong Zhang 
953d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetFromOptions_MPISELL(Mat A, PetscOptionItems *PetscOptionsObject)
954d71ae5a4SJacob Faibussowitsch {
955d4002b98SHong Zhang   PetscFunctionBegin;
956d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "MPISELL options");
957d0609cedSBarry Smith   PetscOptionsHeadEnd();
9583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
959d4002b98SHong Zhang }
960d4002b98SHong Zhang 
961d71ae5a4SJacob Faibussowitsch PetscErrorCode MatShift_MPISELL(Mat Y, PetscScalar a)
962d71ae5a4SJacob Faibussowitsch {
963d4002b98SHong Zhang   Mat_MPISELL *msell = (Mat_MPISELL *)Y->data;
964d4002b98SHong Zhang   Mat_SeqSELL *sell  = (Mat_SeqSELL *)msell->A->data;
965d4002b98SHong Zhang 
966d4002b98SHong Zhang   PetscFunctionBegin;
967d4002b98SHong Zhang   if (!Y->preallocated) {
9689566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(Y, 1, NULL, 0, NULL));
969d4002b98SHong Zhang   } else if (!sell->nz) {
970d4002b98SHong Zhang     PetscInt nonew = sell->nonew;
9719566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(msell->A, 1, NULL));
972d4002b98SHong Zhang     sell->nonew = nonew;
973d4002b98SHong Zhang   }
9749566063dSJacob Faibussowitsch   PetscCall(MatShift_Basic(Y, a));
9753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
976d4002b98SHong Zhang }
977d4002b98SHong Zhang 
978d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMissingDiagonal_MPISELL(Mat A, PetscBool *missing, PetscInt *d)
979d71ae5a4SJacob Faibussowitsch {
980d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
981d4002b98SHong Zhang 
982d4002b98SHong Zhang   PetscFunctionBegin;
98308401ef6SPierre Jolivet   PetscCheck(A->rmap->n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only works for square matrices");
9849566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(a->A, missing, d));
985d4002b98SHong Zhang   if (d) {
986d4002b98SHong Zhang     PetscInt rstart;
9879566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
988d4002b98SHong Zhang     *d += rstart;
989d4002b98SHong Zhang   }
9903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
991d4002b98SHong Zhang }
992d4002b98SHong Zhang 
993d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonalBlock_MPISELL(Mat A, Mat *a)
994d71ae5a4SJacob Faibussowitsch {
995d4002b98SHong Zhang   PetscFunctionBegin;
996d4002b98SHong Zhang   *a = ((Mat_MPISELL *)A->data)->A;
9973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
998d4002b98SHong Zhang }
999d4002b98SHong Zhang 
1000d4002b98SHong Zhang static struct _MatOps MatOps_Values = {MatSetValues_MPISELL,
1001f4259b30SLisandro Dalcin                                        NULL,
1002f4259b30SLisandro Dalcin                                        NULL,
1003d4002b98SHong Zhang                                        MatMult_MPISELL,
1004d4002b98SHong Zhang                                        /* 4*/ MatMultAdd_MPISELL,
1005d4002b98SHong Zhang                                        MatMultTranspose_MPISELL,
1006d4002b98SHong Zhang                                        MatMultTransposeAdd_MPISELL,
1007f4259b30SLisandro Dalcin                                        NULL,
1008f4259b30SLisandro Dalcin                                        NULL,
1009f4259b30SLisandro Dalcin                                        NULL,
1010f4259b30SLisandro Dalcin                                        /*10*/ NULL,
1011f4259b30SLisandro Dalcin                                        NULL,
1012f4259b30SLisandro Dalcin                                        NULL,
1013d4002b98SHong Zhang                                        MatSOR_MPISELL,
1014f4259b30SLisandro Dalcin                                        NULL,
1015d4002b98SHong Zhang                                        /*15*/ MatGetInfo_MPISELL,
1016d4002b98SHong Zhang                                        MatEqual_MPISELL,
1017d4002b98SHong Zhang                                        MatGetDiagonal_MPISELL,
1018d4002b98SHong Zhang                                        MatDiagonalScale_MPISELL,
1019f4259b30SLisandro Dalcin                                        NULL,
1020d4002b98SHong Zhang                                        /*20*/ MatAssemblyBegin_MPISELL,
1021d4002b98SHong Zhang                                        MatAssemblyEnd_MPISELL,
1022d4002b98SHong Zhang                                        MatSetOption_MPISELL,
1023d4002b98SHong Zhang                                        MatZeroEntries_MPISELL,
1024f4259b30SLisandro Dalcin                                        /*24*/ NULL,
1025f4259b30SLisandro Dalcin                                        NULL,
1026f4259b30SLisandro Dalcin                                        NULL,
1027f4259b30SLisandro Dalcin                                        NULL,
1028f4259b30SLisandro Dalcin                                        NULL,
1029d4002b98SHong Zhang                                        /*29*/ MatSetUp_MPISELL,
1030f4259b30SLisandro Dalcin                                        NULL,
1031f4259b30SLisandro Dalcin                                        NULL,
1032d4002b98SHong Zhang                                        MatGetDiagonalBlock_MPISELL,
1033f4259b30SLisandro Dalcin                                        NULL,
1034d4002b98SHong Zhang                                        /*34*/ MatDuplicate_MPISELL,
1035f4259b30SLisandro Dalcin                                        NULL,
1036f4259b30SLisandro Dalcin                                        NULL,
1037f4259b30SLisandro Dalcin                                        NULL,
1038f4259b30SLisandro Dalcin                                        NULL,
1039f4259b30SLisandro Dalcin                                        /*39*/ NULL,
1040f4259b30SLisandro Dalcin                                        NULL,
1041f4259b30SLisandro Dalcin                                        NULL,
1042d4002b98SHong Zhang                                        MatGetValues_MPISELL,
1043d4002b98SHong Zhang                                        MatCopy_MPISELL,
1044f4259b30SLisandro Dalcin                                        /*44*/ NULL,
1045d4002b98SHong Zhang                                        MatScale_MPISELL,
1046d4002b98SHong Zhang                                        MatShift_MPISELL,
1047d4002b98SHong Zhang                                        MatDiagonalSet_MPISELL,
1048f4259b30SLisandro Dalcin                                        NULL,
1049d4002b98SHong Zhang                                        /*49*/ MatSetRandom_MPISELL,
1050f4259b30SLisandro Dalcin                                        NULL,
1051f4259b30SLisandro Dalcin                                        NULL,
1052f4259b30SLisandro Dalcin                                        NULL,
1053f4259b30SLisandro Dalcin                                        NULL,
1054d4002b98SHong Zhang                                        /*54*/ MatFDColoringCreate_MPIXAIJ,
1055f4259b30SLisandro Dalcin                                        NULL,
1056d4002b98SHong Zhang                                        MatSetUnfactored_MPISELL,
1057f4259b30SLisandro Dalcin                                        NULL,
1058f4259b30SLisandro Dalcin                                        NULL,
1059f4259b30SLisandro Dalcin                                        /*59*/ NULL,
1060d4002b98SHong Zhang                                        MatDestroy_MPISELL,
1061d4002b98SHong Zhang                                        MatView_MPISELL,
1062f4259b30SLisandro Dalcin                                        NULL,
1063f4259b30SLisandro Dalcin                                        NULL,
1064f4259b30SLisandro Dalcin                                        /*64*/ NULL,
1065f4259b30SLisandro Dalcin                                        NULL,
1066f4259b30SLisandro Dalcin                                        NULL,
1067f4259b30SLisandro Dalcin                                        NULL,
1068f4259b30SLisandro Dalcin                                        NULL,
1069f4259b30SLisandro Dalcin                                        /*69*/ NULL,
1070f4259b30SLisandro Dalcin                                        NULL,
1071f4259b30SLisandro Dalcin                                        NULL,
1072f4259b30SLisandro Dalcin                                        NULL,
1073f4259b30SLisandro Dalcin                                        NULL,
1074f4259b30SLisandro Dalcin                                        NULL,
1075d4002b98SHong Zhang                                        /*75*/ MatFDColoringApply_AIJ, /* reuse AIJ function */
1076d4002b98SHong Zhang                                        MatSetFromOptions_MPISELL,
1077f4259b30SLisandro Dalcin                                        NULL,
1078f4259b30SLisandro Dalcin                                        NULL,
1079f4259b30SLisandro Dalcin                                        NULL,
1080f4259b30SLisandro Dalcin                                        /*80*/ NULL,
1081f4259b30SLisandro Dalcin                                        NULL,
1082f4259b30SLisandro Dalcin                                        NULL,
1083f4259b30SLisandro Dalcin                                        /*83*/ NULL,
1084f4259b30SLisandro Dalcin                                        NULL,
1085f4259b30SLisandro Dalcin                                        NULL,
1086f4259b30SLisandro Dalcin                                        NULL,
1087f4259b30SLisandro Dalcin                                        NULL,
1088f4259b30SLisandro Dalcin                                        NULL,
1089f4259b30SLisandro Dalcin                                        /*89*/ NULL,
1090f4259b30SLisandro Dalcin                                        NULL,
1091f4259b30SLisandro Dalcin                                        NULL,
1092f4259b30SLisandro Dalcin                                        NULL,
1093f4259b30SLisandro Dalcin                                        NULL,
1094f4259b30SLisandro Dalcin                                        /*94*/ NULL,
1095f4259b30SLisandro Dalcin                                        NULL,
1096f4259b30SLisandro Dalcin                                        NULL,
1097f4259b30SLisandro Dalcin                                        NULL,
1098f4259b30SLisandro Dalcin                                        NULL,
1099f4259b30SLisandro Dalcin                                        /*99*/ NULL,
1100f4259b30SLisandro Dalcin                                        NULL,
1101f4259b30SLisandro Dalcin                                        NULL,
1102d4002b98SHong Zhang                                        MatConjugate_MPISELL,
1103f4259b30SLisandro Dalcin                                        NULL,
1104f4259b30SLisandro Dalcin                                        /*104*/ NULL,
1105d4002b98SHong Zhang                                        MatRealPart_MPISELL,
1106d4002b98SHong Zhang                                        MatImaginaryPart_MPISELL,
1107f4259b30SLisandro Dalcin                                        NULL,
1108f4259b30SLisandro Dalcin                                        NULL,
1109f4259b30SLisandro Dalcin                                        /*109*/ NULL,
1110f4259b30SLisandro Dalcin                                        NULL,
1111f4259b30SLisandro Dalcin                                        NULL,
1112f4259b30SLisandro Dalcin                                        NULL,
1113d4002b98SHong Zhang                                        MatMissingDiagonal_MPISELL,
1114f4259b30SLisandro Dalcin                                        /*114*/ NULL,
1115f4259b30SLisandro Dalcin                                        NULL,
1116d4002b98SHong Zhang                                        MatGetGhosts_MPISELL,
1117f4259b30SLisandro Dalcin                                        NULL,
1118f4259b30SLisandro Dalcin                                        NULL,
1119f4259b30SLisandro Dalcin                                        /*119*/ NULL,
1120f4259b30SLisandro Dalcin                                        NULL,
1121f4259b30SLisandro Dalcin                                        NULL,
1122f4259b30SLisandro Dalcin                                        NULL,
1123f4259b30SLisandro Dalcin                                        NULL,
1124f4259b30SLisandro Dalcin                                        /*124*/ NULL,
1125f4259b30SLisandro Dalcin                                        NULL,
1126d4002b98SHong Zhang                                        MatInvertBlockDiagonal_MPISELL,
1127f4259b30SLisandro Dalcin                                        NULL,
1128f4259b30SLisandro Dalcin                                        NULL,
1129f4259b30SLisandro Dalcin                                        /*129*/ NULL,
1130f4259b30SLisandro Dalcin                                        NULL,
1131f4259b30SLisandro Dalcin                                        NULL,
1132f4259b30SLisandro Dalcin                                        NULL,
1133f4259b30SLisandro Dalcin                                        NULL,
1134f4259b30SLisandro Dalcin                                        /*134*/ NULL,
1135f4259b30SLisandro Dalcin                                        NULL,
1136f4259b30SLisandro Dalcin                                        NULL,
1137f4259b30SLisandro Dalcin                                        NULL,
1138f4259b30SLisandro Dalcin                                        NULL,
1139f4259b30SLisandro Dalcin                                        /*139*/ NULL,
1140f4259b30SLisandro Dalcin                                        NULL,
1141f4259b30SLisandro Dalcin                                        NULL,
1142d4002b98SHong Zhang                                        MatFDColoringSetUp_MPIXAIJ,
1143f4259b30SLisandro Dalcin                                        NULL,
1144d70f29a3SPierre Jolivet                                        /*144*/ NULL,
1145d70f29a3SPierre Jolivet                                        NULL,
1146d70f29a3SPierre Jolivet                                        NULL,
114799a7f59eSMark Adams                                        NULL,
114899a7f59eSMark Adams                                        NULL,
11497fb60732SBarry Smith                                        NULL,
1150dec0b466SHong Zhang                                        /*150*/ NULL,
1151dec0b466SHong Zhang                                        NULL};
1152d4002b98SHong Zhang 
1153d71ae5a4SJacob Faibussowitsch PetscErrorCode MatStoreValues_MPISELL(Mat mat)
1154d71ae5a4SJacob Faibussowitsch {
1155d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
1156d4002b98SHong Zhang 
1157d4002b98SHong Zhang   PetscFunctionBegin;
11589566063dSJacob Faibussowitsch   PetscCall(MatStoreValues(sell->A));
11599566063dSJacob Faibussowitsch   PetscCall(MatStoreValues(sell->B));
11603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1161d4002b98SHong Zhang }
1162d4002b98SHong Zhang 
1163d71ae5a4SJacob Faibussowitsch PetscErrorCode MatRetrieveValues_MPISELL(Mat mat)
1164d71ae5a4SJacob Faibussowitsch {
1165d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
1166d4002b98SHong Zhang 
1167d4002b98SHong Zhang   PetscFunctionBegin;
11689566063dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(sell->A));
11699566063dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(sell->B));
11703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1171d4002b98SHong Zhang }
1172d4002b98SHong Zhang 
1173d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPISELLSetPreallocation_MPISELL(Mat B, PetscInt d_rlenmax, const PetscInt d_rlen[], PetscInt o_rlenmax, const PetscInt o_rlen[])
1174d71ae5a4SJacob Faibussowitsch {
1175d4002b98SHong Zhang   Mat_MPISELL *b;
1176d4002b98SHong Zhang 
1177d4002b98SHong Zhang   PetscFunctionBegin;
11789566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
11799566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
1180d4002b98SHong Zhang   b = (Mat_MPISELL *)B->data;
1181d4002b98SHong Zhang 
1182d4002b98SHong Zhang   if (!B->preallocated) {
1183d4002b98SHong Zhang     /* Explicitly create 2 MATSEQSELL matrices. */
11849566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &b->A));
11859566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(b->A, B->rmap->n, B->cmap->n, B->rmap->n, B->cmap->n));
11869566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(b->A, B, B));
11879566063dSJacob Faibussowitsch     PetscCall(MatSetType(b->A, MATSEQSELL));
11889566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &b->B));
11899566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(b->B, B->rmap->n, B->cmap->N, B->rmap->n, B->cmap->N));
11909566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(b->B, B, B));
11919566063dSJacob Faibussowitsch     PetscCall(MatSetType(b->B, MATSEQSELL));
1192d4002b98SHong Zhang   }
1193d4002b98SHong Zhang 
11949566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(b->A, d_rlenmax, d_rlen));
11959566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(b->B, o_rlenmax, o_rlen));
1196d4002b98SHong Zhang   B->preallocated  = PETSC_TRUE;
1197d4002b98SHong Zhang   B->was_assembled = PETSC_FALSE;
1198d4002b98SHong Zhang   /*
1199d4002b98SHong Zhang     critical for MatAssemblyEnd to work.
1200d4002b98SHong Zhang     MatAssemblyBegin checks it to set up was_assembled
1201d4002b98SHong Zhang     and MatAssemblyEnd checks was_assembled to determine whether to build garray
1202d4002b98SHong Zhang   */
1203d4002b98SHong Zhang   B->assembled = PETSC_FALSE;
12043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1205d4002b98SHong Zhang }
1206d4002b98SHong Zhang 
1207d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDuplicate_MPISELL(Mat matin, MatDuplicateOption cpvalues, Mat *newmat)
1208d71ae5a4SJacob Faibussowitsch {
1209d4002b98SHong Zhang   Mat          mat;
1210d4002b98SHong Zhang   Mat_MPISELL *a, *oldmat = (Mat_MPISELL *)matin->data;
1211d4002b98SHong Zhang 
1212d4002b98SHong Zhang   PetscFunctionBegin;
1213f4259b30SLisandro Dalcin   *newmat = NULL;
12149566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)matin), &mat));
12159566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, matin->rmap->n, matin->cmap->n, matin->rmap->N, matin->cmap->N));
12169566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizesFromMats(mat, matin, matin));
12179566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, ((PetscObject)matin)->type_name));
1218d4002b98SHong Zhang   a = (Mat_MPISELL *)mat->data;
1219d4002b98SHong Zhang 
1220d4002b98SHong Zhang   mat->factortype   = matin->factortype;
1221d4002b98SHong Zhang   mat->assembled    = PETSC_TRUE;
1222d4002b98SHong Zhang   mat->insertmode   = NOT_SET_VALUES;
1223d4002b98SHong Zhang   mat->preallocated = PETSC_TRUE;
1224d4002b98SHong Zhang 
1225d4002b98SHong Zhang   a->size         = oldmat->size;
1226d4002b98SHong Zhang   a->rank         = oldmat->rank;
1227d4002b98SHong Zhang   a->donotstash   = oldmat->donotstash;
1228d4002b98SHong Zhang   a->roworiented  = oldmat->roworiented;
1229f4259b30SLisandro Dalcin   a->rowindices   = NULL;
1230f4259b30SLisandro Dalcin   a->rowvalues    = NULL;
1231d4002b98SHong Zhang   a->getrowactive = PETSC_FALSE;
1232d4002b98SHong Zhang 
12339566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->rmap, &mat->rmap));
12349566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->cmap, &mat->cmap));
1235d4002b98SHong Zhang 
1236d4002b98SHong Zhang   if (oldmat->colmap) {
1237d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
1238eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIDuplicate(oldmat->colmap, &a->colmap));
1239d4002b98SHong Zhang #else
12409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mat->cmap->N, &a->colmap));
12419566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->colmap, oldmat->colmap, mat->cmap->N));
1242d4002b98SHong Zhang #endif
1243f4259b30SLisandro Dalcin   } else a->colmap = NULL;
1244d4002b98SHong Zhang   if (oldmat->garray) {
1245d4002b98SHong Zhang     PetscInt len;
1246d4002b98SHong Zhang     len = oldmat->B->cmap->n;
12479566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len + 1, &a->garray));
12489566063dSJacob Faibussowitsch     if (len) PetscCall(PetscArraycpy(a->garray, oldmat->garray, len));
1249f4259b30SLisandro Dalcin   } else a->garray = NULL;
1250d4002b98SHong Zhang 
12519566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(oldmat->lvec, &a->lvec));
12529566063dSJacob Faibussowitsch   PetscCall(VecScatterCopy(oldmat->Mvctx, &a->Mvctx));
12539566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->A, cpvalues, &a->A));
12549566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->B, cpvalues, &a->B));
12559566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDuplicate(((PetscObject)matin)->qlist, &((PetscObject)mat)->qlist));
1256d4002b98SHong Zhang   *newmat = mat;
12573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1258d4002b98SHong Zhang }
1259d4002b98SHong Zhang 
1260d4002b98SHong Zhang /*@C
126111a5261eSBarry Smith    MatMPISELLSetPreallocation - Preallocates memory for a `MATMPISELL` sparse parallel matrix in sell format.
1262d4002b98SHong Zhang    For good matrix assembly performance the user should preallocate the matrix storage by
126367be906fSBarry Smith    setting the parameters `d_nz` (or `d_nnz`) and `o_nz` (or `o_nnz`).
1264d4002b98SHong Zhang 
1265d083f849SBarry Smith    Collective
1266d4002b98SHong Zhang 
1267d4002b98SHong Zhang    Input Parameters:
1268d4002b98SHong Zhang +  B - the matrix
1269d4002b98SHong Zhang .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1270d4002b98SHong Zhang            (same value is used for all local rows)
1271d4002b98SHong Zhang .  d_nnz - array containing the number of nonzeros in the various rows of the
1272d4002b98SHong Zhang            DIAGONAL portion of the local submatrix (possibly different for each row)
127367be906fSBarry Smith            or NULL (`PETSC_NULL_INTEGER` in Fortran), if `d_nz` is used to specify the nonzero structure.
1274d4002b98SHong Zhang            The size of this array is equal to the number of local rows, i.e 'm'.
1275d4002b98SHong Zhang            For matrices that will be factored, you must leave room for (and set)
1276d4002b98SHong Zhang            the diagonal entry even if it is zero.
1277d4002b98SHong Zhang .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1278d4002b98SHong Zhang            submatrix (same value is used for all local rows).
1279d4002b98SHong Zhang -  o_nnz - array containing the number of nonzeros in the various rows of the
1280d4002b98SHong Zhang            OFF-DIAGONAL portion of the local submatrix (possibly different for
128167be906fSBarry Smith            each row) or NULL (`PETSC_NULL_INTEGER` in Fortran), if `o_nz` is used to specify the nonzero
1282d4002b98SHong Zhang            structure. The size of this array is equal to the number
1283d4002b98SHong Zhang            of local rows, i.e 'm'.
1284d4002b98SHong Zhang 
1285d4002b98SHong Zhang    Example usage:
1286d4002b98SHong Zhang    Consider the following 8x8 matrix with 34 non-zero values, that is
1287d4002b98SHong Zhang    assembled across 3 processors. Lets assume that proc0 owns 3 rows,
1288d4002b98SHong Zhang    proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
128967be906fSBarry Smith    as follows
1290d4002b98SHong Zhang 
1291d4002b98SHong Zhang .vb
1292d4002b98SHong Zhang             1  2  0  |  0  3  0  |  0  4
1293d4002b98SHong Zhang     Proc0   0  5  6  |  7  0  0  |  8  0
1294d4002b98SHong Zhang             9  0 10  | 11  0  0  | 12  0
1295d4002b98SHong Zhang     -------------------------------------
1296d4002b98SHong Zhang            13  0 14  | 15 16 17  |  0  0
1297d4002b98SHong Zhang     Proc1   0 18  0  | 19 20 21  |  0  0
1298d4002b98SHong Zhang             0  0  0  | 22 23  0  | 24  0
1299d4002b98SHong Zhang     -------------------------------------
1300d4002b98SHong Zhang     Proc2  25 26 27  |  0  0 28  | 29  0
1301d4002b98SHong Zhang            30  0  0  | 31 32 33  |  0 34
1302d4002b98SHong Zhang .ve
1303d4002b98SHong Zhang 
130427430b45SBarry Smith    This can be represented as a collection of submatrices as
1305d4002b98SHong Zhang 
1306d4002b98SHong Zhang .vb
1307d4002b98SHong Zhang       A B C
1308d4002b98SHong Zhang       D E F
1309d4002b98SHong Zhang       G H I
1310d4002b98SHong Zhang .ve
1311d4002b98SHong Zhang 
1312d4002b98SHong Zhang    Where the submatrices A,B,C are owned by proc0, D,E,F are
1313d4002b98SHong Zhang    owned by proc1, G,H,I are owned by proc2.
1314d4002b98SHong Zhang 
1315d4002b98SHong Zhang    The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1316d4002b98SHong Zhang    The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1317d4002b98SHong Zhang    The 'M','N' parameters are 8,8, and have the same values on all procs.
1318d4002b98SHong Zhang 
1319d4002b98SHong Zhang    The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
1320d4002b98SHong Zhang    submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
1321d4002b98SHong Zhang    corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
1322d4002b98SHong Zhang    Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
132327430b45SBarry Smith    part as `MATSEQSELL` matrices. For example, proc1 will store [E] as a `MATSEQSELL`
1324d4002b98SHong Zhang    matrix, ans [DF] as another SeqSELL matrix.
1325d4002b98SHong Zhang 
132667be906fSBarry Smith    When `d_nz`, `o_nz` parameters are specified, `d_nz` storage elements are
1327d4002b98SHong Zhang    allocated for every row of the local diagonal submatrix, and o_nz
1328d4002b98SHong Zhang    storage locations are allocated for every row of the OFF-DIAGONAL submat.
132967be906fSBarry Smith    One way to choose `d_nz` and `o_nz` is to use the max nonzerors per local
1330d4002b98SHong Zhang    rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
133127430b45SBarry Smith    In this case, the values of d_nz,o_nz are
1332d4002b98SHong Zhang .vb
133327430b45SBarry Smith      proc0  dnz = 2, o_nz = 2
133427430b45SBarry Smith      proc1  dnz = 3, o_nz = 2
133527430b45SBarry Smith      proc2  dnz = 1, o_nz = 4
1336d4002b98SHong Zhang .ve
1337d4002b98SHong Zhang    We are allocating m*(d_nz+o_nz) storage locations for every proc. This
1338d4002b98SHong Zhang    translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
1339d4002b98SHong Zhang    for proc3. i.e we are using 12+15+10=37 storage locations to store
1340d4002b98SHong Zhang    34 values.
1341d4002b98SHong Zhang 
134267be906fSBarry Smith    When `d_nnz`, `o_nnz` parameters are specified, the storage is specified
1343a5b23f4aSJose E. Roman    for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
134427430b45SBarry Smith    In the above case the values for d_nnz,o_nnz are
1345d4002b98SHong Zhang .vb
134627430b45SBarry Smith      proc0 d_nnz = [2,2,2] and o_nnz = [2,2,2]
134727430b45SBarry Smith      proc1 d_nnz = [3,3,2] and o_nnz = [2,1,1]
134827430b45SBarry Smith      proc2 d_nnz = [1,1]   and o_nnz = [4,4]
1349d4002b98SHong Zhang .ve
1350d4002b98SHong Zhang    Here the space allocated is according to nz (or maximum values in the nnz
1351d4002b98SHong Zhang    if nnz is provided) for DIAGONAL and OFF-DIAGONAL submatrices, i.e (2+2+3+2)*3+(1+4)*2=37
1352d4002b98SHong Zhang 
1353d4002b98SHong Zhang    Level: intermediate
1354d4002b98SHong Zhang 
135567be906fSBarry Smith    Notes:
135667be906fSBarry Smith    If the *_nnz parameter is given then the *_nz parameter is ignored
135767be906fSBarry Smith 
135867be906fSBarry Smith    The stored row and column indices begin with zero.
135967be906fSBarry Smith 
136067be906fSBarry Smith    The parallel matrix is partitioned such that the first m0 rows belong to
136167be906fSBarry Smith    process 0, the next m1 rows belong to process 1, the next m2 rows belong
136267be906fSBarry Smith    to process 2 etc.. where m0,m1,m2... are the input parameter 'm'.
136367be906fSBarry Smith 
136467be906fSBarry Smith    The DIAGONAL portion of the local submatrix of a processor can be defined
136567be906fSBarry Smith    as the submatrix which is obtained by extraction the part corresponding to
136667be906fSBarry Smith    the rows r1-r2 and columns c1-c2 of the global matrix, where r1 is the
136767be906fSBarry Smith    first row that belongs to the processor, r2 is the last row belonging to
136867be906fSBarry Smith    the this processor, and c1-c2 is range of indices of the local part of a
136967be906fSBarry Smith    vector suitable for applying the matrix to.  This is an mxn matrix.  In the
137067be906fSBarry Smith    common case of a square matrix, the row and column ranges are the same and
137167be906fSBarry Smith    the DIAGONAL part is also square. The remaining portion of the local
137267be906fSBarry Smith    submatrix (mxN) constitute the OFF-DIAGONAL portion.
137367be906fSBarry Smith 
137467be906fSBarry Smith    If `o_nnz`, `d_nnz` are specified, then `o_nz`, and `d_nz` are ignored.
137567be906fSBarry Smith 
137667be906fSBarry Smith    You can call `MatGetInfo()` to get information on how effective the preallocation was;
137767be906fSBarry Smith    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
137867be906fSBarry Smith    You can also run with the option -info and look for messages with the string
137967be906fSBarry Smith    malloc in them to see if additional memory allocation was needed.
138067be906fSBarry Smith 
138167be906fSBarry Smith .seealso: `Mat`, `MatCreate()`, `MatCreateSeqSELL()`, `MatSetValues()`, `MatCreatesell()`,
138211a5261eSBarry Smith           `MATMPISELL`, `MatGetInfo()`, `PetscSplitOwnership()`, `MATSELL`
1383d4002b98SHong Zhang @*/
1384d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPISELLSetPreallocation(Mat B, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[])
1385d71ae5a4SJacob Faibussowitsch {
1386d4002b98SHong Zhang   PetscFunctionBegin;
1387d4002b98SHong Zhang   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
1388d4002b98SHong Zhang   PetscValidType(B, 1);
1389cac4c232SBarry Smith   PetscTryMethod(B, "MatMPISELLSetPreallocation_C", (Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]), (B, d_nz, d_nnz, o_nz, o_nnz));
13903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1391d4002b98SHong Zhang }
1392d4002b98SHong Zhang 
1393ed73aabaSBarry Smith /*MC
1394ed73aabaSBarry Smith    MATMPISELL - MATMPISELL = "mpisell" - A matrix type to be used for MPI sparse matrices,
1395ed73aabaSBarry Smith    based on the sliced Ellpack format
1396ed73aabaSBarry Smith 
139727430b45SBarry Smith    Options Database Key:
139820f4b53cSBarry Smith . -mat_type sell - sets the matrix type to `MATSELL` during a call to `MatSetFromOptions()`
1399ed73aabaSBarry Smith 
1400ed73aabaSBarry Smith    Level: beginner
1401ed73aabaSBarry Smith 
140267be906fSBarry Smith .seealso: `Mat`, `MatCreateSell()`, `MATSEQSELL`, `MATSELL`, `MATSEQAIJ`, `MATAIJ`, `MATMPIAIJ`
1403ed73aabaSBarry Smith M*/
1404ed73aabaSBarry Smith 
1405d4002b98SHong Zhang /*@C
140611a5261eSBarry Smith    MatCreateSELL - Creates a sparse parallel matrix in `MATSELL` format.
1407d4002b98SHong Zhang 
1408d083f849SBarry Smith    Collective
1409d4002b98SHong Zhang 
1410d4002b98SHong Zhang    Input Parameters:
1411d4002b98SHong Zhang +  comm - MPI communicator
141211a5261eSBarry Smith .  m - number of local rows (or `PETSC_DECIDE` to have calculated if M is given)
1413d4002b98SHong Zhang            This value should be the same as the local size used in creating the
1414d4002b98SHong Zhang            y vector for the matrix-vector product y = Ax.
1415d4002b98SHong Zhang .  n - This value should be the same as the local size used in creating the
141620f4b53cSBarry Smith        x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
141720f4b53cSBarry Smith        calculated if `N` is given) For square matrices n is almost always `m`.
141820f4b53cSBarry Smith .  M - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
141920f4b53cSBarry Smith .  N - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
1420d4002b98SHong Zhang .  d_rlenmax - max number of nonzeros per row in DIAGONAL portion of local submatrix
1421d4002b98SHong Zhang                (same value is used for all local rows)
1422d4002b98SHong Zhang .  d_rlen - array containing the number of nonzeros in the various rows of the
1423d4002b98SHong Zhang             DIAGONAL portion of the local submatrix (possibly different for each row)
142420f4b53cSBarry Smith             or `NULL`, if d_rlenmax is used to specify the nonzero structure.
142520f4b53cSBarry Smith             The size of this array is equal to the number of local rows, i.e `m`.
1426d4002b98SHong Zhang .  o_rlenmax - max number of nonzeros per row in the OFF-DIAGONAL portion of local
1427d4002b98SHong Zhang                submatrix (same value is used for all local rows).
1428d4002b98SHong Zhang -  o_rlen - array containing the number of nonzeros in the various rows of the
1429d4002b98SHong Zhang             OFF-DIAGONAL portion of the local submatrix (possibly different for
143020f4b53cSBarry Smith             each row) or `NULL`, if `o_rlenmax` is used to specify the nonzero
1431d4002b98SHong Zhang             structure. The size of this array is equal to the number
143220f4b53cSBarry Smith             of local rows, i.e `m`.
1433d4002b98SHong Zhang 
1434d4002b98SHong Zhang    Output Parameter:
1435d4002b98SHong Zhang .  A - the matrix
1436d4002b98SHong Zhang 
143727430b45SBarry Smith    Options Database Key:
143827430b45SBarry Smith -  -mat_sell_oneindex - Internally use indexing starting at 1
143927430b45SBarry Smith         rather than 0.  When calling `MatSetValues()`,
144027430b45SBarry Smith         the user still MUST index entries starting at 0!
144127430b45SBarry Smith 
144227430b45SBarry Smith    Example:
144327430b45SBarry Smith    Consider the following 8x8 matrix with 34 non-zero values, that is
144427430b45SBarry Smith    assembled across 3 processors. Lets assume that proc0 owns 3 rows,
144527430b45SBarry Smith    proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
144627430b45SBarry Smith    as follows
144727430b45SBarry Smith 
144827430b45SBarry Smith .vb
144927430b45SBarry Smith             1  2  0  |  0  3  0  |  0  4
145027430b45SBarry Smith     Proc0   0  5  6  |  7  0  0  |  8  0
145127430b45SBarry Smith             9  0 10  | 11  0  0  | 12  0
145227430b45SBarry Smith     -------------------------------------
145327430b45SBarry Smith            13  0 14  | 15 16 17  |  0  0
145427430b45SBarry Smith     Proc1   0 18  0  | 19 20 21  |  0  0
145527430b45SBarry Smith             0  0  0  | 22 23  0  | 24  0
145627430b45SBarry Smith     -------------------------------------
145727430b45SBarry Smith     Proc2  25 26 27  |  0  0 28  | 29  0
145827430b45SBarry Smith            30  0  0  | 31 32 33  |  0 34
145927430b45SBarry Smith .ve
146027430b45SBarry Smith 
146120f4b53cSBarry Smith    This can be represented as a collection of submatrices as
146227430b45SBarry Smith .vb
146327430b45SBarry Smith       A B C
146427430b45SBarry Smith       D E F
146527430b45SBarry Smith       G H I
146627430b45SBarry Smith .ve
146727430b45SBarry Smith 
146827430b45SBarry Smith    Where the submatrices A,B,C are owned by proc0, D,E,F are
146927430b45SBarry Smith    owned by proc1, G,H,I are owned by proc2.
147027430b45SBarry Smith 
147127430b45SBarry Smith    The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
147227430b45SBarry Smith    The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
147327430b45SBarry Smith    The 'M','N' parameters are 8,8, and have the same values on all procs.
147427430b45SBarry Smith 
147527430b45SBarry Smith    The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
147627430b45SBarry Smith    submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
147727430b45SBarry Smith    corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
147827430b45SBarry Smith    Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
147927430b45SBarry Smith    part as `MATSEQSELL` matrices. For example, proc1 will store [E] as a `MATSEQSELL`
148027430b45SBarry Smith    matrix, ans [DF] as another `MATSEQSELL` matrix.
148127430b45SBarry Smith 
148227430b45SBarry Smith    When d_rlenmax, o_rlenmax parameters are specified, d_rlenmax storage elements are
148327430b45SBarry Smith    allocated for every row of the local diagonal submatrix, and o_rlenmax
148427430b45SBarry Smith    storage locations are allocated for every row of the OFF-DIAGONAL submat.
148527430b45SBarry Smith    One way to choose d_rlenmax and o_rlenmax is to use the max nonzerors per local
148627430b45SBarry Smith    rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
148720f4b53cSBarry Smith    In this case, the values of d_rlenmax,o_rlenmax are
148827430b45SBarry Smith .vb
148920f4b53cSBarry Smith      proc0 - d_rlenmax = 2, o_rlenmax = 2
149020f4b53cSBarry Smith      proc1 - d_rlenmax = 3, o_rlenmax = 2
149120f4b53cSBarry Smith      proc2 - d_rlenmax = 1, o_rlenmax = 4
149227430b45SBarry Smith .ve
149327430b45SBarry Smith    We are allocating m*(d_rlenmax+o_rlenmax) storage locations for every proc. This
149427430b45SBarry Smith    translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
149527430b45SBarry Smith    for proc3. i.e we are using 12+15+10=37 storage locations to store
149627430b45SBarry Smith    34 values.
149727430b45SBarry Smith 
149820f4b53cSBarry Smith    When `d_rlen`, `o_rlen` parameters are specified, the storage is specified
149927430b45SBarry Smith    for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
150020f4b53cSBarry Smith    In the above case the values for `d_nnz`, `o_nnz` are
150127430b45SBarry Smith .vb
150220f4b53cSBarry Smith      proc0 - d_nnz = [2,2,2] and o_nnz = [2,2,2]
150320f4b53cSBarry Smith      proc1 - d_nnz = [3,3,2] and o_nnz = [2,1,1]
150420f4b53cSBarry Smith      proc2 - d_nnz = [1,1]   and o_nnz = [4,4]
150527430b45SBarry Smith .ve
150627430b45SBarry Smith    Here the space allocated is still 37 though there are 34 nonzeros because
150727430b45SBarry Smith    the allocation is always done according to rlenmax.
150827430b45SBarry Smith 
150927430b45SBarry Smith    Level: intermediate
151027430b45SBarry Smith 
151127430b45SBarry Smith    Notes:
151211a5261eSBarry Smith    It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`,
1513f6f02116SRichard Tran Mills    MatXXXXSetPreallocation() paradigm instead of this routine directly.
151411a5261eSBarry Smith    [MatXXXXSetPreallocation() is, for example, `MatSeqSELLSetPreallocation()`]
1515d4002b98SHong Zhang 
1516d4002b98SHong Zhang    If the *_rlen parameter is given then the *_rlenmax parameter is ignored
1517d4002b98SHong Zhang 
151820f4b53cSBarry Smith    `m`, `n`, `M`, `N` parameters specify the size of the matrix, and its partitioning across
151920f4b53cSBarry Smith    processors, while `d_rlenmax`, `d_rlen`, `o_rlenmax` , `o_rlen` parameters specify the approximate
1520d4002b98SHong Zhang    storage requirements for this matrix.
1521d4002b98SHong Zhang 
152211a5261eSBarry Smith    If `PETSC_DECIDE` or  `PETSC_DETERMINE` is used for a particular argument on one
1523d4002b98SHong Zhang    processor than it must be used on all processors that share the object for
1524d4002b98SHong Zhang    that argument.
1525d4002b98SHong Zhang 
1526d4002b98SHong Zhang    The user MUST specify either the local or global matrix dimensions
1527d4002b98SHong Zhang    (possibly both).
1528d4002b98SHong Zhang 
1529d4002b98SHong Zhang    The parallel matrix is partitioned across processors such that the
1530d4002b98SHong Zhang    first m0 rows belong to process 0, the next m1 rows belong to
1531d4002b98SHong Zhang    process 1, the next m2 rows belong to process 2 etc.. where
1532d4002b98SHong Zhang    m0,m1,m2,.. are the input parameter 'm'. i.e each processor stores
153320f4b53cSBarry Smith    values corresponding to [`m` x `N`] submatrix.
1534d4002b98SHong Zhang 
1535d4002b98SHong Zhang    The columns are logically partitioned with the n0 columns belonging
1536d4002b98SHong Zhang    to 0th partition, the next n1 columns belonging to the next
153720f4b53cSBarry Smith    partition etc.. where n0,n1,n2... are the input parameter `n`.
1538d4002b98SHong Zhang 
1539d4002b98SHong Zhang    The DIAGONAL portion of the local submatrix on any given processor
154020f4b53cSBarry Smith    is the submatrix corresponding to the rows and columns `m`, `n`
1541d4002b98SHong Zhang    corresponding to the given processor. i.e diagonal matrix on
1542d4002b98SHong Zhang    process 0 is [m0 x n0], diagonal matrix on process 1 is [m1 x n1]
1543d4002b98SHong Zhang    etc. The remaining portion of the local submatrix [m x (N-n)]
1544d4002b98SHong Zhang    constitute the OFF-DIAGONAL portion. The example below better
1545d4002b98SHong Zhang    illustrates this concept.
1546d4002b98SHong Zhang 
1547d4002b98SHong Zhang    For a square global matrix we define each processor's diagonal portion
1548d4002b98SHong Zhang    to be its local rows and the corresponding columns (a square submatrix);
1549d4002b98SHong Zhang    each processor's off-diagonal portion encompasses the remainder of the
1550d4002b98SHong Zhang    local matrix (a rectangular submatrix).
1551d4002b98SHong Zhang 
155220f4b53cSBarry Smith    If `o_rlen`, `d_rlen` are specified, then `o_rlenmax`, and `d_rlenmax` are ignored.
1553d4002b98SHong Zhang 
1554d4002b98SHong Zhang    When calling this routine with a single process communicator, a matrix of
155511a5261eSBarry Smith    type `MATSEQSELL` is returned.  If a matrix of type `MATMPISELL` is desired for this
155627430b45SBarry Smith    type of communicator, use the construction mechanism
1557d4002b98SHong Zhang .vb
155827430b45SBarry Smith    MatCreate(...,&A);
155927430b45SBarry Smith    MatSetType(A,MATMPISELL);
156027430b45SBarry Smith    MatSetSizes(A, m,n,M,N);
156127430b45SBarry Smith    MatMPISELLSetPreallocation(A,...);
1562d4002b98SHong Zhang .ve
1563d4002b98SHong Zhang 
156467be906fSBarry Smith .seealso: `Mat`, `MATSELL`, `MatCreate()`, `MatCreateSeqSELL()`, `MatSetValues()`, `MatMPISELLSetPreallocation()`, `MatMPISELLSetPreallocationSELL()`,
1565db781477SPatrick Sanan           `MATMPISELL`, `MatCreateMPISELLWithArrays()`
1566d4002b98SHong Zhang @*/
1567d71ae5a4SJacob Faibussowitsch 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)
1568d71ae5a4SJacob Faibussowitsch {
1569d4002b98SHong Zhang   PetscMPIInt size;
1570d4002b98SHong Zhang 
1571d4002b98SHong Zhang   PetscFunctionBegin;
15729566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
15739566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, n, M, N));
15749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1575d4002b98SHong Zhang   if (size > 1) {
15769566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A, MATMPISELL));
15779566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(*A, d_rlenmax, d_rlen, o_rlenmax, o_rlen));
1578d4002b98SHong Zhang   } else {
15799566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A, MATSEQSELL));
15809566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(*A, d_rlenmax, d_rlen));
1581d4002b98SHong Zhang   }
15823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1583d4002b98SHong Zhang }
1584d4002b98SHong Zhang 
1585d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPISELLGetSeqSELL(Mat A, Mat *Ad, Mat *Ao, const PetscInt *colmap[])
1586d71ae5a4SJacob Faibussowitsch {
1587d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1588d4002b98SHong Zhang   PetscBool    flg;
1589d4002b98SHong Zhang 
1590d4002b98SHong Zhang   PetscFunctionBegin;
15919566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISELL, &flg));
159228b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "This function requires a MATMPISELL matrix as input");
1593d4002b98SHong Zhang   if (Ad) *Ad = a->A;
1594d4002b98SHong Zhang   if (Ao) *Ao = a->B;
1595d4002b98SHong Zhang   if (colmap) *colmap = a->garray;
15963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1597d4002b98SHong Zhang }
1598d4002b98SHong Zhang 
1599d4002b98SHong Zhang /*@C
160011a5261eSBarry Smith      MatMPISELLGetLocalMatCondensed - Creates a `MATSEQSELL` matrix from an `MATMPISELL` matrix by taking all its local rows and NON-ZERO columns
1601d4002b98SHong Zhang 
1602d4002b98SHong Zhang     Not Collective
1603d4002b98SHong Zhang 
1604d4002b98SHong Zhang    Input Parameters:
1605d4002b98SHong Zhang +    A - the matrix
160611a5261eSBarry Smith .    scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
160767be906fSBarry Smith .    row - index sets of rows to extract (or `NULL`)
160867be906fSBarry Smith -    col - index sets of columns to extract (or `NULL`)
1609d4002b98SHong Zhang 
1610d4002b98SHong Zhang    Output Parameter:
1611d4002b98SHong Zhang .    A_loc - the local sequential matrix generated
1612d4002b98SHong Zhang 
1613d4002b98SHong Zhang     Level: developer
1614d4002b98SHong Zhang 
161567be906fSBarry Smith .seealso: `Mat`, `MATSEQSELL`, `MATMPISELL`, `MatGetOwnershipRange()`, `MatMPISELLGetLocalMat()`
1616d4002b98SHong Zhang @*/
1617d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPISELLGetLocalMatCondensed(Mat A, MatReuse scall, IS *row, IS *col, Mat *A_loc)
1618d71ae5a4SJacob Faibussowitsch {
1619d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1620d4002b98SHong Zhang   PetscInt     i, start, end, ncols, nzA, nzB, *cmap, imark, *idx;
1621d4002b98SHong Zhang   IS           isrowa, iscola;
1622d4002b98SHong Zhang   Mat         *aloc;
1623d4002b98SHong Zhang   PetscBool    match;
1624d4002b98SHong Zhang 
1625d4002b98SHong Zhang   PetscFunctionBegin;
16269566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISELL, &match));
162728b400f6SJacob Faibussowitsch   PetscCheck(match, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Requires MATMPISELL matrix as input");
16289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_Getlocalmatcondensed, A, 0, 0, 0));
1629d4002b98SHong Zhang   if (!row) {
16309371c9d4SSatish Balay     start = A->rmap->rstart;
16319371c9d4SSatish Balay     end   = A->rmap->rend;
16329566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, end - start, start, 1, &isrowa));
1633d4002b98SHong Zhang   } else {
1634d4002b98SHong Zhang     isrowa = *row;
1635d4002b98SHong Zhang   }
1636d4002b98SHong Zhang   if (!col) {
1637d4002b98SHong Zhang     start = A->cmap->rstart;
1638d4002b98SHong Zhang     cmap  = a->garray;
1639d4002b98SHong Zhang     nzA   = a->A->cmap->n;
1640d4002b98SHong Zhang     nzB   = a->B->cmap->n;
16419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nzA + nzB, &idx));
1642d4002b98SHong Zhang     ncols = 0;
1643d4002b98SHong Zhang     for (i = 0; i < nzB; i++) {
1644d4002b98SHong Zhang       if (cmap[i] < start) idx[ncols++] = cmap[i];
1645d4002b98SHong Zhang       else break;
1646d4002b98SHong Zhang     }
1647d4002b98SHong Zhang     imark = i;
1648d4002b98SHong Zhang     for (i = 0; i < nzA; i++) idx[ncols++] = start + i;
1649d4002b98SHong Zhang     for (i = imark; i < nzB; i++) idx[ncols++] = cmap[i];
16509566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncols, idx, PETSC_OWN_POINTER, &iscola));
1651d4002b98SHong Zhang   } else {
1652d4002b98SHong Zhang     iscola = *col;
1653d4002b98SHong Zhang   }
1654d4002b98SHong Zhang   if (scall != MAT_INITIAL_MATRIX) {
16559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &aloc));
1656d4002b98SHong Zhang     aloc[0] = *A_loc;
1657d4002b98SHong Zhang   }
16589566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, 1, &isrowa, &iscola, scall, &aloc));
1659d4002b98SHong Zhang   *A_loc = aloc[0];
16609566063dSJacob Faibussowitsch   PetscCall(PetscFree(aloc));
166148a46eb9SPierre Jolivet   if (!row) PetscCall(ISDestroy(&isrowa));
166248a46eb9SPierre Jolivet   if (!col) PetscCall(ISDestroy(&iscola));
16639566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_Getlocalmatcondensed, A, 0, 0, 0));
16643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1665d4002b98SHong Zhang }
1666d4002b98SHong Zhang 
1667d4002b98SHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>
1668d4002b98SHong Zhang 
1669d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_MPISELL_MPIAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1670d71ae5a4SJacob Faibussowitsch {
1671d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1672d4002b98SHong Zhang   Mat          B;
1673d4002b98SHong Zhang   Mat_MPIAIJ  *b;
1674d4002b98SHong Zhang 
1675d4002b98SHong Zhang   PetscFunctionBegin;
167628b400f6SJacob Faibussowitsch   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Matrix must be assembled");
1677d4002b98SHong Zhang 
167894a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
167994a8b381SRichard Tran Mills     B = *newmat;
168094a8b381SRichard Tran Mills   } else {
16819566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
16829566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATMPIAIJ));
16839566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
16849566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
16859566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL));
16869566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL));
168794a8b381SRichard Tran Mills   }
1688d4002b98SHong Zhang   b = (Mat_MPIAIJ *)B->data;
168994a8b381SRichard Tran Mills 
169094a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
16919566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_REUSE_MATRIX, &b->A));
16929566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_REUSE_MATRIX, &b->B));
169394a8b381SRichard Tran Mills   } else {
16949566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->A));
16959566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->B));
16969566063dSJacob Faibussowitsch     PetscCall(MatDisAssemble_MPISELL(A));
16979566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->A));
16989566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->B));
16999566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
17009566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
17019566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
17029566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
170394a8b381SRichard Tran Mills   }
1704d4002b98SHong Zhang 
1705d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
17069566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
1707d4002b98SHong Zhang   } else {
1708d4002b98SHong Zhang     *newmat = B;
1709d4002b98SHong Zhang   }
17103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1711d4002b98SHong Zhang }
1712d4002b98SHong Zhang 
1713d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_MPIAIJ_MPISELL(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1714d71ae5a4SJacob Faibussowitsch {
1715d4002b98SHong Zhang   Mat_MPIAIJ  *a = (Mat_MPIAIJ *)A->data;
1716d4002b98SHong Zhang   Mat          B;
1717d4002b98SHong Zhang   Mat_MPISELL *b;
1718d4002b98SHong Zhang 
1719d4002b98SHong Zhang   PetscFunctionBegin;
172028b400f6SJacob Faibussowitsch   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Matrix must be assembled");
1721d4002b98SHong Zhang 
172294a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
172394a8b381SRichard Tran Mills     B = *newmat;
172494a8b381SRichard Tran Mills   } else {
17252d1451d4SHong Zhang     Mat_SeqAIJ *Aa = (Mat_SeqAIJ *)a->A->data, *Ba = (Mat_SeqAIJ *)a->B->data;
17262d1451d4SHong Zhang     PetscInt    i, d_nz = 0, o_nz = 0, m = A->rmap->N, n = A->cmap->N, lm = A->rmap->n, ln = A->cmap->n;
17272d1451d4SHong Zhang     PetscInt   *d_nnz, *o_nnz;
17282d1451d4SHong Zhang     PetscCall(PetscMalloc2(lm, &d_nnz, lm, &o_nnz));
17292d1451d4SHong Zhang     for (i = 0; i < lm; i++) {
17302d1451d4SHong Zhang       d_nnz[i] = Aa->i[i + 1] - Aa->i[i];
17312d1451d4SHong Zhang       o_nnz[i] = Ba->i[i + 1] - Ba->i[i];
17322d1451d4SHong Zhang       if (d_nnz[i] > d_nz) d_nz = d_nnz[i];
17332d1451d4SHong Zhang       if (o_nnz[i] > o_nz) o_nz = o_nnz[i];
17342d1451d4SHong Zhang     }
17359566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
17369566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATMPISELL));
17372d1451d4SHong Zhang     PetscCall(MatSetSizes(B, lm, ln, m, n));
17389566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
17392d1451d4SHong Zhang     PetscCall(MatSeqSELLSetPreallocation(B, d_nz, d_nnz));
17402d1451d4SHong Zhang     PetscCall(MatMPISELLSetPreallocation(B, d_nz, d_nnz, o_nz, o_nnz));
17412d1451d4SHong Zhang     PetscCall(PetscFree2(d_nnz, o_nnz));
174294a8b381SRichard Tran Mills   }
1743d4002b98SHong Zhang   b = (Mat_MPISELL *)B->data;
174494a8b381SRichard Tran Mills 
174594a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
17469566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_REUSE_MATRIX, &b->A));
17479566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_REUSE_MATRIX, &b->B));
174894a8b381SRichard Tran Mills   } else {
17499566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->A));
17509566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->B));
17519566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_INITIAL_MATRIX, &b->A));
17529566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_INITIAL_MATRIX, &b->B));
17539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
17549566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
17552d1451d4SHong Zhang     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
17562d1451d4SHong Zhang     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
175794a8b381SRichard Tran Mills   }
1758d4002b98SHong Zhang 
1759d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
17609566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
1761d4002b98SHong Zhang   } else {
1762d4002b98SHong Zhang     *newmat = B;
1763d4002b98SHong Zhang   }
17643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1765d4002b98SHong Zhang }
1766d4002b98SHong Zhang 
1767d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSOR_MPISELL(Mat matin, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
1768d71ae5a4SJacob Faibussowitsch {
1769d4002b98SHong Zhang   Mat_MPISELL *mat = (Mat_MPISELL *)matin->data;
1770f4259b30SLisandro Dalcin   Vec          bb1 = NULL;
1771d4002b98SHong Zhang 
1772d4002b98SHong Zhang   PetscFunctionBegin;
1773d4002b98SHong Zhang   if (flag == SOR_APPLY_UPPER) {
17749566063dSJacob Faibussowitsch     PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
17753ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1776d4002b98SHong Zhang   }
1777d4002b98SHong Zhang 
177848a46eb9SPierre Jolivet   if (its > 1 || ~flag & SOR_ZERO_INITIAL_GUESS || flag & SOR_EISENSTAT) PetscCall(VecDuplicate(bb, &bb1));
1779d4002b98SHong Zhang 
1780d4002b98SHong Zhang   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
1781d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17829566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1783d4002b98SHong Zhang       its--;
1784d4002b98SHong Zhang     }
1785d4002b98SHong Zhang 
1786d4002b98SHong Zhang     while (its--) {
17879566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
17889566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1789d4002b98SHong Zhang 
1790d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
17919566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec, -1.0));
17929566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));
1793d4002b98SHong Zhang 
1794d4002b98SHong Zhang       /* local sweep */
17959566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_SYMMETRIC_SWEEP, fshift, lits, 1, xx));
1796d4002b98SHong Zhang     }
1797d4002b98SHong Zhang   } else if (flag & SOR_LOCAL_FORWARD_SWEEP) {
1798d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17999566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1800d4002b98SHong Zhang       its--;
1801d4002b98SHong Zhang     }
1802d4002b98SHong Zhang     while (its--) {
18039566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
18049566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1805d4002b98SHong Zhang 
1806d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
18079566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec, -1.0));
18089566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));
1809d4002b98SHong Zhang 
1810d4002b98SHong Zhang       /* local sweep */
18119566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_FORWARD_SWEEP, fshift, lits, 1, xx));
1812d4002b98SHong Zhang     }
1813d4002b98SHong Zhang   } else if (flag & SOR_LOCAL_BACKWARD_SWEEP) {
1814d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
18159566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1816d4002b98SHong Zhang       its--;
1817d4002b98SHong Zhang     }
1818d4002b98SHong Zhang     while (its--) {
18199566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
18209566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1821d4002b98SHong Zhang 
1822d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
18239566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec, -1.0));
18249566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));
1825d4002b98SHong Zhang 
1826d4002b98SHong Zhang       /* local sweep */
18279566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_BACKWARD_SWEEP, fshift, lits, 1, xx));
1828d4002b98SHong Zhang     }
1829d4002b98SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)matin), PETSC_ERR_SUP, "Parallel SOR not supported");
1830d4002b98SHong Zhang 
18319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bb1));
1832d4002b98SHong Zhang 
1833d4002b98SHong Zhang   matin->factorerrortype = mat->A->factorerrortype;
18343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1835d4002b98SHong Zhang }
1836d4002b98SHong Zhang 
1837*b5917f1bSHong Zhang #if defined(PETSC_HAVE_CUDA)
1838*b5917f1bSHong Zhang PETSC_INTERN PetscErrorCode MatConvert_MPISELL_MPISELLCUDA(Mat, MatType, MatReuse, Mat *);
1839*b5917f1bSHong Zhang #endif
1840*b5917f1bSHong Zhang 
1841d4002b98SHong Zhang /*MC
1842d4002b98SHong Zhang    MATMPISELL - MATMPISELL = "MPISELL" - A matrix type to be used for parallel sparse matrices.
1843d4002b98SHong Zhang 
1844d4002b98SHong Zhang    Options Database Keys:
184511a5261eSBarry Smith . -mat_type mpisell - sets the matrix type to `MATMPISELL` during a call to `MatSetFromOptions()`
1846d4002b98SHong Zhang 
1847d4002b98SHong Zhang   Level: beginner
1848d4002b98SHong Zhang 
184967be906fSBarry Smith .seealso: `Mat`, `MATSELL`, `MATSEQSELL` `MatCreateSELL()`
1850d4002b98SHong Zhang M*/
1851d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat B)
1852d71ae5a4SJacob Faibussowitsch {
1853d4002b98SHong Zhang   Mat_MPISELL *b;
1854d4002b98SHong Zhang   PetscMPIInt  size;
1855d4002b98SHong Zhang 
1856d4002b98SHong Zhang   PetscFunctionBegin;
18579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
18584dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
1859d4002b98SHong Zhang   B->data = (void *)b;
18609566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(B->ops, &MatOps_Values, sizeof(struct _MatOps)));
1861d4002b98SHong Zhang   B->assembled  = PETSC_FALSE;
1862d4002b98SHong Zhang   B->insertmode = NOT_SET_VALUES;
1863d4002b98SHong Zhang   b->size       = size;
18649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)B), &b->rank));
1865d4002b98SHong Zhang   /* build cache for off array entries formed */
18669566063dSJacob Faibussowitsch   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)B), 1, &B->stash));
1867d4002b98SHong Zhang 
1868d4002b98SHong Zhang   b->donotstash  = PETSC_FALSE;
1869f4259b30SLisandro Dalcin   b->colmap      = NULL;
1870f4259b30SLisandro Dalcin   b->garray      = NULL;
1871d4002b98SHong Zhang   b->roworiented = PETSC_TRUE;
1872d4002b98SHong Zhang 
1873d4002b98SHong Zhang   /* stuff used for matrix vector multiply */
1874d4002b98SHong Zhang   b->lvec  = NULL;
1875d4002b98SHong Zhang   b->Mvctx = NULL;
1876d4002b98SHong Zhang 
1877d4002b98SHong Zhang   /* stuff for MatGetRow() */
1878f4259b30SLisandro Dalcin   b->rowindices   = NULL;
1879f4259b30SLisandro Dalcin   b->rowvalues    = NULL;
1880d4002b98SHong Zhang   b->getrowactive = PETSC_FALSE;
1881d4002b98SHong Zhang 
18829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_MPISELL));
18839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_MPISELL));
18849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsTranspose_C", MatIsTranspose_MPISELL));
18859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPISELLSetPreallocation_C", MatMPISELLSetPreallocation_MPISELL));
18869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisell_mpiaij_C", MatConvert_MPISELL_MPIAIJ));
1887*b5917f1bSHong Zhang #if defined(PETSC_HAVE_CUDA)
1888*b5917f1bSHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisell_mpisellcuda_C", MatConvert_MPISELL_MPISELLCUDA));
1889*b5917f1bSHong Zhang #endif
18909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatDiagonalScaleLocal_C", MatDiagonalScaleLocal_MPISELL));
18919566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPISELL));
18923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1893d4002b98SHong Zhang }
1894