xref: /petsc/src/mat/impls/sell/mpi/mpisell.c (revision 8bb0f5c6fef004073486bd40bcac1547cd22e816)
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 
25ba38deedSJacob Faibussowitsch static 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++; \
102d4002b98SHong Zhang   a_noinsert:; \
103d4002b98SHong Zhang     a->rlen[row] = nrow1; \
104d4002b98SHong Zhang   }
105d4002b98SHong Zhang 
106d4002b98SHong Zhang #define MatSetValues_SeqSELL_B_Private(row, col, value, addv, orow, ocol) \
107d4002b98SHong Zhang   { \
108d4002b98SHong Zhang     if (col <= lastcol2) low2 = 0; \
109d4002b98SHong Zhang     else high2 = nrow2; \
110d4002b98SHong Zhang     lastcol2 = col; \
111d4002b98SHong Zhang     while (high2 - low2 > 5) { \
112d4002b98SHong Zhang       t = (low2 + high2) / 2; \
1134f8ff0b3SHong Zhang       if (cp2[sliceheight * t] > col) high2 = t; \
114d4002b98SHong Zhang       else low2 = t; \
115d4002b98SHong Zhang     } \
116d4002b98SHong Zhang     for (_i = low2; _i < high2; _i++) { \
1174f8ff0b3SHong Zhang       if (cp2[sliceheight * _i] > col) break; \
1184f8ff0b3SHong Zhang       if (cp2[sliceheight * _i] == col) { \
1194f8ff0b3SHong Zhang         if (addv == ADD_VALUES) vp2[sliceheight * _i] += value; \
1204f8ff0b3SHong Zhang         else vp2[sliceheight * _i] = value; \
1212d1451d4SHong Zhang         inserted = PETSC_TRUE; \
122d4002b98SHong Zhang         goto b_noinsert; \
123d4002b98SHong Zhang       } \
124d4002b98SHong Zhang     } \
1259371c9d4SSatish Balay     if (value == 0.0 && ignorezeroentries) { \
1269371c9d4SSatish Balay       low2  = 0; \
1279371c9d4SSatish Balay       high2 = nrow2; \
1289371c9d4SSatish Balay       goto b_noinsert; \
1299371c9d4SSatish Balay     } \
1309371c9d4SSatish Balay     if (nonew == 1) { \
1319371c9d4SSatish Balay       low2  = 0; \
1329371c9d4SSatish Balay       high2 = nrow2; \
1339371c9d4SSatish Balay       goto b_noinsert; \
1349371c9d4SSatish Balay     } \
13508401ef6SPierre 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); \
13607e43b41SHong Zhang     MatSeqXSELLReallocateSELL(B, bm, 1, nrow2, b->sliidx, b->sliceheight, row / sliceheight, row, col, b->colidx, b->val, cp2, vp2, nonew, MatScalar); \
137d4002b98SHong Zhang     /* shift up all the later entries in this row */ \
138d4002b98SHong Zhang     for (ii = nrow2 - 1; ii >= _i; ii--) { \
1394f8ff0b3SHong Zhang       cp2[sliceheight * (ii + 1)] = cp2[sliceheight * ii]; \
1404f8ff0b3SHong Zhang       vp2[sliceheight * (ii + 1)] = vp2[sliceheight * ii]; \
141d4002b98SHong Zhang     } \
1424f8ff0b3SHong Zhang     cp2[sliceheight * _i] = col; \
1434f8ff0b3SHong Zhang     vp2[sliceheight * _i] = value; \
1449371c9d4SSatish Balay     b->nz++; \
1459371c9d4SSatish Balay     nrow2++; \
146d4002b98SHong Zhang   b_noinsert:; \
147d4002b98SHong Zhang     b->rlen[row] = nrow2; \
148d4002b98SHong Zhang   }
149d4002b98SHong Zhang 
15043b34f9dSJacob Faibussowitsch static PetscErrorCode MatSetValues_MPISELL(Mat mat, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode addv)
151d71ae5a4SJacob Faibussowitsch {
152d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
153d4002b98SHong Zhang   PetscScalar  value;
154d4002b98SHong Zhang   PetscInt     i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend, shift1, shift2;
155d4002b98SHong Zhang   PetscInt     cstart = mat->cmap->rstart, cend = mat->cmap->rend, row, col;
156d4002b98SHong Zhang   PetscBool    roworiented = sell->roworiented;
157d4002b98SHong Zhang 
158d4002b98SHong Zhang   /* Some Variables required in the macro */
159d4002b98SHong Zhang   Mat          A                 = sell->A;
160d4002b98SHong Zhang   Mat_SeqSELL *a                 = (Mat_SeqSELL *)A->data;
161d4002b98SHong Zhang   PetscBool    ignorezeroentries = a->ignorezeroentries, found;
162d4002b98SHong Zhang   Mat          B                 = sell->B;
163d4002b98SHong Zhang   Mat_SeqSELL *b                 = (Mat_SeqSELL *)B->data;
16407e43b41SHong Zhang   PetscInt    *cp1, *cp2, ii, _i, nrow1, nrow2, low1, high1, low2, high2, t, lastcol1, lastcol2, sliceheight = a->sliceheight;
165d4002b98SHong Zhang   MatScalar   *vp1, *vp2;
166d4002b98SHong Zhang 
167d4002b98SHong Zhang   PetscFunctionBegin;
168d4002b98SHong Zhang   for (i = 0; i < m; i++) {
169d4002b98SHong Zhang     if (im[i] < 0) continue;
1706bdcaf15SBarry 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);
171d4002b98SHong Zhang     if (im[i] >= rstart && im[i] < rend) {
172d4002b98SHong Zhang       row      = im[i] - rstart;
173d4002b98SHong Zhang       lastcol1 = -1;
17407e43b41SHong Zhang       shift1   = a->sliidx[row / sliceheight] + (row % sliceheight); /* starting index of the row */
1758e3a54c0SPierre Jolivet       cp1      = PetscSafePointerPlusOffset(a->colidx, shift1);
1768e3a54c0SPierre Jolivet       vp1      = PetscSafePointerPlusOffset(a->val, shift1);
177d4002b98SHong Zhang       nrow1    = a->rlen[row];
178d4002b98SHong Zhang       low1     = 0;
179d4002b98SHong Zhang       high1    = nrow1;
180d4002b98SHong Zhang       lastcol2 = -1;
18107e43b41SHong Zhang       shift2   = b->sliidx[row / sliceheight] + (row % sliceheight); /* starting index of the row */
1828e3a54c0SPierre Jolivet       cp2      = PetscSafePointerPlusOffset(b->colidx, shift2);
1838e3a54c0SPierre Jolivet       vp2      = PetscSafePointerPlusOffset(b->val, shift2);
184d4002b98SHong Zhang       nrow2    = b->rlen[row];
185d4002b98SHong Zhang       low2     = 0;
186d4002b98SHong Zhang       high2    = nrow2;
187d4002b98SHong Zhang 
188d4002b98SHong Zhang       for (j = 0; j < n; j++) {
189d4002b98SHong Zhang         if (roworiented) value = v[i * n + j];
190d4002b98SHong Zhang         else value = v[i + j * m];
191d4002b98SHong Zhang         if (ignorezeroentries && value == 0.0 && (addv == ADD_VALUES)) continue;
192d4002b98SHong Zhang         if (in[j] >= cstart && in[j] < cend) {
193d4002b98SHong Zhang           col = in[j] - cstart;
194d4002b98SHong Zhang           MatSetValue_SeqSELL_Private(A, row, col, value, addv, im[i], in[j], cp1, vp1, lastcol1, low1, high1); /* set one value */
1952d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
1962d1451d4SHong Zhang           if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED && found) A->offloadmask = PETSC_OFFLOAD_CPU;
1972d1451d4SHong Zhang #endif
198f7d195e4SLawrence Mitchell         } else if (in[j] < 0) {
199f7d195e4SLawrence Mitchell           continue;
200f7d195e4SLawrence Mitchell         } else {
201f7d195e4SLawrence 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);
202d4002b98SHong Zhang           if (mat->was_assembled) {
20348a46eb9SPierre Jolivet             if (!sell->colmap) PetscCall(MatCreateColmap_MPISELL_Private(mat));
204d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
205eec179cfSJacob Faibussowitsch             PetscCall(PetscHMapIGetWithDefault(sell->colmap, in[j] + 1, 0, &col));
206d4002b98SHong Zhang             col--;
207d4002b98SHong Zhang #else
208d4002b98SHong Zhang             col = sell->colmap[in[j]] - 1;
209d4002b98SHong Zhang #endif
210f4f49eeaSPierre Jolivet             if (col < 0 && !((Mat_SeqSELL *)sell->B->data)->nonew) {
2119566063dSJacob Faibussowitsch               PetscCall(MatDisAssemble_MPISELL(mat));
212d4002b98SHong Zhang               col = in[j];
213d4002b98SHong Zhang               /* Reinitialize the variables required by MatSetValues_SeqSELL_B_Private() */
214d4002b98SHong Zhang               B      = sell->B;
215d4002b98SHong Zhang               b      = (Mat_SeqSELL *)B->data;
21607e43b41SHong Zhang               shift2 = b->sliidx[row / sliceheight] + (row % sliceheight); /* starting index of the row */
217d4002b98SHong Zhang               cp2    = b->colidx + shift2;
218d4002b98SHong Zhang               vp2    = b->val + shift2;
219d4002b98SHong Zhang               nrow2  = b->rlen[row];
220d4002b98SHong Zhang               low2   = 0;
221d4002b98SHong Zhang               high2  = nrow2;
2222d1451d4SHong Zhang               found  = PETSC_FALSE;
223f7d195e4SLawrence Mitchell             } else {
224f7d195e4SLawrence 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]);
225f7d195e4SLawrence Mitchell             }
226d4002b98SHong Zhang           } else col = in[j];
227d4002b98SHong Zhang           MatSetValue_SeqSELL_Private(B, row, col, value, addv, im[i], in[j], cp2, vp2, lastcol2, low2, high2); /* set one value */
2282d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
2292d1451d4SHong Zhang           if (B->offloadmask != PETSC_OFFLOAD_UNALLOCATED && found) B->offloadmask = PETSC_OFFLOAD_CPU;
2302d1451d4SHong Zhang #endif
231d4002b98SHong Zhang         }
232d4002b98SHong Zhang       }
233d4002b98SHong Zhang     } else {
23428b400f6SJacob 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]);
235d4002b98SHong Zhang       if (!sell->donotstash) {
236d4002b98SHong Zhang         mat->assembled = PETSC_FALSE;
237d4002b98SHong Zhang         if (roworiented) {
2389566063dSJacob Faibussowitsch           PetscCall(MatStashValuesRow_Private(&mat->stash, im[i], n, in, v + i * n, (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
239d4002b98SHong Zhang         } else {
2409566063dSJacob Faibussowitsch           PetscCall(MatStashValuesCol_Private(&mat->stash, im[i], n, in, v + i, m, (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
241d4002b98SHong Zhang         }
242d4002b98SHong Zhang       }
243d4002b98SHong Zhang     }
244d4002b98SHong Zhang   }
2453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
246d4002b98SHong Zhang }
247d4002b98SHong Zhang 
248ba38deedSJacob Faibussowitsch static PetscErrorCode MatGetValues_MPISELL(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[])
249d71ae5a4SJacob Faibussowitsch {
250d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
251d4002b98SHong Zhang   PetscInt     i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend;
252d4002b98SHong Zhang   PetscInt     cstart = mat->cmap->rstart, cend = mat->cmap->rend, row, col;
253d4002b98SHong Zhang 
254d4002b98SHong Zhang   PetscFunctionBegin;
255d4002b98SHong Zhang   for (i = 0; i < m; i++) {
25654c59aa7SJacob Faibussowitsch     if (idxm[i] < 0) continue; /* negative row */
25754c59aa7SJacob 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);
258d4002b98SHong Zhang     if (idxm[i] >= rstart && idxm[i] < rend) {
259d4002b98SHong Zhang       row = idxm[i] - rstart;
260d4002b98SHong Zhang       for (j = 0; j < n; j++) {
26154c59aa7SJacob Faibussowitsch         if (idxn[j] < 0) continue; /* negative column */
26254c59aa7SJacob 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);
263d4002b98SHong Zhang         if (idxn[j] >= cstart && idxn[j] < cend) {
264d4002b98SHong Zhang           col = idxn[j] - cstart;
2659566063dSJacob Faibussowitsch           PetscCall(MatGetValues(sell->A, 1, &row, 1, &col, v + i * n + j));
266d4002b98SHong Zhang         } else {
26748a46eb9SPierre Jolivet           if (!sell->colmap) PetscCall(MatCreateColmap_MPISELL_Private(mat));
268d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
269eec179cfSJacob Faibussowitsch           PetscCall(PetscHMapIGetWithDefault(sell->colmap, idxn[j] + 1, 0, &col));
270d4002b98SHong Zhang           col--;
271d4002b98SHong Zhang #else
272d4002b98SHong Zhang           col = sell->colmap[idxn[j]] - 1;
273d4002b98SHong Zhang #endif
274d4002b98SHong Zhang           if ((col < 0) || (sell->garray[col] != idxn[j])) *(v + i * n + j) = 0.0;
27548a46eb9SPierre Jolivet           else PetscCall(MatGetValues(sell->B, 1, &row, 1, &col, v + i * n + j));
276d4002b98SHong Zhang         }
277d4002b98SHong Zhang       }
278d4002b98SHong Zhang     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only local values currently supported");
279d4002b98SHong Zhang   }
2803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
281d4002b98SHong Zhang }
282d4002b98SHong Zhang 
283ba38deedSJacob Faibussowitsch static PetscErrorCode MatAssemblyBegin_MPISELL(Mat mat, MatAssemblyType mode)
284d71ae5a4SJacob Faibussowitsch {
285d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
286d4002b98SHong Zhang   PetscInt     nstash, reallocs;
287d4002b98SHong Zhang 
288d4002b98SHong Zhang   PetscFunctionBegin;
2893ba16761SJacob Faibussowitsch   if (sell->donotstash || mat->nooffprocentries) PetscFunctionReturn(PETSC_SUCCESS);
290d4002b98SHong Zhang 
2919566063dSJacob Faibussowitsch   PetscCall(MatStashScatterBegin_Private(mat, &mat->stash, mat->rmap->range));
2929566063dSJacob Faibussowitsch   PetscCall(MatStashGetInfo_Private(&mat->stash, &nstash, &reallocs));
2939566063dSJacob Faibussowitsch   PetscCall(PetscInfo(sell->A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
2943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
295d4002b98SHong Zhang }
296d4002b98SHong Zhang 
297d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_MPISELL(Mat mat, MatAssemblyType mode)
298d71ae5a4SJacob Faibussowitsch {
299d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
300d4002b98SHong Zhang   PetscMPIInt  n;
301d4002b98SHong Zhang   PetscInt     i, flg;
302d4002b98SHong Zhang   PetscInt    *row, *col;
303d4002b98SHong Zhang   PetscScalar *val;
304d4002b98SHong Zhang   PetscBool    other_disassembled;
305d4002b98SHong Zhang   /* do not use 'b = (Mat_SeqSELL*)sell->B->data' as B can be reset in disassembly */
306d4002b98SHong Zhang   PetscFunctionBegin;
307d4002b98SHong Zhang   if (!sell->donotstash && !mat->nooffprocentries) {
308d4002b98SHong Zhang     while (1) {
3099566063dSJacob Faibussowitsch       PetscCall(MatStashScatterGetMesg_Private(&mat->stash, &n, &row, &col, &val, &flg));
310d4002b98SHong Zhang       if (!flg) break;
311d4002b98SHong Zhang 
312d4002b98SHong Zhang       for (i = 0; i < n; i++) { /* assemble one by one */
3139566063dSJacob Faibussowitsch         PetscCall(MatSetValues_MPISELL(mat, 1, row + i, 1, col + i, val + i, mat->insertmode));
314d4002b98SHong Zhang       }
315d4002b98SHong Zhang     }
3169566063dSJacob Faibussowitsch     PetscCall(MatStashScatterEnd_Private(&mat->stash));
317d4002b98SHong Zhang   }
3182d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
3192d1451d4SHong Zhang   if (mat->offloadmask == PETSC_OFFLOAD_CPU) sell->A->offloadmask = PETSC_OFFLOAD_CPU;
3202d1451d4SHong Zhang #endif
3219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sell->A, mode));
3229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sell->A, mode));
323d4002b98SHong Zhang 
324d4002b98SHong Zhang   /*
325d4002b98SHong Zhang      determine if any processor has disassembled, if so we must
3266aad120cSJose E. Roman      also disassemble ourselves, in order that we may reassemble.
327d4002b98SHong Zhang   */
328d4002b98SHong Zhang   /*
329d4002b98SHong Zhang      if nonzero structure of submatrix B cannot change then we know that
330d4002b98SHong Zhang      no processor disassembled thus we can skip this stuff
331d4002b98SHong Zhang   */
332d4002b98SHong Zhang   if (!((Mat_SeqSELL *)sell->B->data)->nonew) {
333462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(&mat->was_assembled, &other_disassembled, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)mat)));
3342d1451d4SHong Zhang     if (mat->was_assembled && !other_disassembled) PetscCall(MatDisAssemble_MPISELL(mat));
335d4002b98SHong Zhang   }
33648a46eb9SPierre Jolivet   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) PetscCall(MatSetUpMultiply_MPISELL(mat));
3372d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
3382d1451d4SHong Zhang   if (mat->offloadmask == PETSC_OFFLOAD_CPU && sell->B->offloadmask != PETSC_OFFLOAD_UNALLOCATED) sell->B->offloadmask = PETSC_OFFLOAD_CPU;
3392d1451d4SHong Zhang #endif
3409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sell->B, mode));
3419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sell->B, mode));
3429566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sell->rowvalues, sell->rowindices));
343f4259b30SLisandro Dalcin   sell->rowvalues = NULL;
3449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->diag));
345d4002b98SHong Zhang 
346d4002b98SHong Zhang   /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
347f4f49eeaSPierre Jolivet   if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqSELL *)sell->A->data)->nonew) {
348d4002b98SHong Zhang     PetscObjectState state = sell->A->nonzerostate + sell->B->nonzerostate;
349462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(&state, &mat->nonzerostate, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)mat)));
350d4002b98SHong Zhang   }
3512d1451d4SHong Zhang #if defined(PETSC_HAVE_CUDA)
3522d1451d4SHong Zhang   mat->offloadmask = PETSC_OFFLOAD_BOTH;
3532d1451d4SHong Zhang #endif
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
355d4002b98SHong Zhang }
356d4002b98SHong Zhang 
357ba38deedSJacob Faibussowitsch static PetscErrorCode MatZeroEntries_MPISELL(Mat A)
358d71ae5a4SJacob Faibussowitsch {
359d4002b98SHong Zhang   Mat_MPISELL *l = (Mat_MPISELL *)A->data;
360d4002b98SHong Zhang 
361d4002b98SHong Zhang   PetscFunctionBegin;
3629566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->A));
3639566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(l->B));
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
365d4002b98SHong Zhang }
366d4002b98SHong Zhang 
367ba38deedSJacob Faibussowitsch static PetscErrorCode MatMult_MPISELL(Mat A, Vec xx, Vec yy)
368d71ae5a4SJacob Faibussowitsch {
369d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
370d4002b98SHong Zhang   PetscInt     nt;
371d4002b98SHong Zhang 
372d4002b98SHong Zhang   PetscFunctionBegin;
3739566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xx, &nt));
37408401ef6SPierre 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);
3759566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3769566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->mult)(a->A, xx, yy));
3779566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3789566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B, a->lvec, yy, yy));
3793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
380d4002b98SHong Zhang }
381d4002b98SHong Zhang 
382fa078d78SJacob Faibussowitsch static PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat A, Vec bb, Vec xx)
383d71ae5a4SJacob Faibussowitsch {
384d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
385d4002b98SHong Zhang 
386d4002b98SHong Zhang   PetscFunctionBegin;
3879566063dSJacob Faibussowitsch   PetscCall(MatMultDiagonalBlock(a->A, bb, xx));
3883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
389d4002b98SHong Zhang }
390d4002b98SHong Zhang 
391ba38deedSJacob Faibussowitsch static PetscErrorCode MatMultAdd_MPISELL(Mat A, Vec xx, Vec yy, Vec zz)
392d71ae5a4SJacob Faibussowitsch {
393d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
394d4002b98SHong Zhang 
395d4002b98SHong Zhang   PetscFunctionBegin;
3969566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3979566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multadd)(a->A, xx, yy, zz));
3989566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
3999566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multadd)(a->B, a->lvec, zz, zz));
4003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
401d4002b98SHong Zhang }
402d4002b98SHong Zhang 
403ba38deedSJacob Faibussowitsch static PetscErrorCode MatMultTranspose_MPISELL(Mat A, Vec xx, Vec yy)
404d71ae5a4SJacob Faibussowitsch {
405d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
406d4002b98SHong Zhang 
407d4002b98SHong Zhang   PetscFunctionBegin;
408d4002b98SHong Zhang   /* do nondiagonal part */
4099566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B, xx, a->lvec));
410d4002b98SHong Zhang   /* do local part */
4119566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multtranspose)(a->A, xx, yy));
412a29b4eb7SJunchao Zhang   /* add partial results together */
4139566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, a->lvec, yy, ADD_VALUES, SCATTER_REVERSE));
4149566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, a->lvec, yy, ADD_VALUES, SCATTER_REVERSE));
4153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
416d4002b98SHong Zhang }
417d4002b98SHong Zhang 
418ba38deedSJacob Faibussowitsch static PetscErrorCode MatIsTranspose_MPISELL(Mat Amat, Mat Bmat, PetscReal tol, PetscBool *f)
419d71ae5a4SJacob Faibussowitsch {
420d4002b98SHong Zhang   MPI_Comm     comm;
421d4002b98SHong Zhang   Mat_MPISELL *Asell = (Mat_MPISELL *)Amat->data, *Bsell;
422d4002b98SHong Zhang   Mat          Adia  = Asell->A, Bdia, Aoff, Boff, *Aoffs, *Boffs;
423d4002b98SHong Zhang   IS           Me, Notme;
424d4002b98SHong Zhang   PetscInt     M, N, first, last, *notme, i;
425d4002b98SHong Zhang   PetscMPIInt  size;
426d4002b98SHong Zhang 
427d4002b98SHong Zhang   PetscFunctionBegin;
428d4002b98SHong Zhang   /* Easy test: symmetric diagonal block */
4299371c9d4SSatish Balay   Bsell = (Mat_MPISELL *)Bmat->data;
4309371c9d4SSatish Balay   Bdia  = Bsell->A;
4319566063dSJacob Faibussowitsch   PetscCall(MatIsTranspose(Adia, Bdia, tol, f));
4323ba16761SJacob Faibussowitsch   if (!*f) PetscFunctionReturn(PETSC_SUCCESS);
4339566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
4349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
4353ba16761SJacob Faibussowitsch   if (size == 1) PetscFunctionReturn(PETSC_SUCCESS);
436d4002b98SHong Zhang 
437d4002b98SHong Zhang   /* Hard test: off-diagonal block. This takes a MatCreateSubMatrix. */
4389566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Amat, &M, &N));
4399566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Amat, &first, &last));
4409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N - last + first, &notme));
441d4002b98SHong Zhang   for (i = 0; i < first; i++) notme[i] = i;
442d4002b98SHong Zhang   for (i = last; i < M; i++) notme[i - last + first] = i;
4439566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(MPI_COMM_SELF, N - last + first, notme, PETSC_COPY_VALUES, &Notme));
4449566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(MPI_COMM_SELF, last - first, first, 1, &Me));
4459566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Amat, 1, &Me, &Notme, MAT_INITIAL_MATRIX, &Aoffs));
446d4002b98SHong Zhang   Aoff = Aoffs[0];
4479566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Bmat, 1, &Notme, &Me, MAT_INITIAL_MATRIX, &Boffs));
448d4002b98SHong Zhang   Boff = Boffs[0];
4499566063dSJacob Faibussowitsch   PetscCall(MatIsTranspose(Aoff, Boff, tol, f));
4509566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1, &Aoffs));
4519566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1, &Boffs));
4529566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Me));
4539566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Notme));
4549566063dSJacob Faibussowitsch   PetscCall(PetscFree(notme));
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
456d4002b98SHong Zhang }
457d4002b98SHong Zhang 
458ba38deedSJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_MPISELL(Mat A, Vec xx, Vec yy, Vec zz)
459d71ae5a4SJacob Faibussowitsch {
460d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
461d4002b98SHong Zhang 
462d4002b98SHong Zhang   PetscFunctionBegin;
463d4002b98SHong Zhang   /* do nondiagonal part */
4649566063dSJacob Faibussowitsch   PetscCall((*a->B->ops->multtranspose)(a->B, xx, a->lvec));
465d4002b98SHong Zhang   /* do local part */
4669566063dSJacob Faibussowitsch   PetscCall((*a->A->ops->multtransposeadd)(a->A, xx, yy, zz));
467e4a140f6SJunchao Zhang   /* add partial results together */
4689566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(a->Mvctx, a->lvec, zz, ADD_VALUES, SCATTER_REVERSE));
4699566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(a->Mvctx, a->lvec, zz, ADD_VALUES, SCATTER_REVERSE));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471d4002b98SHong Zhang }
472d4002b98SHong Zhang 
473d4002b98SHong Zhang /*
474d4002b98SHong Zhang   This only works correctly for square matrices where the subblock A->A is the
475d4002b98SHong Zhang    diagonal block
476d4002b98SHong Zhang */
477ba38deedSJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_MPISELL(Mat A, Vec v)
478d71ae5a4SJacob Faibussowitsch {
479d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
480d4002b98SHong Zhang 
481d4002b98SHong Zhang   PetscFunctionBegin;
48208401ef6SPierre Jolivet   PetscCheck(A->rmap->N == A->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Supports only square matrix where A->A is diag block");
483aed4548fSBarry 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");
4849566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(a->A, v));
4853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
486d4002b98SHong Zhang }
487d4002b98SHong Zhang 
488ba38deedSJacob Faibussowitsch static PetscErrorCode MatScale_MPISELL(Mat A, PetscScalar aa)
489d71ae5a4SJacob Faibussowitsch {
490d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
491d4002b98SHong Zhang 
492d4002b98SHong Zhang   PetscFunctionBegin;
4939566063dSJacob Faibussowitsch   PetscCall(MatScale(a->A, aa));
4949566063dSJacob Faibussowitsch   PetscCall(MatScale(a->B, aa));
4953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
496d4002b98SHong Zhang }
497d4002b98SHong Zhang 
498d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_MPISELL(Mat mat)
499d71ae5a4SJacob Faibussowitsch {
500d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
501d4002b98SHong Zhang 
502d4002b98SHong Zhang   PetscFunctionBegin;
5033ba16761SJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT, mat->rmap->N, mat->cmap->N));
5049566063dSJacob Faibussowitsch   PetscCall(MatStashDestroy_Private(&mat->stash));
5059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->diag));
5069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sell->A));
5079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sell->B));
508d4002b98SHong Zhang #if defined(PETSC_USE_CTABLE)
509eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&sell->colmap));
510d4002b98SHong Zhang #else
5119566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->colmap));
512d4002b98SHong Zhang #endif
5139566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->garray));
5149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sell->lvec));
5159566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sell->Mvctx));
5169566063dSJacob Faibussowitsch   PetscCall(PetscFree2(sell->rowvalues, sell->rowindices));
5179566063dSJacob Faibussowitsch   PetscCall(PetscFree(sell->ld));
5189566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
519d4002b98SHong Zhang 
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL));
5219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatStoreValues_C", NULL));
5229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatRetrieveValues_C", NULL));
5239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatIsTranspose_C", NULL));
5249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPISELLSetPreallocation_C", NULL));
5259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisell_mpiaij_C", NULL));
526b5917f1bSHong Zhang #if defined(PETSC_HAVE_CUDA)
527b5917f1bSHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisell_mpisellcuda_C", NULL));
528b5917f1bSHong Zhang #endif
5299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalScaleLocal_C", NULL));
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
531d4002b98SHong Zhang }
532d4002b98SHong Zhang 
533d4002b98SHong Zhang #include <petscdraw.h>
534ba38deedSJacob Faibussowitsch static PetscErrorCode MatView_MPISELL_ASCIIorDraworSocket(Mat mat, PetscViewer viewer)
535d71ae5a4SJacob Faibussowitsch {
536d4002b98SHong Zhang   Mat_MPISELL      *sell = (Mat_MPISELL *)mat->data;
537d4002b98SHong Zhang   PetscMPIInt       rank = sell->rank, size = sell->size;
538d4002b98SHong Zhang   PetscBool         isdraw, iascii, isbinary;
539d4002b98SHong Zhang   PetscViewer       sviewer;
540d4002b98SHong Zhang   PetscViewerFormat format;
541d4002b98SHong Zhang 
542d4002b98SHong Zhang   PetscFunctionBegin;
5439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
546d4002b98SHong Zhang   if (iascii) {
5479566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
548d4002b98SHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
549d4002b98SHong Zhang       MatInfo   info;
5506335e310SSatish Balay       PetscInt *inodes;
551d4002b98SHong Zhang 
5529566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat), &rank));
5539566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(mat, MAT_LOCAL, &info));
5549566063dSJacob Faibussowitsch       PetscCall(MatInodeGetInodeSizes(sell->A, NULL, &inodes, NULL));
5559566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
556d4002b98SHong Zhang       if (!inodes) {
5579371c9d4SSatish 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,
5589371c9d4SSatish Balay                                                      (PetscInt)info.nz_allocated, (PetscInt)info.memory));
559d4002b98SHong Zhang       } else {
5609371c9d4SSatish 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,
5619371c9d4SSatish Balay                                                      (PetscInt)info.nz_allocated, (PetscInt)info.memory));
562d4002b98SHong Zhang       }
5639566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(sell->A, MAT_LOCAL, &info));
5649566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] on-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
5659566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(sell->B, MAT_LOCAL, &info));
5669566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] off-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
5679566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
5689566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
5699566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Information on VecScatter used in matrix-vector product: \n"));
5709566063dSJacob Faibussowitsch       PetscCall(VecScatterView(sell->Mvctx, viewer));
5713ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
572d4002b98SHong Zhang     } else if (format == PETSC_VIEWER_ASCII_INFO) {
573d4002b98SHong Zhang       PetscInt inodecount, inodelimit, *inodes;
5749566063dSJacob Faibussowitsch       PetscCall(MatInodeGetInodeSizes(sell->A, &inodecount, &inodes, &inodelimit));
575d4002b98SHong Zhang       if (inodes) {
5769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "using I-node (on process 0) routines: found %" PetscInt_FMT " nodes, limit used is %" PetscInt_FMT "\n", inodecount, inodelimit));
577d4002b98SHong Zhang       } else {
5789566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "not using I-node (on process 0) routines\n"));
579d4002b98SHong Zhang       }
5803ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
581d4002b98SHong Zhang     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
5823ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
583d4002b98SHong Zhang     }
584d4002b98SHong Zhang   } else if (isbinary) {
585d4002b98SHong Zhang     if (size == 1) {
5869566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)sell->A, ((PetscObject)mat)->name));
5879566063dSJacob Faibussowitsch       PetscCall(MatView(sell->A, viewer));
588d4002b98SHong Zhang     } else {
5899566063dSJacob Faibussowitsch       /* PetscCall(MatView_MPISELL_Binary(mat,viewer)); */
590d4002b98SHong Zhang     }
5913ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
592d4002b98SHong Zhang   } else if (isdraw) {
593d4002b98SHong Zhang     PetscDraw draw;
594d4002b98SHong Zhang     PetscBool isnull;
5959566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
5969566063dSJacob Faibussowitsch     PetscCall(PetscDrawIsNull(draw, &isnull));
5973ba16761SJacob Faibussowitsch     if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
598d4002b98SHong Zhang   }
599d4002b98SHong Zhang 
600d4002b98SHong Zhang   {
601d4002b98SHong Zhang     /* assemble the entire matrix onto first processor. */
602d4002b98SHong Zhang     Mat          A;
603d4002b98SHong Zhang     Mat_SeqSELL *Aloc;
604d4002b98SHong Zhang     PetscInt     M = mat->rmap->N, N = mat->cmap->N, *acolidx, row, col, i, j;
605d4002b98SHong Zhang     MatScalar   *aval;
606d4002b98SHong Zhang     PetscBool    isnonzero;
607d4002b98SHong Zhang 
6089566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &A));
609dd400576SPatrick Sanan     if (rank == 0) {
6109566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A, M, N, M, N));
611d4002b98SHong Zhang     } else {
6129566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(A, 0, 0, M, N));
613d4002b98SHong Zhang     }
614d4002b98SHong Zhang     /* This is just a temporary matrix, so explicitly using MATMPISELL is probably best */
6159566063dSJacob Faibussowitsch     PetscCall(MatSetType(A, MATMPISELL));
6169566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(A, 0, NULL, 0, NULL));
6179566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
618d4002b98SHong Zhang 
619d4002b98SHong Zhang     /* copy over the A part */
620d4002b98SHong Zhang     Aloc    = (Mat_SeqSELL *)sell->A->data;
6219371c9d4SSatish Balay     acolidx = Aloc->colidx;
6229371c9d4SSatish Balay     aval    = Aloc->val;
623d4002b98SHong Zhang     for (i = 0; i < Aloc->totalslices; i++) { /* loop over slices */
624d4002b98SHong Zhang       for (j = Aloc->sliidx[i]; j < Aloc->sliidx[i + 1]; j++) {
62507e43b41SHong Zhang         isnonzero = (PetscBool)((j - Aloc->sliidx[i]) / Aloc->sliceheight < Aloc->rlen[i * Aloc->sliceheight + j % Aloc->sliceheight]);
626d4002b98SHong Zhang         if (isnonzero) { /* check the mask bit */
62707e43b41SHong Zhang           row = i * Aloc->sliceheight + j % Aloc->sliceheight + mat->rmap->rstart;
628d4002b98SHong Zhang           col = *acolidx + mat->rmap->rstart;
6299566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A, 1, &row, 1, &col, aval, INSERT_VALUES));
630d4002b98SHong Zhang         }
6319371c9d4SSatish Balay         aval++;
6329371c9d4SSatish Balay         acolidx++;
633d4002b98SHong Zhang       }
634d4002b98SHong Zhang     }
635d4002b98SHong Zhang 
636d4002b98SHong Zhang     /* copy over the B part */
637d4002b98SHong Zhang     Aloc    = (Mat_SeqSELL *)sell->B->data;
6389371c9d4SSatish Balay     acolidx = Aloc->colidx;
6399371c9d4SSatish Balay     aval    = Aloc->val;
640d4002b98SHong Zhang     for (i = 0; i < Aloc->totalslices; i++) {
641d4002b98SHong Zhang       for (j = Aloc->sliidx[i]; j < Aloc->sliidx[i + 1]; j++) {
64207e43b41SHong Zhang         isnonzero = (PetscBool)((j - Aloc->sliidx[i]) / Aloc->sliceheight < Aloc->rlen[i * Aloc->sliceheight + j % Aloc->sliceheight]);
643d4002b98SHong Zhang         if (isnonzero) {
64407e43b41SHong Zhang           row = i * Aloc->sliceheight + j % Aloc->sliceheight + mat->rmap->rstart;
645d4002b98SHong Zhang           col = sell->garray[*acolidx];
6469566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A, 1, &row, 1, &col, aval, INSERT_VALUES));
647d4002b98SHong Zhang         }
6489371c9d4SSatish Balay         aval++;
6499371c9d4SSatish Balay         acolidx++;
650d4002b98SHong Zhang       }
651d4002b98SHong Zhang     }
652d4002b98SHong Zhang 
6539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
6549566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
655d4002b98SHong Zhang     /*
656d4002b98SHong Zhang        Everyone has to call to draw the matrix since the graphics waits are
657d4002b98SHong Zhang        synchronized across all processors that share the PetscDraw object
658d4002b98SHong Zhang     */
6599566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
660dd400576SPatrick Sanan     if (rank == 0) {
661f4f49eeaSPierre Jolivet       PetscCall(PetscObjectSetName((PetscObject)((Mat_MPISELL *)A->data)->A, ((PetscObject)mat)->name));
662f4f49eeaSPierre Jolivet       PetscCall(MatView_SeqSELL(((Mat_MPISELL *)A->data)->A, sviewer));
663d4002b98SHong Zhang     }
6649566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
6659566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A));
666d4002b98SHong Zhang   }
6673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
668d4002b98SHong Zhang }
669d4002b98SHong Zhang 
670ba38deedSJacob Faibussowitsch static PetscErrorCode MatView_MPISELL(Mat mat, PetscViewer viewer)
671d71ae5a4SJacob Faibussowitsch {
672d4002b98SHong Zhang   PetscBool iascii, isdraw, issocket, isbinary;
673d4002b98SHong Zhang 
674d4002b98SHong Zhang   PetscFunctionBegin;
6759566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
6769566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
6779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
6789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSOCKET, &issocket));
67948a46eb9SPierre Jolivet   if (iascii || isdraw || isbinary || issocket) PetscCall(MatView_MPISELL_ASCIIorDraworSocket(mat, viewer));
6803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
681d4002b98SHong Zhang }
682d4002b98SHong Zhang 
683ba38deedSJacob Faibussowitsch static PetscErrorCode MatGetGhosts_MPISELL(Mat mat, PetscInt *nghosts, const PetscInt *ghosts[])
684d71ae5a4SJacob Faibussowitsch {
685d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
686d4002b98SHong Zhang 
687d4002b98SHong Zhang   PetscFunctionBegin;
6889566063dSJacob Faibussowitsch   PetscCall(MatGetSize(sell->B, NULL, nghosts));
689d4002b98SHong Zhang   if (ghosts) *ghosts = sell->garray;
6903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
691d4002b98SHong Zhang }
692d4002b98SHong Zhang 
693ba38deedSJacob Faibussowitsch static PetscErrorCode MatGetInfo_MPISELL(Mat matin, MatInfoType flag, MatInfo *info)
694d71ae5a4SJacob Faibussowitsch {
695d4002b98SHong Zhang   Mat_MPISELL   *mat = (Mat_MPISELL *)matin->data;
696d4002b98SHong Zhang   Mat            A = mat->A, B = mat->B;
6973966268fSBarry Smith   PetscLogDouble isend[5], irecv[5];
698d4002b98SHong Zhang 
699d4002b98SHong Zhang   PetscFunctionBegin;
700d4002b98SHong Zhang   info->block_size = 1.0;
7019566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_LOCAL, info));
702d4002b98SHong Zhang 
7039371c9d4SSatish Balay   isend[0] = info->nz_used;
7049371c9d4SSatish Balay   isend[1] = info->nz_allocated;
7059371c9d4SSatish Balay   isend[2] = info->nz_unneeded;
7069371c9d4SSatish Balay   isend[3] = info->memory;
7079371c9d4SSatish Balay   isend[4] = info->mallocs;
708d4002b98SHong Zhang 
7099566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(B, MAT_LOCAL, info));
710d4002b98SHong Zhang 
7119371c9d4SSatish Balay   isend[0] += info->nz_used;
7129371c9d4SSatish Balay   isend[1] += info->nz_allocated;
7139371c9d4SSatish Balay   isend[2] += info->nz_unneeded;
7149371c9d4SSatish Balay   isend[3] += info->memory;
7159371c9d4SSatish Balay   isend[4] += info->mallocs;
716d4002b98SHong Zhang   if (flag == MAT_LOCAL) {
717d4002b98SHong Zhang     info->nz_used      = isend[0];
718d4002b98SHong Zhang     info->nz_allocated = isend[1];
719d4002b98SHong Zhang     info->nz_unneeded  = isend[2];
720d4002b98SHong Zhang     info->memory       = isend[3];
721d4002b98SHong Zhang     info->mallocs      = isend[4];
722d4002b98SHong Zhang   } else if (flag == MAT_GLOBAL_MAX) {
723462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_MAX, PetscObjectComm((PetscObject)matin)));
724d4002b98SHong Zhang 
725d4002b98SHong Zhang     info->nz_used      = irecv[0];
726d4002b98SHong Zhang     info->nz_allocated = irecv[1];
727d4002b98SHong Zhang     info->nz_unneeded  = irecv[2];
728d4002b98SHong Zhang     info->memory       = irecv[3];
729d4002b98SHong Zhang     info->mallocs      = irecv[4];
730d4002b98SHong Zhang   } else if (flag == MAT_GLOBAL_SUM) {
731462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_SUM, PetscObjectComm((PetscObject)matin)));
732d4002b98SHong Zhang 
733d4002b98SHong Zhang     info->nz_used      = irecv[0];
734d4002b98SHong Zhang     info->nz_allocated = irecv[1];
735d4002b98SHong Zhang     info->nz_unneeded  = irecv[2];
736d4002b98SHong Zhang     info->memory       = irecv[3];
737d4002b98SHong Zhang     info->mallocs      = irecv[4];
738d4002b98SHong Zhang   }
739d4002b98SHong Zhang   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
740d4002b98SHong Zhang   info->fill_ratio_needed = 0;
741d4002b98SHong Zhang   info->factor_mallocs    = 0;
7423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
743d4002b98SHong Zhang }
744d4002b98SHong Zhang 
74543b34f9dSJacob Faibussowitsch static PetscErrorCode MatSetOption_MPISELL(Mat A, MatOption op, PetscBool flg)
746d71ae5a4SJacob Faibussowitsch {
747d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
748d4002b98SHong Zhang 
749d4002b98SHong Zhang   PetscFunctionBegin;
750d4002b98SHong Zhang   switch (op) {
751d4002b98SHong Zhang   case MAT_NEW_NONZERO_LOCATIONS:
752d4002b98SHong Zhang   case MAT_NEW_NONZERO_ALLOCATION_ERR:
753d4002b98SHong Zhang   case MAT_UNUSED_NONZERO_LOCATION_ERR:
754d4002b98SHong Zhang   case MAT_KEEP_NONZERO_PATTERN:
755d4002b98SHong Zhang   case MAT_NEW_NONZERO_LOCATION_ERR:
756d4002b98SHong Zhang   case MAT_USE_INODES:
757d4002b98SHong Zhang   case MAT_IGNORE_ZERO_ENTRIES:
758d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7599566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
7609566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B, op, flg));
761d4002b98SHong Zhang     break;
762d4002b98SHong Zhang   case MAT_ROW_ORIENTED:
763d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
764d4002b98SHong Zhang     a->roworiented = flg;
765d4002b98SHong Zhang 
7669566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
7679566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->B, op, flg));
768d4002b98SHong Zhang     break;
769d71ae5a4SJacob Faibussowitsch   case MAT_IGNORE_OFF_PROC_ENTRIES:
770d71ae5a4SJacob Faibussowitsch     a->donotstash = flg;
771d71ae5a4SJacob Faibussowitsch     break;
772d4002b98SHong Zhang   case MAT_SYMMETRIC:
773d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7749566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
775d4002b98SHong Zhang     break;
776d4002b98SHong Zhang   case MAT_STRUCTURALLY_SYMMETRIC:
777d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7789566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
779d4002b98SHong Zhang     break;
780d4002b98SHong Zhang   case MAT_HERMITIAN:
781d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7829566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
783d4002b98SHong Zhang     break;
784d4002b98SHong Zhang   case MAT_SYMMETRY_ETERNAL:
785d4002b98SHong Zhang     MatCheckPreallocated(A, 1);
7869566063dSJacob Faibussowitsch     PetscCall(MatSetOption(a->A, op, flg));
787d4002b98SHong Zhang     break;
788b94d7dedSBarry Smith   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
789b94d7dedSBarry Smith     MatCheckPreallocated(A, 1);
790b94d7dedSBarry Smith     PetscCall(MatSetOption(a->A, op, flg));
791b94d7dedSBarry Smith     break;
792d71ae5a4SJacob Faibussowitsch   default:
793888c827cSStefano Zampini     break;
794d4002b98SHong Zhang   }
7953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
796d4002b98SHong Zhang }
797d4002b98SHong Zhang 
798ba38deedSJacob Faibussowitsch static PetscErrorCode MatDiagonalScale_MPISELL(Mat mat, Vec ll, Vec rr)
799d71ae5a4SJacob Faibussowitsch {
800d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
801d4002b98SHong Zhang   Mat          a = sell->A, b = sell->B;
802d4002b98SHong Zhang   PetscInt     s1, s2, s3;
803d4002b98SHong Zhang 
804d4002b98SHong Zhang   PetscFunctionBegin;
8059566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mat, &s2, &s3));
806d4002b98SHong Zhang   if (rr) {
8079566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(rr, &s1));
80808401ef6SPierre Jolivet     PetscCheck(s1 == s3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "right vector non-conforming local size");
809d4002b98SHong Zhang     /* Overlap communication with computation. */
8109566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(sell->Mvctx, rr, sell->lvec, INSERT_VALUES, SCATTER_FORWARD));
811d4002b98SHong Zhang   }
812d4002b98SHong Zhang   if (ll) {
8139566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(ll, &s1));
81408401ef6SPierre Jolivet     PetscCheck(s1 == s2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "left vector non-conforming local size");
815dbbe0bcdSBarry Smith     PetscUseTypeMethod(b, diagonalscale, ll, NULL);
816d4002b98SHong Zhang   }
817d4002b98SHong Zhang   /* scale  the diagonal block */
818dbbe0bcdSBarry Smith   PetscUseTypeMethod(a, diagonalscale, ll, rr);
819d4002b98SHong Zhang 
820d4002b98SHong Zhang   if (rr) {
821d4002b98SHong Zhang     /* Do a scatter end and then right scale the off-diagonal block */
8229566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(sell->Mvctx, rr, sell->lvec, INSERT_VALUES, SCATTER_FORWARD));
823dbbe0bcdSBarry Smith     PetscUseTypeMethod(b, diagonalscale, NULL, sell->lvec);
824d4002b98SHong Zhang   }
8253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
826d4002b98SHong Zhang }
827d4002b98SHong Zhang 
828ba38deedSJacob Faibussowitsch static PetscErrorCode MatSetUnfactored_MPISELL(Mat A)
829d71ae5a4SJacob Faibussowitsch {
830d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
831d4002b98SHong Zhang 
832d4002b98SHong Zhang   PetscFunctionBegin;
8339566063dSJacob Faibussowitsch   PetscCall(MatSetUnfactored(a->A));
8343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
835d4002b98SHong Zhang }
836d4002b98SHong Zhang 
837ba38deedSJacob Faibussowitsch static PetscErrorCode MatEqual_MPISELL(Mat A, Mat B, PetscBool *flag)
838d71ae5a4SJacob Faibussowitsch {
839d4002b98SHong Zhang   Mat_MPISELL *matB = (Mat_MPISELL *)B->data, *matA = (Mat_MPISELL *)A->data;
840d4002b98SHong Zhang   Mat          a, b, c, d;
841d4002b98SHong Zhang   PetscBool    flg;
842d4002b98SHong Zhang 
843d4002b98SHong Zhang   PetscFunctionBegin;
8449371c9d4SSatish Balay   a = matA->A;
8459371c9d4SSatish Balay   b = matA->B;
8469371c9d4SSatish Balay   c = matB->A;
8479371c9d4SSatish Balay   d = matB->B;
848d4002b98SHong Zhang 
8499566063dSJacob Faibussowitsch   PetscCall(MatEqual(a, c, &flg));
85048a46eb9SPierre Jolivet   if (flg) PetscCall(MatEqual(b, d, &flg));
851462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&flg, flag, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
8523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
853d4002b98SHong Zhang }
854d4002b98SHong Zhang 
855ba38deedSJacob Faibussowitsch static PetscErrorCode MatCopy_MPISELL(Mat A, Mat B, MatStructure str)
856d71ae5a4SJacob Faibussowitsch {
857d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
858d4002b98SHong Zhang   Mat_MPISELL *b = (Mat_MPISELL *)B->data;
859d4002b98SHong Zhang 
860d4002b98SHong Zhang   PetscFunctionBegin;
861d4002b98SHong Zhang   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
862d4002b98SHong Zhang   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
863d4002b98SHong Zhang     /* because of the column compression in the off-processor part of the matrix a->B,
864d4002b98SHong Zhang        the number of columns in a->B and b->B may be different, hence we cannot call
865d4002b98SHong Zhang        the MatCopy() directly on the two parts. If need be, we can provide a more
866d4002b98SHong Zhang        efficient copy than the MatCopy_Basic() by first uncompressing the a->B matrices
867d4002b98SHong Zhang        then copying the submatrices */
8689566063dSJacob Faibussowitsch     PetscCall(MatCopy_Basic(A, B, str));
869d4002b98SHong Zhang   } else {
8709566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->A, b->A, str));
8719566063dSJacob Faibussowitsch     PetscCall(MatCopy(a->B, b->B, str));
872d4002b98SHong Zhang   }
8733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
874d4002b98SHong Zhang }
875d4002b98SHong Zhang 
876ba38deedSJacob Faibussowitsch static PetscErrorCode MatSetUp_MPISELL(Mat A)
877d71ae5a4SJacob Faibussowitsch {
878d4002b98SHong Zhang   PetscFunctionBegin;
8799566063dSJacob Faibussowitsch   PetscCall(MatMPISELLSetPreallocation(A, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
8803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
881d4002b98SHong Zhang }
882d4002b98SHong Zhang 
883ba38deedSJacob Faibussowitsch static PetscErrorCode MatConjugate_MPISELL(Mat mat)
884d71ae5a4SJacob Faibussowitsch {
8855f80ce2aSJacob Faibussowitsch   PetscFunctionBegin;
8865f80ce2aSJacob Faibussowitsch   if (PetscDefined(USE_COMPLEX)) {
887d4002b98SHong Zhang     Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
888d4002b98SHong Zhang 
8899566063dSJacob Faibussowitsch     PetscCall(MatConjugate_SeqSELL(sell->A));
8909566063dSJacob Faibussowitsch     PetscCall(MatConjugate_SeqSELL(sell->B));
8915f80ce2aSJacob Faibussowitsch   }
8923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
893d4002b98SHong Zhang }
894d4002b98SHong Zhang 
895ba38deedSJacob Faibussowitsch static PetscErrorCode MatRealPart_MPISELL(Mat A)
896d71ae5a4SJacob Faibussowitsch {
897d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
898d4002b98SHong Zhang 
899d4002b98SHong Zhang   PetscFunctionBegin;
9009566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->A));
9019566063dSJacob Faibussowitsch   PetscCall(MatRealPart(a->B));
9023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
903d4002b98SHong Zhang }
904d4002b98SHong Zhang 
905ba38deedSJacob Faibussowitsch static PetscErrorCode MatImaginaryPart_MPISELL(Mat A)
906d71ae5a4SJacob Faibussowitsch {
907d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
908d4002b98SHong Zhang 
909d4002b98SHong Zhang   PetscFunctionBegin;
9109566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->A));
9119566063dSJacob Faibussowitsch   PetscCall(MatImaginaryPart(a->B));
9123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
913d4002b98SHong Zhang }
914d4002b98SHong Zhang 
915ba38deedSJacob Faibussowitsch static PetscErrorCode MatInvertBlockDiagonal_MPISELL(Mat A, const PetscScalar **values)
916d71ae5a4SJacob Faibussowitsch {
917d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
918d4002b98SHong Zhang 
919d4002b98SHong Zhang   PetscFunctionBegin;
9209566063dSJacob Faibussowitsch   PetscCall(MatInvertBlockDiagonal(a->A, values));
921d4002b98SHong Zhang   A->factorerrortype = a->A->factorerrortype;
9223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
923d4002b98SHong Zhang }
924d4002b98SHong Zhang 
925d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetRandom_MPISELL(Mat x, PetscRandom rctx)
926d71ae5a4SJacob Faibussowitsch {
927d4002b98SHong Zhang   Mat_MPISELL *sell = (Mat_MPISELL *)x->data;
928d4002b98SHong Zhang 
929d4002b98SHong Zhang   PetscFunctionBegin;
9309566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(sell->A, rctx));
9319566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(sell->B, rctx));
9329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
9339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
9343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
935d4002b98SHong Zhang }
936d4002b98SHong Zhang 
937ce78bad3SBarry Smith static PetscErrorCode MatSetFromOptions_MPISELL(Mat A, PetscOptionItems PetscOptionsObject)
938d71ae5a4SJacob Faibussowitsch {
939d4002b98SHong Zhang   PetscFunctionBegin;
940d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "MPISELL options");
941d0609cedSBarry Smith   PetscOptionsHeadEnd();
9423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
943d4002b98SHong Zhang }
944d4002b98SHong Zhang 
945ba38deedSJacob Faibussowitsch static PetscErrorCode MatShift_MPISELL(Mat Y, PetscScalar a)
946d71ae5a4SJacob Faibussowitsch {
947d4002b98SHong Zhang   Mat_MPISELL *msell = (Mat_MPISELL *)Y->data;
948d4002b98SHong Zhang   Mat_SeqSELL *sell  = (Mat_SeqSELL *)msell->A->data;
949d4002b98SHong Zhang 
950d4002b98SHong Zhang   PetscFunctionBegin;
951d4002b98SHong Zhang   if (!Y->preallocated) {
9529566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(Y, 1, NULL, 0, NULL));
953d4002b98SHong Zhang   } else if (!sell->nz) {
954d4002b98SHong Zhang     PetscInt nonew = sell->nonew;
9559566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(msell->A, 1, NULL));
956d4002b98SHong Zhang     sell->nonew = nonew;
957d4002b98SHong Zhang   }
9589566063dSJacob Faibussowitsch   PetscCall(MatShift_Basic(Y, a));
9593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
960d4002b98SHong Zhang }
961d4002b98SHong Zhang 
962ba38deedSJacob Faibussowitsch static PetscErrorCode MatMissingDiagonal_MPISELL(Mat A, PetscBool *missing, PetscInt *d)
963d71ae5a4SJacob Faibussowitsch {
964d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
965d4002b98SHong Zhang 
966d4002b98SHong Zhang   PetscFunctionBegin;
96708401ef6SPierre Jolivet   PetscCheck(A->rmap->n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only works for square matrices");
9689566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(a->A, missing, d));
969d4002b98SHong Zhang   if (d) {
970d4002b98SHong Zhang     PetscInt rstart;
9719566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
972d4002b98SHong Zhang     *d += rstart;
973d4002b98SHong Zhang   }
9743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
975d4002b98SHong Zhang }
976d4002b98SHong Zhang 
977ba38deedSJacob Faibussowitsch static PetscErrorCode MatGetDiagonalBlock_MPISELL(Mat A, Mat *a)
978d71ae5a4SJacob Faibussowitsch {
979d4002b98SHong Zhang   PetscFunctionBegin;
980d4002b98SHong Zhang   *a = ((Mat_MPISELL *)A->data)->A;
9813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
982d4002b98SHong Zhang }
983d4002b98SHong Zhang 
98443b34f9dSJacob Faibussowitsch static PetscErrorCode MatStoreValues_MPISELL(Mat mat)
98543b34f9dSJacob Faibussowitsch {
98643b34f9dSJacob Faibussowitsch   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
98743b34f9dSJacob Faibussowitsch 
98843b34f9dSJacob Faibussowitsch   PetscFunctionBegin;
98943b34f9dSJacob Faibussowitsch   PetscCall(MatStoreValues(sell->A));
99043b34f9dSJacob Faibussowitsch   PetscCall(MatStoreValues(sell->B));
99143b34f9dSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99243b34f9dSJacob Faibussowitsch }
99343b34f9dSJacob Faibussowitsch 
99443b34f9dSJacob Faibussowitsch static PetscErrorCode MatRetrieveValues_MPISELL(Mat mat)
99543b34f9dSJacob Faibussowitsch {
99643b34f9dSJacob Faibussowitsch   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
99743b34f9dSJacob Faibussowitsch 
99843b34f9dSJacob Faibussowitsch   PetscFunctionBegin;
99943b34f9dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(sell->A));
100043b34f9dSJacob Faibussowitsch   PetscCall(MatRetrieveValues(sell->B));
100143b34f9dSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100243b34f9dSJacob Faibussowitsch }
100343b34f9dSJacob Faibussowitsch 
100443b34f9dSJacob Faibussowitsch static PetscErrorCode MatMPISELLSetPreallocation_MPISELL(Mat B, PetscInt d_rlenmax, const PetscInt d_rlen[], PetscInt o_rlenmax, const PetscInt o_rlen[])
100543b34f9dSJacob Faibussowitsch {
100643b34f9dSJacob Faibussowitsch   Mat_MPISELL *b;
100743b34f9dSJacob Faibussowitsch 
100843b34f9dSJacob Faibussowitsch   PetscFunctionBegin;
100943b34f9dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
101043b34f9dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
101143b34f9dSJacob Faibussowitsch   b = (Mat_MPISELL *)B->data;
101243b34f9dSJacob Faibussowitsch 
101343b34f9dSJacob Faibussowitsch   if (!B->preallocated) {
101443b34f9dSJacob Faibussowitsch     /* Explicitly create 2 MATSEQSELL matrices. */
101543b34f9dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &b->A));
101643b34f9dSJacob Faibussowitsch     PetscCall(MatSetSizes(b->A, B->rmap->n, B->cmap->n, B->rmap->n, B->cmap->n));
101743b34f9dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(b->A, B, B));
101843b34f9dSJacob Faibussowitsch     PetscCall(MatSetType(b->A, MATSEQSELL));
101943b34f9dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &b->B));
102043b34f9dSJacob Faibussowitsch     PetscCall(MatSetSizes(b->B, B->rmap->n, B->cmap->N, B->rmap->n, B->cmap->N));
102143b34f9dSJacob Faibussowitsch     PetscCall(MatSetBlockSizesFromMats(b->B, B, B));
102243b34f9dSJacob Faibussowitsch     PetscCall(MatSetType(b->B, MATSEQSELL));
102343b34f9dSJacob Faibussowitsch   }
102443b34f9dSJacob Faibussowitsch 
102543b34f9dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(b->A, d_rlenmax, d_rlen));
102643b34f9dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(b->B, o_rlenmax, o_rlen));
102743b34f9dSJacob Faibussowitsch   B->preallocated  = PETSC_TRUE;
102843b34f9dSJacob Faibussowitsch   B->was_assembled = PETSC_FALSE;
102943b34f9dSJacob Faibussowitsch   /*
103043b34f9dSJacob Faibussowitsch     critical for MatAssemblyEnd to work.
103143b34f9dSJacob Faibussowitsch     MatAssemblyBegin checks it to set up was_assembled
103243b34f9dSJacob Faibussowitsch     and MatAssemblyEnd checks was_assembled to determine whether to build garray
103343b34f9dSJacob Faibussowitsch   */
103443b34f9dSJacob Faibussowitsch   B->assembled = PETSC_FALSE;
103543b34f9dSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
103643b34f9dSJacob Faibussowitsch }
103743b34f9dSJacob Faibussowitsch 
103843b34f9dSJacob Faibussowitsch static PetscErrorCode MatDuplicate_MPISELL(Mat matin, MatDuplicateOption cpvalues, Mat *newmat)
103943b34f9dSJacob Faibussowitsch {
104043b34f9dSJacob Faibussowitsch   Mat          mat;
104143b34f9dSJacob Faibussowitsch   Mat_MPISELL *a, *oldmat = (Mat_MPISELL *)matin->data;
104243b34f9dSJacob Faibussowitsch 
104343b34f9dSJacob Faibussowitsch   PetscFunctionBegin;
104443b34f9dSJacob Faibussowitsch   *newmat = NULL;
104543b34f9dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)matin), &mat));
104643b34f9dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, matin->rmap->n, matin->cmap->n, matin->rmap->N, matin->cmap->N));
104743b34f9dSJacob Faibussowitsch   PetscCall(MatSetBlockSizesFromMats(mat, matin, matin));
104843b34f9dSJacob Faibussowitsch   PetscCall(MatSetType(mat, ((PetscObject)matin)->type_name));
104943b34f9dSJacob Faibussowitsch   a = (Mat_MPISELL *)mat->data;
105043b34f9dSJacob Faibussowitsch 
105143b34f9dSJacob Faibussowitsch   mat->factortype   = matin->factortype;
105243b34f9dSJacob Faibussowitsch   mat->assembled    = PETSC_TRUE;
105343b34f9dSJacob Faibussowitsch   mat->insertmode   = NOT_SET_VALUES;
105443b34f9dSJacob Faibussowitsch   mat->preallocated = PETSC_TRUE;
105543b34f9dSJacob Faibussowitsch 
105643b34f9dSJacob Faibussowitsch   a->size         = oldmat->size;
105743b34f9dSJacob Faibussowitsch   a->rank         = oldmat->rank;
105843b34f9dSJacob Faibussowitsch   a->donotstash   = oldmat->donotstash;
105943b34f9dSJacob Faibussowitsch   a->roworiented  = oldmat->roworiented;
106043b34f9dSJacob Faibussowitsch   a->rowindices   = NULL;
106143b34f9dSJacob Faibussowitsch   a->rowvalues    = NULL;
106243b34f9dSJacob Faibussowitsch   a->getrowactive = PETSC_FALSE;
106343b34f9dSJacob Faibussowitsch 
106443b34f9dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->rmap, &mat->rmap));
106543b34f9dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(matin->cmap, &mat->cmap));
106643b34f9dSJacob Faibussowitsch 
106743b34f9dSJacob Faibussowitsch   if (oldmat->colmap) {
106843b34f9dSJacob Faibussowitsch #if defined(PETSC_USE_CTABLE)
106943b34f9dSJacob Faibussowitsch     PetscCall(PetscHMapIDuplicate(oldmat->colmap, &a->colmap));
107043b34f9dSJacob Faibussowitsch #else
107143b34f9dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mat->cmap->N, &a->colmap));
107243b34f9dSJacob Faibussowitsch     PetscCall(PetscArraycpy(a->colmap, oldmat->colmap, mat->cmap->N));
107343b34f9dSJacob Faibussowitsch #endif
107443b34f9dSJacob Faibussowitsch   } else a->colmap = NULL;
107543b34f9dSJacob Faibussowitsch   if (oldmat->garray) {
107643b34f9dSJacob Faibussowitsch     PetscInt len;
107743b34f9dSJacob Faibussowitsch     len = oldmat->B->cmap->n;
107843b34f9dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len + 1, &a->garray));
107943b34f9dSJacob Faibussowitsch     if (len) PetscCall(PetscArraycpy(a->garray, oldmat->garray, len));
108043b34f9dSJacob Faibussowitsch   } else a->garray = NULL;
108143b34f9dSJacob Faibussowitsch 
108243b34f9dSJacob Faibussowitsch   PetscCall(VecDuplicate(oldmat->lvec, &a->lvec));
108343b34f9dSJacob Faibussowitsch   PetscCall(VecScatterCopy(oldmat->Mvctx, &a->Mvctx));
108443b34f9dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->A, cpvalues, &a->A));
108543b34f9dSJacob Faibussowitsch   PetscCall(MatDuplicate(oldmat->B, cpvalues, &a->B));
108643b34f9dSJacob Faibussowitsch   PetscCall(PetscFunctionListDuplicate(((PetscObject)matin)->qlist, &((PetscObject)mat)->qlist));
108743b34f9dSJacob Faibussowitsch   *newmat = mat;
108843b34f9dSJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108943b34f9dSJacob Faibussowitsch }
109043b34f9dSJacob Faibussowitsch 
109143b34f9dSJacob Faibussowitsch static const struct _MatOps MatOps_Values = {MatSetValues_MPISELL,
1092f4259b30SLisandro Dalcin                                              NULL,
1093f4259b30SLisandro Dalcin                                              NULL,
1094d4002b98SHong Zhang                                              MatMult_MPISELL,
1095d4002b98SHong Zhang                                              /* 4*/ MatMultAdd_MPISELL,
1096d4002b98SHong Zhang                                              MatMultTranspose_MPISELL,
1097d4002b98SHong Zhang                                              MatMultTransposeAdd_MPISELL,
1098f4259b30SLisandro Dalcin                                              NULL,
1099f4259b30SLisandro Dalcin                                              NULL,
1100f4259b30SLisandro Dalcin                                              NULL,
1101f4259b30SLisandro Dalcin                                              /*10*/ NULL,
1102f4259b30SLisandro Dalcin                                              NULL,
1103f4259b30SLisandro Dalcin                                              NULL,
1104d4002b98SHong Zhang                                              MatSOR_MPISELL,
1105f4259b30SLisandro Dalcin                                              NULL,
1106d4002b98SHong Zhang                                              /*15*/ MatGetInfo_MPISELL,
1107d4002b98SHong Zhang                                              MatEqual_MPISELL,
1108d4002b98SHong Zhang                                              MatGetDiagonal_MPISELL,
1109d4002b98SHong Zhang                                              MatDiagonalScale_MPISELL,
1110f4259b30SLisandro Dalcin                                              NULL,
1111d4002b98SHong Zhang                                              /*20*/ MatAssemblyBegin_MPISELL,
1112d4002b98SHong Zhang                                              MatAssemblyEnd_MPISELL,
1113d4002b98SHong Zhang                                              MatSetOption_MPISELL,
1114d4002b98SHong Zhang                                              MatZeroEntries_MPISELL,
1115f4259b30SLisandro Dalcin                                              /*24*/ NULL,
1116f4259b30SLisandro Dalcin                                              NULL,
1117f4259b30SLisandro Dalcin                                              NULL,
1118f4259b30SLisandro Dalcin                                              NULL,
1119f4259b30SLisandro Dalcin                                              NULL,
1120d4002b98SHong Zhang                                              /*29*/ MatSetUp_MPISELL,
1121f4259b30SLisandro Dalcin                                              NULL,
1122f4259b30SLisandro Dalcin                                              NULL,
1123d4002b98SHong Zhang                                              MatGetDiagonalBlock_MPISELL,
1124f4259b30SLisandro Dalcin                                              NULL,
1125d4002b98SHong Zhang                                              /*34*/ MatDuplicate_MPISELL,
1126f4259b30SLisandro Dalcin                                              NULL,
1127f4259b30SLisandro Dalcin                                              NULL,
1128f4259b30SLisandro Dalcin                                              NULL,
1129f4259b30SLisandro Dalcin                                              NULL,
1130f4259b30SLisandro Dalcin                                              /*39*/ NULL,
1131f4259b30SLisandro Dalcin                                              NULL,
1132f4259b30SLisandro Dalcin                                              NULL,
1133d4002b98SHong Zhang                                              MatGetValues_MPISELL,
1134d4002b98SHong Zhang                                              MatCopy_MPISELL,
1135f4259b30SLisandro Dalcin                                              /*44*/ NULL,
1136d4002b98SHong Zhang                                              MatScale_MPISELL,
1137d4002b98SHong Zhang                                              MatShift_MPISELL,
1138d4002b98SHong Zhang                                              MatDiagonalSet_MPISELL,
1139f4259b30SLisandro Dalcin                                              NULL,
1140d4002b98SHong Zhang                                              /*49*/ MatSetRandom_MPISELL,
1141f4259b30SLisandro Dalcin                                              NULL,
1142f4259b30SLisandro Dalcin                                              NULL,
1143f4259b30SLisandro Dalcin                                              NULL,
1144f4259b30SLisandro Dalcin                                              NULL,
1145d4002b98SHong Zhang                                              /*54*/ MatFDColoringCreate_MPIXAIJ,
1146f4259b30SLisandro Dalcin                                              NULL,
1147d4002b98SHong Zhang                                              MatSetUnfactored_MPISELL,
1148f4259b30SLisandro Dalcin                                              NULL,
1149f4259b30SLisandro Dalcin                                              NULL,
1150f4259b30SLisandro Dalcin                                              /*59*/ NULL,
1151d4002b98SHong Zhang                                              MatDestroy_MPISELL,
1152d4002b98SHong Zhang                                              MatView_MPISELL,
1153f4259b30SLisandro Dalcin                                              NULL,
1154f4259b30SLisandro Dalcin                                              NULL,
1155f4259b30SLisandro Dalcin                                              /*64*/ NULL,
1156f4259b30SLisandro Dalcin                                              NULL,
1157f4259b30SLisandro Dalcin                                              NULL,
1158f4259b30SLisandro Dalcin                                              NULL,
1159f4259b30SLisandro Dalcin                                              NULL,
1160f4259b30SLisandro Dalcin                                              /*69*/ NULL,
1161f4259b30SLisandro Dalcin                                              NULL,
1162f4259b30SLisandro Dalcin                                              NULL,
1163*8bb0f5c6SPierre Jolivet                                              MatFDColoringApply_AIJ, /* reuse AIJ function */
1164d4002b98SHong Zhang                                              MatSetFromOptions_MPISELL,
1165f4259b30SLisandro Dalcin                                              NULL,
1166*8bb0f5c6SPierre Jolivet                                              /*75*/ NULL,
1167*8bb0f5c6SPierre Jolivet                                              NULL,
1168*8bb0f5c6SPierre Jolivet                                              NULL,
1169f4259b30SLisandro Dalcin                                              NULL,
1170f4259b30SLisandro Dalcin                                              NULL,
1171f4259b30SLisandro Dalcin                                              /*80*/ NULL,
1172f4259b30SLisandro Dalcin                                              NULL,
1173f4259b30SLisandro Dalcin                                              NULL,
1174f4259b30SLisandro Dalcin                                              /*83*/ NULL,
1175f4259b30SLisandro Dalcin                                              NULL,
1176f4259b30SLisandro Dalcin                                              NULL,
1177f4259b30SLisandro Dalcin                                              NULL,
1178f4259b30SLisandro Dalcin                                              NULL,
1179f4259b30SLisandro Dalcin                                              NULL,
1180f4259b30SLisandro Dalcin                                              /*89*/ NULL,
1181f4259b30SLisandro Dalcin                                              NULL,
1182f4259b30SLisandro Dalcin                                              NULL,
1183f4259b30SLisandro Dalcin                                              NULL,
1184*8bb0f5c6SPierre Jolivet                                              MatConjugate_MPISELL,
1185f4259b30SLisandro Dalcin                                              /*94*/ NULL,
1186f4259b30SLisandro Dalcin                                              NULL,
1187*8bb0f5c6SPierre Jolivet                                              MatRealPart_MPISELL,
1188*8bb0f5c6SPierre Jolivet                                              MatImaginaryPart_MPISELL,
1189f4259b30SLisandro Dalcin                                              NULL,
1190f4259b30SLisandro Dalcin                                              /*99*/ NULL,
1191f4259b30SLisandro Dalcin                                              NULL,
1192f4259b30SLisandro Dalcin                                              NULL,
1193f4259b30SLisandro Dalcin                                              NULL,
1194f4259b30SLisandro Dalcin                                              NULL,
1195*8bb0f5c6SPierre Jolivet                                              /*104*/ MatMissingDiagonal_MPISELL,
1196f4259b30SLisandro Dalcin                                              NULL,
1197f4259b30SLisandro Dalcin                                              NULL,
1198d4002b98SHong Zhang                                              MatGetGhosts_MPISELL,
1199f4259b30SLisandro Dalcin                                              NULL,
1200*8bb0f5c6SPierre Jolivet                                              /*109*/ NULL,
1201*8bb0f5c6SPierre Jolivet                                              MatMultDiagonalBlock_MPISELL,
1202f4259b30SLisandro Dalcin                                              NULL,
1203*8bb0f5c6SPierre Jolivet                                              NULL,
1204*8bb0f5c6SPierre Jolivet                                              NULL,
1205*8bb0f5c6SPierre Jolivet                                              /*114*/ NULL,
1206*8bb0f5c6SPierre Jolivet                                              NULL,
1207*8bb0f5c6SPierre Jolivet                                              NULL,
1208*8bb0f5c6SPierre Jolivet                                              MatInvertBlockDiagonal_MPISELL,
1209*8bb0f5c6SPierre Jolivet                                              NULL,
1210*8bb0f5c6SPierre Jolivet                                              /*119*/ NULL,
1211f4259b30SLisandro Dalcin                                              NULL,
1212f4259b30SLisandro Dalcin                                              NULL,
1213f4259b30SLisandro Dalcin                                              NULL,
1214f4259b30SLisandro Dalcin                                              NULL,
1215f4259b30SLisandro Dalcin                                              /*124*/ NULL,
1216f4259b30SLisandro Dalcin                                              NULL,
1217f4259b30SLisandro Dalcin                                              NULL,
1218f4259b30SLisandro Dalcin                                              NULL,
1219*8bb0f5c6SPierre Jolivet                                              NULL,
1220*8bb0f5c6SPierre Jolivet                                              /*129*/ MatFDColoringSetUp_MPIXAIJ,
1221f4259b30SLisandro Dalcin                                              NULL,
1222f4259b30SLisandro Dalcin                                              NULL,
1223f4259b30SLisandro Dalcin                                              NULL,
1224f4259b30SLisandro Dalcin                                              NULL,
1225f4259b30SLisandro Dalcin                                              /*134*/ NULL,
1226f4259b30SLisandro Dalcin                                              NULL,
1227f4259b30SLisandro Dalcin                                              NULL,
1228f4259b30SLisandro Dalcin                                              NULL,
1229f4259b30SLisandro Dalcin                                              NULL,
1230f4259b30SLisandro Dalcin                                              /*139*/ NULL,
1231f4259b30SLisandro Dalcin                                              NULL,
1232f4259b30SLisandro Dalcin                                              NULL,
1233dec0b466SHong Zhang                                              NULL};
1234d4002b98SHong Zhang 
1235d4002b98SHong Zhang /*@C
123611a5261eSBarry Smith   MatMPISELLSetPreallocation - Preallocates memory for a `MATMPISELL` sparse parallel matrix in sell format.
1237d4002b98SHong Zhang   For good matrix assembly performance the user should preallocate the matrix storage by
123867be906fSBarry Smith   setting the parameters `d_nz` (or `d_nnz`) and `o_nz` (or `o_nnz`).
1239d4002b98SHong Zhang 
1240d083f849SBarry Smith   Collective
1241d4002b98SHong Zhang 
1242d4002b98SHong Zhang   Input Parameters:
1243d4002b98SHong Zhang + B     - the matrix
1244d4002b98SHong Zhang . d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1245d4002b98SHong Zhang            (same value is used for all local rows)
1246d4002b98SHong Zhang . d_nnz - array containing the number of nonzeros in the various rows of the
1247d4002b98SHong Zhang            DIAGONAL portion of the local submatrix (possibly different for each row)
124867be906fSBarry Smith            or NULL (`PETSC_NULL_INTEGER` in Fortran), if `d_nz` is used to specify the nonzero structure.
1249d4002b98SHong Zhang            The size of this array is equal to the number of local rows, i.e 'm'.
1250d4002b98SHong Zhang            For matrices that will be factored, you must leave room for (and set)
1251d4002b98SHong Zhang            the diagonal entry even if it is zero.
1252d4002b98SHong Zhang . o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1253d4002b98SHong Zhang            submatrix (same value is used for all local rows).
1254d4002b98SHong Zhang - o_nnz - array containing the number of nonzeros in the various rows of the
1255d4002b98SHong Zhang            OFF-DIAGONAL portion of the local submatrix (possibly different for
125667be906fSBarry Smith            each row) or NULL (`PETSC_NULL_INTEGER` in Fortran), if `o_nz` is used to specify the nonzero
1257d4002b98SHong Zhang            structure. The size of this array is equal to the number
1258d4002b98SHong Zhang            of local rows, i.e 'm'.
1259d4002b98SHong Zhang 
1260d4002b98SHong Zhang   Example usage:
1261d4002b98SHong Zhang   Consider the following 8x8 matrix with 34 non-zero values, that is
1262d4002b98SHong Zhang   assembled across 3 processors. Lets assume that proc0 owns 3 rows,
1263d4002b98SHong Zhang   proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
126467be906fSBarry Smith   as follows
1265d4002b98SHong Zhang 
1266d4002b98SHong Zhang .vb
1267d4002b98SHong Zhang             1  2  0  |  0  3  0  |  0  4
1268d4002b98SHong Zhang     Proc0   0  5  6  |  7  0  0  |  8  0
1269d4002b98SHong Zhang             9  0 10  | 11  0  0  | 12  0
1270d4002b98SHong Zhang     -------------------------------------
1271d4002b98SHong Zhang            13  0 14  | 15 16 17  |  0  0
1272d4002b98SHong Zhang     Proc1   0 18  0  | 19 20 21  |  0  0
1273d4002b98SHong Zhang             0  0  0  | 22 23  0  | 24  0
1274d4002b98SHong Zhang     -------------------------------------
1275d4002b98SHong Zhang     Proc2  25 26 27  |  0  0 28  | 29  0
1276d4002b98SHong Zhang            30  0  0  | 31 32 33  |  0 34
1277d4002b98SHong Zhang .ve
1278d4002b98SHong Zhang 
127927430b45SBarry Smith   This can be represented as a collection of submatrices as
1280d4002b98SHong Zhang 
1281d4002b98SHong Zhang .vb
1282d4002b98SHong Zhang       A B C
1283d4002b98SHong Zhang       D E F
1284d4002b98SHong Zhang       G H I
1285d4002b98SHong Zhang .ve
1286d4002b98SHong Zhang 
1287d4002b98SHong Zhang   Where the submatrices A,B,C are owned by proc0, D,E,F are
1288d4002b98SHong Zhang   owned by proc1, G,H,I are owned by proc2.
1289d4002b98SHong Zhang 
1290d4002b98SHong Zhang   The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1291d4002b98SHong Zhang   The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1292d4002b98SHong Zhang   The 'M','N' parameters are 8,8, and have the same values on all procs.
1293d4002b98SHong Zhang 
1294d4002b98SHong Zhang   The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
1295d4002b98SHong Zhang   submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
1296d4002b98SHong Zhang   corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
1297d4002b98SHong Zhang   Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
129827430b45SBarry Smith   part as `MATSEQSELL` matrices. For example, proc1 will store [E] as a `MATSEQSELL`
12995163b989SNuno Nobre   matrix, and [DF] as another SeqSELL matrix.
1300d4002b98SHong Zhang 
130167be906fSBarry Smith   When `d_nz`, `o_nz` parameters are specified, `d_nz` storage elements are
13025163b989SNuno Nobre   allocated for every row of the local DIAGONAL submatrix, and o_nz
13035163b989SNuno Nobre   storage locations are allocated for every row of the OFF-DIAGONAL submatrix.
13045163b989SNuno Nobre   One way to choose `d_nz` and `o_nz` is to use the maximum number of nonzeros over
13055163b989SNuno Nobre   the local rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
130627430b45SBarry Smith   In this case, the values of d_nz,o_nz are
1307d4002b98SHong Zhang .vb
130827430b45SBarry Smith      proc0  dnz = 2, o_nz = 2
130927430b45SBarry Smith      proc1  dnz = 3, o_nz = 2
131027430b45SBarry Smith      proc2  dnz = 1, o_nz = 4
1311d4002b98SHong Zhang .ve
1312d4002b98SHong Zhang   We are allocating m*(d_nz+o_nz) storage locations for every proc. This
1313d4002b98SHong Zhang   translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
1314d4002b98SHong Zhang   for proc3. i.e we are using 12+15+10=37 storage locations to store
1315d4002b98SHong Zhang   34 values.
1316d4002b98SHong Zhang 
131767be906fSBarry Smith   When `d_nnz`, `o_nnz` parameters are specified, the storage is specified
1318a5b23f4aSJose E. Roman   for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
131927430b45SBarry Smith   In the above case the values for d_nnz,o_nnz are
1320d4002b98SHong Zhang .vb
132127430b45SBarry Smith      proc0 d_nnz = [2,2,2] and o_nnz = [2,2,2]
132227430b45SBarry Smith      proc1 d_nnz = [3,3,2] and o_nnz = [2,1,1]
132327430b45SBarry Smith      proc2 d_nnz = [1,1]   and o_nnz = [4,4]
1324d4002b98SHong Zhang .ve
1325d4002b98SHong Zhang   Here the space allocated is according to nz (or maximum values in the nnz
1326d4002b98SHong Zhang   if nnz is provided) for DIAGONAL and OFF-DIAGONAL submatrices, i.e (2+2+3+2)*3+(1+4)*2=37
1327d4002b98SHong Zhang 
1328d4002b98SHong Zhang   Level: intermediate
1329d4002b98SHong Zhang 
133067be906fSBarry Smith   Notes:
133167be906fSBarry Smith   If the *_nnz parameter is given then the *_nz parameter is ignored
133267be906fSBarry Smith 
133367be906fSBarry Smith   The stored row and column indices begin with zero.
133467be906fSBarry Smith 
133567be906fSBarry Smith   The parallel matrix is partitioned such that the first m0 rows belong to
133667be906fSBarry Smith   process 0, the next m1 rows belong to process 1, the next m2 rows belong
133767be906fSBarry Smith   to process 2 etc.. where m0,m1,m2... are the input parameter 'm'.
133867be906fSBarry Smith 
133967be906fSBarry Smith   The DIAGONAL portion of the local submatrix of a processor can be defined
134067be906fSBarry Smith   as the submatrix which is obtained by extraction the part corresponding to
134167be906fSBarry Smith   the rows r1-r2 and columns c1-c2 of the global matrix, where r1 is the
134267be906fSBarry Smith   first row that belongs to the processor, r2 is the last row belonging to
134367be906fSBarry Smith   the this processor, and c1-c2 is range of indices of the local part of a
134467be906fSBarry Smith   vector suitable for applying the matrix to.  This is an mxn matrix.  In the
134567be906fSBarry Smith   common case of a square matrix, the row and column ranges are the same and
134667be906fSBarry Smith   the DIAGONAL part is also square. The remaining portion of the local
134767be906fSBarry Smith   submatrix (mxN) constitute the OFF-DIAGONAL portion.
134867be906fSBarry Smith 
134967be906fSBarry Smith   If `o_nnz`, `d_nnz` are specified, then `o_nz`, and `d_nz` are ignored.
135067be906fSBarry Smith 
135167be906fSBarry Smith   You can call `MatGetInfo()` to get information on how effective the preallocation was;
135267be906fSBarry Smith   for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
135367be906fSBarry Smith   You can also run with the option -info and look for messages with the string
135467be906fSBarry Smith   malloc in them to see if additional memory allocation was needed.
135567be906fSBarry Smith 
135694764886SPierre Jolivet .seealso: `Mat`, `MatCreate()`, `MatCreateSeqSELL()`, `MatSetValues()`, `MatCreateSELL()`,
135711a5261eSBarry Smith           `MATMPISELL`, `MatGetInfo()`, `PetscSplitOwnership()`, `MATSELL`
1358d4002b98SHong Zhang @*/
1359d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPISELLSetPreallocation(Mat B, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[])
1360d71ae5a4SJacob Faibussowitsch {
1361d4002b98SHong Zhang   PetscFunctionBegin;
1362d4002b98SHong Zhang   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
1363d4002b98SHong Zhang   PetscValidType(B, 1);
1364cac4c232SBarry Smith   PetscTryMethod(B, "MatMPISELLSetPreallocation_C", (Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]), (B, d_nz, d_nnz, o_nz, o_nnz));
13653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1366d4002b98SHong Zhang }
1367d4002b98SHong Zhang 
1368ed73aabaSBarry Smith /*MC
1369ed73aabaSBarry Smith    MATMPISELL - MATMPISELL = "mpisell" - A matrix type to be used for MPI sparse matrices,
1370ed73aabaSBarry Smith    based on the sliced Ellpack format
1371ed73aabaSBarry Smith 
137227430b45SBarry Smith    Options Database Key:
137320f4b53cSBarry Smith . -mat_type sell - sets the matrix type to `MATSELL` during a call to `MatSetFromOptions()`
1374ed73aabaSBarry Smith 
1375ed73aabaSBarry Smith    Level: beginner
1376ed73aabaSBarry Smith 
137794764886SPierre Jolivet .seealso: `Mat`, `MatCreateSELL()`, `MATSEQSELL`, `MATSELL`, `MATSEQAIJ`, `MATAIJ`, `MATMPIAIJ`
1378ed73aabaSBarry Smith M*/
1379ed73aabaSBarry Smith 
1380d4002b98SHong Zhang /*@C
138111a5261eSBarry Smith   MatCreateSELL - Creates a sparse parallel matrix in `MATSELL` format.
1382d4002b98SHong Zhang 
1383d083f849SBarry Smith   Collective
1384d4002b98SHong Zhang 
1385d4002b98SHong Zhang   Input Parameters:
1386d4002b98SHong Zhang + comm      - MPI communicator
138711a5261eSBarry Smith . m         - number of local rows (or `PETSC_DECIDE` to have calculated if M is given)
1388d4002b98SHong Zhang               This value should be the same as the local size used in creating the
1389d4002b98SHong Zhang               y vector for the matrix-vector product y = Ax.
1390d4002b98SHong Zhang . n         - This value should be the same as the local size used in creating the
139120f4b53cSBarry Smith               x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
139220f4b53cSBarry Smith               calculated if `N` is given) For square matrices n is almost always `m`.
139320f4b53cSBarry Smith . M         - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
139420f4b53cSBarry Smith . N         - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
1395d4002b98SHong Zhang . d_rlenmax - max number of nonzeros per row in DIAGONAL portion of local submatrix
1396d4002b98SHong Zhang              (same value is used for all local rows)
1397d4002b98SHong Zhang . d_rlen    - array containing the number of nonzeros in the various rows of the
1398d4002b98SHong Zhang               DIAGONAL portion of the local submatrix (possibly different for each row)
139920f4b53cSBarry Smith               or `NULL`, if d_rlenmax is used to specify the nonzero structure.
140020f4b53cSBarry Smith               The size of this array is equal to the number of local rows, i.e `m`.
1401d4002b98SHong Zhang . o_rlenmax - max number of nonzeros per row in the OFF-DIAGONAL portion of local
1402d4002b98SHong Zhang               submatrix (same value is used for all local rows).
1403d4002b98SHong Zhang - o_rlen    - array containing the number of nonzeros in the various rows of the
1404d4002b98SHong Zhang               OFF-DIAGONAL portion of the local submatrix (possibly different for
140520f4b53cSBarry Smith               each row) or `NULL`, if `o_rlenmax` is used to specify the nonzero
1406d4002b98SHong Zhang               structure. The size of this array is equal to the number
140720f4b53cSBarry Smith               of local rows, i.e `m`.
1408d4002b98SHong Zhang 
1409d4002b98SHong Zhang   Output Parameter:
1410d4002b98SHong Zhang . A - the matrix
1411d4002b98SHong Zhang 
141227430b45SBarry Smith   Options Database Key:
1413fe59aa6dSJacob Faibussowitsch . -mat_sell_oneindex - Internally use indexing starting at 1
141427430b45SBarry Smith         rather than 0.  When calling `MatSetValues()`,
141527430b45SBarry Smith         the user still MUST index entries starting at 0!
141627430b45SBarry Smith 
141727430b45SBarry Smith   Example:
141827430b45SBarry Smith   Consider the following 8x8 matrix with 34 non-zero values, that is
141927430b45SBarry Smith   assembled across 3 processors. Lets assume that proc0 owns 3 rows,
142027430b45SBarry Smith   proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
142127430b45SBarry Smith   as follows
142227430b45SBarry Smith 
142327430b45SBarry Smith .vb
142427430b45SBarry Smith             1  2  0  |  0  3  0  |  0  4
142527430b45SBarry Smith     Proc0   0  5  6  |  7  0  0  |  8  0
142627430b45SBarry Smith             9  0 10  | 11  0  0  | 12  0
142727430b45SBarry Smith     -------------------------------------
142827430b45SBarry Smith            13  0 14  | 15 16 17  |  0  0
142927430b45SBarry Smith     Proc1   0 18  0  | 19 20 21  |  0  0
143027430b45SBarry Smith             0  0  0  | 22 23  0  | 24  0
143127430b45SBarry Smith     -------------------------------------
143227430b45SBarry Smith     Proc2  25 26 27  |  0  0 28  | 29  0
143327430b45SBarry Smith            30  0  0  | 31 32 33  |  0 34
143427430b45SBarry Smith .ve
143527430b45SBarry Smith 
143620f4b53cSBarry Smith   This can be represented as a collection of submatrices as
143727430b45SBarry Smith .vb
143827430b45SBarry Smith       A B C
143927430b45SBarry Smith       D E F
144027430b45SBarry Smith       G H I
144127430b45SBarry Smith .ve
144227430b45SBarry Smith 
144327430b45SBarry Smith   Where the submatrices A,B,C are owned by proc0, D,E,F are
144427430b45SBarry Smith   owned by proc1, G,H,I are owned by proc2.
144527430b45SBarry Smith 
144627430b45SBarry Smith   The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
144727430b45SBarry Smith   The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
144827430b45SBarry Smith   The 'M','N' parameters are 8,8, and have the same values on all procs.
144927430b45SBarry Smith 
145027430b45SBarry Smith   The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
145127430b45SBarry Smith   submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
145227430b45SBarry Smith   corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
145327430b45SBarry Smith   Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
145427430b45SBarry Smith   part as `MATSEQSELL` matrices. For example, proc1 will store [E] as a `MATSEQSELL`
14555163b989SNuno Nobre   matrix, and [DF] as another `MATSEQSELL` matrix.
145627430b45SBarry Smith 
145727430b45SBarry Smith   When d_rlenmax, o_rlenmax parameters are specified, d_rlenmax storage elements are
14585163b989SNuno Nobre   allocated for every row of the local DIAGONAL submatrix, and o_rlenmax
14595163b989SNuno Nobre   storage locations are allocated for every row of the OFF-DIAGONAL submatrix.
14605163b989SNuno Nobre   One way to choose `d_rlenmax` and `o_rlenmax` is to use the maximum number of nonzeros over
14615163b989SNuno Nobre   the local rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
146220f4b53cSBarry Smith   In this case, the values of d_rlenmax,o_rlenmax are
146327430b45SBarry Smith .vb
146420f4b53cSBarry Smith      proc0 - d_rlenmax = 2, o_rlenmax = 2
146520f4b53cSBarry Smith      proc1 - d_rlenmax = 3, o_rlenmax = 2
146620f4b53cSBarry Smith      proc2 - d_rlenmax = 1, o_rlenmax = 4
146727430b45SBarry Smith .ve
146827430b45SBarry Smith   We are allocating m*(d_rlenmax+o_rlenmax) storage locations for every proc. This
146927430b45SBarry Smith   translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
147027430b45SBarry Smith   for proc3. i.e we are using 12+15+10=37 storage locations to store
147127430b45SBarry Smith   34 values.
147227430b45SBarry Smith 
147320f4b53cSBarry Smith   When `d_rlen`, `o_rlen` parameters are specified, the storage is specified
147427430b45SBarry Smith   for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
147520f4b53cSBarry Smith   In the above case the values for `d_nnz`, `o_nnz` are
147627430b45SBarry Smith .vb
147720f4b53cSBarry Smith      proc0 - d_nnz = [2,2,2] and o_nnz = [2,2,2]
147820f4b53cSBarry Smith      proc1 - d_nnz = [3,3,2] and o_nnz = [2,1,1]
147920f4b53cSBarry Smith      proc2 - d_nnz = [1,1]   and o_nnz = [4,4]
148027430b45SBarry Smith .ve
148127430b45SBarry Smith   Here the space allocated is still 37 though there are 34 nonzeros because
148227430b45SBarry Smith   the allocation is always done according to rlenmax.
148327430b45SBarry Smith 
148427430b45SBarry Smith   Level: intermediate
148527430b45SBarry Smith 
148627430b45SBarry Smith   Notes:
148711a5261eSBarry Smith   It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`,
1488f6f02116SRichard Tran Mills   MatXXXXSetPreallocation() paradigm instead of this routine directly.
148911a5261eSBarry Smith   [MatXXXXSetPreallocation() is, for example, `MatSeqSELLSetPreallocation()`]
1490d4002b98SHong Zhang 
1491d4002b98SHong Zhang   If the *_rlen parameter is given then the *_rlenmax parameter is ignored
1492d4002b98SHong Zhang 
149320f4b53cSBarry Smith   `m`, `n`, `M`, `N` parameters specify the size of the matrix, and its partitioning across
149420f4b53cSBarry Smith   processors, while `d_rlenmax`, `d_rlen`, `o_rlenmax` , `o_rlen` parameters specify the approximate
1495d4002b98SHong Zhang   storage requirements for this matrix.
1496d4002b98SHong Zhang 
149711a5261eSBarry Smith   If `PETSC_DECIDE` or  `PETSC_DETERMINE` is used for a particular argument on one
1498d4002b98SHong Zhang   processor than it must be used on all processors that share the object for
1499d4002b98SHong Zhang   that argument.
1500d4002b98SHong Zhang 
1501d4002b98SHong Zhang   The user MUST specify either the local or global matrix dimensions
1502d4002b98SHong Zhang   (possibly both).
1503d4002b98SHong Zhang 
1504d4002b98SHong Zhang   The parallel matrix is partitioned across processors such that the
1505d4002b98SHong Zhang   first m0 rows belong to process 0, the next m1 rows belong to
1506d4002b98SHong Zhang   process 1, the next m2 rows belong to process 2 etc.. where
1507d4002b98SHong Zhang   m0,m1,m2,.. are the input parameter 'm'. i.e each processor stores
150820f4b53cSBarry Smith   values corresponding to [`m` x `N`] submatrix.
1509d4002b98SHong Zhang 
1510d4002b98SHong Zhang   The columns are logically partitioned with the n0 columns belonging
1511d4002b98SHong Zhang   to 0th partition, the next n1 columns belonging to the next
151220f4b53cSBarry Smith   partition etc.. where n0,n1,n2... are the input parameter `n`.
1513d4002b98SHong Zhang 
1514d4002b98SHong Zhang   The DIAGONAL portion of the local submatrix on any given processor
151520f4b53cSBarry Smith   is the submatrix corresponding to the rows and columns `m`, `n`
1516d4002b98SHong Zhang   corresponding to the given processor. i.e diagonal matrix on
1517d4002b98SHong Zhang   process 0 is [m0 x n0], diagonal matrix on process 1 is [m1 x n1]
1518d4002b98SHong Zhang   etc. The remaining portion of the local submatrix [m x (N-n)]
1519d4002b98SHong Zhang   constitute the OFF-DIAGONAL portion. The example below better
1520d4002b98SHong Zhang   illustrates this concept.
1521d4002b98SHong Zhang 
1522d4002b98SHong Zhang   For a square global matrix we define each processor's diagonal portion
1523d4002b98SHong Zhang   to be its local rows and the corresponding columns (a square submatrix);
1524d4002b98SHong Zhang   each processor's off-diagonal portion encompasses the remainder of the
1525d4002b98SHong Zhang   local matrix (a rectangular submatrix).
1526d4002b98SHong Zhang 
152720f4b53cSBarry Smith   If `o_rlen`, `d_rlen` are specified, then `o_rlenmax`, and `d_rlenmax` are ignored.
1528d4002b98SHong Zhang 
1529d4002b98SHong Zhang   When calling this routine with a single process communicator, a matrix of
153011a5261eSBarry Smith   type `MATSEQSELL` is returned.  If a matrix of type `MATMPISELL` is desired for this
153127430b45SBarry Smith   type of communicator, use the construction mechanism
1532d4002b98SHong Zhang .vb
153327430b45SBarry Smith    MatCreate(...,&A);
153427430b45SBarry Smith    MatSetType(A,MATMPISELL);
153527430b45SBarry Smith    MatSetSizes(A, m,n,M,N);
153627430b45SBarry Smith    MatMPISELLSetPreallocation(A,...);
1537d4002b98SHong Zhang .ve
1538d4002b98SHong Zhang 
153994764886SPierre Jolivet .seealso: `Mat`, `MATSELL`, `MatCreate()`, `MatCreateSeqSELL()`, `MatSetValues()`, `MatMPISELLSetPreallocation()`, `MATMPISELL`
1540d4002b98SHong Zhang @*/
1541d71ae5a4SJacob 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)
1542d71ae5a4SJacob Faibussowitsch {
1543d4002b98SHong Zhang   PetscMPIInt size;
1544d4002b98SHong Zhang 
1545d4002b98SHong Zhang   PetscFunctionBegin;
15469566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
15479566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, n, M, N));
15489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1549d4002b98SHong Zhang   if (size > 1) {
15509566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A, MATMPISELL));
15519566063dSJacob Faibussowitsch     PetscCall(MatMPISELLSetPreallocation(*A, d_rlenmax, d_rlen, o_rlenmax, o_rlen));
1552d4002b98SHong Zhang   } else {
15539566063dSJacob Faibussowitsch     PetscCall(MatSetType(*A, MATSEQSELL));
15549566063dSJacob Faibussowitsch     PetscCall(MatSeqSELLSetPreallocation(*A, d_rlenmax, d_rlen));
1555d4002b98SHong Zhang   }
15563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1557d4002b98SHong Zhang }
1558d4002b98SHong Zhang 
1559fa078d78SJacob Faibussowitsch /*@C
1560fa078d78SJacob Faibussowitsch   MatMPISELLGetSeqSELL - Returns the local pieces of this distributed matrix
1561fa078d78SJacob Faibussowitsch 
1562fa078d78SJacob Faibussowitsch   Not Collective
1563fa078d78SJacob Faibussowitsch 
1564fa078d78SJacob Faibussowitsch   Input Parameter:
1565fa078d78SJacob Faibussowitsch . A - the `MATMPISELL` matrix
1566fa078d78SJacob Faibussowitsch 
1567fa078d78SJacob Faibussowitsch   Output Parameters:
1568fa078d78SJacob Faibussowitsch + Ad     - The diagonal portion of `A`
15694cf0e950SBarry Smith . Ao     - The off-diagonal portion of `A`
1570fa078d78SJacob Faibussowitsch - colmap - An array mapping local column numbers of `Ao` to global column numbers of the parallel matrix
1571fa078d78SJacob Faibussowitsch 
1572fa078d78SJacob Faibussowitsch   Level: advanced
1573fa078d78SJacob Faibussowitsch 
1574fa078d78SJacob Faibussowitsch .seealso: `Mat`, `MATSEQSELL`, `MATMPISELL`
1575fa078d78SJacob Faibussowitsch @*/
1576fa078d78SJacob Faibussowitsch PetscErrorCode MatMPISELLGetSeqSELL(Mat A, Mat *Ad, Mat *Ao, const PetscInt *colmap[])
1577fa078d78SJacob Faibussowitsch {
1578fa078d78SJacob Faibussowitsch   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1579fa078d78SJacob Faibussowitsch   PetscBool    flg;
1580fa078d78SJacob Faibussowitsch 
1581fa078d78SJacob Faibussowitsch   PetscFunctionBegin;
1582fa078d78SJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISELL, &flg));
1583fa078d78SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "This function requires a MATMPISELL matrix as input");
1584fa078d78SJacob Faibussowitsch   if (Ad) *Ad = a->A;
1585fa078d78SJacob Faibussowitsch   if (Ao) *Ao = a->B;
1586fa078d78SJacob Faibussowitsch   if (colmap) *colmap = a->garray;
1587fa078d78SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1588fa078d78SJacob Faibussowitsch }
1589fa078d78SJacob Faibussowitsch 
1590fa078d78SJacob Faibussowitsch /*@C
1591fa078d78SJacob Faibussowitsch   MatMPISELLGetLocalMatCondensed - Creates a `MATSEQSELL` matrix from an `MATMPISELL` matrix by
1592fa078d78SJacob Faibussowitsch   taking all its local rows and NON-ZERO columns
1593fa078d78SJacob Faibussowitsch 
1594fa078d78SJacob Faibussowitsch   Not Collective
1595fa078d78SJacob Faibussowitsch 
1596fa078d78SJacob Faibussowitsch   Input Parameters:
1597fa078d78SJacob Faibussowitsch + A     - the matrix
1598fa078d78SJacob Faibussowitsch . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
1599fa078d78SJacob Faibussowitsch . row   - index sets of rows to extract (or `NULL`)
1600fa078d78SJacob Faibussowitsch - col   - index sets of columns to extract (or `NULL`)
1601fa078d78SJacob Faibussowitsch 
1602fa078d78SJacob Faibussowitsch   Output Parameter:
1603fa078d78SJacob Faibussowitsch . A_loc - the local sequential matrix generated
1604fa078d78SJacob Faibussowitsch 
1605fa078d78SJacob Faibussowitsch   Level: advanced
1606fa078d78SJacob Faibussowitsch 
1607fa078d78SJacob Faibussowitsch .seealso: `Mat`, `MATSEQSELL`, `MATMPISELL`, `MatGetOwnershipRange()`, `MatMPISELLGetLocalMat()`
1608fa078d78SJacob Faibussowitsch @*/
1609fa078d78SJacob Faibussowitsch PetscErrorCode MatMPISELLGetLocalMatCondensed(Mat A, MatReuse scall, IS *row, IS *col, Mat *A_loc)
1610fa078d78SJacob Faibussowitsch {
1611fa078d78SJacob Faibussowitsch   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1612fa078d78SJacob Faibussowitsch   PetscInt     i, start, end, ncols, nzA, nzB, *cmap, imark, *idx;
1613fa078d78SJacob Faibussowitsch   IS           isrowa, iscola;
1614fa078d78SJacob Faibussowitsch   Mat         *aloc;
1615fa078d78SJacob Faibussowitsch   PetscBool    match;
1616fa078d78SJacob Faibussowitsch 
1617fa078d78SJacob Faibussowitsch   PetscFunctionBegin;
1618fa078d78SJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISELL, &match));
1619fa078d78SJacob Faibussowitsch   PetscCheck(match, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Requires MATMPISELL matrix as input");
1620fa078d78SJacob Faibussowitsch   PetscCall(PetscLogEventBegin(MAT_Getlocalmatcondensed, A, 0, 0, 0));
1621fa078d78SJacob Faibussowitsch   if (!row) {
1622fa078d78SJacob Faibussowitsch     start = A->rmap->rstart;
1623fa078d78SJacob Faibussowitsch     end   = A->rmap->rend;
1624fa078d78SJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, end - start, start, 1, &isrowa));
1625fa078d78SJacob Faibussowitsch   } else {
1626fa078d78SJacob Faibussowitsch     isrowa = *row;
1627fa078d78SJacob Faibussowitsch   }
1628fa078d78SJacob Faibussowitsch   if (!col) {
1629fa078d78SJacob Faibussowitsch     start = A->cmap->rstart;
1630fa078d78SJacob Faibussowitsch     cmap  = a->garray;
1631fa078d78SJacob Faibussowitsch     nzA   = a->A->cmap->n;
1632fa078d78SJacob Faibussowitsch     nzB   = a->B->cmap->n;
1633fa078d78SJacob Faibussowitsch     PetscCall(PetscMalloc1(nzA + nzB, &idx));
1634fa078d78SJacob Faibussowitsch     ncols = 0;
1635fa078d78SJacob Faibussowitsch     for (i = 0; i < nzB; i++) {
1636fa078d78SJacob Faibussowitsch       if (cmap[i] < start) idx[ncols++] = cmap[i];
1637fa078d78SJacob Faibussowitsch       else break;
1638fa078d78SJacob Faibussowitsch     }
1639fa078d78SJacob Faibussowitsch     imark = i;
1640fa078d78SJacob Faibussowitsch     for (i = 0; i < nzA; i++) idx[ncols++] = start + i;
1641fa078d78SJacob Faibussowitsch     for (i = imark; i < nzB; i++) idx[ncols++] = cmap[i];
1642fa078d78SJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncols, idx, PETSC_OWN_POINTER, &iscola));
1643fa078d78SJacob Faibussowitsch   } else {
1644fa078d78SJacob Faibussowitsch     iscola = *col;
1645fa078d78SJacob Faibussowitsch   }
1646fa078d78SJacob Faibussowitsch   if (scall != MAT_INITIAL_MATRIX) {
1647fa078d78SJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &aloc));
1648fa078d78SJacob Faibussowitsch     aloc[0] = *A_loc;
1649fa078d78SJacob Faibussowitsch   }
1650fa078d78SJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, 1, &isrowa, &iscola, scall, &aloc));
1651fa078d78SJacob Faibussowitsch   *A_loc = aloc[0];
1652fa078d78SJacob Faibussowitsch   PetscCall(PetscFree(aloc));
1653fa078d78SJacob Faibussowitsch   if (!row) PetscCall(ISDestroy(&isrowa));
1654fa078d78SJacob Faibussowitsch   if (!col) PetscCall(ISDestroy(&iscola));
1655fa078d78SJacob Faibussowitsch   PetscCall(PetscLogEventEnd(MAT_Getlocalmatcondensed, A, 0, 0, 0));
1656fa078d78SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1657fa078d78SJacob Faibussowitsch }
1658fa078d78SJacob Faibussowitsch 
1659d4002b98SHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>
1660d4002b98SHong Zhang 
1661d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_MPISELL_MPIAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1662d71ae5a4SJacob Faibussowitsch {
1663d4002b98SHong Zhang   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1664d4002b98SHong Zhang   Mat          B;
1665d4002b98SHong Zhang   Mat_MPIAIJ  *b;
1666d4002b98SHong Zhang 
1667d4002b98SHong Zhang   PetscFunctionBegin;
166828b400f6SJacob Faibussowitsch   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Matrix must be assembled");
1669d4002b98SHong Zhang 
167094a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
167194a8b381SRichard Tran Mills     B = *newmat;
167294a8b381SRichard Tran Mills   } else {
16739566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
16749566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATMPIAIJ));
16759566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
16769566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
16779566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL));
16789566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL));
167994a8b381SRichard Tran Mills   }
1680d4002b98SHong Zhang   b = (Mat_MPIAIJ *)B->data;
168194a8b381SRichard Tran Mills 
168294a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
16839566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_REUSE_MATRIX, &b->A));
16849566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_REUSE_MATRIX, &b->B));
168594a8b381SRichard Tran Mills   } else {
16869566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->A));
16879566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->B));
16889566063dSJacob Faibussowitsch     PetscCall(MatDisAssemble_MPISELL(A));
16899566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->A));
16909566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->B));
16919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
16929566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
16939566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
16949566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
169594a8b381SRichard Tran Mills   }
1696d4002b98SHong Zhang 
1697d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
16989566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
1699d4002b98SHong Zhang   } else {
1700d4002b98SHong Zhang     *newmat = B;
1701d4002b98SHong Zhang   }
17023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1703d4002b98SHong Zhang }
1704d4002b98SHong Zhang 
1705d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvert_MPIAIJ_MPISELL(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1706d71ae5a4SJacob Faibussowitsch {
1707d4002b98SHong Zhang   Mat_MPIAIJ  *a = (Mat_MPIAIJ *)A->data;
1708d4002b98SHong Zhang   Mat          B;
1709d4002b98SHong Zhang   Mat_MPISELL *b;
1710d4002b98SHong Zhang 
1711d4002b98SHong Zhang   PetscFunctionBegin;
171228b400f6SJacob Faibussowitsch   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Matrix must be assembled");
1713d4002b98SHong Zhang 
171494a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
171594a8b381SRichard Tran Mills     B = *newmat;
171694a8b381SRichard Tran Mills   } else {
17172d1451d4SHong Zhang     Mat_SeqAIJ *Aa = (Mat_SeqAIJ *)a->A->data, *Ba = (Mat_SeqAIJ *)a->B->data;
17182d1451d4SHong 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;
17192d1451d4SHong Zhang     PetscInt   *d_nnz, *o_nnz;
17202d1451d4SHong Zhang     PetscCall(PetscMalloc2(lm, &d_nnz, lm, &o_nnz));
17212d1451d4SHong Zhang     for (i = 0; i < lm; i++) {
17222d1451d4SHong Zhang       d_nnz[i] = Aa->i[i + 1] - Aa->i[i];
17232d1451d4SHong Zhang       o_nnz[i] = Ba->i[i + 1] - Ba->i[i];
17242d1451d4SHong Zhang       if (d_nnz[i] > d_nz) d_nz = d_nnz[i];
17252d1451d4SHong Zhang       if (o_nnz[i] > o_nz) o_nz = o_nnz[i];
17262d1451d4SHong Zhang     }
17279566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
17289566063dSJacob Faibussowitsch     PetscCall(MatSetType(B, MATMPISELL));
17292d1451d4SHong Zhang     PetscCall(MatSetSizes(B, lm, ln, m, n));
17309566063dSJacob Faibussowitsch     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
17312d1451d4SHong Zhang     PetscCall(MatSeqSELLSetPreallocation(B, d_nz, d_nnz));
17322d1451d4SHong Zhang     PetscCall(MatMPISELLSetPreallocation(B, d_nz, d_nnz, o_nz, o_nnz));
17332d1451d4SHong Zhang     PetscCall(PetscFree2(d_nnz, o_nnz));
173494a8b381SRichard Tran Mills   }
1735d4002b98SHong Zhang   b = (Mat_MPISELL *)B->data;
173694a8b381SRichard Tran Mills 
173794a8b381SRichard Tran Mills   if (reuse == MAT_REUSE_MATRIX) {
17389566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_REUSE_MATRIX, &b->A));
17399566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_REUSE_MATRIX, &b->B));
174094a8b381SRichard Tran Mills   } else {
17419566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->A));
17429566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&b->B));
17439566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_INITIAL_MATRIX, &b->A));
17449566063dSJacob Faibussowitsch     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_INITIAL_MATRIX, &b->B));
17459566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
17469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
17472d1451d4SHong Zhang     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
17482d1451d4SHong Zhang     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
174994a8b381SRichard Tran Mills   }
1750d4002b98SHong Zhang 
1751d4002b98SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
17529566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
1753d4002b98SHong Zhang   } else {
1754d4002b98SHong Zhang     *newmat = B;
1755d4002b98SHong Zhang   }
17563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1757d4002b98SHong Zhang }
1758d4002b98SHong Zhang 
1759d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSOR_MPISELL(Mat matin, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
1760d71ae5a4SJacob Faibussowitsch {
1761d4002b98SHong Zhang   Mat_MPISELL *mat = (Mat_MPISELL *)matin->data;
1762f4259b30SLisandro Dalcin   Vec          bb1 = NULL;
1763d4002b98SHong Zhang 
1764d4002b98SHong Zhang   PetscFunctionBegin;
1765d4002b98SHong Zhang   if (flag == SOR_APPLY_UPPER) {
17669566063dSJacob Faibussowitsch     PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
17673ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1768d4002b98SHong Zhang   }
1769d4002b98SHong Zhang 
177048a46eb9SPierre Jolivet   if (its > 1 || ~flag & SOR_ZERO_INITIAL_GUESS || flag & SOR_EISENSTAT) PetscCall(VecDuplicate(bb, &bb1));
1771d4002b98SHong Zhang 
1772d4002b98SHong Zhang   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
1773d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17749566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1775d4002b98SHong Zhang       its--;
1776d4002b98SHong Zhang     }
1777d4002b98SHong Zhang 
1778d4002b98SHong Zhang     while (its--) {
17799566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
17809566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1781d4002b98SHong Zhang 
1782d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
17839566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec, -1.0));
17849566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));
1785d4002b98SHong Zhang 
1786d4002b98SHong Zhang       /* local sweep */
17879566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_SYMMETRIC_SWEEP, fshift, lits, 1, xx));
1788d4002b98SHong Zhang     }
1789d4002b98SHong Zhang   } else if (flag & SOR_LOCAL_FORWARD_SWEEP) {
1790d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
17919566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1792d4002b98SHong Zhang       its--;
1793d4002b98SHong Zhang     }
1794d4002b98SHong Zhang     while (its--) {
17959566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
17969566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1797d4002b98SHong Zhang 
1798d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
17999566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec, -1.0));
18009566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));
1801d4002b98SHong Zhang 
1802d4002b98SHong Zhang       /* local sweep */
18039566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_FORWARD_SWEEP, fshift, lits, 1, xx));
1804d4002b98SHong Zhang     }
1805d4002b98SHong Zhang   } else if (flag & SOR_LOCAL_BACKWARD_SWEEP) {
1806d4002b98SHong Zhang     if (flag & SOR_ZERO_INITIAL_GUESS) {
18079566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1808d4002b98SHong Zhang       its--;
1809d4002b98SHong Zhang     }
1810d4002b98SHong Zhang     while (its--) {
18119566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
18129566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1813d4002b98SHong Zhang 
1814d4002b98SHong Zhang       /* update rhs: bb1 = bb - B*x */
18159566063dSJacob Faibussowitsch       PetscCall(VecScale(mat->lvec, -1.0));
18169566063dSJacob Faibussowitsch       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));
1817d4002b98SHong Zhang 
1818d4002b98SHong Zhang       /* local sweep */
18199566063dSJacob Faibussowitsch       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_BACKWARD_SWEEP, fshift, lits, 1, xx));
1820d4002b98SHong Zhang     }
1821d4002b98SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)matin), PETSC_ERR_SUP, "Parallel SOR not supported");
1822d4002b98SHong Zhang 
18239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bb1));
1824d4002b98SHong Zhang 
1825d4002b98SHong Zhang   matin->factorerrortype = mat->A->factorerrortype;
18263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1827d4002b98SHong Zhang }
1828d4002b98SHong Zhang 
1829b5917f1bSHong Zhang #if defined(PETSC_HAVE_CUDA)
1830b5917f1bSHong Zhang PETSC_INTERN PetscErrorCode MatConvert_MPISELL_MPISELLCUDA(Mat, MatType, MatReuse, Mat *);
1831b5917f1bSHong Zhang #endif
1832b5917f1bSHong Zhang 
1833d4002b98SHong Zhang /*MC
1834d4002b98SHong Zhang    MATMPISELL - MATMPISELL = "MPISELL" - A matrix type to be used for parallel sparse matrices.
1835d4002b98SHong Zhang 
1836d4002b98SHong Zhang    Options Database Keys:
183711a5261eSBarry Smith . -mat_type mpisell - sets the matrix type to `MATMPISELL` during a call to `MatSetFromOptions()`
1838d4002b98SHong Zhang 
1839d4002b98SHong Zhang   Level: beginner
1840d4002b98SHong Zhang 
184167be906fSBarry Smith .seealso: `Mat`, `MATSELL`, `MATSEQSELL` `MatCreateSELL()`
1842d4002b98SHong Zhang M*/
1843d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat B)
1844d71ae5a4SJacob Faibussowitsch {
1845d4002b98SHong Zhang   Mat_MPISELL *b;
1846d4002b98SHong Zhang   PetscMPIInt  size;
1847d4002b98SHong Zhang 
1848d4002b98SHong Zhang   PetscFunctionBegin;
18499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
18504dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
1851d4002b98SHong Zhang   B->data       = (void *)b;
1852aea10558SJacob Faibussowitsch   B->ops[0]     = MatOps_Values;
1853d4002b98SHong Zhang   B->assembled  = PETSC_FALSE;
1854d4002b98SHong Zhang   B->insertmode = NOT_SET_VALUES;
1855d4002b98SHong Zhang   b->size       = size;
18569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)B), &b->rank));
1857d4002b98SHong Zhang   /* build cache for off array entries formed */
18589566063dSJacob Faibussowitsch   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)B), 1, &B->stash));
1859d4002b98SHong Zhang 
1860d4002b98SHong Zhang   b->donotstash  = PETSC_FALSE;
1861f4259b30SLisandro Dalcin   b->colmap      = NULL;
1862f4259b30SLisandro Dalcin   b->garray      = NULL;
1863d4002b98SHong Zhang   b->roworiented = PETSC_TRUE;
1864d4002b98SHong Zhang 
1865d4002b98SHong Zhang   /* stuff used for matrix vector multiply */
1866d4002b98SHong Zhang   b->lvec  = NULL;
1867d4002b98SHong Zhang   b->Mvctx = NULL;
1868d4002b98SHong Zhang 
1869d4002b98SHong Zhang   /* stuff for MatGetRow() */
1870f4259b30SLisandro Dalcin   b->rowindices   = NULL;
1871f4259b30SLisandro Dalcin   b->rowvalues    = NULL;
1872d4002b98SHong Zhang   b->getrowactive = PETSC_FALSE;
1873d4002b98SHong Zhang 
18749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_MPISELL));
18759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_MPISELL));
18769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsTranspose_C", MatIsTranspose_MPISELL));
18779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPISELLSetPreallocation_C", MatMPISELLSetPreallocation_MPISELL));
18789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisell_mpiaij_C", MatConvert_MPISELL_MPIAIJ));
1879b5917f1bSHong Zhang #if defined(PETSC_HAVE_CUDA)
1880b5917f1bSHong Zhang   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisell_mpisellcuda_C", MatConvert_MPISELL_MPISELLCUDA));
1881b5917f1bSHong Zhang #endif
18829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatDiagonalScaleLocal_C", MatDiagonalScaleLocal_MPISELL));
18839566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPISELL));
18843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1885d4002b98SHong Zhang }
1886