Actual source code: mpisell.c

  1: #include <../src/mat/impls/aij/mpi/mpiaij.h>
  2: #include <../src/mat/impls/sell/mpi/mpisell.h>
  3: #include <petsc/private/vecimpl.h>
  4: #include <petsc/private/isimpl.h>
  5: #include <petscblaslapack.h>
  6: #include <petscsf.h>

  8: /*MC
  9:    MATSELL - MATSELL = "sell" - A matrix type to be used for sparse matrices.

 11:    This matrix type is identical to `MATSEQSELL` when constructed with a single process communicator,
 12:    and `MATMPISELL` otherwise.  As a result, for single process communicators,
 13:   `MatSeqSELLSetPreallocation()` is supported, and similarly `MatMPISELLSetPreallocation()` is supported
 14:   for communicators controlling multiple processes.  It is recommended that you call both of
 15:   the above preallocation routines for simplicity.

 17:    Options Database Keys:
 18: . -mat_type sell - sets the matrix type to `MATSELL` during a call to `MatSetFromOptions()`

 20:   Level: beginner

 22: .seealso: `Mat`, `MATAIJ`, `MATBAIJ`, `MATSBAIJ`, `MatCreateSELL()`, `MatCreateSeqSELL()`, `MATSEQSELL`, `MATMPISELL`
 23: M*/

 25: static PetscErrorCode MatDiagonalSet_MPISELL(Mat Y, Vec D, InsertMode is)
 26: {
 27:   Mat_MPISELL *sell = (Mat_MPISELL *)Y->data;

 29:   PetscFunctionBegin;
 30:   if (Y->assembled && Y->rmap->rstart == Y->cmap->rstart && Y->rmap->rend == Y->cmap->rend) {
 31:     PetscCall(MatDiagonalSet(sell->A, D, is));
 32:   } else {
 33:     PetscCall(MatDiagonalSet_Default(Y, D, is));
 34:   }
 35:   PetscFunctionReturn(PETSC_SUCCESS);
 36: }

 38: /*
 39:   Local utility routine that creates a mapping from the global column
 40: number to the local number in the off-diagonal part of the local
 41: storage of the matrix.  When PETSC_USE_CTABLE is used this is scalable at
 42: a slightly higher hash table cost; without it it is not scalable (each processor
 43: has an order N integer array but is fast to access.
 44: */
 45: PetscErrorCode MatCreateColmap_MPISELL_Private(Mat mat)
 46: {
 47:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
 48:   PetscInt     n    = sell->B->cmap->n, i;

 50:   PetscFunctionBegin;
 51:   PetscCheck(sell->garray, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MPISELL Matrix was assembled but is missing garray");
 52: #if defined(PETSC_USE_CTABLE)
 53:   PetscCall(PetscHMapICreateWithSize(n, &sell->colmap));
 54:   for (i = 0; i < n; i++) PetscCall(PetscHMapISet(sell->colmap, sell->garray[i] + 1, i + 1));
 55: #else
 56:   PetscCall(PetscCalloc1(mat->cmap->N + 1, &sell->colmap));
 57:   for (i = 0; i < n; i++) sell->colmap[sell->garray[i]] = i + 1;
 58: #endif
 59:   PetscFunctionReturn(PETSC_SUCCESS);
 60: }

 62: #define MatSetValues_SeqSELL_A_Private(row, col, value, addv, orow, ocol) \
 63:   { \
 64:     if (col <= lastcol1) low1 = 0; \
 65:     else high1 = nrow1; \
 66:     lastcol1 = col; \
 67:     while (high1 - low1 > 5) { \
 68:       t = (low1 + high1) / 2; \
 69:       if (cp1[sliceheight * t] > col) high1 = t; \
 70:       else low1 = t; \
 71:     } \
 72:     for (_i = low1; _i < high1; _i++) { \
 73:       if (cp1[sliceheight * _i] > col) break; \
 74:       if (cp1[sliceheight * _i] == col) { \
 75:         if (addv == ADD_VALUES) vp1[sliceheight * _i] += value; \
 76:         else vp1[sliceheight * _i] = value; \
 77:         inserted = PETSC_TRUE; \
 78:         goto a_noinsert; \
 79:       } \
 80:     } \
 81:     if (value == 0.0 && ignorezeroentries) { \
 82:       low1  = 0; \
 83:       high1 = nrow1; \
 84:       goto a_noinsert; \
 85:     } \
 86:     if (nonew == 1) { \
 87:       low1  = 0; \
 88:       high1 = nrow1; \
 89:       goto a_noinsert; \
 90:     } \
 91:     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); \
 92:     MatSeqXSELLReallocateSELL(A, am, 1, nrow1, a->sliidx, a->sliceheight, row / sliceheight, row, col, a->colidx, a->val, cp1, vp1, nonew, MatScalar); \
 93:     /* shift up all the later entries in this row */ \
 94:     for (ii = nrow1 - 1; ii >= _i; ii--) { \
 95:       cp1[sliceheight * (ii + 1)] = cp1[sliceheight * ii]; \
 96:       vp1[sliceheight * (ii + 1)] = vp1[sliceheight * ii]; \
 97:     } \
 98:     cp1[sliceheight * _i] = col; \
 99:     vp1[sliceheight * _i] = value; \
100:     a->nz++; \
101:     nrow1++; \
102:   a_noinsert:; \
103:     a->rlen[row] = nrow1; \
104:   }

106: #define MatSetValues_SeqSELL_B_Private(row, col, value, addv, orow, ocol) \
107:   { \
108:     if (col <= lastcol2) low2 = 0; \
109:     else high2 = nrow2; \
110:     lastcol2 = col; \
111:     while (high2 - low2 > 5) { \
112:       t = (low2 + high2) / 2; \
113:       if (cp2[sliceheight * t] > col) high2 = t; \
114:       else low2 = t; \
115:     } \
116:     for (_i = low2; _i < high2; _i++) { \
117:       if (cp2[sliceheight * _i] > col) break; \
118:       if (cp2[sliceheight * _i] == col) { \
119:         if (addv == ADD_VALUES) vp2[sliceheight * _i] += value; \
120:         else vp2[sliceheight * _i] = value; \
121:         inserted = PETSC_TRUE; \
122:         goto b_noinsert; \
123:       } \
124:     } \
125:     if (value == 0.0 && ignorezeroentries) { \
126:       low2  = 0; \
127:       high2 = nrow2; \
128:       goto b_noinsert; \
129:     } \
130:     if (nonew == 1) { \
131:       low2  = 0; \
132:       high2 = nrow2; \
133:       goto b_noinsert; \
134:     } \
135:     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); \
136:     MatSeqXSELLReallocateSELL(B, bm, 1, nrow2, b->sliidx, b->sliceheight, row / sliceheight, row, col, b->colidx, b->val, cp2, vp2, nonew, MatScalar); \
137:     /* shift up all the later entries in this row */ \
138:     for (ii = nrow2 - 1; ii >= _i; ii--) { \
139:       cp2[sliceheight * (ii + 1)] = cp2[sliceheight * ii]; \
140:       vp2[sliceheight * (ii + 1)] = vp2[sliceheight * ii]; \
141:     } \
142:     cp2[sliceheight * _i] = col; \
143:     vp2[sliceheight * _i] = value; \
144:     b->nz++; \
145:     nrow2++; \
146:   b_noinsert:; \
147:     b->rlen[row] = nrow2; \
148:   }

150: static PetscErrorCode MatSetValues_MPISELL(Mat mat, PetscInt m, const PetscInt im[], PetscInt n, const PetscInt in[], const PetscScalar v[], InsertMode addv)
151: {
152:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
153:   PetscScalar  value;
154:   PetscInt     i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend, shift1, shift2;
155:   PetscInt     cstart = mat->cmap->rstart, cend = mat->cmap->rend, row, col;
156:   PetscBool    roworiented = sell->roworiented;

158:   /* Some Variables required in the macro */
159:   Mat          A                 = sell->A;
160:   Mat_SeqSELL *a                 = (Mat_SeqSELL *)A->data;
161:   PetscBool    ignorezeroentries = a->ignorezeroentries, found;
162:   Mat          B                 = sell->B;
163:   Mat_SeqSELL *b                 = (Mat_SeqSELL *)B->data;
164:   PetscInt    *cp1, *cp2, ii, _i, nrow1, nrow2, low1, high1, low2, high2, t, lastcol1, lastcol2, sliceheight = a->sliceheight;
165:   MatScalar   *vp1, *vp2;

167:   PetscFunctionBegin;
168:   for (i = 0; i < m; i++) {
169:     if (im[i] < 0) continue;
170:     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);
171:     if (im[i] >= rstart && im[i] < rend) {
172:       row      = im[i] - rstart;
173:       lastcol1 = -1;
174:       shift1   = a->sliidx[row / sliceheight] + (row % sliceheight); /* starting index of the row */
175:       cp1      = PetscSafePointerPlusOffset(a->colidx, shift1);
176:       vp1      = PetscSafePointerPlusOffset(a->val, shift1);
177:       nrow1    = a->rlen[row];
178:       low1     = 0;
179:       high1    = nrow1;
180:       lastcol2 = -1;
181:       shift2   = b->sliidx[row / sliceheight] + (row % sliceheight); /* starting index of the row */
182:       cp2      = PetscSafePointerPlusOffset(b->colidx, shift2);
183:       vp2      = PetscSafePointerPlusOffset(b->val, shift2);
184:       nrow2    = b->rlen[row];
185:       low2     = 0;
186:       high2    = nrow2;

188:       for (j = 0; j < n; j++) {
189:         if (roworiented) value = v[i * n + j];
190:         else value = v[i + j * m];
191:         if (ignorezeroentries && value == 0.0 && (addv == ADD_VALUES)) continue;
192:         if (in[j] >= cstart && in[j] < cend) {
193:           col = in[j] - cstart;
194:           MatSetValue_SeqSELL_Private(A, row, col, value, addv, im[i], in[j], cp1, vp1, lastcol1, low1, high1); /* set one value */
195: #if defined(PETSC_HAVE_CUDA)
196:           if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED && found) A->offloadmask = PETSC_OFFLOAD_CPU;
197: #endif
198:         } else if (in[j] < 0) {
199:           continue;
200:         } else {
201:           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);
202:           if (mat->was_assembled) {
203:             if (!sell->colmap) PetscCall(MatCreateColmap_MPISELL_Private(mat));
204: #if defined(PETSC_USE_CTABLE)
205:             PetscCall(PetscHMapIGetWithDefault(sell->colmap, in[j] + 1, 0, &col));
206:             col--;
207: #else
208:             col = sell->colmap[in[j]] - 1;
209: #endif
210:             if (col < 0 && !((Mat_SeqSELL *)sell->B->data)->nonew) {
211:               PetscCall(MatDisAssemble_MPISELL(mat));
212:               col = in[j];
213:               /* Reinitialize the variables required by MatSetValues_SeqSELL_B_Private() */
214:               B      = sell->B;
215:               b      = (Mat_SeqSELL *)B->data;
216:               shift2 = b->sliidx[row / sliceheight] + (row % sliceheight); /* starting index of the row */
217:               cp2    = b->colidx + shift2;
218:               vp2    = b->val + shift2;
219:               nrow2  = b->rlen[row];
220:               low2   = 0;
221:               high2  = nrow2;
222:               found  = PETSC_FALSE;
223:             } else {
224:               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]);
225:             }
226:           } else col = in[j];
227:           MatSetValue_SeqSELL_Private(B, row, col, value, addv, im[i], in[j], cp2, vp2, lastcol2, low2, high2); /* set one value */
228: #if defined(PETSC_HAVE_CUDA)
229:           if (B->offloadmask != PETSC_OFFLOAD_UNALLOCATED && found) B->offloadmask = PETSC_OFFLOAD_CPU;
230: #endif
231:         }
232:       }
233:     } else {
234:       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]);
235:       if (!sell->donotstash) {
236:         mat->assembled = PETSC_FALSE;
237:         if (roworiented) {
238:           PetscCall(MatStashValuesRow_Private(&mat->stash, im[i], n, in, v + i * n, (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
239:         } else {
240:           PetscCall(MatStashValuesCol_Private(&mat->stash, im[i], n, in, v + i, m, (PetscBool)(ignorezeroentries && (addv == ADD_VALUES))));
241:         }
242:       }
243:     }
244:   }
245:   PetscFunctionReturn(PETSC_SUCCESS);
246: }

248: static PetscErrorCode MatGetValues_MPISELL(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[])
249: {
250:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
251:   PetscInt     i, j, rstart = mat->rmap->rstart, rend = mat->rmap->rend;
252:   PetscInt     cstart = mat->cmap->rstart, cend = mat->cmap->rend, row, col;

254:   PetscFunctionBegin;
255:   for (i = 0; i < m; i++) {
256:     if (idxm[i] < 0) continue; /* negative row */
257:     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);
258:     if (idxm[i] >= rstart && idxm[i] < rend) {
259:       row = idxm[i] - rstart;
260:       for (j = 0; j < n; j++) {
261:         if (idxn[j] < 0) continue; /* negative column */
262:         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);
263:         if (idxn[j] >= cstart && idxn[j] < cend) {
264:           col = idxn[j] - cstart;
265:           PetscCall(MatGetValues(sell->A, 1, &row, 1, &col, v + i * n + j));
266:         } else {
267:           if (!sell->colmap) PetscCall(MatCreateColmap_MPISELL_Private(mat));
268: #if defined(PETSC_USE_CTABLE)
269:           PetscCall(PetscHMapIGetWithDefault(sell->colmap, idxn[j] + 1, 0, &col));
270:           col--;
271: #else
272:           col = sell->colmap[idxn[j]] - 1;
273: #endif
274:           if ((col < 0) || (sell->garray[col] != idxn[j])) *(v + i * n + j) = 0.0;
275:           else PetscCall(MatGetValues(sell->B, 1, &row, 1, &col, v + i * n + j));
276:         }
277:       }
278:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only local values currently supported");
279:   }
280:   PetscFunctionReturn(PETSC_SUCCESS);
281: }

283: static PetscErrorCode MatAssemblyBegin_MPISELL(Mat mat, MatAssemblyType mode)
284: {
285:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
286:   PetscInt     nstash, reallocs;

288:   PetscFunctionBegin;
289:   if (sell->donotstash || mat->nooffprocentries) PetscFunctionReturn(PETSC_SUCCESS);

291:   PetscCall(MatStashScatterBegin_Private(mat, &mat->stash, mat->rmap->range));
292:   PetscCall(MatStashGetInfo_Private(&mat->stash, &nstash, &reallocs));
293:   PetscCall(PetscInfo(sell->A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
294:   PetscFunctionReturn(PETSC_SUCCESS);
295: }

297: PetscErrorCode MatAssemblyEnd_MPISELL(Mat mat, MatAssemblyType mode)
298: {
299:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
300:   PetscMPIInt  n;
301:   PetscInt     i, flg;
302:   PetscInt    *row, *col;
303:   PetscScalar *val;
304:   PetscBool    other_disassembled;
305:   /* do not use 'b = (Mat_SeqSELL*)sell->B->data' as B can be reset in disassembly */
306:   PetscFunctionBegin;
307:   if (!sell->donotstash && !mat->nooffprocentries) {
308:     while (1) {
309:       PetscCall(MatStashScatterGetMesg_Private(&mat->stash, &n, &row, &col, &val, &flg));
310:       if (!flg) break;

312:       for (i = 0; i < n; i++) { /* assemble one by one */
313:         PetscCall(MatSetValues_MPISELL(mat, 1, row + i, 1, col + i, val + i, mat->insertmode));
314:       }
315:     }
316:     PetscCall(MatStashScatterEnd_Private(&mat->stash));
317:   }
318: #if defined(PETSC_HAVE_CUDA)
319:   if (mat->offloadmask == PETSC_OFFLOAD_CPU) sell->A->offloadmask = PETSC_OFFLOAD_CPU;
320: #endif
321:   PetscCall(MatAssemblyBegin(sell->A, mode));
322:   PetscCall(MatAssemblyEnd(sell->A, mode));

324:   /*
325:      determine if any processor has disassembled, if so we must
326:      also disassemble ourselves, in order that we may reassemble.
327:   */
328:   /*
329:      if nonzero structure of submatrix B cannot change then we know that
330:      no processor disassembled thus we can skip this stuff
331:   */
332:   if (!((Mat_SeqSELL *)sell->B->data)->nonew) {
333:     PetscCallMPI(MPIU_Allreduce(&mat->was_assembled, &other_disassembled, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)mat)));
334:     if (mat->was_assembled && !other_disassembled) PetscCall(MatDisAssemble_MPISELL(mat));
335:   }
336:   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) PetscCall(MatSetUpMultiply_MPISELL(mat));
337: #if defined(PETSC_HAVE_CUDA)
338:   if (mat->offloadmask == PETSC_OFFLOAD_CPU && sell->B->offloadmask != PETSC_OFFLOAD_UNALLOCATED) sell->B->offloadmask = PETSC_OFFLOAD_CPU;
339: #endif
340:   PetscCall(MatAssemblyBegin(sell->B, mode));
341:   PetscCall(MatAssemblyEnd(sell->B, mode));
342:   PetscCall(PetscFree2(sell->rowvalues, sell->rowindices));
343:   sell->rowvalues = NULL;
344:   PetscCall(VecDestroy(&sell->diag));

346:   /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
347:   if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqSELL *)sell->A->data)->nonew) {
348:     PetscObjectState state = sell->A->nonzerostate + sell->B->nonzerostate;
349:     PetscCallMPI(MPIU_Allreduce(&state, &mat->nonzerostate, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)mat)));
350:   }
351: #if defined(PETSC_HAVE_CUDA)
352:   mat->offloadmask = PETSC_OFFLOAD_BOTH;
353: #endif
354:   PetscFunctionReturn(PETSC_SUCCESS);
355: }

357: static PetscErrorCode MatZeroEntries_MPISELL(Mat A)
358: {
359:   Mat_MPISELL *l = (Mat_MPISELL *)A->data;

361:   PetscFunctionBegin;
362:   PetscCall(MatZeroEntries(l->A));
363:   PetscCall(MatZeroEntries(l->B));
364:   PetscFunctionReturn(PETSC_SUCCESS);
365: }

367: static PetscErrorCode MatMult_MPISELL(Mat A, Vec xx, Vec yy)
368: {
369:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
370:   PetscInt     nt;

372:   PetscFunctionBegin;
373:   PetscCall(VecGetLocalSize(xx, &nt));
374:   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);
375:   PetscCall(VecScatterBegin(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
376:   PetscCall((*a->A->ops->mult)(a->A, xx, yy));
377:   PetscCall(VecScatterEnd(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
378:   PetscCall((*a->B->ops->multadd)(a->B, a->lvec, yy, yy));
379:   PetscFunctionReturn(PETSC_SUCCESS);
380: }

382: static PetscErrorCode MatMultDiagonalBlock_MPISELL(Mat A, Vec bb, Vec xx)
383: {
384:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

386:   PetscFunctionBegin;
387:   PetscCall(MatMultDiagonalBlock(a->A, bb, xx));
388:   PetscFunctionReturn(PETSC_SUCCESS);
389: }

391: static PetscErrorCode MatMultAdd_MPISELL(Mat A, Vec xx, Vec yy, Vec zz)
392: {
393:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

395:   PetscFunctionBegin;
396:   PetscCall(VecScatterBegin(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
397:   PetscCall((*a->A->ops->multadd)(a->A, xx, yy, zz));
398:   PetscCall(VecScatterEnd(a->Mvctx, xx, a->lvec, INSERT_VALUES, SCATTER_FORWARD));
399:   PetscCall((*a->B->ops->multadd)(a->B, a->lvec, zz, zz));
400:   PetscFunctionReturn(PETSC_SUCCESS);
401: }

403: static PetscErrorCode MatMultTranspose_MPISELL(Mat A, Vec xx, Vec yy)
404: {
405:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

407:   PetscFunctionBegin;
408:   /* do nondiagonal part */
409:   PetscCall((*a->B->ops->multtranspose)(a->B, xx, a->lvec));
410:   /* do local part */
411:   PetscCall((*a->A->ops->multtranspose)(a->A, xx, yy));
412:   /* add partial results together */
413:   PetscCall(VecScatterBegin(a->Mvctx, a->lvec, yy, ADD_VALUES, SCATTER_REVERSE));
414:   PetscCall(VecScatterEnd(a->Mvctx, a->lvec, yy, ADD_VALUES, SCATTER_REVERSE));
415:   PetscFunctionReturn(PETSC_SUCCESS);
416: }

418: static PetscErrorCode MatIsTranspose_MPISELL(Mat Amat, Mat Bmat, PetscReal tol, PetscBool *f)
419: {
420:   MPI_Comm     comm;
421:   Mat_MPISELL *Asell = (Mat_MPISELL *)Amat->data, *Bsell;
422:   Mat          Adia  = Asell->A, Bdia, Aoff, Boff, *Aoffs, *Boffs;
423:   IS           Me, Notme;
424:   PetscInt     M, N, first, last, *notme, i;
425:   PetscMPIInt  size;

427:   PetscFunctionBegin;
428:   /* Easy test: symmetric diagonal block */
429:   Bsell = (Mat_MPISELL *)Bmat->data;
430:   Bdia  = Bsell->A;
431:   PetscCall(MatIsTranspose(Adia, Bdia, tol, f));
432:   if (!*f) PetscFunctionReturn(PETSC_SUCCESS);
433:   PetscCall(PetscObjectGetComm((PetscObject)Amat, &comm));
434:   PetscCallMPI(MPI_Comm_size(comm, &size));
435:   if (size == 1) PetscFunctionReturn(PETSC_SUCCESS);

437:   /* Hard test: off-diagonal block. This takes a MatCreateSubMatrix. */
438:   PetscCall(MatGetSize(Amat, &M, &N));
439:   PetscCall(MatGetOwnershipRange(Amat, &first, &last));
440:   PetscCall(PetscMalloc1(N - last + first, &notme));
441:   for (i = 0; i < first; i++) notme[i] = i;
442:   for (i = last; i < M; i++) notme[i - last + first] = i;
443:   PetscCall(ISCreateGeneral(MPI_COMM_SELF, N - last + first, notme, PETSC_COPY_VALUES, &Notme));
444:   PetscCall(ISCreateStride(MPI_COMM_SELF, last - first, first, 1, &Me));
445:   PetscCall(MatCreateSubMatrices(Amat, 1, &Me, &Notme, MAT_INITIAL_MATRIX, &Aoffs));
446:   Aoff = Aoffs[0];
447:   PetscCall(MatCreateSubMatrices(Bmat, 1, &Notme, &Me, MAT_INITIAL_MATRIX, &Boffs));
448:   Boff = Boffs[0];
449:   PetscCall(MatIsTranspose(Aoff, Boff, tol, f));
450:   PetscCall(MatDestroyMatrices(1, &Aoffs));
451:   PetscCall(MatDestroyMatrices(1, &Boffs));
452:   PetscCall(ISDestroy(&Me));
453:   PetscCall(ISDestroy(&Notme));
454:   PetscCall(PetscFree(notme));
455:   PetscFunctionReturn(PETSC_SUCCESS);
456: }

458: static PetscErrorCode MatMultTransposeAdd_MPISELL(Mat A, Vec xx, Vec yy, Vec zz)
459: {
460:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

462:   PetscFunctionBegin;
463:   /* do nondiagonal part */
464:   PetscCall((*a->B->ops->multtranspose)(a->B, xx, a->lvec));
465:   /* do local part */
466:   PetscCall((*a->A->ops->multtransposeadd)(a->A, xx, yy, zz));
467:   /* add partial results together */
468:   PetscCall(VecScatterBegin(a->Mvctx, a->lvec, zz, ADD_VALUES, SCATTER_REVERSE));
469:   PetscCall(VecScatterEnd(a->Mvctx, a->lvec, zz, ADD_VALUES, SCATTER_REVERSE));
470:   PetscFunctionReturn(PETSC_SUCCESS);
471: }

473: /*
474:   This only works correctly for square matrices where the subblock A->A is the
475:    diagonal block
476: */
477: static PetscErrorCode MatGetDiagonal_MPISELL(Mat A, Vec v)
478: {
479:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

481:   PetscFunctionBegin;
482:   PetscCheck(A->rmap->N == A->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Supports only square matrix where A->A is diag block");
483:   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");
484:   PetscCall(MatGetDiagonal(a->A, v));
485:   PetscFunctionReturn(PETSC_SUCCESS);
486: }

488: static PetscErrorCode MatScale_MPISELL(Mat A, PetscScalar aa)
489: {
490:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

492:   PetscFunctionBegin;
493:   PetscCall(MatScale(a->A, aa));
494:   PetscCall(MatScale(a->B, aa));
495:   PetscFunctionReturn(PETSC_SUCCESS);
496: }

498: PetscErrorCode MatDestroy_MPISELL(Mat mat)
499: {
500:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;

502:   PetscFunctionBegin;
503:   PetscCall(PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT, mat->rmap->N, mat->cmap->N));
504:   PetscCall(MatStashDestroy_Private(&mat->stash));
505:   PetscCall(VecDestroy(&sell->diag));
506:   PetscCall(MatDestroy(&sell->A));
507:   PetscCall(MatDestroy(&sell->B));
508: #if defined(PETSC_USE_CTABLE)
509:   PetscCall(PetscHMapIDestroy(&sell->colmap));
510: #else
511:   PetscCall(PetscFree(sell->colmap));
512: #endif
513:   PetscCall(PetscFree(sell->garray));
514:   PetscCall(VecDestroy(&sell->lvec));
515:   PetscCall(VecScatterDestroy(&sell->Mvctx));
516:   PetscCall(PetscFree2(sell->rowvalues, sell->rowindices));
517:   PetscCall(PetscFree(sell->ld));
518:   PetscCall(PetscFree(mat->data));

520:   PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL));
521:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatStoreValues_C", NULL));
522:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatRetrieveValues_C", NULL));
523:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatIsTranspose_C", NULL));
524:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPISELLSetPreallocation_C", NULL));
525:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisell_mpiaij_C", NULL));
526: #if defined(PETSC_HAVE_CUDA)
527:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConvert_mpisell_mpisellcuda_C", NULL));
528: #endif
529:   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatDiagonalScaleLocal_C", NULL));
530:   PetscFunctionReturn(PETSC_SUCCESS);
531: }

533: #include <petscdraw.h>
534: static PetscErrorCode MatView_MPISELL_ASCIIorDraworSocket(Mat mat, PetscViewer viewer)
535: {
536:   Mat_MPISELL      *sell = (Mat_MPISELL *)mat->data;
537:   PetscMPIInt       rank = sell->rank, size = sell->size;
538:   PetscBool         isdraw, iascii, isbinary;
539:   PetscViewer       sviewer;
540:   PetscViewerFormat format;

542:   PetscFunctionBegin;
543:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
544:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
545:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
546:   if (iascii) {
547:     PetscCall(PetscViewerGetFormat(viewer, &format));
548:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
549:       MatInfo   info;
550:       PetscInt *inodes;

552:       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat), &rank));
553:       PetscCall(MatGetInfo(mat, MAT_LOCAL, &info));
554:       PetscCall(MatInodeGetInodeSizes(sell->A, NULL, &inodes, NULL));
555:       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
556:       if (!inodes) {
557:         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,
558:                                                      (PetscInt)info.nz_allocated, (PetscInt)info.memory));
559:       } else {
560:         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,
561:                                                      (PetscInt)info.nz_allocated, (PetscInt)info.memory));
562:       }
563:       PetscCall(MatGetInfo(sell->A, MAT_LOCAL, &info));
564:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] on-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
565:       PetscCall(MatGetInfo(sell->B, MAT_LOCAL, &info));
566:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] off-diagonal part: nz %" PetscInt_FMT " \n", rank, (PetscInt)info.nz_used));
567:       PetscCall(PetscViewerFlush(viewer));
568:       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
569:       PetscCall(PetscViewerASCIIPrintf(viewer, "Information on VecScatter used in matrix-vector product: \n"));
570:       PetscCall(VecScatterView(sell->Mvctx, viewer));
571:       PetscFunctionReturn(PETSC_SUCCESS);
572:     } else if (format == PETSC_VIEWER_ASCII_INFO) {
573:       PetscInt inodecount, inodelimit, *inodes;
574:       PetscCall(MatInodeGetInodeSizes(sell->A, &inodecount, &inodes, &inodelimit));
575:       if (inodes) {
576:         PetscCall(PetscViewerASCIIPrintf(viewer, "using I-node (on process 0) routines: found %" PetscInt_FMT " nodes, limit used is %" PetscInt_FMT "\n", inodecount, inodelimit));
577:       } else {
578:         PetscCall(PetscViewerASCIIPrintf(viewer, "not using I-node (on process 0) routines\n"));
579:       }
580:       PetscFunctionReturn(PETSC_SUCCESS);
581:     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
582:       PetscFunctionReturn(PETSC_SUCCESS);
583:     }
584:   } else if (isbinary) {
585:     if (size == 1) {
586:       PetscCall(PetscObjectSetName((PetscObject)sell->A, ((PetscObject)mat)->name));
587:       PetscCall(MatView(sell->A, viewer));
588:     } else {
589:       /* PetscCall(MatView_MPISELL_Binary(mat,viewer)); */
590:     }
591:     PetscFunctionReturn(PETSC_SUCCESS);
592:   } else if (isdraw) {
593:     PetscDraw draw;
594:     PetscBool isnull;
595:     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
596:     PetscCall(PetscDrawIsNull(draw, &isnull));
597:     if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
598:   }

600:   {
601:     /* assemble the entire matrix onto first processor. */
602:     Mat          A;
603:     Mat_SeqSELL *Aloc;
604:     PetscInt     M = mat->rmap->N, N = mat->cmap->N, *acolidx, row, col, i, j;
605:     MatScalar   *aval;
606:     PetscBool    isnonzero;

608:     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &A));
609:     if (rank == 0) {
610:       PetscCall(MatSetSizes(A, M, N, M, N));
611:     } else {
612:       PetscCall(MatSetSizes(A, 0, 0, M, N));
613:     }
614:     /* This is just a temporary matrix, so explicitly using MATMPISELL is probably best */
615:     PetscCall(MatSetType(A, MATMPISELL));
616:     PetscCall(MatMPISELLSetPreallocation(A, 0, NULL, 0, NULL));
617:     PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));

619:     /* copy over the A part */
620:     Aloc    = (Mat_SeqSELL *)sell->A->data;
621:     acolidx = Aloc->colidx;
622:     aval    = Aloc->val;
623:     for (i = 0; i < Aloc->totalslices; i++) { /* loop over slices */
624:       for (j = Aloc->sliidx[i]; j < Aloc->sliidx[i + 1]; j++) {
625:         isnonzero = (PetscBool)((j - Aloc->sliidx[i]) / Aloc->sliceheight < Aloc->rlen[i * Aloc->sliceheight + j % Aloc->sliceheight]);
626:         if (isnonzero) { /* check the mask bit */
627:           row = i * Aloc->sliceheight + j % Aloc->sliceheight + mat->rmap->rstart;
628:           col = *acolidx + mat->rmap->rstart;
629:           PetscCall(MatSetValues(A, 1, &row, 1, &col, aval, INSERT_VALUES));
630:         }
631:         aval++;
632:         acolidx++;
633:       }
634:     }

636:     /* copy over the B part */
637:     Aloc    = (Mat_SeqSELL *)sell->B->data;
638:     acolidx = Aloc->colidx;
639:     aval    = Aloc->val;
640:     for (i = 0; i < Aloc->totalslices; i++) {
641:       for (j = Aloc->sliidx[i]; j < Aloc->sliidx[i + 1]; j++) {
642:         isnonzero = (PetscBool)((j - Aloc->sliidx[i]) / Aloc->sliceheight < Aloc->rlen[i * Aloc->sliceheight + j % Aloc->sliceheight]);
643:         if (isnonzero) {
644:           row = i * Aloc->sliceheight + j % Aloc->sliceheight + mat->rmap->rstart;
645:           col = sell->garray[*acolidx];
646:           PetscCall(MatSetValues(A, 1, &row, 1, &col, aval, INSERT_VALUES));
647:         }
648:         aval++;
649:         acolidx++;
650:       }
651:     }

653:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
654:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
655:     /*
656:        Everyone has to call to draw the matrix since the graphics waits are
657:        synchronized across all processors that share the PetscDraw object
658:     */
659:     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
660:     if (rank == 0) {
661:       PetscCall(PetscObjectSetName((PetscObject)((Mat_MPISELL *)A->data)->A, ((PetscObject)mat)->name));
662:       PetscCall(MatView_SeqSELL(((Mat_MPISELL *)A->data)->A, sviewer));
663:     }
664:     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
665:     PetscCall(MatDestroy(&A));
666:   }
667:   PetscFunctionReturn(PETSC_SUCCESS);
668: }

670: static PetscErrorCode MatView_MPISELL(Mat mat, PetscViewer viewer)
671: {
672:   PetscBool iascii, isdraw, issocket, isbinary;

674:   PetscFunctionBegin;
675:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
676:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
677:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
678:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSOCKET, &issocket));
679:   if (iascii || isdraw || isbinary || issocket) PetscCall(MatView_MPISELL_ASCIIorDraworSocket(mat, viewer));
680:   PetscFunctionReturn(PETSC_SUCCESS);
681: }

683: static PetscErrorCode MatGetGhosts_MPISELL(Mat mat, PetscInt *nghosts, const PetscInt *ghosts[])
684: {
685:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;

687:   PetscFunctionBegin;
688:   PetscCall(MatGetSize(sell->B, NULL, nghosts));
689:   if (ghosts) *ghosts = sell->garray;
690:   PetscFunctionReturn(PETSC_SUCCESS);
691: }

693: static PetscErrorCode MatGetInfo_MPISELL(Mat matin, MatInfoType flag, MatInfo *info)
694: {
695:   Mat_MPISELL   *mat = (Mat_MPISELL *)matin->data;
696:   Mat            A = mat->A, B = mat->B;
697:   PetscLogDouble isend[5], irecv[5];

699:   PetscFunctionBegin;
700:   info->block_size = 1.0;
701:   PetscCall(MatGetInfo(A, MAT_LOCAL, info));

703:   isend[0] = info->nz_used;
704:   isend[1] = info->nz_allocated;
705:   isend[2] = info->nz_unneeded;
706:   isend[3] = info->memory;
707:   isend[4] = info->mallocs;

709:   PetscCall(MatGetInfo(B, MAT_LOCAL, info));

711:   isend[0] += info->nz_used;
712:   isend[1] += info->nz_allocated;
713:   isend[2] += info->nz_unneeded;
714:   isend[3] += info->memory;
715:   isend[4] += info->mallocs;
716:   if (flag == MAT_LOCAL) {
717:     info->nz_used      = isend[0];
718:     info->nz_allocated = isend[1];
719:     info->nz_unneeded  = isend[2];
720:     info->memory       = isend[3];
721:     info->mallocs      = isend[4];
722:   } else if (flag == MAT_GLOBAL_MAX) {
723:     PetscCallMPI(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_MAX, PetscObjectComm((PetscObject)matin)));

725:     info->nz_used      = irecv[0];
726:     info->nz_allocated = irecv[1];
727:     info->nz_unneeded  = irecv[2];
728:     info->memory       = irecv[3];
729:     info->mallocs      = irecv[4];
730:   } else if (flag == MAT_GLOBAL_SUM) {
731:     PetscCallMPI(MPIU_Allreduce(isend, irecv, 5, MPIU_PETSCLOGDOUBLE, MPI_SUM, PetscObjectComm((PetscObject)matin)));

733:     info->nz_used      = irecv[0];
734:     info->nz_allocated = irecv[1];
735:     info->nz_unneeded  = irecv[2];
736:     info->memory       = irecv[3];
737:     info->mallocs      = irecv[4];
738:   }
739:   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
740:   info->fill_ratio_needed = 0;
741:   info->factor_mallocs    = 0;
742:   PetscFunctionReturn(PETSC_SUCCESS);
743: }

745: static PetscErrorCode MatSetOption_MPISELL(Mat A, MatOption op, PetscBool flg)
746: {
747:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

749:   PetscFunctionBegin;
750:   switch (op) {
751:   case MAT_NEW_NONZERO_LOCATIONS:
752:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
753:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
754:   case MAT_KEEP_NONZERO_PATTERN:
755:   case MAT_NEW_NONZERO_LOCATION_ERR:
756:   case MAT_USE_INODES:
757:   case MAT_IGNORE_ZERO_ENTRIES:
758:     MatCheckPreallocated(A, 1);
759:     PetscCall(MatSetOption(a->A, op, flg));
760:     PetscCall(MatSetOption(a->B, op, flg));
761:     break;
762:   case MAT_ROW_ORIENTED:
763:     MatCheckPreallocated(A, 1);
764:     a->roworiented = flg;

766:     PetscCall(MatSetOption(a->A, op, flg));
767:     PetscCall(MatSetOption(a->B, op, flg));
768:     break;
769:   case MAT_IGNORE_OFF_PROC_ENTRIES:
770:     a->donotstash = flg;
771:     break;
772:   case MAT_SYMMETRIC:
773:     MatCheckPreallocated(A, 1);
774:     PetscCall(MatSetOption(a->A, op, flg));
775:     break;
776:   case MAT_STRUCTURALLY_SYMMETRIC:
777:     MatCheckPreallocated(A, 1);
778:     PetscCall(MatSetOption(a->A, op, flg));
779:     break;
780:   case MAT_HERMITIAN:
781:     MatCheckPreallocated(A, 1);
782:     PetscCall(MatSetOption(a->A, op, flg));
783:     break;
784:   case MAT_SYMMETRY_ETERNAL:
785:     MatCheckPreallocated(A, 1);
786:     PetscCall(MatSetOption(a->A, op, flg));
787:     break;
788:   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
789:     MatCheckPreallocated(A, 1);
790:     PetscCall(MatSetOption(a->A, op, flg));
791:     break;
792:   default:
793:     break;
794:   }
795:   PetscFunctionReturn(PETSC_SUCCESS);
796: }

798: static PetscErrorCode MatDiagonalScale_MPISELL(Mat mat, Vec ll, Vec rr)
799: {
800:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;
801:   Mat          a = sell->A, b = sell->B;
802:   PetscInt     s1, s2, s3;

804:   PetscFunctionBegin;
805:   PetscCall(MatGetLocalSize(mat, &s2, &s3));
806:   if (rr) {
807:     PetscCall(VecGetLocalSize(rr, &s1));
808:     PetscCheck(s1 == s3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "right vector non-conforming local size");
809:     /* Overlap communication with computation. */
810:     PetscCall(VecScatterBegin(sell->Mvctx, rr, sell->lvec, INSERT_VALUES, SCATTER_FORWARD));
811:   }
812:   if (ll) {
813:     PetscCall(VecGetLocalSize(ll, &s1));
814:     PetscCheck(s1 == s2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "left vector non-conforming local size");
815:     PetscUseTypeMethod(b, diagonalscale, ll, NULL);
816:   }
817:   /* scale  the diagonal block */
818:   PetscUseTypeMethod(a, diagonalscale, ll, rr);

820:   if (rr) {
821:     /* Do a scatter end and then right scale the off-diagonal block */
822:     PetscCall(VecScatterEnd(sell->Mvctx, rr, sell->lvec, INSERT_VALUES, SCATTER_FORWARD));
823:     PetscUseTypeMethod(b, diagonalscale, NULL, sell->lvec);
824:   }
825:   PetscFunctionReturn(PETSC_SUCCESS);
826: }

828: static PetscErrorCode MatSetUnfactored_MPISELL(Mat A)
829: {
830:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

832:   PetscFunctionBegin;
833:   PetscCall(MatSetUnfactored(a->A));
834:   PetscFunctionReturn(PETSC_SUCCESS);
835: }

837: static PetscErrorCode MatEqual_MPISELL(Mat A, Mat B, PetscBool *flag)
838: {
839:   Mat_MPISELL *matB = (Mat_MPISELL *)B->data, *matA = (Mat_MPISELL *)A->data;
840:   Mat          a, b, c, d;
841:   PetscBool    flg;

843:   PetscFunctionBegin;
844:   a = matA->A;
845:   b = matA->B;
846:   c = matB->A;
847:   d = matB->B;

849:   PetscCall(MatEqual(a, c, &flg));
850:   if (flg) PetscCall(MatEqual(b, d, &flg));
851:   PetscCallMPI(MPIU_Allreduce(&flg, flag, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
852:   PetscFunctionReturn(PETSC_SUCCESS);
853: }

855: static PetscErrorCode MatCopy_MPISELL(Mat A, Mat B, MatStructure str)
856: {
857:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
858:   Mat_MPISELL *b = (Mat_MPISELL *)B->data;

860:   PetscFunctionBegin;
861:   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
862:   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
863:     /* because of the column compression in the off-processor part of the matrix a->B,
864:        the number of columns in a->B and b->B may be different, hence we cannot call
865:        the MatCopy() directly on the two parts. If need be, we can provide a more
866:        efficient copy than the MatCopy_Basic() by first uncompressing the a->B matrices
867:        then copying the submatrices */
868:     PetscCall(MatCopy_Basic(A, B, str));
869:   } else {
870:     PetscCall(MatCopy(a->A, b->A, str));
871:     PetscCall(MatCopy(a->B, b->B, str));
872:   }
873:   PetscFunctionReturn(PETSC_SUCCESS);
874: }

876: static PetscErrorCode MatSetUp_MPISELL(Mat A)
877: {
878:   PetscFunctionBegin;
879:   PetscCall(MatMPISELLSetPreallocation(A, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
880:   PetscFunctionReturn(PETSC_SUCCESS);
881: }

883: static PetscErrorCode MatConjugate_MPISELL(Mat mat)
884: {
885:   PetscFunctionBegin;
886:   if (PetscDefined(USE_COMPLEX)) {
887:     Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;

889:     PetscCall(MatConjugate_SeqSELL(sell->A));
890:     PetscCall(MatConjugate_SeqSELL(sell->B));
891:   }
892:   PetscFunctionReturn(PETSC_SUCCESS);
893: }

895: static PetscErrorCode MatRealPart_MPISELL(Mat A)
896: {
897:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

899:   PetscFunctionBegin;
900:   PetscCall(MatRealPart(a->A));
901:   PetscCall(MatRealPart(a->B));
902:   PetscFunctionReturn(PETSC_SUCCESS);
903: }

905: static PetscErrorCode MatImaginaryPart_MPISELL(Mat A)
906: {
907:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

909:   PetscFunctionBegin;
910:   PetscCall(MatImaginaryPart(a->A));
911:   PetscCall(MatImaginaryPart(a->B));
912:   PetscFunctionReturn(PETSC_SUCCESS);
913: }

915: static PetscErrorCode MatInvertBlockDiagonal_MPISELL(Mat A, const PetscScalar **values)
916: {
917:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

919:   PetscFunctionBegin;
920:   PetscCall(MatInvertBlockDiagonal(a->A, values));
921:   A->factorerrortype = a->A->factorerrortype;
922:   PetscFunctionReturn(PETSC_SUCCESS);
923: }

925: static PetscErrorCode MatSetRandom_MPISELL(Mat x, PetscRandom rctx)
926: {
927:   Mat_MPISELL *sell = (Mat_MPISELL *)x->data;

929:   PetscFunctionBegin;
930:   PetscCall(MatSetRandom(sell->A, rctx));
931:   PetscCall(MatSetRandom(sell->B, rctx));
932:   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
933:   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
934:   PetscFunctionReturn(PETSC_SUCCESS);
935: }

937: static PetscErrorCode MatSetFromOptions_MPISELL(Mat A, PetscOptionItems PetscOptionsObject)
938: {
939:   PetscFunctionBegin;
940:   PetscOptionsHeadBegin(PetscOptionsObject, "MPISELL options");
941:   PetscOptionsHeadEnd();
942:   PetscFunctionReturn(PETSC_SUCCESS);
943: }

945: static PetscErrorCode MatShift_MPISELL(Mat Y, PetscScalar a)
946: {
947:   Mat_MPISELL *msell = (Mat_MPISELL *)Y->data;
948:   Mat_SeqSELL *sell  = (Mat_SeqSELL *)msell->A->data;

950:   PetscFunctionBegin;
951:   if (!Y->preallocated) {
952:     PetscCall(MatMPISELLSetPreallocation(Y, 1, NULL, 0, NULL));
953:   } else if (!sell->nz) {
954:     PetscInt nonew = sell->nonew;
955:     PetscCall(MatSeqSELLSetPreallocation(msell->A, 1, NULL));
956:     sell->nonew = nonew;
957:   }
958:   PetscCall(MatShift_Basic(Y, a));
959:   PetscFunctionReturn(PETSC_SUCCESS);
960: }

962: static PetscErrorCode MatMissingDiagonal_MPISELL(Mat A, PetscBool *missing, PetscInt *d)
963: {
964:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;

966:   PetscFunctionBegin;
967:   PetscCheck(A->rmap->n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only works for square matrices");
968:   PetscCall(MatMissingDiagonal(a->A, missing, d));
969:   if (d) {
970:     PetscInt rstart;
971:     PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
972:     *d += rstart;
973:   }
974:   PetscFunctionReturn(PETSC_SUCCESS);
975: }

977: static PetscErrorCode MatGetDiagonalBlock_MPISELL(Mat A, Mat *a)
978: {
979:   PetscFunctionBegin;
980:   *a = ((Mat_MPISELL *)A->data)->A;
981:   PetscFunctionReturn(PETSC_SUCCESS);
982: }

984: static PetscErrorCode MatStoreValues_MPISELL(Mat mat)
985: {
986:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;

988:   PetscFunctionBegin;
989:   PetscCall(MatStoreValues(sell->A));
990:   PetscCall(MatStoreValues(sell->B));
991:   PetscFunctionReturn(PETSC_SUCCESS);
992: }

994: static PetscErrorCode MatRetrieveValues_MPISELL(Mat mat)
995: {
996:   Mat_MPISELL *sell = (Mat_MPISELL *)mat->data;

998:   PetscFunctionBegin;
999:   PetscCall(MatRetrieveValues(sell->A));
1000:   PetscCall(MatRetrieveValues(sell->B));
1001:   PetscFunctionReturn(PETSC_SUCCESS);
1002: }

1004: static PetscErrorCode MatMPISELLSetPreallocation_MPISELL(Mat B, PetscInt d_rlenmax, const PetscInt d_rlen[], PetscInt o_rlenmax, const PetscInt o_rlen[])
1005: {
1006:   Mat_MPISELL *b;

1008:   PetscFunctionBegin;
1009:   PetscCall(PetscLayoutSetUp(B->rmap));
1010:   PetscCall(PetscLayoutSetUp(B->cmap));
1011:   b = (Mat_MPISELL *)B->data;

1013:   if (!B->preallocated) {
1014:     /* Explicitly create 2 MATSEQSELL matrices. */
1015:     PetscCall(MatCreate(PETSC_COMM_SELF, &b->A));
1016:     PetscCall(MatSetSizes(b->A, B->rmap->n, B->cmap->n, B->rmap->n, B->cmap->n));
1017:     PetscCall(MatSetBlockSizesFromMats(b->A, B, B));
1018:     PetscCall(MatSetType(b->A, MATSEQSELL));
1019:     PetscCall(MatCreate(PETSC_COMM_SELF, &b->B));
1020:     PetscCall(MatSetSizes(b->B, B->rmap->n, B->cmap->N, B->rmap->n, B->cmap->N));
1021:     PetscCall(MatSetBlockSizesFromMats(b->B, B, B));
1022:     PetscCall(MatSetType(b->B, MATSEQSELL));
1023:   }

1025:   PetscCall(MatSeqSELLSetPreallocation(b->A, d_rlenmax, d_rlen));
1026:   PetscCall(MatSeqSELLSetPreallocation(b->B, o_rlenmax, o_rlen));
1027:   B->preallocated  = PETSC_TRUE;
1028:   B->was_assembled = PETSC_FALSE;
1029:   /*
1030:     critical for MatAssemblyEnd to work.
1031:     MatAssemblyBegin checks it to set up was_assembled
1032:     and MatAssemblyEnd checks was_assembled to determine whether to build garray
1033:   */
1034:   B->assembled = PETSC_FALSE;
1035:   PetscFunctionReturn(PETSC_SUCCESS);
1036: }

1038: static PetscErrorCode MatDuplicate_MPISELL(Mat matin, MatDuplicateOption cpvalues, Mat *newmat)
1039: {
1040:   Mat          mat;
1041:   Mat_MPISELL *a, *oldmat = (Mat_MPISELL *)matin->data;

1043:   PetscFunctionBegin;
1044:   *newmat = NULL;
1045:   PetscCall(MatCreate(PetscObjectComm((PetscObject)matin), &mat));
1046:   PetscCall(MatSetSizes(mat, matin->rmap->n, matin->cmap->n, matin->rmap->N, matin->cmap->N));
1047:   PetscCall(MatSetBlockSizesFromMats(mat, matin, matin));
1048:   PetscCall(MatSetType(mat, ((PetscObject)matin)->type_name));
1049:   a = (Mat_MPISELL *)mat->data;

1051:   mat->factortype   = matin->factortype;
1052:   mat->assembled    = PETSC_TRUE;
1053:   mat->insertmode   = NOT_SET_VALUES;
1054:   mat->preallocated = PETSC_TRUE;

1056:   a->size         = oldmat->size;
1057:   a->rank         = oldmat->rank;
1058:   a->donotstash   = oldmat->donotstash;
1059:   a->roworiented  = oldmat->roworiented;
1060:   a->rowindices   = NULL;
1061:   a->rowvalues    = NULL;
1062:   a->getrowactive = PETSC_FALSE;

1064:   PetscCall(PetscLayoutReference(matin->rmap, &mat->rmap));
1065:   PetscCall(PetscLayoutReference(matin->cmap, &mat->cmap));

1067:   if (oldmat->colmap) {
1068: #if defined(PETSC_USE_CTABLE)
1069:     PetscCall(PetscHMapIDuplicate(oldmat->colmap, &a->colmap));
1070: #else
1071:     PetscCall(PetscMalloc1(mat->cmap->N, &a->colmap));
1072:     PetscCall(PetscArraycpy(a->colmap, oldmat->colmap, mat->cmap->N));
1073: #endif
1074:   } else a->colmap = NULL;
1075:   if (oldmat->garray) {
1076:     PetscInt len;
1077:     len = oldmat->B->cmap->n;
1078:     PetscCall(PetscMalloc1(len + 1, &a->garray));
1079:     if (len) PetscCall(PetscArraycpy(a->garray, oldmat->garray, len));
1080:   } else a->garray = NULL;

1082:   PetscCall(VecDuplicate(oldmat->lvec, &a->lvec));
1083:   PetscCall(VecScatterCopy(oldmat->Mvctx, &a->Mvctx));
1084:   PetscCall(MatDuplicate(oldmat->A, cpvalues, &a->A));
1085:   PetscCall(MatDuplicate(oldmat->B, cpvalues, &a->B));
1086:   PetscCall(PetscFunctionListDuplicate(((PetscObject)matin)->qlist, &((PetscObject)mat)->qlist));
1087:   *newmat = mat;
1088:   PetscFunctionReturn(PETSC_SUCCESS);
1089: }

1091: static const struct _MatOps MatOps_Values = {MatSetValues_MPISELL,
1092:                                              NULL,
1093:                                              NULL,
1094:                                              MatMult_MPISELL,
1095:                                              /* 4*/ MatMultAdd_MPISELL,
1096:                                              MatMultTranspose_MPISELL,
1097:                                              MatMultTransposeAdd_MPISELL,
1098:                                              NULL,
1099:                                              NULL,
1100:                                              NULL,
1101:                                              /*10*/ NULL,
1102:                                              NULL,
1103:                                              NULL,
1104:                                              MatSOR_MPISELL,
1105:                                              NULL,
1106:                                              /*15*/ MatGetInfo_MPISELL,
1107:                                              MatEqual_MPISELL,
1108:                                              MatGetDiagonal_MPISELL,
1109:                                              MatDiagonalScale_MPISELL,
1110:                                              NULL,
1111:                                              /*20*/ MatAssemblyBegin_MPISELL,
1112:                                              MatAssemblyEnd_MPISELL,
1113:                                              MatSetOption_MPISELL,
1114:                                              MatZeroEntries_MPISELL,
1115:                                              /*24*/ NULL,
1116:                                              NULL,
1117:                                              NULL,
1118:                                              NULL,
1119:                                              NULL,
1120:                                              /*29*/ MatSetUp_MPISELL,
1121:                                              NULL,
1122:                                              NULL,
1123:                                              MatGetDiagonalBlock_MPISELL,
1124:                                              NULL,
1125:                                              /*34*/ MatDuplicate_MPISELL,
1126:                                              NULL,
1127:                                              NULL,
1128:                                              NULL,
1129:                                              NULL,
1130:                                              /*39*/ NULL,
1131:                                              NULL,
1132:                                              NULL,
1133:                                              MatGetValues_MPISELL,
1134:                                              MatCopy_MPISELL,
1135:                                              /*44*/ NULL,
1136:                                              MatScale_MPISELL,
1137:                                              MatShift_MPISELL,
1138:                                              MatDiagonalSet_MPISELL,
1139:                                              NULL,
1140:                                              /*49*/ MatSetRandom_MPISELL,
1141:                                              NULL,
1142:                                              NULL,
1143:                                              NULL,
1144:                                              NULL,
1145:                                              /*54*/ MatFDColoringCreate_MPIXAIJ,
1146:                                              NULL,
1147:                                              MatSetUnfactored_MPISELL,
1148:                                              NULL,
1149:                                              NULL,
1150:                                              /*59*/ NULL,
1151:                                              MatDestroy_MPISELL,
1152:                                              MatView_MPISELL,
1153:                                              NULL,
1154:                                              NULL,
1155:                                              /*64*/ NULL,
1156:                                              NULL,
1157:                                              NULL,
1158:                                              NULL,
1159:                                              NULL,
1160:                                              /*69*/ NULL,
1161:                                              NULL,
1162:                                              NULL,
1163:                                              NULL,
1164:                                              NULL,
1165:                                              NULL,
1166:                                              /*75*/ MatFDColoringApply_AIJ, /* reuse AIJ function */
1167:                                              MatSetFromOptions_MPISELL,
1168:                                              NULL,
1169:                                              NULL,
1170:                                              NULL,
1171:                                              /*80*/ NULL,
1172:                                              NULL,
1173:                                              NULL,
1174:                                              /*83*/ NULL,
1175:                                              NULL,
1176:                                              NULL,
1177:                                              NULL,
1178:                                              NULL,
1179:                                              NULL,
1180:                                              /*89*/ NULL,
1181:                                              NULL,
1182:                                              NULL,
1183:                                              NULL,
1184:                                              NULL,
1185:                                              /*94*/ NULL,
1186:                                              NULL,
1187:                                              NULL,
1188:                                              NULL,
1189:                                              NULL,
1190:                                              /*99*/ NULL,
1191:                                              NULL,
1192:                                              NULL,
1193:                                              MatConjugate_MPISELL,
1194:                                              NULL,
1195:                                              /*104*/ NULL,
1196:                                              MatRealPart_MPISELL,
1197:                                              MatImaginaryPart_MPISELL,
1198:                                              NULL,
1199:                                              NULL,
1200:                                              /*109*/ NULL,
1201:                                              NULL,
1202:                                              NULL,
1203:                                              NULL,
1204:                                              MatMissingDiagonal_MPISELL,
1205:                                              /*114*/ NULL,
1206:                                              NULL,
1207:                                              MatGetGhosts_MPISELL,
1208:                                              NULL,
1209:                                              NULL,
1210:                                              /*119*/ MatMultDiagonalBlock_MPISELL,
1211:                                              NULL,
1212:                                              NULL,
1213:                                              NULL,
1214:                                              NULL,
1215:                                              /*124*/ NULL,
1216:                                              NULL,
1217:                                              MatInvertBlockDiagonal_MPISELL,
1218:                                              NULL,
1219:                                              NULL,
1220:                                              /*129*/ NULL,
1221:                                              NULL,
1222:                                              NULL,
1223:                                              NULL,
1224:                                              NULL,
1225:                                              /*134*/ NULL,
1226:                                              NULL,
1227:                                              NULL,
1228:                                              NULL,
1229:                                              NULL,
1230:                                              /*139*/ NULL,
1231:                                              NULL,
1232:                                              NULL,
1233:                                              MatFDColoringSetUp_MPIXAIJ,
1234:                                              NULL,
1235:                                              /*144*/ NULL,
1236:                                              NULL,
1237:                                              NULL,
1238:                                              NULL,
1239:                                              NULL,
1240:                                              NULL,
1241:                                              /*150*/ NULL,
1242:                                              NULL,
1243:                                              NULL,
1244:                                              NULL,
1245:                                              NULL,
1246:                                              /*155*/ NULL,
1247:                                              NULL};

1249: /*@C
1250:   MatMPISELLSetPreallocation - Preallocates memory for a `MATMPISELL` sparse parallel matrix in sell format.
1251:   For good matrix assembly performance the user should preallocate the matrix storage by
1252:   setting the parameters `d_nz` (or `d_nnz`) and `o_nz` (or `o_nnz`).

1254:   Collective

1256:   Input Parameters:
1257: + B     - the matrix
1258: . d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1259:            (same value is used for all local rows)
1260: . d_nnz - array containing the number of nonzeros in the various rows of the
1261:            DIAGONAL portion of the local submatrix (possibly different for each row)
1262:            or NULL (`PETSC_NULL_INTEGER` in Fortran), if `d_nz` is used to specify the nonzero structure.
1263:            The size of this array is equal to the number of local rows, i.e 'm'.
1264:            For matrices that will be factored, you must leave room for (and set)
1265:            the diagonal entry even if it is zero.
1266: . o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1267:            submatrix (same value is used for all local rows).
1268: - o_nnz - array containing the number of nonzeros in the various rows of the
1269:            OFF-DIAGONAL portion of the local submatrix (possibly different for
1270:            each row) or NULL (`PETSC_NULL_INTEGER` in Fortran), if `o_nz` is used to specify the nonzero
1271:            structure. The size of this array is equal to the number
1272:            of local rows, i.e 'm'.

1274:   Example usage:
1275:   Consider the following 8x8 matrix with 34 non-zero values, that is
1276:   assembled across 3 processors. Lets assume that proc0 owns 3 rows,
1277:   proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
1278:   as follows

1280: .vb
1281:             1  2  0  |  0  3  0  |  0  4
1282:     Proc0   0  5  6  |  7  0  0  |  8  0
1283:             9  0 10  | 11  0  0  | 12  0
1284:     -------------------------------------
1285:            13  0 14  | 15 16 17  |  0  0
1286:     Proc1   0 18  0  | 19 20 21  |  0  0
1287:             0  0  0  | 22 23  0  | 24  0
1288:     -------------------------------------
1289:     Proc2  25 26 27  |  0  0 28  | 29  0
1290:            30  0  0  | 31 32 33  |  0 34
1291: .ve

1293:   This can be represented as a collection of submatrices as

1295: .vb
1296:       A B C
1297:       D E F
1298:       G H I
1299: .ve

1301:   Where the submatrices A,B,C are owned by proc0, D,E,F are
1302:   owned by proc1, G,H,I are owned by proc2.

1304:   The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1305:   The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1306:   The 'M','N' parameters are 8,8, and have the same values on all procs.

1308:   The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
1309:   submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
1310:   corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
1311:   Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
1312:   part as `MATSEQSELL` matrices. For example, proc1 will store [E] as a `MATSEQSELL`
1313:   matrix, and [DF] as another SeqSELL matrix.

1315:   When `d_nz`, `o_nz` parameters are specified, `d_nz` storage elements are
1316:   allocated for every row of the local DIAGONAL submatrix, and o_nz
1317:   storage locations are allocated for every row of the OFF-DIAGONAL submatrix.
1318:   One way to choose `d_nz` and `o_nz` is to use the maximum number of nonzeros over
1319:   the local rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
1320:   In this case, the values of d_nz,o_nz are
1321: .vb
1322:      proc0  dnz = 2, o_nz = 2
1323:      proc1  dnz = 3, o_nz = 2
1324:      proc2  dnz = 1, o_nz = 4
1325: .ve
1326:   We are allocating m*(d_nz+o_nz) storage locations for every proc. This
1327:   translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
1328:   for proc3. i.e we are using 12+15+10=37 storage locations to store
1329:   34 values.

1331:   When `d_nnz`, `o_nnz` parameters are specified, the storage is specified
1332:   for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
1333:   In the above case the values for d_nnz,o_nnz are
1334: .vb
1335:      proc0 d_nnz = [2,2,2] and o_nnz = [2,2,2]
1336:      proc1 d_nnz = [3,3,2] and o_nnz = [2,1,1]
1337:      proc2 d_nnz = [1,1]   and o_nnz = [4,4]
1338: .ve
1339:   Here the space allocated is according to nz (or maximum values in the nnz
1340:   if nnz is provided) for DIAGONAL and OFF-DIAGONAL submatrices, i.e (2+2+3+2)*3+(1+4)*2=37

1342:   Level: intermediate

1344:   Notes:
1345:   If the *_nnz parameter is given then the *_nz parameter is ignored

1347:   The stored row and column indices begin with zero.

1349:   The parallel matrix is partitioned such that the first m0 rows belong to
1350:   process 0, the next m1 rows belong to process 1, the next m2 rows belong
1351:   to process 2 etc.. where m0,m1,m2... are the input parameter 'm'.

1353:   The DIAGONAL portion of the local submatrix of a processor can be defined
1354:   as the submatrix which is obtained by extraction the part corresponding to
1355:   the rows r1-r2 and columns c1-c2 of the global matrix, where r1 is the
1356:   first row that belongs to the processor, r2 is the last row belonging to
1357:   the this processor, and c1-c2 is range of indices of the local part of a
1358:   vector suitable for applying the matrix to.  This is an mxn matrix.  In the
1359:   common case of a square matrix, the row and column ranges are the same and
1360:   the DIAGONAL part is also square. The remaining portion of the local
1361:   submatrix (mxN) constitute the OFF-DIAGONAL portion.

1363:   If `o_nnz`, `d_nnz` are specified, then `o_nz`, and `d_nz` are ignored.

1365:   You can call `MatGetInfo()` to get information on how effective the preallocation was;
1366:   for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
1367:   You can also run with the option -info and look for messages with the string
1368:   malloc in them to see if additional memory allocation was needed.

1370: .seealso: `Mat`, `MatCreate()`, `MatCreateSeqSELL()`, `MatSetValues()`, `MatCreateSELL()`,
1371:           `MATMPISELL`, `MatGetInfo()`, `PetscSplitOwnership()`, `MATSELL`
1372: @*/
1373: PetscErrorCode MatMPISELLSetPreallocation(Mat B, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[])
1374: {
1375:   PetscFunctionBegin;
1378:   PetscTryMethod(B, "MatMPISELLSetPreallocation_C", (Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]), (B, d_nz, d_nnz, o_nz, o_nnz));
1379:   PetscFunctionReturn(PETSC_SUCCESS);
1380: }

1382: /*MC
1383:    MATMPISELL - MATMPISELL = "mpisell" - A matrix type to be used for MPI sparse matrices,
1384:    based on the sliced Ellpack format

1386:    Options Database Key:
1387: . -mat_type sell - sets the matrix type to `MATSELL` during a call to `MatSetFromOptions()`

1389:    Level: beginner

1391: .seealso: `Mat`, `MatCreateSELL()`, `MATSEQSELL`, `MATSELL`, `MATSEQAIJ`, `MATAIJ`, `MATMPIAIJ`
1392: M*/

1394: /*@C
1395:   MatCreateSELL - Creates a sparse parallel matrix in `MATSELL` format.

1397:   Collective

1399:   Input Parameters:
1400: + comm      - MPI communicator
1401: . m         - number of local rows (or `PETSC_DECIDE` to have calculated if M is given)
1402:               This value should be the same as the local size used in creating the
1403:               y vector for the matrix-vector product y = Ax.
1404: . n         - This value should be the same as the local size used in creating the
1405:               x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
1406:               calculated if `N` is given) For square matrices n is almost always `m`.
1407: . M         - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
1408: . N         - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
1409: . d_rlenmax - max number of nonzeros per row in DIAGONAL portion of local submatrix
1410:              (same value is used for all local rows)
1411: . d_rlen    - array containing the number of nonzeros in the various rows of the
1412:               DIAGONAL portion of the local submatrix (possibly different for each row)
1413:               or `NULL`, if d_rlenmax is used to specify the nonzero structure.
1414:               The size of this array is equal to the number of local rows, i.e `m`.
1415: . o_rlenmax - max number of nonzeros per row in the OFF-DIAGONAL portion of local
1416:               submatrix (same value is used for all local rows).
1417: - o_rlen    - array containing the number of nonzeros in the various rows of the
1418:               OFF-DIAGONAL portion of the local submatrix (possibly different for
1419:               each row) or `NULL`, if `o_rlenmax` is used to specify the nonzero
1420:               structure. The size of this array is equal to the number
1421:               of local rows, i.e `m`.

1423:   Output Parameter:
1424: . A - the matrix

1426:   Options Database Key:
1427: . -mat_sell_oneindex - Internally use indexing starting at 1
1428:         rather than 0.  When calling `MatSetValues()`,
1429:         the user still MUST index entries starting at 0!

1431:   Example:
1432:   Consider the following 8x8 matrix with 34 non-zero values, that is
1433:   assembled across 3 processors. Lets assume that proc0 owns 3 rows,
1434:   proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
1435:   as follows

1437: .vb
1438:             1  2  0  |  0  3  0  |  0  4
1439:     Proc0   0  5  6  |  7  0  0  |  8  0
1440:             9  0 10  | 11  0  0  | 12  0
1441:     -------------------------------------
1442:            13  0 14  | 15 16 17  |  0  0
1443:     Proc1   0 18  0  | 19 20 21  |  0  0
1444:             0  0  0  | 22 23  0  | 24  0
1445:     -------------------------------------
1446:     Proc2  25 26 27  |  0  0 28  | 29  0
1447:            30  0  0  | 31 32 33  |  0 34
1448: .ve

1450:   This can be represented as a collection of submatrices as
1451: .vb
1452:       A B C
1453:       D E F
1454:       G H I
1455: .ve

1457:   Where the submatrices A,B,C are owned by proc0, D,E,F are
1458:   owned by proc1, G,H,I are owned by proc2.

1460:   The 'm' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1461:   The 'n' parameters for proc0,proc1,proc2 are 3,3,2 respectively.
1462:   The 'M','N' parameters are 8,8, and have the same values on all procs.

1464:   The DIAGONAL submatrices corresponding to proc0,proc1,proc2 are
1465:   submatrices [A], [E], [I] respectively. The OFF-DIAGONAL submatrices
1466:   corresponding to proc0,proc1,proc2 are [BC], [DF], [GH] respectively.
1467:   Internally, each processor stores the DIAGONAL part, and the OFF-DIAGONAL
1468:   part as `MATSEQSELL` matrices. For example, proc1 will store [E] as a `MATSEQSELL`
1469:   matrix, and [DF] as another `MATSEQSELL` matrix.

1471:   When d_rlenmax, o_rlenmax parameters are specified, d_rlenmax storage elements are
1472:   allocated for every row of the local DIAGONAL submatrix, and o_rlenmax
1473:   storage locations are allocated for every row of the OFF-DIAGONAL submatrix.
1474:   One way to choose `d_rlenmax` and `o_rlenmax` is to use the maximum number of nonzeros over
1475:   the local rows for each of the local DIAGONAL, and the OFF-DIAGONAL submatrices.
1476:   In this case, the values of d_rlenmax,o_rlenmax are
1477: .vb
1478:      proc0 - d_rlenmax = 2, o_rlenmax = 2
1479:      proc1 - d_rlenmax = 3, o_rlenmax = 2
1480:      proc2 - d_rlenmax = 1, o_rlenmax = 4
1481: .ve
1482:   We are allocating m*(d_rlenmax+o_rlenmax) storage locations for every proc. This
1483:   translates to 3*(2+2)=12 for proc0, 3*(3+2)=15 for proc1, 2*(1+4)=10
1484:   for proc3. i.e we are using 12+15+10=37 storage locations to store
1485:   34 values.

1487:   When `d_rlen`, `o_rlen` parameters are specified, the storage is specified
1488:   for every row, corresponding to both DIAGONAL and OFF-DIAGONAL submatrices.
1489:   In the above case the values for `d_nnz`, `o_nnz` are
1490: .vb
1491:      proc0 - d_nnz = [2,2,2] and o_nnz = [2,2,2]
1492:      proc1 - d_nnz = [3,3,2] and o_nnz = [2,1,1]
1493:      proc2 - d_nnz = [1,1]   and o_nnz = [4,4]
1494: .ve
1495:   Here the space allocated is still 37 though there are 34 nonzeros because
1496:   the allocation is always done according to rlenmax.

1498:   Level: intermediate

1500:   Notes:
1501:   It is recommended that one use the `MatCreate()`, `MatSetType()` and/or `MatSetFromOptions()`,
1502:   MatXXXXSetPreallocation() paradigm instead of this routine directly.
1503:   [MatXXXXSetPreallocation() is, for example, `MatSeqSELLSetPreallocation()`]

1505:   If the *_rlen parameter is given then the *_rlenmax parameter is ignored

1507:   `m`, `n`, `M`, `N` parameters specify the size of the matrix, and its partitioning across
1508:   processors, while `d_rlenmax`, `d_rlen`, `o_rlenmax` , `o_rlen` parameters specify the approximate
1509:   storage requirements for this matrix.

1511:   If `PETSC_DECIDE` or  `PETSC_DETERMINE` is used for a particular argument on one
1512:   processor than it must be used on all processors that share the object for
1513:   that argument.

1515:   The user MUST specify either the local or global matrix dimensions
1516:   (possibly both).

1518:   The parallel matrix is partitioned across processors such that the
1519:   first m0 rows belong to process 0, the next m1 rows belong to
1520:   process 1, the next m2 rows belong to process 2 etc.. where
1521:   m0,m1,m2,.. are the input parameter 'm'. i.e each processor stores
1522:   values corresponding to [`m` x `N`] submatrix.

1524:   The columns are logically partitioned with the n0 columns belonging
1525:   to 0th partition, the next n1 columns belonging to the next
1526:   partition etc.. where n0,n1,n2... are the input parameter `n`.

1528:   The DIAGONAL portion of the local submatrix on any given processor
1529:   is the submatrix corresponding to the rows and columns `m`, `n`
1530:   corresponding to the given processor. i.e diagonal matrix on
1531:   process 0 is [m0 x n0], diagonal matrix on process 1 is [m1 x n1]
1532:   etc. The remaining portion of the local submatrix [m x (N-n)]
1533:   constitute the OFF-DIAGONAL portion. The example below better
1534:   illustrates this concept.

1536:   For a square global matrix we define each processor's diagonal portion
1537:   to be its local rows and the corresponding columns (a square submatrix);
1538:   each processor's off-diagonal portion encompasses the remainder of the
1539:   local matrix (a rectangular submatrix).

1541:   If `o_rlen`, `d_rlen` are specified, then `o_rlenmax`, and `d_rlenmax` are ignored.

1543:   When calling this routine with a single process communicator, a matrix of
1544:   type `MATSEQSELL` is returned.  If a matrix of type `MATMPISELL` is desired for this
1545:   type of communicator, use the construction mechanism
1546: .vb
1547:    MatCreate(...,&A);
1548:    MatSetType(A,MATMPISELL);
1549:    MatSetSizes(A, m,n,M,N);
1550:    MatMPISELLSetPreallocation(A,...);
1551: .ve

1553: .seealso: `Mat`, `MATSELL`, `MatCreate()`, `MatCreateSeqSELL()`, `MatSetValues()`, `MatMPISELLSetPreallocation()`, `MATMPISELL`
1554: @*/
1555: 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)
1556: {
1557:   PetscMPIInt size;

1559:   PetscFunctionBegin;
1560:   PetscCall(MatCreate(comm, A));
1561:   PetscCall(MatSetSizes(*A, m, n, M, N));
1562:   PetscCallMPI(MPI_Comm_size(comm, &size));
1563:   if (size > 1) {
1564:     PetscCall(MatSetType(*A, MATMPISELL));
1565:     PetscCall(MatMPISELLSetPreallocation(*A, d_rlenmax, d_rlen, o_rlenmax, o_rlen));
1566:   } else {
1567:     PetscCall(MatSetType(*A, MATSEQSELL));
1568:     PetscCall(MatSeqSELLSetPreallocation(*A, d_rlenmax, d_rlen));
1569:   }
1570:   PetscFunctionReturn(PETSC_SUCCESS);
1571: }

1573: /*@C
1574:   MatMPISELLGetSeqSELL - Returns the local pieces of this distributed matrix

1576:   Not Collective

1578:   Input Parameter:
1579: . A - the `MATMPISELL` matrix

1581:   Output Parameters:
1582: + Ad     - The diagonal portion of `A`
1583: . Ao     - The off-diagonal portion of `A`
1584: - colmap - An array mapping local column numbers of `Ao` to global column numbers of the parallel matrix

1586:   Level: advanced

1588: .seealso: `Mat`, `MATSEQSELL`, `MATMPISELL`
1589: @*/
1590: PetscErrorCode MatMPISELLGetSeqSELL(Mat A, Mat *Ad, Mat *Ao, const PetscInt *colmap[])
1591: {
1592:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1593:   PetscBool    flg;

1595:   PetscFunctionBegin;
1596:   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISELL, &flg));
1597:   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "This function requires a MATMPISELL matrix as input");
1598:   if (Ad) *Ad = a->A;
1599:   if (Ao) *Ao = a->B;
1600:   if (colmap) *colmap = a->garray;
1601:   PetscFunctionReturn(PETSC_SUCCESS);
1602: }

1604: /*@C
1605:   MatMPISELLGetLocalMatCondensed - Creates a `MATSEQSELL` matrix from an `MATMPISELL` matrix by
1606:   taking all its local rows and NON-ZERO columns

1608:   Not Collective

1610:   Input Parameters:
1611: + A     - the matrix
1612: . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
1613: . row   - index sets of rows to extract (or `NULL`)
1614: - col   - index sets of columns to extract (or `NULL`)

1616:   Output Parameter:
1617: . A_loc - the local sequential matrix generated

1619:   Level: advanced

1621: .seealso: `Mat`, `MATSEQSELL`, `MATMPISELL`, `MatGetOwnershipRange()`, `MatMPISELLGetLocalMat()`
1622: @*/
1623: PetscErrorCode MatMPISELLGetLocalMatCondensed(Mat A, MatReuse scall, IS *row, IS *col, Mat *A_loc)
1624: {
1625:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1626:   PetscInt     i, start, end, ncols, nzA, nzB, *cmap, imark, *idx;
1627:   IS           isrowa, iscola;
1628:   Mat         *aloc;
1629:   PetscBool    match;

1631:   PetscFunctionBegin;
1632:   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISELL, &match));
1633:   PetscCheck(match, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Requires MATMPISELL matrix as input");
1634:   PetscCall(PetscLogEventBegin(MAT_Getlocalmatcondensed, A, 0, 0, 0));
1635:   if (!row) {
1636:     start = A->rmap->rstart;
1637:     end   = A->rmap->rend;
1638:     PetscCall(ISCreateStride(PETSC_COMM_SELF, end - start, start, 1, &isrowa));
1639:   } else {
1640:     isrowa = *row;
1641:   }
1642:   if (!col) {
1643:     start = A->cmap->rstart;
1644:     cmap  = a->garray;
1645:     nzA   = a->A->cmap->n;
1646:     nzB   = a->B->cmap->n;
1647:     PetscCall(PetscMalloc1(nzA + nzB, &idx));
1648:     ncols = 0;
1649:     for (i = 0; i < nzB; i++) {
1650:       if (cmap[i] < start) idx[ncols++] = cmap[i];
1651:       else break;
1652:     }
1653:     imark = i;
1654:     for (i = 0; i < nzA; i++) idx[ncols++] = start + i;
1655:     for (i = imark; i < nzB; i++) idx[ncols++] = cmap[i];
1656:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncols, idx, PETSC_OWN_POINTER, &iscola));
1657:   } else {
1658:     iscola = *col;
1659:   }
1660:   if (scall != MAT_INITIAL_MATRIX) {
1661:     PetscCall(PetscMalloc1(1, &aloc));
1662:     aloc[0] = *A_loc;
1663:   }
1664:   PetscCall(MatCreateSubMatrices(A, 1, &isrowa, &iscola, scall, &aloc));
1665:   *A_loc = aloc[0];
1666:   PetscCall(PetscFree(aloc));
1667:   if (!row) PetscCall(ISDestroy(&isrowa));
1668:   if (!col) PetscCall(ISDestroy(&iscola));
1669:   PetscCall(PetscLogEventEnd(MAT_Getlocalmatcondensed, A, 0, 0, 0));
1670:   PetscFunctionReturn(PETSC_SUCCESS);
1671: }

1673: #include <../src/mat/impls/aij/mpi/mpiaij.h>

1675: PetscErrorCode MatConvert_MPISELL_MPIAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1676: {
1677:   Mat_MPISELL *a = (Mat_MPISELL *)A->data;
1678:   Mat          B;
1679:   Mat_MPIAIJ  *b;

1681:   PetscFunctionBegin;
1682:   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Matrix must be assembled");

1684:   if (reuse == MAT_REUSE_MATRIX) {
1685:     B = *newmat;
1686:   } else {
1687:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
1688:     PetscCall(MatSetType(B, MATMPIAIJ));
1689:     PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
1690:     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
1691:     PetscCall(MatSeqAIJSetPreallocation(B, 0, NULL));
1692:     PetscCall(MatMPIAIJSetPreallocation(B, 0, NULL, 0, NULL));
1693:   }
1694:   b = (Mat_MPIAIJ *)B->data;

1696:   if (reuse == MAT_REUSE_MATRIX) {
1697:     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_REUSE_MATRIX, &b->A));
1698:     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_REUSE_MATRIX, &b->B));
1699:   } else {
1700:     PetscCall(MatDestroy(&b->A));
1701:     PetscCall(MatDestroy(&b->B));
1702:     PetscCall(MatDisAssemble_MPISELL(A));
1703:     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->A, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->A));
1704:     PetscCall(MatConvert_SeqSELL_SeqAIJ(a->B, MATSEQAIJ, MAT_INITIAL_MATRIX, &b->B));
1705:     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1706:     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1707:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1708:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1709:   }

1711:   if (reuse == MAT_INPLACE_MATRIX) {
1712:     PetscCall(MatHeaderReplace(A, &B));
1713:   } else {
1714:     *newmat = B;
1715:   }
1716:   PetscFunctionReturn(PETSC_SUCCESS);
1717: }

1719: PetscErrorCode MatConvert_MPIAIJ_MPISELL(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1720: {
1721:   Mat_MPIAIJ  *a = (Mat_MPIAIJ *)A->data;
1722:   Mat          B;
1723:   Mat_MPISELL *b;

1725:   PetscFunctionBegin;
1726:   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Matrix must be assembled");

1728:   if (reuse == MAT_REUSE_MATRIX) {
1729:     B = *newmat;
1730:   } else {
1731:     Mat_SeqAIJ *Aa = (Mat_SeqAIJ *)a->A->data, *Ba = (Mat_SeqAIJ *)a->B->data;
1732:     PetscInt    i, d_nz = 0, o_nz = 0, m = A->rmap->N, n = A->cmap->N, lm = A->rmap->n, ln = A->cmap->n;
1733:     PetscInt   *d_nnz, *o_nnz;
1734:     PetscCall(PetscMalloc2(lm, &d_nnz, lm, &o_nnz));
1735:     for (i = 0; i < lm; i++) {
1736:       d_nnz[i] = Aa->i[i + 1] - Aa->i[i];
1737:       o_nnz[i] = Ba->i[i + 1] - Ba->i[i];
1738:       if (d_nnz[i] > d_nz) d_nz = d_nnz[i];
1739:       if (o_nnz[i] > o_nz) o_nz = o_nnz[i];
1740:     }
1741:     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
1742:     PetscCall(MatSetType(B, MATMPISELL));
1743:     PetscCall(MatSetSizes(B, lm, ln, m, n));
1744:     PetscCall(MatSetBlockSizes(B, A->rmap->bs, A->cmap->bs));
1745:     PetscCall(MatSeqSELLSetPreallocation(B, d_nz, d_nnz));
1746:     PetscCall(MatMPISELLSetPreallocation(B, d_nz, d_nnz, o_nz, o_nnz));
1747:     PetscCall(PetscFree2(d_nnz, o_nnz));
1748:   }
1749:   b = (Mat_MPISELL *)B->data;

1751:   if (reuse == MAT_REUSE_MATRIX) {
1752:     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_REUSE_MATRIX, &b->A));
1753:     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_REUSE_MATRIX, &b->B));
1754:   } else {
1755:     PetscCall(MatDestroy(&b->A));
1756:     PetscCall(MatDestroy(&b->B));
1757:     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->A, MATSEQSELL, MAT_INITIAL_MATRIX, &b->A));
1758:     PetscCall(MatConvert_SeqAIJ_SeqSELL(a->B, MATSEQSELL, MAT_INITIAL_MATRIX, &b->B));
1759:     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1760:     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1761:     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1762:     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1763:   }

1765:   if (reuse == MAT_INPLACE_MATRIX) {
1766:     PetscCall(MatHeaderReplace(A, &B));
1767:   } else {
1768:     *newmat = B;
1769:   }
1770:   PetscFunctionReturn(PETSC_SUCCESS);
1771: }

1773: PetscErrorCode MatSOR_MPISELL(Mat matin, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx)
1774: {
1775:   Mat_MPISELL *mat = (Mat_MPISELL *)matin->data;
1776:   Vec          bb1 = NULL;

1778:   PetscFunctionBegin;
1779:   if (flag == SOR_APPLY_UPPER) {
1780:     PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1781:     PetscFunctionReturn(PETSC_SUCCESS);
1782:   }

1784:   if (its > 1 || ~flag & SOR_ZERO_INITIAL_GUESS || flag & SOR_EISENSTAT) PetscCall(VecDuplicate(bb, &bb1));

1786:   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
1787:     if (flag & SOR_ZERO_INITIAL_GUESS) {
1788:       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1789:       its--;
1790:     }

1792:     while (its--) {
1793:       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1794:       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));

1796:       /* update rhs: bb1 = bb - B*x */
1797:       PetscCall(VecScale(mat->lvec, -1.0));
1798:       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));

1800:       /* local sweep */
1801:       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_SYMMETRIC_SWEEP, fshift, lits, 1, xx));
1802:     }
1803:   } else if (flag & SOR_LOCAL_FORWARD_SWEEP) {
1804:     if (flag & SOR_ZERO_INITIAL_GUESS) {
1805:       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1806:       its--;
1807:     }
1808:     while (its--) {
1809:       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1810:       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));

1812:       /* update rhs: bb1 = bb - B*x */
1813:       PetscCall(VecScale(mat->lvec, -1.0));
1814:       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));

1816:       /* local sweep */
1817:       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_FORWARD_SWEEP, fshift, lits, 1, xx));
1818:     }
1819:   } else if (flag & SOR_LOCAL_BACKWARD_SWEEP) {
1820:     if (flag & SOR_ZERO_INITIAL_GUESS) {
1821:       PetscCall((*mat->A->ops->sor)(mat->A, bb, omega, flag, fshift, lits, 1, xx));
1822:       its--;
1823:     }
1824:     while (its--) {
1825:       PetscCall(VecScatterBegin(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));
1826:       PetscCall(VecScatterEnd(mat->Mvctx, xx, mat->lvec, INSERT_VALUES, SCATTER_FORWARD));

1828:       /* update rhs: bb1 = bb - B*x */
1829:       PetscCall(VecScale(mat->lvec, -1.0));
1830:       PetscCall((*mat->B->ops->multadd)(mat->B, mat->lvec, bb, bb1));

1832:       /* local sweep */
1833:       PetscCall((*mat->A->ops->sor)(mat->A, bb1, omega, SOR_BACKWARD_SWEEP, fshift, lits, 1, xx));
1834:     }
1835:   } else SETERRQ(PetscObjectComm((PetscObject)matin), PETSC_ERR_SUP, "Parallel SOR not supported");

1837:   PetscCall(VecDestroy(&bb1));

1839:   matin->factorerrortype = mat->A->factorerrortype;
1840:   PetscFunctionReturn(PETSC_SUCCESS);
1841: }

1843: #if defined(PETSC_HAVE_CUDA)
1844: PETSC_INTERN PetscErrorCode MatConvert_MPISELL_MPISELLCUDA(Mat, MatType, MatReuse, Mat *);
1845: #endif

1847: /*MC
1848:    MATMPISELL - MATMPISELL = "MPISELL" - A matrix type to be used for parallel sparse matrices.

1850:    Options Database Keys:
1851: . -mat_type mpisell - sets the matrix type to `MATMPISELL` during a call to `MatSetFromOptions()`

1853:   Level: beginner

1855: .seealso: `Mat`, `MATSELL`, `MATSEQSELL` `MatCreateSELL()`
1856: M*/
1857: PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat B)
1858: {
1859:   Mat_MPISELL *b;
1860:   PetscMPIInt  size;

1862:   PetscFunctionBegin;
1863:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B), &size));
1864:   PetscCall(PetscNew(&b));
1865:   B->data       = (void *)b;
1866:   B->ops[0]     = MatOps_Values;
1867:   B->assembled  = PETSC_FALSE;
1868:   B->insertmode = NOT_SET_VALUES;
1869:   b->size       = size;
1870:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)B), &b->rank));
1871:   /* build cache for off array entries formed */
1872:   PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)B), 1, &B->stash));

1874:   b->donotstash  = PETSC_FALSE;
1875:   b->colmap      = NULL;
1876:   b->garray      = NULL;
1877:   b->roworiented = PETSC_TRUE;

1879:   /* stuff used for matrix vector multiply */
1880:   b->lvec  = NULL;
1881:   b->Mvctx = NULL;

1883:   /* stuff for MatGetRow() */
1884:   b->rowindices   = NULL;
1885:   b->rowvalues    = NULL;
1886:   b->getrowactive = PETSC_FALSE;

1888:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatStoreValues_C", MatStoreValues_MPISELL));
1889:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatRetrieveValues_C", MatRetrieveValues_MPISELL));
1890:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatIsTranspose_C", MatIsTranspose_MPISELL));
1891:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPISELLSetPreallocation_C", MatMPISELLSetPreallocation_MPISELL));
1892:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisell_mpiaij_C", MatConvert_MPISELL_MPIAIJ));
1893: #if defined(PETSC_HAVE_CUDA)
1894:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpisell_mpisellcuda_C", MatConvert_MPISELL_MPISELLCUDA));
1895: #endif
1896:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatDiagonalScaleLocal_C", MatDiagonalScaleLocal_MPISELL));
1897:   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPISELL));
1898:   PetscFunctionReturn(PETSC_SUCCESS);
1899: }