Actual source code: matrix.c

  1: /*
  2:    This is where the abstract matrix operations are defined
  3:    Portions of this code are under:
  4:    Copyright (c) 2022 Advanced Micro Devices, Inc. All rights reserved.
  5: */

  7: #include <petsc/private/matimpl.h>
  8: #include <petsc/private/isimpl.h>
  9: #include <petsc/private/vecimpl.h>

 11: /* Logging support */
 12: PetscClassId MAT_CLASSID;
 13: PetscClassId MAT_COLORING_CLASSID;
 14: PetscClassId MAT_FDCOLORING_CLASSID;
 15: PetscClassId MAT_TRANSPOSECOLORING_CLASSID;

 17: PetscLogEvent MAT_Mult, MAT_MultAdd, MAT_MultTranspose;
 18: PetscLogEvent MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose, MAT_MatSolve, MAT_MatTrSolve;
 19: PetscLogEvent MAT_SolveTransposeAdd, MAT_SOR, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
 20: PetscLogEvent MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
 21: PetscLogEvent MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
 22: PetscLogEvent MAT_QRFactorNumeric, MAT_QRFactorSymbolic, MAT_QRFactor;
 23: PetscLogEvent MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetRowIJ, MAT_CreateSubMats, MAT_GetOrdering, MAT_RedundantMat, MAT_GetSeqNonzeroStructure;
 24: PetscLogEvent MAT_IncreaseOverlap, MAT_Partitioning, MAT_PartitioningND, MAT_Coarsen, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate;
 25: PetscLogEvent MAT_FDColoringSetUp, MAT_FDColoringApply, MAT_Transpose, MAT_FDColoringFunction, MAT_CreateSubMat;
 26: PetscLogEvent MAT_TransposeColoringCreate;
 27: PetscLogEvent MAT_MatMult, MAT_MatMultSymbolic, MAT_MatMultNumeric;
 28: PetscLogEvent MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric, MAT_RARt, MAT_RARtSymbolic, MAT_RARtNumeric;
 29: PetscLogEvent MAT_MatTransposeMult, MAT_MatTransposeMultSymbolic, MAT_MatTransposeMultNumeric;
 30: PetscLogEvent MAT_TransposeMatMult, MAT_TransposeMatMultSymbolic, MAT_TransposeMatMultNumeric;
 31: PetscLogEvent MAT_MatMatMult, MAT_MatMatMultSymbolic, MAT_MatMatMultNumeric;
 32: PetscLogEvent MAT_MultHermitianTranspose, MAT_MultHermitianTransposeAdd;
 33: PetscLogEvent MAT_Getsymtransreduced, MAT_GetBrowsOfAcols;
 34: PetscLogEvent MAT_GetBrowsOfAocols, MAT_Getlocalmat, MAT_Getlocalmatcondensed, MAT_Seqstompi, MAT_Seqstompinum, MAT_Seqstompisym;
 35: PetscLogEvent MAT_GetMultiProcBlock;
 36: PetscLogEvent MAT_CUSPARSECopyToGPU, MAT_CUSPARSECopyFromGPU, MAT_CUSPARSEGenerateTranspose, MAT_CUSPARSESolveAnalysis;
 37: PetscLogEvent MAT_HIPSPARSECopyToGPU, MAT_HIPSPARSECopyFromGPU, MAT_HIPSPARSEGenerateTranspose, MAT_HIPSPARSESolveAnalysis;
 38: PetscLogEvent MAT_PreallCOO, MAT_SetVCOO;
 39: PetscLogEvent MAT_CreateGraph;
 40: PetscLogEvent MAT_SetValuesBatch;
 41: PetscLogEvent MAT_ViennaCLCopyToGPU;
 42: PetscLogEvent MAT_CUDACopyToGPU, MAT_HIPCopyToGPU;
 43: PetscLogEvent MAT_DenseCopyToGPU, MAT_DenseCopyFromGPU;
 44: PetscLogEvent MAT_Merge, MAT_Residual, MAT_SetRandom;
 45: PetscLogEvent MAT_FactorFactS, MAT_FactorInvS;
 46: PetscLogEvent MATCOLORING_Apply, MATCOLORING_Comm, MATCOLORING_Local, MATCOLORING_ISCreate, MATCOLORING_SetUp, MATCOLORING_Weights;
 47: PetscLogEvent MAT_H2Opus_Build, MAT_H2Opus_Compress, MAT_H2Opus_Orthog, MAT_H2Opus_LR;

 49: const char *const MatFactorTypes[] = {"NONE", "LU", "CHOLESKY", "ILU", "ICC", "ILUDT", "QR", "MatFactorType", "MAT_FACTOR_", NULL};

 51: /*@
 52:   MatSetRandom - Sets all components of a matrix to random numbers.

 54:   Logically Collective

 56:   Input Parameters:
 57: + x    - the matrix
 58: - rctx - the `PetscRandom` object, formed by `PetscRandomCreate()`, or `NULL` and
 59:           it will create one internally.

 61:   Example:
 62: .vb
 63:      PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 64:      MatSetRandom(x,rctx);
 65:      PetscRandomDestroy(rctx);
 66: .ve

 68:   Level: intermediate

 70:   Notes:
 71:   For sparse matrices that have been preallocated but not been assembled, it randomly selects appropriate locations,

 73:   for sparse matrices that already have nonzero locations, it fills the locations with random numbers.

 75:   It generates an error if used on unassembled sparse matrices that have not been preallocated.

 77: .seealso: [](ch_matrices), `Mat`, `PetscRandom`, `PetscRandomCreate()`, `MatZeroEntries()`, `MatSetValues()`, `PetscRandomDestroy()`
 78: @*/
 79: PetscErrorCode MatSetRandom(Mat x, PetscRandom rctx)
 80: {
 81:   PetscRandom randObj = NULL;

 83:   PetscFunctionBegin;
 87:   MatCheckPreallocated(x, 1);

 89:   if (!rctx) {
 90:     MPI_Comm comm;
 91:     PetscCall(PetscObjectGetComm((PetscObject)x, &comm));
 92:     PetscCall(PetscRandomCreate(comm, &randObj));
 93:     PetscCall(PetscRandomSetType(randObj, x->defaultrandtype));
 94:     PetscCall(PetscRandomSetFromOptions(randObj));
 95:     rctx = randObj;
 96:   }
 97:   PetscCall(PetscLogEventBegin(MAT_SetRandom, x, rctx, 0, 0));
 98:   PetscUseTypeMethod(x, setrandom, rctx);
 99:   PetscCall(PetscLogEventEnd(MAT_SetRandom, x, rctx, 0, 0));

101:   PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY));
102:   PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY));
103:   PetscCall(PetscRandomDestroy(&randObj));
104:   PetscFunctionReturn(PETSC_SUCCESS);
105: }

107: /*@
108:   MatFactorGetErrorZeroPivot - returns the pivot value that was determined to be zero and the row it occurred in

110:   Logically Collective

112:   Input Parameter:
113: . mat - the factored matrix

115:   Output Parameters:
116: + pivot - the pivot value computed
117: - row   - the row that the zero pivot occurred. This row value must be interpreted carefully due to row reorderings and which processes
118:          the share the matrix

120:   Level: advanced

122:   Notes:
123:   This routine does not work for factorizations done with external packages.

125:   This routine should only be called if `MatGetFactorError()` returns a value of `MAT_FACTOR_NUMERIC_ZEROPIVOT`

127:   This can also be called on non-factored matrices that come from, for example, matrices used in SOR.

129: .seealso: [](ch_matrices), `Mat`, `MatZeroEntries()`, `MatFactor()`, `MatGetFactor()`,
130: `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`, `MatFactorClearError()`,
131: `MAT_FACTOR_NUMERIC_ZEROPIVOT`
132: @*/
133: PetscErrorCode MatFactorGetErrorZeroPivot(Mat mat, PetscReal *pivot, PetscInt *row)
134: {
135:   PetscFunctionBegin;
137:   PetscAssertPointer(pivot, 2);
138:   PetscAssertPointer(row, 3);
139:   *pivot = mat->factorerror_zeropivot_value;
140:   *row   = mat->factorerror_zeropivot_row;
141:   PetscFunctionReturn(PETSC_SUCCESS);
142: }

144: /*@
145:   MatFactorGetError - gets the error code from a factorization

147:   Logically Collective

149:   Input Parameter:
150: . mat - the factored matrix

152:   Output Parameter:
153: . err - the error code

155:   Level: advanced

157:   Note:
158:   This can also be called on non-factored matrices that come from, for example, matrices used in SOR.

160: .seealso: [](ch_matrices), `Mat`, `MatZeroEntries()`, `MatFactor()`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`,
161:           `MatFactorClearError()`, `MatFactorGetErrorZeroPivot()`, `MatFactorError`
162: @*/
163: PetscErrorCode MatFactorGetError(Mat mat, MatFactorError *err)
164: {
165:   PetscFunctionBegin;
167:   PetscAssertPointer(err, 2);
168:   *err = mat->factorerrortype;
169:   PetscFunctionReturn(PETSC_SUCCESS);
170: }

172: /*@
173:   MatFactorClearError - clears the error code in a factorization

175:   Logically Collective

177:   Input Parameter:
178: . mat - the factored matrix

180:   Level: developer

182:   Note:
183:   This can also be called on non-factored matrices that come from, for example, matrices used in SOR.

185: .seealso: [](ch_matrices), `Mat`, `MatZeroEntries()`, `MatFactor()`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`, `MatFactorGetError()`, `MatFactorGetErrorZeroPivot()`,
186:           `MatGetErrorCode()`, `MatFactorError`
187: @*/
188: PetscErrorCode MatFactorClearError(Mat mat)
189: {
190:   PetscFunctionBegin;
192:   mat->factorerrortype             = MAT_FACTOR_NOERROR;
193:   mat->factorerror_zeropivot_value = 0.0;
194:   mat->factorerror_zeropivot_row   = 0;
195:   PetscFunctionReturn(PETSC_SUCCESS);
196: }

198: PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat mat, PetscBool cols, PetscReal tol, IS *nonzero)
199: {
200:   Vec                r, l;
201:   const PetscScalar *al;
202:   PetscInt           i, nz, gnz, N, n, st;

204:   PetscFunctionBegin;
205:   PetscCall(MatCreateVecs(mat, &r, &l));
206:   if (!cols) { /* nonzero rows */
207:     PetscCall(MatGetOwnershipRange(mat, &st, NULL));
208:     PetscCall(MatGetSize(mat, &N, NULL));
209:     PetscCall(MatGetLocalSize(mat, &n, NULL));
210:     PetscCall(VecSet(l, 0.0));
211:     PetscCall(VecSetRandom(r, NULL));
212:     PetscCall(MatMult(mat, r, l));
213:     PetscCall(VecGetArrayRead(l, &al));
214:   } else { /* nonzero columns */
215:     PetscCall(MatGetOwnershipRangeColumn(mat, &st, NULL));
216:     PetscCall(MatGetSize(mat, NULL, &N));
217:     PetscCall(MatGetLocalSize(mat, NULL, &n));
218:     PetscCall(VecSet(r, 0.0));
219:     PetscCall(VecSetRandom(l, NULL));
220:     PetscCall(MatMultTranspose(mat, l, r));
221:     PetscCall(VecGetArrayRead(r, &al));
222:   }
223:   if (tol <= 0.0) {
224:     for (i = 0, nz = 0; i < n; i++)
225:       if (al[i] != 0.0) nz++;
226:   } else {
227:     for (i = 0, nz = 0; i < n; i++)
228:       if (PetscAbsScalar(al[i]) > tol) nz++;
229:   }
230:   PetscCallMPI(MPIU_Allreduce(&nz, &gnz, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)mat)));
231:   if (gnz != N) {
232:     PetscInt *nzr;
233:     PetscCall(PetscMalloc1(nz, &nzr));
234:     if (nz) {
235:       if (tol < 0) {
236:         for (i = 0, nz = 0; i < n; i++)
237:           if (al[i] != 0.0) nzr[nz++] = i + st;
238:       } else {
239:         for (i = 0, nz = 0; i < n; i++)
240:           if (PetscAbsScalar(al[i]) > tol) nzr[nz++] = i + st;
241:       }
242:     }
243:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)mat), nz, nzr, PETSC_OWN_POINTER, nonzero));
244:   } else *nonzero = NULL;
245:   if (!cols) { /* nonzero rows */
246:     PetscCall(VecRestoreArrayRead(l, &al));
247:   } else {
248:     PetscCall(VecRestoreArrayRead(r, &al));
249:   }
250:   PetscCall(VecDestroy(&l));
251:   PetscCall(VecDestroy(&r));
252:   PetscFunctionReturn(PETSC_SUCCESS);
253: }

255: /*@
256:   MatFindNonzeroRows - Locate all rows that are not completely zero in the matrix

258:   Input Parameter:
259: . mat - the matrix

261:   Output Parameter:
262: . keptrows - the rows that are not completely zero

264:   Level: intermediate

266:   Note:
267:   `keptrows` is set to `NULL` if all rows are nonzero.

269:   Developer Note:
270:   If `keptrows` is not `NULL`, it must be sorted.

272: .seealso: [](ch_matrices), `Mat`, `MatFindZeroRows()`
273:  @*/
274: PetscErrorCode MatFindNonzeroRows(Mat mat, IS *keptrows)
275: {
276:   PetscFunctionBegin;
279:   PetscAssertPointer(keptrows, 2);
280:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
281:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
282:   if (mat->ops->findnonzerorows) PetscUseTypeMethod(mat, findnonzerorows, keptrows);
283:   else PetscCall(MatFindNonzeroRowsOrCols_Basic(mat, PETSC_FALSE, 0.0, keptrows));
284:   if (keptrows && *keptrows) PetscCall(ISSetInfo(*keptrows, IS_SORTED, IS_GLOBAL, PETSC_FALSE, PETSC_TRUE));
285:   PetscFunctionReturn(PETSC_SUCCESS);
286: }

288: /*@
289:   MatFindZeroRows - Locate all rows that are completely zero in the matrix

291:   Input Parameter:
292: . mat - the matrix

294:   Output Parameter:
295: . zerorows - the rows that are completely zero

297:   Level: intermediate

299:   Note:
300:   `zerorows` is set to `NULL` if no rows are zero.

302: .seealso: [](ch_matrices), `Mat`, `MatFindNonzeroRows()`
303:  @*/
304: PetscErrorCode MatFindZeroRows(Mat mat, IS *zerorows)
305: {
306:   IS       keptrows;
307:   PetscInt m, n;

309:   PetscFunctionBegin;
312:   PetscAssertPointer(zerorows, 2);
313:   PetscCall(MatFindNonzeroRows(mat, &keptrows));
314:   /* MatFindNonzeroRows sets keptrows to NULL if there are no zero rows.
315:      In keeping with this convention, we set zerorows to NULL if there are no zero
316:      rows. */
317:   if (keptrows == NULL) {
318:     *zerorows = NULL;
319:   } else {
320:     PetscCall(MatGetOwnershipRange(mat, &m, &n));
321:     PetscCall(ISComplement(keptrows, m, n, zerorows));
322:     PetscCall(ISDestroy(&keptrows));
323:   }
324:   PetscFunctionReturn(PETSC_SUCCESS);
325: }

327: /*@
328:   MatGetDiagonalBlock - Returns the part of the matrix associated with the on-process coupling

330:   Not Collective

332:   Input Parameter:
333: . A - the matrix

335:   Output Parameter:
336: . a - the diagonal part (which is a SEQUENTIAL matrix)

338:   Level: advanced

340:   Notes:
341:   See `MatCreateAIJ()` for more information on the "diagonal part" of the matrix.

343:   Use caution, as the reference count on the returned matrix is not incremented and it is used as part of `A`'s normal operation.

345: .seealso: [](ch_matrices), `Mat`, `MatCreateAIJ()`, `MATAIJ`, `MATBAIJ`, `MATSBAIJ`
346: @*/
347: PetscErrorCode MatGetDiagonalBlock(Mat A, Mat *a)
348: {
349:   PetscFunctionBegin;
352:   PetscAssertPointer(a, 2);
353:   PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
354:   if (A->ops->getdiagonalblock) PetscUseTypeMethod(A, getdiagonalblock, a);
355:   else {
356:     PetscMPIInt size;

358:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
359:     PetscCheck(size == 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for parallel matrix type %s", ((PetscObject)A)->type_name);
360:     *a = A;
361:   }
362:   PetscFunctionReturn(PETSC_SUCCESS);
363: }

365: /*@
366:   MatGetTrace - Gets the trace of a matrix. The sum of the diagonal entries.

368:   Collective

370:   Input Parameter:
371: . mat - the matrix

373:   Output Parameter:
374: . trace - the sum of the diagonal entries

376:   Level: advanced

378: .seealso: [](ch_matrices), `Mat`
379: @*/
380: PetscErrorCode MatGetTrace(Mat mat, PetscScalar *trace)
381: {
382:   Vec diag;

384:   PetscFunctionBegin;
386:   PetscAssertPointer(trace, 2);
387:   PetscCall(MatCreateVecs(mat, &diag, NULL));
388:   PetscCall(MatGetDiagonal(mat, diag));
389:   PetscCall(VecSum(diag, trace));
390:   PetscCall(VecDestroy(&diag));
391:   PetscFunctionReturn(PETSC_SUCCESS);
392: }

394: /*@
395:   MatRealPart - Zeros out the imaginary part of the matrix

397:   Logically Collective

399:   Input Parameter:
400: . mat - the matrix

402:   Level: advanced

404: .seealso: [](ch_matrices), `Mat`, `MatImaginaryPart()`
405: @*/
406: PetscErrorCode MatRealPart(Mat mat)
407: {
408:   PetscFunctionBegin;
411:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
412:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
413:   MatCheckPreallocated(mat, 1);
414:   PetscUseTypeMethod(mat, realpart);
415:   PetscFunctionReturn(PETSC_SUCCESS);
416: }

418: /*@C
419:   MatGetGhosts - Get the global indices of all ghost nodes defined by the sparse matrix

421:   Collective

423:   Input Parameter:
424: . mat - the matrix

426:   Output Parameters:
427: + nghosts - number of ghosts (for `MATBAIJ` and `MATSBAIJ` matrices there is one ghost for each matrix block)
428: - ghosts  - the global indices of the ghost points

430:   Level: advanced

432:   Note:
433:   `nghosts` and `ghosts` are suitable to pass into `VecCreateGhost()` or `VecCreateGhostBlock()`

435: .seealso: [](ch_matrices), `Mat`, `VecCreateGhost()`, `VecCreateGhostBlock()`
436: @*/
437: PetscErrorCode MatGetGhosts(Mat mat, PetscInt *nghosts, const PetscInt *ghosts[])
438: {
439:   PetscFunctionBegin;
442:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
443:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
444:   if (mat->ops->getghosts) PetscUseTypeMethod(mat, getghosts, nghosts, ghosts);
445:   else {
446:     if (nghosts) *nghosts = 0;
447:     if (ghosts) *ghosts = NULL;
448:   }
449:   PetscFunctionReturn(PETSC_SUCCESS);
450: }

452: /*@
453:   MatImaginaryPart - Moves the imaginary part of the matrix to the real part and zeros the imaginary part

455:   Logically Collective

457:   Input Parameter:
458: . mat - the matrix

460:   Level: advanced

462: .seealso: [](ch_matrices), `Mat`, `MatRealPart()`
463: @*/
464: PetscErrorCode MatImaginaryPart(Mat mat)
465: {
466:   PetscFunctionBegin;
469:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
470:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
471:   MatCheckPreallocated(mat, 1);
472:   PetscUseTypeMethod(mat, imaginarypart);
473:   PetscFunctionReturn(PETSC_SUCCESS);
474: }

476: /*@
477:   MatMissingDiagonal - Determine if sparse matrix is missing a diagonal entry (or block entry for `MATBAIJ` and `MATSBAIJ` matrices) in the nonzero structure

479:   Not Collective

481:   Input Parameter:
482: . mat - the matrix

484:   Output Parameters:
485: + missing - is any diagonal entry missing
486: - dd      - first diagonal entry that is missing (optional) on this process

488:   Level: advanced

490:   Note:
491:   This does not return diagonal entries that are in the nonzero structure but happen to have a zero numerical value

493: .seealso: [](ch_matrices), `Mat`
494: @*/
495: PetscErrorCode MatMissingDiagonal(Mat mat, PetscBool *missing, PetscInt *dd)
496: {
497:   PetscFunctionBegin;
500:   PetscAssertPointer(missing, 2);
501:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix %s", ((PetscObject)mat)->type_name);
502:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
503:   PetscUseTypeMethod(mat, missingdiagonal, missing, dd);
504:   PetscFunctionReturn(PETSC_SUCCESS);
505: }

507: // PetscClangLinter pragma disable: -fdoc-section-header-unknown
508: /*@C
509:   MatGetRow - Gets a row of a matrix.  You MUST call `MatRestoreRow()`
510:   for each row that you get to ensure that your application does
511:   not bleed memory.

513:   Not Collective

515:   Input Parameters:
516: + mat - the matrix
517: - row - the row to get

519:   Output Parameters:
520: + ncols - if not `NULL`, the number of nonzeros in `row`
521: . cols  - if not `NULL`, the column numbers
522: - vals  - if not `NULL`, the numerical values

524:   Level: advanced

526:   Notes:
527:   This routine is provided for people who need to have direct access
528:   to the structure of a matrix.  We hope that we provide enough
529:   high-level matrix routines that few users will need it.

531:   `MatGetRow()` always returns 0-based column indices, regardless of
532:   whether the internal representation is 0-based (default) or 1-based.

534:   For better efficiency, set `cols` and/or `vals` to `NULL` if you do
535:   not wish to extract these quantities.

537:   The user can only examine the values extracted with `MatGetRow()`;
538:   the values CANNOT be altered.  To change the matrix entries, one
539:   must use `MatSetValues()`.

541:   You can only have one call to `MatGetRow()` outstanding for a particular
542:   matrix at a time, per processor. `MatGetRow()` can only obtain rows
543:   associated with the given processor, it cannot get rows from the
544:   other processors; for that we suggest using `MatCreateSubMatrices()`, then
545:   `MatGetRow()` on the submatrix. The row index passed to `MatGetRow()`
546:   is in the global number of rows.

548:   Use `MatGetRowIJ()` and `MatRestoreRowIJ()` to access all the local indices of the sparse matrix.

550:   Use `MatSeqAIJGetArray()` and similar functions to access the numerical values for certain matrix types directly.

552:   Fortran Note:
553:   The calling sequence is
554: .vb
555:    MatGetRow(matrix,row,ncols,cols,values,ierr)
556:          Mat         matrix (input)
557:          PetscInt    row    (input)
558:          PetscInt    ncols  (output)
559:          PetscInt    cols(maxcols) (output)
560:          PetscScalar values(maxcols) output
561: .ve
562:   where maxcols >= maximum nonzeros in any row of the matrix.

564: .seealso: [](ch_matrices), `Mat`, `MatRestoreRow()`, `MatSetValues()`, `MatGetValues()`, `MatCreateSubMatrices()`, `MatGetDiagonal()`, `MatGetRowIJ()`, `MatRestoreRowIJ()`
565: @*/
566: PetscErrorCode MatGetRow(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt *cols[], const PetscScalar *vals[])
567: {
568:   PetscInt incols;

570:   PetscFunctionBegin;
573:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
574:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
575:   MatCheckPreallocated(mat, 1);
576:   PetscCheck(row >= mat->rmap->rstart && row < mat->rmap->rend, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Only for local rows, %" PetscInt_FMT " not in [%" PetscInt_FMT ",%" PetscInt_FMT ")", row, mat->rmap->rstart, mat->rmap->rend);
577:   PetscCall(PetscLogEventBegin(MAT_GetRow, mat, 0, 0, 0));
578:   PetscUseTypeMethod(mat, getrow, row, &incols, (PetscInt **)cols, (PetscScalar **)vals);
579:   if (ncols) *ncols = incols;
580:   PetscCall(PetscLogEventEnd(MAT_GetRow, mat, 0, 0, 0));
581:   PetscFunctionReturn(PETSC_SUCCESS);
582: }

584: /*@
585:   MatConjugate - replaces the matrix values with their complex conjugates

587:   Logically Collective

589:   Input Parameter:
590: . mat - the matrix

592:   Level: advanced

594: .seealso: [](ch_matrices), `Mat`, `MatRealPart()`, `MatImaginaryPart()`, `VecConjugate()`, `MatTranspose()`
595: @*/
596: PetscErrorCode MatConjugate(Mat mat)
597: {
598:   PetscFunctionBegin;
600:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
601:   if (PetscDefined(USE_COMPLEX) && mat->hermitian != PETSC_BOOL3_TRUE) {
602:     PetscUseTypeMethod(mat, conjugate);
603:     PetscCall(PetscObjectStateIncrease((PetscObject)mat));
604:   }
605:   PetscFunctionReturn(PETSC_SUCCESS);
606: }

608: /*@C
609:   MatRestoreRow - Frees any temporary space allocated by `MatGetRow()`.

611:   Not Collective

613:   Input Parameters:
614: + mat   - the matrix
615: . row   - the row to get
616: . ncols - the number of nonzeros
617: . cols  - the columns of the nonzeros
618: - vals  - if nonzero the column values

620:   Level: advanced

622:   Notes:
623:   This routine should be called after you have finished examining the entries.

625:   This routine zeros out `ncols`, `cols`, and `vals`. This is to prevent accidental
626:   us of the array after it has been restored. If you pass `NULL`, it will
627:   not zero the pointers.  Use of `cols` or `vals` after `MatRestoreRow()` is invalid.

629:   Fortran Note:
630:   `MatRestoreRow()` MUST be called after `MatGetRow()`
631:   before another call to `MatGetRow()` can be made.

633: .seealso: [](ch_matrices), `Mat`, `MatGetRow()`
634: @*/
635: PetscErrorCode MatRestoreRow(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt *cols[], const PetscScalar *vals[])
636: {
637:   PetscFunctionBegin;
639:   if (ncols) PetscAssertPointer(ncols, 3);
640:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
641:   if (!mat->ops->restorerow) PetscFunctionReturn(PETSC_SUCCESS);
642:   PetscUseTypeMethod(mat, restorerow, row, ncols, (PetscInt **)cols, (PetscScalar **)vals);
643:   if (ncols) *ncols = 0;
644:   if (cols) *cols = NULL;
645:   if (vals) *vals = NULL;
646:   PetscFunctionReturn(PETSC_SUCCESS);
647: }

649: /*@
650:   MatGetRowUpperTriangular - Sets a flag to enable calls to `MatGetRow()` for matrix in `MATSBAIJ` format.
651:   You should call `MatRestoreRowUpperTriangular()` after calling` MatGetRow()` and `MatRestoreRow()` to disable the flag.

653:   Not Collective

655:   Input Parameter:
656: . mat - the matrix

658:   Level: advanced

660:   Note:
661:   The flag is to ensure that users are aware that `MatGetRow()` only provides the upper triangular part of the row for the matrices in `MATSBAIJ` format.

663: .seealso: [](ch_matrices), `Mat`, `MATSBAIJ`, `MatRestoreRowUpperTriangular()`
664: @*/
665: PetscErrorCode MatGetRowUpperTriangular(Mat mat)
666: {
667:   PetscFunctionBegin;
670:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
671:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
672:   MatCheckPreallocated(mat, 1);
673:   if (!mat->ops->getrowuppertriangular) PetscFunctionReturn(PETSC_SUCCESS);
674:   PetscUseTypeMethod(mat, getrowuppertriangular);
675:   PetscFunctionReturn(PETSC_SUCCESS);
676: }

678: /*@
679:   MatRestoreRowUpperTriangular - Disable calls to `MatGetRow()` for matrix in `MATSBAIJ` format.

681:   Not Collective

683:   Input Parameter:
684: . mat - the matrix

686:   Level: advanced

688:   Note:
689:   This routine should be called after you have finished calls to `MatGetRow()` and `MatRestoreRow()`.

691: .seealso: [](ch_matrices), `Mat`, `MATSBAIJ`, `MatGetRowUpperTriangular()`
692: @*/
693: PetscErrorCode MatRestoreRowUpperTriangular(Mat mat)
694: {
695:   PetscFunctionBegin;
698:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
699:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
700:   MatCheckPreallocated(mat, 1);
701:   if (!mat->ops->restorerowuppertriangular) PetscFunctionReturn(PETSC_SUCCESS);
702:   PetscUseTypeMethod(mat, restorerowuppertriangular);
703:   PetscFunctionReturn(PETSC_SUCCESS);
704: }

706: /*@
707:   MatSetOptionsPrefix - Sets the prefix used for searching for all
708:   `Mat` options in the database.

710:   Logically Collective

712:   Input Parameters:
713: + A      - the matrix
714: - prefix - the prefix to prepend to all option names

716:   Level: advanced

718:   Notes:
719:   A hyphen (-) must NOT be given at the beginning of the prefix name.
720:   The first character of all runtime options is AUTOMATICALLY the hyphen.

722:   This is NOT used for options for the factorization of the matrix. Normally the
723:   prefix is automatically passed in from the PC calling the factorization. To set
724:   it directly use  `MatSetOptionsPrefixFactor()`

726: .seealso: [](ch_matrices), `Mat`, `MatSetFromOptions()`, `MatSetOptionsPrefixFactor()`
727: @*/
728: PetscErrorCode MatSetOptionsPrefix(Mat A, const char prefix[])
729: {
730:   PetscFunctionBegin;
732:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)A, prefix));
733:   PetscFunctionReturn(PETSC_SUCCESS);
734: }

736: /*@
737:   MatSetOptionsPrefixFactor - Sets the prefix used for searching for all matrix factor options in the database for
738:   for matrices created with `MatGetFactor()`

740:   Logically Collective

742:   Input Parameters:
743: + A      - the matrix
744: - prefix - the prefix to prepend to all option names for the factored matrix

746:   Level: developer

748:   Notes:
749:   A hyphen (-) must NOT be given at the beginning of the prefix name.
750:   The first character of all runtime options is AUTOMATICALLY the hyphen.

752:   Normally the prefix is automatically passed in from the `PC` calling the factorization. To set
753:   it directly when not using `KSP`/`PC` use  `MatSetOptionsPrefixFactor()`

755: .seealso: [](ch_matrices), `Mat`,   [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatSetFromOptions()`, `MatSetOptionsPrefix()`, `MatAppendOptionsPrefixFactor()`
756: @*/
757: PetscErrorCode MatSetOptionsPrefixFactor(Mat A, const char prefix[])
758: {
759:   PetscFunctionBegin;
761:   if (prefix) {
762:     PetscAssertPointer(prefix, 2);
763:     PetscCheck(prefix[0] != '-', PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Options prefix should not begin with a hyphen");
764:     if (prefix != A->factorprefix) {
765:       PetscCall(PetscFree(A->factorprefix));
766:       PetscCall(PetscStrallocpy(prefix, &A->factorprefix));
767:     }
768:   } else PetscCall(PetscFree(A->factorprefix));
769:   PetscFunctionReturn(PETSC_SUCCESS);
770: }

772: /*@
773:   MatAppendOptionsPrefixFactor - Appends to the prefix used for searching for all matrix factor options in the database for
774:   for matrices created with `MatGetFactor()`

776:   Logically Collective

778:   Input Parameters:
779: + A      - the matrix
780: - prefix - the prefix to prepend to all option names for the factored matrix

782:   Level: developer

784:   Notes:
785:   A hyphen (-) must NOT be given at the beginning of the prefix name.
786:   The first character of all runtime options is AUTOMATICALLY the hyphen.

788:   Normally the prefix is automatically passed in from the `PC` calling the factorization. To set
789:   it directly when not using `KSP`/`PC` use  `MatAppendOptionsPrefixFactor()`

791: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `PetscOptionsCreate()`, `PetscOptionsDestroy()`, `PetscObjectSetOptionsPrefix()`, `PetscObjectPrependOptionsPrefix()`,
792:           `PetscObjectGetOptionsPrefix()`, `TSAppendOptionsPrefix()`, `SNESAppendOptionsPrefix()`, `KSPAppendOptionsPrefix()`, `MatSetOptionsPrefixFactor()`,
793:           `MatSetOptionsPrefix()`
794: @*/
795: PetscErrorCode MatAppendOptionsPrefixFactor(Mat A, const char prefix[])
796: {
797:   size_t len1, len2, new_len;

799:   PetscFunctionBegin;
801:   if (!prefix) PetscFunctionReturn(PETSC_SUCCESS);
802:   if (!A->factorprefix) {
803:     PetscCall(MatSetOptionsPrefixFactor(A, prefix));
804:     PetscFunctionReturn(PETSC_SUCCESS);
805:   }
806:   PetscCheck(prefix[0] != '-', PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Options prefix should not begin with a hyphen");

808:   PetscCall(PetscStrlen(A->factorprefix, &len1));
809:   PetscCall(PetscStrlen(prefix, &len2));
810:   new_len = len1 + len2 + 1;
811:   PetscCall(PetscRealloc(new_len * sizeof(*A->factorprefix), &A->factorprefix));
812:   PetscCall(PetscStrncpy(A->factorprefix + len1, prefix, len2 + 1));
813:   PetscFunctionReturn(PETSC_SUCCESS);
814: }

816: /*@
817:   MatAppendOptionsPrefix - Appends to the prefix used for searching for all
818:   matrix options in the database.

820:   Logically Collective

822:   Input Parameters:
823: + A      - the matrix
824: - prefix - the prefix to prepend to all option names

826:   Level: advanced

828:   Note:
829:   A hyphen (-) must NOT be given at the beginning of the prefix name.
830:   The first character of all runtime options is AUTOMATICALLY the hyphen.

832: .seealso: [](ch_matrices), `Mat`, `MatGetOptionsPrefix()`, `MatAppendOptionsPrefixFactor()`, `MatSetOptionsPrefix()`
833: @*/
834: PetscErrorCode MatAppendOptionsPrefix(Mat A, const char prefix[])
835: {
836:   PetscFunctionBegin;
838:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)A, prefix));
839:   PetscFunctionReturn(PETSC_SUCCESS);
840: }

842: /*@
843:   MatGetOptionsPrefix - Gets the prefix used for searching for all
844:   matrix options in the database.

846:   Not Collective

848:   Input Parameter:
849: . A - the matrix

851:   Output Parameter:
852: . prefix - pointer to the prefix string used

854:   Level: advanced

856:   Fortran Note:
857:   The user should pass in a string `prefix` of
858:   sufficient length to hold the prefix.

860: .seealso: [](ch_matrices), `Mat`, `MatAppendOptionsPrefix()`, `MatSetOptionsPrefix()`, `MatAppendOptionsPrefixFactor()`, `MatSetOptionsPrefixFactor()`
861: @*/
862: PetscErrorCode MatGetOptionsPrefix(Mat A, const char *prefix[])
863: {
864:   PetscFunctionBegin;
866:   PetscAssertPointer(prefix, 2);
867:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)A, prefix));
868:   PetscFunctionReturn(PETSC_SUCCESS);
869: }

871: /*@
872:   MatGetState - Gets the state of a `Mat`. Same value as returned by `PetscObjectStateGet()`

874:   Not Collective

876:   Input Parameter:
877: . A - the matrix

879:   Output Parameter:
880: . state - the object state

882:   Level: advanced

884:   Note:
885:   Object state is an integer which gets increased every time
886:   the object is changed. By saving and later querying the object state
887:   one can determine whether information about the object is still current.

889:   See `MatGetNonzeroState()` to determine if the nonzero structure of the matrix has changed.

891: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `PetscObjectStateGet()`, `MatGetNonzeroState()`
892: @*/
893: PetscErrorCode MatGetState(Mat A, PetscObjectState *state)
894: {
895:   PetscFunctionBegin;
897:   PetscAssertPointer(state, 2);
898:   PetscCall(PetscObjectStateGet((PetscObject)A, state));
899:   PetscFunctionReturn(PETSC_SUCCESS);
900: }

902: /*@
903:   MatResetPreallocation - Reset matrix to use the original preallocation values provided by the user, for example with `MatXAIJSetPreallocation()`

905:   Collective

907:   Input Parameter:
908: . A - the matrix

910:   Level: beginner

912:   Notes:
913:   After calling `MatAssemblyBegin()` and `MatAssemblyEnd()` with `MAT_FINAL_ASSEMBLY` the matrix data structures represent the nonzeros assigned to the
914:   matrix. If that space is less than the preallocated space that extra preallocated space is no longer available to take on new values. `MatResetPreallocation()`
915:   makes all of the preallocation space available

917:   Current values in the matrix are lost in this call.

919:   Currently only supported for  `MATAIJ` matrices.

921: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatXAIJSetPreallocation()`
922: @*/
923: PetscErrorCode MatResetPreallocation(Mat A)
924: {
925:   PetscFunctionBegin;
928:   PetscCheck(A->insertmode == NOT_SET_VALUES, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot reset preallocation after setting some values but not yet calling MatAssemblyBegin()/MatAssemblyEnd()");
929:   if (A->num_ass == 0) PetscFunctionReturn(PETSC_SUCCESS);
930:   PetscUseMethod(A, "MatResetPreallocation_C", (Mat), (A));
931:   PetscFunctionReturn(PETSC_SUCCESS);
932: }

934: /*@
935:   MatSetUp - Sets up the internal matrix data structures for later use by the matrix

937:   Collective

939:   Input Parameter:
940: . A - the matrix

942:   Level: advanced

944:   Notes:
945:   If the user has not set preallocation for this matrix then an efficient algorithm will be used for the first round of
946:   setting values in the matrix.

948:   This routine is called internally by other `Mat` functions when needed so rarely needs to be called by users

950: .seealso: [](ch_matrices), `Mat`, `MatMult()`, `MatCreate()`, `MatDestroy()`, `MatXAIJSetPreallocation()`
951: @*/
952: PetscErrorCode MatSetUp(Mat A)
953: {
954:   PetscFunctionBegin;
956:   if (!((PetscObject)A)->type_name) {
957:     PetscMPIInt size;

959:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
960:     PetscCall(MatSetType(A, size == 1 ? MATSEQAIJ : MATMPIAIJ));
961:   }
962:   if (!A->preallocated) PetscTryTypeMethod(A, setup);
963:   PetscCall(PetscLayoutSetUp(A->rmap));
964:   PetscCall(PetscLayoutSetUp(A->cmap));
965:   A->preallocated = PETSC_TRUE;
966:   PetscFunctionReturn(PETSC_SUCCESS);
967: }

969: #if defined(PETSC_HAVE_SAWS)
970: #include <petscviewersaws.h>
971: #endif

973: /*
974:    If threadsafety is on extraneous matrices may be printed

976:    This flag cannot be stored in the matrix because the original matrix in MatView() may assemble a new matrix which is passed into MatViewFromOptions()
977: */
978: #if !defined(PETSC_HAVE_THREADSAFETY)
979: static PetscInt insidematview = 0;
980: #endif

982: /*@
983:   MatViewFromOptions - View properties of the matrix based on options set in the options database

985:   Collective

987:   Input Parameters:
988: + A    - the matrix
989: . obj  - optional additional object that provides the options prefix to use
990: - name - command line option

992:   Options Database Key:
993: . -mat_view [viewertype]:... - the viewer and its options

995:   Level: intermediate

997:   Note:
998: .vb
999:     If no value is provided ascii:stdout is used
1000:        ascii[:[filename][:[format][:append]]]    defaults to stdout - format can be one of ascii_info, ascii_info_detail, or ascii_matlab,
1001:                                                   for example ascii::ascii_info prints just the information about the object not all details
1002:                                                   unless :append is given filename opens in write mode, overwriting what was already there
1003:        binary[:[filename][:[format][:append]]]   defaults to the file binaryoutput
1004:        draw[:drawtype[:filename]]                for example, draw:tikz, draw:tikz:figure.tex  or draw:x
1005:        socket[:port]                             defaults to the standard output port
1006:        saws[:communicatorname]                    publishes object to the Scientific Application Webserver (SAWs)
1007: .ve

1009: .seealso: [](ch_matrices), `Mat`, `MatView()`, `PetscObjectViewFromOptions()`, `MatCreate()`
1010: @*/
1011: PetscErrorCode MatViewFromOptions(Mat A, PetscObject obj, const char name[])
1012: {
1013:   PetscFunctionBegin;
1015: #if !defined(PETSC_HAVE_THREADSAFETY)
1016:   if (insidematview) PetscFunctionReturn(PETSC_SUCCESS);
1017: #endif
1018:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
1019:   PetscFunctionReturn(PETSC_SUCCESS);
1020: }

1022: /*@
1023:   MatView - display information about a matrix in a variety ways

1025:   Collective on viewer

1027:   Input Parameters:
1028: + mat    - the matrix
1029: - viewer - visualization context

1031:   Options Database Keys:
1032: + -mat_view ::ascii_info           - Prints info on matrix at conclusion of `MatAssemblyEnd()`
1033: . -mat_view ::ascii_info_detail    - Prints more detailed info
1034: . -mat_view                        - Prints matrix in ASCII format
1035: . -mat_view ::ascii_matlab         - Prints matrix in MATLAB format
1036: . -mat_view draw                   - PetscDraws nonzero structure of matrix, using `MatView()` and `PetscDrawOpenX()`.
1037: . -display <name>                  - Sets display name (default is host)
1038: . -draw_pause <sec>                - Sets number of seconds to pause after display
1039: . -mat_view socket                 - Sends matrix to socket, can be accessed from MATLAB (see Users-Manual: ch_matlab for details)
1040: . -viewer_socket_machine <machine> - -
1041: . -viewer_socket_port <port>       - -
1042: . -mat_view binary                 - save matrix to file in binary format
1043: - -viewer_binary_filename <name>   - -

1045:   Level: beginner

1047:   Notes:
1048:   The available visualization contexts include
1049: +    `PETSC_VIEWER_STDOUT_SELF`   - for sequential matrices
1050: .    `PETSC_VIEWER_STDOUT_WORLD`  - for parallel matrices created on `PETSC_COMM_WORLD`
1051: .    `PETSC_VIEWER_STDOUT_`(comm) - for matrices created on MPI communicator comm
1052: -     `PETSC_VIEWER_DRAW_WORLD`   - graphical display of nonzero structure

1054:   The user can open alternative visualization contexts with
1055: +    `PetscViewerASCIIOpen()`  - Outputs matrix to a specified file
1056: .    `PetscViewerBinaryOpen()` - Outputs matrix in binary to a  specified file; corresponding input uses `MatLoad()`
1057: .    `PetscViewerDrawOpen()`   - Outputs nonzero matrix nonzero structure to an X window display
1058: -    `PetscViewerSocketOpen()` - Outputs matrix to Socket viewer, `PETSCVIEWERSOCKET`. Only the `MATSEQDENSE` and `MATAIJ` types support this viewer.

1060:   The user can call `PetscViewerPushFormat()` to specify the output
1061:   format of ASCII printed objects (when using `PETSC_VIEWER_STDOUT_SELF`,
1062:   `PETSC_VIEWER_STDOUT_WORLD` and `PetscViewerASCIIOpen()`).  Available formats include
1063: +    `PETSC_VIEWER_DEFAULT`           - default, prints matrix contents
1064: .    `PETSC_VIEWER_ASCII_MATLAB`      - prints matrix contents in MATLAB format
1065: .    `PETSC_VIEWER_ASCII_DENSE`       - prints entire matrix including zeros
1066: .    `PETSC_VIEWER_ASCII_COMMON`      - prints matrix contents, using a sparse  format common among all matrix types
1067: .    `PETSC_VIEWER_ASCII_IMPL`        - prints matrix contents, using an implementation-specific format (which is in many cases the same as the default)
1068: .    `PETSC_VIEWER_ASCII_INFO`        - prints basic information about the matrix size and structure (not the matrix entries)
1069: -    `PETSC_VIEWER_ASCII_INFO_DETAIL` - prints more detailed information about the matrix nonzero structure (still not vector or matrix entries)

1071:   The ASCII viewers are only recommended for small matrices on at most a moderate number of processes,
1072:   the program will seemingly hang and take hours for larger matrices, for larger matrices one should use the binary format.

1074:   In the debugger you can do "call MatView(mat,0)" to display the matrix. (The same holds for any PETSc object viewer).

1076:   See the manual page for `MatLoad()` for the exact format of the binary file when the binary
1077:   viewer is used.

1079:   See share/petsc/matlab/PetscBinaryRead.m for a MATLAB code that can read in the binary file when the binary
1080:   viewer is used and lib/petsc/bin/PetscBinaryIO.py for loading them into Python.

1082:   One can use '-mat_view draw -draw_pause -1' to pause the graphical display of matrix nonzero structure,
1083:   and then use the following mouse functions.
1084: .vb
1085:   left mouse: zoom in
1086:   middle mouse: zoom out
1087:   right mouse: continue with the simulation
1088: .ve

1090: .seealso: [](ch_matrices), `Mat`, `PetscViewerPushFormat()`, `PetscViewerASCIIOpen()`, `PetscViewerDrawOpen()`, `PetscViewer`,
1091:           `PetscViewerSocketOpen()`, `PetscViewerBinaryOpen()`, `MatLoad()`, `MatViewFromOptions()`
1092: @*/
1093: PetscErrorCode MatView(Mat mat, PetscViewer viewer)
1094: {
1095:   PetscInt          rows, cols, rbs, cbs;
1096:   PetscBool         isascii, isstring, issaws;
1097:   PetscViewerFormat format;
1098:   PetscMPIInt       size;

1100:   PetscFunctionBegin;
1103:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mat), &viewer));

1106:   PetscCall(PetscViewerGetFormat(viewer, &format));
1107:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
1108:   if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) PetscFunctionReturn(PETSC_SUCCESS);

1110: #if !defined(PETSC_HAVE_THREADSAFETY)
1111:   insidematview++;
1112: #endif
1113:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1114:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1115:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws));
1116:   PetscCheck((isascii && (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL)) || !mat->factortype, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "No viewers for factored matrix except ASCII, info, or info_detail");

1118:   PetscCall(PetscLogEventBegin(MAT_View, mat, viewer, 0, 0));
1119:   if (isascii) {
1120:     if (!mat->preallocated) {
1121:       PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix has not been preallocated yet\n"));
1122: #if !defined(PETSC_HAVE_THREADSAFETY)
1123:       insidematview--;
1124: #endif
1125:       PetscCall(PetscLogEventEnd(MAT_View, mat, viewer, 0, 0));
1126:       PetscFunctionReturn(PETSC_SUCCESS);
1127:     }
1128:     if (!mat->assembled) {
1129:       PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix has not been assembled yet\n"));
1130: #if !defined(PETSC_HAVE_THREADSAFETY)
1131:       insidematview--;
1132: #endif
1133:       PetscCall(PetscLogEventEnd(MAT_View, mat, viewer, 0, 0));
1134:       PetscFunctionReturn(PETSC_SUCCESS);
1135:     }
1136:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mat, viewer));
1137:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1138:       MatNullSpace nullsp, transnullsp;

1140:       PetscCall(PetscViewerASCIIPushTab(viewer));
1141:       PetscCall(MatGetSize(mat, &rows, &cols));
1142:       PetscCall(MatGetBlockSizes(mat, &rbs, &cbs));
1143:       if (rbs != 1 || cbs != 1) {
1144:         if (rbs != cbs) PetscCall(PetscViewerASCIIPrintf(viewer, "rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT ", rbs=%" PetscInt_FMT ", cbs=%" PetscInt_FMT "%s\n", rows, cols, rbs, cbs, mat->bsizes ? " variable blocks set" : ""));
1145:         else PetscCall(PetscViewerASCIIPrintf(viewer, "rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT ", bs=%" PetscInt_FMT "%s\n", rows, cols, rbs, mat->bsizes ? " variable blocks set" : ""));
1146:       } else PetscCall(PetscViewerASCIIPrintf(viewer, "rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT "\n", rows, cols));
1147:       if (mat->factortype) {
1148:         MatSolverType solver;
1149:         PetscCall(MatFactorGetSolverType(mat, &solver));
1150:         PetscCall(PetscViewerASCIIPrintf(viewer, "package used to perform factorization: %s\n", solver));
1151:       }
1152:       if (mat->ops->getinfo) {
1153:         MatInfo info;
1154:         PetscCall(MatGetInfo(mat, MAT_GLOBAL_SUM, &info));
1155:         PetscCall(PetscViewerASCIIPrintf(viewer, "total: nonzeros=%.f, allocated nonzeros=%.f\n", info.nz_used, info.nz_allocated));
1156:         if (!mat->factortype) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of mallocs used during MatSetValues calls=%" PetscInt_FMT "\n", (PetscInt)info.mallocs));
1157:       }
1158:       PetscCall(MatGetNullSpace(mat, &nullsp));
1159:       PetscCall(MatGetTransposeNullSpace(mat, &transnullsp));
1160:       if (nullsp) PetscCall(PetscViewerASCIIPrintf(viewer, "  has attached null space\n"));
1161:       if (transnullsp && transnullsp != nullsp) PetscCall(PetscViewerASCIIPrintf(viewer, "  has attached transposed null space\n"));
1162:       PetscCall(MatGetNearNullSpace(mat, &nullsp));
1163:       if (nullsp) PetscCall(PetscViewerASCIIPrintf(viewer, "  has attached near null space\n"));
1164:       PetscCall(PetscViewerASCIIPushTab(viewer));
1165:       PetscCall(MatProductView(mat, viewer));
1166:       PetscCall(PetscViewerASCIIPopTab(viewer));
1167:       if (mat->bsizes && format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1168:         IS tmp;

1170:         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)viewer), mat->nblocks, mat->bsizes, PETSC_USE_POINTER, &tmp));
1171:         PetscCall(PetscObjectSetName((PetscObject)tmp, "Block Sizes"));
1172:         PetscCall(PetscViewerASCIIPushTab(viewer));
1173:         PetscCall(ISView(tmp, viewer));
1174:         PetscCall(PetscViewerASCIIPopTab(viewer));
1175:         PetscCall(ISDestroy(&tmp));
1176:       }
1177:     }
1178:   } else if (issaws) {
1179: #if defined(PETSC_HAVE_SAWS)
1180:     PetscMPIInt rank;

1182:     PetscCall(PetscObjectName((PetscObject)mat));
1183:     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1184:     if (!((PetscObject)mat)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)mat, viewer));
1185: #endif
1186:   } else if (isstring) {
1187:     const char *type;
1188:     PetscCall(MatGetType(mat, &type));
1189:     PetscCall(PetscViewerStringSPrintf(viewer, " MatType: %-7.7s", type));
1190:     PetscTryTypeMethod(mat, view, viewer);
1191:   }
1192:   if ((format == PETSC_VIEWER_NATIVE || format == PETSC_VIEWER_LOAD_BALANCE) && mat->ops->viewnative) {
1193:     PetscCall(PetscViewerASCIIPushTab(viewer));
1194:     PetscUseTypeMethod(mat, viewnative, viewer);
1195:     PetscCall(PetscViewerASCIIPopTab(viewer));
1196:   } else if (mat->ops->view) {
1197:     PetscCall(PetscViewerASCIIPushTab(viewer));
1198:     PetscUseTypeMethod(mat, view, viewer);
1199:     PetscCall(PetscViewerASCIIPopTab(viewer));
1200:   }
1201:   if (isascii) {
1202:     PetscCall(PetscViewerGetFormat(viewer, &format));
1203:     if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscCall(PetscViewerASCIIPopTab(viewer));
1204:   }
1205:   PetscCall(PetscLogEventEnd(MAT_View, mat, viewer, 0, 0));
1206: #if !defined(PETSC_HAVE_THREADSAFETY)
1207:   insidematview--;
1208: #endif
1209:   PetscFunctionReturn(PETSC_SUCCESS);
1210: }

1212: #if defined(PETSC_USE_DEBUG)
1213: #include <../src/sys/totalview/tv_data_display.h>
1214: PETSC_UNUSED static int TV_display_type(const struct _p_Mat *mat)
1215: {
1216:   TV_add_row("Local rows", "int", &mat->rmap->n);
1217:   TV_add_row("Local columns", "int", &mat->cmap->n);
1218:   TV_add_row("Global rows", "int", &mat->rmap->N);
1219:   TV_add_row("Global columns", "int", &mat->cmap->N);
1220:   TV_add_row("Typename", TV_ascii_string_type, ((PetscObject)mat)->type_name);
1221:   return TV_format_OK;
1222: }
1223: #endif

1225: /*@
1226:   MatLoad - Loads a matrix that has been stored in binary/HDF5 format
1227:   with `MatView()`.  The matrix format is determined from the options database.
1228:   Generates a parallel MPI matrix if the communicator has more than one
1229:   processor.  The default matrix type is `MATAIJ`.

1231:   Collective

1233:   Input Parameters:
1234: + mat    - the newly loaded matrix, this needs to have been created with `MatCreate()`
1235:             or some related function before a call to `MatLoad()`
1236: - viewer - `PETSCVIEWERBINARY`/`PETSCVIEWERHDF5` file viewer

1238:   Options Database Key:
1239: . -matload_block_size <bs> - set block size

1241:   Level: beginner

1243:   Notes:
1244:   If the `Mat` type has not yet been given then `MATAIJ` is used, call `MatSetFromOptions()` on the
1245:   `Mat` before calling this routine if you wish to set it from the options database.

1247:   `MatLoad()` automatically loads into the options database any options
1248:   given in the file filename.info where filename is the name of the file
1249:   that was passed to the `PetscViewerBinaryOpen()`. The options in the info
1250:   file will be ignored if you use the -viewer_binary_skip_info option.

1252:   If the type or size of mat is not set before a call to `MatLoad()`, PETSc
1253:   sets the default matrix type AIJ and sets the local and global sizes.
1254:   If type and/or size is already set, then the same are used.

1256:   In parallel, each processor can load a subset of rows (or the
1257:   entire matrix).  This routine is especially useful when a large
1258:   matrix is stored on disk and only part of it is desired on each
1259:   processor.  For example, a parallel solver may access only some of
1260:   the rows from each processor.  The algorithm used here reads
1261:   relatively small blocks of data rather than reading the entire
1262:   matrix and then subsetting it.

1264:   Viewer's `PetscViewerType` must be either `PETSCVIEWERBINARY` or `PETSCVIEWERHDF5`.
1265:   Such viewer can be created using `PetscViewerBinaryOpen()` or `PetscViewerHDF5Open()`,
1266:   or the sequence like
1267: .vb
1268:     `PetscViewer` v;
1269:     `PetscViewerCreate`(`PETSC_COMM_WORLD`,&v);
1270:     `PetscViewerSetType`(v,`PETSCVIEWERBINARY`);
1271:     `PetscViewerSetFromOptions`(v);
1272:     `PetscViewerFileSetMode`(v,`FILE_MODE_READ`);
1273:     `PetscViewerFileSetName`(v,"datafile");
1274: .ve
1275:   The optional `PetscViewerSetFromOptions()` call allows overriding `PetscViewerSetType()` using the option
1276: $ -viewer_type {binary, hdf5}

1278:   See the example src/ksp/ksp/tutorials/ex27.c with the first approach,
1279:   and src/mat/tutorials/ex10.c with the second approach.

1281:   In case of `PETSCVIEWERBINARY`, a native PETSc binary format is used. Each of the blocks
1282:   is read onto MPI rank 0 and then shipped to its destination MPI rank, one after another.
1283:   Multiple objects, both matrices and vectors, can be stored within the same file.
1284:   Their `PetscObject` name is ignored; they are loaded in the order of their storage.

1286:   Most users should not need to know the details of the binary storage
1287:   format, since `MatLoad()` and `MatView()` completely hide these details.
1288:   But for anyone who is interested, the standard binary matrix storage
1289:   format is

1291: .vb
1292:     PetscInt    MAT_FILE_CLASSID
1293:     PetscInt    number of rows
1294:     PetscInt    number of columns
1295:     PetscInt    total number of nonzeros
1296:     PetscInt    *number nonzeros in each row
1297:     PetscInt    *column indices of all nonzeros (starting index is zero)
1298:     PetscScalar *values of all nonzeros
1299: .ve
1300:   If PETSc was not configured with `--with-64-bit-indices` then only `MATMPIAIJ` matrices with more than `PETSC_INT_MAX` non-zeros can be
1301:   stored or loaded (each MPI process part of the matrix must have less than `PETSC_INT_MAX` nonzeros). Since the total nonzero count in this
1302:   case will not fit in a (32-bit) `PetscInt` the value `PETSC_INT_MAX` is used for the header entry `total number of nonzeros`.

1304:   PETSc automatically does the byte swapping for
1305:   machines that store the bytes reversed. Thus if you write your own binary
1306:   read/write routines you have to swap the bytes; see `PetscBinaryRead()`
1307:   and `PetscBinaryWrite()` to see how this may be done.

1309:   In case of `PETSCVIEWERHDF5`, a parallel HDF5 reader is used.
1310:   Each processor's chunk is loaded independently by its owning MPI process.
1311:   Multiple objects, both matrices and vectors, can be stored within the same file.
1312:   They are looked up by their PetscObject name.

1314:   As the MATLAB MAT-File Version 7.3 format is also a HDF5 flavor, we decided to use
1315:   by default the same structure and naming of the AIJ arrays and column count
1316:   within the HDF5 file. This means that a MAT file saved with -v7.3 flag, e.g.
1317: $    save example.mat A b -v7.3
1318:   can be directly read by this routine (see Reference 1 for details).

1320:   Depending on your MATLAB version, this format might be a default,
1321:   otherwise you can set it as default in Preferences.

1323:   Unless -nocompression flag is used to save the file in MATLAB,
1324:   PETSc must be configured with ZLIB package.

1326:   See also examples src/mat/tutorials/ex10.c and src/ksp/ksp/tutorials/ex27.c

1328:   This reader currently supports only real `MATSEQAIJ`, `MATMPIAIJ`, `MATSEQDENSE` and `MATMPIDENSE` matrices for `PETSCVIEWERHDF5`

1330:   Corresponding `MatView()` is not yet implemented.

1332:   The loaded matrix is actually a transpose of the original one in MATLAB,
1333:   unless you push `PETSC_VIEWER_HDF5_MAT` format (see examples above).
1334:   With this format, matrix is automatically transposed by PETSc,
1335:   unless the matrix is marked as SPD or symmetric
1336:   (see `MatSetOption()`, `MAT_SPD`, `MAT_SYMMETRIC`).

1338:   See MATLAB Documentation on `save()`, <https://www.mathworks.com/help/matlab/ref/save.html#btox10b-1-version>

1340: .seealso: [](ch_matrices), `Mat`, `PetscViewerBinaryOpen()`, `PetscViewerSetType()`, `MatView()`, `VecLoad()`
1341:  @*/
1342: PetscErrorCode MatLoad(Mat mat, PetscViewer viewer)
1343: {
1344:   PetscBool flg;

1346:   PetscFunctionBegin;

1350:   if (!((PetscObject)mat)->type_name) PetscCall(MatSetType(mat, MATAIJ));

1352:   flg = PETSC_FALSE;
1353:   PetscCall(PetscOptionsGetBool(((PetscObject)mat)->options, ((PetscObject)mat)->prefix, "-matload_symmetric", &flg, NULL));
1354:   if (flg) {
1355:     PetscCall(MatSetOption(mat, MAT_SYMMETRIC, PETSC_TRUE));
1356:     PetscCall(MatSetOption(mat, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
1357:   }
1358:   flg = PETSC_FALSE;
1359:   PetscCall(PetscOptionsGetBool(((PetscObject)mat)->options, ((PetscObject)mat)->prefix, "-matload_spd", &flg, NULL));
1360:   if (flg) PetscCall(MatSetOption(mat, MAT_SPD, PETSC_TRUE));

1362:   PetscCall(PetscLogEventBegin(MAT_Load, mat, viewer, 0, 0));
1363:   PetscUseTypeMethod(mat, load, viewer);
1364:   PetscCall(PetscLogEventEnd(MAT_Load, mat, viewer, 0, 0));
1365:   PetscFunctionReturn(PETSC_SUCCESS);
1366: }

1368: static PetscErrorCode MatDestroy_Redundant(Mat_Redundant **redundant)
1369: {
1370:   Mat_Redundant *redund = *redundant;

1372:   PetscFunctionBegin;
1373:   if (redund) {
1374:     if (redund->matseq) { /* via MatCreateSubMatrices()  */
1375:       PetscCall(ISDestroy(&redund->isrow));
1376:       PetscCall(ISDestroy(&redund->iscol));
1377:       PetscCall(MatDestroySubMatrices(1, &redund->matseq));
1378:     } else {
1379:       PetscCall(PetscFree2(redund->send_rank, redund->recv_rank));
1380:       PetscCall(PetscFree(redund->sbuf_j));
1381:       PetscCall(PetscFree(redund->sbuf_a));
1382:       for (PetscInt i = 0; i < redund->nrecvs; i++) {
1383:         PetscCall(PetscFree(redund->rbuf_j[i]));
1384:         PetscCall(PetscFree(redund->rbuf_a[i]));
1385:       }
1386:       PetscCall(PetscFree4(redund->sbuf_nz, redund->rbuf_nz, redund->rbuf_j, redund->rbuf_a));
1387:     }

1389:     if (redund->subcomm) PetscCall(PetscCommDestroy(&redund->subcomm));
1390:     PetscCall(PetscFree(redund));
1391:   }
1392:   PetscFunctionReturn(PETSC_SUCCESS);
1393: }

1395: /*@
1396:   MatDestroy - Frees space taken by a matrix.

1398:   Collective

1400:   Input Parameter:
1401: . A - the matrix

1403:   Level: beginner

1405:   Developer Note:
1406:   Some special arrays of matrices are not destroyed in this routine but instead by the routines called by
1407:   `MatDestroySubMatrices()`. Thus one must be sure that any changes here must also be made in those routines.
1408:   `MatHeaderMerge()` and `MatHeaderReplace()` also manipulate the data in the `Mat` object and likely need changes
1409:   if changes are needed here.

1411: .seealso: [](ch_matrices), `Mat`, `MatCreate()`
1412: @*/
1413: PetscErrorCode MatDestroy(Mat *A)
1414: {
1415:   PetscFunctionBegin;
1416:   if (!*A) PetscFunctionReturn(PETSC_SUCCESS);
1418:   if (--((PetscObject)*A)->refct > 0) {
1419:     *A = NULL;
1420:     PetscFunctionReturn(PETSC_SUCCESS);
1421:   }

1423:   /* if memory was published with SAWs then destroy it */
1424:   PetscCall(PetscObjectSAWsViewOff((PetscObject)*A));
1425:   PetscTryTypeMethod(*A, destroy);

1427:   PetscCall(PetscFree((*A)->factorprefix));
1428:   PetscCall(PetscFree((*A)->defaultvectype));
1429:   PetscCall(PetscFree((*A)->defaultrandtype));
1430:   PetscCall(PetscFree((*A)->bsizes));
1431:   PetscCall(PetscFree((*A)->solvertype));
1432:   for (PetscInt i = 0; i < MAT_FACTOR_NUM_TYPES; i++) PetscCall(PetscFree((*A)->preferredordering[i]));
1433:   if ((*A)->redundant && (*A)->redundant->matseq[0] == *A) (*A)->redundant->matseq[0] = NULL;
1434:   PetscCall(MatDestroy_Redundant(&(*A)->redundant));
1435:   PetscCall(MatProductClear(*A));
1436:   PetscCall(MatNullSpaceDestroy(&(*A)->nullsp));
1437:   PetscCall(MatNullSpaceDestroy(&(*A)->transnullsp));
1438:   PetscCall(MatNullSpaceDestroy(&(*A)->nearnullsp));
1439:   PetscCall(MatDestroy(&(*A)->schur));
1440:   PetscCall(PetscLayoutDestroy(&(*A)->rmap));
1441:   PetscCall(PetscLayoutDestroy(&(*A)->cmap));
1442:   PetscCall(PetscHeaderDestroy(A));
1443:   PetscFunctionReturn(PETSC_SUCCESS);
1444: }

1446: // PetscClangLinter pragma disable: -fdoc-section-header-unknown
1447: /*@
1448:   MatSetValues - Inserts or adds a block of values into a matrix.
1449:   These values may be cached, so `MatAssemblyBegin()` and `MatAssemblyEnd()`
1450:   MUST be called after all calls to `MatSetValues()` have been completed.

1452:   Not Collective

1454:   Input Parameters:
1455: + mat  - the matrix
1456: . v    - a logically two-dimensional array of values
1457: . m    - the number of rows
1458: . idxm - the global indices of the rows
1459: . n    - the number of columns
1460: . idxn - the global indices of the columns
1461: - addv - either `ADD_VALUES` to add values to any existing entries, or `INSERT_VALUES` to replace existing entries with new values

1463:   Level: beginner

1465:   Notes:
1466:   By default the values, `v`, are stored row-oriented. See `MatSetOption()` for other options.

1468:   Calls to `MatSetValues()` with the `INSERT_VALUES` and `ADD_VALUES`
1469:   options cannot be mixed without intervening calls to the assembly
1470:   routines.

1472:   `MatSetValues()` uses 0-based row and column numbers in Fortran
1473:   as well as in C.

1475:   Negative indices may be passed in `idxm` and `idxn`, these rows and columns are
1476:   simply ignored. This allows easily inserting element stiffness matrices
1477:   with homogeneous Dirichlet boundary conditions that you don't want represented
1478:   in the matrix.

1480:   Efficiency Alert:
1481:   The routine `MatSetValuesBlocked()` may offer much better efficiency
1482:   for users of block sparse formats (`MATSEQBAIJ` and `MATMPIBAIJ`).

1484:   Fortran Notes:
1485:   If any of `idxm`, `idxn`, and `v` are scalars pass them using, for example,
1486: .vb
1487:   MatSetValues(mat, one, [idxm], one, [idxn], [v], INSERT_VALUES)
1488: .ve

1490:   If `v` is a two-dimensional array use `reshape()` to pass it as a one dimensional array

1492:   Developer Note:
1493:   This is labeled with C so does not automatically generate Fortran stubs and interfaces
1494:   because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.

1496: .seealso: [](ch_matrices), `Mat`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`,
1497:           `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`
1498: @*/
1499: PetscErrorCode MatSetValues(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], const PetscScalar v[], InsertMode addv)
1500: {
1501:   PetscFunctionBeginHot;
1504:   if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */
1505:   PetscAssertPointer(idxm, 3);
1506:   PetscAssertPointer(idxn, 5);
1507:   MatCheckPreallocated(mat, 1);

1509:   if (mat->insertmode == NOT_SET_VALUES) mat->insertmode = addv;
1510:   else PetscCheck(mat->insertmode == addv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add values and insert values");

1512:   if (PetscDefined(USE_DEBUG)) {
1513:     PetscInt i, j;

1515:     PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
1516:     if (v) {
1517:       for (i = 0; i < m; i++) {
1518:         for (j = 0; j < n; j++) {
1519:           if (mat->erroriffailure && PetscIsInfOrNanScalar(v[i * n + j]))
1520: #if defined(PETSC_USE_COMPLEX)
1521:             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Inserting %g+i%g at matrix entry (%" PetscInt_FMT ",%" PetscInt_FMT ")", (double)PetscRealPart(v[i * n + j]), (double)PetscImaginaryPart(v[i * n + j]), idxm[i], idxn[j]);
1522: #else
1523:             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Inserting %g at matrix entry (%" PetscInt_FMT ",%" PetscInt_FMT ")", (double)v[i * n + j], idxm[i], idxn[j]);
1524: #endif
1525:         }
1526:       }
1527:     }
1528:     for (i = 0; i < m; i++) PetscCheck(idxm[i] < mat->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot insert in row %" PetscInt_FMT ", maximum is %" PetscInt_FMT, idxm[i], mat->rmap->N - 1);
1529:     for (i = 0; i < n; i++) PetscCheck(idxn[i] < mat->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot insert in column %" PetscInt_FMT ", maximum is %" PetscInt_FMT, idxn[i], mat->cmap->N - 1);
1530:   }

1532:   if (mat->assembled) {
1533:     mat->was_assembled = PETSC_TRUE;
1534:     mat->assembled     = PETSC_FALSE;
1535:   }
1536:   PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0));
1537:   PetscUseTypeMethod(mat, setvalues, m, idxm, n, idxn, v, addv);
1538:   PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0));
1539:   PetscFunctionReturn(PETSC_SUCCESS);
1540: }

1542: // PetscClangLinter pragma disable: -fdoc-section-header-unknown
1543: /*@
1544:   MatSetValuesIS - Inserts or adds a block of values into a matrix using an `IS` to indicate the rows and columns
1545:   These values may be cached, so `MatAssemblyBegin()` and `MatAssemblyEnd()`
1546:   MUST be called after all calls to `MatSetValues()` have been completed.

1548:   Not Collective

1550:   Input Parameters:
1551: + mat  - the matrix
1552: . v    - a logically two-dimensional array of values
1553: . ism  - the rows to provide
1554: . isn  - the columns to provide
1555: - addv - either `ADD_VALUES` to add values to any existing entries, or `INSERT_VALUES` to replace existing entries with new values

1557:   Level: beginner

1559:   Notes:
1560:   By default the values, `v`, are stored row-oriented. See `MatSetOption()` for other options.

1562:   Calls to `MatSetValues()` with the `INSERT_VALUES` and `ADD_VALUES`
1563:   options cannot be mixed without intervening calls to the assembly
1564:   routines.

1566:   `MatSetValues()` uses 0-based row and column numbers in Fortran
1567:   as well as in C.

1569:   Negative indices may be passed in `ism` and `isn`, these rows and columns are
1570:   simply ignored. This allows easily inserting element stiffness matrices
1571:   with homogeneous Dirichlet boundary conditions that you don't want represented
1572:   in the matrix.

1574:   Efficiency Alert:
1575:   The routine `MatSetValuesBlocked()` may offer much better efficiency
1576:   for users of block sparse formats (`MATSEQBAIJ` and `MATMPIBAIJ`).

1578:   This is currently not optimized for any particular `ISType`

1580:   Developer Note:
1581:   This is labeled with C so does not automatically generate Fortran stubs and interfaces
1582:   because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.

1584: .seealso: [](ch_matrices), `Mat`, `MatSetOption()`, `MatSetValues()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`,
1585:           `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`
1586: @*/
1587: PetscErrorCode MatSetValuesIS(Mat mat, IS ism, IS isn, const PetscScalar v[], InsertMode addv)
1588: {
1589:   PetscInt        m, n;
1590:   const PetscInt *rows, *cols;

1592:   PetscFunctionBeginHot;
1594:   PetscCall(ISGetIndices(ism, &rows));
1595:   PetscCall(ISGetIndices(isn, &cols));
1596:   PetscCall(ISGetLocalSize(ism, &m));
1597:   PetscCall(ISGetLocalSize(isn, &n));
1598:   PetscCall(MatSetValues(mat, m, rows, n, cols, v, addv));
1599:   PetscCall(ISRestoreIndices(ism, &rows));
1600:   PetscCall(ISRestoreIndices(isn, &cols));
1601:   PetscFunctionReturn(PETSC_SUCCESS);
1602: }

1604: /*@
1605:   MatSetValuesRowLocal - Inserts a row (block row for `MATBAIJ` matrices) of nonzero
1606:   values into a matrix

1608:   Not Collective

1610:   Input Parameters:
1611: + mat - the matrix
1612: . row - the (block) row to set
1613: - v   - a logically two-dimensional array of values

1615:   Level: intermediate

1617:   Notes:
1618:   The values, `v`, are column-oriented (for the block version) and sorted

1620:   All the nonzero values in `row` must be provided

1622:   The matrix must have previously had its column indices set, likely by having been assembled.

1624:   `row` must belong to this MPI process

1626: .seealso: [](ch_matrices), `Mat`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`,
1627:           `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`, `MatSetValues()`, `MatSetValuesRow()`, `MatSetLocalToGlobalMapping()`
1628: @*/
1629: PetscErrorCode MatSetValuesRowLocal(Mat mat, PetscInt row, const PetscScalar v[])
1630: {
1631:   PetscInt globalrow;

1633:   PetscFunctionBegin;
1636:   PetscAssertPointer(v, 3);
1637:   PetscCall(ISLocalToGlobalMappingApply(mat->rmap->mapping, 1, &row, &globalrow));
1638:   PetscCall(MatSetValuesRow(mat, globalrow, v));
1639:   PetscFunctionReturn(PETSC_SUCCESS);
1640: }

1642: /*@
1643:   MatSetValuesRow - Inserts a row (block row for `MATBAIJ` matrices) of nonzero
1644:   values into a matrix

1646:   Not Collective

1648:   Input Parameters:
1649: + mat - the matrix
1650: . row - the (block) row to set
1651: - v   - a logically two-dimensional (column major) array of values for  block matrices with blocksize larger than one, otherwise a one dimensional array of values

1653:   Level: advanced

1655:   Notes:
1656:   The values, `v`, are column-oriented for the block version.

1658:   All the nonzeros in `row` must be provided

1660:   THE MATRIX MUST HAVE PREVIOUSLY HAD ITS COLUMN INDICES SET. IT IS RARE THAT THIS ROUTINE IS USED, usually `MatSetValues()` is used.

1662:   `row` must belong to this process

1664: .seealso: [](ch_matrices), `Mat`, `MatSetValues()`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`,
1665:           `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`
1666: @*/
1667: PetscErrorCode MatSetValuesRow(Mat mat, PetscInt row, const PetscScalar v[])
1668: {
1669:   PetscFunctionBeginHot;
1672:   MatCheckPreallocated(mat, 1);
1673:   PetscAssertPointer(v, 3);
1674:   PetscCheck(mat->insertmode != ADD_VALUES, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add and insert values");
1675:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
1676:   mat->insertmode = INSERT_VALUES;

1678:   if (mat->assembled) {
1679:     mat->was_assembled = PETSC_TRUE;
1680:     mat->assembled     = PETSC_FALSE;
1681:   }
1682:   PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0));
1683:   PetscUseTypeMethod(mat, setvaluesrow, row, v);
1684:   PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0));
1685:   PetscFunctionReturn(PETSC_SUCCESS);
1686: }

1688: // PetscClangLinter pragma disable: -fdoc-section-header-unknown
1689: /*@
1690:   MatSetValuesStencil - Inserts or adds a block of values into a matrix.
1691:   Using structured grid indexing

1693:   Not Collective

1695:   Input Parameters:
1696: + mat  - the matrix
1697: . m    - number of rows being entered
1698: . idxm - grid coordinates (and component number when dof > 1) for matrix rows being entered
1699: . n    - number of columns being entered
1700: . idxn - grid coordinates (and component number when dof > 1) for matrix columns being entered
1701: . v    - a logically two-dimensional array of values
1702: - addv - either `ADD_VALUES` to add to existing entries at that location or `INSERT_VALUES` to replace existing entries with new values

1704:   Level: beginner

1706:   Notes:
1707:   By default the values, `v`, are row-oriented.  See `MatSetOption()` for other options.

1709:   Calls to `MatSetValuesStencil()` with the `INSERT_VALUES` and `ADD_VALUES`
1710:   options cannot be mixed without intervening calls to the assembly
1711:   routines.

1713:   The grid coordinates are across the entire grid, not just the local portion

1715:   `MatSetValuesStencil()` uses 0-based row and column numbers in Fortran
1716:   as well as in C.

1718:   For setting/accessing vector values via array coordinates you can use the `DMDAVecGetArray()` routine

1720:   In order to use this routine you must either obtain the matrix with `DMCreateMatrix()`
1721:   or call `MatSetLocalToGlobalMapping()` and `MatSetStencil()` first.

1723:   The columns and rows in the stencil passed in MUST be contained within the
1724:   ghost region of the given process as set with DMDACreateXXX() or `MatSetStencil()`. For example,
1725:   if you create a `DMDA` with an overlap of one grid level and on a particular process its first
1726:   local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the
1727:   first i index you can use in your column and row indices in `MatSetStencil()` is 5.

1729:   For periodic boundary conditions use negative indices for values to the left (below 0; that are to be
1730:   obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one
1731:   etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the
1732:   `DM_BOUNDARY_PERIODIC` boundary type.

1734:   For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have
1735:   a single value per point) you can skip filling those indices.

1737:   Inspired by the structured grid interface to the HYPRE package
1738:   (https://computation.llnl.gov/projects/hypre-scalable-linear-solvers-multigrid-methods)

1740:   Efficiency Alert:
1741:   The routine `MatSetValuesBlockedStencil()` may offer much better efficiency
1742:   for users of block sparse formats (`MATSEQBAIJ` and `MATMPIBAIJ`).

1744:   Fortran Note:
1745:   `idxm` and `idxn` should be declared as
1746: $     MatStencil idxm(4,m),idxn(4,n)
1747:   and the values inserted using
1748: .vb
1749:     idxm(MatStencil_i,1) = i
1750:     idxm(MatStencil_j,1) = j
1751:     idxm(MatStencil_k,1) = k
1752:     idxm(MatStencil_c,1) = c
1753:     etc
1754: .ve

1756: .seealso: [](ch_matrices), `Mat`, `DMDA`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`
1757:           `MatSetValues()`, `MatSetValuesBlockedStencil()`, `MatSetStencil()`, `DMCreateMatrix()`, `DMDAVecGetArray()`, `MatStencil`
1758: @*/
1759: PetscErrorCode MatSetValuesStencil(Mat mat, PetscInt m, const MatStencil idxm[], PetscInt n, const MatStencil idxn[], const PetscScalar v[], InsertMode addv)
1760: {
1761:   PetscInt  buf[8192], *bufm = NULL, *bufn = NULL, *jdxm, *jdxn;
1762:   PetscInt  j, i, dim = mat->stencil.dim, *dims = mat->stencil.dims + 1, tmp;
1763:   PetscInt *starts = mat->stencil.starts, *dxm = (PetscInt *)idxm, *dxn = (PetscInt *)idxn, sdim = dim - (1 - (PetscInt)mat->stencil.noc);

1765:   PetscFunctionBegin;
1766:   if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */
1769:   PetscAssertPointer(idxm, 3);
1770:   PetscAssertPointer(idxn, 5);

1772:   if ((m + n) <= (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) {
1773:     jdxm = buf;
1774:     jdxn = buf + m;
1775:   } else {
1776:     PetscCall(PetscMalloc2(m, &bufm, n, &bufn));
1777:     jdxm = bufm;
1778:     jdxn = bufn;
1779:   }
1780:   for (i = 0; i < m; i++) {
1781:     for (j = 0; j < 3 - sdim; j++) dxm++;
1782:     tmp = *dxm++ - starts[0];
1783:     for (j = 0; j < dim - 1; j++) {
1784:       if ((*dxm++ - starts[j + 1]) < 0 || tmp < 0) tmp = -1;
1785:       else tmp = tmp * dims[j] + *(dxm - 1) - starts[j + 1];
1786:     }
1787:     if (mat->stencil.noc) dxm++;
1788:     jdxm[i] = tmp;
1789:   }
1790:   for (i = 0; i < n; i++) {
1791:     for (j = 0; j < 3 - sdim; j++) dxn++;
1792:     tmp = *dxn++ - starts[0];
1793:     for (j = 0; j < dim - 1; j++) {
1794:       if ((*dxn++ - starts[j + 1]) < 0 || tmp < 0) tmp = -1;
1795:       else tmp = tmp * dims[j] + *(dxn - 1) - starts[j + 1];
1796:     }
1797:     if (mat->stencil.noc) dxn++;
1798:     jdxn[i] = tmp;
1799:   }
1800:   PetscCall(MatSetValuesLocal(mat, m, jdxm, n, jdxn, v, addv));
1801:   PetscCall(PetscFree2(bufm, bufn));
1802:   PetscFunctionReturn(PETSC_SUCCESS);
1803: }

1805: /*@
1806:   MatSetValuesBlockedStencil - Inserts or adds a block of values into a matrix.
1807:   Using structured grid indexing

1809:   Not Collective

1811:   Input Parameters:
1812: + mat  - the matrix
1813: . m    - number of rows being entered
1814: . idxm - grid coordinates for matrix rows being entered
1815: . n    - number of columns being entered
1816: . idxn - grid coordinates for matrix columns being entered
1817: . v    - a logically two-dimensional array of values
1818: - addv - either `ADD_VALUES` to add to existing entries or `INSERT_VALUES` to replace existing entries with new values

1820:   Level: beginner

1822:   Notes:
1823:   By default the values, `v`, are row-oriented and unsorted.
1824:   See `MatSetOption()` for other options.

1826:   Calls to `MatSetValuesBlockedStencil()` with the `INSERT_VALUES` and `ADD_VALUES`
1827:   options cannot be mixed without intervening calls to the assembly
1828:   routines.

1830:   The grid coordinates are across the entire grid, not just the local portion

1832:   `MatSetValuesBlockedStencil()` uses 0-based row and column numbers in Fortran
1833:   as well as in C.

1835:   For setting/accessing vector values via array coordinates you can use the `DMDAVecGetArray()` routine

1837:   In order to use this routine you must either obtain the matrix with `DMCreateMatrix()`
1838:   or call `MatSetBlockSize()`, `MatSetLocalToGlobalMapping()` and `MatSetStencil()` first.

1840:   The columns and rows in the stencil passed in MUST be contained within the
1841:   ghost region of the given process as set with DMDACreateXXX() or `MatSetStencil()`. For example,
1842:   if you create a `DMDA` with an overlap of one grid level and on a particular process its first
1843:   local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the
1844:   first i index you can use in your column and row indices in `MatSetStencil()` is 5.

1846:   Negative indices may be passed in idxm and idxn, these rows and columns are
1847:   simply ignored. This allows easily inserting element stiffness matrices
1848:   with homogeneous Dirichlet boundary conditions that you don't want represented
1849:   in the matrix.

1851:   Inspired by the structured grid interface to the HYPRE package
1852:   (https://computation.llnl.gov/projects/hypre-scalable-linear-solvers-multigrid-methods)

1854:   Fortran Note:
1855:   `idxm` and `idxn` should be declared as
1856: $     MatStencil idxm(4,m),idxn(4,n)
1857:   and the values inserted using
1858: .vb
1859:     idxm(MatStencil_i,1) = i
1860:     idxm(MatStencil_j,1) = j
1861:     idxm(MatStencil_k,1) = k
1862:    etc
1863: .ve

1865: .seealso: [](ch_matrices), `Mat`, `DMDA`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`
1866:           `MatSetValues()`, `MatSetValuesStencil()`, `MatSetStencil()`, `DMCreateMatrix()`, `DMDAVecGetArray()`, `MatStencil`,
1867:           `MatSetBlockSize()`, `MatSetLocalToGlobalMapping()`
1868: @*/
1869: PetscErrorCode MatSetValuesBlockedStencil(Mat mat, PetscInt m, const MatStencil idxm[], PetscInt n, const MatStencil idxn[], const PetscScalar v[], InsertMode addv)
1870: {
1871:   PetscInt  buf[8192], *bufm = NULL, *bufn = NULL, *jdxm, *jdxn;
1872:   PetscInt  j, i, dim = mat->stencil.dim, *dims = mat->stencil.dims + 1, tmp;
1873:   PetscInt *starts = mat->stencil.starts, *dxm = (PetscInt *)idxm, *dxn = (PetscInt *)idxn, sdim = dim - (1 - (PetscInt)mat->stencil.noc);

1875:   PetscFunctionBegin;
1876:   if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */
1879:   PetscAssertPointer(idxm, 3);
1880:   PetscAssertPointer(idxn, 5);
1881:   PetscAssertPointer(v, 6);

1883:   if ((m + n) <= (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) {
1884:     jdxm = buf;
1885:     jdxn = buf + m;
1886:   } else {
1887:     PetscCall(PetscMalloc2(m, &bufm, n, &bufn));
1888:     jdxm = bufm;
1889:     jdxn = bufn;
1890:   }
1891:   for (i = 0; i < m; i++) {
1892:     for (j = 0; j < 3 - sdim; j++) dxm++;
1893:     tmp = *dxm++ - starts[0];
1894:     for (j = 0; j < sdim - 1; j++) {
1895:       if ((*dxm++ - starts[j + 1]) < 0 || tmp < 0) tmp = -1;
1896:       else tmp = tmp * dims[j] + *(dxm - 1) - starts[j + 1];
1897:     }
1898:     dxm++;
1899:     jdxm[i] = tmp;
1900:   }
1901:   for (i = 0; i < n; i++) {
1902:     for (j = 0; j < 3 - sdim; j++) dxn++;
1903:     tmp = *dxn++ - starts[0];
1904:     for (j = 0; j < sdim - 1; j++) {
1905:       if ((*dxn++ - starts[j + 1]) < 0 || tmp < 0) tmp = -1;
1906:       else tmp = tmp * dims[j] + *(dxn - 1) - starts[j + 1];
1907:     }
1908:     dxn++;
1909:     jdxn[i] = tmp;
1910:   }
1911:   PetscCall(MatSetValuesBlockedLocal(mat, m, jdxm, n, jdxn, v, addv));
1912:   PetscCall(PetscFree2(bufm, bufn));
1913:   PetscFunctionReturn(PETSC_SUCCESS);
1914: }

1916: /*@
1917:   MatSetStencil - Sets the grid information for setting values into a matrix via
1918:   `MatSetValuesStencil()`

1920:   Not Collective

1922:   Input Parameters:
1923: + mat    - the matrix
1924: . dim    - dimension of the grid 1, 2, or 3
1925: . dims   - number of grid points in x, y, and z direction, including ghost points on your processor
1926: . starts - starting point of ghost nodes on your processor in x, y, and z direction
1927: - dof    - number of degrees of freedom per node

1929:   Level: beginner

1931:   Notes:
1932:   Inspired by the structured grid interface to the HYPRE package
1933:   (www.llnl.gov/CASC/hyper)

1935:   For matrices generated with `DMCreateMatrix()` this routine is automatically called and so not needed by the
1936:   user.

1938: .seealso: [](ch_matrices), `Mat`, `MatStencil`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`
1939:           `MatSetValues()`, `MatSetValuesBlockedStencil()`, `MatSetValuesStencil()`
1940: @*/
1941: PetscErrorCode MatSetStencil(Mat mat, PetscInt dim, const PetscInt dims[], const PetscInt starts[], PetscInt dof)
1942: {
1943:   PetscFunctionBegin;
1945:   PetscAssertPointer(dims, 3);
1946:   PetscAssertPointer(starts, 4);

1948:   mat->stencil.dim = dim + (dof > 1);
1949:   for (PetscInt i = 0; i < dim; i++) {
1950:     mat->stencil.dims[i]   = dims[dim - i - 1]; /* copy the values in backwards */
1951:     mat->stencil.starts[i] = starts[dim - i - 1];
1952:   }
1953:   mat->stencil.dims[dim]   = dof;
1954:   mat->stencil.starts[dim] = 0;
1955:   mat->stencil.noc         = (PetscBool)(dof == 1);
1956:   PetscFunctionReturn(PETSC_SUCCESS);
1957: }

1959: /*@
1960:   MatSetValuesBlocked - Inserts or adds a block of values into a matrix.

1962:   Not Collective

1964:   Input Parameters:
1965: + mat  - the matrix
1966: . v    - a logically two-dimensional array of values
1967: . m    - the number of block rows
1968: . idxm - the global block indices
1969: . n    - the number of block columns
1970: . idxn - the global block indices
1971: - addv - either `ADD_VALUES` to add values to any existing entries, or `INSERT_VALUES` replaces existing entries with new values

1973:   Level: intermediate

1975:   Notes:
1976:   If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call
1977:   MatXXXXSetPreallocation() or `MatSetUp()` before using this routine.

1979:   The `m` and `n` count the NUMBER of blocks in the row direction and column direction,
1980:   NOT the total number of rows/columns; for example, if the block size is 2 and
1981:   you are passing in values for rows 2,3,4,5  then `m` would be 2 (not 4).
1982:   The values in `idxm` would be 1 2; that is the first index for each block divided by
1983:   the block size.

1985:   You must call `MatSetBlockSize()` when constructing this matrix (before
1986:   preallocating it).

1988:   By default the values, `v`, are row-oriented, so the layout of
1989:   `v` is the same as for `MatSetValues()`. See `MatSetOption()` for other options.

1991:   Calls to `MatSetValuesBlocked()` with the `INSERT_VALUES` and `ADD_VALUES`
1992:   options cannot be mixed without intervening calls to the assembly
1993:   routines.

1995:   `MatSetValuesBlocked()` uses 0-based row and column numbers in Fortran
1996:   as well as in C.

1998:   Negative indices may be passed in `idxm` and `idxn`, these rows and columns are
1999:   simply ignored. This allows easily inserting element stiffness matrices
2000:   with homogeneous Dirichlet boundary conditions that you don't want represented
2001:   in the matrix.

2003:   Each time an entry is set within a sparse matrix via `MatSetValues()`,
2004:   internal searching must be done to determine where to place the
2005:   data in the matrix storage space.  By instead inserting blocks of
2006:   entries via `MatSetValuesBlocked()`, the overhead of matrix assembly is
2007:   reduced.

2009:   Example:
2010: .vb
2011:    Suppose m=n=2 and block size(bs) = 2 The array is

2013:    1  2  | 3  4
2014:    5  6  | 7  8
2015:    - - - | - - -
2016:    9  10 | 11 12
2017:    13 14 | 15 16

2019:    v[] should be passed in like
2020:    v[] = [1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16]

2022:   If you are not using row-oriented storage of v (that is you called MatSetOption(mat,MAT_ROW_ORIENTED,PETSC_FALSE)) then
2023:    v[] = [1,5,9,13,2,6,10,14,3,7,11,15,4,8,12,16]
2024: .ve

2026:   Fortran Notes:
2027:   If any of `idmx`, `idxn`, and `v` are scalars pass them using, for example,
2028: .vb
2029:   MatSetValuesBlocked(mat, one, [idxm], one, [idxn], [v], INSERT_VALUES)
2030: .ve

2032:   If `v` is a two-dimensional array use `reshape()` to pass it as a one dimensional array

2034: .seealso: [](ch_matrices), `Mat`, `MatSetBlockSize()`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValues()`, `MatSetValuesBlockedLocal()`
2035: @*/
2036: PetscErrorCode MatSetValuesBlocked(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], const PetscScalar v[], InsertMode addv)
2037: {
2038:   PetscFunctionBeginHot;
2041:   if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */
2042:   PetscAssertPointer(idxm, 3);
2043:   PetscAssertPointer(idxn, 5);
2044:   MatCheckPreallocated(mat, 1);
2045:   if (mat->insertmode == NOT_SET_VALUES) mat->insertmode = addv;
2046:   else PetscCheck(mat->insertmode == addv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add values and insert values");
2047:   if (PetscDefined(USE_DEBUG)) {
2048:     PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2049:     PetscCheck(mat->ops->setvaluesblocked || mat->ops->setvalues, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name);
2050:   }
2051:   if (PetscDefined(USE_DEBUG)) {
2052:     PetscInt rbs, cbs, M, N, i;
2053:     PetscCall(MatGetBlockSizes(mat, &rbs, &cbs));
2054:     PetscCall(MatGetSize(mat, &M, &N));
2055:     for (i = 0; i < m; i++) PetscCheck(idxm[i] * rbs < M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row block %" PetscInt_FMT " contains an index %" PetscInt_FMT "*%" PetscInt_FMT " greater than row length %" PetscInt_FMT, i, idxm[i], rbs, M);
2056:     for (i = 0; i < n; i++)
2057:       PetscCheck(idxn[i] * cbs < N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column block %" PetscInt_FMT " contains an index %" PetscInt_FMT "*%" PetscInt_FMT " greater than column length %" PetscInt_FMT, i, idxn[i], cbs, N);
2058:   }
2059:   if (mat->assembled) {
2060:     mat->was_assembled = PETSC_TRUE;
2061:     mat->assembled     = PETSC_FALSE;
2062:   }
2063:   PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0));
2064:   if (mat->ops->setvaluesblocked) {
2065:     PetscUseTypeMethod(mat, setvaluesblocked, m, idxm, n, idxn, v, addv);
2066:   } else {
2067:     PetscInt buf[8192], *bufr = NULL, *bufc = NULL, *iidxm, *iidxn;
2068:     PetscInt i, j, bs, cbs;

2070:     PetscCall(MatGetBlockSizes(mat, &bs, &cbs));
2071:     if ((m * bs + n * cbs) <= (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) {
2072:       iidxm = buf;
2073:       iidxn = buf + m * bs;
2074:     } else {
2075:       PetscCall(PetscMalloc2(m * bs, &bufr, n * cbs, &bufc));
2076:       iidxm = bufr;
2077:       iidxn = bufc;
2078:     }
2079:     for (i = 0; i < m; i++) {
2080:       for (j = 0; j < bs; j++) iidxm[i * bs + j] = bs * idxm[i] + j;
2081:     }
2082:     if (m != n || bs != cbs || idxm != idxn) {
2083:       for (i = 0; i < n; i++) {
2084:         for (j = 0; j < cbs; j++) iidxn[i * cbs + j] = cbs * idxn[i] + j;
2085:       }
2086:     } else iidxn = iidxm;
2087:     PetscCall(MatSetValues(mat, m * bs, iidxm, n * cbs, iidxn, v, addv));
2088:     PetscCall(PetscFree2(bufr, bufc));
2089:   }
2090:   PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0));
2091:   PetscFunctionReturn(PETSC_SUCCESS);
2092: }

2094: /*@
2095:   MatGetValues - Gets a block of local values from a matrix.

2097:   Not Collective; can only return values that are owned by the give process

2099:   Input Parameters:
2100: + mat  - the matrix
2101: . v    - a logically two-dimensional array for storing the values
2102: . m    - the number of rows
2103: . idxm - the  global indices of the rows
2104: . n    - the number of columns
2105: - idxn - the global indices of the columns

2107:   Level: advanced

2109:   Notes:
2110:   The user must allocate space (m*n `PetscScalar`s) for the values, `v`.
2111:   The values, `v`, are then returned in a row-oriented format,
2112:   analogous to that used by default in `MatSetValues()`.

2114:   `MatGetValues()` uses 0-based row and column numbers in
2115:   Fortran as well as in C.

2117:   `MatGetValues()` requires that the matrix has been assembled
2118:   with `MatAssemblyBegin()`/`MatAssemblyEnd()`.  Thus, calls to
2119:   `MatSetValues()` and `MatGetValues()` CANNOT be made in succession
2120:   without intermediate matrix assembly.

2122:   Negative row or column indices will be ignored and those locations in `v` will be
2123:   left unchanged.

2125:   For the standard row-based matrix formats, `idxm` can only contain rows owned by the requesting MPI process.
2126:   That is, rows with global index greater than or equal to rstart and less than rend where rstart and rend are obtainable
2127:   from `MatGetOwnershipRange`(mat,&rstart,&rend).

2129: .seealso: [](ch_matrices), `Mat`, `MatGetRow()`, `MatCreateSubMatrices()`, `MatSetValues()`, `MatGetOwnershipRange()`, `MatGetValuesLocal()`, `MatGetValue()`
2130: @*/
2131: PetscErrorCode MatGetValues(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[])
2132: {
2133:   PetscFunctionBegin;
2136:   if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS);
2137:   PetscAssertPointer(idxm, 3);
2138:   PetscAssertPointer(idxn, 5);
2139:   PetscAssertPointer(v, 6);
2140:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2141:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2142:   MatCheckPreallocated(mat, 1);

2144:   PetscCall(PetscLogEventBegin(MAT_GetValues, mat, 0, 0, 0));
2145:   PetscUseTypeMethod(mat, getvalues, m, idxm, n, idxn, v);
2146:   PetscCall(PetscLogEventEnd(MAT_GetValues, mat, 0, 0, 0));
2147:   PetscFunctionReturn(PETSC_SUCCESS);
2148: }

2150: /*@
2151:   MatGetValuesLocal - retrieves values from certain locations in a matrix using the local numbering of the indices
2152:   defined previously by `MatSetLocalToGlobalMapping()`

2154:   Not Collective

2156:   Input Parameters:
2157: + mat  - the matrix
2158: . nrow - number of rows
2159: . irow - the row local indices
2160: . ncol - number of columns
2161: - icol - the column local indices

2163:   Output Parameter:
2164: . y - a logically two-dimensional array of values

2166:   Level: advanced

2168:   Notes:
2169:   If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call `MatSetLocalToGlobalMapping()` before using this routine.

2171:   This routine can only return values that are owned by the requesting MPI process. That is, for standard matrix formats, rows that, in the global numbering,
2172:   are greater than or equal to rstart and less than rend where rstart and rend are obtainable from `MatGetOwnershipRange`(mat,&rstart,&rend). One can
2173:   determine if the resulting global row associated with the local row r is owned by the requesting MPI process by applying the `ISLocalToGlobalMapping` set
2174:   with `MatSetLocalToGlobalMapping()`.

2176:   Developer Note:
2177:   This is labelled with C so does not automatically generate Fortran stubs and interfaces
2178:   because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.

2180: .seealso: [](ch_matrices), `Mat`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValues()`, `MatSetLocalToGlobalMapping()`,
2181:           `MatSetValuesLocal()`, `MatGetValues()`
2182: @*/
2183: PetscErrorCode MatGetValuesLocal(Mat mat, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], PetscScalar y[])
2184: {
2185:   PetscFunctionBeginHot;
2188:   MatCheckPreallocated(mat, 1);
2189:   if (!nrow || !ncol) PetscFunctionReturn(PETSC_SUCCESS); /* no values to retrieve */
2190:   PetscAssertPointer(irow, 3);
2191:   PetscAssertPointer(icol, 5);
2192:   if (PetscDefined(USE_DEBUG)) {
2193:     PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2194:     PetscCheck(mat->ops->getvalueslocal || mat->ops->getvalues, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name);
2195:   }
2196:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2197:   PetscCall(PetscLogEventBegin(MAT_GetValues, mat, 0, 0, 0));
2198:   if (mat->ops->getvalueslocal) PetscUseTypeMethod(mat, getvalueslocal, nrow, irow, ncol, icol, y);
2199:   else {
2200:     PetscInt buf[8192], *bufr = NULL, *bufc = NULL, *irowm, *icolm;
2201:     if ((nrow + ncol) <= (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) {
2202:       irowm = buf;
2203:       icolm = buf + nrow;
2204:     } else {
2205:       PetscCall(PetscMalloc2(nrow, &bufr, ncol, &bufc));
2206:       irowm = bufr;
2207:       icolm = bufc;
2208:     }
2209:     PetscCheck(mat->rmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MatGetValuesLocal() cannot proceed without local-to-global row mapping (See MatSetLocalToGlobalMapping()).");
2210:     PetscCheck(mat->cmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MatGetValuesLocal() cannot proceed without local-to-global column mapping (See MatSetLocalToGlobalMapping()).");
2211:     PetscCall(ISLocalToGlobalMappingApply(mat->rmap->mapping, nrow, irow, irowm));
2212:     PetscCall(ISLocalToGlobalMappingApply(mat->cmap->mapping, ncol, icol, icolm));
2213:     PetscCall(MatGetValues(mat, nrow, irowm, ncol, icolm, y));
2214:     PetscCall(PetscFree2(bufr, bufc));
2215:   }
2216:   PetscCall(PetscLogEventEnd(MAT_GetValues, mat, 0, 0, 0));
2217:   PetscFunctionReturn(PETSC_SUCCESS);
2218: }

2220: /*@
2221:   MatSetValuesBatch - Adds (`ADD_VALUES`) many blocks of values into a matrix at once. The blocks must all be square and
2222:   the same size. Currently, this can only be called once and creates the given matrix.

2224:   Not Collective

2226:   Input Parameters:
2227: + mat  - the matrix
2228: . nb   - the number of blocks
2229: . bs   - the number of rows (and columns) in each block
2230: . rows - a concatenation of the rows for each block
2231: - v    - a concatenation of logically two-dimensional arrays of values

2233:   Level: advanced

2235:   Notes:
2236:   `MatSetPreallocationCOO()` and `MatSetValuesCOO()` may be a better way to provide the values

2238:   In the future, we will extend this routine to handle rectangular blocks, and to allow multiple calls for a given matrix.

2240: .seealso: [](ch_matrices), `Mat`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`,
2241:           `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`, `MatSetValues()`, `MatSetPreallocationCOO()`, `MatSetValuesCOO()`
2242: @*/
2243: PetscErrorCode MatSetValuesBatch(Mat mat, PetscInt nb, PetscInt bs, PetscInt rows[], const PetscScalar v[])
2244: {
2245:   PetscFunctionBegin;
2248:   PetscAssertPointer(rows, 4);
2249:   PetscAssertPointer(v, 5);
2250:   PetscAssert(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");

2252:   PetscCall(PetscLogEventBegin(MAT_SetValuesBatch, mat, 0, 0, 0));
2253:   if (mat->ops->setvaluesbatch) PetscUseTypeMethod(mat, setvaluesbatch, nb, bs, rows, v);
2254:   else {
2255:     for (PetscInt b = 0; b < nb; ++b) PetscCall(MatSetValues(mat, bs, &rows[b * bs], bs, &rows[b * bs], &v[b * bs * bs], ADD_VALUES));
2256:   }
2257:   PetscCall(PetscLogEventEnd(MAT_SetValuesBatch, mat, 0, 0, 0));
2258:   PetscFunctionReturn(PETSC_SUCCESS);
2259: }

2261: /*@
2262:   MatSetLocalToGlobalMapping - Sets a local-to-global numbering for use by
2263:   the routine `MatSetValuesLocal()` to allow users to insert matrix entries
2264:   using a local (per-processor) numbering.

2266:   Not Collective

2268:   Input Parameters:
2269: + x        - the matrix
2270: . rmapping - row mapping created with `ISLocalToGlobalMappingCreate()` or `ISLocalToGlobalMappingCreateIS()`
2271: - cmapping - column mapping

2273:   Level: intermediate

2275:   Note:
2276:   If the matrix is obtained with `DMCreateMatrix()` then this may already have been called on the matrix

2278: .seealso: [](ch_matrices), `Mat`, `DM`, `DMCreateMatrix()`, `MatGetLocalToGlobalMapping()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValues()`, `MatSetValuesLocal()`, `MatGetValuesLocal()`
2279: @*/
2280: PetscErrorCode MatSetLocalToGlobalMapping(Mat x, ISLocalToGlobalMapping rmapping, ISLocalToGlobalMapping cmapping)
2281: {
2282:   PetscFunctionBegin;
2287:   if (x->ops->setlocaltoglobalmapping) PetscUseTypeMethod(x, setlocaltoglobalmapping, rmapping, cmapping);
2288:   else {
2289:     PetscCall(PetscLayoutSetISLocalToGlobalMapping(x->rmap, rmapping));
2290:     PetscCall(PetscLayoutSetISLocalToGlobalMapping(x->cmap, cmapping));
2291:   }
2292:   PetscFunctionReturn(PETSC_SUCCESS);
2293: }

2295: /*@
2296:   MatGetLocalToGlobalMapping - Gets the local-to-global numbering set by `MatSetLocalToGlobalMapping()`

2298:   Not Collective

2300:   Input Parameter:
2301: . A - the matrix

2303:   Output Parameters:
2304: + rmapping - row mapping
2305: - cmapping - column mapping

2307:   Level: advanced

2309: .seealso: [](ch_matrices), `Mat`, `MatSetLocalToGlobalMapping()`, `MatSetValuesLocal()`
2310: @*/
2311: PetscErrorCode MatGetLocalToGlobalMapping(Mat A, ISLocalToGlobalMapping *rmapping, ISLocalToGlobalMapping *cmapping)
2312: {
2313:   PetscFunctionBegin;
2316:   if (rmapping) {
2317:     PetscAssertPointer(rmapping, 2);
2318:     *rmapping = A->rmap->mapping;
2319:   }
2320:   if (cmapping) {
2321:     PetscAssertPointer(cmapping, 3);
2322:     *cmapping = A->cmap->mapping;
2323:   }
2324:   PetscFunctionReturn(PETSC_SUCCESS);
2325: }

2327: /*@
2328:   MatSetLayouts - Sets the `PetscLayout` objects for rows and columns of a matrix

2330:   Logically Collective

2332:   Input Parameters:
2333: + A    - the matrix
2334: . rmap - row layout
2335: - cmap - column layout

2337:   Level: advanced

2339:   Note:
2340:   The `PetscLayout` objects are usually created automatically for the matrix so this routine rarely needs to be called.

2342: .seealso: [](ch_matrices), `Mat`, `PetscLayout`, `MatCreateVecs()`, `MatGetLocalToGlobalMapping()`, `MatGetLayouts()`
2343: @*/
2344: PetscErrorCode MatSetLayouts(Mat A, PetscLayout rmap, PetscLayout cmap)
2345: {
2346:   PetscFunctionBegin;
2348:   PetscCall(PetscLayoutReference(rmap, &A->rmap));
2349:   PetscCall(PetscLayoutReference(cmap, &A->cmap));
2350:   PetscFunctionReturn(PETSC_SUCCESS);
2351: }

2353: /*@
2354:   MatGetLayouts - Gets the `PetscLayout` objects for rows and columns

2356:   Not Collective

2358:   Input Parameter:
2359: . A - the matrix

2361:   Output Parameters:
2362: + rmap - row layout
2363: - cmap - column layout

2365:   Level: advanced

2367: .seealso: [](ch_matrices), `Mat`, [Matrix Layouts](sec_matlayout), `PetscLayout`, `MatCreateVecs()`, `MatGetLocalToGlobalMapping()`, `MatSetLayouts()`
2368: @*/
2369: PetscErrorCode MatGetLayouts(Mat A, PetscLayout *rmap, PetscLayout *cmap)
2370: {
2371:   PetscFunctionBegin;
2374:   if (rmap) {
2375:     PetscAssertPointer(rmap, 2);
2376:     *rmap = A->rmap;
2377:   }
2378:   if (cmap) {
2379:     PetscAssertPointer(cmap, 3);
2380:     *cmap = A->cmap;
2381:   }
2382:   PetscFunctionReturn(PETSC_SUCCESS);
2383: }

2385: /*@
2386:   MatSetValuesLocal - Inserts or adds values into certain locations of a matrix,
2387:   using a local numbering of the rows and columns.

2389:   Not Collective

2391:   Input Parameters:
2392: + mat  - the matrix
2393: . nrow - number of rows
2394: . irow - the row local indices
2395: . ncol - number of columns
2396: . icol - the column local indices
2397: . y    - a logically two-dimensional array of values
2398: - addv - either `INSERT_VALUES` to add values to any existing entries, or `INSERT_VALUES` to replace existing entries with new values

2400:   Level: intermediate

2402:   Notes:
2403:   If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call `MatSetLocalToGlobalMapping()` before using this routine

2405:   Calls to `MatSetValuesLocal()` with the `INSERT_VALUES` and `ADD_VALUES`
2406:   options cannot be mixed without intervening calls to the assembly
2407:   routines.

2409:   These values may be cached, so `MatAssemblyBegin()` and `MatAssemblyEnd()`
2410:   MUST be called after all calls to `MatSetValuesLocal()` have been completed.

2412:   Fortran Notes:
2413:   If any of `irow`, `icol`, and `y` are scalars pass them using, for example,
2414: .vb
2415:   MatSetValuesLocal(mat, one, [irow], one, [icol], [y], INSERT_VALUES)
2416: .ve

2418:   If `y` is a two-dimensional array use `reshape()` to pass it as a one dimensional array

2420:   Developer Note:
2421:   This is labeled with C so does not automatically generate Fortran stubs and interfaces
2422:   because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.

2424: .seealso: [](ch_matrices), `Mat`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValues()`, `MatSetLocalToGlobalMapping()`,
2425:           `MatGetValuesLocal()`
2426: @*/
2427: PetscErrorCode MatSetValuesLocal(Mat mat, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv)
2428: {
2429:   PetscFunctionBeginHot;
2432:   MatCheckPreallocated(mat, 1);
2433:   if (!nrow || !ncol) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */
2434:   PetscAssertPointer(irow, 3);
2435:   PetscAssertPointer(icol, 5);
2436:   if (mat->insertmode == NOT_SET_VALUES) mat->insertmode = addv;
2437:   else PetscCheck(mat->insertmode == addv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add values and insert values");
2438:   if (PetscDefined(USE_DEBUG)) {
2439:     PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2440:     PetscCheck(mat->ops->setvalueslocal || mat->ops->setvalues, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name);
2441:   }

2443:   if (mat->assembled) {
2444:     mat->was_assembled = PETSC_TRUE;
2445:     mat->assembled     = PETSC_FALSE;
2446:   }
2447:   PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0));
2448:   if (mat->ops->setvalueslocal) PetscUseTypeMethod(mat, setvalueslocal, nrow, irow, ncol, icol, y, addv);
2449:   else {
2450:     PetscInt        buf[8192], *bufr = NULL, *bufc = NULL;
2451:     const PetscInt *irowm, *icolm;

2453:     if ((!mat->rmap->mapping && !mat->cmap->mapping) || (nrow + ncol) <= (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) {
2454:       bufr  = buf;
2455:       bufc  = buf + nrow;
2456:       irowm = bufr;
2457:       icolm = bufc;
2458:     } else {
2459:       PetscCall(PetscMalloc2(nrow, &bufr, ncol, &bufc));
2460:       irowm = bufr;
2461:       icolm = bufc;
2462:     }
2463:     if (mat->rmap->mapping) PetscCall(ISLocalToGlobalMappingApply(mat->rmap->mapping, nrow, irow, bufr));
2464:     else irowm = irow;
2465:     if (mat->cmap->mapping) {
2466:       if (mat->cmap->mapping != mat->rmap->mapping || ncol != nrow || icol != irow) {
2467:         PetscCall(ISLocalToGlobalMappingApply(mat->cmap->mapping, ncol, icol, bufc));
2468:       } else icolm = irowm;
2469:     } else icolm = icol;
2470:     PetscCall(MatSetValues(mat, nrow, irowm, ncol, icolm, y, addv));
2471:     if (bufr != buf) PetscCall(PetscFree2(bufr, bufc));
2472:   }
2473:   PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0));
2474:   PetscFunctionReturn(PETSC_SUCCESS);
2475: }

2477: /*@
2478:   MatSetValuesBlockedLocal - Inserts or adds values into certain locations of a matrix,
2479:   using a local ordering of the nodes a block at a time.

2481:   Not Collective

2483:   Input Parameters:
2484: + mat  - the matrix
2485: . nrow - number of rows
2486: . irow - the row local indices
2487: . ncol - number of columns
2488: . icol - the column local indices
2489: . y    - a logically two-dimensional array of values
2490: - addv - either `ADD_VALUES` to add values to any existing entries, or `INSERT_VALUES` to replace existing entries with new values

2492:   Level: intermediate

2494:   Notes:
2495:   If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call `MatSetBlockSize()` and `MatSetLocalToGlobalMapping()`
2496:   before using this routineBefore calling `MatSetValuesLocal()`, the user must first set the

2498:   Calls to `MatSetValuesBlockedLocal()` with the `INSERT_VALUES` and `ADD_VALUES`
2499:   options cannot be mixed without intervening calls to the assembly
2500:   routines.

2502:   These values may be cached, so `MatAssemblyBegin()` and `MatAssemblyEnd()`
2503:   MUST be called after all calls to `MatSetValuesBlockedLocal()` have been completed.

2505:   Fortran Notes:
2506:   If any of `irow`, `icol`, and `y` are scalars pass them using, for example,
2507: .vb
2508:   MatSetValuesBlockedLocal(mat, one, [irow], one, [icol], [y], INSERT_VALUES)
2509: .ve

2511:   If `y` is a two-dimensional array use `reshape()` to pass it as a one dimensional array

2513:   Developer Note:
2514:   This is labeled with C so does not automatically generate Fortran stubs and interfaces
2515:   because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays.

2517: .seealso: [](ch_matrices), `Mat`, `MatSetBlockSize()`, `MatSetLocalToGlobalMapping()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`,
2518:           `MatSetValuesLocal()`, `MatSetValuesBlocked()`
2519: @*/
2520: PetscErrorCode MatSetValuesBlockedLocal(Mat mat, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv)
2521: {
2522:   PetscFunctionBeginHot;
2525:   MatCheckPreallocated(mat, 1);
2526:   if (!nrow || !ncol) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */
2527:   PetscAssertPointer(irow, 3);
2528:   PetscAssertPointer(icol, 5);
2529:   if (mat->insertmode == NOT_SET_VALUES) mat->insertmode = addv;
2530:   else PetscCheck(mat->insertmode == addv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add values and insert values");
2531:   if (PetscDefined(USE_DEBUG)) {
2532:     PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2533:     PetscCheck(mat->ops->setvaluesblockedlocal || mat->ops->setvaluesblocked || mat->ops->setvalueslocal || mat->ops->setvalues, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name);
2534:   }

2536:   if (mat->assembled) {
2537:     mat->was_assembled = PETSC_TRUE;
2538:     mat->assembled     = PETSC_FALSE;
2539:   }
2540:   if (PetscUnlikelyDebug(mat->rmap->mapping)) { /* Condition on the mapping existing, because MatSetValuesBlockedLocal_IS does not require it to be set. */
2541:     PetscInt irbs, rbs;
2542:     PetscCall(MatGetBlockSizes(mat, &rbs, NULL));
2543:     PetscCall(ISLocalToGlobalMappingGetBlockSize(mat->rmap->mapping, &irbs));
2544:     PetscCheck(rbs == irbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Different row block sizes! mat %" PetscInt_FMT ", row l2g map %" PetscInt_FMT, rbs, irbs);
2545:   }
2546:   if (PetscUnlikelyDebug(mat->cmap->mapping)) {
2547:     PetscInt icbs, cbs;
2548:     PetscCall(MatGetBlockSizes(mat, NULL, &cbs));
2549:     PetscCall(ISLocalToGlobalMappingGetBlockSize(mat->cmap->mapping, &icbs));
2550:     PetscCheck(cbs == icbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Different col block sizes! mat %" PetscInt_FMT ", col l2g map %" PetscInt_FMT, cbs, icbs);
2551:   }
2552:   PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0));
2553:   if (mat->ops->setvaluesblockedlocal) PetscUseTypeMethod(mat, setvaluesblockedlocal, nrow, irow, ncol, icol, y, addv);
2554:   else {
2555:     PetscInt        buf[8192], *bufr = NULL, *bufc = NULL;
2556:     const PetscInt *irowm, *icolm;

2558:     if ((!mat->rmap->mapping && !mat->cmap->mapping) || (nrow + ncol) <= ((PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf))) {
2559:       bufr  = buf;
2560:       bufc  = buf + nrow;
2561:       irowm = bufr;
2562:       icolm = bufc;
2563:     } else {
2564:       PetscCall(PetscMalloc2(nrow, &bufr, ncol, &bufc));
2565:       irowm = bufr;
2566:       icolm = bufc;
2567:     }
2568:     if (mat->rmap->mapping) PetscCall(ISLocalToGlobalMappingApplyBlock(mat->rmap->mapping, nrow, irow, bufr));
2569:     else irowm = irow;
2570:     if (mat->cmap->mapping) {
2571:       if (mat->cmap->mapping != mat->rmap->mapping || ncol != nrow || icol != irow) {
2572:         PetscCall(ISLocalToGlobalMappingApplyBlock(mat->cmap->mapping, ncol, icol, bufc));
2573:       } else icolm = irowm;
2574:     } else icolm = icol;
2575:     PetscCall(MatSetValuesBlocked(mat, nrow, irowm, ncol, icolm, y, addv));
2576:     if (bufr != buf) PetscCall(PetscFree2(bufr, bufc));
2577:   }
2578:   PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0));
2579:   PetscFunctionReturn(PETSC_SUCCESS);
2580: }

2582: /*@
2583:   MatMultDiagonalBlock - Computes the matrix-vector product, $y = Dx$. Where `D` is defined by the inode or block structure of the diagonal

2585:   Collective

2587:   Input Parameters:
2588: + mat - the matrix
2589: - x   - the vector to be multiplied

2591:   Output Parameter:
2592: . y - the result

2594:   Level: developer

2596:   Note:
2597:   The vectors `x` and `y` cannot be the same.  I.e., one cannot
2598:   call `MatMultDiagonalBlock`(A,y,y).

2600: .seealso: [](ch_matrices), `Mat`, `MatMult()`, `MatMultTranspose()`, `MatMultAdd()`, `MatMultTransposeAdd()`
2601: @*/
2602: PetscErrorCode MatMultDiagonalBlock(Mat mat, Vec x, Vec y)
2603: {
2604:   PetscFunctionBegin;

2610:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2611:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2612:   PetscCheck(x != y, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "x and y must be different vectors");
2613:   MatCheckPreallocated(mat, 1);

2615:   PetscUseTypeMethod(mat, multdiagonalblock, x, y);
2616:   PetscCall(PetscObjectStateIncrease((PetscObject)y));
2617:   PetscFunctionReturn(PETSC_SUCCESS);
2618: }

2620: /*@
2621:   MatMult - Computes the matrix-vector product, $y = Ax$.

2623:   Neighbor-wise Collective

2625:   Input Parameters:
2626: + mat - the matrix
2627: - x   - the vector to be multiplied

2629:   Output Parameter:
2630: . y - the result

2632:   Level: beginner

2634:   Note:
2635:   The vectors `x` and `y` cannot be the same.  I.e., one cannot
2636:   call `MatMult`(A,y,y).

2638: .seealso: [](ch_matrices), `Mat`, `MatMultTranspose()`, `MatMultAdd()`, `MatMultTransposeAdd()`
2639: @*/
2640: PetscErrorCode MatMult(Mat mat, Vec x, Vec y)
2641: {
2642:   PetscFunctionBegin;
2646:   VecCheckAssembled(x);
2648:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2649:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2650:   PetscCheck(x != y, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "x and y must be different vectors");
2651:   PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N);
2652:   PetscCheck(mat->rmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, y->map->N);
2653:   PetscCheck(mat->cmap->n == x->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->n, x->map->n);
2654:   PetscCheck(mat->rmap->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, y->map->n);
2655:   PetscCall(VecSetErrorIfLocked(y, 3));
2656:   if (mat->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
2657:   MatCheckPreallocated(mat, 1);

2659:   PetscCall(VecLockReadPush(x));
2660:   PetscCall(PetscLogEventBegin(MAT_Mult, mat, x, y, 0));
2661:   PetscUseTypeMethod(mat, mult, x, y);
2662:   PetscCall(PetscLogEventEnd(MAT_Mult, mat, x, y, 0));
2663:   if (mat->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
2664:   PetscCall(VecLockReadPop(x));
2665:   PetscFunctionReturn(PETSC_SUCCESS);
2666: }

2668: /*@
2669:   MatMultTranspose - Computes matrix transpose times a vector $y = A^T * x$.

2671:   Neighbor-wise Collective

2673:   Input Parameters:
2674: + mat - the matrix
2675: - x   - the vector to be multiplied

2677:   Output Parameter:
2678: . y - the result

2680:   Level: beginner

2682:   Notes:
2683:   The vectors `x` and `y` cannot be the same.  I.e., one cannot
2684:   call `MatMultTranspose`(A,y,y).

2686:   For complex numbers this does NOT compute the Hermitian (complex conjugate) transpose multiple,
2687:   use `MatMultHermitianTranspose()`

2689: .seealso: [](ch_matrices), `Mat`, `MatMult()`, `MatMultAdd()`, `MatMultTransposeAdd()`, `MatMultHermitianTranspose()`, `MatTranspose()`
2690: @*/
2691: PetscErrorCode MatMultTranspose(Mat mat, Vec x, Vec y)
2692: {
2693:   PetscErrorCode (*op)(Mat, Vec, Vec) = NULL;

2695:   PetscFunctionBegin;
2699:   VecCheckAssembled(x);

2702:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2703:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2704:   PetscCheck(x != y, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "x and y must be different vectors");
2705:   PetscCheck(mat->cmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, y->map->N);
2706:   PetscCheck(mat->rmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, x->map->N);
2707:   PetscCheck(mat->cmap->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->n, y->map->n);
2708:   PetscCheck(mat->rmap->n == x->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, x->map->n);
2709:   if (mat->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE));
2710:   MatCheckPreallocated(mat, 1);

2712:   if (!mat->ops->multtranspose) {
2713:     if (mat->symmetric == PETSC_BOOL3_TRUE && mat->ops->mult) op = mat->ops->mult;
2714:     PetscCheck(op, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Matrix type %s does not have a multiply transpose defined or is symmetric and does not have a multiply defined", ((PetscObject)mat)->type_name);
2715:   } else op = mat->ops->multtranspose;
2716:   PetscCall(PetscLogEventBegin(MAT_MultTranspose, mat, x, y, 0));
2717:   PetscCall(VecLockReadPush(x));
2718:   PetscCall((*op)(mat, x, y));
2719:   PetscCall(VecLockReadPop(x));
2720:   PetscCall(PetscLogEventEnd(MAT_MultTranspose, mat, x, y, 0));
2721:   PetscCall(PetscObjectStateIncrease((PetscObject)y));
2722:   if (mat->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE));
2723:   PetscFunctionReturn(PETSC_SUCCESS);
2724: }

2726: /*@
2727:   MatMultHermitianTranspose - Computes matrix Hermitian-transpose times a vector $y = A^H * x$.

2729:   Neighbor-wise Collective

2731:   Input Parameters:
2732: + mat - the matrix
2733: - x   - the vector to be multiplied

2735:   Output Parameter:
2736: . y - the result

2738:   Level: beginner

2740:   Notes:
2741:   The vectors `x` and `y` cannot be the same.  I.e., one cannot
2742:   call `MatMultHermitianTranspose`(A,y,y).

2744:   Also called the conjugate transpose, complex conjugate transpose, or adjoint.

2746:   For real numbers `MatMultTranspose()` and `MatMultHermitianTranspose()` are identical.

2748: .seealso: [](ch_matrices), `Mat`, `MatMult()`, `MatMultAdd()`, `MatMultHermitianTransposeAdd()`, `MatMultTranspose()`
2749: @*/
2750: PetscErrorCode MatMultHermitianTranspose(Mat mat, Vec x, Vec y)
2751: {
2752:   PetscFunctionBegin;

2758:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2759:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2760:   PetscCheck(x != y, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "x and y must be different vectors");
2761:   PetscCheck(mat->cmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, y->map->N);
2762:   PetscCheck(mat->rmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, x->map->N);
2763:   PetscCheck(mat->cmap->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->n, y->map->n);
2764:   PetscCheck(mat->rmap->n == x->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, x->map->n);
2765:   MatCheckPreallocated(mat, 1);

2767:   PetscCall(PetscLogEventBegin(MAT_MultHermitianTranspose, mat, x, y, 0));
2768: #if defined(PETSC_USE_COMPLEX)
2769:   if (mat->ops->multhermitiantranspose || (mat->hermitian == PETSC_BOOL3_TRUE && mat->ops->mult)) {
2770:     PetscCall(VecLockReadPush(x));
2771:     if (mat->ops->multhermitiantranspose) PetscUseTypeMethod(mat, multhermitiantranspose, x, y);
2772:     else PetscUseTypeMethod(mat, mult, x, y);
2773:     PetscCall(VecLockReadPop(x));
2774:   } else {
2775:     Vec w;
2776:     PetscCall(VecDuplicate(x, &w));
2777:     PetscCall(VecCopy(x, w));
2778:     PetscCall(VecConjugate(w));
2779:     PetscCall(MatMultTranspose(mat, w, y));
2780:     PetscCall(VecDestroy(&w));
2781:     PetscCall(VecConjugate(y));
2782:   }
2783:   PetscCall(PetscObjectStateIncrease((PetscObject)y));
2784: #else
2785:   PetscCall(MatMultTranspose(mat, x, y));
2786: #endif
2787:   PetscCall(PetscLogEventEnd(MAT_MultHermitianTranspose, mat, x, y, 0));
2788:   PetscFunctionReturn(PETSC_SUCCESS);
2789: }

2791: /*@
2792:   MatMultAdd -  Computes $v3 = v2 + A * v1$.

2794:   Neighbor-wise Collective

2796:   Input Parameters:
2797: + mat - the matrix
2798: . v1  - the vector to be multiplied by `mat`
2799: - v2  - the vector to be added to the result

2801:   Output Parameter:
2802: . v3 - the result

2804:   Level: beginner

2806:   Note:
2807:   The vectors `v1` and `v3` cannot be the same.  I.e., one cannot
2808:   call `MatMultAdd`(A,v1,v2,v1).

2810: .seealso: [](ch_matrices), `Mat`, `MatMultTranspose()`, `MatMult()`, `MatMultTransposeAdd()`
2811: @*/
2812: PetscErrorCode MatMultAdd(Mat mat, Vec v1, Vec v2, Vec v3)
2813: {
2814:   PetscFunctionBegin;

2821:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2822:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2823:   PetscCheck(mat->cmap->N == v1->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v1: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v1->map->N);
2824:   /* PetscCheck(mat->rmap->N == v2->map->N,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: global dim %" PetscInt_FMT " %" PetscInt_FMT,mat->rmap->N,v2->map->N);
2825:      PetscCheck(mat->rmap->N == v3->map->N,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: global dim %" PetscInt_FMT " %" PetscInt_FMT,mat->rmap->N,v3->map->N); */
2826:   PetscCheck(mat->rmap->n == v3->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec v3: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, v3->map->n);
2827:   PetscCheck(mat->rmap->n == v2->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec v2: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, v2->map->n);
2828:   PetscCheck(v1 != v3, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "v1 and v3 must be different vectors");
2829:   MatCheckPreallocated(mat, 1);

2831:   PetscCall(PetscLogEventBegin(MAT_MultAdd, mat, v1, v2, v3));
2832:   PetscCall(VecLockReadPush(v1));
2833:   PetscUseTypeMethod(mat, multadd, v1, v2, v3);
2834:   PetscCall(VecLockReadPop(v1));
2835:   PetscCall(PetscLogEventEnd(MAT_MultAdd, mat, v1, v2, v3));
2836:   PetscCall(PetscObjectStateIncrease((PetscObject)v3));
2837:   PetscFunctionReturn(PETSC_SUCCESS);
2838: }

2840: /*@
2841:   MatMultTransposeAdd - Computes $v3 = v2 + A^T * v1$.

2843:   Neighbor-wise Collective

2845:   Input Parameters:
2846: + mat - the matrix
2847: . v1  - the vector to be multiplied by the transpose of the matrix
2848: - v2  - the vector to be added to the result

2850:   Output Parameter:
2851: . v3 - the result

2853:   Level: beginner

2855:   Note:
2856:   The vectors `v1` and `v3` cannot be the same.  I.e., one cannot
2857:   call `MatMultTransposeAdd`(A,v1,v2,v1).

2859: .seealso: [](ch_matrices), `Mat`, `MatMultTranspose()`, `MatMultAdd()`, `MatMult()`
2860: @*/
2861: PetscErrorCode MatMultTransposeAdd(Mat mat, Vec v1, Vec v2, Vec v3)
2862: {
2863:   PetscErrorCode (*op)(Mat, Vec, Vec, Vec) = (!mat->ops->multtransposeadd && mat->symmetric) ? mat->ops->multadd : mat->ops->multtransposeadd;

2865:   PetscFunctionBegin;

2872:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2873:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2874:   PetscCheck(mat->rmap->N == v1->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v1: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, v1->map->N);
2875:   PetscCheck(mat->cmap->N == v2->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v2: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v2->map->N);
2876:   PetscCheck(mat->cmap->N == v3->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v3: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v3->map->N);
2877:   PetscCheck(v1 != v3, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "v1 and v3 must be different vectors");
2878:   PetscCheck(op, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name);
2879:   MatCheckPreallocated(mat, 1);

2881:   PetscCall(PetscLogEventBegin(MAT_MultTransposeAdd, mat, v1, v2, v3));
2882:   PetscCall(VecLockReadPush(v1));
2883:   PetscCall((*op)(mat, v1, v2, v3));
2884:   PetscCall(VecLockReadPop(v1));
2885:   PetscCall(PetscLogEventEnd(MAT_MultTransposeAdd, mat, v1, v2, v3));
2886:   PetscCall(PetscObjectStateIncrease((PetscObject)v3));
2887:   PetscFunctionReturn(PETSC_SUCCESS);
2888: }

2890: /*@
2891:   MatMultHermitianTransposeAdd - Computes $v3 = v2 + A^H * v1$.

2893:   Neighbor-wise Collective

2895:   Input Parameters:
2896: + mat - the matrix
2897: . v1  - the vector to be multiplied by the Hermitian transpose
2898: - v2  - the vector to be added to the result

2900:   Output Parameter:
2901: . v3 - the result

2903:   Level: beginner

2905:   Note:
2906:   The vectors `v1` and `v3` cannot be the same.  I.e., one cannot
2907:   call `MatMultHermitianTransposeAdd`(A,v1,v2,v1).

2909: .seealso: [](ch_matrices), `Mat`, `MatMultHermitianTranspose()`, `MatMultTranspose()`, `MatMultAdd()`, `MatMult()`
2910: @*/
2911: PetscErrorCode MatMultHermitianTransposeAdd(Mat mat, Vec v1, Vec v2, Vec v3)
2912: {
2913:   PetscFunctionBegin;

2920:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
2921:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
2922:   PetscCheck(v1 != v3, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "v1 and v3 must be different vectors");
2923:   PetscCheck(mat->rmap->N == v1->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v1: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, v1->map->N);
2924:   PetscCheck(mat->cmap->N == v2->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v2: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v2->map->N);
2925:   PetscCheck(mat->cmap->N == v3->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v3: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v3->map->N);
2926:   MatCheckPreallocated(mat, 1);

2928:   PetscCall(PetscLogEventBegin(MAT_MultHermitianTransposeAdd, mat, v1, v2, v3));
2929:   PetscCall(VecLockReadPush(v1));
2930:   if (mat->ops->multhermitiantransposeadd) PetscUseTypeMethod(mat, multhermitiantransposeadd, v1, v2, v3);
2931:   else {
2932:     Vec w, z;
2933:     PetscCall(VecDuplicate(v1, &w));
2934:     PetscCall(VecCopy(v1, w));
2935:     PetscCall(VecConjugate(w));
2936:     PetscCall(VecDuplicate(v3, &z));
2937:     PetscCall(MatMultTranspose(mat, w, z));
2938:     PetscCall(VecDestroy(&w));
2939:     PetscCall(VecConjugate(z));
2940:     if (v2 != v3) {
2941:       PetscCall(VecWAXPY(v3, 1.0, v2, z));
2942:     } else {
2943:       PetscCall(VecAXPY(v3, 1.0, z));
2944:     }
2945:     PetscCall(VecDestroy(&z));
2946:   }
2947:   PetscCall(VecLockReadPop(v1));
2948:   PetscCall(PetscLogEventEnd(MAT_MultHermitianTransposeAdd, mat, v1, v2, v3));
2949:   PetscCall(PetscObjectStateIncrease((PetscObject)v3));
2950:   PetscFunctionReturn(PETSC_SUCCESS);
2951: }

2953: /*@
2954:   MatGetFactorType - gets the type of factorization a matrix is

2956:   Not Collective

2958:   Input Parameter:
2959: . mat - the matrix

2961:   Output Parameter:
2962: . t - the type, one of `MAT_FACTOR_NONE`, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ILU`, `MAT_FACTOR_ICC,MAT_FACTOR_ILUDT`, `MAT_FACTOR_QR`

2964:   Level: intermediate

2966: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorType`, `MatGetFactor()`, `MatSetFactorType()`, `MAT_FACTOR_NONE`, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ILU`,
2967:           `MAT_FACTOR_ICC`,`MAT_FACTOR_ILUDT`, `MAT_FACTOR_QR`
2968: @*/
2969: PetscErrorCode MatGetFactorType(Mat mat, MatFactorType *t)
2970: {
2971:   PetscFunctionBegin;
2974:   PetscAssertPointer(t, 2);
2975:   *t = mat->factortype;
2976:   PetscFunctionReturn(PETSC_SUCCESS);
2977: }

2979: /*@
2980:   MatSetFactorType - sets the type of factorization a matrix is

2982:   Logically Collective

2984:   Input Parameters:
2985: + mat - the matrix
2986: - t   - the type, one of `MAT_FACTOR_NONE`, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ILU`, `MAT_FACTOR_ICC,MAT_FACTOR_ILUDT`, `MAT_FACTOR_QR`

2988:   Level: intermediate

2990: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorType`, `MatGetFactor()`, `MatGetFactorType()`, `MAT_FACTOR_NONE`, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ILU`,
2991:           `MAT_FACTOR_ICC`,`MAT_FACTOR_ILUDT`, `MAT_FACTOR_QR`
2992: @*/
2993: PetscErrorCode MatSetFactorType(Mat mat, MatFactorType t)
2994: {
2995:   PetscFunctionBegin;
2998:   mat->factortype = t;
2999:   PetscFunctionReturn(PETSC_SUCCESS);
3000: }

3002: /*@
3003:   MatGetInfo - Returns information about matrix storage (number of
3004:   nonzeros, memory, etc.).

3006:   Collective if `MAT_GLOBAL_MAX` or `MAT_GLOBAL_SUM` is used as the flag

3008:   Input Parameters:
3009: + mat  - the matrix
3010: - flag - flag indicating the type of parameters to be returned (`MAT_LOCAL` - local matrix, `MAT_GLOBAL_MAX` - maximum over all processors, `MAT_GLOBAL_SUM` - sum over all processors)

3012:   Output Parameter:
3013: . info - matrix information context

3015:   Options Database Key:
3016: . -mat_view ::ascii_info - print matrix info to `PETSC_STDOUT`

3018:   Level: intermediate

3020:   Notes:
3021:   The `MatInfo` context contains a variety of matrix data, including
3022:   number of nonzeros allocated and used, number of mallocs during
3023:   matrix assembly, etc.  Additional information for factored matrices
3024:   is provided (such as the fill ratio, number of mallocs during
3025:   factorization, etc.).

3027:   Example:
3028:   See the file ${PETSC_DIR}/include/petscmat.h for a complete list of
3029:   data within the `MatInfo` context.  For example,
3030: .vb
3031:       MatInfo info;
3032:       Mat     A;
3033:       double  mal, nz_a, nz_u;

3035:       MatGetInfo(A, MAT_LOCAL, &info);
3036:       mal  = info.mallocs;
3037:       nz_a = info.nz_allocated;
3038: .ve

3040:   Fortran Note:
3041:   Declare info as a `MatInfo` array of dimension `MAT_INFO_SIZE`, and then extract the parameters
3042:   of interest.  See the file ${PETSC_DIR}/include/petsc/finclude/petscmat.h
3043:   a complete list of parameter names.
3044: .vb
3045:       MatInfo info(MAT_INFO_SIZE)
3046:       double  precision mal, nz_a
3047:       Mat     A
3048:       integer ierr

3050:       call MatGetInfo(A, MAT_LOCAL, info, ierr)
3051:       mal = info(MAT_INFO_MALLOCS)
3052:       nz_a = info(MAT_INFO_NZ_ALLOCATED)
3053: .ve

3055: .seealso: [](ch_matrices), `Mat`, `MatInfo`, `MatStashGetInfo()`
3056: @*/
3057: PetscErrorCode MatGetInfo(Mat mat, MatInfoType flag, MatInfo *info)
3058: {
3059:   PetscFunctionBegin;
3062:   PetscAssertPointer(info, 3);
3063:   MatCheckPreallocated(mat, 1);
3064:   PetscUseTypeMethod(mat, getinfo, flag, info);
3065:   PetscFunctionReturn(PETSC_SUCCESS);
3066: }

3068: /*
3069:    This is used by external packages where it is not easy to get the info from the actual
3070:    matrix factorization.
3071: */
3072: PetscErrorCode MatGetInfo_External(Mat A, MatInfoType flag, MatInfo *info)
3073: {
3074:   PetscFunctionBegin;
3075:   PetscCall(PetscMemzero(info, sizeof(MatInfo)));
3076:   PetscFunctionReturn(PETSC_SUCCESS);
3077: }

3079: /*@
3080:   MatLUFactor - Performs in-place LU factorization of matrix.

3082:   Collective

3084:   Input Parameters:
3085: + mat  - the matrix
3086: . row  - row permutation
3087: . col  - column permutation
3088: - info - options for factorization, includes
3089: .vb
3090:           fill - expected fill as ratio of original fill.
3091:           dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3092:                    Run with the option -info to determine an optimal value to use
3093: .ve

3095:   Level: developer

3097:   Notes:
3098:   Most users should employ the `KSP` interface for linear solvers
3099:   instead of working directly with matrix algebra routines such as this.
3100:   See, e.g., `KSPCreate()`.

3102:   This changes the state of the matrix to a factored matrix; it cannot be used
3103:   for example with `MatSetValues()` unless one first calls `MatSetUnfactored()`.

3105:   This is really in-place only for dense matrices, the preferred approach is to use `MatGetFactor()`, `MatLUFactorSymbolic()`, and `MatLUFactorNumeric()`
3106:   when not using `KSP`.

3108:   Developer Note:
3109:   The Fortran interface is not autogenerated as the
3110:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

3112: .seealso: [](ch_matrices), [Matrix Factorization](sec_matfactor), `Mat`, `MatFactorType`, `MatLUFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`,
3113:           `MatGetOrdering()`, `MatSetUnfactored()`, `MatFactorInfo`, `MatGetFactor()`
3114: @*/
3115: PetscErrorCode MatLUFactor(Mat mat, IS row, IS col, const MatFactorInfo *info)
3116: {
3117:   MatFactorInfo tinfo;

3119:   PetscFunctionBegin;
3123:   if (info) PetscAssertPointer(info, 4);
3125:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3126:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3127:   MatCheckPreallocated(mat, 1);
3128:   if (!info) {
3129:     PetscCall(MatFactorInfoInitialize(&tinfo));
3130:     info = &tinfo;
3131:   }

3133:   PetscCall(PetscLogEventBegin(MAT_LUFactor, mat, row, col, 0));
3134:   PetscUseTypeMethod(mat, lufactor, row, col, info);
3135:   PetscCall(PetscLogEventEnd(MAT_LUFactor, mat, row, col, 0));
3136:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
3137:   PetscFunctionReturn(PETSC_SUCCESS);
3138: }

3140: /*@
3141:   MatILUFactor - Performs in-place ILU factorization of matrix.

3143:   Collective

3145:   Input Parameters:
3146: + mat  - the matrix
3147: . row  - row permutation
3148: . col  - column permutation
3149: - info - structure containing
3150: .vb
3151:       levels - number of levels of fill.
3152:       expected fill - as ratio of original fill.
3153:       1 or 0 - indicating force fill on diagonal (improves robustness for matrices
3154:                 missing diagonal entries)
3155: .ve

3157:   Level: developer

3159:   Notes:
3160:   Most users should employ the `KSP` interface for linear solvers
3161:   instead of working directly with matrix algebra routines such as this.
3162:   See, e.g., `KSPCreate()`.

3164:   Probably really in-place only when level of fill is zero, otherwise allocates
3165:   new space to store factored matrix and deletes previous memory. The preferred approach is to use `MatGetFactor()`, `MatILUFactorSymbolic()`, and `MatILUFactorNumeric()`
3166:   when not using `KSP`.

3168:   Developer Note:
3169:   The Fortran interface is not autogenerated as the
3170:   interface definition cannot be generated correctly [due to MatFactorInfo]

3172: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatILUFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`, `MatFactorInfo`
3173: @*/
3174: PetscErrorCode MatILUFactor(Mat mat, IS row, IS col, const MatFactorInfo *info)
3175: {
3176:   PetscFunctionBegin;
3180:   PetscAssertPointer(info, 4);
3182:   PetscCheck(mat->rmap->N == mat->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "matrix must be square");
3183:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3184:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3185:   MatCheckPreallocated(mat, 1);

3187:   PetscCall(PetscLogEventBegin(MAT_ILUFactor, mat, row, col, 0));
3188:   PetscUseTypeMethod(mat, ilufactor, row, col, info);
3189:   PetscCall(PetscLogEventEnd(MAT_ILUFactor, mat, row, col, 0));
3190:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
3191:   PetscFunctionReturn(PETSC_SUCCESS);
3192: }

3194: /*@
3195:   MatLUFactorSymbolic - Performs symbolic LU factorization of matrix.
3196:   Call this routine before calling `MatLUFactorNumeric()` and after `MatGetFactor()`.

3198:   Collective

3200:   Input Parameters:
3201: + fact - the factor matrix obtained with `MatGetFactor()`
3202: . mat  - the matrix
3203: . row  - the row permutation
3204: . col  - the column permutation
3205: - info - options for factorization, includes
3206: .vb
3207:           fill - expected fill as ratio of original fill. Run with the option -info to determine an optimal value to use
3208:           dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3209: .ve

3211:   Level: developer

3213:   Notes:
3214:   See [Matrix Factorization](sec_matfactor) for additional information about factorizations

3216:   Most users should employ the simplified `KSP` interface for linear solvers
3217:   instead of working directly with matrix algebra routines such as this.
3218:   See, e.g., `KSPCreate()`.

3220:   Developer Note:
3221:   The Fortran interface is not autogenerated as the
3222:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

3224: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatLUFactor()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`, `MatFactorInfo`, `MatFactorInfoInitialize()`
3225: @*/
3226: PetscErrorCode MatLUFactorSymbolic(Mat fact, Mat mat, IS row, IS col, const MatFactorInfo *info)
3227: {
3228:   MatFactorInfo tinfo;

3230:   PetscFunctionBegin;
3235:   if (info) PetscAssertPointer(info, 5);
3238:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3239:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3240:   MatCheckPreallocated(mat, 2);
3241:   if (!info) {
3242:     PetscCall(MatFactorInfoInitialize(&tinfo));
3243:     info = &tinfo;
3244:   }

3246:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_LUFactorSymbolic, mat, row, col, 0));
3247:   PetscUseTypeMethod(fact, lufactorsymbolic, mat, row, col, info);
3248:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_LUFactorSymbolic, mat, row, col, 0));
3249:   PetscCall(PetscObjectStateIncrease((PetscObject)fact));
3250:   PetscFunctionReturn(PETSC_SUCCESS);
3251: }

3253: /*@
3254:   MatLUFactorNumeric - Performs numeric LU factorization of a matrix.
3255:   Call this routine after first calling `MatLUFactorSymbolic()` and `MatGetFactor()`.

3257:   Collective

3259:   Input Parameters:
3260: + fact - the factor matrix obtained with `MatGetFactor()`
3261: . mat  - the matrix
3262: - info - options for factorization

3264:   Level: developer

3266:   Notes:
3267:   See `MatLUFactor()` for in-place factorization.  See
3268:   `MatCholeskyFactorNumeric()` for the symmetric, positive definite case.

3270:   Most users should employ the `KSP` interface for linear solvers
3271:   instead of working directly with matrix algebra routines such as this.
3272:   See, e.g., `KSPCreate()`.

3274:   Developer Note:
3275:   The Fortran interface is not autogenerated as the
3276:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

3278: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatFactorInfo`, `MatLUFactorSymbolic()`, `MatLUFactor()`, `MatCholeskyFactor()`
3279: @*/
3280: PetscErrorCode MatLUFactorNumeric(Mat fact, Mat mat, const MatFactorInfo *info)
3281: {
3282:   MatFactorInfo tinfo;

3284:   PetscFunctionBegin;
3289:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3290:   PetscCheck(mat->rmap->N == (fact)->rmap->N && mat->cmap->N == (fact)->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Mat fact: global dimensions are different %" PetscInt_FMT " should = %" PetscInt_FMT " %" PetscInt_FMT " should = %" PetscInt_FMT,
3291:              mat->rmap->N, (fact)->rmap->N, mat->cmap->N, (fact)->cmap->N);

3293:   MatCheckPreallocated(mat, 2);
3294:   if (!info) {
3295:     PetscCall(MatFactorInfoInitialize(&tinfo));
3296:     info = &tinfo;
3297:   }

3299:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_LUFactorNumeric, mat, fact, 0, 0));
3300:   else PetscCall(PetscLogEventBegin(MAT_LUFactor, mat, fact, 0, 0));
3301:   PetscUseTypeMethod(fact, lufactornumeric, mat, info);
3302:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_LUFactorNumeric, mat, fact, 0, 0));
3303:   else PetscCall(PetscLogEventEnd(MAT_LUFactor, mat, fact, 0, 0));
3304:   PetscCall(MatViewFromOptions(fact, NULL, "-mat_factor_view"));
3305:   PetscCall(PetscObjectStateIncrease((PetscObject)fact));
3306:   PetscFunctionReturn(PETSC_SUCCESS);
3307: }

3309: /*@
3310:   MatCholeskyFactor - Performs in-place Cholesky factorization of a
3311:   symmetric matrix.

3313:   Collective

3315:   Input Parameters:
3316: + mat  - the matrix
3317: . perm - row and column permutations
3318: - info - expected fill as ratio of original fill

3320:   Level: developer

3322:   Notes:
3323:   See `MatLUFactor()` for the nonsymmetric case.  See also `MatGetFactor()`,
3324:   `MatCholeskyFactorSymbolic()`, and `MatCholeskyFactorNumeric()`.

3326:   Most users should employ the `KSP` interface for linear solvers
3327:   instead of working directly with matrix algebra routines such as this.
3328:   See, e.g., `KSPCreate()`.

3330:   Developer Note:
3331:   The Fortran interface is not autogenerated as the
3332:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

3334: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatFactorInfo`, `MatLUFactor()`, `MatCholeskyFactorSymbolic()`, `MatCholeskyFactorNumeric()`
3335:           `MatGetOrdering()`
3336: @*/
3337: PetscErrorCode MatCholeskyFactor(Mat mat, IS perm, const MatFactorInfo *info)
3338: {
3339:   MatFactorInfo tinfo;

3341:   PetscFunctionBegin;
3344:   if (info) PetscAssertPointer(info, 3);
3346:   PetscCheck(mat->rmap->N == mat->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix must be square");
3347:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3348:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3349:   MatCheckPreallocated(mat, 1);
3350:   if (!info) {
3351:     PetscCall(MatFactorInfoInitialize(&tinfo));
3352:     info = &tinfo;
3353:   }

3355:   PetscCall(PetscLogEventBegin(MAT_CholeskyFactor, mat, perm, 0, 0));
3356:   PetscUseTypeMethod(mat, choleskyfactor, perm, info);
3357:   PetscCall(PetscLogEventEnd(MAT_CholeskyFactor, mat, perm, 0, 0));
3358:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
3359:   PetscFunctionReturn(PETSC_SUCCESS);
3360: }

3362: /*@
3363:   MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization
3364:   of a symmetric matrix.

3366:   Collective

3368:   Input Parameters:
3369: + fact - the factor matrix obtained with `MatGetFactor()`
3370: . mat  - the matrix
3371: . perm - row and column permutations
3372: - info - options for factorization, includes
3373: .vb
3374:           fill - expected fill as ratio of original fill.
3375:           dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3376:                    Run with the option -info to determine an optimal value to use
3377: .ve

3379:   Level: developer

3381:   Notes:
3382:   See `MatLUFactorSymbolic()` for the nonsymmetric case.  See also
3383:   `MatCholeskyFactor()` and `MatCholeskyFactorNumeric()`.

3385:   Most users should employ the `KSP` interface for linear solvers
3386:   instead of working directly with matrix algebra routines such as this.
3387:   See, e.g., `KSPCreate()`.

3389:   Developer Note:
3390:   The Fortran interface is not autogenerated as the
3391:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

3393: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorInfo`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactor()`, `MatCholeskyFactorNumeric()`
3394:           `MatGetOrdering()`
3395: @*/
3396: PetscErrorCode MatCholeskyFactorSymbolic(Mat fact, Mat mat, IS perm, const MatFactorInfo *info)
3397: {
3398:   MatFactorInfo tinfo;

3400:   PetscFunctionBegin;
3404:   if (info) PetscAssertPointer(info, 4);
3407:   PetscCheck(mat->rmap->N == mat->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix must be square");
3408:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3409:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3410:   MatCheckPreallocated(mat, 2);
3411:   if (!info) {
3412:     PetscCall(MatFactorInfoInitialize(&tinfo));
3413:     info = &tinfo;
3414:   }

3416:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_CholeskyFactorSymbolic, mat, perm, 0, 0));
3417:   PetscUseTypeMethod(fact, choleskyfactorsymbolic, mat, perm, info);
3418:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_CholeskyFactorSymbolic, mat, perm, 0, 0));
3419:   PetscCall(PetscObjectStateIncrease((PetscObject)fact));
3420:   PetscFunctionReturn(PETSC_SUCCESS);
3421: }

3423: /*@
3424:   MatCholeskyFactorNumeric - Performs numeric Cholesky factorization
3425:   of a symmetric matrix. Call this routine after first calling `MatGetFactor()` and
3426:   `MatCholeskyFactorSymbolic()`.

3428:   Collective

3430:   Input Parameters:
3431: + fact - the factor matrix obtained with `MatGetFactor()`, where the factored values are stored
3432: . mat  - the initial matrix that is to be factored
3433: - info - options for factorization

3435:   Level: developer

3437:   Note:
3438:   Most users should employ the `KSP` interface for linear solvers
3439:   instead of working directly with matrix algebra routines such as this.
3440:   See, e.g., `KSPCreate()`.

3442:   Developer Note:
3443:   The Fortran interface is not autogenerated as the
3444:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

3446: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorInfo`, `MatGetFactor()`, `MatCholeskyFactorSymbolic()`, `MatCholeskyFactor()`, `MatLUFactorNumeric()`
3447: @*/
3448: PetscErrorCode MatCholeskyFactorNumeric(Mat fact, Mat mat, const MatFactorInfo *info)
3449: {
3450:   MatFactorInfo tinfo;

3452:   PetscFunctionBegin;
3457:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3458:   PetscCheck(mat->rmap->N == (fact)->rmap->N && mat->cmap->N == (fact)->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Mat fact: global dim %" PetscInt_FMT " should = %" PetscInt_FMT " %" PetscInt_FMT " should = %" PetscInt_FMT,
3459:              mat->rmap->N, (fact)->rmap->N, mat->cmap->N, (fact)->cmap->N);
3460:   MatCheckPreallocated(mat, 2);
3461:   if (!info) {
3462:     PetscCall(MatFactorInfoInitialize(&tinfo));
3463:     info = &tinfo;
3464:   }

3466:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_CholeskyFactorNumeric, mat, fact, 0, 0));
3467:   else PetscCall(PetscLogEventBegin(MAT_CholeskyFactor, mat, fact, 0, 0));
3468:   PetscUseTypeMethod(fact, choleskyfactornumeric, mat, info);
3469:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_CholeskyFactorNumeric, mat, fact, 0, 0));
3470:   else PetscCall(PetscLogEventEnd(MAT_CholeskyFactor, mat, fact, 0, 0));
3471:   PetscCall(MatViewFromOptions(fact, NULL, "-mat_factor_view"));
3472:   PetscCall(PetscObjectStateIncrease((PetscObject)fact));
3473:   PetscFunctionReturn(PETSC_SUCCESS);
3474: }

3476: /*@
3477:   MatQRFactor - Performs in-place QR factorization of matrix.

3479:   Collective

3481:   Input Parameters:
3482: + mat  - the matrix
3483: . col  - column permutation
3484: - info - options for factorization, includes
3485: .vb
3486:           fill - expected fill as ratio of original fill.
3487:           dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3488:                    Run with the option -info to determine an optimal value to use
3489: .ve

3491:   Level: developer

3493:   Notes:
3494:   Most users should employ the `KSP` interface for linear solvers
3495:   instead of working directly with matrix algebra routines such as this.
3496:   See, e.g., `KSPCreate()`.

3498:   This changes the state of the matrix to a factored matrix; it cannot be used
3499:   for example with `MatSetValues()` unless one first calls `MatSetUnfactored()`.

3501:   Developer Note:
3502:   The Fortran interface is not autogenerated as the
3503:   interface definition cannot be generated correctly [due to MatFactorInfo]

3505: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorInfo`, `MatGetFactor()`, `MatQRFactorSymbolic()`, `MatQRFactorNumeric()`, `MatLUFactor()`,
3506:           `MatSetUnfactored()`
3507: @*/
3508: PetscErrorCode MatQRFactor(Mat mat, IS col, const MatFactorInfo *info)
3509: {
3510:   PetscFunctionBegin;
3513:   if (info) PetscAssertPointer(info, 3);
3515:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3516:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3517:   MatCheckPreallocated(mat, 1);
3518:   PetscCall(PetscLogEventBegin(MAT_QRFactor, mat, col, 0, 0));
3519:   PetscUseMethod(mat, "MatQRFactor_C", (Mat, IS, const MatFactorInfo *), (mat, col, info));
3520:   PetscCall(PetscLogEventEnd(MAT_QRFactor, mat, col, 0, 0));
3521:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
3522:   PetscFunctionReturn(PETSC_SUCCESS);
3523: }

3525: /*@
3526:   MatQRFactorSymbolic - Performs symbolic QR factorization of matrix.
3527:   Call this routine after `MatGetFactor()` but before calling `MatQRFactorNumeric()`.

3529:   Collective

3531:   Input Parameters:
3532: + fact - the factor matrix obtained with `MatGetFactor()`
3533: . mat  - the matrix
3534: . col  - column permutation
3535: - info - options for factorization, includes
3536: .vb
3537:           fill - expected fill as ratio of original fill.
3538:           dtcol - pivot tolerance (0 no pivot, 1 full column pivoting)
3539:                    Run with the option -info to determine an optimal value to use
3540: .ve

3542:   Level: developer

3544:   Note:
3545:   Most users should employ the `KSP` interface for linear solvers
3546:   instead of working directly with matrix algebra routines such as this.
3547:   See, e.g., `KSPCreate()`.

3549:   Developer Note:
3550:   The Fortran interface is not autogenerated as the
3551:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

3553: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatFactorInfo`, `MatQRFactor()`, `MatQRFactorNumeric()`, `MatLUFactor()`, `MatFactorInfoInitialize()`
3554: @*/
3555: PetscErrorCode MatQRFactorSymbolic(Mat fact, Mat mat, IS col, const MatFactorInfo *info)
3556: {
3557:   MatFactorInfo tinfo;

3559:   PetscFunctionBegin;
3563:   if (info) PetscAssertPointer(info, 4);
3566:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3567:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
3568:   MatCheckPreallocated(mat, 2);
3569:   if (!info) {
3570:     PetscCall(MatFactorInfoInitialize(&tinfo));
3571:     info = &tinfo;
3572:   }

3574:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_QRFactorSymbolic, fact, mat, col, 0));
3575:   PetscUseMethod(fact, "MatQRFactorSymbolic_C", (Mat, Mat, IS, const MatFactorInfo *), (fact, mat, col, info));
3576:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_QRFactorSymbolic, fact, mat, col, 0));
3577:   PetscCall(PetscObjectStateIncrease((PetscObject)fact));
3578:   PetscFunctionReturn(PETSC_SUCCESS);
3579: }

3581: /*@
3582:   MatQRFactorNumeric - Performs numeric QR factorization of a matrix.
3583:   Call this routine after first calling `MatGetFactor()`, and `MatQRFactorSymbolic()`.

3585:   Collective

3587:   Input Parameters:
3588: + fact - the factor matrix obtained with `MatGetFactor()`
3589: . mat  - the matrix
3590: - info - options for factorization

3592:   Level: developer

3594:   Notes:
3595:   See `MatQRFactor()` for in-place factorization.

3597:   Most users should employ the `KSP` interface for linear solvers
3598:   instead of working directly with matrix algebra routines such as this.
3599:   See, e.g., `KSPCreate()`.

3601:   Developer Note:
3602:   The Fortran interface is not autogenerated as the
3603:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

3605: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorInfo`, `MatGetFactor()`, `MatQRFactor()`, `MatQRFactorSymbolic()`, `MatLUFactor()`
3606: @*/
3607: PetscErrorCode MatQRFactorNumeric(Mat fact, Mat mat, const MatFactorInfo *info)
3608: {
3609:   MatFactorInfo tinfo;

3611:   PetscFunctionBegin;
3616:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
3617:   PetscCheck(mat->rmap->N == fact->rmap->N && mat->cmap->N == fact->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Mat fact: global dimensions are different %" PetscInt_FMT " should = %" PetscInt_FMT " %" PetscInt_FMT " should = %" PetscInt_FMT,
3618:              mat->rmap->N, (fact)->rmap->N, mat->cmap->N, (fact)->cmap->N);

3620:   MatCheckPreallocated(mat, 2);
3621:   if (!info) {
3622:     PetscCall(MatFactorInfoInitialize(&tinfo));
3623:     info = &tinfo;
3624:   }

3626:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_QRFactorNumeric, mat, fact, 0, 0));
3627:   else PetscCall(PetscLogEventBegin(MAT_QRFactor, mat, fact, 0, 0));
3628:   PetscUseMethod(fact, "MatQRFactorNumeric_C", (Mat, Mat, const MatFactorInfo *), (fact, mat, info));
3629:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_QRFactorNumeric, mat, fact, 0, 0));
3630:   else PetscCall(PetscLogEventEnd(MAT_QRFactor, mat, fact, 0, 0));
3631:   PetscCall(MatViewFromOptions(fact, NULL, "-mat_factor_view"));
3632:   PetscCall(PetscObjectStateIncrease((PetscObject)fact));
3633:   PetscFunctionReturn(PETSC_SUCCESS);
3634: }

3636: /*@
3637:   MatSolve - Solves $A x = b$, given a factored matrix.

3639:   Neighbor-wise Collective

3641:   Input Parameters:
3642: + mat - the factored matrix
3643: - b   - the right-hand-side vector

3645:   Output Parameter:
3646: . x - the result vector

3648:   Level: developer

3650:   Notes:
3651:   The vectors `b` and `x` cannot be the same.  I.e., one cannot
3652:   call `MatSolve`(A,x,x).

3654:   Most users should employ the `KSP` interface for linear solvers
3655:   instead of working directly with matrix algebra routines such as this.
3656:   See, e.g., `KSPCreate()`.

3658: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatLUFactor()`, `MatSolveAdd()`, `MatSolveTranspose()`, `MatSolveTransposeAdd()`
3659: @*/
3660: PetscErrorCode MatSolve(Mat mat, Vec b, Vec x)
3661: {
3662:   PetscFunctionBegin;
3667:   PetscCheckSameComm(mat, 1, b, 2);
3668:   PetscCheckSameComm(mat, 1, x, 3);
3669:   PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors");
3670:   PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N);
3671:   PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N);
3672:   PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n);
3673:   if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
3674:   MatCheckPreallocated(mat, 1);

3676:   PetscCall(PetscLogEventBegin(MAT_Solve, mat, b, x, 0));
3677:   PetscCall(VecFlag(x, mat->factorerrortype));
3678:   if (mat->factorerrortype) {
3679:     PetscCall(PetscInfo(mat, "MatFactorError %d\n", mat->factorerrortype));
3680:   } else PetscUseTypeMethod(mat, solve, b, x);
3681:   PetscCall(PetscLogEventEnd(MAT_Solve, mat, b, x, 0));
3682:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
3683:   PetscFunctionReturn(PETSC_SUCCESS);
3684: }

3686: static PetscErrorCode MatMatSolve_Basic(Mat A, Mat B, Mat X, PetscBool trans)
3687: {
3688:   Vec      b, x;
3689:   PetscInt N, i;
3690:   PetscErrorCode (*f)(Mat, Vec, Vec);
3691:   PetscBool Abound, Bneedconv = PETSC_FALSE, Xneedconv = PETSC_FALSE;

3693:   PetscFunctionBegin;
3694:   if (A->factorerrortype) {
3695:     PetscCall(PetscInfo(A, "MatFactorError %d\n", A->factorerrortype));
3696:     PetscCall(MatSetInf(X));
3697:     PetscFunctionReturn(PETSC_SUCCESS);
3698:   }
3699:   f = (!trans || (!A->ops->solvetranspose && A->symmetric)) ? A->ops->solve : A->ops->solvetranspose;
3700:   PetscCheck(f, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Mat type %s", ((PetscObject)A)->type_name);
3701:   PetscCall(MatBoundToCPU(A, &Abound));
3702:   if (!Abound) {
3703:     PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &Bneedconv, MATSEQDENSE, MATMPIDENSE, ""));
3704:     PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &Xneedconv, MATSEQDENSE, MATMPIDENSE, ""));
3705:   }
3706: #if PetscDefined(HAVE_CUDA)
3707:   if (Bneedconv) PetscCall(MatConvert(B, MATDENSECUDA, MAT_INPLACE_MATRIX, &B));
3708:   if (Xneedconv) PetscCall(MatConvert(X, MATDENSECUDA, MAT_INPLACE_MATRIX, &X));
3709: #elif PetscDefined(HAVE_HIP)
3710:   if (Bneedconv) PetscCall(MatConvert(B, MATDENSEHIP, MAT_INPLACE_MATRIX, &B));
3711:   if (Xneedconv) PetscCall(MatConvert(X, MATDENSEHIP, MAT_INPLACE_MATRIX, &X));
3712: #endif
3713:   PetscCall(MatGetSize(B, NULL, &N));
3714:   for (i = 0; i < N; i++) {
3715:     PetscCall(MatDenseGetColumnVecRead(B, i, &b));
3716:     PetscCall(MatDenseGetColumnVecWrite(X, i, &x));
3717:     PetscCall((*f)(A, b, x));
3718:     PetscCall(MatDenseRestoreColumnVecWrite(X, i, &x));
3719:     PetscCall(MatDenseRestoreColumnVecRead(B, i, &b));
3720:   }
3721:   if (Bneedconv) PetscCall(MatConvert(B, MATDENSE, MAT_INPLACE_MATRIX, &B));
3722:   if (Xneedconv) PetscCall(MatConvert(X, MATDENSE, MAT_INPLACE_MATRIX, &X));
3723:   PetscFunctionReturn(PETSC_SUCCESS);
3724: }

3726: /*@
3727:   MatMatSolve - Solves $A X = B$, given a factored matrix.

3729:   Neighbor-wise Collective

3731:   Input Parameters:
3732: + A - the factored matrix
3733: - B - the right-hand-side matrix `MATDENSE` (or sparse `MATAIJ`-- when using MUMPS)

3735:   Output Parameter:
3736: . X - the result matrix (dense matrix)

3738:   Level: developer

3740:   Note:
3741:   If `B` is a `MATDENSE` matrix then one can call `MatMatSolve`(A,B,B) except with `MATSOLVERMKL_CPARDISO`;
3742:   otherwise, `B` and `X` cannot be the same.

3744: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatSolve()`, `MatMatSolveTranspose()`, `MatLUFactor()`, `MatCholeskyFactor()`
3745: @*/
3746: PetscErrorCode MatMatSolve(Mat A, Mat B, Mat X)
3747: {
3748:   PetscFunctionBegin;
3753:   PetscCheckSameComm(A, 1, B, 2);
3754:   PetscCheckSameComm(A, 1, X, 3);
3755:   PetscCheck(A->cmap->N == X->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat X: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->cmap->N, X->rmap->N);
3756:   PetscCheck(A->rmap->N == B->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat B: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->N, B->rmap->N);
3757:   PetscCheck(X->cmap->N == B->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Solution matrix must have same number of columns as rhs matrix");
3758:   if (!A->rmap->N && !A->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
3759:   MatCheckPreallocated(A, 1);

3761:   PetscCall(PetscLogEventBegin(MAT_MatSolve, A, B, X, 0));
3762:   if (!A->ops->matsolve) {
3763:     PetscCall(PetscInfo(A, "Mat type %s using basic MatMatSolve\n", ((PetscObject)A)->type_name));
3764:     PetscCall(MatMatSolve_Basic(A, B, X, PETSC_FALSE));
3765:   } else PetscUseTypeMethod(A, matsolve, B, X);
3766:   PetscCall(PetscLogEventEnd(MAT_MatSolve, A, B, X, 0));
3767:   PetscCall(PetscObjectStateIncrease((PetscObject)X));
3768:   PetscFunctionReturn(PETSC_SUCCESS);
3769: }

3771: /*@
3772:   MatMatSolveTranspose - Solves $A^T X = B $, given a factored matrix.

3774:   Neighbor-wise Collective

3776:   Input Parameters:
3777: + A - the factored matrix
3778: - B - the right-hand-side matrix  (`MATDENSE` matrix)

3780:   Output Parameter:
3781: . X - the result matrix (dense matrix)

3783:   Level: developer

3785:   Note:
3786:   The matrices `B` and `X` cannot be the same.  I.e., one cannot
3787:   call `MatMatSolveTranspose`(A,X,X).

3789: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatSolveTranspose()`, `MatMatSolve()`, `MatLUFactor()`, `MatCholeskyFactor()`
3790: @*/
3791: PetscErrorCode MatMatSolveTranspose(Mat A, Mat B, Mat X)
3792: {
3793:   PetscFunctionBegin;
3798:   PetscCheckSameComm(A, 1, B, 2);
3799:   PetscCheckSameComm(A, 1, X, 3);
3800:   PetscCheck(X != B, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_IDN, "X and B must be different matrices");
3801:   PetscCheck(A->cmap->N == X->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat X: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->cmap->N, X->rmap->N);
3802:   PetscCheck(A->rmap->N == B->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat B: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->N, B->rmap->N);
3803:   PetscCheck(A->rmap->n == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat A,Mat B: local dim %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->n, B->rmap->n);
3804:   PetscCheck(X->cmap->N >= B->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Solution matrix must have same number of columns as rhs matrix");
3805:   if (!A->rmap->N && !A->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
3806:   MatCheckPreallocated(A, 1);

3808:   PetscCall(PetscLogEventBegin(MAT_MatSolve, A, B, X, 0));
3809:   if (!A->ops->matsolvetranspose) {
3810:     PetscCall(PetscInfo(A, "Mat type %s using basic MatMatSolveTranspose\n", ((PetscObject)A)->type_name));
3811:     PetscCall(MatMatSolve_Basic(A, B, X, PETSC_TRUE));
3812:   } else PetscUseTypeMethod(A, matsolvetranspose, B, X);
3813:   PetscCall(PetscLogEventEnd(MAT_MatSolve, A, B, X, 0));
3814:   PetscCall(PetscObjectStateIncrease((PetscObject)X));
3815:   PetscFunctionReturn(PETSC_SUCCESS);
3816: }

3818: /*@
3819:   MatMatTransposeSolve - Solves $A X = B^T$, given a factored matrix.

3821:   Neighbor-wise Collective

3823:   Input Parameters:
3824: + A  - the factored matrix
3825: - Bt - the transpose of right-hand-side matrix as a `MATDENSE`

3827:   Output Parameter:
3828: . X - the result matrix (dense matrix)

3830:   Level: developer

3832:   Note:
3833:   For MUMPS, it only supports centralized sparse compressed column format on the host processor for right-hand side matrix. User must create `Bt` in sparse compressed row
3834:   format on the host processor and call `MatMatTransposeSolve()` to implement MUMPS' `MatMatSolve()`.

3836: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatMatSolve()`, `MatMatSolveTranspose()`, `MatLUFactor()`, `MatCholeskyFactor()`
3837: @*/
3838: PetscErrorCode MatMatTransposeSolve(Mat A, Mat Bt, Mat X)
3839: {
3840:   PetscFunctionBegin;
3845:   PetscCheckSameComm(A, 1, Bt, 2);
3846:   PetscCheckSameComm(A, 1, X, 3);

3848:   PetscCheck(X != Bt, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_IDN, "X and B must be different matrices");
3849:   PetscCheck(A->cmap->N == X->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat X: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->cmap->N, X->rmap->N);
3850:   PetscCheck(A->rmap->N == Bt->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat Bt: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->N, Bt->cmap->N);
3851:   PetscCheck(X->cmap->N >= Bt->rmap->N, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_SIZ, "Solution matrix must have same number of columns as row number of the rhs matrix");
3852:   if (!A->rmap->N && !A->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
3853:   PetscCheck(A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Unfactored matrix");
3854:   MatCheckPreallocated(A, 1);

3856:   PetscCall(PetscLogEventBegin(MAT_MatTrSolve, A, Bt, X, 0));
3857:   PetscUseTypeMethod(A, mattransposesolve, Bt, X);
3858:   PetscCall(PetscLogEventEnd(MAT_MatTrSolve, A, Bt, X, 0));
3859:   PetscCall(PetscObjectStateIncrease((PetscObject)X));
3860:   PetscFunctionReturn(PETSC_SUCCESS);
3861: }

3863: /*@
3864:   MatForwardSolve - Solves $ L x = b $, given a factored matrix, $A = LU $, or
3865:   $U^T*D^(1/2) x = b$, given a factored symmetric matrix, $A = U^T*D*U$,

3867:   Neighbor-wise Collective

3869:   Input Parameters:
3870: + mat - the factored matrix
3871: - b   - the right-hand-side vector

3873:   Output Parameter:
3874: . x - the result vector

3876:   Level: developer

3878:   Notes:
3879:   `MatSolve()` should be used for most applications, as it performs
3880:   a forward solve followed by a backward solve.

3882:   The vectors `b` and `x` cannot be the same,  i.e., one cannot
3883:   call `MatForwardSolve`(A,x,x).

3885:   For matrix in `MATSEQBAIJ` format with block size larger than 1,
3886:   the diagonal blocks are not implemented as $D = D^(1/2) * D^(1/2)$ yet.
3887:   `MatForwardSolve()` solves $U^T*D y = b$, and
3888:   `MatBackwardSolve()` solves $U x = y$.
3889:   Thus they do not provide a symmetric preconditioner.

3891: .seealso: [](ch_matrices), `Mat`, `MatBackwardSolve()`, `MatGetFactor()`, `MatSolve()`
3892: @*/
3893: PetscErrorCode MatForwardSolve(Mat mat, Vec b, Vec x)
3894: {
3895:   PetscFunctionBegin;
3900:   PetscCheckSameComm(mat, 1, b, 2);
3901:   PetscCheckSameComm(mat, 1, x, 3);
3902:   PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors");
3903:   PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N);
3904:   PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N);
3905:   PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n);
3906:   if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
3907:   MatCheckPreallocated(mat, 1);

3909:   PetscCall(PetscLogEventBegin(MAT_ForwardSolve, mat, b, x, 0));
3910:   PetscUseTypeMethod(mat, forwardsolve, b, x);
3911:   PetscCall(PetscLogEventEnd(MAT_ForwardSolve, mat, b, x, 0));
3912:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
3913:   PetscFunctionReturn(PETSC_SUCCESS);
3914: }

3916: /*@
3917:   MatBackwardSolve - Solves $U x = b$, given a factored matrix, $A = LU$.
3918:   $D^(1/2) U x = b$, given a factored symmetric matrix, $A = U^T*D*U$,

3920:   Neighbor-wise Collective

3922:   Input Parameters:
3923: + mat - the factored matrix
3924: - b   - the right-hand-side vector

3926:   Output Parameter:
3927: . x - the result vector

3929:   Level: developer

3931:   Notes:
3932:   `MatSolve()` should be used for most applications, as it performs
3933:   a forward solve followed by a backward solve.

3935:   The vectors `b` and `x` cannot be the same.  I.e., one cannot
3936:   call `MatBackwardSolve`(A,x,x).

3938:   For matrix in `MATSEQBAIJ` format with block size larger than 1,
3939:   the diagonal blocks are not implemented as $D = D^(1/2) * D^(1/2)$ yet.
3940:   `MatForwardSolve()` solves $U^T*D y = b$, and
3941:   `MatBackwardSolve()` solves $U x = y$.
3942:   Thus they do not provide a symmetric preconditioner.

3944: .seealso: [](ch_matrices), `Mat`, `MatForwardSolve()`, `MatGetFactor()`, `MatSolve()`
3945: @*/
3946: PetscErrorCode MatBackwardSolve(Mat mat, Vec b, Vec x)
3947: {
3948:   PetscFunctionBegin;
3953:   PetscCheckSameComm(mat, 1, b, 2);
3954:   PetscCheckSameComm(mat, 1, x, 3);
3955:   PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors");
3956:   PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N);
3957:   PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N);
3958:   PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n);
3959:   if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
3960:   MatCheckPreallocated(mat, 1);

3962:   PetscCall(PetscLogEventBegin(MAT_BackwardSolve, mat, b, x, 0));
3963:   PetscUseTypeMethod(mat, backwardsolve, b, x);
3964:   PetscCall(PetscLogEventEnd(MAT_BackwardSolve, mat, b, x, 0));
3965:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
3966:   PetscFunctionReturn(PETSC_SUCCESS);
3967: }

3969: /*@
3970:   MatSolveAdd - Computes $x = y + A^{-1}*b$, given a factored matrix.

3972:   Neighbor-wise Collective

3974:   Input Parameters:
3975: + mat - the factored matrix
3976: . b   - the right-hand-side vector
3977: - y   - the vector to be added to

3979:   Output Parameter:
3980: . x - the result vector

3982:   Level: developer

3984:   Note:
3985:   The vectors `b` and `x` cannot be the same.  I.e., one cannot
3986:   call `MatSolveAdd`(A,x,y,x).

3988: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatSolve()`, `MatGetFactor()`, `MatSolveTranspose()`, `MatSolveTransposeAdd()`
3989: @*/
3990: PetscErrorCode MatSolveAdd(Mat mat, Vec b, Vec y, Vec x)
3991: {
3992:   PetscScalar one = 1.0;
3993:   Vec         tmp;

3995:   PetscFunctionBegin;
4001:   PetscCheckSameComm(mat, 1, b, 2);
4002:   PetscCheckSameComm(mat, 1, y, 3);
4003:   PetscCheckSameComm(mat, 1, x, 4);
4004:   PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors");
4005:   PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N);
4006:   PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N);
4007:   PetscCheck(mat->rmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, y->map->N);
4008:   PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n);
4009:   PetscCheck(x->map->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Vec x,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, x->map->n, y->map->n);
4010:   if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
4011:   MatCheckPreallocated(mat, 1);

4013:   PetscCall(PetscLogEventBegin(MAT_SolveAdd, mat, b, x, y));
4014:   PetscCall(VecFlag(x, mat->factorerrortype));
4015:   if (mat->factorerrortype) {
4016:     PetscCall(PetscInfo(mat, "MatFactorError %d\n", mat->factorerrortype));
4017:   } else if (mat->ops->solveadd) {
4018:     PetscUseTypeMethod(mat, solveadd, b, y, x);
4019:   } else {
4020:     /* do the solve then the add manually */
4021:     if (x != y) {
4022:       PetscCall(MatSolve(mat, b, x));
4023:       PetscCall(VecAXPY(x, one, y));
4024:     } else {
4025:       PetscCall(VecDuplicate(x, &tmp));
4026:       PetscCall(VecCopy(x, tmp));
4027:       PetscCall(MatSolve(mat, b, x));
4028:       PetscCall(VecAXPY(x, one, tmp));
4029:       PetscCall(VecDestroy(&tmp));
4030:     }
4031:   }
4032:   PetscCall(PetscLogEventEnd(MAT_SolveAdd, mat, b, x, y));
4033:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
4034:   PetscFunctionReturn(PETSC_SUCCESS);
4035: }

4037: /*@
4038:   MatSolveTranspose - Solves $A^T x = b$, given a factored matrix.

4040:   Neighbor-wise Collective

4042:   Input Parameters:
4043: + mat - the factored matrix
4044: - b   - the right-hand-side vector

4046:   Output Parameter:
4047: . x - the result vector

4049:   Level: developer

4051:   Notes:
4052:   The vectors `b` and `x` cannot be the same.  I.e., one cannot
4053:   call `MatSolveTranspose`(A,x,x).

4055:   Most users should employ the `KSP` interface for linear solvers
4056:   instead of working directly with matrix algebra routines such as this.
4057:   See, e.g., `KSPCreate()`.

4059: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `KSP`, `MatSolve()`, `MatSolveAdd()`, `MatSolveTransposeAdd()`
4060: @*/
4061: PetscErrorCode MatSolveTranspose(Mat mat, Vec b, Vec x)
4062: {
4063:   PetscErrorCode (*f)(Mat, Vec, Vec) = (!mat->ops->solvetranspose && mat->symmetric) ? mat->ops->solve : mat->ops->solvetranspose;

4065:   PetscFunctionBegin;
4070:   PetscCheckSameComm(mat, 1, b, 2);
4071:   PetscCheckSameComm(mat, 1, x, 3);
4072:   PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors");
4073:   PetscCheck(mat->rmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, x->map->N);
4074:   PetscCheck(mat->cmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, b->map->N);
4075:   if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
4076:   MatCheckPreallocated(mat, 1);
4077:   PetscCall(PetscLogEventBegin(MAT_SolveTranspose, mat, b, x, 0));
4078:   PetscCall(VecFlag(x, mat->factorerrortype));
4079:   if (mat->factorerrortype) {
4080:     PetscCall(PetscInfo(mat, "MatFactorError %d\n", mat->factorerrortype));
4081:   } else {
4082:     PetscCheck(f, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Matrix type %s", ((PetscObject)mat)->type_name);
4083:     PetscCall((*f)(mat, b, x));
4084:   }
4085:   PetscCall(PetscLogEventEnd(MAT_SolveTranspose, mat, b, x, 0));
4086:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
4087:   PetscFunctionReturn(PETSC_SUCCESS);
4088: }

4090: /*@
4091:   MatSolveTransposeAdd - Computes $x = y + A^{-T} b$
4092:   factored matrix.

4094:   Neighbor-wise Collective

4096:   Input Parameters:
4097: + mat - the factored matrix
4098: . b   - the right-hand-side vector
4099: - y   - the vector to be added to

4101:   Output Parameter:
4102: . x - the result vector

4104:   Level: developer

4106:   Note:
4107:   The vectors `b` and `x` cannot be the same.  I.e., one cannot
4108:   call `MatSolveTransposeAdd`(A,x,y,x).

4110: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatSolve()`, `MatSolveAdd()`, `MatSolveTranspose()`
4111: @*/
4112: PetscErrorCode MatSolveTransposeAdd(Mat mat, Vec b, Vec y, Vec x)
4113: {
4114:   PetscScalar one = 1.0;
4115:   Vec         tmp;
4116:   PetscErrorCode (*f)(Mat, Vec, Vec, Vec) = (!mat->ops->solvetransposeadd && mat->symmetric) ? mat->ops->solveadd : mat->ops->solvetransposeadd;

4118:   PetscFunctionBegin;
4124:   PetscCheckSameComm(mat, 1, b, 2);
4125:   PetscCheckSameComm(mat, 1, y, 3);
4126:   PetscCheckSameComm(mat, 1, x, 4);
4127:   PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors");
4128:   PetscCheck(mat->rmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, x->map->N);
4129:   PetscCheck(mat->cmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, b->map->N);
4130:   PetscCheck(mat->cmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, y->map->N);
4131:   PetscCheck(x->map->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Vec x,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, x->map->n, y->map->n);
4132:   if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
4133:   MatCheckPreallocated(mat, 1);

4135:   PetscCall(PetscLogEventBegin(MAT_SolveTransposeAdd, mat, b, x, y));
4136:   PetscCall(VecFlag(x, mat->factorerrortype));
4137:   if (mat->factorerrortype) {
4138:     PetscCall(PetscInfo(mat, "MatFactorError %d\n", mat->factorerrortype));
4139:   } else if (f) {
4140:     PetscCall((*f)(mat, b, y, x));
4141:   } else {
4142:     /* do the solve then the add manually */
4143:     if (x != y) {
4144:       PetscCall(MatSolveTranspose(mat, b, x));
4145:       PetscCall(VecAXPY(x, one, y));
4146:     } else {
4147:       PetscCall(VecDuplicate(x, &tmp));
4148:       PetscCall(VecCopy(x, tmp));
4149:       PetscCall(MatSolveTranspose(mat, b, x));
4150:       PetscCall(VecAXPY(x, one, tmp));
4151:       PetscCall(VecDestroy(&tmp));
4152:     }
4153:   }
4154:   PetscCall(PetscLogEventEnd(MAT_SolveTransposeAdd, mat, b, x, y));
4155:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
4156:   PetscFunctionReturn(PETSC_SUCCESS);
4157: }

4159: // PetscClangLinter pragma disable: -fdoc-section-header-unknown
4160: /*@
4161:   MatSOR - Computes relaxation (SOR, Gauss-Seidel) sweeps.

4163:   Neighbor-wise Collective

4165:   Input Parameters:
4166: + mat   - the matrix
4167: . b     - the right-hand side
4168: . omega - the relaxation factor
4169: . flag  - flag indicating the type of SOR (see below)
4170: . shift - diagonal shift
4171: . its   - the number of iterations
4172: - lits  - the number of local iterations

4174:   Output Parameter:
4175: . x - the solution (can contain an initial guess, use option `SOR_ZERO_INITIAL_GUESS` to indicate no guess)

4177:   SOR Flags:
4178: +     `SOR_FORWARD_SWEEP` - forward SOR
4179: .     `SOR_BACKWARD_SWEEP` - backward SOR
4180: .     `SOR_SYMMETRIC_SWEEP` - SSOR (symmetric SOR)
4181: .     `SOR_LOCAL_FORWARD_SWEEP` - local forward SOR
4182: .     `SOR_LOCAL_BACKWARD_SWEEP` - local forward SOR
4183: .     `SOR_LOCAL_SYMMETRIC_SWEEP` - local SSOR
4184: .     `SOR_EISENSTAT` - SOR with Eisenstat trick
4185: .     `SOR_APPLY_UPPER`, `SOR_APPLY_LOWER` - applies
4186:   upper/lower triangular part of matrix to
4187:   vector (with omega)
4188: -     `SOR_ZERO_INITIAL_GUESS` - zero initial guess

4190:   Level: developer

4192:   Notes:
4193:   `SOR_LOCAL_FORWARD_SWEEP`, `SOR_LOCAL_BACKWARD_SWEEP`, and
4194:   `SOR_LOCAL_SYMMETRIC_SWEEP` perform separate independent smoothings
4195:   on each processor.

4197:   Application programmers will not generally use `MatSOR()` directly,
4198:   but instead will employ the `KSP`/`PC` interface.

4200:   For `MATBAIJ`, `MATSBAIJ`, and `MATAIJ` matrices with Inodes this does a block SOR smoothing, otherwise it does a pointwise smoothing

4202:   Most users should employ the `KSP` interface for linear solvers
4203:   instead of working directly with matrix algebra routines such as this.
4204:   See, e.g., `KSPCreate()`.

4206:   Vectors `x` and `b` CANNOT be the same

4208:   The flags are implemented as bitwise inclusive or operations.
4209:   For example, use (`SOR_ZERO_INITIAL_GUESS` | `SOR_SYMMETRIC_SWEEP`)
4210:   to specify a zero initial guess for SSOR.

4212:   Developer Note:
4213:   We should add block SOR support for `MATAIJ` matrices with block size set to great than one and no inodes

4215: .seealso: [](ch_matrices), `Mat`, `MatMult()`, `KSP`, `PC`, `MatGetFactor()`
4216: @*/
4217: PetscErrorCode MatSOR(Mat mat, Vec b, PetscReal omega, MatSORType flag, PetscReal shift, PetscInt its, PetscInt lits, Vec x)
4218: {
4219:   PetscFunctionBegin;
4224:   PetscCheckSameComm(mat, 1, b, 2);
4225:   PetscCheckSameComm(mat, 1, x, 8);
4226:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
4227:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4228:   PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N);
4229:   PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N);
4230:   PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n);
4231:   PetscCheck(its > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Relaxation requires global its %" PetscInt_FMT " positive", its);
4232:   PetscCheck(lits > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Relaxation requires local its %" PetscInt_FMT " positive", lits);
4233:   PetscCheck(b != x, PETSC_COMM_SELF, PETSC_ERR_ARG_IDN, "b and x vector cannot be the same");

4235:   MatCheckPreallocated(mat, 1);
4236:   PetscCall(PetscLogEventBegin(MAT_SOR, mat, b, x, 0));
4237:   PetscUseTypeMethod(mat, sor, b, omega, flag, shift, its, lits, x);
4238:   PetscCall(PetscLogEventEnd(MAT_SOR, mat, b, x, 0));
4239:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
4240:   PetscFunctionReturn(PETSC_SUCCESS);
4241: }

4243: /*
4244:       Default matrix copy routine.
4245: */
4246: PetscErrorCode MatCopy_Basic(Mat A, Mat B, MatStructure str)
4247: {
4248:   PetscInt           i, rstart = 0, rend = 0, nz;
4249:   const PetscInt    *cwork;
4250:   const PetscScalar *vwork;

4252:   PetscFunctionBegin;
4253:   if (B->assembled) PetscCall(MatZeroEntries(B));
4254:   if (str == SAME_NONZERO_PATTERN) {
4255:     PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
4256:     for (i = rstart; i < rend; i++) {
4257:       PetscCall(MatGetRow(A, i, &nz, &cwork, &vwork));
4258:       PetscCall(MatSetValues(B, 1, &i, nz, cwork, vwork, INSERT_VALUES));
4259:       PetscCall(MatRestoreRow(A, i, &nz, &cwork, &vwork));
4260:     }
4261:   } else {
4262:     PetscCall(MatAYPX(B, 0.0, A, str));
4263:   }
4264:   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
4265:   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
4266:   PetscFunctionReturn(PETSC_SUCCESS);
4267: }

4269: /*@
4270:   MatCopy - Copies a matrix to another matrix.

4272:   Collective

4274:   Input Parameters:
4275: + A   - the matrix
4276: - str - `SAME_NONZERO_PATTERN` or `DIFFERENT_NONZERO_PATTERN`

4278:   Output Parameter:
4279: . B - where the copy is put

4281:   Level: intermediate

4283:   Notes:
4284:   If you use `SAME_NONZERO_PATTERN`, then the two matrices must have the same nonzero pattern or the routine will crash.

4286:   `MatCopy()` copies the matrix entries of a matrix to another existing
4287:   matrix (after first zeroing the second matrix).  A related routine is
4288:   `MatConvert()`, which first creates a new matrix and then copies the data.

4290: .seealso: [](ch_matrices), `Mat`, `MatConvert()`, `MatDuplicate()`
4291: @*/
4292: PetscErrorCode MatCopy(Mat A, Mat B, MatStructure str)
4293: {
4294:   PetscInt i;

4296:   PetscFunctionBegin;
4301:   PetscCheckSameComm(A, 1, B, 2);
4302:   MatCheckPreallocated(B, 2);
4303:   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
4304:   PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4305:   PetscCheck(A->rmap->N == B->rmap->N && A->cmap->N == B->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat B: global dim (%" PetscInt_FMT ",%" PetscInt_FMT ") (%" PetscInt_FMT ",%" PetscInt_FMT ")", A->rmap->N, B->rmap->N,
4306:              A->cmap->N, B->cmap->N);
4307:   MatCheckPreallocated(A, 1);
4308:   if (A == B) PetscFunctionReturn(PETSC_SUCCESS);

4310:   PetscCall(PetscLogEventBegin(MAT_Copy, A, B, 0, 0));
4311:   if (A->ops->copy) PetscUseTypeMethod(A, copy, B, str);
4312:   else PetscCall(MatCopy_Basic(A, B, str));

4314:   B->stencil.dim = A->stencil.dim;
4315:   B->stencil.noc = A->stencil.noc;
4316:   for (i = 0; i <= A->stencil.dim + (A->stencil.noc ? 0 : -1); i++) {
4317:     B->stencil.dims[i]   = A->stencil.dims[i];
4318:     B->stencil.starts[i] = A->stencil.starts[i];
4319:   }

4321:   PetscCall(PetscLogEventEnd(MAT_Copy, A, B, 0, 0));
4322:   PetscCall(PetscObjectStateIncrease((PetscObject)B));
4323:   PetscFunctionReturn(PETSC_SUCCESS);
4324: }

4326: /*@
4327:   MatConvert - Converts a matrix to another matrix, either of the same
4328:   or different type.

4330:   Collective

4332:   Input Parameters:
4333: + mat     - the matrix
4334: . newtype - new matrix type.  Use `MATSAME` to create a new matrix of the
4335:             same type as the original matrix.
4336: - reuse   - denotes if the destination matrix is to be created or reused.
4337:             Use `MAT_INPLACE_MATRIX` for inplace conversion (that is when you want the input `Mat` to be changed to contain the matrix in the new format), otherwise use
4338:             `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` (can only be used after the first call was made with `MAT_INITIAL_MATRIX`, causes the matrix space in M to be reused).

4340:   Output Parameter:
4341: . M - pointer to place new matrix

4343:   Level: intermediate

4345:   Notes:
4346:   `MatConvert()` first creates a new matrix and then copies the data from
4347:   the first matrix.  A related routine is `MatCopy()`, which copies the matrix
4348:   entries of one matrix to another already existing matrix context.

4350:   Cannot be used to convert a sequential matrix to parallel or parallel to sequential,
4351:   the MPI communicator of the generated matrix is always the same as the communicator
4352:   of the input matrix.

4354: .seealso: [](ch_matrices), `Mat`, `MatCopy()`, `MatDuplicate()`, `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, `MAT_INPLACE_MATRIX`
4355: @*/
4356: PetscErrorCode MatConvert(Mat mat, MatType newtype, MatReuse reuse, Mat *M)
4357: {
4358:   PetscBool  sametype, issame, flg;
4359:   PetscBool3 issymmetric, ishermitian;
4360:   char       convname[256], mtype[256];
4361:   Mat        B;

4363:   PetscFunctionBegin;
4366:   PetscAssertPointer(M, 4);
4367:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
4368:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4369:   MatCheckPreallocated(mat, 1);

4371:   PetscCall(PetscOptionsGetString(((PetscObject)mat)->options, ((PetscObject)mat)->prefix, "-matconvert_type", mtype, sizeof(mtype), &flg));
4372:   if (flg) newtype = mtype;

4374:   PetscCall(PetscObjectTypeCompare((PetscObject)mat, newtype, &sametype));
4375:   PetscCall(PetscStrcmp(newtype, "same", &issame));
4376:   PetscCheck(!(reuse == MAT_INPLACE_MATRIX) || !(mat != *M), PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "MAT_INPLACE_MATRIX requires same input and output matrix");
4377:   if (reuse == MAT_REUSE_MATRIX) {
4379:     PetscCheck(mat != *M, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "MAT_REUSE_MATRIX means reuse matrix in final argument, perhaps you mean MAT_INPLACE_MATRIX");
4380:   }

4382:   if ((reuse == MAT_INPLACE_MATRIX) && (issame || sametype)) {
4383:     PetscCall(PetscInfo(mat, "Early return for inplace %s %d %d\n", ((PetscObject)mat)->type_name, sametype, issame));
4384:     PetscFunctionReturn(PETSC_SUCCESS);
4385:   }

4387:   /* Cache Mat options because some converters use MatHeaderReplace  */
4388:   issymmetric = mat->symmetric;
4389:   ishermitian = mat->hermitian;

4391:   if ((sametype || issame) && (reuse == MAT_INITIAL_MATRIX) && mat->ops->duplicate) {
4392:     PetscCall(PetscInfo(mat, "Calling duplicate for initial matrix %s %d %d\n", ((PetscObject)mat)->type_name, sametype, issame));
4393:     PetscUseTypeMethod(mat, duplicate, MAT_COPY_VALUES, M);
4394:   } else {
4395:     PetscErrorCode (*conv)(Mat, MatType, MatReuse, Mat *) = NULL;
4396:     const char *prefix[3]                                 = {"seq", "mpi", ""};
4397:     PetscInt    i;
4398:     /*
4399:        Order of precedence:
4400:        0) See if newtype is a superclass of the current matrix.
4401:        1) See if a specialized converter is known to the current matrix.
4402:        2) See if a specialized converter is known to the desired matrix class.
4403:        3) See if a good general converter is registered for the desired class
4404:           (as of 6/27/03 only MATMPIADJ falls into this category).
4405:        4) See if a good general converter is known for the current matrix.
4406:        5) Use a really basic converter.
4407:     */

4409:     /* 0) See if newtype is a superclass of the current matrix.
4410:           i.e mat is mpiaij and newtype is aij */
4411:     for (i = 0; i < 2; i++) {
4412:       PetscCall(PetscStrncpy(convname, prefix[i], sizeof(convname)));
4413:       PetscCall(PetscStrlcat(convname, newtype, sizeof(convname)));
4414:       PetscCall(PetscStrcmp(convname, ((PetscObject)mat)->type_name, &flg));
4415:       PetscCall(PetscInfo(mat, "Check superclass %s %s -> %d\n", convname, ((PetscObject)mat)->type_name, flg));
4416:       if (flg) {
4417:         if (reuse == MAT_INPLACE_MATRIX) {
4418:           PetscCall(PetscInfo(mat, "Early return\n"));
4419:           PetscFunctionReturn(PETSC_SUCCESS);
4420:         } else if (reuse == MAT_INITIAL_MATRIX && mat->ops->duplicate) {
4421:           PetscCall(PetscInfo(mat, "Calling MatDuplicate\n"));
4422:           PetscUseTypeMethod(mat, duplicate, MAT_COPY_VALUES, M);
4423:           PetscFunctionReturn(PETSC_SUCCESS);
4424:         } else if (reuse == MAT_REUSE_MATRIX && mat->ops->copy) {
4425:           PetscCall(PetscInfo(mat, "Calling MatCopy\n"));
4426:           PetscCall(MatCopy(mat, *M, SAME_NONZERO_PATTERN));
4427:           PetscFunctionReturn(PETSC_SUCCESS);
4428:         }
4429:       }
4430:     }
4431:     /* 1) See if a specialized converter is known to the current matrix and the desired class */
4432:     for (i = 0; i < 3; i++) {
4433:       PetscCall(PetscStrncpy(convname, "MatConvert_", sizeof(convname)));
4434:       PetscCall(PetscStrlcat(convname, ((PetscObject)mat)->type_name, sizeof(convname)));
4435:       PetscCall(PetscStrlcat(convname, "_", sizeof(convname)));
4436:       PetscCall(PetscStrlcat(convname, prefix[i], sizeof(convname)));
4437:       PetscCall(PetscStrlcat(convname, issame ? ((PetscObject)mat)->type_name : newtype, sizeof(convname)));
4438:       PetscCall(PetscStrlcat(convname, "_C", sizeof(convname)));
4439:       PetscCall(PetscObjectQueryFunction((PetscObject)mat, convname, &conv));
4440:       PetscCall(PetscInfo(mat, "Check specialized (1) %s (%s) -> %d\n", convname, ((PetscObject)mat)->type_name, !!conv));
4441:       if (conv) goto foundconv;
4442:     }

4444:     /* 2)  See if a specialized converter is known to the desired matrix class. */
4445:     PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &B));
4446:     PetscCall(MatSetSizes(B, mat->rmap->n, mat->cmap->n, mat->rmap->N, mat->cmap->N));
4447:     PetscCall(MatSetType(B, newtype));
4448:     for (i = 0; i < 3; i++) {
4449:       PetscCall(PetscStrncpy(convname, "MatConvert_", sizeof(convname)));
4450:       PetscCall(PetscStrlcat(convname, ((PetscObject)mat)->type_name, sizeof(convname)));
4451:       PetscCall(PetscStrlcat(convname, "_", sizeof(convname)));
4452:       PetscCall(PetscStrlcat(convname, prefix[i], sizeof(convname)));
4453:       PetscCall(PetscStrlcat(convname, newtype, sizeof(convname)));
4454:       PetscCall(PetscStrlcat(convname, "_C", sizeof(convname)));
4455:       PetscCall(PetscObjectQueryFunction((PetscObject)B, convname, &conv));
4456:       PetscCall(PetscInfo(mat, "Check specialized (2) %s (%s) -> %d\n", convname, ((PetscObject)B)->type_name, !!conv));
4457:       if (conv) {
4458:         PetscCall(MatDestroy(&B));
4459:         goto foundconv;
4460:       }
4461:     }

4463:     /* 3) See if a good general converter is registered for the desired class */
4464:     conv = B->ops->convertfrom;
4465:     PetscCall(PetscInfo(mat, "Check convertfrom (%s) -> %d\n", ((PetscObject)B)->type_name, !!conv));
4466:     PetscCall(MatDestroy(&B));
4467:     if (conv) goto foundconv;

4469:     /* 4) See if a good general converter is known for the current matrix */
4470:     if (mat->ops->convert) conv = mat->ops->convert;
4471:     PetscCall(PetscInfo(mat, "Check general convert (%s) -> %d\n", ((PetscObject)mat)->type_name, !!conv));
4472:     if (conv) goto foundconv;

4474:     /* 5) Use a really basic converter. */
4475:     PetscCall(PetscInfo(mat, "Using MatConvert_Basic\n"));
4476:     conv = MatConvert_Basic;

4478:   foundconv:
4479:     PetscCall(PetscLogEventBegin(MAT_Convert, mat, 0, 0, 0));
4480:     PetscCall((*conv)(mat, newtype, reuse, M));
4481:     if (mat->rmap->mapping && mat->cmap->mapping && !(*M)->rmap->mapping && !(*M)->cmap->mapping) {
4482:       /* the block sizes must be same if the mappings are copied over */
4483:       (*M)->rmap->bs = mat->rmap->bs;
4484:       (*M)->cmap->bs = mat->cmap->bs;
4485:       PetscCall(PetscObjectReference((PetscObject)mat->rmap->mapping));
4486:       PetscCall(PetscObjectReference((PetscObject)mat->cmap->mapping));
4487:       (*M)->rmap->mapping = mat->rmap->mapping;
4488:       (*M)->cmap->mapping = mat->cmap->mapping;
4489:     }
4490:     (*M)->stencil.dim = mat->stencil.dim;
4491:     (*M)->stencil.noc = mat->stencil.noc;
4492:     for (i = 0; i <= mat->stencil.dim + (mat->stencil.noc ? 0 : -1); i++) {
4493:       (*M)->stencil.dims[i]   = mat->stencil.dims[i];
4494:       (*M)->stencil.starts[i] = mat->stencil.starts[i];
4495:     }
4496:     PetscCall(PetscLogEventEnd(MAT_Convert, mat, 0, 0, 0));
4497:   }
4498:   PetscCall(PetscObjectStateIncrease((PetscObject)*M));

4500:   /* Copy Mat options */
4501:   if (issymmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetOption(*M, MAT_SYMMETRIC, PETSC_TRUE));
4502:   else if (issymmetric == PETSC_BOOL3_FALSE) PetscCall(MatSetOption(*M, MAT_SYMMETRIC, PETSC_FALSE));
4503:   if (ishermitian == PETSC_BOOL3_TRUE) PetscCall(MatSetOption(*M, MAT_HERMITIAN, PETSC_TRUE));
4504:   else if (ishermitian == PETSC_BOOL3_FALSE) PetscCall(MatSetOption(*M, MAT_HERMITIAN, PETSC_FALSE));
4505:   PetscFunctionReturn(PETSC_SUCCESS);
4506: }

4508: /*@
4509:   MatFactorGetSolverType - Returns name of the package providing the factorization routines

4511:   Not Collective

4513:   Input Parameter:
4514: . mat - the matrix, must be a factored matrix

4516:   Output Parameter:
4517: . type - the string name of the package (do not free this string)

4519:   Level: intermediate

4521:   Fortran Note:
4522:   Pass in an empty string that is long enough and the package name will be copied into it.

4524: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatSolverType`, `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`
4525: @*/
4526: PetscErrorCode MatFactorGetSolverType(Mat mat, MatSolverType *type)
4527: {
4528:   PetscErrorCode (*conv)(Mat, MatSolverType *);

4530:   PetscFunctionBegin;
4533:   PetscAssertPointer(type, 2);
4534:   PetscCheck(mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
4535:   PetscCall(PetscObjectQueryFunction((PetscObject)mat, "MatFactorGetSolverType_C", &conv));
4536:   if (conv) PetscCall((*conv)(mat, type));
4537:   else *type = MATSOLVERPETSC;
4538:   PetscFunctionReturn(PETSC_SUCCESS);
4539: }

4541: typedef struct _MatSolverTypeForSpecifcType *MatSolverTypeForSpecifcType;
4542: struct _MatSolverTypeForSpecifcType {
4543:   MatType mtype;
4544:   /* no entry for MAT_FACTOR_NONE */
4545:   PetscErrorCode (*createfactor[MAT_FACTOR_NUM_TYPES - 1])(Mat, MatFactorType, Mat *);
4546:   MatSolverTypeForSpecifcType next;
4547: };

4549: typedef struct _MatSolverTypeHolder *MatSolverTypeHolder;
4550: struct _MatSolverTypeHolder {
4551:   char                       *name;
4552:   MatSolverTypeForSpecifcType handlers;
4553:   MatSolverTypeHolder         next;
4554: };

4556: static MatSolverTypeHolder MatSolverTypeHolders = NULL;

4558: /*@C
4559:   MatSolverTypeRegister - Registers a `MatSolverType` that works for a particular matrix type

4561:   Logically Collective, No Fortran Support

4563:   Input Parameters:
4564: + package      - name of the package, for example petsc or superlu
4565: . mtype        - the matrix type that works with this package
4566: . ftype        - the type of factorization supported by the package
4567: - createfactor - routine that will create the factored matrix ready to be used

4569:   Level: developer

4571: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorGetSolverType()`, `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`,
4572:   `MatGetFactor()`
4573: @*/
4574: PetscErrorCode MatSolverTypeRegister(MatSolverType package, MatType mtype, MatFactorType ftype, PetscErrorCode (*createfactor)(Mat, MatFactorType, Mat *))
4575: {
4576:   MatSolverTypeHolder         next = MatSolverTypeHolders, prev = NULL;
4577:   PetscBool                   flg;
4578:   MatSolverTypeForSpecifcType inext, iprev = NULL;

4580:   PetscFunctionBegin;
4581:   PetscCall(MatInitializePackage());
4582:   if (!next) {
4583:     PetscCall(PetscNew(&MatSolverTypeHolders));
4584:     PetscCall(PetscStrallocpy(package, &MatSolverTypeHolders->name));
4585:     PetscCall(PetscNew(&MatSolverTypeHolders->handlers));
4586:     PetscCall(PetscStrallocpy(mtype, (char **)&MatSolverTypeHolders->handlers->mtype));
4587:     MatSolverTypeHolders->handlers->createfactor[(int)ftype - 1] = createfactor;
4588:     PetscFunctionReturn(PETSC_SUCCESS);
4589:   }
4590:   while (next) {
4591:     PetscCall(PetscStrcasecmp(package, next->name, &flg));
4592:     if (flg) {
4593:       PetscCheck(next->handlers, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatSolverTypeHolder is missing handlers");
4594:       inext = next->handlers;
4595:       while (inext) {
4596:         PetscCall(PetscStrcasecmp(mtype, inext->mtype, &flg));
4597:         if (flg) {
4598:           inext->createfactor[(int)ftype - 1] = createfactor;
4599:           PetscFunctionReturn(PETSC_SUCCESS);
4600:         }
4601:         iprev = inext;
4602:         inext = inext->next;
4603:       }
4604:       PetscCall(PetscNew(&iprev->next));
4605:       PetscCall(PetscStrallocpy(mtype, (char **)&iprev->next->mtype));
4606:       iprev->next->createfactor[(int)ftype - 1] = createfactor;
4607:       PetscFunctionReturn(PETSC_SUCCESS);
4608:     }
4609:     prev = next;
4610:     next = next->next;
4611:   }
4612:   PetscCall(PetscNew(&prev->next));
4613:   PetscCall(PetscStrallocpy(package, &prev->next->name));
4614:   PetscCall(PetscNew(&prev->next->handlers));
4615:   PetscCall(PetscStrallocpy(mtype, (char **)&prev->next->handlers->mtype));
4616:   prev->next->handlers->createfactor[(int)ftype - 1] = createfactor;
4617:   PetscFunctionReturn(PETSC_SUCCESS);
4618: }

4620: /*@C
4621:   MatSolverTypeGet - Gets the function that creates the factor matrix if it exist

4623:   Input Parameters:
4624: + type  - name of the package, for example petsc or superlu, if this is 'NULL', then the first result that satisfies the other criteria is returned
4625: . ftype - the type of factorization supported by the type
4626: - mtype - the matrix type that works with this type

4628:   Output Parameters:
4629: + foundtype    - `PETSC_TRUE` if the type was registered
4630: . foundmtype   - `PETSC_TRUE` if the type supports the requested mtype
4631: - createfactor - routine that will create the factored matrix ready to be used or `NULL` if not found

4633:   Calling sequence of `createfactor`:
4634: + A     - the matrix providing the factor matrix
4635: . ftype - the `MatFactorType` of the factor requested
4636: - B     - the new factor matrix that responds to MatXXFactorSymbolic,Numeric() functions, such as `MatLUFactorSymbolic()`

4638:   Level: developer

4640:   Note:
4641:   When `type` is `NULL` the available functions are searched for based on the order of the calls to `MatSolverTypeRegister()` in `MatInitializePackage()`.
4642:   Since different PETSc configurations may have different external solvers, seemingly identical runs with different PETSc configurations may use a different solver.
4643:   For example if one configuration had `--download-mumps` while a different one had `--download-superlu_dist`.

4645: .seealso: [](ch_matrices), `Mat`, `MatFactorType`, `MatType`, `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`, `MatSolverTypeRegister()`, `MatGetFactor()`,
4646:           `MatInitializePackage()`
4647: @*/
4648: PetscErrorCode MatSolverTypeGet(MatSolverType type, MatType mtype, MatFactorType ftype, PetscBool *foundtype, PetscBool *foundmtype, PetscErrorCode (**createfactor)(Mat A, MatFactorType ftype, Mat *B))
4649: {
4650:   MatSolverTypeHolder         next = MatSolverTypeHolders;
4651:   PetscBool                   flg;
4652:   MatSolverTypeForSpecifcType inext;

4654:   PetscFunctionBegin;
4655:   if (foundtype) *foundtype = PETSC_FALSE;
4656:   if (foundmtype) *foundmtype = PETSC_FALSE;
4657:   if (createfactor) *createfactor = NULL;

4659:   if (type) {
4660:     while (next) {
4661:       PetscCall(PetscStrcasecmp(type, next->name, &flg));
4662:       if (flg) {
4663:         if (foundtype) *foundtype = PETSC_TRUE;
4664:         inext = next->handlers;
4665:         while (inext) {
4666:           PetscCall(PetscStrbeginswith(mtype, inext->mtype, &flg));
4667:           if (flg) {
4668:             if (foundmtype) *foundmtype = PETSC_TRUE;
4669:             if (createfactor) *createfactor = inext->createfactor[(int)ftype - 1];
4670:             PetscFunctionReturn(PETSC_SUCCESS);
4671:           }
4672:           inext = inext->next;
4673:         }
4674:       }
4675:       next = next->next;
4676:     }
4677:   } else {
4678:     while (next) {
4679:       inext = next->handlers;
4680:       while (inext) {
4681:         PetscCall(PetscStrcmp(mtype, inext->mtype, &flg));
4682:         if (flg && inext->createfactor[(int)ftype - 1]) {
4683:           if (foundtype) *foundtype = PETSC_TRUE;
4684:           if (foundmtype) *foundmtype = PETSC_TRUE;
4685:           if (createfactor) *createfactor = inext->createfactor[(int)ftype - 1];
4686:           PetscFunctionReturn(PETSC_SUCCESS);
4687:         }
4688:         inext = inext->next;
4689:       }
4690:       next = next->next;
4691:     }
4692:     /* try with base classes inext->mtype */
4693:     next = MatSolverTypeHolders;
4694:     while (next) {
4695:       inext = next->handlers;
4696:       while (inext) {
4697:         PetscCall(PetscStrbeginswith(mtype, inext->mtype, &flg));
4698:         if (flg && inext->createfactor[(int)ftype - 1]) {
4699:           if (foundtype) *foundtype = PETSC_TRUE;
4700:           if (foundmtype) *foundmtype = PETSC_TRUE;
4701:           if (createfactor) *createfactor = inext->createfactor[(int)ftype - 1];
4702:           PetscFunctionReturn(PETSC_SUCCESS);
4703:         }
4704:         inext = inext->next;
4705:       }
4706:       next = next->next;
4707:     }
4708:   }
4709:   PetscFunctionReturn(PETSC_SUCCESS);
4710: }

4712: PetscErrorCode MatSolverTypeDestroy(void)
4713: {
4714:   MatSolverTypeHolder         next = MatSolverTypeHolders, prev;
4715:   MatSolverTypeForSpecifcType inext, iprev;

4717:   PetscFunctionBegin;
4718:   while (next) {
4719:     PetscCall(PetscFree(next->name));
4720:     inext = next->handlers;
4721:     while (inext) {
4722:       PetscCall(PetscFree(inext->mtype));
4723:       iprev = inext;
4724:       inext = inext->next;
4725:       PetscCall(PetscFree(iprev));
4726:     }
4727:     prev = next;
4728:     next = next->next;
4729:     PetscCall(PetscFree(prev));
4730:   }
4731:   MatSolverTypeHolders = NULL;
4732:   PetscFunctionReturn(PETSC_SUCCESS);
4733: }

4735: /*@
4736:   MatFactorGetCanUseOrdering - Indicates if the factorization can use the ordering provided in `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`

4738:   Logically Collective

4740:   Input Parameter:
4741: . mat - the matrix

4743:   Output Parameter:
4744: . flg - `PETSC_TRUE` if uses the ordering

4746:   Level: developer

4748:   Note:
4749:   Most internal PETSc factorizations use the ordering passed to the factorization routine but external
4750:   packages do not, thus we want to skip generating the ordering when it is not needed or used.

4752: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`
4753: @*/
4754: PetscErrorCode MatFactorGetCanUseOrdering(Mat mat, PetscBool *flg)
4755: {
4756:   PetscFunctionBegin;
4757:   *flg = mat->canuseordering;
4758:   PetscFunctionReturn(PETSC_SUCCESS);
4759: }

4761: /*@
4762:   MatFactorGetPreferredOrdering - The preferred ordering for a particular matrix factor object

4764:   Logically Collective

4766:   Input Parameters:
4767: + mat   - the matrix obtained with `MatGetFactor()`
4768: - ftype - the factorization type to be used

4770:   Output Parameter:
4771: . otype - the preferred ordering type

4773:   Level: developer

4775: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorType`, `MatOrderingType`, `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`
4776: @*/
4777: PetscErrorCode MatFactorGetPreferredOrdering(Mat mat, MatFactorType ftype, MatOrderingType *otype)
4778: {
4779:   PetscFunctionBegin;
4780:   *otype = mat->preferredordering[ftype];
4781:   PetscCheck(*otype, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatFactor did not have a preferred ordering");
4782:   PetscFunctionReturn(PETSC_SUCCESS);
4783: }

4785: /*@
4786:   MatGetFactor - Returns a matrix suitable to calls to MatXXFactorSymbolic,Numeric()

4788:   Collective

4790:   Input Parameters:
4791: + mat   - the matrix
4792: . type  - name of solver type, for example, superlu, petsc (to use PETSc's solver if it is available), if this is 'NULL', then the first result that satisfies
4793:           the other criteria is returned
4794: - ftype - factor type, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ICC`, `MAT_FACTOR_ILU`, `MAT_FACTOR_QR`

4796:   Output Parameter:
4797: . f - the factor matrix used with MatXXFactorSymbolic,Numeric() calls. Can be `NULL` in some cases, see notes below.

4799:   Options Database Keys:
4800: + -pc_factor_mat_solver_type <type>             - choose the type at run time. When using `KSP` solvers
4801: - -mat_factor_bind_factorization <host, device> - Where to do matrix factorization? Default is device (might consume more device memory.
4802:                                                   One can choose host to save device memory). Currently only supported with `MATSEQAIJCUSPARSE` matrices.

4804:   Level: intermediate

4806:   Notes:
4807:   The return matrix can be `NULL` if the requested factorization is not available, since some combinations of matrix types and factorization
4808:   types registered with `MatSolverTypeRegister()` cannot be fully tested if not at runtime.

4810:   Users usually access the factorization solvers via `KSP`

4812:   Some PETSc matrix formats have alternative solvers available that are contained in alternative packages
4813:   such as pastix, superlu, mumps etc. PETSc must have been ./configure to use the external solver, using the option --download-package or --with-package-dir

4815:   When `type` is `NULL` the available results are searched for based on the order of the calls to `MatSolverTypeRegister()` in `MatInitializePackage()`.
4816:   Since different PETSc configurations may have different external solvers, seemingly identical runs with different PETSc configurations may use a different solver.
4817:   For example if one configuration had --download-mumps while a different one had --download-superlu_dist.

4819:   Some of the packages have options for controlling the factorization, these are in the form -prefix_mat_packagename_packageoption
4820:   where prefix is normally obtained from the calling `KSP`/`PC`. If `MatGetFactor()` is called directly one can set
4821:   call `MatSetOptionsPrefixFactor()` on the originating matrix or  `MatSetOptionsPrefix()` on the resulting factor matrix.

4823:   Developer Note:
4824:   This should actually be called `MatCreateFactor()` since it creates a new factor object

4826: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `KSP`, `MatSolverType`, `MatFactorType`, `MatCopy()`, `MatDuplicate()`,
4827:           `MatGetFactorAvailable()`, `MatFactorGetCanUseOrdering()`, `MatSolverTypeRegister()`, `MatSolverTypeGet()`
4828:           `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ICC`, `MAT_FACTOR_ILU`, `MAT_FACTOR_QR`, `MatInitializePackage()`
4829: @*/
4830: PetscErrorCode MatGetFactor(Mat mat, MatSolverType type, MatFactorType ftype, Mat *f)
4831: {
4832:   PetscBool foundtype, foundmtype, shell, hasop = PETSC_FALSE;
4833:   PetscErrorCode (*conv)(Mat, MatFactorType, Mat *);

4835:   PetscFunctionBegin;

4839:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4840:   MatCheckPreallocated(mat, 1);

4842:   PetscCall(MatIsShell(mat, &shell));
4843:   if (shell) PetscCall(MatHasOperation(mat, MATOP_GET_FACTOR, &hasop));
4844:   if (hasop) {
4845:     PetscUseTypeMethod(mat, getfactor, type, ftype, f);
4846:     PetscFunctionReturn(PETSC_SUCCESS);
4847:   }

4849:   PetscCall(MatSolverTypeGet(type, ((PetscObject)mat)->type_name, ftype, &foundtype, &foundmtype, &conv));
4850:   if (!foundtype) {
4851:     if (type) {
4852:       SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_MISSING_FACTOR, "Could not locate solver type %s for factorization type %s and matrix type %s. Perhaps you must ./configure with --download-%s", type, MatFactorTypes[ftype],
4853:               ((PetscObject)mat)->type_name, type);
4854:     } else {
4855:       SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_MISSING_FACTOR, "Could not locate a solver type for factorization type %s and matrix type %s.", MatFactorTypes[ftype], ((PetscObject)mat)->type_name);
4856:     }
4857:   }
4858:   PetscCheck(foundmtype, PetscObjectComm((PetscObject)mat), PETSC_ERR_MISSING_FACTOR, "MatSolverType %s does not support matrix type %s", type, ((PetscObject)mat)->type_name);
4859:   PetscCheck(conv, PetscObjectComm((PetscObject)mat), PETSC_ERR_MISSING_FACTOR, "MatSolverType %s does not support factorization type %s for matrix type %s", type, MatFactorTypes[ftype], ((PetscObject)mat)->type_name);

4861:   PetscCall((*conv)(mat, ftype, f));
4862:   if (mat->factorprefix) PetscCall(MatSetOptionsPrefix(*f, mat->factorprefix));
4863:   PetscFunctionReturn(PETSC_SUCCESS);
4864: }

4866: /*@
4867:   MatGetFactorAvailable - Returns a flag if matrix supports particular type and factor type

4869:   Not Collective

4871:   Input Parameters:
4872: + mat   - the matrix
4873: . type  - name of solver type, for example, superlu, petsc (to use PETSc's default)
4874: - ftype - factor type, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ICC`, `MAT_FACTOR_ILU`, `MAT_FACTOR_QR`

4876:   Output Parameter:
4877: . flg - PETSC_TRUE if the factorization is available

4879:   Level: intermediate

4881:   Notes:
4882:   Some PETSc matrix formats have alternative solvers available that are contained in alternative packages
4883:   such as pastix, superlu, mumps etc.

4885:   PETSc must have been ./configure to use the external solver, using the option --download-package

4887:   Developer Note:
4888:   This should actually be called `MatCreateFactorAvailable()` since `MatGetFactor()` creates a new factor object

4890: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatSolverType`, `MatFactorType`, `MatGetFactor()`, `MatCopy()`, `MatDuplicate()`, `MatSolverTypeRegister()`,
4891:           `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ICC`, `MAT_FACTOR_ILU`, `MAT_FACTOR_QR`, `MatSolverTypeGet()`
4892: @*/
4893: PetscErrorCode MatGetFactorAvailable(Mat mat, MatSolverType type, MatFactorType ftype, PetscBool *flg)
4894: {
4895:   PetscErrorCode (*gconv)(Mat, MatFactorType, Mat *);

4897:   PetscFunctionBegin;
4899:   PetscAssertPointer(flg, 4);

4901:   *flg = PETSC_FALSE;
4902:   if (!((PetscObject)mat)->type_name) PetscFunctionReturn(PETSC_SUCCESS);

4904:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4905:   MatCheckPreallocated(mat, 1);

4907:   PetscCall(MatSolverTypeGet(type, ((PetscObject)mat)->type_name, ftype, NULL, NULL, &gconv));
4908:   *flg = gconv ? PETSC_TRUE : PETSC_FALSE;
4909:   PetscFunctionReturn(PETSC_SUCCESS);
4910: }

4912: /*@
4913:   MatDuplicate - Duplicates a matrix including the non-zero structure.

4915:   Collective

4917:   Input Parameters:
4918: + mat - the matrix
4919: - op  - One of `MAT_DO_NOT_COPY_VALUES`, `MAT_COPY_VALUES`, or `MAT_SHARE_NONZERO_PATTERN`.
4920:         See the manual page for `MatDuplicateOption()` for an explanation of these options.

4922:   Output Parameter:
4923: . M - pointer to place new matrix

4925:   Level: intermediate

4927:   Notes:
4928:   You cannot change the nonzero pattern for the parent or child matrix later if you use `MAT_SHARE_NONZERO_PATTERN`.

4930:   If `op` is not `MAT_COPY_VALUES` the numerical values in the new matrix are zeroed.

4932:   May be called with an unassembled input `Mat` if `MAT_DO_NOT_COPY_VALUES` is used, in which case the output `Mat` is unassembled as well.

4934:   When original mat is a product of matrix operation, e.g., an output of `MatMatMult()` or `MatCreateSubMatrix()`, only the matrix data structure of `mat`
4935:   is duplicated and the internal data structures created for the reuse of previous matrix operations are not duplicated.
4936:   User should not use `MatDuplicate()` to create new matrix `M` if `M` is intended to be reused as the product of matrix operation.

4938: .seealso: [](ch_matrices), `Mat`, `MatCopy()`, `MatConvert()`, `MatDuplicateOption`
4939: @*/
4940: PetscErrorCode MatDuplicate(Mat mat, MatDuplicateOption op, Mat *M)
4941: {
4942:   Mat         B;
4943:   VecType     vtype;
4944:   PetscInt    i;
4945:   PetscObject dm, container_h, container_d;
4946:   void (*viewf)(void);

4948:   PetscFunctionBegin;
4951:   PetscAssertPointer(M, 3);
4952:   PetscCheck(op != MAT_COPY_VALUES || mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MAT_COPY_VALUES not allowed for unassembled matrix");
4953:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
4954:   MatCheckPreallocated(mat, 1);

4956:   *M = NULL;
4957:   PetscCall(PetscLogEventBegin(MAT_Convert, mat, 0, 0, 0));
4958:   PetscUseTypeMethod(mat, duplicate, op, M);
4959:   PetscCall(PetscLogEventEnd(MAT_Convert, mat, 0, 0, 0));
4960:   B = *M;

4962:   PetscCall(MatGetOperation(mat, MATOP_VIEW, &viewf));
4963:   if (viewf) PetscCall(MatSetOperation(B, MATOP_VIEW, viewf));
4964:   PetscCall(MatGetVecType(mat, &vtype));
4965:   PetscCall(MatSetVecType(B, vtype));

4967:   B->stencil.dim = mat->stencil.dim;
4968:   B->stencil.noc = mat->stencil.noc;
4969:   for (i = 0; i <= mat->stencil.dim + (mat->stencil.noc ? 0 : -1); i++) {
4970:     B->stencil.dims[i]   = mat->stencil.dims[i];
4971:     B->stencil.starts[i] = mat->stencil.starts[i];
4972:   }

4974:   B->nooffproczerorows = mat->nooffproczerorows;
4975:   B->nooffprocentries  = mat->nooffprocentries;

4977:   PetscCall(PetscObjectQuery((PetscObject)mat, "__PETSc_dm", &dm));
4978:   if (dm) PetscCall(PetscObjectCompose((PetscObject)B, "__PETSc_dm", dm));
4979:   PetscCall(PetscObjectQuery((PetscObject)mat, "__PETSc_MatCOOStruct_Host", &container_h));
4980:   if (container_h) PetscCall(PetscObjectCompose((PetscObject)B, "__PETSc_MatCOOStruct_Host", container_h));
4981:   PetscCall(PetscObjectQuery((PetscObject)mat, "__PETSc_MatCOOStruct_Device", &container_d));
4982:   if (container_d) PetscCall(PetscObjectCompose((PetscObject)B, "__PETSc_MatCOOStruct_Device", container_d));
4983:   if (op == MAT_COPY_VALUES) PetscCall(MatPropagateSymmetryOptions(mat, B));
4984:   PetscCall(PetscObjectStateIncrease((PetscObject)B));
4985:   PetscFunctionReturn(PETSC_SUCCESS);
4986: }

4988: /*@
4989:   MatGetDiagonal - Gets the diagonal of a matrix as a `Vec`

4991:   Logically Collective

4993:   Input Parameter:
4994: . mat - the matrix

4996:   Output Parameter:
4997: . v - the diagonal of the matrix

4999:   Level: intermediate

5001:   Note:
5002:   If `mat` has local sizes `n` x `m`, this routine fills the first `ndiag = min(n, m)` entries
5003:   of `v` with the diagonal values. Thus `v` must have local size of at least `ndiag`. If `v`
5004:   is larger than `ndiag`, the values of the remaining entries are unspecified.

5006:   Currently only correct in parallel for square matrices.

5008: .seealso: [](ch_matrices), `Mat`, `Vec`, `MatGetRow()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMaxAbs()`
5009: @*/
5010: PetscErrorCode MatGetDiagonal(Mat mat, Vec v)
5011: {
5012:   PetscFunctionBegin;
5016:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5017:   MatCheckPreallocated(mat, 1);
5018:   if (PetscDefined(USE_DEBUG)) {
5019:     PetscInt nv, row, col, ndiag;

5021:     PetscCall(VecGetLocalSize(v, &nv));
5022:     PetscCall(MatGetLocalSize(mat, &row, &col));
5023:     ndiag = PetscMin(row, col);
5024:     PetscCheck(nv >= ndiag, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Nonconforming Mat and Vec. Vec local size %" PetscInt_FMT " < Mat local diagonal length %" PetscInt_FMT, nv, ndiag);
5025:   }

5027:   PetscUseTypeMethod(mat, getdiagonal, v);
5028:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
5029:   PetscFunctionReturn(PETSC_SUCCESS);
5030: }

5032: /*@
5033:   MatGetRowMin - Gets the minimum value (of the real part) of each
5034:   row of the matrix

5036:   Logically Collective

5038:   Input Parameter:
5039: . mat - the matrix

5041:   Output Parameters:
5042: + v   - the vector for storing the maximums
5043: - idx - the indices of the column found for each row (optional, pass `NULL` if not needed)

5045:   Level: intermediate

5047:   Note:
5048:   The result of this call are the same as if one converted the matrix to dense format
5049:   and found the minimum value in each row (i.e. the implicit zeros are counted as zeros).

5051:   This code is only implemented for a couple of matrix formats.

5053: .seealso: [](ch_matrices), `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMaxAbs()`, `MatGetRowMinAbs()`,
5054:           `MatGetRowMax()`
5055: @*/
5056: PetscErrorCode MatGetRowMin(Mat mat, Vec v, PetscInt idx[])
5057: {
5058:   PetscFunctionBegin;
5062:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");

5064:   if (!mat->cmap->N) {
5065:     PetscCall(VecSet(v, PETSC_MAX_REAL));
5066:     if (idx) {
5067:       PetscInt i, m = mat->rmap->n;
5068:       for (i = 0; i < m; i++) idx[i] = -1;
5069:     }
5070:   } else {
5071:     MatCheckPreallocated(mat, 1);
5072:   }
5073:   PetscUseTypeMethod(mat, getrowmin, v, idx);
5074:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
5075:   PetscFunctionReturn(PETSC_SUCCESS);
5076: }

5078: /*@
5079:   MatGetRowMinAbs - Gets the minimum value (in absolute value) of each
5080:   row of the matrix

5082:   Logically Collective

5084:   Input Parameter:
5085: . mat - the matrix

5087:   Output Parameters:
5088: + v   - the vector for storing the minimums
5089: - idx - the indices of the column found for each row (or `NULL` if not needed)

5091:   Level: intermediate

5093:   Notes:
5094:   if a row is completely empty or has only 0.0 values, then the `idx` value for that
5095:   row is 0 (the first column).

5097:   This code is only implemented for a couple of matrix formats.

5099: .seealso: [](ch_matrices), `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMax()`, `MatGetRowMaxAbs()`, `MatGetRowMin()`
5100: @*/
5101: PetscErrorCode MatGetRowMinAbs(Mat mat, Vec v, PetscInt idx[])
5102: {
5103:   PetscFunctionBegin;
5107:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5108:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");

5110:   if (!mat->cmap->N) {
5111:     PetscCall(VecSet(v, 0.0));
5112:     if (idx) {
5113:       PetscInt i, m = mat->rmap->n;
5114:       for (i = 0; i < m; i++) idx[i] = -1;
5115:     }
5116:   } else {
5117:     MatCheckPreallocated(mat, 1);
5118:     if (idx) PetscCall(PetscArrayzero(idx, mat->rmap->n));
5119:     PetscUseTypeMethod(mat, getrowminabs, v, idx);
5120:   }
5121:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
5122:   PetscFunctionReturn(PETSC_SUCCESS);
5123: }

5125: /*@
5126:   MatGetRowMax - Gets the maximum value (of the real part) of each
5127:   row of the matrix

5129:   Logically Collective

5131:   Input Parameter:
5132: . mat - the matrix

5134:   Output Parameters:
5135: + v   - the vector for storing the maximums
5136: - idx - the indices of the column found for each row (optional, otherwise pass `NULL`)

5138:   Level: intermediate

5140:   Notes:
5141:   The result of this call are the same as if one converted the matrix to dense format
5142:   and found the minimum value in each row (i.e. the implicit zeros are counted as zeros).

5144:   This code is only implemented for a couple of matrix formats.

5146: .seealso: [](ch_matrices), `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMaxAbs()`, `MatGetRowMin()`, `MatGetRowMinAbs()`
5147: @*/
5148: PetscErrorCode MatGetRowMax(Mat mat, Vec v, PetscInt idx[])
5149: {
5150:   PetscFunctionBegin;
5154:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");

5156:   if (!mat->cmap->N) {
5157:     PetscCall(VecSet(v, PETSC_MIN_REAL));
5158:     if (idx) {
5159:       PetscInt i, m = mat->rmap->n;
5160:       for (i = 0; i < m; i++) idx[i] = -1;
5161:     }
5162:   } else {
5163:     MatCheckPreallocated(mat, 1);
5164:     PetscUseTypeMethod(mat, getrowmax, v, idx);
5165:   }
5166:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
5167:   PetscFunctionReturn(PETSC_SUCCESS);
5168: }

5170: /*@
5171:   MatGetRowMaxAbs - Gets the maximum value (in absolute value) of each
5172:   row of the matrix

5174:   Logically Collective

5176:   Input Parameter:
5177: . mat - the matrix

5179:   Output Parameters:
5180: + v   - the vector for storing the maximums
5181: - idx - the indices of the column found for each row (or `NULL` if not needed)

5183:   Level: intermediate

5185:   Notes:
5186:   if a row is completely empty or has only 0.0 values, then the `idx` value for that
5187:   row is 0 (the first column).

5189:   This code is only implemented for a couple of matrix formats.

5191: .seealso: [](ch_matrices), `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowSum()`, `MatGetRowMin()`, `MatGetRowMinAbs()`
5192: @*/
5193: PetscErrorCode MatGetRowMaxAbs(Mat mat, Vec v, PetscInt idx[])
5194: {
5195:   PetscFunctionBegin;
5199:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");

5201:   if (!mat->cmap->N) {
5202:     PetscCall(VecSet(v, 0.0));
5203:     if (idx) {
5204:       PetscInt i, m = mat->rmap->n;
5205:       for (i = 0; i < m; i++) idx[i] = -1;
5206:     }
5207:   } else {
5208:     MatCheckPreallocated(mat, 1);
5209:     if (idx) PetscCall(PetscArrayzero(idx, mat->rmap->n));
5210:     PetscUseTypeMethod(mat, getrowmaxabs, v, idx);
5211:   }
5212:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
5213:   PetscFunctionReturn(PETSC_SUCCESS);
5214: }

5216: /*@
5217:   MatGetRowSumAbs - Gets the sum value (in absolute value) of each row of the matrix

5219:   Logically Collective

5221:   Input Parameter:
5222: . mat - the matrix

5224:   Output Parameter:
5225: . v - the vector for storing the sum

5227:   Level: intermediate

5229:   This code is only implemented for a couple of matrix formats.

5231: .seealso: [](ch_matrices), `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMax()`, `MatGetRowMin()`, `MatGetRowMinAbs()`
5232: @*/
5233: PetscErrorCode MatGetRowSumAbs(Mat mat, Vec v)
5234: {
5235:   PetscFunctionBegin;
5239:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");

5241:   if (!mat->cmap->N) {
5242:     PetscCall(VecSet(v, 0.0));
5243:   } else {
5244:     MatCheckPreallocated(mat, 1);
5245:     PetscUseTypeMethod(mat, getrowsumabs, v);
5246:   }
5247:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
5248:   PetscFunctionReturn(PETSC_SUCCESS);
5249: }

5251: /*@
5252:   MatGetRowSum - Gets the sum of each row of the matrix

5254:   Logically or Neighborhood Collective

5256:   Input Parameter:
5257: . mat - the matrix

5259:   Output Parameter:
5260: . v - the vector for storing the sum of rows

5262:   Level: intermediate

5264:   Note:
5265:   This code is slow since it is not currently specialized for different formats

5267: .seealso: [](ch_matrices), `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMax()`, `MatGetRowMin()`, `MatGetRowMaxAbs()`, `MatGetRowMinAbs()`, `MatGetRowSumAbs()`
5268: @*/
5269: PetscErrorCode MatGetRowSum(Mat mat, Vec v)
5270: {
5271:   Vec ones;

5273:   PetscFunctionBegin;
5277:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5278:   MatCheckPreallocated(mat, 1);
5279:   PetscCall(MatCreateVecs(mat, &ones, NULL));
5280:   PetscCall(VecSet(ones, 1.));
5281:   PetscCall(MatMult(mat, ones, v));
5282:   PetscCall(VecDestroy(&ones));
5283:   PetscFunctionReturn(PETSC_SUCCESS);
5284: }

5286: /*@
5287:   MatTransposeSetPrecursor - Set the matrix from which the second matrix will receive numerical transpose data with a call to `MatTranspose`(A,`MAT_REUSE_MATRIX`,&B)
5288:   when B was not obtained with `MatTranspose`(A,`MAT_INITIAL_MATRIX`,&B)

5290:   Collective

5292:   Input Parameter:
5293: . mat - the matrix to provide the transpose

5295:   Output Parameter:
5296: . B - the matrix to contain the transpose; it MUST have the nonzero structure of the transpose of A or the code will crash or generate incorrect results

5298:   Level: advanced

5300:   Note:
5301:   Normally the use of `MatTranspose`(A, `MAT_REUSE_MATRIX`, &B) requires that `B` was obtained with a call to `MatTranspose`(A, `MAT_INITIAL_MATRIX`, &B). This
5302:   routine allows bypassing that call.

5304: .seealso: [](ch_matrices), `Mat`, `MatTransposeSymbolic()`, `MatTranspose()`, `MatMultTranspose()`, `MatMultTransposeAdd()`, `MatIsTranspose()`, `MatReuse`, `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, `MAT_INPLACE_MATRIX`
5305: @*/
5306: PetscErrorCode MatTransposeSetPrecursor(Mat mat, Mat B)
5307: {
5308:   MatParentState *rb = NULL;

5310:   PetscFunctionBegin;
5311:   PetscCall(PetscNew(&rb));
5312:   rb->id    = ((PetscObject)mat)->id;
5313:   rb->state = 0;
5314:   PetscCall(MatGetNonzeroState(mat, &rb->nonzerostate));
5315:   PetscCall(PetscObjectContainerCompose((PetscObject)B, "MatTransposeParent", rb, PetscContainerUserDestroyDefault));
5316:   PetscFunctionReturn(PETSC_SUCCESS);
5317: }

5319: /*@
5320:   MatTranspose - Computes the transpose of a matrix, either in-place or out-of-place.

5322:   Collective

5324:   Input Parameters:
5325: + mat   - the matrix to transpose
5326: - reuse - either `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, or `MAT_INPLACE_MATRIX`

5328:   Output Parameter:
5329: . B - the transpose of the matrix

5331:   Level: intermediate

5333:   Notes:
5334:   If you use `MAT_INPLACE_MATRIX` then you must pass in `&mat` for `B`

5336:   `MAT_REUSE_MATRIX` uses the `B` matrix obtained from a previous call to this function with `MAT_INITIAL_MATRIX` to store the transpose. If you already have a matrix to contain the
5337:   transpose, call `MatTransposeSetPrecursor(mat, B)` before calling this routine.

5339:   If the nonzero structure of `mat` changed from the previous call to this function with the same matrices an error will be generated for some matrix types.

5341:   Consider using `MatCreateTranspose()` instead if you only need a matrix that behaves like the transpose but don't need the storage to be changed.
5342:   For example, the result of `MatCreateTranspose()` will compute the transpose of the given matrix times a vector for matrix-vector products computed with `MatMult()`.

5344:   If `mat` is unchanged from the last call this function returns immediately without recomputing the result

5346:   If you only need the symbolic transpose of a matrix, and not the numerical values, use `MatTransposeSymbolic()`

5348: .seealso: [](ch_matrices), `Mat`, `MatTransposeSetPrecursor()`, `MatMultTranspose()`, `MatMultTransposeAdd()`, `MatIsTranspose()`, `MatReuse`, `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, `MAT_INPLACE_MATRIX`,
5349:           `MatTransposeSymbolic()`, `MatCreateTranspose()`
5350: @*/
5351: PetscErrorCode MatTranspose(Mat mat, MatReuse reuse, Mat *B)
5352: {
5353:   PetscContainer  rB = NULL;
5354:   MatParentState *rb = NULL;

5356:   PetscFunctionBegin;
5359:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5360:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
5361:   PetscCheck(reuse != MAT_INPLACE_MATRIX || mat == *B, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "MAT_INPLACE_MATRIX requires last matrix to match first");
5362:   PetscCheck(reuse != MAT_REUSE_MATRIX || mat != *B, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Perhaps you mean MAT_INPLACE_MATRIX");
5363:   MatCheckPreallocated(mat, 1);
5364:   if (reuse == MAT_REUSE_MATRIX) {
5365:     PetscCall(PetscObjectQuery((PetscObject)*B, "MatTransposeParent", (PetscObject *)&rB));
5366:     PetscCheck(rB, PetscObjectComm((PetscObject)*B), PETSC_ERR_ARG_WRONG, "Reuse matrix used was not generated from call to MatTranspose(). Suggest MatTransposeSetPrecursor().");
5367:     PetscCall(PetscContainerGetPointer(rB, (void **)&rb));
5368:     PetscCheck(rb->id == ((PetscObject)mat)->id, PetscObjectComm((PetscObject)*B), PETSC_ERR_ARG_WRONG, "Reuse matrix used was not generated from input matrix");
5369:     if (rb->state == ((PetscObject)mat)->state) PetscFunctionReturn(PETSC_SUCCESS);
5370:   }

5372:   PetscCall(PetscLogEventBegin(MAT_Transpose, mat, 0, 0, 0));
5373:   if (reuse != MAT_INPLACE_MATRIX || mat->symmetric != PETSC_BOOL3_TRUE) {
5374:     PetscUseTypeMethod(mat, transpose, reuse, B);
5375:     PetscCall(PetscObjectStateIncrease((PetscObject)*B));
5376:   }
5377:   PetscCall(PetscLogEventEnd(MAT_Transpose, mat, 0, 0, 0));

5379:   if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatTransposeSetPrecursor(mat, *B));
5380:   if (reuse != MAT_INPLACE_MATRIX) {
5381:     PetscCall(PetscObjectQuery((PetscObject)*B, "MatTransposeParent", (PetscObject *)&rB));
5382:     PetscCall(PetscContainerGetPointer(rB, (void **)&rb));
5383:     rb->state        = ((PetscObject)mat)->state;
5384:     rb->nonzerostate = mat->nonzerostate;
5385:   }
5386:   PetscFunctionReturn(PETSC_SUCCESS);
5387: }

5389: /*@
5390:   MatTransposeSymbolic - Computes the symbolic part of the transpose of a matrix.

5392:   Collective

5394:   Input Parameter:
5395: . A - the matrix to transpose

5397:   Output Parameter:
5398: . B - the transpose. This is a complete matrix but the numerical portion is invalid. One can call `MatTranspose`(A,`MAT_REUSE_MATRIX`,&B) to compute the
5399:       numerical portion.

5401:   Level: intermediate

5403:   Note:
5404:   This is not supported for many matrix types, use `MatTranspose()` in those cases

5406: .seealso: [](ch_matrices), `Mat`, `MatTransposeSetPrecursor()`, `MatTranspose()`, `MatMultTranspose()`, `MatMultTransposeAdd()`, `MatIsTranspose()`, `MatReuse`, `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, `MAT_INPLACE_MATRIX`
5407: @*/
5408: PetscErrorCode MatTransposeSymbolic(Mat A, Mat *B)
5409: {
5410:   PetscFunctionBegin;
5413:   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5414:   PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
5415:   PetscCall(PetscLogEventBegin(MAT_Transpose, A, 0, 0, 0));
5416:   PetscUseTypeMethod(A, transposesymbolic, B);
5417:   PetscCall(PetscLogEventEnd(MAT_Transpose, A, 0, 0, 0));

5419:   PetscCall(MatTransposeSetPrecursor(A, *B));
5420:   PetscFunctionReturn(PETSC_SUCCESS);
5421: }

5423: PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat A, Mat B)
5424: {
5425:   PetscContainer  rB;
5426:   MatParentState *rb;

5428:   PetscFunctionBegin;
5431:   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5432:   PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
5433:   PetscCall(PetscObjectQuery((PetscObject)B, "MatTransposeParent", (PetscObject *)&rB));
5434:   PetscCheck(rB, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONG, "Reuse matrix used was not generated from call to MatTranspose()");
5435:   PetscCall(PetscContainerGetPointer(rB, (void **)&rb));
5436:   PetscCheck(rb->id == ((PetscObject)A)->id, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONG, "Reuse matrix used was not generated from input matrix");
5437:   PetscCheck(rb->nonzerostate == A->nonzerostate, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONGSTATE, "Reuse matrix has changed nonzero structure");
5438:   PetscFunctionReturn(PETSC_SUCCESS);
5439: }

5441: /*@
5442:   MatIsTranspose - Test whether a matrix is another one's transpose,
5443:   or its own, in which case it tests symmetry.

5445:   Collective

5447:   Input Parameters:
5448: + A   - the matrix to test
5449: . B   - the matrix to test against, this can equal the first parameter
5450: - tol - tolerance, differences between entries smaller than this are counted as zero

5452:   Output Parameter:
5453: . flg - the result

5455:   Level: intermediate

5457:   Notes:
5458:   The sequential algorithm has a running time of the order of the number of nonzeros; the parallel
5459:   test involves parallel copies of the block off-diagonal parts of the matrix.

5461: .seealso: [](ch_matrices), `Mat`, `MatTranspose()`, `MatIsSymmetric()`, `MatIsHermitian()`
5462: @*/
5463: PetscErrorCode MatIsTranspose(Mat A, Mat B, PetscReal tol, PetscBool *flg)
5464: {
5465:   PetscErrorCode (*f)(Mat, Mat, PetscReal, PetscBool *), (*g)(Mat, Mat, PetscReal, PetscBool *);

5467:   PetscFunctionBegin;
5470:   PetscAssertPointer(flg, 4);
5471:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatIsTranspose_C", &f));
5472:   PetscCall(PetscObjectQueryFunction((PetscObject)B, "MatIsTranspose_C", &g));
5473:   *flg = PETSC_FALSE;
5474:   if (f && g) {
5475:     PetscCheck(f == g, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_NOTSAMETYPE, "Matrices do not have the same comparator for symmetry test");
5476:     PetscCall((*f)(A, B, tol, flg));
5477:   } else {
5478:     MatType mattype;

5480:     PetscCall(MatGetType(f ? B : A, &mattype));
5481:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix of type %s does not support checking for transpose", mattype);
5482:   }
5483:   PetscFunctionReturn(PETSC_SUCCESS);
5484: }

5486: /*@
5487:   MatHermitianTranspose - Computes an in-place or out-of-place Hermitian transpose of a matrix in complex conjugate.

5489:   Collective

5491:   Input Parameters:
5492: + mat   - the matrix to transpose and complex conjugate
5493: - reuse - either `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, or `MAT_INPLACE_MATRIX`

5495:   Output Parameter:
5496: . B - the Hermitian transpose

5498:   Level: intermediate

5500: .seealso: [](ch_matrices), `Mat`, `MatTranspose()`, `MatMultTranspose()`, `MatMultTransposeAdd()`, `MatIsTranspose()`, `MatReuse`
5501: @*/
5502: PetscErrorCode MatHermitianTranspose(Mat mat, MatReuse reuse, Mat *B)
5503: {
5504:   PetscFunctionBegin;
5505:   PetscCall(MatTranspose(mat, reuse, B));
5506: #if defined(PETSC_USE_COMPLEX)
5507:   PetscCall(MatConjugate(*B));
5508: #endif
5509:   PetscFunctionReturn(PETSC_SUCCESS);
5510: }

5512: /*@
5513:   MatIsHermitianTranspose - Test whether a matrix is another one's Hermitian transpose,

5515:   Collective

5517:   Input Parameters:
5518: + A   - the matrix to test
5519: . B   - the matrix to test against, this can equal the first parameter
5520: - tol - tolerance, differences between entries smaller than this are counted as zero

5522:   Output Parameter:
5523: . flg - the result

5525:   Level: intermediate

5527:   Notes:
5528:   Only available for `MATAIJ` matrices.

5530:   The sequential algorithm
5531:   has a running time of the order of the number of nonzeros; the parallel
5532:   test involves parallel copies of the block off-diagonal parts of the matrix.

5534: .seealso: [](ch_matrices), `Mat`, `MatTranspose()`, `MatIsSymmetric()`, `MatIsHermitian()`, `MatIsTranspose()`
5535: @*/
5536: PetscErrorCode MatIsHermitianTranspose(Mat A, Mat B, PetscReal tol, PetscBool *flg)
5537: {
5538:   PetscErrorCode (*f)(Mat, Mat, PetscReal, PetscBool *), (*g)(Mat, Mat, PetscReal, PetscBool *);

5540:   PetscFunctionBegin;
5543:   PetscAssertPointer(flg, 4);
5544:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatIsHermitianTranspose_C", &f));
5545:   PetscCall(PetscObjectQueryFunction((PetscObject)B, "MatIsHermitianTranspose_C", &g));
5546:   if (f && g) {
5547:     PetscCheck(f != g, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_NOTSAMETYPE, "Matrices do not have the same comparator for Hermitian test");
5548:     PetscCall((*f)(A, B, tol, flg));
5549:   }
5550:   PetscFunctionReturn(PETSC_SUCCESS);
5551: }

5553: /*@
5554:   MatPermute - Creates a new matrix with rows and columns permuted from the
5555:   original.

5557:   Collective

5559:   Input Parameters:
5560: + mat - the matrix to permute
5561: . row - row permutation, each processor supplies only the permutation for its rows
5562: - col - column permutation, each processor supplies only the permutation for its columns

5564:   Output Parameter:
5565: . B - the permuted matrix

5567:   Level: advanced

5569:   Note:
5570:   The index sets map from row/col of permuted matrix to row/col of original matrix.
5571:   The index sets should be on the same communicator as mat and have the same local sizes.

5573:   Developer Note:
5574:   If you want to implement `MatPermute()` for a matrix type, and your approach doesn't
5575:   exploit the fact that row and col are permutations, consider implementing the
5576:   more general `MatCreateSubMatrix()` instead.

5578: .seealso: [](ch_matrices), `Mat`, `MatGetOrdering()`, `ISAllGather()`, `MatCreateSubMatrix()`
5579: @*/
5580: PetscErrorCode MatPermute(Mat mat, IS row, IS col, Mat *B)
5581: {
5582:   PetscFunctionBegin;
5587:   PetscAssertPointer(B, 4);
5588:   PetscCheckSameComm(mat, 1, row, 2);
5589:   if (row != col) PetscCheckSameComm(row, 2, col, 3);
5590:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5591:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
5592:   PetscCheck(mat->ops->permute || mat->ops->createsubmatrix, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatPermute not available for Mat type %s", ((PetscObject)mat)->type_name);
5593:   MatCheckPreallocated(mat, 1);

5595:   if (mat->ops->permute) {
5596:     PetscUseTypeMethod(mat, permute, row, col, B);
5597:     PetscCall(PetscObjectStateIncrease((PetscObject)*B));
5598:   } else {
5599:     PetscCall(MatCreateSubMatrix(mat, row, col, MAT_INITIAL_MATRIX, B));
5600:   }
5601:   PetscFunctionReturn(PETSC_SUCCESS);
5602: }

5604: /*@
5605:   MatEqual - Compares two matrices.

5607:   Collective

5609:   Input Parameters:
5610: + A - the first matrix
5611: - B - the second matrix

5613:   Output Parameter:
5614: . flg - `PETSC_TRUE` if the matrices are equal; `PETSC_FALSE` otherwise.

5616:   Level: intermediate

5618:   Note:
5619:   If either of the matrix is "matrix-free", meaning the matrix entries are not stored explicitly then equality is determined by comparing the results of several matrix-vector product
5620:   using several randomly created vectors, see `MatMultEqual()`.

5622: .seealso: [](ch_matrices), `Mat`, `MatMultEqual()`
5623: @*/
5624: PetscErrorCode MatEqual(Mat A, Mat B, PetscBool *flg)
5625: {
5626:   PetscFunctionBegin;
5631:   PetscAssertPointer(flg, 3);
5632:   PetscCheckSameComm(A, 1, B, 2);
5633:   MatCheckPreallocated(A, 1);
5634:   MatCheckPreallocated(B, 2);
5635:   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5636:   PetscCheck(B->assembled, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5637:   PetscCheck(A->rmap->N == B->rmap->N && A->cmap->N == B->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat B: global dim %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->N, B->rmap->N, A->cmap->N,
5638:              B->cmap->N);
5639:   if (A->ops->equal && A->ops->equal == B->ops->equal) {
5640:     PetscUseTypeMethod(A, equal, B, flg);
5641:   } else {
5642:     PetscCall(MatMultEqual(A, B, 10, flg));
5643:   }
5644:   PetscFunctionReturn(PETSC_SUCCESS);
5645: }

5647: /*@
5648:   MatDiagonalScale - Scales a matrix on the left and right by diagonal
5649:   matrices that are stored as vectors.  Either of the two scaling
5650:   matrices can be `NULL`.

5652:   Collective

5654:   Input Parameters:
5655: + mat - the matrix to be scaled
5656: . l   - the left scaling vector (or `NULL`)
5657: - r   - the right scaling vector (or `NULL`)

5659:   Level: intermediate

5661:   Note:
5662:   `MatDiagonalScale()` computes $A = LAR$, where
5663:   L = a diagonal matrix (stored as a vector), R = a diagonal matrix (stored as a vector)
5664:   The L scales the rows of the matrix, the R scales the columns of the matrix.

5666: .seealso: [](ch_matrices), `Mat`, `MatScale()`, `MatShift()`, `MatDiagonalSet()`
5667: @*/
5668: PetscErrorCode MatDiagonalScale(Mat mat, Vec l, Vec r)
5669: {
5670:   PetscFunctionBegin;
5673:   if (l) {
5675:     PetscCheckSameComm(mat, 1, l, 2);
5676:   }
5677:   if (r) {
5679:     PetscCheckSameComm(mat, 1, r, 3);
5680:   }
5681:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5682:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
5683:   MatCheckPreallocated(mat, 1);
5684:   if (!l && !r) PetscFunctionReturn(PETSC_SUCCESS);

5686:   PetscCall(PetscLogEventBegin(MAT_Scale, mat, 0, 0, 0));
5687:   PetscUseTypeMethod(mat, diagonalscale, l, r);
5688:   PetscCall(PetscLogEventEnd(MAT_Scale, mat, 0, 0, 0));
5689:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
5690:   if (l != r) mat->symmetric = PETSC_BOOL3_FALSE;
5691:   PetscFunctionReturn(PETSC_SUCCESS);
5692: }

5694: /*@
5695:   MatScale - Scales all elements of a matrix by a given number.

5697:   Logically Collective

5699:   Input Parameters:
5700: + mat - the matrix to be scaled
5701: - a   - the scaling value

5703:   Level: intermediate

5705: .seealso: [](ch_matrices), `Mat`, `MatDiagonalScale()`
5706: @*/
5707: PetscErrorCode MatScale(Mat mat, PetscScalar a)
5708: {
5709:   PetscFunctionBegin;
5712:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5713:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
5715:   MatCheckPreallocated(mat, 1);

5717:   PetscCall(PetscLogEventBegin(MAT_Scale, mat, 0, 0, 0));
5718:   if (a != (PetscScalar)1.0) {
5719:     PetscUseTypeMethod(mat, scale, a);
5720:     PetscCall(PetscObjectStateIncrease((PetscObject)mat));
5721:   }
5722:   PetscCall(PetscLogEventEnd(MAT_Scale, mat, 0, 0, 0));
5723:   PetscFunctionReturn(PETSC_SUCCESS);
5724: }

5726: /*@
5727:   MatNorm - Calculates various norms of a matrix.

5729:   Collective

5731:   Input Parameters:
5732: + mat  - the matrix
5733: - type - the type of norm, `NORM_1`, `NORM_FROBENIUS`, `NORM_INFINITY`

5735:   Output Parameter:
5736: . nrm - the resulting norm

5738:   Level: intermediate

5740: .seealso: [](ch_matrices), `Mat`
5741: @*/
5742: PetscErrorCode MatNorm(Mat mat, NormType type, PetscReal *nrm)
5743: {
5744:   PetscFunctionBegin;
5747:   PetscAssertPointer(nrm, 3);

5749:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
5750:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
5751:   MatCheckPreallocated(mat, 1);

5753:   PetscUseTypeMethod(mat, norm, type, nrm);
5754:   PetscFunctionReturn(PETSC_SUCCESS);
5755: }

5757: /*
5758:      This variable is used to prevent counting of MatAssemblyBegin() that
5759:    are called from within a MatAssemblyEnd().
5760: */
5761: static PetscInt MatAssemblyEnd_InUse = 0;
5762: /*@
5763:   MatAssemblyBegin - Begins assembling the matrix.  This routine should
5764:   be called after completing all calls to `MatSetValues()`.

5766:   Collective

5768:   Input Parameters:
5769: + mat  - the matrix
5770: - type - type of assembly, either `MAT_FLUSH_ASSEMBLY` or `MAT_FINAL_ASSEMBLY`

5772:   Level: beginner

5774:   Notes:
5775:   `MatSetValues()` generally caches the values that belong to other MPI processes.  The matrix is ready to
5776:   use only after `MatAssemblyBegin()` and `MatAssemblyEnd()` have been called.

5778:   Use `MAT_FLUSH_ASSEMBLY` when switching between `ADD_VALUES` and `INSERT_VALUES`
5779:   in `MatSetValues()`; use `MAT_FINAL_ASSEMBLY` for the final assembly before
5780:   using the matrix.

5782:   ALL processes that share a matrix MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` the SAME NUMBER of times, and each time with the
5783:   same flag of `MAT_FLUSH_ASSEMBLY` or `MAT_FINAL_ASSEMBLY` for all processes. Thus you CANNOT locally change from `ADD_VALUES` to `INSERT_VALUES`, that is
5784:   a global collective operation requiring all processes that share the matrix.

5786:   Space for preallocated nonzeros that is not filled by a call to `MatSetValues()` or a related routine are compressed
5787:   out by assembly. If you intend to use that extra space on a subsequent assembly, be sure to insert explicit zeros
5788:   before `MAT_FINAL_ASSEMBLY` so the space is not compressed out.

5790: .seealso: [](ch_matrices), `Mat`, `MatAssemblyEnd()`, `MatSetValues()`, `MatAssembled()`
5791: @*/
5792: PetscErrorCode MatAssemblyBegin(Mat mat, MatAssemblyType type)
5793: {
5794:   PetscFunctionBegin;
5797:   MatCheckPreallocated(mat, 1);
5798:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix. Did you forget to call MatSetUnfactored()?");
5799:   if (mat->assembled) {
5800:     mat->was_assembled = PETSC_TRUE;
5801:     mat->assembled     = PETSC_FALSE;
5802:   }

5804:   if (!MatAssemblyEnd_InUse) {
5805:     PetscCall(PetscLogEventBegin(MAT_AssemblyBegin, mat, 0, 0, 0));
5806:     PetscTryTypeMethod(mat, assemblybegin, type);
5807:     PetscCall(PetscLogEventEnd(MAT_AssemblyBegin, mat, 0, 0, 0));
5808:   } else PetscTryTypeMethod(mat, assemblybegin, type);
5809:   PetscFunctionReturn(PETSC_SUCCESS);
5810: }

5812: /*@
5813:   MatAssembled - Indicates if a matrix has been assembled and is ready for
5814:   use; for example, in matrix-vector product.

5816:   Not Collective

5818:   Input Parameter:
5819: . mat - the matrix

5821:   Output Parameter:
5822: . assembled - `PETSC_TRUE` or `PETSC_FALSE`

5824:   Level: advanced

5826: .seealso: [](ch_matrices), `Mat`, `MatAssemblyEnd()`, `MatSetValues()`, `MatAssemblyBegin()`
5827: @*/
5828: PetscErrorCode MatAssembled(Mat mat, PetscBool *assembled)
5829: {
5830:   PetscFunctionBegin;
5832:   PetscAssertPointer(assembled, 2);
5833:   *assembled = mat->assembled;
5834:   PetscFunctionReturn(PETSC_SUCCESS);
5835: }

5837: /*@
5838:   MatAssemblyEnd - Completes assembling the matrix.  This routine should
5839:   be called after `MatAssemblyBegin()`.

5841:   Collective

5843:   Input Parameters:
5844: + mat  - the matrix
5845: - type - type of assembly, either `MAT_FLUSH_ASSEMBLY` or `MAT_FINAL_ASSEMBLY`

5847:   Options Database Keys:
5848: + -mat_view ::ascii_info             - Prints info on matrix at conclusion of `MatAssemblyEnd()`
5849: . -mat_view ::ascii_info_detail      - Prints more detailed info
5850: . -mat_view                          - Prints matrix in ASCII format
5851: . -mat_view ::ascii_matlab           - Prints matrix in MATLAB format
5852: . -mat_view draw                     - draws nonzero structure of matrix, using `MatView()` and `PetscDrawOpenX()`.
5853: . -display <name>                    - Sets display name (default is host)
5854: . -draw_pause <sec>                  - Sets number of seconds to pause after display
5855: . -mat_view socket                   - Sends matrix to socket, can be accessed from MATLAB (See [Using MATLAB with PETSc](ch_matlab))
5856: . -viewer_socket_machine <machine>   - Machine to use for socket
5857: . -viewer_socket_port <port>         - Port number to use for socket
5858: - -mat_view binary:filename[:append] - Save matrix to file in binary format

5860:   Level: beginner

5862: .seealso: [](ch_matrices), `Mat`, `MatAssemblyBegin()`, `MatSetValues()`, `PetscDrawOpenX()`, `PetscDrawCreate()`, `MatView()`, `MatAssembled()`, `PetscViewerSocketOpen()`
5863: @*/
5864: PetscErrorCode MatAssemblyEnd(Mat mat, MatAssemblyType type)
5865: {
5866:   static PetscInt inassm = 0;
5867:   PetscBool       flg    = PETSC_FALSE;

5869:   PetscFunctionBegin;

5873:   inassm++;
5874:   MatAssemblyEnd_InUse++;
5875:   if (MatAssemblyEnd_InUse == 1) { /* Do the logging only the first time through */
5876:     PetscCall(PetscLogEventBegin(MAT_AssemblyEnd, mat, 0, 0, 0));
5877:     PetscTryTypeMethod(mat, assemblyend, type);
5878:     PetscCall(PetscLogEventEnd(MAT_AssemblyEnd, mat, 0, 0, 0));
5879:   } else PetscTryTypeMethod(mat, assemblyend, type);

5881:   /* Flush assembly is not a true assembly */
5882:   if (type != MAT_FLUSH_ASSEMBLY) {
5883:     if (mat->num_ass) {
5884:       if (!mat->symmetry_eternal) {
5885:         mat->symmetric = PETSC_BOOL3_UNKNOWN;
5886:         mat->hermitian = PETSC_BOOL3_UNKNOWN;
5887:       }
5888:       if (!mat->structural_symmetry_eternal && mat->ass_nonzerostate != mat->nonzerostate) mat->structurally_symmetric = PETSC_BOOL3_UNKNOWN;
5889:       if (!mat->spd_eternal) mat->spd = PETSC_BOOL3_UNKNOWN;
5890:     }
5891:     mat->num_ass++;
5892:     mat->assembled        = PETSC_TRUE;
5893:     mat->ass_nonzerostate = mat->nonzerostate;
5894:   }

5896:   mat->insertmode = NOT_SET_VALUES;
5897:   MatAssemblyEnd_InUse--;
5898:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
5899:   if (inassm == 1 && type != MAT_FLUSH_ASSEMBLY) {
5900:     PetscCall(MatViewFromOptions(mat, NULL, "-mat_view"));

5902:     if (mat->checksymmetryonassembly) {
5903:       PetscCall(MatIsSymmetric(mat, mat->checksymmetrytol, &flg));
5904:       if (flg) {
5905:         PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat), "Matrix is symmetric (tolerance %g)\n", (double)mat->checksymmetrytol));
5906:       } else {
5907:         PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat), "Matrix is not symmetric (tolerance %g)\n", (double)mat->checksymmetrytol));
5908:       }
5909:     }
5910:     if (mat->nullsp && mat->checknullspaceonassembly) PetscCall(MatNullSpaceTest(mat->nullsp, mat, NULL));
5911:   }
5912:   inassm--;
5913:   PetscFunctionReturn(PETSC_SUCCESS);
5914: }

5916: // PetscClangLinter pragma disable: -fdoc-section-header-unknown
5917: /*@
5918:   MatSetOption - Sets a parameter option for a matrix. Some options
5919:   may be specific to certain storage formats.  Some options
5920:   determine how values will be inserted (or added). Sorted,
5921:   row-oriented input will generally assemble the fastest. The default
5922:   is row-oriented.

5924:   Logically Collective for certain operations, such as `MAT_SPD`, not collective for `MAT_ROW_ORIENTED`, see `MatOption`

5926:   Input Parameters:
5927: + mat - the matrix
5928: . op  - the option, one of those listed below (and possibly others),
5929: - flg - turn the option on (`PETSC_TRUE`) or off (`PETSC_FALSE`)

5931:   Options Describing Matrix Structure:
5932: + `MAT_SPD`                         - symmetric positive definite
5933: . `MAT_SYMMETRIC`                   - symmetric in terms of both structure and value
5934: . `MAT_HERMITIAN`                   - transpose is the complex conjugation
5935: . `MAT_STRUCTURALLY_SYMMETRIC`      - symmetric nonzero structure
5936: . `MAT_SYMMETRY_ETERNAL`            - indicates the symmetry (or Hermitian structure) or its absence will persist through any changes to the matrix
5937: . `MAT_STRUCTURAL_SYMMETRY_ETERNAL` - indicates the structural symmetry or its absence will persist through any changes to the matrix
5938: . `MAT_SPD_ETERNAL`                 - indicates the value of `MAT_SPD` (true or false) will persist through any changes to the matrix

5940:    These are not really options of the matrix, they are knowledge about the structure of the matrix that users may provide so that they
5941:    do not need to be computed (usually at a high cost)

5943:    Options For Use with `MatSetValues()`:
5944:    Insert a logically dense subblock, which can be
5945: . `MAT_ROW_ORIENTED`                - row-oriented (default)

5947:    These options reflect the data you pass in with `MatSetValues()`; it has
5948:    nothing to do with how the data is stored internally in the matrix
5949:    data structure.

5951:    When (re)assembling a matrix, we can restrict the input for
5952:    efficiency/debugging purposes.  These options include
5953: . `MAT_NEW_NONZERO_LOCATIONS`       - additional insertions will be allowed if they generate a new nonzero (slow)
5954: . `MAT_FORCE_DIAGONAL_ENTRIES`      - forces diagonal entries to be allocated
5955: . `MAT_IGNORE_OFF_PROC_ENTRIES`     - drops off-processor entries
5956: . `MAT_NEW_NONZERO_LOCATION_ERR`    - generates an error for new matrix entry
5957: . `MAT_USE_HASH_TABLE`              - uses a hash table to speed up matrix assembly
5958: . `MAT_NO_OFF_PROC_ENTRIES`         - you know each process will only set values for its own rows, will generate an error if
5959:         any process sets values for another process. This avoids all reductions in the MatAssembly routines and thus improves
5960:         performance for very large process counts.
5961: - `MAT_SUBSET_OFF_PROC_ENTRIES`     - you know that the first assembly after setting this flag will set a superset
5962:         of the off-process entries required for all subsequent assemblies. This avoids a rendezvous step in the MatAssembly
5963:         functions, instead sending only neighbor messages.

5965:   Level: intermediate

5967:   Notes:
5968:   Except for `MAT_UNUSED_NONZERO_LOCATION_ERR` and  `MAT_ROW_ORIENTED` all processes that share the matrix must pass the same value in flg!

5970:   Some options are relevant only for particular matrix types and
5971:   are thus ignored by others.  Other options are not supported by
5972:   certain matrix types and will generate an error message if set.

5974:   If using Fortran to compute a matrix, one may need to
5975:   use the column-oriented option (or convert to the row-oriented
5976:   format).

5978:   `MAT_NEW_NONZERO_LOCATIONS` set to `PETSC_FALSE` indicates that any add or insertion
5979:   that would generate a new entry in the nonzero structure is instead
5980:   ignored.  Thus, if memory has not already been allocated for this particular
5981:   data, then the insertion is ignored. For dense matrices, in which
5982:   the entire array is allocated, no entries are ever ignored.
5983:   Set after the first `MatAssemblyEnd()`. If this option is set, then the `MatAssemblyBegin()`/`MatAssemblyEnd()` processes has one less global reduction

5985:   `MAT_NEW_NONZERO_LOCATION_ERR` set to PETSC_TRUE indicates that any add or insertion
5986:   that would generate a new entry in the nonzero structure instead produces
5987:   an error. (Currently supported for `MATAIJ` and `MATBAIJ` formats only.) If this option is set, then the `MatAssemblyBegin()`/`MatAssemblyEnd()` processes has one less global reduction

5989:   `MAT_NEW_NONZERO_ALLOCATION_ERR` set to `PETSC_TRUE` indicates that any add or insertion
5990:   that would generate a new entry that has not been preallocated will
5991:   instead produce an error. (Currently supported for `MATAIJ` and `MATBAIJ` formats
5992:   only.) This is a useful flag when debugging matrix memory preallocation.
5993:   If this option is set, then the `MatAssemblyBegin()`/`MatAssemblyEnd()` processes has one less global reduction

5995:   `MAT_IGNORE_OFF_PROC_ENTRIES` set to `PETSC_TRUE` indicates entries destined for
5996:   other processors should be dropped, rather than stashed.
5997:   This is useful if you know that the "owning" processor is also
5998:   always generating the correct matrix entries, so that PETSc need
5999:   not transfer duplicate entries generated on another processor.

6001:   `MAT_USE_HASH_TABLE` indicates that a hash table be used to improve the
6002:   searches during matrix assembly. When this flag is set, the hash table
6003:   is created during the first matrix assembly. This hash table is
6004:   used the next time through, during `MatSetValues()`/`MatSetValuesBlocked()`
6005:   to improve the searching of indices. `MAT_NEW_NONZERO_LOCATIONS` flag
6006:   should be used with `MAT_USE_HASH_TABLE` flag. This option is currently
6007:   supported by `MATMPIBAIJ` format only.

6009:   `MAT_KEEP_NONZERO_PATTERN` indicates when `MatZeroRows()` is called the zeroed entries
6010:   are kept in the nonzero structure. This flag is not used for `MatZeroRowsColumns()`

6012:   `MAT_IGNORE_ZERO_ENTRIES` - for `MATAIJ` and `MATIS` matrices this will stop zero values from creating
6013:   a zero location in the matrix

6015:   `MAT_USE_INODES` - indicates using inode version of the code - works with `MATAIJ` matrix types

6017:   `MAT_NO_OFF_PROC_ZERO_ROWS` - you know each process will only zero its own rows. This avoids all reductions in the
6018:   zero row routines and thus improves performance for very large process counts.

6020:   `MAT_IGNORE_LOWER_TRIANGULAR` - For `MATSBAIJ` matrices will ignore any insertions you make in the lower triangular
6021:   part of the matrix (since they should match the upper triangular part).

6023:   `MAT_SORTED_FULL` - each process provides exactly its local rows; all column indices for a given row are passed in a
6024:   single call to `MatSetValues()`, preallocation is perfect, row-oriented, `INSERT_VALUES` is used. Common
6025:   with finite difference schemes with non-periodic boundary conditions.

6027:   Developer Note:
6028:   `MAT_SYMMETRY_ETERNAL`, `MAT_STRUCTURAL_SYMMETRY_ETERNAL`, and `MAT_SPD_ETERNAL` are used by `MatAssemblyEnd()` and in other
6029:   places where otherwise the value of `MAT_SYMMETRIC`, `MAT_STRUCTURALLY_SYMMETRIC` or `MAT_SPD` would need to be changed back
6030:   to `PETSC_BOOL3_UNKNOWN` because the matrix values had changed so the code cannot be certain that the related property had
6031:   not changed.

6033: .seealso: [](ch_matrices), `MatOption`, `Mat`, `MatGetOption()`
6034: @*/
6035: PetscErrorCode MatSetOption(Mat mat, MatOption op, PetscBool flg)
6036: {
6037:   PetscFunctionBegin;
6039:   if (op > 0) {
6042:   }

6044:   PetscCheck(((int)op) > MAT_OPTION_MIN && ((int)op) < MAT_OPTION_MAX, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Options %d is out of range", (int)op);

6046:   switch (op) {
6047:   case MAT_FORCE_DIAGONAL_ENTRIES:
6048:     mat->force_diagonals = flg;
6049:     PetscFunctionReturn(PETSC_SUCCESS);
6050:   case MAT_NO_OFF_PROC_ENTRIES:
6051:     mat->nooffprocentries = flg;
6052:     PetscFunctionReturn(PETSC_SUCCESS);
6053:   case MAT_SUBSET_OFF_PROC_ENTRIES:
6054:     mat->assembly_subset = flg;
6055:     if (!mat->assembly_subset) { /* See the same logic in VecAssembly wrt VEC_SUBSET_OFF_PROC_ENTRIES */
6056: #if !defined(PETSC_HAVE_MPIUNI)
6057:       PetscCall(MatStashScatterDestroy_BTS(&mat->stash));
6058: #endif
6059:       mat->stash.first_assembly_done = PETSC_FALSE;
6060:     }
6061:     PetscFunctionReturn(PETSC_SUCCESS);
6062:   case MAT_NO_OFF_PROC_ZERO_ROWS:
6063:     mat->nooffproczerorows = flg;
6064:     PetscFunctionReturn(PETSC_SUCCESS);
6065:   case MAT_SPD:
6066:     if (flg) {
6067:       mat->spd                    = PETSC_BOOL3_TRUE;
6068:       mat->symmetric              = PETSC_BOOL3_TRUE;
6069:       mat->structurally_symmetric = PETSC_BOOL3_TRUE;
6070:     } else {
6071:       mat->spd = PETSC_BOOL3_FALSE;
6072:     }
6073:     break;
6074:   case MAT_SYMMETRIC:
6075:     mat->symmetric = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE;
6076:     if (flg) mat->structurally_symmetric = PETSC_BOOL3_TRUE;
6077: #if !defined(PETSC_USE_COMPLEX)
6078:     mat->hermitian = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE;
6079: #endif
6080:     break;
6081:   case MAT_HERMITIAN:
6082:     mat->hermitian = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE;
6083:     if (flg) mat->structurally_symmetric = PETSC_BOOL3_TRUE;
6084: #if !defined(PETSC_USE_COMPLEX)
6085:     mat->symmetric = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE;
6086: #endif
6087:     break;
6088:   case MAT_STRUCTURALLY_SYMMETRIC:
6089:     mat->structurally_symmetric = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE;
6090:     break;
6091:   case MAT_SYMMETRY_ETERNAL:
6092:     PetscCheck(mat->symmetric != PETSC_BOOL3_UNKNOWN, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Cannot set MAT_SYMMETRY_ETERNAL without first setting MAT_SYMMETRIC to true or false");
6093:     mat->symmetry_eternal = flg;
6094:     if (flg) mat->structural_symmetry_eternal = PETSC_TRUE;
6095:     break;
6096:   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
6097:     PetscCheck(mat->structurally_symmetric != PETSC_BOOL3_UNKNOWN, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Cannot set MAT_STRUCTURAL_SYMMETRY_ETERNAL without first setting MAT_STRUCTURALLY_SYMMETRIC to true or false");
6098:     mat->structural_symmetry_eternal = flg;
6099:     break;
6100:   case MAT_SPD_ETERNAL:
6101:     PetscCheck(mat->spd != PETSC_BOOL3_UNKNOWN, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Cannot set MAT_SPD_ETERNAL without first setting MAT_SPD to true or false");
6102:     mat->spd_eternal = flg;
6103:     if (flg) {
6104:       mat->structural_symmetry_eternal = PETSC_TRUE;
6105:       mat->symmetry_eternal            = PETSC_TRUE;
6106:     }
6107:     break;
6108:   case MAT_STRUCTURE_ONLY:
6109:     mat->structure_only = flg;
6110:     break;
6111:   case MAT_SORTED_FULL:
6112:     mat->sortedfull = flg;
6113:     break;
6114:   default:
6115:     break;
6116:   }
6117:   PetscTryTypeMethod(mat, setoption, op, flg);
6118:   PetscFunctionReturn(PETSC_SUCCESS);
6119: }

6121: /*@
6122:   MatGetOption - Gets a parameter option that has been set for a matrix.

6124:   Logically Collective

6126:   Input Parameters:
6127: + mat - the matrix
6128: - op  - the option, this only responds to certain options, check the code for which ones

6130:   Output Parameter:
6131: . flg - turn the option on (`PETSC_TRUE`) or off (`PETSC_FALSE`)

6133:   Level: intermediate

6135:   Notes:
6136:   Can only be called after `MatSetSizes()` and `MatSetType()` have been set.

6138:   Certain option values may be unknown, for those use the routines `MatIsSymmetric()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, or
6139:   `MatIsSymmetricKnown()`, `MatIsHermitianKnown()`, `MatIsStructurallySymmetricKnown()`

6141: .seealso: [](ch_matrices), `Mat`, `MatOption`, `MatSetOption()`, `MatIsSymmetric()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`,
6142:     `MatIsSymmetricKnown()`, `MatIsHermitianKnown()`, `MatIsStructurallySymmetricKnown()`
6143: @*/
6144: PetscErrorCode MatGetOption(Mat mat, MatOption op, PetscBool *flg)
6145: {
6146:   PetscFunctionBegin;

6150:   PetscCheck(((int)op) > MAT_OPTION_MIN && ((int)op) < MAT_OPTION_MAX, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Options %d is out of range", (int)op);
6151:   PetscCheck(((PetscObject)mat)->type_name, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_TYPENOTSET, "Cannot get options until type and size have been set, see MatSetType() and MatSetSizes()");

6153:   switch (op) {
6154:   case MAT_NO_OFF_PROC_ENTRIES:
6155:     *flg = mat->nooffprocentries;
6156:     break;
6157:   case MAT_NO_OFF_PROC_ZERO_ROWS:
6158:     *flg = mat->nooffproczerorows;
6159:     break;
6160:   case MAT_SYMMETRIC:
6161:     SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Use MatIsSymmetric() or MatIsSymmetricKnown()");
6162:     break;
6163:   case MAT_HERMITIAN:
6164:     SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Use MatIsHermitian() or MatIsHermitianKnown()");
6165:     break;
6166:   case MAT_STRUCTURALLY_SYMMETRIC:
6167:     SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Use MatIsStructurallySymmetric() or MatIsStructurallySymmetricKnown()");
6168:     break;
6169:   case MAT_SPD:
6170:     SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Use MatIsSPDKnown()");
6171:     break;
6172:   case MAT_SYMMETRY_ETERNAL:
6173:     *flg = mat->symmetry_eternal;
6174:     break;
6175:   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
6176:     *flg = mat->symmetry_eternal;
6177:     break;
6178:   default:
6179:     break;
6180:   }
6181:   PetscFunctionReturn(PETSC_SUCCESS);
6182: }

6184: /*@
6185:   MatZeroEntries - Zeros all entries of a matrix.  For sparse matrices
6186:   this routine retains the old nonzero structure.

6188:   Logically Collective

6190:   Input Parameter:
6191: . mat - the matrix

6193:   Level: intermediate

6195:   Note:
6196:   If the matrix was not preallocated then a default, likely poor preallocation will be set in the matrix, so this should be called after the preallocation phase.
6197:   See the Performance chapter of the users manual for information on preallocating matrices.

6199: .seealso: [](ch_matrices), `Mat`, `MatZeroRows()`, `MatZeroRowsColumns()`
6200: @*/
6201: PetscErrorCode MatZeroEntries(Mat mat)
6202: {
6203:   PetscFunctionBegin;
6206:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
6207:   PetscCheck(mat->insertmode == NOT_SET_VALUES, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for matrices where you have set values but not yet assembled");
6208:   MatCheckPreallocated(mat, 1);

6210:   PetscCall(PetscLogEventBegin(MAT_ZeroEntries, mat, 0, 0, 0));
6211:   PetscUseTypeMethod(mat, zeroentries);
6212:   PetscCall(PetscLogEventEnd(MAT_ZeroEntries, mat, 0, 0, 0));
6213:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
6214:   PetscFunctionReturn(PETSC_SUCCESS);
6215: }

6217: /*@
6218:   MatZeroRowsColumns - Zeros all entries (except possibly the main diagonal)
6219:   of a set of rows and columns of a matrix.

6221:   Collective

6223:   Input Parameters:
6224: + mat     - the matrix
6225: . numRows - the number of rows/columns to zero
6226: . rows    - the global row indices
6227: . diag    - value put in the diagonal of the eliminated rows
6228: . x       - optional vector of the solution for zeroed rows (other entries in vector are not used), these must be set before this call
6229: - b       - optional vector of the right-hand side, that will be adjusted by provided solution entries

6231:   Level: intermediate

6233:   Notes:
6234:   This routine, along with `MatZeroRows()`, is typically used to eliminate known Dirichlet boundary conditions from a linear system.

6236:   For each zeroed row, the value of the corresponding `b` is set to diag times the value of the corresponding `x`.
6237:   The other entries of `b` will be adjusted by the known values of `x` times the corresponding matrix entries in the columns that are being eliminated

6239:   If the resulting linear system is to be solved with `KSP` then one can (but does not have to) call `KSPSetInitialGuessNonzero()` to allow the
6240:   Krylov method to take advantage of the known solution on the zeroed rows.

6242:   For the parallel case, all processes that share the matrix (i.e.,
6243:   those in the communicator used for matrix creation) MUST call this
6244:   routine, regardless of whether any rows being zeroed are owned by
6245:   them.

6247:   Unlike `MatZeroRows()`, this ignores the `MAT_KEEP_NONZERO_PATTERN` option value set with `MatSetOption()`, it merely zeros those entries in the matrix, but never
6248:   removes them from the nonzero pattern. The nonzero pattern of the matrix can still change if a nonzero needs to be inserted on a diagonal entry that was previously
6249:   missing.

6251:   Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
6252:   list only rows local to itself).

6254:   The option `MAT_NO_OFF_PROC_ZERO_ROWS` does not apply to this routine.

6256: .seealso: [](ch_matrices), `Mat`, `MatZeroRowsIS()`, `MatZeroRows()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`,
6257:           `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()`
6258: @*/
6259: PetscErrorCode MatZeroRowsColumns(Mat mat, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
6260: {
6261:   PetscFunctionBegin;
6264:   if (numRows) PetscAssertPointer(rows, 3);
6265:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
6266:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
6267:   MatCheckPreallocated(mat, 1);

6269:   PetscUseTypeMethod(mat, zerorowscolumns, numRows, rows, diag, x, b);
6270:   PetscCall(MatViewFromOptions(mat, NULL, "-mat_view"));
6271:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
6272:   PetscFunctionReturn(PETSC_SUCCESS);
6273: }

6275: /*@
6276:   MatZeroRowsColumnsIS - Zeros all entries (except possibly the main diagonal)
6277:   of a set of rows and columns of a matrix.

6279:   Collective

6281:   Input Parameters:
6282: + mat  - the matrix
6283: . is   - the rows to zero
6284: . diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
6285: . x    - optional vector of solutions for zeroed rows (other entries in vector are not used)
6286: - b    - optional vector of right-hand side, that will be adjusted by provided solution

6288:   Level: intermediate

6290:   Note:
6291:   See `MatZeroRowsColumns()` for details on how this routine operates.

6293: .seealso: [](ch_matrices), `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`,
6294:           `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRows()`, `MatZeroRowsColumnsStencil()`
6295: @*/
6296: PetscErrorCode MatZeroRowsColumnsIS(Mat mat, IS is, PetscScalar diag, Vec x, Vec b)
6297: {
6298:   PetscInt        numRows;
6299:   const PetscInt *rows;

6301:   PetscFunctionBegin;
6306:   PetscCall(ISGetLocalSize(is, &numRows));
6307:   PetscCall(ISGetIndices(is, &rows));
6308:   PetscCall(MatZeroRowsColumns(mat, numRows, rows, diag, x, b));
6309:   PetscCall(ISRestoreIndices(is, &rows));
6310:   PetscFunctionReturn(PETSC_SUCCESS);
6311: }

6313: /*@
6314:   MatZeroRows - Zeros all entries (except possibly the main diagonal)
6315:   of a set of rows of a matrix.

6317:   Collective

6319:   Input Parameters:
6320: + mat     - the matrix
6321: . numRows - the number of rows to zero
6322: . rows    - the global row indices
6323: . diag    - value put in the diagonal of the zeroed rows
6324: . x       - optional vector of solutions for zeroed rows (other entries in vector are not used), these must be set before this call
6325: - b       - optional vector of right-hand side, that will be adjusted by provided solution entries

6327:   Level: intermediate

6329:   Notes:
6330:   This routine, along with `MatZeroRowsColumns()`, is typically used to eliminate known Dirichlet boundary conditions from a linear system.

6332:   For each zeroed row, the value of the corresponding `b` is set to `diag` times the value of the corresponding `x`.

6334:   If the resulting linear system is to be solved with `KSP` then one can (but does not have to) call `KSPSetInitialGuessNonzero()` to allow the
6335:   Krylov method to take advantage of the known solution on the zeroed rows.

6337:   May be followed by using a `PC` of type `PCREDISTRIBUTE` to solve the reduced problem (`PCDISTRIBUTE` completely eliminates the zeroed rows and their corresponding columns)
6338:   from the matrix.

6340:   Unlike `MatZeroRowsColumns()` for the `MATAIJ` and `MATBAIJ` matrix formats this removes the old nonzero structure, from the eliminated rows of the matrix
6341:   but does not release memory.  Because of this removal matrix-vector products with the adjusted matrix will be a bit faster. For the dense
6342:   formats this does not alter the nonzero structure.

6344:   If the option `MatSetOption`(mat,`MAT_KEEP_NONZERO_PATTERN`,`PETSC_TRUE`) the nonzero structure
6345:   of the matrix is not changed the values are
6346:   merely zeroed.

6348:   The user can set a value in the diagonal entry (or for the `MATAIJ` format
6349:   formats can optionally remove the main diagonal entry from the
6350:   nonzero structure as well, by passing 0.0 as the final argument).

6352:   For the parallel case, all processes that share the matrix (i.e.,
6353:   those in the communicator used for matrix creation) MUST call this
6354:   routine, regardless of whether any rows being zeroed are owned by
6355:   them.

6357:   Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to
6358:   list only rows local to itself).

6360:   You can call `MatSetOption`(mat,`MAT_NO_OFF_PROC_ZERO_ROWS`,`PETSC_TRUE`) if each process indicates only rows it
6361:   owns that are to be zeroed. This saves a global synchronization in the implementation.

6363: .seealso: [](ch_matrices), `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`,
6364:           `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()`, `PCREDISTRIBUTE`, `MAT_KEEP_NONZERO_PATTERN`
6365: @*/
6366: PetscErrorCode MatZeroRows(Mat mat, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
6367: {
6368:   PetscFunctionBegin;
6371:   if (numRows) PetscAssertPointer(rows, 3);
6372:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
6373:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
6374:   MatCheckPreallocated(mat, 1);

6376:   PetscUseTypeMethod(mat, zerorows, numRows, rows, diag, x, b);
6377:   PetscCall(MatViewFromOptions(mat, NULL, "-mat_view"));
6378:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
6379:   PetscFunctionReturn(PETSC_SUCCESS);
6380: }

6382: /*@
6383:   MatZeroRowsIS - Zeros all entries (except possibly the main diagonal)
6384:   of a set of rows of a matrix indicated by an `IS`

6386:   Collective

6388:   Input Parameters:
6389: + mat  - the matrix
6390: . is   - index set, `IS`, of rows to remove (if `NULL` then no row is removed)
6391: . diag - value put in all diagonals of eliminated rows
6392: . x    - optional vector of solutions for zeroed rows (other entries in vector are not used)
6393: - b    - optional vector of right-hand side, that will be adjusted by provided solution

6395:   Level: intermediate

6397:   Note:
6398:   See `MatZeroRows()` for details on how this routine operates.

6400: .seealso: [](ch_matrices), `Mat`, `MatZeroRows()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`,
6401:           `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()`, `IS`
6402: @*/
6403: PetscErrorCode MatZeroRowsIS(Mat mat, IS is, PetscScalar diag, Vec x, Vec b)
6404: {
6405:   PetscInt        numRows = 0;
6406:   const PetscInt *rows    = NULL;

6408:   PetscFunctionBegin;
6411:   if (is) {
6413:     PetscCall(ISGetLocalSize(is, &numRows));
6414:     PetscCall(ISGetIndices(is, &rows));
6415:   }
6416:   PetscCall(MatZeroRows(mat, numRows, rows, diag, x, b));
6417:   if (is) PetscCall(ISRestoreIndices(is, &rows));
6418:   PetscFunctionReturn(PETSC_SUCCESS);
6419: }

6421: /*@
6422:   MatZeroRowsStencil - Zeros all entries (except possibly the main diagonal)
6423:   of a set of rows of a matrix indicated by a `MatStencil`. These rows must be local to the process.

6425:   Collective

6427:   Input Parameters:
6428: + mat     - the matrix
6429: . numRows - the number of rows to remove
6430: . rows    - the grid coordinates (and component number when dof > 1) for matrix rows indicated by an array of `MatStencil`
6431: . diag    - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
6432: . x       - optional vector of solutions for zeroed rows (other entries in vector are not used)
6433: - b       - optional vector of right-hand side, that will be adjusted by provided solution

6435:   Level: intermediate

6437:   Notes:
6438:   See `MatZeroRows()` for details on how this routine operates.

6440:   The grid coordinates are across the entire grid, not just the local portion

6442:   For periodic boundary conditions use negative indices for values to the left (below 0; that are to be
6443:   obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one
6444:   etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the
6445:   `DM_BOUNDARY_PERIODIC` boundary type.

6447:   For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have
6448:   a single value per point) you can skip filling those indices.

6450:   Fortran Note:
6451:   `idxm` and `idxn` should be declared as
6452: $     MatStencil idxm(4, m)
6453:   and the values inserted using
6454: .vb
6455:     idxm(MatStencil_i, 1) = i
6456:     idxm(MatStencil_j, 1) = j
6457:     idxm(MatStencil_k, 1) = k
6458:     idxm(MatStencil_c, 1) = c
6459:    etc
6460: .ve

6462: .seealso: [](ch_matrices), `Mat`, `MatStencil`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRows()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`,
6463:           `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()`
6464: @*/
6465: PetscErrorCode MatZeroRowsStencil(Mat mat, PetscInt numRows, const MatStencil rows[], PetscScalar diag, Vec x, Vec b)
6466: {
6467:   PetscInt  dim    = mat->stencil.dim;
6468:   PetscInt  sdim   = dim - (1 - (PetscInt)mat->stencil.noc);
6469:   PetscInt *dims   = mat->stencil.dims + 1;
6470:   PetscInt *starts = mat->stencil.starts;
6471:   PetscInt *dxm    = (PetscInt *)rows;
6472:   PetscInt *jdxm, i, j, tmp, numNewRows = 0;

6474:   PetscFunctionBegin;
6477:   if (numRows) PetscAssertPointer(rows, 3);

6479:   PetscCall(PetscMalloc1(numRows, &jdxm));
6480:   for (i = 0; i < numRows; ++i) {
6481:     /* Skip unused dimensions (they are ordered k, j, i, c) */
6482:     for (j = 0; j < 3 - sdim; ++j) dxm++;
6483:     /* Local index in X dir */
6484:     tmp = *dxm++ - starts[0];
6485:     /* Loop over remaining dimensions */
6486:     for (j = 0; j < dim - 1; ++j) {
6487:       /* If nonlocal, set index to be negative */
6488:       if ((*dxm++ - starts[j + 1]) < 0 || tmp < 0) tmp = PETSC_INT_MIN;
6489:       /* Update local index */
6490:       else tmp = tmp * dims[j] + *(dxm - 1) - starts[j + 1];
6491:     }
6492:     /* Skip component slot if necessary */
6493:     if (mat->stencil.noc) dxm++;
6494:     /* Local row number */
6495:     if (tmp >= 0) jdxm[numNewRows++] = tmp;
6496:   }
6497:   PetscCall(MatZeroRowsLocal(mat, numNewRows, jdxm, diag, x, b));
6498:   PetscCall(PetscFree(jdxm));
6499:   PetscFunctionReturn(PETSC_SUCCESS);
6500: }

6502: /*@
6503:   MatZeroRowsColumnsStencil - Zeros all row and column entries (except possibly the main diagonal)
6504:   of a set of rows and columns of a matrix.

6506:   Collective

6508:   Input Parameters:
6509: + mat     - the matrix
6510: . numRows - the number of rows/columns to remove
6511: . rows    - the grid coordinates (and component number when dof > 1) for matrix rows
6512: . diag    - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry)
6513: . x       - optional vector of solutions for zeroed rows (other entries in vector are not used)
6514: - b       - optional vector of right-hand side, that will be adjusted by provided solution

6516:   Level: intermediate

6518:   Notes:
6519:   See `MatZeroRowsColumns()` for details on how this routine operates.

6521:   The grid coordinates are across the entire grid, not just the local portion

6523:   For periodic boundary conditions use negative indices for values to the left (below 0; that are to be
6524:   obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one
6525:   etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the
6526:   `DM_BOUNDARY_PERIODIC` boundary type.

6528:   For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have
6529:   a single value per point) you can skip filling those indices.

6531:   Fortran Note:
6532:   `idxm` and `idxn` should be declared as
6533: $     MatStencil idxm(4, m)
6534:   and the values inserted using
6535: .vb
6536:     idxm(MatStencil_i, 1) = i
6537:     idxm(MatStencil_j, 1) = j
6538:     idxm(MatStencil_k, 1) = k
6539:     idxm(MatStencil_c, 1) = c
6540:     etc
6541: .ve

6543: .seealso: [](ch_matrices), `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`,
6544:           `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRows()`
6545: @*/
6546: PetscErrorCode MatZeroRowsColumnsStencil(Mat mat, PetscInt numRows, const MatStencil rows[], PetscScalar diag, Vec x, Vec b)
6547: {
6548:   PetscInt  dim    = mat->stencil.dim;
6549:   PetscInt  sdim   = dim - (1 - (PetscInt)mat->stencil.noc);
6550:   PetscInt *dims   = mat->stencil.dims + 1;
6551:   PetscInt *starts = mat->stencil.starts;
6552:   PetscInt *dxm    = (PetscInt *)rows;
6553:   PetscInt *jdxm, i, j, tmp, numNewRows = 0;

6555:   PetscFunctionBegin;
6558:   if (numRows) PetscAssertPointer(rows, 3);

6560:   PetscCall(PetscMalloc1(numRows, &jdxm));
6561:   for (i = 0; i < numRows; ++i) {
6562:     /* Skip unused dimensions (they are ordered k, j, i, c) */
6563:     for (j = 0; j < 3 - sdim; ++j) dxm++;
6564:     /* Local index in X dir */
6565:     tmp = *dxm++ - starts[0];
6566:     /* Loop over remaining dimensions */
6567:     for (j = 0; j < dim - 1; ++j) {
6568:       /* If nonlocal, set index to be negative */
6569:       if ((*dxm++ - starts[j + 1]) < 0 || tmp < 0) tmp = PETSC_INT_MIN;
6570:       /* Update local index */
6571:       else tmp = tmp * dims[j] + *(dxm - 1) - starts[j + 1];
6572:     }
6573:     /* Skip component slot if necessary */
6574:     if (mat->stencil.noc) dxm++;
6575:     /* Local row number */
6576:     if (tmp >= 0) jdxm[numNewRows++] = tmp;
6577:   }
6578:   PetscCall(MatZeroRowsColumnsLocal(mat, numNewRows, jdxm, diag, x, b));
6579:   PetscCall(PetscFree(jdxm));
6580:   PetscFunctionReturn(PETSC_SUCCESS);
6581: }

6583: /*@
6584:   MatZeroRowsLocal - Zeros all entries (except possibly the main diagonal)
6585:   of a set of rows of a matrix; using local numbering of rows.

6587:   Collective

6589:   Input Parameters:
6590: + mat     - the matrix
6591: . numRows - the number of rows to remove
6592: . rows    - the local row indices
6593: . diag    - value put in all diagonals of eliminated rows
6594: . x       - optional vector of solutions for zeroed rows (other entries in vector are not used)
6595: - b       - optional vector of right-hand side, that will be adjusted by provided solution

6597:   Level: intermediate

6599:   Notes:
6600:   Before calling `MatZeroRowsLocal()`, the user must first set the
6601:   local-to-global mapping by calling MatSetLocalToGlobalMapping(), this is often already set for matrices obtained with `DMCreateMatrix()`.

6603:   See `MatZeroRows()` for details on how this routine operates.

6605: .seealso: [](ch_matrices), `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRows()`, `MatSetOption()`,
6606:           `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()`
6607: @*/
6608: PetscErrorCode MatZeroRowsLocal(Mat mat, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
6609: {
6610:   PetscFunctionBegin;
6613:   if (numRows) PetscAssertPointer(rows, 3);
6614:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
6615:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
6616:   MatCheckPreallocated(mat, 1);

6618:   if (mat->ops->zerorowslocal) {
6619:     PetscUseTypeMethod(mat, zerorowslocal, numRows, rows, diag, x, b);
6620:   } else {
6621:     IS              is, newis;
6622:     const PetscInt *newRows;

6624:     PetscCheck(mat->rmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Need to provide local to global mapping to matrix first");
6625:     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numRows, rows, PETSC_COPY_VALUES, &is));
6626:     PetscCall(ISLocalToGlobalMappingApplyIS(mat->rmap->mapping, is, &newis));
6627:     PetscCall(ISGetIndices(newis, &newRows));
6628:     PetscUseTypeMethod(mat, zerorows, numRows, newRows, diag, x, b);
6629:     PetscCall(ISRestoreIndices(newis, &newRows));
6630:     PetscCall(ISDestroy(&newis));
6631:     PetscCall(ISDestroy(&is));
6632:   }
6633:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
6634:   PetscFunctionReturn(PETSC_SUCCESS);
6635: }

6637: /*@
6638:   MatZeroRowsLocalIS - Zeros all entries (except possibly the main diagonal)
6639:   of a set of rows of a matrix; using local numbering of rows.

6641:   Collective

6643:   Input Parameters:
6644: + mat  - the matrix
6645: . is   - index set of rows to remove
6646: . diag - value put in all diagonals of eliminated rows
6647: . x    - optional vector of solutions for zeroed rows (other entries in vector are not used)
6648: - b    - optional vector of right-hand side, that will be adjusted by provided solution

6650:   Level: intermediate

6652:   Notes:
6653:   Before calling `MatZeroRowsLocalIS()`, the user must first set the
6654:   local-to-global mapping by calling `MatSetLocalToGlobalMapping()`, this is often already set for matrices obtained with `DMCreateMatrix()`.

6656:   See `MatZeroRows()` for details on how this routine operates.

6658: .seealso: [](ch_matrices), `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRows()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`,
6659:           `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()`
6660: @*/
6661: PetscErrorCode MatZeroRowsLocalIS(Mat mat, IS is, PetscScalar diag, Vec x, Vec b)
6662: {
6663:   PetscInt        numRows;
6664:   const PetscInt *rows;

6666:   PetscFunctionBegin;
6670:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
6671:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
6672:   MatCheckPreallocated(mat, 1);

6674:   PetscCall(ISGetLocalSize(is, &numRows));
6675:   PetscCall(ISGetIndices(is, &rows));
6676:   PetscCall(MatZeroRowsLocal(mat, numRows, rows, diag, x, b));
6677:   PetscCall(ISRestoreIndices(is, &rows));
6678:   PetscFunctionReturn(PETSC_SUCCESS);
6679: }

6681: /*@
6682:   MatZeroRowsColumnsLocal - Zeros all entries (except possibly the main diagonal)
6683:   of a set of rows and columns of a matrix; using local numbering of rows.

6685:   Collective

6687:   Input Parameters:
6688: + mat     - the matrix
6689: . numRows - the number of rows to remove
6690: . rows    - the global row indices
6691: . diag    - value put in all diagonals of eliminated rows
6692: . x       - optional vector of solutions for zeroed rows (other entries in vector are not used)
6693: - b       - optional vector of right-hand side, that will be adjusted by provided solution

6695:   Level: intermediate

6697:   Notes:
6698:   Before calling `MatZeroRowsColumnsLocal()`, the user must first set the
6699:   local-to-global mapping by calling `MatSetLocalToGlobalMapping()`, this is often already set for matrices obtained with `DMCreateMatrix()`.

6701:   See `MatZeroRowsColumns()` for details on how this routine operates.

6703: .seealso: [](ch_matrices), `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`,
6704:           `MatZeroRows()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()`
6705: @*/
6706: PetscErrorCode MatZeroRowsColumnsLocal(Mat mat, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
6707: {
6708:   IS              is, newis;
6709:   const PetscInt *newRows;

6711:   PetscFunctionBegin;
6714:   if (numRows) PetscAssertPointer(rows, 3);
6715:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
6716:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
6717:   MatCheckPreallocated(mat, 1);

6719:   PetscCheck(mat->cmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Need to provide local to global mapping to matrix first");
6720:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numRows, rows, PETSC_COPY_VALUES, &is));
6721:   PetscCall(ISLocalToGlobalMappingApplyIS(mat->cmap->mapping, is, &newis));
6722:   PetscCall(ISGetIndices(newis, &newRows));
6723:   PetscUseTypeMethod(mat, zerorowscolumns, numRows, newRows, diag, x, b);
6724:   PetscCall(ISRestoreIndices(newis, &newRows));
6725:   PetscCall(ISDestroy(&newis));
6726:   PetscCall(ISDestroy(&is));
6727:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
6728:   PetscFunctionReturn(PETSC_SUCCESS);
6729: }

6731: /*@
6732:   MatZeroRowsColumnsLocalIS - Zeros all entries (except possibly the main diagonal)
6733:   of a set of rows and columns of a matrix; using local numbering of rows.

6735:   Collective

6737:   Input Parameters:
6738: + mat  - the matrix
6739: . is   - index set of rows to remove
6740: . diag - value put in all diagonals of eliminated rows
6741: . x    - optional vector of solutions for zeroed rows (other entries in vector are not used)
6742: - b    - optional vector of right-hand side, that will be adjusted by provided solution

6744:   Level: intermediate

6746:   Notes:
6747:   Before calling `MatZeroRowsColumnsLocalIS()`, the user must first set the
6748:   local-to-global mapping by calling `MatSetLocalToGlobalMapping()`, this is often already set for matrices obtained with `DMCreateMatrix()`.

6750:   See `MatZeroRowsColumns()` for details on how this routine operates.

6752: .seealso: [](ch_matrices), `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`,
6753:           `MatZeroRowsColumnsLocal()`, `MatZeroRows()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()`
6754: @*/
6755: PetscErrorCode MatZeroRowsColumnsLocalIS(Mat mat, IS is, PetscScalar diag, Vec x, Vec b)
6756: {
6757:   PetscInt        numRows;
6758:   const PetscInt *rows;

6760:   PetscFunctionBegin;
6764:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
6765:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
6766:   MatCheckPreallocated(mat, 1);

6768:   PetscCall(ISGetLocalSize(is, &numRows));
6769:   PetscCall(ISGetIndices(is, &rows));
6770:   PetscCall(MatZeroRowsColumnsLocal(mat, numRows, rows, diag, x, b));
6771:   PetscCall(ISRestoreIndices(is, &rows));
6772:   PetscFunctionReturn(PETSC_SUCCESS);
6773: }

6775: /*@
6776:   MatGetSize - Returns the numbers of rows and columns in a matrix.

6778:   Not Collective

6780:   Input Parameter:
6781: . mat - the matrix

6783:   Output Parameters:
6784: + m - the number of global rows
6785: - n - the number of global columns

6787:   Level: beginner

6789:   Note:
6790:   Both output parameters can be `NULL` on input.

6792: .seealso: [](ch_matrices), `Mat`, `MatSetSizes()`, `MatGetLocalSize()`
6793: @*/
6794: PetscErrorCode MatGetSize(Mat mat, PetscInt *m, PetscInt *n)
6795: {
6796:   PetscFunctionBegin;
6798:   if (m) *m = mat->rmap->N;
6799:   if (n) *n = mat->cmap->N;
6800:   PetscFunctionReturn(PETSC_SUCCESS);
6801: }

6803: /*@
6804:   MatGetLocalSize - For most matrix formats, excluding `MATELEMENTAL` and `MATSCALAPACK`, Returns the number of local rows and local columns
6805:   of a matrix. For all matrices this is the local size of the left and right vectors as returned by `MatCreateVecs()`.

6807:   Not Collective

6809:   Input Parameter:
6810: . mat - the matrix

6812:   Output Parameters:
6813: + m - the number of local rows, use `NULL` to not obtain this value
6814: - n - the number of local columns, use `NULL` to not obtain this value

6816:   Level: beginner

6818: .seealso: [](ch_matrices), `Mat`, `MatSetSizes()`, `MatGetSize()`
6819: @*/
6820: PetscErrorCode MatGetLocalSize(Mat mat, PetscInt *m, PetscInt *n)
6821: {
6822:   PetscFunctionBegin;
6824:   if (m) PetscAssertPointer(m, 2);
6825:   if (n) PetscAssertPointer(n, 3);
6826:   if (m) *m = mat->rmap->n;
6827:   if (n) *n = mat->cmap->n;
6828:   PetscFunctionReturn(PETSC_SUCCESS);
6829: }

6831: /*@
6832:   MatGetOwnershipRangeColumn - Returns the range of matrix columns associated with rows of a
6833:   vector one multiplies this matrix by that are owned by this processor.

6835:   Not Collective, unless matrix has not been allocated, then collective

6837:   Input Parameter:
6838: . mat - the matrix

6840:   Output Parameters:
6841: + m - the global index of the first local column, use `NULL` to not obtain this value
6842: - n - one more than the global index of the last local column, use `NULL` to not obtain this value

6844:   Level: developer

6846:   Notes:
6847:   If the `Mat` was obtained from a `DM` with `DMCreateMatrix()`, then the range values are determined by the specific `DM`.

6849:   If the `Mat` was created directly the range values are determined by the local size passed to `MatSetSizes()` or `MatCreateAIJ()`.
6850:   If `PETSC_DECIDE` was passed as the local size, then the vector uses default values for the range using `PetscSplitOwnership()`.

6852:   For certain `DM`, such as `DMDA`, it is better to use `DM` specific routines, such as `DMDAGetGhostCorners()`, to determine
6853:   the local values in the matrix.

6855:   Returns the columns of the "diagonal block" for most sparse matrix formats. See [Matrix
6856:   Layouts](sec_matlayout) for details on matrix layouts.

6858: .seealso: [](ch_matrices), `Mat`, `MatGetOwnershipRange()`, `MatGetOwnershipRanges()`, `MatGetOwnershipRangesColumn()`, `PetscLayout`,
6859:           `MatSetSizes()`, `MatCreateAIJ()`, `DMDAGetGhostCorners()`, `DM`
6860: @*/
6861: PetscErrorCode MatGetOwnershipRangeColumn(Mat mat, PetscInt *m, PetscInt *n)
6862: {
6863:   PetscFunctionBegin;
6866:   if (m) PetscAssertPointer(m, 2);
6867:   if (n) PetscAssertPointer(n, 3);
6868:   MatCheckPreallocated(mat, 1);
6869:   if (m) *m = mat->cmap->rstart;
6870:   if (n) *n = mat->cmap->rend;
6871:   PetscFunctionReturn(PETSC_SUCCESS);
6872: }

6874: /*@
6875:   MatGetOwnershipRange - For matrices that own values by row, excludes `MATELEMENTAL` and `MATSCALAPACK`, returns the range of matrix rows owned by
6876:   this MPI process.

6878:   Not Collective

6880:   Input Parameter:
6881: . mat - the matrix

6883:   Output Parameters:
6884: + m - the global index of the first local row, use `NULL` to not obtain this value
6885: - n - one more than the global index of the last local row, use `NULL` to not obtain this value

6887:   Level: beginner

6889:   Notes:
6890:   If the `Mat` was obtained from a `DM` with `DMCreateMatrix()`, then the range values are determined by the specific `DM`.

6892:   If the `Mat` was created directly the range values are determined by the local size passed to `MatSetSizes()` or `MatCreateAIJ()`.
6893:   If `PETSC_DECIDE` was passed as the local size, then the vector uses default values for the range using `PetscSplitOwnership()`.

6895:   For certain `DM`, such as `DMDA`, it is better to use `DM` specific routines, such as `DMDAGetGhostCorners()`, to determine
6896:   the local values in the matrix.

6898:   The high argument is one more than the last element stored locally.

6900:   For all matrices  it returns the range of matrix rows associated with rows of a vector that
6901:   would contain the result of a matrix vector product with this matrix. See [Matrix
6902:   Layouts](sec_matlayout) for details on matrix layouts.

6904: .seealso: [](ch_matrices), `Mat`, `MatGetOwnershipRanges()`, `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRangesColumn()`, `PetscSplitOwnership()`,
6905:           `PetscSplitOwnershipBlock()`, `PetscLayout`, `MatSetSizes()`, `MatCreateAIJ()`, `DMDAGetGhostCorners()`, `DM`
6906: @*/
6907: PetscErrorCode MatGetOwnershipRange(Mat mat, PetscInt *m, PetscInt *n)
6908: {
6909:   PetscFunctionBegin;
6912:   if (m) PetscAssertPointer(m, 2);
6913:   if (n) PetscAssertPointer(n, 3);
6914:   MatCheckPreallocated(mat, 1);
6915:   if (m) *m = mat->rmap->rstart;
6916:   if (n) *n = mat->rmap->rend;
6917:   PetscFunctionReturn(PETSC_SUCCESS);
6918: }

6920: /*@C
6921:   MatGetOwnershipRanges - For matrices that own values by row, excludes `MATELEMENTAL` and
6922:   `MATSCALAPACK`, returns the range of matrix rows owned by each process.

6924:   Not Collective, unless matrix has not been allocated

6926:   Input Parameter:
6927: . mat - the matrix

6929:   Output Parameter:
6930: . ranges - start of each processors portion plus one more than the total length at the end, of length `size` + 1
6931:            where `size` is the number of MPI processes used by `mat`

6933:   Level: beginner

6935:   Notes:
6936:   If the `Mat` was obtained from a `DM` with `DMCreateMatrix()`, then the range values are determined by the specific `DM`.

6938:   If the `Mat` was created directly the range values are determined by the local size passed to `MatSetSizes()` or `MatCreateAIJ()`.
6939:   If `PETSC_DECIDE` was passed as the local size, then the vector uses default values for the range using `PetscSplitOwnership()`.

6941:   For certain `DM`, such as `DMDA`, it is better to use `DM` specific routines, such as `DMDAGetGhostCorners()`, to determine
6942:   the local values in the matrix.

6944:   For all matrices  it returns the ranges of matrix rows associated with rows of a vector that
6945:   would contain the result of a matrix vector product with this matrix. See [Matrix
6946:   Layouts](sec_matlayout) for details on matrix layouts.

6948: .seealso: [](ch_matrices), `Mat`, `MatGetOwnershipRange()`, `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRangesColumn()`, `PetscLayout`,
6949:           `PetscSplitOwnership()`, `PetscSplitOwnershipBlock()`, `MatSetSizes()`, `MatCreateAIJ()`,
6950:           `DMDAGetGhostCorners()`, `DM`
6951: @*/
6952: PetscErrorCode MatGetOwnershipRanges(Mat mat, const PetscInt *ranges[])
6953: {
6954:   PetscFunctionBegin;
6957:   MatCheckPreallocated(mat, 1);
6958:   PetscCall(PetscLayoutGetRanges(mat->rmap, ranges));
6959:   PetscFunctionReturn(PETSC_SUCCESS);
6960: }

6962: /*@C
6963:   MatGetOwnershipRangesColumn - Returns the ranges of matrix columns associated with rows of a
6964:   vector one multiplies this vector by that are owned by each processor.

6966:   Not Collective, unless matrix has not been allocated

6968:   Input Parameter:
6969: . mat - the matrix

6971:   Output Parameter:
6972: . ranges - start of each processors portion plus one more than the total length at the end

6974:   Level: beginner

6976:   Notes:
6977:   If the `Mat` was obtained from a `DM` with `DMCreateMatrix()`, then the range values are determined by the specific `DM`.

6979:   If the `Mat` was created directly the range values are determined by the local size passed to `MatSetSizes()` or `MatCreateAIJ()`.
6980:   If `PETSC_DECIDE` was passed as the local size, then the vector uses default values for the range using `PetscSplitOwnership()`.

6982:   For certain `DM`, such as `DMDA`, it is better to use `DM` specific routines, such as `DMDAGetGhostCorners()`, to determine
6983:   the local values in the matrix.

6985:   Returns the columns of the "diagonal blocks", for most sparse matrix formats. See [Matrix
6986:   Layouts](sec_matlayout) for details on matrix layouts.

6988: .seealso: [](ch_matrices), `Mat`, `MatGetOwnershipRange()`, `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRanges()`,
6989:           `PetscSplitOwnership()`, `PetscSplitOwnershipBlock()`, `PetscLayout`, `MatSetSizes()`, `MatCreateAIJ()`,
6990:           `DMDAGetGhostCorners()`, `DM`
6991: @*/
6992: PetscErrorCode MatGetOwnershipRangesColumn(Mat mat, const PetscInt *ranges[])
6993: {
6994:   PetscFunctionBegin;
6997:   MatCheckPreallocated(mat, 1);
6998:   PetscCall(PetscLayoutGetRanges(mat->cmap, ranges));
6999:   PetscFunctionReturn(PETSC_SUCCESS);
7000: }

7002: /*@
7003:   MatGetOwnershipIS - Get row and column ownership of a matrices' values as index sets.

7005:   Not Collective

7007:   Input Parameter:
7008: . A - matrix

7010:   Output Parameters:
7011: + rows - rows in which this process owns elements, , use `NULL` to not obtain this value
7012: - cols - columns in which this process owns elements, use `NULL` to not obtain this value

7014:   Level: intermediate

7016:   Note:
7017:   You should call `ISDestroy()` on the returned `IS`

7019:   For most matrices, excluding `MATELEMENTAL` and `MATSCALAPACK`, this corresponds to values
7020:   returned by `MatGetOwnershipRange()`, `MatGetOwnershipRangeColumn()`. For `MATELEMENTAL` and
7021:   `MATSCALAPACK` the ownership is more complicated. See [Matrix Layouts](sec_matlayout) for
7022:   details on matrix layouts.

7024: .seealso: [](ch_matrices), `IS`, `Mat`, `MatGetOwnershipRanges()`, `MatSetValues()`, `MATELEMENTAL`, `MATSCALAPACK`
7025: @*/
7026: PetscErrorCode MatGetOwnershipIS(Mat A, IS *rows, IS *cols)
7027: {
7028:   PetscErrorCode (*f)(Mat, IS *, IS *);

7030:   PetscFunctionBegin;
7033:   MatCheckPreallocated(A, 1);
7034:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatGetOwnershipIS_C", &f));
7035:   if (f) {
7036:     PetscCall((*f)(A, rows, cols));
7037:   } else { /* Create a standard row-based partition, each process is responsible for ALL columns in their row block */
7038:     if (rows) PetscCall(ISCreateStride(PETSC_COMM_SELF, A->rmap->n, A->rmap->rstart, 1, rows));
7039:     if (cols) PetscCall(ISCreateStride(PETSC_COMM_SELF, A->cmap->N, 0, 1, cols));
7040:   }
7041:   PetscFunctionReturn(PETSC_SUCCESS);
7042: }

7044: /*@
7045:   MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix obtained with `MatGetFactor()`
7046:   Uses levels of fill only, not drop tolerance. Use `MatLUFactorNumeric()`
7047:   to complete the factorization.

7049:   Collective

7051:   Input Parameters:
7052: + fact - the factorized matrix obtained with `MatGetFactor()`
7053: . mat  - the matrix
7054: . row  - row permutation
7055: . col  - column permutation
7056: - info - structure containing
7057: .vb
7058:       levels - number of levels of fill.
7059:       expected fill - as ratio of original fill.
7060:       1 or 0 - indicating force fill on diagonal (improves robustness for matrices
7061:                 missing diagonal entries)
7062: .ve

7064:   Level: developer

7066:   Notes:
7067:   See [Matrix Factorization](sec_matfactor) for additional information.

7069:   Most users should employ the `KSP` interface for linear solvers
7070:   instead of working directly with matrix algebra routines such as this.
7071:   See, e.g., `KSPCreate()`.

7073:   Uses the definition of level of fill as in Y. Saad, {cite}`saad2003`

7075:   Developer Note:
7076:   The Fortran interface is not autogenerated as the
7077:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

7079: .seealso: [](ch_matrices), `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`
7080:           `MatGetOrdering()`, `MatFactorInfo`
7081: @*/
7082: PetscErrorCode MatILUFactorSymbolic(Mat fact, Mat mat, IS row, IS col, const MatFactorInfo *info)
7083: {
7084:   PetscFunctionBegin;
7089:   PetscAssertPointer(info, 5);
7090:   PetscAssertPointer(fact, 1);
7091:   PetscCheck(info->levels >= 0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Levels of fill negative %" PetscInt_FMT, (PetscInt)info->levels);
7092:   PetscCheck(info->fill >= 1.0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Expected fill less than 1.0 %g", (double)info->fill);
7093:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
7094:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
7095:   MatCheckPreallocated(mat, 2);

7097:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_ILUFactorSymbolic, mat, row, col, 0));
7098:   PetscUseTypeMethod(fact, ilufactorsymbolic, mat, row, col, info);
7099:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_ILUFactorSymbolic, mat, row, col, 0));
7100:   PetscFunctionReturn(PETSC_SUCCESS);
7101: }

7103: /*@
7104:   MatICCFactorSymbolic - Performs symbolic incomplete
7105:   Cholesky factorization for a symmetric matrix.  Use
7106:   `MatCholeskyFactorNumeric()` to complete the factorization.

7108:   Collective

7110:   Input Parameters:
7111: + fact - the factorized matrix obtained with `MatGetFactor()`
7112: . mat  - the matrix to be factored
7113: . perm - row and column permutation
7114: - info - structure containing
7115: .vb
7116:       levels - number of levels of fill.
7117:       expected fill - as ratio of original fill.
7118: .ve

7120:   Level: developer

7122:   Notes:
7123:   Most users should employ the `KSP` interface for linear solvers
7124:   instead of working directly with matrix algebra routines such as this.
7125:   See, e.g., `KSPCreate()`.

7127:   This uses the definition of level of fill as in Y. Saad {cite}`saad2003`

7129:   Developer Note:
7130:   The Fortran interface is not autogenerated as the
7131:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

7133: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatCholeskyFactorNumeric()`, `MatCholeskyFactor()`, `MatFactorInfo`
7134: @*/
7135: PetscErrorCode MatICCFactorSymbolic(Mat fact, Mat mat, IS perm, const MatFactorInfo *info)
7136: {
7137:   PetscFunctionBegin;
7141:   PetscAssertPointer(info, 4);
7142:   PetscAssertPointer(fact, 1);
7143:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
7144:   PetscCheck(info->levels >= 0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Levels negative %" PetscInt_FMT, (PetscInt)info->levels);
7145:   PetscCheck(info->fill >= 1.0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Expected fill less than 1.0 %g", (double)info->fill);
7146:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
7147:   MatCheckPreallocated(mat, 2);

7149:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_ICCFactorSymbolic, mat, perm, 0, 0));
7150:   PetscUseTypeMethod(fact, iccfactorsymbolic, mat, perm, info);
7151:   if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_ICCFactorSymbolic, mat, perm, 0, 0));
7152:   PetscFunctionReturn(PETSC_SUCCESS);
7153: }

7155: /*@C
7156:   MatCreateSubMatrices - Extracts several submatrices from a matrix. If submat
7157:   points to an array of valid matrices, they may be reused to store the new
7158:   submatrices.

7160:   Collective

7162:   Input Parameters:
7163: + mat   - the matrix
7164: . n     - the number of submatrixes to be extracted (on this processor, may be zero)
7165: . irow  - index set of rows to extract
7166: . icol  - index set of columns to extract
7167: - scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`

7169:   Output Parameter:
7170: . submat - the array of submatrices

7172:   Level: advanced

7174:   Notes:
7175:   `MatCreateSubMatrices()` can extract ONLY sequential submatrices
7176:   (from both sequential and parallel matrices). Use `MatCreateSubMatrix()`
7177:   to extract a parallel submatrix.

7179:   Some matrix types place restrictions on the row and column
7180:   indices, such as that they be sorted or that they be equal to each other.

7182:   The index sets may not have duplicate entries.

7184:   When extracting submatrices from a parallel matrix, each processor can
7185:   form a different submatrix by setting the rows and columns of its
7186:   individual index sets according to the local submatrix desired.

7188:   When finished using the submatrices, the user should destroy
7189:   them with `MatDestroySubMatrices()`.

7191:   `MAT_REUSE_MATRIX` can only be used when the nonzero structure of the
7192:   original matrix has not changed from that last call to `MatCreateSubMatrices()`.

7194:   This routine creates the matrices in submat; you should NOT create them before
7195:   calling it. It also allocates the array of matrix pointers submat.

7197:   For `MATBAIJ` matrices the index sets must respect the block structure, that is if they
7198:   request one row/column in a block, they must request all rows/columns that are in
7199:   that block. For example, if the block size is 2 you cannot request just row 0 and
7200:   column 0.

7202:   Fortran Note:
7203:   One must pass in as `submat` a `Mat` array of size at least `n`+1.

7205: .seealso: [](ch_matrices), `Mat`, `MatDestroySubMatrices()`, `MatCreateSubMatrix()`, `MatGetRow()`, `MatGetDiagonal()`, `MatReuse`
7206: @*/
7207: PetscErrorCode MatCreateSubMatrices(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
7208: {
7209:   PetscInt  i;
7210:   PetscBool eq;

7212:   PetscFunctionBegin;
7215:   if (n) {
7216:     PetscAssertPointer(irow, 3);
7218:     PetscAssertPointer(icol, 4);
7220:   }
7221:   PetscAssertPointer(submat, 6);
7222:   if (n && scall == MAT_REUSE_MATRIX) {
7223:     PetscAssertPointer(*submat, 6);
7225:   }
7226:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
7227:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
7228:   MatCheckPreallocated(mat, 1);
7229:   PetscCall(PetscLogEventBegin(MAT_CreateSubMats, mat, 0, 0, 0));
7230:   PetscUseTypeMethod(mat, createsubmatrices, n, irow, icol, scall, submat);
7231:   PetscCall(PetscLogEventEnd(MAT_CreateSubMats, mat, 0, 0, 0));
7232:   for (i = 0; i < n; i++) {
7233:     (*submat)[i]->factortype = MAT_FACTOR_NONE; /* in case in place factorization was previously done on submatrix */
7234:     PetscCall(ISEqualUnsorted(irow[i], icol[i], &eq));
7235:     if (eq) PetscCall(MatPropagateSymmetryOptions(mat, (*submat)[i]));
7236: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
7237:     if (mat->boundtocpu && mat->bindingpropagates) {
7238:       PetscCall(MatBindToCPU((*submat)[i], PETSC_TRUE));
7239:       PetscCall(MatSetBindingPropagates((*submat)[i], PETSC_TRUE));
7240:     }
7241: #endif
7242:   }
7243:   PetscFunctionReturn(PETSC_SUCCESS);
7244: }

7246: /*@C
7247:   MatCreateSubMatricesMPI - Extracts MPI submatrices across a sub communicator of mat (by pairs of `IS` that may live on subcomms).

7249:   Collective

7251:   Input Parameters:
7252: + mat   - the matrix
7253: . n     - the number of submatrixes to be extracted
7254: . irow  - index set of rows to extract
7255: . icol  - index set of columns to extract
7256: - scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`

7258:   Output Parameter:
7259: . submat - the array of submatrices

7261:   Level: advanced

7263:   Note:
7264:   This is used by `PCGASM`

7266: .seealso: [](ch_matrices), `Mat`, `PCGASM`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRow()`, `MatGetDiagonal()`, `MatReuse`
7267: @*/
7268: PetscErrorCode MatCreateSubMatricesMPI(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
7269: {
7270:   PetscInt  i;
7271:   PetscBool eq;

7273:   PetscFunctionBegin;
7276:   if (n) {
7277:     PetscAssertPointer(irow, 3);
7279:     PetscAssertPointer(icol, 4);
7281:   }
7282:   PetscAssertPointer(submat, 6);
7283:   if (n && scall == MAT_REUSE_MATRIX) {
7284:     PetscAssertPointer(*submat, 6);
7286:   }
7287:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
7288:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
7289:   MatCheckPreallocated(mat, 1);

7291:   PetscCall(PetscLogEventBegin(MAT_CreateSubMats, mat, 0, 0, 0));
7292:   PetscUseTypeMethod(mat, createsubmatricesmpi, n, irow, icol, scall, submat);
7293:   PetscCall(PetscLogEventEnd(MAT_CreateSubMats, mat, 0, 0, 0));
7294:   for (i = 0; i < n; i++) {
7295:     PetscCall(ISEqualUnsorted(irow[i], icol[i], &eq));
7296:     if (eq) PetscCall(MatPropagateSymmetryOptions(mat, (*submat)[i]));
7297:   }
7298:   PetscFunctionReturn(PETSC_SUCCESS);
7299: }

7301: /*@C
7302:   MatDestroyMatrices - Destroys an array of matrices.

7304:   Collective

7306:   Input Parameters:
7307: + n   - the number of local matrices
7308: - mat - the matrices (this is a pointer to the array of matrices)

7310:   Level: advanced

7312:   Notes:
7313:   Frees not only the matrices, but also the array that contains the matrices

7315:   For matrices obtained with  `MatCreateSubMatrices()` use `MatDestroySubMatrices()`

7317:   Fortran Note:
7318:   Does not free the `mat` array.

7320: .seealso: [](ch_matrices), `Mat`, `MatCreateSubMatrices()`, `MatDestroySubMatrices()`
7321: @*/
7322: PetscErrorCode MatDestroyMatrices(PetscInt n, Mat *mat[])
7323: {
7324:   PetscInt i;

7326:   PetscFunctionBegin;
7327:   if (!*mat) PetscFunctionReturn(PETSC_SUCCESS);
7328:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Trying to destroy negative number of matrices %" PetscInt_FMT, n);
7329:   PetscAssertPointer(mat, 2);

7331:   for (i = 0; i < n; i++) PetscCall(MatDestroy(&(*mat)[i]));

7333:   /* memory is allocated even if n = 0 */
7334:   PetscCall(PetscFree(*mat));
7335:   PetscFunctionReturn(PETSC_SUCCESS);
7336: }

7338: /*@C
7339:   MatDestroySubMatrices - Destroys a set of matrices obtained with `MatCreateSubMatrices()`.

7341:   Collective

7343:   Input Parameters:
7344: + n   - the number of local matrices
7345: - mat - the matrices (this is a pointer to the array of matrices, just to match the calling
7346:                        sequence of `MatCreateSubMatrices()`)

7348:   Level: advanced

7350:   Note:
7351:   Frees not only the matrices, but also the array that contains the matrices

7353:   Fortran Note:
7354:   Does not free the `mat` array.

7356: .seealso: [](ch_matrices), `Mat`, `MatCreateSubMatrices()`, `MatDestroyMatrices()`
7357: @*/
7358: PetscErrorCode MatDestroySubMatrices(PetscInt n, Mat *mat[])
7359: {
7360:   Mat mat0;

7362:   PetscFunctionBegin;
7363:   if (!*mat) PetscFunctionReturn(PETSC_SUCCESS);
7364:   /* mat[] is an array of length n+1, see MatCreateSubMatrices_xxx() */
7365:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Trying to destroy negative number of matrices %" PetscInt_FMT, n);
7366:   PetscAssertPointer(mat, 2);

7368:   mat0 = (*mat)[0];
7369:   if (mat0 && mat0->ops->destroysubmatrices) {
7370:     PetscCall((*mat0->ops->destroysubmatrices)(n, mat));
7371:   } else {
7372:     PetscCall(MatDestroyMatrices(n, mat));
7373:   }
7374:   PetscFunctionReturn(PETSC_SUCCESS);
7375: }

7377: /*@
7378:   MatGetSeqNonzeroStructure - Extracts the nonzero structure from a matrix and stores it, in its entirety, on each process

7380:   Collective

7382:   Input Parameter:
7383: . mat - the matrix

7385:   Output Parameter:
7386: . matstruct - the sequential matrix with the nonzero structure of `mat`

7388:   Level: developer

7390: .seealso: [](ch_matrices), `Mat`, `MatDestroySeqNonzeroStructure()`, `MatCreateSubMatrices()`, `MatDestroyMatrices()`
7391: @*/
7392: PetscErrorCode MatGetSeqNonzeroStructure(Mat mat, Mat *matstruct)
7393: {
7394:   PetscFunctionBegin;
7396:   PetscAssertPointer(matstruct, 2);

7399:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
7400:   MatCheckPreallocated(mat, 1);

7402:   PetscCall(PetscLogEventBegin(MAT_GetSeqNonzeroStructure, mat, 0, 0, 0));
7403:   PetscUseTypeMethod(mat, getseqnonzerostructure, matstruct);
7404:   PetscCall(PetscLogEventEnd(MAT_GetSeqNonzeroStructure, mat, 0, 0, 0));
7405:   PetscFunctionReturn(PETSC_SUCCESS);
7406: }

7408: /*@C
7409:   MatDestroySeqNonzeroStructure - Destroys matrix obtained with `MatGetSeqNonzeroStructure()`.

7411:   Collective

7413:   Input Parameter:
7414: . mat - the matrix

7416:   Level: advanced

7418:   Note:
7419:   This is not needed, one can just call `MatDestroy()`

7421: .seealso: [](ch_matrices), `Mat`, `MatGetSeqNonzeroStructure()`
7422: @*/
7423: PetscErrorCode MatDestroySeqNonzeroStructure(Mat *mat)
7424: {
7425:   PetscFunctionBegin;
7426:   PetscAssertPointer(mat, 1);
7427:   PetscCall(MatDestroy(mat));
7428:   PetscFunctionReturn(PETSC_SUCCESS);
7429: }

7431: /*@
7432:   MatIncreaseOverlap - Given a set of submatrices indicated by index sets,
7433:   replaces the index sets by larger ones that represent submatrices with
7434:   additional overlap.

7436:   Collective

7438:   Input Parameters:
7439: + mat - the matrix
7440: . n   - the number of index sets
7441: . is  - the array of index sets (these index sets will changed during the call)
7442: - ov  - the additional overlap requested

7444:   Options Database Key:
7445: . -mat_increase_overlap_scalable - use a scalable algorithm to compute the overlap (supported by MPIAIJ matrix)

7447:   Level: developer

7449:   Note:
7450:   The computed overlap preserves the matrix block sizes when the blocks are square.
7451:   That is: if a matrix nonzero for a given block would increase the overlap all columns associated with
7452:   that block are included in the overlap regardless of whether each specific column would increase the overlap.

7454: .seealso: [](ch_matrices), `Mat`, `PCASM`, `MatSetBlockSize()`, `MatIncreaseOverlapSplit()`, `MatCreateSubMatrices()`
7455: @*/
7456: PetscErrorCode MatIncreaseOverlap(Mat mat, PetscInt n, IS is[], PetscInt ov)
7457: {
7458:   PetscInt i, bs, cbs;

7460:   PetscFunctionBegin;
7464:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have one or more domains, you have %" PetscInt_FMT, n);
7465:   if (n) {
7466:     PetscAssertPointer(is, 3);
7468:   }
7469:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
7470:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
7471:   MatCheckPreallocated(mat, 1);

7473:   if (!ov || !n) PetscFunctionReturn(PETSC_SUCCESS);
7474:   PetscCall(PetscLogEventBegin(MAT_IncreaseOverlap, mat, 0, 0, 0));
7475:   PetscUseTypeMethod(mat, increaseoverlap, n, is, ov);
7476:   PetscCall(PetscLogEventEnd(MAT_IncreaseOverlap, mat, 0, 0, 0));
7477:   PetscCall(MatGetBlockSizes(mat, &bs, &cbs));
7478:   if (bs == cbs) {
7479:     for (i = 0; i < n; i++) PetscCall(ISSetBlockSize(is[i], bs));
7480:   }
7481:   PetscFunctionReturn(PETSC_SUCCESS);
7482: }

7484: PetscErrorCode MatIncreaseOverlapSplit_Single(Mat, IS *, PetscInt);

7486: /*@
7487:   MatIncreaseOverlapSplit - Given a set of submatrices indicated by index sets across
7488:   a sub communicator, replaces the index sets by larger ones that represent submatrices with
7489:   additional overlap.

7491:   Collective

7493:   Input Parameters:
7494: + mat - the matrix
7495: . n   - the number of index sets
7496: . is  - the array of index sets (these index sets will changed during the call)
7497: - ov  - the additional overlap requested

7499:   `   Options Database Key:
7500: . -mat_increase_overlap_scalable - use a scalable algorithm to compute the overlap (supported by MPIAIJ matrix)

7502:   Level: developer

7504: .seealso: [](ch_matrices), `Mat`, `MatCreateSubMatrices()`, `MatIncreaseOverlap()`
7505: @*/
7506: PetscErrorCode MatIncreaseOverlapSplit(Mat mat, PetscInt n, IS is[], PetscInt ov)
7507: {
7508:   PetscInt i;

7510:   PetscFunctionBegin;
7513:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have one or more domains, you have %" PetscInt_FMT, n);
7514:   if (n) {
7515:     PetscAssertPointer(is, 3);
7517:   }
7518:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
7519:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
7520:   MatCheckPreallocated(mat, 1);
7521:   if (!ov) PetscFunctionReturn(PETSC_SUCCESS);
7522:   PetscCall(PetscLogEventBegin(MAT_IncreaseOverlap, mat, 0, 0, 0));
7523:   for (i = 0; i < n; i++) PetscCall(MatIncreaseOverlapSplit_Single(mat, &is[i], ov));
7524:   PetscCall(PetscLogEventEnd(MAT_IncreaseOverlap, mat, 0, 0, 0));
7525:   PetscFunctionReturn(PETSC_SUCCESS);
7526: }

7528: /*@
7529:   MatGetBlockSize - Returns the matrix block size.

7531:   Not Collective

7533:   Input Parameter:
7534: . mat - the matrix

7536:   Output Parameter:
7537: . bs - block size

7539:   Level: intermediate

7541:   Notes:
7542:   Block row formats are `MATBAIJ` and `MATSBAIJ` ALWAYS have square block storage in the matrix.

7544:   If the block size has not been set yet this routine returns 1.

7546: .seealso: [](ch_matrices), `Mat`, `MATBAIJ`, `MATSBAIJ`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSizes()`
7547: @*/
7548: PetscErrorCode MatGetBlockSize(Mat mat, PetscInt *bs)
7549: {
7550:   PetscFunctionBegin;
7552:   PetscAssertPointer(bs, 2);
7553:   *bs = PetscAbs(mat->rmap->bs);
7554:   PetscFunctionReturn(PETSC_SUCCESS);
7555: }

7557: /*@
7558:   MatGetBlockSizes - Returns the matrix block row and column sizes.

7560:   Not Collective

7562:   Input Parameter:
7563: . mat - the matrix

7565:   Output Parameters:
7566: + rbs - row block size
7567: - cbs - column block size

7569:   Level: intermediate

7571:   Notes:
7572:   Block row formats are `MATBAIJ` and `MATSBAIJ` ALWAYS have square block storage in the matrix.
7573:   If you pass a different block size for the columns than the rows, the row block size determines the square block storage.

7575:   If a block size has not been set yet this routine returns 1.

7577: .seealso: [](ch_matrices), `Mat`, `MATBAIJ`, `MATSBAIJ`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSize()`, `MatSetBlockSizes()`
7578: @*/
7579: PetscErrorCode MatGetBlockSizes(Mat mat, PetscInt *rbs, PetscInt *cbs)
7580: {
7581:   PetscFunctionBegin;
7583:   if (rbs) PetscAssertPointer(rbs, 2);
7584:   if (cbs) PetscAssertPointer(cbs, 3);
7585:   if (rbs) *rbs = PetscAbs(mat->rmap->bs);
7586:   if (cbs) *cbs = PetscAbs(mat->cmap->bs);
7587:   PetscFunctionReturn(PETSC_SUCCESS);
7588: }

7590: /*@
7591:   MatSetBlockSize - Sets the matrix block size.

7593:   Logically Collective

7595:   Input Parameters:
7596: + mat - the matrix
7597: - bs  - block size

7599:   Level: intermediate

7601:   Notes:
7602:   Block row formats are `MATBAIJ` and `MATSBAIJ` formats ALWAYS have square block storage in the matrix.
7603:   This must be called before `MatSetUp()` or MatXXXSetPreallocation() (or will default to 1) and the block size cannot be changed later.

7605:   For `MATAIJ` matrix format, this function can be called at a later stage, provided that the specified block size
7606:   is compatible with the matrix local sizes.

7608: .seealso: [](ch_matrices), `Mat`, `MATBAIJ`, `MATSBAIJ`, `MATAIJ`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSizes()`, `MatGetBlockSizes()`
7609: @*/
7610: PetscErrorCode MatSetBlockSize(Mat mat, PetscInt bs)
7611: {
7612:   PetscFunctionBegin;
7615:   PetscCall(MatSetBlockSizes(mat, bs, bs));
7616:   PetscFunctionReturn(PETSC_SUCCESS);
7617: }

7619: typedef struct {
7620:   PetscInt         n;
7621:   IS              *is;
7622:   Mat             *mat;
7623:   PetscObjectState nonzerostate;
7624:   Mat              C;
7625: } EnvelopeData;

7627: static PetscErrorCode EnvelopeDataDestroy(void *ptr)
7628: {
7629:   EnvelopeData *edata = (EnvelopeData *)ptr;

7631:   PetscFunctionBegin;
7632:   for (PetscInt i = 0; i < edata->n; i++) PetscCall(ISDestroy(&edata->is[i]));
7633:   PetscCall(PetscFree(edata->is));
7634:   PetscCall(PetscFree(edata));
7635:   PetscFunctionReturn(PETSC_SUCCESS);
7636: }

7638: /*@
7639:   MatComputeVariableBlockEnvelope - Given a matrix whose nonzeros are in blocks along the diagonal this computes and stores
7640:   the sizes of these blocks in the matrix. An individual block may lie over several processes.

7642:   Collective

7644:   Input Parameter:
7645: . mat - the matrix

7647:   Level: intermediate

7649:   Notes:
7650:   There can be zeros within the blocks

7652:   The blocks can overlap between processes, including laying on more than two processes

7654: .seealso: [](ch_matrices), `Mat`, `MatInvertVariableBlockEnvelope()`, `MatSetVariableBlockSizes()`
7655: @*/
7656: PetscErrorCode MatComputeVariableBlockEnvelope(Mat mat)
7657: {
7658:   PetscInt           n, *sizes, *starts, i = 0, env = 0, tbs = 0, lblocks = 0, rstart, II, ln = 0, cnt = 0, cstart, cend;
7659:   PetscInt          *diag, *odiag, sc;
7660:   VecScatter         scatter;
7661:   PetscScalar       *seqv;
7662:   const PetscScalar *parv;
7663:   const PetscInt    *ia, *ja;
7664:   PetscBool          set, flag, done;
7665:   Mat                AA = mat, A;
7666:   MPI_Comm           comm;
7667:   PetscMPIInt        rank, size, tag;
7668:   MPI_Status         status;
7669:   PetscContainer     container;
7670:   EnvelopeData      *edata;
7671:   Vec                seq, par;
7672:   IS                 isglobal;

7674:   PetscFunctionBegin;
7676:   PetscCall(MatIsSymmetricKnown(mat, &set, &flag));
7677:   if (!set || !flag) {
7678:     /* TODO: only needs nonzero structure of transpose */
7679:     PetscCall(MatTranspose(mat, MAT_INITIAL_MATRIX, &AA));
7680:     PetscCall(MatAXPY(AA, 1.0, mat, DIFFERENT_NONZERO_PATTERN));
7681:   }
7682:   PetscCall(MatAIJGetLocalMat(AA, &A));
7683:   PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &n, &ia, &ja, &done));
7684:   PetscCheck(done, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unable to get IJ structure from matrix");

7686:   PetscCall(MatGetLocalSize(mat, &n, NULL));
7687:   PetscCall(PetscObjectGetNewTag((PetscObject)mat, &tag));
7688:   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
7689:   PetscCallMPI(MPI_Comm_size(comm, &size));
7690:   PetscCallMPI(MPI_Comm_rank(comm, &rank));

7692:   PetscCall(PetscMalloc2(n, &sizes, n, &starts));

7694:   if (rank > 0) {
7695:     PetscCallMPI(MPI_Recv(&env, 1, MPIU_INT, rank - 1, tag, comm, &status));
7696:     PetscCallMPI(MPI_Recv(&tbs, 1, MPIU_INT, rank - 1, tag, comm, &status));
7697:   }
7698:   PetscCall(MatGetOwnershipRange(mat, &rstart, NULL));
7699:   for (i = 0; i < n; i++) {
7700:     env = PetscMax(env, ja[ia[i + 1] - 1]);
7701:     II  = rstart + i;
7702:     if (env == II) {
7703:       starts[lblocks]  = tbs;
7704:       sizes[lblocks++] = 1 + II - tbs;
7705:       tbs              = 1 + II;
7706:     }
7707:   }
7708:   if (rank < size - 1) {
7709:     PetscCallMPI(MPI_Send(&env, 1, MPIU_INT, rank + 1, tag, comm));
7710:     PetscCallMPI(MPI_Send(&tbs, 1, MPIU_INT, rank + 1, tag, comm));
7711:   }

7713:   PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &n, &ia, &ja, &done));
7714:   if (!set || !flag) PetscCall(MatDestroy(&AA));
7715:   PetscCall(MatDestroy(&A));

7717:   PetscCall(PetscNew(&edata));
7718:   PetscCall(MatGetNonzeroState(mat, &edata->nonzerostate));
7719:   edata->n = lblocks;
7720:   /* create IS needed for extracting blocks from the original matrix */
7721:   PetscCall(PetscMalloc1(lblocks, &edata->is));
7722:   for (PetscInt i = 0; i < lblocks; i++) PetscCall(ISCreateStride(PETSC_COMM_SELF, sizes[i], starts[i], 1, &edata->is[i]));

7724:   /* Create the resulting inverse matrix nonzero structure with preallocation information */
7725:   PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &edata->C));
7726:   PetscCall(MatSetSizes(edata->C, mat->rmap->n, mat->cmap->n, mat->rmap->N, mat->cmap->N));
7727:   PetscCall(MatSetBlockSizesFromMats(edata->C, mat, mat));
7728:   PetscCall(MatSetType(edata->C, MATAIJ));

7730:   /* Communicate the start and end of each row, from each block to the correct rank */
7731:   /* TODO: Use PetscSF instead of VecScatter */
7732:   for (PetscInt i = 0; i < lblocks; i++) ln += sizes[i];
7733:   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 2 * ln, &seq));
7734:   PetscCall(VecGetArrayWrite(seq, &seqv));
7735:   for (PetscInt i = 0; i < lblocks; i++) {
7736:     for (PetscInt j = 0; j < sizes[i]; j++) {
7737:       seqv[cnt]     = starts[i];
7738:       seqv[cnt + 1] = starts[i] + sizes[i];
7739:       cnt += 2;
7740:     }
7741:   }
7742:   PetscCall(VecRestoreArrayWrite(seq, &seqv));
7743:   PetscCallMPI(MPI_Scan(&cnt, &sc, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)mat)));
7744:   sc -= cnt;
7745:   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)mat), 2 * mat->rmap->n, 2 * mat->rmap->N, &par));
7746:   PetscCall(ISCreateStride(PETSC_COMM_SELF, cnt, sc, 1, &isglobal));
7747:   PetscCall(VecScatterCreate(seq, NULL, par, isglobal, &scatter));
7748:   PetscCall(ISDestroy(&isglobal));
7749:   PetscCall(VecScatterBegin(scatter, seq, par, INSERT_VALUES, SCATTER_FORWARD));
7750:   PetscCall(VecScatterEnd(scatter, seq, par, INSERT_VALUES, SCATTER_FORWARD));
7751:   PetscCall(VecScatterDestroy(&scatter));
7752:   PetscCall(VecDestroy(&seq));
7753:   PetscCall(MatGetOwnershipRangeColumn(mat, &cstart, &cend));
7754:   PetscCall(PetscMalloc2(mat->rmap->n, &diag, mat->rmap->n, &odiag));
7755:   PetscCall(VecGetArrayRead(par, &parv));
7756:   cnt = 0;
7757:   PetscCall(MatGetSize(mat, NULL, &n));
7758:   for (PetscInt i = 0; i < mat->rmap->n; i++) {
7759:     PetscInt start, end, d = 0, od = 0;

7761:     start = (PetscInt)PetscRealPart(parv[cnt]);
7762:     end   = (PetscInt)PetscRealPart(parv[cnt + 1]);
7763:     cnt += 2;

7765:     if (start < cstart) {
7766:       od += cstart - start + n - cend;
7767:       d += cend - cstart;
7768:     } else if (start < cend) {
7769:       od += n - cend;
7770:       d += cend - start;
7771:     } else od += n - start;
7772:     if (end <= cstart) {
7773:       od -= cstart - end + n - cend;
7774:       d -= cend - cstart;
7775:     } else if (end < cend) {
7776:       od -= n - cend;
7777:       d -= cend - end;
7778:     } else od -= n - end;

7780:     odiag[i] = od;
7781:     diag[i]  = d;
7782:   }
7783:   PetscCall(VecRestoreArrayRead(par, &parv));
7784:   PetscCall(VecDestroy(&par));
7785:   PetscCall(MatXAIJSetPreallocation(edata->C, mat->rmap->bs, diag, odiag, NULL, NULL));
7786:   PetscCall(PetscFree2(diag, odiag));
7787:   PetscCall(PetscFree2(sizes, starts));

7789:   PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
7790:   PetscCall(PetscContainerSetPointer(container, edata));
7791:   PetscCall(PetscContainerSetUserDestroy(container, (PetscErrorCode (*)(void *))EnvelopeDataDestroy));
7792:   PetscCall(PetscObjectCompose((PetscObject)mat, "EnvelopeData", (PetscObject)container));
7793:   PetscCall(PetscObjectDereference((PetscObject)container));
7794:   PetscFunctionReturn(PETSC_SUCCESS);
7795: }

7797: /*@
7798:   MatInvertVariableBlockEnvelope - set matrix C to be the inverted block diagonal of matrix A

7800:   Collective

7802:   Input Parameters:
7803: + A     - the matrix
7804: - reuse - indicates if the `C` matrix was obtained from a previous call to this routine

7806:   Output Parameter:
7807: . C - matrix with inverted block diagonal of `A`

7809:   Level: advanced

7811:   Note:
7812:   For efficiency the matrix `A` should have all the nonzero entries clustered in smallish blocks along the diagonal.

7814: .seealso: [](ch_matrices), `Mat`, `MatInvertBlockDiagonal()`, `MatComputeBlockDiagonal()`
7815: @*/
7816: PetscErrorCode MatInvertVariableBlockEnvelope(Mat A, MatReuse reuse, Mat *C)
7817: {
7818:   PetscContainer   container;
7819:   EnvelopeData    *edata;
7820:   PetscObjectState nonzerostate;

7822:   PetscFunctionBegin;
7823:   PetscCall(PetscObjectQuery((PetscObject)A, "EnvelopeData", (PetscObject *)&container));
7824:   if (!container) {
7825:     PetscCall(MatComputeVariableBlockEnvelope(A));
7826:     PetscCall(PetscObjectQuery((PetscObject)A, "EnvelopeData", (PetscObject *)&container));
7827:   }
7828:   PetscCall(PetscContainerGetPointer(container, (void **)&edata));
7829:   PetscCall(MatGetNonzeroState(A, &nonzerostate));
7830:   PetscCheck(nonzerostate <= edata->nonzerostate, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot handle changes to matrix nonzero structure");
7831:   PetscCheck(reuse != MAT_REUSE_MATRIX || *C == edata->C, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "C matrix must be the same as previously output");

7833:   PetscCall(MatCreateSubMatrices(A, edata->n, edata->is, edata->is, MAT_INITIAL_MATRIX, &edata->mat));
7834:   *C = edata->C;

7836:   for (PetscInt i = 0; i < edata->n; i++) {
7837:     Mat          D;
7838:     PetscScalar *dvalues;

7840:     PetscCall(MatConvert(edata->mat[i], MATSEQDENSE, MAT_INITIAL_MATRIX, &D));
7841:     PetscCall(MatSetOption(*C, MAT_ROW_ORIENTED, PETSC_FALSE));
7842:     PetscCall(MatSeqDenseInvert(D));
7843:     PetscCall(MatDenseGetArray(D, &dvalues));
7844:     PetscCall(MatSetValuesIS(*C, edata->is[i], edata->is[i], dvalues, INSERT_VALUES));
7845:     PetscCall(MatDestroy(&D));
7846:   }
7847:   PetscCall(MatDestroySubMatrices(edata->n, &edata->mat));
7848:   PetscCall(MatAssemblyBegin(*C, MAT_FINAL_ASSEMBLY));
7849:   PetscCall(MatAssemblyEnd(*C, MAT_FINAL_ASSEMBLY));
7850:   PetscFunctionReturn(PETSC_SUCCESS);
7851: }

7853: /*@
7854:   MatSetVariableBlockSizes - Sets diagonal point-blocks of the matrix that need not be of the same size

7856:   Not Collective

7858:   Input Parameters:
7859: + mat     - the matrix
7860: . nblocks - the number of blocks on this process, each block can only exist on a single process
7861: - bsizes  - the block sizes

7863:   Level: intermediate

7865:   Notes:
7866:   Currently used by `PCVPBJACOBI` for `MATAIJ` matrices

7868:   Each variable point-block set of degrees of freedom must live on a single MPI process. That is a point block cannot straddle two MPI processes.

7870: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSizes()`, `MatGetBlockSizes()`, `MatGetVariableBlockSizes()`,
7871:           `MatComputeVariableBlockEnvelope()`, `PCVPBJACOBI`
7872: @*/
7873: PetscErrorCode MatSetVariableBlockSizes(Mat mat, PetscInt nblocks, const PetscInt bsizes[])
7874: {
7875:   PetscInt ncnt = 0, nlocal;

7877:   PetscFunctionBegin;
7879:   PetscCall(MatGetLocalSize(mat, &nlocal, NULL));
7880:   PetscCheck(nblocks >= 0 && nblocks <= nlocal, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of local blocks %" PetscInt_FMT " is not in [0, %" PetscInt_FMT "]", nblocks, nlocal);
7881:   for (PetscInt i = 0; i < nblocks; i++) ncnt += bsizes[i];
7882:   PetscCheck(ncnt == nlocal, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Sum of local block sizes %" PetscInt_FMT " does not equal local size of matrix %" PetscInt_FMT, ncnt, nlocal);
7883:   PetscCall(PetscFree(mat->bsizes));
7884:   mat->nblocks = nblocks;
7885:   PetscCall(PetscMalloc1(nblocks, &mat->bsizes));
7886:   PetscCall(PetscArraycpy(mat->bsizes, bsizes, nblocks));
7887:   PetscFunctionReturn(PETSC_SUCCESS);
7888: }

7890: /*@C
7891:   MatGetVariableBlockSizes - Gets a diagonal blocks of the matrix that need not be of the same size

7893:   Not Collective; No Fortran Support

7895:   Input Parameter:
7896: . mat - the matrix

7898:   Output Parameters:
7899: + nblocks - the number of blocks on this process
7900: - bsizes  - the block sizes

7902:   Level: intermediate

7904: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSizes()`, `MatGetBlockSizes()`, `MatSetVariableBlockSizes()`, `MatComputeVariableBlockEnvelope()`
7905: @*/
7906: PetscErrorCode MatGetVariableBlockSizes(Mat mat, PetscInt *nblocks, const PetscInt *bsizes[])
7907: {
7908:   PetscFunctionBegin;
7910:   if (nblocks) *nblocks = mat->nblocks;
7911:   if (bsizes) *bsizes = mat->bsizes;
7912:   PetscFunctionReturn(PETSC_SUCCESS);
7913: }

7915: /*@
7916:   MatSetBlockSizes - Sets the matrix block row and column sizes.

7918:   Logically Collective

7920:   Input Parameters:
7921: + mat - the matrix
7922: . rbs - row block size
7923: - cbs - column block size

7925:   Level: intermediate

7927:   Notes:
7928:   Block row formats are `MATBAIJ` and  `MATSBAIJ`. These formats ALWAYS have square block storage in the matrix.
7929:   If you pass a different block size for the columns than the rows, the row block size determines the square block storage.
7930:   This must be called before `MatSetUp()` or MatXXXSetPreallocation() (or will default to 1) and the block size cannot be changed later.

7932:   For `MATAIJ` matrix this function can be called at a later stage, provided that the specified block sizes
7933:   are compatible with the matrix local sizes.

7935:   The row and column block size determine the blocksize of the "row" and "column" vectors returned by `MatCreateVecs()`.

7937: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSize()`, `MatGetBlockSizes()`
7938: @*/
7939: PetscErrorCode MatSetBlockSizes(Mat mat, PetscInt rbs, PetscInt cbs)
7940: {
7941:   PetscFunctionBegin;
7945:   PetscTryTypeMethod(mat, setblocksizes, rbs, cbs);
7946:   if (mat->rmap->refcnt) {
7947:     ISLocalToGlobalMapping l2g  = NULL;
7948:     PetscLayout            nmap = NULL;

7950:     PetscCall(PetscLayoutDuplicate(mat->rmap, &nmap));
7951:     if (mat->rmap->mapping) PetscCall(ISLocalToGlobalMappingDuplicate(mat->rmap->mapping, &l2g));
7952:     PetscCall(PetscLayoutDestroy(&mat->rmap));
7953:     mat->rmap          = nmap;
7954:     mat->rmap->mapping = l2g;
7955:   }
7956:   if (mat->cmap->refcnt) {
7957:     ISLocalToGlobalMapping l2g  = NULL;
7958:     PetscLayout            nmap = NULL;

7960:     PetscCall(PetscLayoutDuplicate(mat->cmap, &nmap));
7961:     if (mat->cmap->mapping) PetscCall(ISLocalToGlobalMappingDuplicate(mat->cmap->mapping, &l2g));
7962:     PetscCall(PetscLayoutDestroy(&mat->cmap));
7963:     mat->cmap          = nmap;
7964:     mat->cmap->mapping = l2g;
7965:   }
7966:   PetscCall(PetscLayoutSetBlockSize(mat->rmap, rbs));
7967:   PetscCall(PetscLayoutSetBlockSize(mat->cmap, cbs));
7968:   PetscFunctionReturn(PETSC_SUCCESS);
7969: }

7971: /*@
7972:   MatSetBlockSizesFromMats - Sets the matrix block row and column sizes to match a pair of matrices

7974:   Logically Collective

7976:   Input Parameters:
7977: + mat     - the matrix
7978: . fromRow - matrix from which to copy row block size
7979: - fromCol - matrix from which to copy column block size (can be same as fromRow)

7981:   Level: developer

7983: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSizes()`
7984: @*/
7985: PetscErrorCode MatSetBlockSizesFromMats(Mat mat, Mat fromRow, Mat fromCol)
7986: {
7987:   PetscFunctionBegin;
7991:   if (fromRow->rmap->bs > 0) PetscCall(PetscLayoutSetBlockSize(mat->rmap, fromRow->rmap->bs));
7992:   if (fromCol->cmap->bs > 0) PetscCall(PetscLayoutSetBlockSize(mat->cmap, fromCol->cmap->bs));
7993:   PetscFunctionReturn(PETSC_SUCCESS);
7994: }

7996: /*@
7997:   MatResidual - Default routine to calculate the residual r = b - Ax

7999:   Collective

8001:   Input Parameters:
8002: + mat - the matrix
8003: . b   - the right-hand-side
8004: - x   - the approximate solution

8006:   Output Parameter:
8007: . r - location to store the residual

8009:   Level: developer

8011: .seealso: [](ch_matrices), `Mat`, `MatMult()`, `MatMultAdd()`, `PCMGSetResidual()`
8012: @*/
8013: PetscErrorCode MatResidual(Mat mat, Vec b, Vec x, Vec r)
8014: {
8015:   PetscFunctionBegin;
8021:   MatCheckPreallocated(mat, 1);
8022:   PetscCall(PetscLogEventBegin(MAT_Residual, mat, 0, 0, 0));
8023:   if (!mat->ops->residual) {
8024:     PetscCall(MatMult(mat, x, r));
8025:     PetscCall(VecAYPX(r, -1.0, b));
8026:   } else {
8027:     PetscUseTypeMethod(mat, residual, b, x, r);
8028:   }
8029:   PetscCall(PetscLogEventEnd(MAT_Residual, mat, 0, 0, 0));
8030:   PetscFunctionReturn(PETSC_SUCCESS);
8031: }

8033: /*MC
8034:     MatGetRowIJF90 - Obtains the compressed row storage i and j indices for the local rows of a sparse matrix

8036:     Synopsis:
8037:     MatGetRowIJF90(Mat A, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt n, {PetscInt, pointer :: ia(:)}, {PetscInt, pointer :: ja(:)}, PetscBool done,integer ierr)

8039:     Not Collective

8041:     Input Parameters:
8042: +   A - the matrix
8043: .   shift -  0 or 1 indicating we want the indices starting at 0 or 1
8044: .   symmetric - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized
8045: -   inodecompressed - `PETSC_TRUE` or `PETSC_FALSE`  indicating if the nonzero structure of the
8046:                  inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is
8047:                  always used.

8049:     Output Parameters:
8050: +   n - number of local rows in the (possibly compressed) matrix
8051: .   ia - the row pointers; that is ia[0] = 0, ia[row] = ia[row-1] + number of elements in that row of the matrix
8052: .   ja - the column indices
8053: -   done - indicates if the routine actually worked and returned appropriate ia[] and ja[] arrays; callers
8054:            are responsible for handling the case when done == `PETSC_FALSE` and ia and ja are not set

8056:     Level: developer

8058:     Note:
8059:     Use  `MatRestoreRowIJF90()` when you no longer need access to the data

8061: .seealso: [](ch_matrices), [](sec_fortranarrays), `Mat`, `MATMPIAIJ`, `MatGetRowIJ()`, `MatRestoreRowIJ()`, `MatRestoreRowIJF90()`
8062: M*/

8064: /*MC
8065:     MatRestoreRowIJF90 - restores the compressed row storage i and j indices for the local rows of a sparse matrix obtained with `MatGetRowIJF90()`

8067:     Synopsis:
8068:     MatRestoreRowIJF90(Mat A, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt n, {PetscInt, pointer :: ia(:)}, {PetscInt, pointer :: ja(:)}, PetscBool done,integer ierr)

8070:     Not Collective

8072:     Input Parameters:
8073: +   A - the  matrix
8074: .   shift -  0 or 1 indicating we want the indices starting at 0 or 1
8075: .   symmetric - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized
8076:     inodecompressed - `PETSC_TRUE` or `PETSC_FALSE`  indicating if the nonzero structure of the
8077:                  inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is
8078:                  always used.
8079: .   n - number of local rows in the (possibly compressed) matrix
8080: .   ia - the row pointers; that is ia[0] = 0, ia[row] = ia[row-1] + number of elements in that row of the matrix
8081: .   ja - the column indices
8082: -   done - indicates if the routine actually worked and returned appropriate ia[] and ja[] arrays; callers
8083:            are responsible for handling the case when done == `PETSC_FALSE` and ia and ja are not set

8085:     Level: developer

8087: .seealso: [](ch_matrices), [](sec_fortranarrays), `Mat`, `MATMPIAIJ`, `MatGetRowIJ()`, `MatRestoreRowIJ()`, `MatGetRowIJF90()`
8088: M*/

8090: /*@C
8091:   MatGetRowIJ - Returns the compressed row storage i and j indices for the local rows of a sparse matrix

8093:   Collective

8095:   Input Parameters:
8096: + mat             - the matrix
8097: . shift           - 0 or 1 indicating we want the indices starting at 0 or 1
8098: . symmetric       - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized
8099: - inodecompressed - `PETSC_TRUE` or `PETSC_FALSE`  indicating if the nonzero structure of the
8100:                  inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is
8101:                  always used.

8103:   Output Parameters:
8104: + n    - number of local rows in the (possibly compressed) matrix, use `NULL` if not needed
8105: . ia   - the row pointers; that is ia[0] = 0, ia[row] = ia[row-1] + number of elements in that row of the matrix, use `NULL` if not needed
8106: . ja   - the column indices, use `NULL` if not needed
8107: - done - indicates if the routine actually worked and returned appropriate ia[] and ja[] arrays; callers
8108:            are responsible for handling the case when done == `PETSC_FALSE` and ia and ja are not set

8110:   Level: developer

8112:   Notes:
8113:   You CANNOT change any of the ia[] or ja[] values.

8115:   Use `MatRestoreRowIJ()` when you are finished accessing the ia[] and ja[] values.

8117:   Fortran Notes:
8118:   Use
8119: .vb
8120:     PetscInt, pointer :: ia(:),ja(:)
8121:     call MatGetRowIJF90(mat,shift,symmetric,inodecompressed,n,ia,ja,done,ierr)
8122:     ! Access the ith and jth entries via ia(i) and ja(j)
8123: .ve

8125:   `MatGetRowIJ()` Fortran binding is deprecated (since PETSc 3.19), use `MatGetRowIJF90()`

8127: .seealso: [](ch_matrices), `Mat`, `MATAIJ`, `MatGetRowIJF90()`, `MatGetColumnIJ()`, `MatRestoreRowIJ()`, `MatSeqAIJGetArray()`
8128: @*/
8129: PetscErrorCode MatGetRowIJ(Mat mat, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
8130: {
8131:   PetscFunctionBegin;
8134:   if (n) PetscAssertPointer(n, 5);
8135:   if (ia) PetscAssertPointer(ia, 6);
8136:   if (ja) PetscAssertPointer(ja, 7);
8137:   if (done) PetscAssertPointer(done, 8);
8138:   MatCheckPreallocated(mat, 1);
8139:   if (!mat->ops->getrowij && done) *done = PETSC_FALSE;
8140:   else {
8141:     if (done) *done = PETSC_TRUE;
8142:     PetscCall(PetscLogEventBegin(MAT_GetRowIJ, mat, 0, 0, 0));
8143:     PetscUseTypeMethod(mat, getrowij, shift, symmetric, inodecompressed, n, ia, ja, done);
8144:     PetscCall(PetscLogEventEnd(MAT_GetRowIJ, mat, 0, 0, 0));
8145:   }
8146:   PetscFunctionReturn(PETSC_SUCCESS);
8147: }

8149: /*@C
8150:   MatGetColumnIJ - Returns the compressed column storage i and j indices for sequential matrices.

8152:   Collective

8154:   Input Parameters:
8155: + mat             - the matrix
8156: . shift           - 1 or zero indicating we want the indices starting at 0 or 1
8157: . symmetric       - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be
8158:                 symmetrized
8159: . inodecompressed - `PETSC_TRUE` or `PETSC_FALSE` indicating if the nonzero structure of the
8160:                  inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is
8161:                  always used.
8162: . n               - number of columns in the (possibly compressed) matrix
8163: . ia              - the column pointers; that is ia[0] = 0, ia[col] = i[col-1] + number of elements in that col of the matrix
8164: - ja              - the row indices

8166:   Output Parameter:
8167: . done - `PETSC_TRUE` or `PETSC_FALSE`, indicating whether the values have been returned

8169:   Level: developer

8171: .seealso: [](ch_matrices), `Mat`, `MatGetRowIJ()`, `MatRestoreColumnIJ()`
8172: @*/
8173: PetscErrorCode MatGetColumnIJ(Mat mat, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
8174: {
8175:   PetscFunctionBegin;
8178:   PetscAssertPointer(n, 5);
8179:   if (ia) PetscAssertPointer(ia, 6);
8180:   if (ja) PetscAssertPointer(ja, 7);
8181:   PetscAssertPointer(done, 8);
8182:   MatCheckPreallocated(mat, 1);
8183:   if (!mat->ops->getcolumnij) *done = PETSC_FALSE;
8184:   else {
8185:     *done = PETSC_TRUE;
8186:     PetscUseTypeMethod(mat, getcolumnij, shift, symmetric, inodecompressed, n, ia, ja, done);
8187:   }
8188:   PetscFunctionReturn(PETSC_SUCCESS);
8189: }

8191: /*@C
8192:   MatRestoreRowIJ - Call after you are completed with the ia,ja indices obtained with `MatGetRowIJ()`.

8194:   Collective

8196:   Input Parameters:
8197: + mat             - the matrix
8198: . shift           - 1 or zero indicating we want the indices starting at 0 or 1
8199: . symmetric       - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized
8200: . inodecompressed - `PETSC_TRUE` or `PETSC_FALSE` indicating if the nonzero structure of the
8201:                     inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is
8202:                     always used.
8203: . n               - size of (possibly compressed) matrix
8204: . ia              - the row pointers
8205: - ja              - the column indices

8207:   Output Parameter:
8208: . done - `PETSC_TRUE` or `PETSC_FALSE` indicated that the values have been returned

8210:   Level: developer

8212:   Note:
8213:   This routine zeros out `n`, `ia`, and `ja`. This is to prevent accidental
8214:   us of the array after it has been restored. If you pass `NULL`, it will
8215:   not zero the pointers.  Use of ia or ja after `MatRestoreRowIJ()` is invalid.

8217:   Fortran Note:
8218:   `MatRestoreRowIJ()` Fortran binding is deprecated (since PETSc 3.19), use `MatRestoreRowIJF90()`

8220: .seealso: [](ch_matrices), `Mat`, `MatGetRowIJ()`, `MatRestoreRowIJF90()`, `MatRestoreColumnIJ()`
8221: @*/
8222: PetscErrorCode MatRestoreRowIJ(Mat mat, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
8223: {
8224:   PetscFunctionBegin;
8227:   if (ia) PetscAssertPointer(ia, 6);
8228:   if (ja) PetscAssertPointer(ja, 7);
8229:   if (done) PetscAssertPointer(done, 8);
8230:   MatCheckPreallocated(mat, 1);

8232:   if (!mat->ops->restorerowij && done) *done = PETSC_FALSE;
8233:   else {
8234:     if (done) *done = PETSC_TRUE;
8235:     PetscUseTypeMethod(mat, restorerowij, shift, symmetric, inodecompressed, n, ia, ja, done);
8236:     if (n) *n = 0;
8237:     if (ia) *ia = NULL;
8238:     if (ja) *ja = NULL;
8239:   }
8240:   PetscFunctionReturn(PETSC_SUCCESS);
8241: }

8243: /*@C
8244:   MatRestoreColumnIJ - Call after you are completed with the ia,ja indices obtained with `MatGetColumnIJ()`.

8246:   Collective

8248:   Input Parameters:
8249: + mat             - the matrix
8250: . shift           - 1 or zero indicating we want the indices starting at 0 or 1
8251: . symmetric       - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized
8252: - inodecompressed - `PETSC_TRUE` or `PETSC_FALSE` indicating if the nonzero structure of the
8253:                     inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is
8254:                     always used.

8256:   Output Parameters:
8257: + n    - size of (possibly compressed) matrix
8258: . ia   - the column pointers
8259: . ja   - the row indices
8260: - done - `PETSC_TRUE` or `PETSC_FALSE` indicated that the values have been returned

8262:   Level: developer

8264: .seealso: [](ch_matrices), `Mat`, `MatGetColumnIJ()`, `MatRestoreRowIJ()`
8265: @*/
8266: PetscErrorCode MatRestoreColumnIJ(Mat mat, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done)
8267: {
8268:   PetscFunctionBegin;
8271:   if (ia) PetscAssertPointer(ia, 6);
8272:   if (ja) PetscAssertPointer(ja, 7);
8273:   PetscAssertPointer(done, 8);
8274:   MatCheckPreallocated(mat, 1);

8276:   if (!mat->ops->restorecolumnij) *done = PETSC_FALSE;
8277:   else {
8278:     *done = PETSC_TRUE;
8279:     PetscUseTypeMethod(mat, restorecolumnij, shift, symmetric, inodecompressed, n, ia, ja, done);
8280:     if (n) *n = 0;
8281:     if (ia) *ia = NULL;
8282:     if (ja) *ja = NULL;
8283:   }
8284:   PetscFunctionReturn(PETSC_SUCCESS);
8285: }

8287: /*@
8288:   MatColoringPatch - Used inside matrix coloring routines that use `MatGetRowIJ()` and/or
8289:   `MatGetColumnIJ()`.

8291:   Collective

8293:   Input Parameters:
8294: + mat        - the matrix
8295: . ncolors    - maximum color value
8296: . n          - number of entries in colorarray
8297: - colorarray - array indicating color for each column

8299:   Output Parameter:
8300: . iscoloring - coloring generated using colorarray information

8302:   Level: developer

8304: .seealso: [](ch_matrices), `Mat`, `MatGetRowIJ()`, `MatGetColumnIJ()`
8305: @*/
8306: PetscErrorCode MatColoringPatch(Mat mat, PetscInt ncolors, PetscInt n, ISColoringValue colorarray[], ISColoring *iscoloring)
8307: {
8308:   PetscFunctionBegin;
8311:   PetscAssertPointer(colorarray, 4);
8312:   PetscAssertPointer(iscoloring, 5);
8313:   MatCheckPreallocated(mat, 1);

8315:   if (!mat->ops->coloringpatch) {
8316:     PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)mat), ncolors, n, colorarray, PETSC_OWN_POINTER, iscoloring));
8317:   } else {
8318:     PetscUseTypeMethod(mat, coloringpatch, ncolors, n, colorarray, iscoloring);
8319:   }
8320:   PetscFunctionReturn(PETSC_SUCCESS);
8321: }

8323: /*@
8324:   MatSetUnfactored - Resets a factored matrix to be treated as unfactored.

8326:   Logically Collective

8328:   Input Parameter:
8329: . mat - the factored matrix to be reset

8331:   Level: developer

8333:   Notes:
8334:   This routine should be used only with factored matrices formed by in-place
8335:   factorization via ILU(0) (or by in-place LU factorization for the `MATSEQDENSE`
8336:   format).  This option can save memory, for example, when solving nonlinear
8337:   systems with a matrix-free Newton-Krylov method and a matrix-based, in-place
8338:   ILU(0) preconditioner.

8340:   One can specify in-place ILU(0) factorization by calling
8341: .vb
8342:      PCType(pc,PCILU);
8343:      PCFactorSeUseInPlace(pc);
8344: .ve
8345:   or by using the options -pc_type ilu -pc_factor_in_place

8347:   In-place factorization ILU(0) can also be used as a local
8348:   solver for the blocks within the block Jacobi or additive Schwarz
8349:   methods (runtime option: -sub_pc_factor_in_place).  See Users-Manual: ch_pc
8350:   for details on setting local solver options.

8352:   Most users should employ the `KSP` interface for linear solvers
8353:   instead of working directly with matrix algebra routines such as this.
8354:   See, e.g., `KSPCreate()`.

8356: .seealso: [](ch_matrices), `Mat`, `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`
8357: @*/
8358: PetscErrorCode MatSetUnfactored(Mat mat)
8359: {
8360:   PetscFunctionBegin;
8363:   MatCheckPreallocated(mat, 1);
8364:   mat->factortype = MAT_FACTOR_NONE;
8365:   if (!mat->ops->setunfactored) PetscFunctionReturn(PETSC_SUCCESS);
8366:   PetscUseTypeMethod(mat, setunfactored);
8367:   PetscFunctionReturn(PETSC_SUCCESS);
8368: }

8370: /*MC
8371:     MatDenseGetArrayF90 - Accesses a matrix array from Fortran

8373:     Synopsis:
8374:     MatDenseGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:,:)},integer ierr)

8376:     Not Collective

8378:     Input Parameter:
8379: .   x - matrix

8381:     Output Parameters:
8382: +   xx_v - the Fortran pointer to the array
8383: -   ierr - error code

8385:     Example of Usage:
8386: .vb
8387:       PetscScalar, pointer xx_v(:,:)
8388:       ....
8389:       call MatDenseGetArrayF90(x,xx_v,ierr)
8390:       a = xx_v(3)
8391:       call MatDenseRestoreArrayF90(x,xx_v,ierr)
8392: .ve

8394:     Level: advanced

8396: .seealso: [](ch_matrices), `Mat`, `MatDenseRestoreArrayF90()`, `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatSeqAIJGetArrayF90()`
8397: M*/

8399: /*MC
8400:     MatDenseRestoreArrayF90 - Restores a matrix array that has been
8401:     accessed with `MatDenseGetArrayF90()`.

8403:     Synopsis:
8404:     MatDenseRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:,:)},integer ierr)

8406:     Not Collective

8408:     Input Parameters:
8409: +   x - matrix
8410: -   xx_v - the Fortran90 pointer to the array

8412:     Output Parameter:
8413: .   ierr - error code

8415:     Example of Usage:
8416: .vb
8417:        PetscScalar, pointer xx_v(:,:)
8418:        ....
8419:        call MatDenseGetArrayF90(x,xx_v,ierr)
8420:        a = xx_v(3)
8421:        call MatDenseRestoreArrayF90(x,xx_v,ierr)
8422: .ve

8424:     Level: advanced

8426: .seealso: [](ch_matrices), `Mat`, `MatDenseGetArrayF90()`, `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatSeqAIJRestoreArrayF90()`
8427: M*/

8429: /*MC
8430:     MatSeqAIJGetArrayF90 - Accesses a matrix array from Fortran.

8432:     Synopsis:
8433:     MatSeqAIJGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)

8435:     Not Collective

8437:     Input Parameter:
8438: .   x - matrix

8440:     Output Parameters:
8441: +   xx_v - the Fortran pointer to the array
8442: -   ierr - error code

8444:     Example of Usage:
8445: .vb
8446:       PetscScalar, pointer xx_v(:)
8447:       ....
8448:       call MatSeqAIJGetArrayF90(x,xx_v,ierr)
8449:       a = xx_v(3)
8450:       call MatSeqAIJRestoreArrayF90(x,xx_v,ierr)
8451: .ve

8453:     Level: advanced

8455: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJRestoreArrayF90()`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArray()`, `MatDenseGetArrayF90()`
8456: M*/

8458: /*MC
8459:     MatSeqAIJRestoreArrayF90 - Restores a matrix array that has been
8460:     accessed with `MatSeqAIJGetArrayF90()`.

8462:     Synopsis:
8463:     MatSeqAIJRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr)

8465:     Not Collective

8467:     Input Parameters:
8468: +   x - matrix
8469: -   xx_v - the Fortran90 pointer to the array

8471:     Output Parameter:
8472: .   ierr - error code

8474:     Example of Usage:
8475: .vb
8476:        PetscScalar, pointer xx_v(:)
8477:        ....
8478:        call MatSeqAIJGetArrayF90(x,xx_v,ierr)
8479:        a = xx_v(3)
8480:        call MatSeqAIJRestoreArrayF90(x,xx_v,ierr)
8481: .ve

8483:     Level: advanced

8485: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJGetArrayF90()`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArray()`, `MatDenseRestoreArrayF90()`
8486: M*/

8488: /*@
8489:   MatCreateSubMatrix - Gets a single submatrix on the same number of processors
8490:   as the original matrix.

8492:   Collective

8494:   Input Parameters:
8495: + mat   - the original matrix
8496: . isrow - parallel `IS` containing the rows this processor should obtain
8497: . iscol - parallel `IS` containing all columns you wish to keep. Each process should list the columns that will be in IT's "diagonal part" in the new matrix.
8498: - cll   - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`

8500:   Output Parameter:
8501: . newmat - the new submatrix, of the same type as the original matrix

8503:   Level: advanced

8505:   Notes:
8506:   The submatrix will be able to be multiplied with vectors using the same layout as `iscol`.

8508:   Some matrix types place restrictions on the row and column indices, such
8509:   as that they be sorted or that they be equal to each other. For `MATBAIJ` and `MATSBAIJ` matrices the indices must include all rows/columns of a block;
8510:   for example, if the block size is 3 one cannot select the 0 and 2 rows without selecting the 1 row.

8512:   The index sets may not have duplicate entries.

8514:   The first time this is called you should use a cll of `MAT_INITIAL_MATRIX`,
8515:   the `MatCreateSubMatrix()` routine will create the newmat for you. Any additional calls
8516:   to this routine with a mat of the same nonzero structure and with a call of `MAT_REUSE_MATRIX`
8517:   will reuse the matrix generated the first time.  You should call `MatDestroy()` on `newmat` when
8518:   you are finished using it.

8520:   The communicator of the newly obtained matrix is ALWAYS the same as the communicator of
8521:   the input matrix.

8523:   If `iscol` is `NULL` then all columns are obtained (not supported in Fortran).

8525:   If `isrow` and `iscol` have a nontrivial block-size, then the resulting matrix has this block-size as well. This feature
8526:   is used by `PCFIELDSPLIT` to allow easy nesting of its use.

8528:   Example usage:
8529:   Consider the following 8x8 matrix with 34 non-zero values, that is
8530:   assembled across 3 processors. Let's assume that proc0 owns 3 rows,
8531:   proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown
8532:   as follows
8533: .vb
8534:             1  2  0  |  0  3  0  |  0  4
8535:     Proc0   0  5  6  |  7  0  0  |  8  0
8536:             9  0 10  | 11  0  0  | 12  0
8537:     -------------------------------------
8538:            13  0 14  | 15 16 17  |  0  0
8539:     Proc1   0 18  0  | 19 20 21  |  0  0
8540:             0  0  0  | 22 23  0  | 24  0
8541:     -------------------------------------
8542:     Proc2  25 26 27  |  0  0 28  | 29  0
8543:            30  0  0  | 31 32 33  |  0 34
8544: .ve

8546:   Suppose `isrow` = [0 1 | 4 | 6 7] and `iscol` = [1 2 | 3 4 5 | 6].  The resulting submatrix is

8548: .vb
8549:             2  0  |  0  3  0  |  0
8550:     Proc0   5  6  |  7  0  0  |  8
8551:     -------------------------------
8552:     Proc1  18  0  | 19 20 21  |  0
8553:     -------------------------------
8554:     Proc2  26 27  |  0  0 28  | 29
8555:             0  0  | 31 32 33  |  0
8556: .ve

8558: .seealso: [](ch_matrices), `Mat`, `MatCreateSubMatrices()`, `MatCreateSubMatricesMPI()`, `MatCreateSubMatrixVirtual()`, `MatSubMatrixVirtualUpdate()`
8559: @*/
8560: PetscErrorCode MatCreateSubMatrix(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat)
8561: {
8562:   PetscMPIInt size;
8563:   Mat        *local;
8564:   IS          iscoltmp;
8565:   PetscBool   flg;

8567:   PetscFunctionBegin;
8571:   PetscAssertPointer(newmat, 5);
8574:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
8575:   PetscCheck(cll != MAT_IGNORE_MATRIX, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Cannot use MAT_IGNORE_MATRIX");

8577:   MatCheckPreallocated(mat, 1);
8578:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));

8580:   if (!iscol || isrow == iscol) {
8581:     PetscBool   stride;
8582:     PetscMPIInt grabentirematrix = 0, grab;
8583:     PetscCall(PetscObjectTypeCompare((PetscObject)isrow, ISSTRIDE, &stride));
8584:     if (stride) {
8585:       PetscInt first, step, n, rstart, rend;
8586:       PetscCall(ISStrideGetInfo(isrow, &first, &step));
8587:       if (step == 1) {
8588:         PetscCall(MatGetOwnershipRange(mat, &rstart, &rend));
8589:         if (rstart == first) {
8590:           PetscCall(ISGetLocalSize(isrow, &n));
8591:           if (n == rend - rstart) grabentirematrix = 1;
8592:         }
8593:       }
8594:     }
8595:     PetscCallMPI(MPIU_Allreduce(&grabentirematrix, &grab, 1, MPI_INT, MPI_MIN, PetscObjectComm((PetscObject)mat)));
8596:     if (grab) {
8597:       PetscCall(PetscInfo(mat, "Getting entire matrix as submatrix\n"));
8598:       if (cll == MAT_INITIAL_MATRIX) {
8599:         *newmat = mat;
8600:         PetscCall(PetscObjectReference((PetscObject)mat));
8601:       }
8602:       PetscFunctionReturn(PETSC_SUCCESS);
8603:     }
8604:   }

8606:   if (!iscol) {
8607:     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)mat), mat->cmap->n, mat->cmap->rstart, 1, &iscoltmp));
8608:   } else {
8609:     iscoltmp = iscol;
8610:   }

8612:   /* if original matrix is on just one processor then use submatrix generated */
8613:   if (mat->ops->createsubmatrices && !mat->ops->createsubmatrix && size == 1 && cll == MAT_REUSE_MATRIX) {
8614:     PetscCall(MatCreateSubMatrices(mat, 1, &isrow, &iscoltmp, MAT_REUSE_MATRIX, &newmat));
8615:     goto setproperties;
8616:   } else if (mat->ops->createsubmatrices && !mat->ops->createsubmatrix && size == 1) {
8617:     PetscCall(MatCreateSubMatrices(mat, 1, &isrow, &iscoltmp, MAT_INITIAL_MATRIX, &local));
8618:     *newmat = *local;
8619:     PetscCall(PetscFree(local));
8620:     goto setproperties;
8621:   } else if (!mat->ops->createsubmatrix) {
8622:     /* Create a new matrix type that implements the operation using the full matrix */
8623:     PetscCall(PetscLogEventBegin(MAT_CreateSubMat, mat, 0, 0, 0));
8624:     switch (cll) {
8625:     case MAT_INITIAL_MATRIX:
8626:       PetscCall(MatCreateSubMatrixVirtual(mat, isrow, iscoltmp, newmat));
8627:       break;
8628:     case MAT_REUSE_MATRIX:
8629:       PetscCall(MatSubMatrixVirtualUpdate(*newmat, mat, isrow, iscoltmp));
8630:       break;
8631:     default:
8632:       SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Invalid MatReuse, must be either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX");
8633:     }
8634:     PetscCall(PetscLogEventEnd(MAT_CreateSubMat, mat, 0, 0, 0));
8635:     goto setproperties;
8636:   }

8638:   PetscCall(PetscLogEventBegin(MAT_CreateSubMat, mat, 0, 0, 0));
8639:   PetscUseTypeMethod(mat, createsubmatrix, isrow, iscoltmp, cll, newmat);
8640:   PetscCall(PetscLogEventEnd(MAT_CreateSubMat, mat, 0, 0, 0));

8642: setproperties:
8643:   if ((*newmat)->symmetric == PETSC_BOOL3_UNKNOWN && (*newmat)->structurally_symmetric == PETSC_BOOL3_UNKNOWN && (*newmat)->spd == PETSC_BOOL3_UNKNOWN && (*newmat)->hermitian == PETSC_BOOL3_UNKNOWN) {
8644:     PetscCall(ISEqualUnsorted(isrow, iscoltmp, &flg));
8645:     if (flg) PetscCall(MatPropagateSymmetryOptions(mat, *newmat));
8646:   }
8647:   if (!iscol) PetscCall(ISDestroy(&iscoltmp));
8648:   if (*newmat && cll == MAT_INITIAL_MATRIX) PetscCall(PetscObjectStateIncrease((PetscObject)*newmat));
8649:   PetscFunctionReturn(PETSC_SUCCESS);
8650: }

8652: /*@
8653:   MatPropagateSymmetryOptions - Propagates symmetry options set on a matrix to another matrix

8655:   Not Collective

8657:   Input Parameters:
8658: + A - the matrix we wish to propagate options from
8659: - B - the matrix we wish to propagate options to

8661:   Level: beginner

8663:   Note:
8664:   Propagates the options associated to `MAT_SYMMETRY_ETERNAL`, `MAT_STRUCTURALLY_SYMMETRIC`, `MAT_HERMITIAN`, `MAT_SPD`, `MAT_SYMMETRIC`, and `MAT_STRUCTURAL_SYMMETRY_ETERNAL`

8666: .seealso: [](ch_matrices), `Mat`, `MatSetOption()`, `MatIsSymmetricKnown()`, `MatIsSPDKnown()`, `MatIsHermitianKnown()`, `MatIsStructurallySymmetricKnown()`
8667: @*/
8668: PetscErrorCode MatPropagateSymmetryOptions(Mat A, Mat B)
8669: {
8670:   PetscFunctionBegin;
8673:   B->symmetry_eternal            = A->symmetry_eternal;
8674:   B->structural_symmetry_eternal = A->structural_symmetry_eternal;
8675:   B->symmetric                   = A->symmetric;
8676:   B->structurally_symmetric      = A->structurally_symmetric;
8677:   B->spd                         = A->spd;
8678:   B->hermitian                   = A->hermitian;
8679:   PetscFunctionReturn(PETSC_SUCCESS);
8680: }

8682: /*@
8683:   MatStashSetInitialSize - sets the sizes of the matrix stash, that is
8684:   used during the assembly process to store values that belong to
8685:   other processors.

8687:   Not Collective

8689:   Input Parameters:
8690: + mat   - the matrix
8691: . size  - the initial size of the stash.
8692: - bsize - the initial size of the block-stash(if used).

8694:   Options Database Keys:
8695: + -matstash_initial_size <size> or <size0,size1,...sizep-1>            - set initial size
8696: - -matstash_block_initial_size <bsize>  or <bsize0,bsize1,...bsizep-1> - set initial block size

8698:   Level: intermediate

8700:   Notes:
8701:   The block-stash is used for values set with `MatSetValuesBlocked()` while
8702:   the stash is used for values set with `MatSetValues()`

8704:   Run with the option -info and look for output of the form
8705:   MatAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs.
8706:   to determine the appropriate value, MM, to use for size and
8707:   MatAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs.
8708:   to determine the value, BMM to use for bsize

8710: .seealso: [](ch_matrices), `MatAssemblyBegin()`, `MatAssemblyEnd()`, `Mat`, `MatStashGetInfo()`
8711: @*/
8712: PetscErrorCode MatStashSetInitialSize(Mat mat, PetscInt size, PetscInt bsize)
8713: {
8714:   PetscFunctionBegin;
8717:   PetscCall(MatStashSetInitialSize_Private(&mat->stash, size));
8718:   PetscCall(MatStashSetInitialSize_Private(&mat->bstash, bsize));
8719:   PetscFunctionReturn(PETSC_SUCCESS);
8720: }

8722: /*@
8723:   MatInterpolateAdd - $w = y + A*x$ or $A^T*x$ depending on the shape of
8724:   the matrix

8726:   Neighbor-wise Collective

8728:   Input Parameters:
8729: + A - the matrix
8730: . x - the vector to be multiplied by the interpolation operator
8731: - y - the vector to be added to the result

8733:   Output Parameter:
8734: . w - the resulting vector

8736:   Level: intermediate

8738:   Notes:
8739:   `w` may be the same vector as `y`.

8741:   This allows one to use either the restriction or interpolation (its transpose)
8742:   matrix to do the interpolation

8744: .seealso: [](ch_matrices), `Mat`, `MatMultAdd()`, `MatMultTransposeAdd()`, `MatRestrict()`, `PCMG`
8745: @*/
8746: PetscErrorCode MatInterpolateAdd(Mat A, Vec x, Vec y, Vec w)
8747: {
8748:   PetscInt M, N, Ny;

8750:   PetscFunctionBegin;
8755:   PetscCall(MatGetSize(A, &M, &N));
8756:   PetscCall(VecGetSize(y, &Ny));
8757:   if (M == Ny) {
8758:     PetscCall(MatMultAdd(A, x, y, w));
8759:   } else {
8760:     PetscCall(MatMultTransposeAdd(A, x, y, w));
8761:   }
8762:   PetscFunctionReturn(PETSC_SUCCESS);
8763: }

8765: /*@
8766:   MatInterpolate - $y = A*x$ or $A^T*x$ depending on the shape of
8767:   the matrix

8769:   Neighbor-wise Collective

8771:   Input Parameters:
8772: + A - the matrix
8773: - x - the vector to be interpolated

8775:   Output Parameter:
8776: . y - the resulting vector

8778:   Level: intermediate

8780:   Note:
8781:   This allows one to use either the restriction or interpolation (its transpose)
8782:   matrix to do the interpolation

8784: .seealso: [](ch_matrices), `Mat`, `MatMultAdd()`, `MatMultTransposeAdd()`, `MatRestrict()`, `PCMG`
8785: @*/
8786: PetscErrorCode MatInterpolate(Mat A, Vec x, Vec y)
8787: {
8788:   PetscInt M, N, Ny;

8790:   PetscFunctionBegin;
8794:   PetscCall(MatGetSize(A, &M, &N));
8795:   PetscCall(VecGetSize(y, &Ny));
8796:   if (M == Ny) {
8797:     PetscCall(MatMult(A, x, y));
8798:   } else {
8799:     PetscCall(MatMultTranspose(A, x, y));
8800:   }
8801:   PetscFunctionReturn(PETSC_SUCCESS);
8802: }

8804: /*@
8805:   MatRestrict - $y = A*x$ or $A^T*x$

8807:   Neighbor-wise Collective

8809:   Input Parameters:
8810: + A - the matrix
8811: - x - the vector to be restricted

8813:   Output Parameter:
8814: . y - the resulting vector

8816:   Level: intermediate

8818:   Note:
8819:   This allows one to use either the restriction or interpolation (its transpose)
8820:   matrix to do the restriction

8822: .seealso: [](ch_matrices), `Mat`, `MatMultAdd()`, `MatMultTransposeAdd()`, `MatInterpolate()`, `PCMG`
8823: @*/
8824: PetscErrorCode MatRestrict(Mat A, Vec x, Vec y)
8825: {
8826:   PetscInt M, N, Nx;

8828:   PetscFunctionBegin;
8832:   PetscCall(MatGetSize(A, &M, &N));
8833:   PetscCall(VecGetSize(x, &Nx));
8834:   if (M == Nx) {
8835:     PetscCall(MatMultTranspose(A, x, y));
8836:   } else {
8837:     PetscCall(MatMult(A, x, y));
8838:   }
8839:   PetscFunctionReturn(PETSC_SUCCESS);
8840: }

8842: /*@
8843:   MatMatInterpolateAdd - $Y = W + A*X$ or $W + A^T*X$ depending on the shape of `A`

8845:   Neighbor-wise Collective

8847:   Input Parameters:
8848: + A - the matrix
8849: . x - the input dense matrix to be multiplied
8850: - w - the input dense matrix to be added to the result

8852:   Output Parameter:
8853: . y - the output dense matrix

8855:   Level: intermediate

8857:   Note:
8858:   This allows one to use either the restriction or interpolation (its transpose)
8859:   matrix to do the interpolation. `y` matrix can be reused if already created with the proper sizes,
8860:   otherwise it will be recreated. `y` must be initialized to `NULL` if not supplied.

8862: .seealso: [](ch_matrices), `Mat`, `MatInterpolateAdd()`, `MatMatInterpolate()`, `MatMatRestrict()`, `PCMG`
8863: @*/
8864: PetscErrorCode MatMatInterpolateAdd(Mat A, Mat x, Mat w, Mat *y)
8865: {
8866:   PetscInt  M, N, Mx, Nx, Mo, My = 0, Ny = 0;
8867:   PetscBool trans = PETSC_TRUE;
8868:   MatReuse  reuse = MAT_INITIAL_MATRIX;

8870:   PetscFunctionBegin;
8876:   PetscCall(MatGetSize(A, &M, &N));
8877:   PetscCall(MatGetSize(x, &Mx, &Nx));
8878:   if (N == Mx) trans = PETSC_FALSE;
8879:   else PetscCheck(M == Mx, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Size mismatch: A %" PetscInt_FMT "x%" PetscInt_FMT ", X %" PetscInt_FMT "x%" PetscInt_FMT, M, N, Mx, Nx);
8880:   Mo = trans ? N : M;
8881:   if (*y) {
8882:     PetscCall(MatGetSize(*y, &My, &Ny));
8883:     if (Mo == My && Nx == Ny) {
8884:       reuse = MAT_REUSE_MATRIX;
8885:     } else {
8886:       PetscCheck(w || *y != w, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot reuse y and w, size mismatch: A %" PetscInt_FMT "x%" PetscInt_FMT ", X %" PetscInt_FMT "x%" PetscInt_FMT ", Y %" PetscInt_FMT "x%" PetscInt_FMT, M, N, Mx, Nx, My, Ny);
8887:       PetscCall(MatDestroy(y));
8888:     }
8889:   }

8891:   if (w && *y == w) { /* this is to minimize changes in PCMG */
8892:     PetscBool flg;

8894:     PetscCall(PetscObjectQuery((PetscObject)*y, "__MatMatIntAdd_w", (PetscObject *)&w));
8895:     if (w) {
8896:       PetscInt My, Ny, Mw, Nw;

8898:       PetscCall(PetscObjectTypeCompare((PetscObject)*y, ((PetscObject)w)->type_name, &flg));
8899:       PetscCall(MatGetSize(*y, &My, &Ny));
8900:       PetscCall(MatGetSize(w, &Mw, &Nw));
8901:       if (!flg || My != Mw || Ny != Nw) w = NULL;
8902:     }
8903:     if (!w) {
8904:       PetscCall(MatDuplicate(*y, MAT_COPY_VALUES, &w));
8905:       PetscCall(PetscObjectCompose((PetscObject)*y, "__MatMatIntAdd_w", (PetscObject)w));
8906:       PetscCall(PetscObjectDereference((PetscObject)w));
8907:     } else {
8908:       PetscCall(MatCopy(*y, w, UNKNOWN_NONZERO_PATTERN));
8909:     }
8910:   }
8911:   if (!trans) {
8912:     PetscCall(MatMatMult(A, x, reuse, PETSC_DETERMINE, y));
8913:   } else {
8914:     PetscCall(MatTransposeMatMult(A, x, reuse, PETSC_DETERMINE, y));
8915:   }
8916:   if (w) PetscCall(MatAXPY(*y, 1.0, w, UNKNOWN_NONZERO_PATTERN));
8917:   PetscFunctionReturn(PETSC_SUCCESS);
8918: }

8920: /*@
8921:   MatMatInterpolate - $Y = A*X$ or $A^T*X$ depending on the shape of `A`

8923:   Neighbor-wise Collective

8925:   Input Parameters:
8926: + A - the matrix
8927: - x - the input dense matrix

8929:   Output Parameter:
8930: . y - the output dense matrix

8932:   Level: intermediate

8934:   Note:
8935:   This allows one to use either the restriction or interpolation (its transpose)
8936:   matrix to do the interpolation. `y` matrix can be reused if already created with the proper sizes,
8937:   otherwise it will be recreated. `y` must be initialized to `NULL` if not supplied.

8939: .seealso: [](ch_matrices), `Mat`, `MatInterpolate()`, `MatRestrict()`, `MatMatRestrict()`, `PCMG`
8940: @*/
8941: PetscErrorCode MatMatInterpolate(Mat A, Mat x, Mat *y)
8942: {
8943:   PetscFunctionBegin;
8944:   PetscCall(MatMatInterpolateAdd(A, x, NULL, y));
8945:   PetscFunctionReturn(PETSC_SUCCESS);
8946: }

8948: /*@
8949:   MatMatRestrict - $Y = A*X$ or $A^T*X$ depending on the shape of `A`

8951:   Neighbor-wise Collective

8953:   Input Parameters:
8954: + A - the matrix
8955: - x - the input dense matrix

8957:   Output Parameter:
8958: . y - the output dense matrix

8960:   Level: intermediate

8962:   Note:
8963:   This allows one to use either the restriction or interpolation (its transpose)
8964:   matrix to do the restriction. `y` matrix can be reused if already created with the proper sizes,
8965:   otherwise it will be recreated. `y` must be initialized to `NULL` if not supplied.

8967: .seealso: [](ch_matrices), `Mat`, `MatRestrict()`, `MatInterpolate()`, `MatMatInterpolate()`, `PCMG`
8968: @*/
8969: PetscErrorCode MatMatRestrict(Mat A, Mat x, Mat *y)
8970: {
8971:   PetscFunctionBegin;
8972:   PetscCall(MatMatInterpolateAdd(A, x, NULL, y));
8973:   PetscFunctionReturn(PETSC_SUCCESS);
8974: }

8976: /*@
8977:   MatGetNullSpace - retrieves the null space of a matrix.

8979:   Logically Collective

8981:   Input Parameters:
8982: + mat    - the matrix
8983: - nullsp - the null space object

8985:   Level: developer

8987: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatSetNullSpace()`, `MatNullSpace`
8988: @*/
8989: PetscErrorCode MatGetNullSpace(Mat mat, MatNullSpace *nullsp)
8990: {
8991:   PetscFunctionBegin;
8993:   PetscAssertPointer(nullsp, 2);
8994:   *nullsp = (mat->symmetric == PETSC_BOOL3_TRUE && !mat->nullsp) ? mat->transnullsp : mat->nullsp;
8995:   PetscFunctionReturn(PETSC_SUCCESS);
8996: }

8998: /*@C
8999:   MatGetNullSpaces - gets the null spaces, transpose null spaces, and near null spaces from an array of matrices

9001:   Logically Collective

9003:   Input Parameters:
9004: + n   - the number of matrices
9005: - mat - the array of matrices

9007:   Output Parameters:
9008: . nullsp - an array of null spaces, `NULL` for each matrix that does not have a null space, length 3 * `n`

9010:   Level: developer

9012:   Note:
9013:   Call `MatRestoreNullspaces()` to provide these to another array of matrices

9015: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatGetNullSpace()`, `MatSetTransposeNullSpace()`, `MatGetTransposeNullSpace()`,
9016:           `MatNullSpaceRemove()`, `MatRestoreNullSpaces()`
9017: @*/
9018: PetscErrorCode MatGetNullSpaces(PetscInt n, Mat mat[], MatNullSpace *nullsp[])
9019: {
9020:   PetscFunctionBegin;
9021:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of matrices %" PetscInt_FMT " must be non-negative", n);
9022:   PetscAssertPointer(mat, 2);
9023:   PetscAssertPointer(nullsp, 3);

9025:   PetscCall(PetscCalloc1(3 * n, nullsp));
9026:   for (PetscInt i = 0; i < n; i++) {
9028:     (*nullsp)[i] = mat[i]->nullsp;
9029:     PetscCall(PetscObjectReference((PetscObject)(*nullsp)[i]));
9030:     (*nullsp)[n + i] = mat[i]->nearnullsp;
9031:     PetscCall(PetscObjectReference((PetscObject)(*nullsp)[n + i]));
9032:     (*nullsp)[2 * n + i] = mat[i]->transnullsp;
9033:     PetscCall(PetscObjectReference((PetscObject)(*nullsp)[2 * n + i]));
9034:   }
9035:   PetscFunctionReturn(PETSC_SUCCESS);
9036: }

9038: /*@C
9039:   MatRestoreNullSpaces - sets the null spaces, transpose null spaces, and near null spaces obtained with `MatGetNullSpaces()` for an array of matrices

9041:   Logically Collective

9043:   Input Parameters:
9044: + n      - the number of matrices
9045: . mat    - the array of matrices
9046: - nullsp - an array of null spaces

9048:   Level: developer

9050:   Note:
9051:   Call `MatGetNullSpaces()` to create `nullsp`

9053: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatGetNullSpace()`, `MatSetTransposeNullSpace()`, `MatGetTransposeNullSpace()`,
9054:           `MatNullSpaceRemove()`, `MatGetNullSpaces()`
9055: @*/
9056: PetscErrorCode MatRestoreNullSpaces(PetscInt n, Mat mat[], MatNullSpace *nullsp[])
9057: {
9058:   PetscFunctionBegin;
9059:   PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of matrices %" PetscInt_FMT " must be non-negative", n);
9060:   PetscAssertPointer(mat, 2);
9061:   PetscAssertPointer(nullsp, 3);
9062:   PetscAssertPointer(*nullsp, 3);

9064:   for (PetscInt i = 0; i < n; i++) {
9066:     PetscCall(MatSetNullSpace(mat[i], (*nullsp)[i]));
9067:     PetscCall(PetscObjectDereference((PetscObject)(*nullsp)[i]));
9068:     PetscCall(MatSetNearNullSpace(mat[i], (*nullsp)[n + i]));
9069:     PetscCall(PetscObjectDereference((PetscObject)(*nullsp)[n + i]));
9070:     PetscCall(MatSetTransposeNullSpace(mat[i], (*nullsp)[2 * n + i]));
9071:     PetscCall(PetscObjectDereference((PetscObject)(*nullsp)[2 * n + i]));
9072:   }
9073:   PetscCall(PetscFree(*nullsp));
9074:   PetscFunctionReturn(PETSC_SUCCESS);
9075: }

9077: /*@
9078:   MatSetNullSpace - attaches a null space to a matrix.

9080:   Logically Collective

9082:   Input Parameters:
9083: + mat    - the matrix
9084: - nullsp - the null space object

9086:   Level: advanced

9088:   Notes:
9089:   This null space is used by the `KSP` linear solvers to solve singular systems.

9091:   Overwrites any previous null space that may have been attached. You can remove the null space from the matrix object by calling this routine with an nullsp of `NULL`

9093:   For inconsistent singular systems (linear systems where the right-hand side is not in the range of the operator) the `KSP` residuals will not converge
9094:   to zero but the linear system will still be solved in a least squares sense.

9096:   The fundamental theorem of linear algebra (Gilbert Strang, Introduction to Applied Mathematics, page 72) states that
9097:   the domain of a matrix A (from $R^n$ to $R^m$ (m rows, n columns) $R^n$ = the direct sum of the null space of A, n(A), + the range of $A^T$, $R(A^T)$.
9098:   Similarly $R^m$ = direct sum n($A^T$) + R(A).  Hence the linear system $A x = b$ has a solution only if b in R(A) (or correspondingly b is orthogonal to
9099:   n($A^T$)) and if x is a solution then x + alpha n(A) is a solution for any alpha. The minimum norm solution is orthogonal to n(A). For problems without a solution
9100:   the solution that minimizes the norm of the residual (the least squares solution) can be obtained by solving A x = \hat{b} where \hat{b} is b orthogonalized to the n($A^T$).
9101:   This  \hat{b} can be obtained by calling `MatNullSpaceRemove()` with the null space of the transpose of the matrix.

9103:   If the matrix is known to be symmetric because it is an `MATSBAIJ` matrix or one as called
9104:   `MatSetOption`(mat,`MAT_SYMMETRIC` or possibly `MAT_SYMMETRY_ETERNAL`,`PETSC_TRUE`); this
9105:   routine also automatically calls `MatSetTransposeNullSpace()`.

9107:   The user should call `MatNullSpaceDestroy()`.

9109: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatGetNullSpace()`, `MatSetTransposeNullSpace()`, `MatGetTransposeNullSpace()`, `MatNullSpaceRemove()`,
9110:           `KSPSetPCSide()`
9111: @*/
9112: PetscErrorCode MatSetNullSpace(Mat mat, MatNullSpace nullsp)
9113: {
9114:   PetscFunctionBegin;
9117:   if (nullsp) PetscCall(PetscObjectReference((PetscObject)nullsp));
9118:   PetscCall(MatNullSpaceDestroy(&mat->nullsp));
9119:   mat->nullsp = nullsp;
9120:   if (mat->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetTransposeNullSpace(mat, nullsp));
9121:   PetscFunctionReturn(PETSC_SUCCESS);
9122: }

9124: /*@
9125:   MatGetTransposeNullSpace - retrieves the null space of the transpose of a matrix.

9127:   Logically Collective

9129:   Input Parameters:
9130: + mat    - the matrix
9131: - nullsp - the null space object

9133:   Level: developer

9135: .seealso: [](ch_matrices), `Mat`, `MatNullSpace`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatSetTransposeNullSpace()`, `MatSetNullSpace()`, `MatGetNullSpace()`
9136: @*/
9137: PetscErrorCode MatGetTransposeNullSpace(Mat mat, MatNullSpace *nullsp)
9138: {
9139:   PetscFunctionBegin;
9142:   PetscAssertPointer(nullsp, 2);
9143:   *nullsp = (mat->symmetric == PETSC_BOOL3_TRUE && !mat->transnullsp) ? mat->nullsp : mat->transnullsp;
9144:   PetscFunctionReturn(PETSC_SUCCESS);
9145: }

9147: /*@
9148:   MatSetTransposeNullSpace - attaches the null space of a transpose of a matrix to the matrix

9150:   Logically Collective

9152:   Input Parameters:
9153: + mat    - the matrix
9154: - nullsp - the null space object

9156:   Level: advanced

9158:   Notes:
9159:   This allows solving singular linear systems defined by the transpose of the matrix using `KSP` solvers with left preconditioning.

9161:   See `MatSetNullSpace()`

9163: .seealso: [](ch_matrices), `Mat`, `MatNullSpace`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatGetNullSpace()`, `MatSetNullSpace()`, `MatGetTransposeNullSpace()`, `MatNullSpaceRemove()`, `KSPSetPCSide()`
9164: @*/
9165: PetscErrorCode MatSetTransposeNullSpace(Mat mat, MatNullSpace nullsp)
9166: {
9167:   PetscFunctionBegin;
9170:   if (nullsp) PetscCall(PetscObjectReference((PetscObject)nullsp));
9171:   PetscCall(MatNullSpaceDestroy(&mat->transnullsp));
9172:   mat->transnullsp = nullsp;
9173:   PetscFunctionReturn(PETSC_SUCCESS);
9174: }

9176: /*@
9177:   MatSetNearNullSpace - attaches a null space to a matrix, which is often the null space (rigid body modes) of the operator without boundary conditions
9178:   This null space will be used to provide near null space vectors to a multigrid preconditioner built from this matrix.

9180:   Logically Collective

9182:   Input Parameters:
9183: + mat    - the matrix
9184: - nullsp - the null space object

9186:   Level: advanced

9188:   Notes:
9189:   Overwrites any previous near null space that may have been attached

9191:   You can remove the null space by calling this routine with an `nullsp` of `NULL`

9193: .seealso: [](ch_matrices), `Mat`, `MatNullSpace`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNullSpace()`, `MatNullSpaceCreateRigidBody()`, `MatGetNearNullSpace()`
9194: @*/
9195: PetscErrorCode MatSetNearNullSpace(Mat mat, MatNullSpace nullsp)
9196: {
9197:   PetscFunctionBegin;
9201:   MatCheckPreallocated(mat, 1);
9202:   if (nullsp) PetscCall(PetscObjectReference((PetscObject)nullsp));
9203:   PetscCall(MatNullSpaceDestroy(&mat->nearnullsp));
9204:   mat->nearnullsp = nullsp;
9205:   PetscFunctionReturn(PETSC_SUCCESS);
9206: }

9208: /*@
9209:   MatGetNearNullSpace - Get null space attached with `MatSetNearNullSpace()`

9211:   Not Collective

9213:   Input Parameter:
9214: . mat - the matrix

9216:   Output Parameter:
9217: . nullsp - the null space object, `NULL` if not set

9219:   Level: advanced

9221: .seealso: [](ch_matrices), `Mat`, `MatNullSpace`, `MatSetNearNullSpace()`, `MatGetNullSpace()`, `MatNullSpaceCreate()`
9222: @*/
9223: PetscErrorCode MatGetNearNullSpace(Mat mat, MatNullSpace *nullsp)
9224: {
9225:   PetscFunctionBegin;
9228:   PetscAssertPointer(nullsp, 2);
9229:   MatCheckPreallocated(mat, 1);
9230:   *nullsp = mat->nearnullsp;
9231:   PetscFunctionReturn(PETSC_SUCCESS);
9232: }

9234: /*@
9235:   MatICCFactor - Performs in-place incomplete Cholesky factorization of matrix.

9237:   Collective

9239:   Input Parameters:
9240: + mat  - the matrix
9241: . row  - row/column permutation
9242: - info - information on desired factorization process

9244:   Level: developer

9246:   Notes:
9247:   Probably really in-place only when level of fill is zero, otherwise allocates
9248:   new space to store factored matrix and deletes previous memory.

9250:   Most users should employ the `KSP` interface for linear solvers
9251:   instead of working directly with matrix algebra routines such as this.
9252:   See, e.g., `KSPCreate()`.

9254:   Developer Note:
9255:   The Fortran interface is not autogenerated as the
9256:   interface definition cannot be generated correctly [due to `MatFactorInfo`]

9258: .seealso: [](ch_matrices), `Mat`, `MatFactorInfo`, `MatGetFactor()`, `MatICCFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`
9259: @*/
9260: PetscErrorCode MatICCFactor(Mat mat, IS row, const MatFactorInfo *info)
9261: {
9262:   PetscFunctionBegin;
9266:   PetscAssertPointer(info, 3);
9267:   PetscCheck(mat->rmap->N == mat->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "matrix must be square");
9268:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
9269:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
9270:   MatCheckPreallocated(mat, 1);
9271:   PetscUseTypeMethod(mat, iccfactor, row, info);
9272:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
9273:   PetscFunctionReturn(PETSC_SUCCESS);
9274: }

9276: /*@
9277:   MatDiagonalScaleLocal - Scales columns of a matrix given the scaling values including the
9278:   ghosted ones.

9280:   Not Collective

9282:   Input Parameters:
9283: + mat  - the matrix
9284: - diag - the diagonal values, including ghost ones

9286:   Level: developer

9288:   Notes:
9289:   Works only for `MATMPIAIJ` and `MATMPIBAIJ` matrices

9291:   This allows one to avoid during communication to perform the scaling that must be done with `MatDiagonalScale()`

9293: .seealso: [](ch_matrices), `Mat`, `MatDiagonalScale()`
9294: @*/
9295: PetscErrorCode MatDiagonalScaleLocal(Mat mat, Vec diag)
9296: {
9297:   PetscMPIInt size;

9299:   PetscFunctionBegin;

9304:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be already assembled");
9305:   PetscCall(PetscLogEventBegin(MAT_Scale, mat, 0, 0, 0));
9306:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
9307:   if (size == 1) {
9308:     PetscInt n, m;
9309:     PetscCall(VecGetSize(diag, &n));
9310:     PetscCall(MatGetSize(mat, NULL, &m));
9311:     PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supported for sequential matrices when no ghost points/periodic conditions");
9312:     PetscCall(MatDiagonalScale(mat, NULL, diag));
9313:   } else {
9314:     PetscUseMethod(mat, "MatDiagonalScaleLocal_C", (Mat, Vec), (mat, diag));
9315:   }
9316:   PetscCall(PetscLogEventEnd(MAT_Scale, mat, 0, 0, 0));
9317:   PetscCall(PetscObjectStateIncrease((PetscObject)mat));
9318:   PetscFunctionReturn(PETSC_SUCCESS);
9319: }

9321: /*@
9322:   MatGetInertia - Gets the inertia from a factored matrix

9324:   Collective

9326:   Input Parameter:
9327: . mat - the matrix

9329:   Output Parameters:
9330: + nneg  - number of negative eigenvalues
9331: . nzero - number of zero eigenvalues
9332: - npos  - number of positive eigenvalues

9334:   Level: advanced

9336:   Note:
9337:   Matrix must have been factored by `MatCholeskyFactor()`

9339: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatCholeskyFactor()`
9340: @*/
9341: PetscErrorCode MatGetInertia(Mat mat, PetscInt *nneg, PetscInt *nzero, PetscInt *npos)
9342: {
9343:   PetscFunctionBegin;
9346:   PetscCheck(mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Unfactored matrix");
9347:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Numeric factor mat is not assembled");
9348:   PetscUseTypeMethod(mat, getinertia, nneg, nzero, npos);
9349:   PetscFunctionReturn(PETSC_SUCCESS);
9350: }

9352: /*@C
9353:   MatSolves - Solves $A x = b$, given a factored matrix, for a collection of vectors

9355:   Neighbor-wise Collective

9357:   Input Parameters:
9358: + mat - the factored matrix obtained with `MatGetFactor()`
9359: - b   - the right-hand-side vectors

9361:   Output Parameter:
9362: . x - the result vectors

9364:   Level: developer

9366:   Note:
9367:   The vectors `b` and `x` cannot be the same.  I.e., one cannot
9368:   call `MatSolves`(A,x,x).

9370: .seealso: [](ch_matrices), `Mat`, `Vecs`, `MatSolveAdd()`, `MatSolveTranspose()`, `MatSolveTransposeAdd()`, `MatSolve()`
9371: @*/
9372: PetscErrorCode MatSolves(Mat mat, Vecs b, Vecs x)
9373: {
9374:   PetscFunctionBegin;
9377:   PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors");
9378:   PetscCheck(mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Unfactored matrix");
9379:   if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);

9381:   MatCheckPreallocated(mat, 1);
9382:   PetscCall(PetscLogEventBegin(MAT_Solves, mat, 0, 0, 0));
9383:   PetscUseTypeMethod(mat, solves, b, x);
9384:   PetscCall(PetscLogEventEnd(MAT_Solves, mat, 0, 0, 0));
9385:   PetscFunctionReturn(PETSC_SUCCESS);
9386: }

9388: /*@
9389:   MatIsSymmetric - Test whether a matrix is symmetric

9391:   Collective

9393:   Input Parameters:
9394: + A   - the matrix to test
9395: - tol - difference between value and its transpose less than this amount counts as equal (use 0.0 for exact transpose)

9397:   Output Parameter:
9398: . flg - the result

9400:   Level: intermediate

9402:   Notes:
9403:   For real numbers `MatIsSymmetric()` and `MatIsHermitian()` return identical results

9405:   If the matrix does not yet know if it is symmetric or not this can be an expensive operation, also available `MatIsSymmetricKnown()`

9407:   One can declare that a matrix is symmetric with `MatSetOption`(mat,`MAT_SYMMETRIC`,`PETSC_TRUE`) and if it is known to remain symmetric
9408:   after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SYMMETRY_ETERNAL`,`PETSC_TRUE`)

9410: .seealso: [](ch_matrices), `Mat`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetricKnown()`,
9411:           `MAT_SYMMETRIC`, `MAT_SYMMETRY_ETERNAL`
9412: @*/
9413: PetscErrorCode MatIsSymmetric(Mat A, PetscReal tol, PetscBool *flg)
9414: {
9415:   PetscFunctionBegin;
9417:   PetscAssertPointer(flg, 3);
9418:   if (A->symmetric != PETSC_BOOL3_UNKNOWN && !tol) *flg = PetscBool3ToBool(A->symmetric);
9419:   else {
9420:     if (A->ops->issymmetric) PetscUseTypeMethod(A, issymmetric, tol, flg);
9421:     else PetscCall(MatIsTranspose(A, A, tol, flg));
9422:     if (!tol) PetscCall(MatSetOption(A, MAT_SYMMETRIC, *flg));
9423:   }
9424:   PetscFunctionReturn(PETSC_SUCCESS);
9425: }

9427: /*@
9428:   MatIsHermitian - Test whether a matrix is Hermitian

9430:   Collective

9432:   Input Parameters:
9433: + A   - the matrix to test
9434: - tol - difference between value and its transpose less than this amount counts as equal (use 0.0 for exact Hermitian)

9436:   Output Parameter:
9437: . flg - the result

9439:   Level: intermediate

9441:   Notes:
9442:   For real numbers `MatIsSymmetric()` and `MatIsHermitian()` return identical results

9444:   If the matrix does not yet know if it is Hermitian or not this can be an expensive operation, also available `MatIsHermitianKnown()`

9446:   One can declare that a matrix is Hermitian with `MatSetOption`(mat,`MAT_HERMITIAN`,`PETSC_TRUE`) and if it is known to remain Hermitian
9447:   after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SYMEMTRY_ETERNAL`,`PETSC_TRUE`)

9449: .seealso: [](ch_matrices), `Mat`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitianKnown()`, `MatIsStructurallySymmetric()`, `MatSetOption()`,
9450:           `MatIsSymmetricKnown()`, `MatIsSymmetric()`, `MAT_HERMITIAN`, `MAT_SYMMETRY_ETERNAL`
9451: @*/
9452: PetscErrorCode MatIsHermitian(Mat A, PetscReal tol, PetscBool *flg)
9453: {
9454:   PetscFunctionBegin;
9456:   PetscAssertPointer(flg, 3);
9457:   if (A->hermitian != PETSC_BOOL3_UNKNOWN && !tol) *flg = PetscBool3ToBool(A->hermitian);
9458:   else {
9459:     if (A->ops->ishermitian) PetscUseTypeMethod(A, ishermitian, tol, flg);
9460:     else PetscCall(MatIsHermitianTranspose(A, A, tol, flg));
9461:     if (!tol) PetscCall(MatSetOption(A, MAT_HERMITIAN, *flg));
9462:   }
9463:   PetscFunctionReturn(PETSC_SUCCESS);
9464: }

9466: /*@
9467:   MatIsSymmetricKnown - Checks if a matrix knows if it is symmetric or not and its symmetric state

9469:   Not Collective

9471:   Input Parameter:
9472: . A - the matrix to check

9474:   Output Parameters:
9475: + set - `PETSC_TRUE` if the matrix knows its symmetry state (this tells you if the next flag is valid)
9476: - flg - the result (only valid if set is `PETSC_TRUE`)

9478:   Level: advanced

9480:   Notes:
9481:   Does not check the matrix values directly, so this may return unknown (set = `PETSC_FALSE`). Use `MatIsSymmetric()`
9482:   if you want it explicitly checked

9484:   One can declare that a matrix is symmetric with `MatSetOption`(mat,`MAT_SYMMETRIC`,`PETSC_TRUE`) and if it is known to remain symmetric
9485:   after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SYMMETRY_ETERNAL`,`PETSC_TRUE`)

9487: .seealso: [](ch_matrices), `Mat`, `MAT_SYMMETRY_ETERNAL`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetric()`, `MatIsHermitianKnown()`
9488: @*/
9489: PetscErrorCode MatIsSymmetricKnown(Mat A, PetscBool *set, PetscBool *flg)
9490: {
9491:   PetscFunctionBegin;
9493:   PetscAssertPointer(set, 2);
9494:   PetscAssertPointer(flg, 3);
9495:   if (A->symmetric != PETSC_BOOL3_UNKNOWN) {
9496:     *set = PETSC_TRUE;
9497:     *flg = PetscBool3ToBool(A->symmetric);
9498:   } else {
9499:     *set = PETSC_FALSE;
9500:   }
9501:   PetscFunctionReturn(PETSC_SUCCESS);
9502: }

9504: /*@
9505:   MatIsSPDKnown - Checks if a matrix knows if it is symmetric positive definite or not and its symmetric positive definite state

9507:   Not Collective

9509:   Input Parameter:
9510: . A - the matrix to check

9512:   Output Parameters:
9513: + set - `PETSC_TRUE` if the matrix knows its symmetric positive definite state (this tells you if the next flag is valid)
9514: - flg - the result (only valid if set is `PETSC_TRUE`)

9516:   Level: advanced

9518:   Notes:
9519:   Does not check the matrix values directly, so this may return unknown (set = `PETSC_FALSE`).

9521:   One can declare that a matrix is SPD with `MatSetOption`(mat,`MAT_SPD`,`PETSC_TRUE`) and if it is known to remain SPD
9522:   after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SPD_ETERNAL`,`PETSC_TRUE`)

9524: .seealso: [](ch_matrices), `Mat`, `MAT_SPD_ETERNAL`, `MAT_SPD`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetric()`, `MatIsHermitianKnown()`
9525: @*/
9526: PetscErrorCode MatIsSPDKnown(Mat A, PetscBool *set, PetscBool *flg)
9527: {
9528:   PetscFunctionBegin;
9530:   PetscAssertPointer(set, 2);
9531:   PetscAssertPointer(flg, 3);
9532:   if (A->spd != PETSC_BOOL3_UNKNOWN) {
9533:     *set = PETSC_TRUE;
9534:     *flg = PetscBool3ToBool(A->spd);
9535:   } else {
9536:     *set = PETSC_FALSE;
9537:   }
9538:   PetscFunctionReturn(PETSC_SUCCESS);
9539: }

9541: /*@
9542:   MatIsHermitianKnown - Checks if a matrix knows if it is Hermitian or not and its Hermitian state

9544:   Not Collective

9546:   Input Parameter:
9547: . A - the matrix to check

9549:   Output Parameters:
9550: + set - `PETSC_TRUE` if the matrix knows its Hermitian state (this tells you if the next flag is valid)
9551: - flg - the result (only valid if set is `PETSC_TRUE`)

9553:   Level: advanced

9555:   Notes:
9556:   Does not check the matrix values directly, so this may return unknown (set = `PETSC_FALSE`). Use `MatIsHermitian()`
9557:   if you want it explicitly checked

9559:   One can declare that a matrix is Hermitian with `MatSetOption`(mat,`MAT_HERMITIAN`,`PETSC_TRUE`) and if it is known to remain Hermitian
9560:   after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SYMMETRY_ETERNAL`,`PETSC_TRUE`)

9562: .seealso: [](ch_matrices), `Mat`, `MAT_SYMMETRY_ETERNAL`, `MAT_HERMITIAN`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetric()`
9563: @*/
9564: PetscErrorCode MatIsHermitianKnown(Mat A, PetscBool *set, PetscBool *flg)
9565: {
9566:   PetscFunctionBegin;
9568:   PetscAssertPointer(set, 2);
9569:   PetscAssertPointer(flg, 3);
9570:   if (A->hermitian != PETSC_BOOL3_UNKNOWN) {
9571:     *set = PETSC_TRUE;
9572:     *flg = PetscBool3ToBool(A->hermitian);
9573:   } else {
9574:     *set = PETSC_FALSE;
9575:   }
9576:   PetscFunctionReturn(PETSC_SUCCESS);
9577: }

9579: /*@
9580:   MatIsStructurallySymmetric - Test whether a matrix is structurally symmetric

9582:   Collective

9584:   Input Parameter:
9585: . A - the matrix to test

9587:   Output Parameter:
9588: . flg - the result

9590:   Level: intermediate

9592:   Notes:
9593:   If the matrix does yet know it is structurally symmetric this can be an expensive operation, also available `MatIsStructurallySymmetricKnown()`

9595:   One can declare that a matrix is structurally symmetric with `MatSetOption`(mat,`MAT_STRUCTURALLY_SYMMETRIC`,`PETSC_TRUE`) and if it is known to remain structurally
9596:   symmetric after changes to the matrices values one can call `MatSetOption`(mat,`MAT_STRUCTURAL_SYMMETRY_ETERNAL`,`PETSC_TRUE`)

9598: .seealso: [](ch_matrices), `Mat`, `MAT_STRUCTURALLY_SYMMETRIC`, `MAT_STRUCTURAL_SYMMETRY_ETERNAL`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsSymmetric()`, `MatSetOption()`, `MatIsStructurallySymmetricKnown()`
9599: @*/
9600: PetscErrorCode MatIsStructurallySymmetric(Mat A, PetscBool *flg)
9601: {
9602:   PetscFunctionBegin;
9604:   PetscAssertPointer(flg, 2);
9605:   if (A->structurally_symmetric != PETSC_BOOL3_UNKNOWN) {
9606:     *flg = PetscBool3ToBool(A->structurally_symmetric);
9607:   } else {
9608:     PetscUseTypeMethod(A, isstructurallysymmetric, flg);
9609:     PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, *flg));
9610:   }
9611:   PetscFunctionReturn(PETSC_SUCCESS);
9612: }

9614: /*@
9615:   MatIsStructurallySymmetricKnown - Checks if a matrix knows if it is structurally symmetric or not and its structurally symmetric state

9617:   Not Collective

9619:   Input Parameter:
9620: . A - the matrix to check

9622:   Output Parameters:
9623: + set - PETSC_TRUE if the matrix knows its structurally symmetric state (this tells you if the next flag is valid)
9624: - flg - the result (only valid if set is PETSC_TRUE)

9626:   Level: advanced

9628:   Notes:
9629:   One can declare that a matrix is structurally symmetric with `MatSetOption`(mat,`MAT_STRUCTURALLY_SYMMETRIC`,`PETSC_TRUE`) and if it is known to remain structurally
9630:   symmetric after changes to the matrices values one can call `MatSetOption`(mat,`MAT_STRUCTURAL_SYMMETRY_ETERNAL`,`PETSC_TRUE`)

9632:   Use `MatIsStructurallySymmetric()` to explicitly check if a matrix is structurally symmetric (this is an expensive operation)

9634: .seealso: [](ch_matrices), `Mat`, `MAT_STRUCTURALLY_SYMMETRIC`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetric()`, `MatIsHermitianKnown()`
9635: @*/
9636: PetscErrorCode MatIsStructurallySymmetricKnown(Mat A, PetscBool *set, PetscBool *flg)
9637: {
9638:   PetscFunctionBegin;
9640:   PetscAssertPointer(set, 2);
9641:   PetscAssertPointer(flg, 3);
9642:   if (A->structurally_symmetric != PETSC_BOOL3_UNKNOWN) {
9643:     *set = PETSC_TRUE;
9644:     *flg = PetscBool3ToBool(A->structurally_symmetric);
9645:   } else {
9646:     *set = PETSC_FALSE;
9647:   }
9648:   PetscFunctionReturn(PETSC_SUCCESS);
9649: }

9651: /*@
9652:   MatStashGetInfo - Gets how many values are currently in the matrix stash, i.e. need
9653:   to be communicated to other processors during the `MatAssemblyBegin()`/`MatAssemblyEnd()` process

9655:   Not Collective

9657:   Input Parameter:
9658: . mat - the matrix

9660:   Output Parameters:
9661: + nstash    - the size of the stash
9662: . reallocs  - the number of additional mallocs incurred.
9663: . bnstash   - the size of the block stash
9664: - breallocs - the number of additional mallocs incurred.in the block stash

9666:   Level: advanced

9668: .seealso: [](ch_matrices), `MatAssemblyBegin()`, `MatAssemblyEnd()`, `Mat`, `MatStashSetInitialSize()`
9669: @*/
9670: PetscErrorCode MatStashGetInfo(Mat mat, PetscInt *nstash, PetscInt *reallocs, PetscInt *bnstash, PetscInt *breallocs)
9671: {
9672:   PetscFunctionBegin;
9673:   PetscCall(MatStashGetInfo_Private(&mat->stash, nstash, reallocs));
9674:   PetscCall(MatStashGetInfo_Private(&mat->bstash, bnstash, breallocs));
9675:   PetscFunctionReturn(PETSC_SUCCESS);
9676: }

9678: /*@
9679:   MatCreateVecs - Get vector(s) compatible with the matrix, i.e. with the same
9680:   parallel layout, `PetscLayout` for rows and columns

9682:   Collective

9684:   Input Parameter:
9685: . mat - the matrix

9687:   Output Parameters:
9688: + right - (optional) vector that the matrix can be multiplied against
9689: - left  - (optional) vector that the matrix vector product can be stored in

9691:   Level: advanced

9693:   Notes:
9694:   The blocksize of the returned vectors is determined by the row and column block sizes set with `MatSetBlockSizes()` or the single blocksize (same for both) set by `MatSetBlockSize()`.

9696:   These are new vectors which are not owned by the mat, they should be destroyed in `VecDestroy()` when no longer needed

9698: .seealso: [](ch_matrices), `Mat`, `Vec`, `VecCreate()`, `VecDestroy()`, `DMCreateGlobalVector()`
9699: @*/
9700: PetscErrorCode MatCreateVecs(Mat mat, Vec *right, Vec *left)
9701: {
9702:   PetscFunctionBegin;
9705:   if (mat->ops->getvecs) {
9706:     PetscUseTypeMethod(mat, getvecs, right, left);
9707:   } else {
9708:     if (right) {
9709:       PetscCheck(mat->cmap->n >= 0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "PetscLayout for columns not yet setup");
9710:       PetscCall(VecCreateWithLayout_Private(mat->cmap, right));
9711:       PetscCall(VecSetType(*right, mat->defaultvectype));
9712: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
9713:       if (mat->boundtocpu && mat->bindingpropagates) {
9714:         PetscCall(VecSetBindingPropagates(*right, PETSC_TRUE));
9715:         PetscCall(VecBindToCPU(*right, PETSC_TRUE));
9716:       }
9717: #endif
9718:     }
9719:     if (left) {
9720:       PetscCheck(mat->rmap->n >= 0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "PetscLayout for rows not yet setup");
9721:       PetscCall(VecCreateWithLayout_Private(mat->rmap, left));
9722:       PetscCall(VecSetType(*left, mat->defaultvectype));
9723: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
9724:       if (mat->boundtocpu && mat->bindingpropagates) {
9725:         PetscCall(VecSetBindingPropagates(*left, PETSC_TRUE));
9726:         PetscCall(VecBindToCPU(*left, PETSC_TRUE));
9727:       }
9728: #endif
9729:     }
9730:   }
9731:   PetscFunctionReturn(PETSC_SUCCESS);
9732: }

9734: /*@
9735:   MatFactorInfoInitialize - Initializes a `MatFactorInfo` data structure
9736:   with default values.

9738:   Not Collective

9740:   Input Parameter:
9741: . info - the `MatFactorInfo` data structure

9743:   Level: developer

9745:   Notes:
9746:   The solvers are generally used through the `KSP` and `PC` objects, for example
9747:   `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC`

9749:   Once the data structure is initialized one may change certain entries as desired for the particular factorization to be performed

9751: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatFactorInfo`
9752: @*/
9753: PetscErrorCode MatFactorInfoInitialize(MatFactorInfo *info)
9754: {
9755:   PetscFunctionBegin;
9756:   PetscCall(PetscMemzero(info, sizeof(MatFactorInfo)));
9757:   PetscFunctionReturn(PETSC_SUCCESS);
9758: }

9760: /*@
9761:   MatFactorSetSchurIS - Set indices corresponding to the Schur complement you wish to have computed

9763:   Collective

9765:   Input Parameters:
9766: + mat - the factored matrix
9767: - is  - the index set defining the Schur indices (0-based)

9769:   Level: advanced

9771:   Notes:
9772:   Call `MatFactorSolveSchurComplement()` or `MatFactorSolveSchurComplementTranspose()` after this call to solve a Schur complement system.

9774:   You can call `MatFactorGetSchurComplement()` or `MatFactorCreateSchurComplement()` after this call.

9776:   This functionality is only supported for `MATSOLVERMUMPS` and `MATSOLVERMKL_PARDISO`

9778: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatFactorGetSchurComplement()`, `MatFactorRestoreSchurComplement()`, `MatFactorCreateSchurComplement()`, `MatFactorSolveSchurComplement()`,
9779:           `MatFactorSolveSchurComplementTranspose()`, `MATSOLVERMUMPS`, `MATSOLVERMKL_PARDISO`
9780: @*/
9781: PetscErrorCode MatFactorSetSchurIS(Mat mat, IS is)
9782: {
9783:   PetscErrorCode (*f)(Mat, IS);

9785:   PetscFunctionBegin;
9790:   PetscCheckSameComm(mat, 1, is, 2);
9791:   PetscCheck(mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix");
9792:   PetscCall(PetscObjectQueryFunction((PetscObject)mat, "MatFactorSetSchurIS_C", &f));
9793:   PetscCheck(f, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "The selected MatSolverType does not support Schur complement computation. You should use MATSOLVERMUMPS or MATSOLVERMKL_PARDISO");
9794:   PetscCall(MatDestroy(&mat->schur));
9795:   PetscCall((*f)(mat, is));
9796:   PetscCheck(mat->schur, PetscObjectComm((PetscObject)mat), PETSC_ERR_PLIB, "Schur complement has not been created");
9797:   PetscFunctionReturn(PETSC_SUCCESS);
9798: }

9800: /*@
9801:   MatFactorCreateSchurComplement - Create a Schur complement matrix object using Schur data computed during the factorization step

9803:   Logically Collective

9805:   Input Parameters:
9806: + F      - the factored matrix obtained by calling `MatGetFactor()`
9807: . S      - location where to return the Schur complement, can be `NULL`
9808: - status - the status of the Schur complement matrix, can be `NULL`

9810:   Level: advanced

9812:   Notes:
9813:   You must call `MatFactorSetSchurIS()` before calling this routine.

9815:   This functionality is only supported for `MATSOLVERMUMPS` and `MATSOLVERMKL_PARDISO`

9817:   The routine provides a copy of the Schur matrix stored within the solver data structures.
9818:   The caller must destroy the object when it is no longer needed.
9819:   If `MatFactorInvertSchurComplement()` has been called, the routine gets back the inverse.

9821:   Use `MatFactorGetSchurComplement()` to get access to the Schur complement matrix inside the factored matrix instead of making a copy of it (which this function does)

9823:   See `MatCreateSchurComplement()` or `MatGetSchurComplement()` for ways to create virtual or approximate Schur complements.

9825:   Developer Note:
9826:   The reason this routine exists is because the representation of the Schur complement within the factor matrix may be different than a standard PETSc
9827:   matrix representation and we normally do not want to use the time or memory to make a copy as a regular PETSc matrix.

9829: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorGetSchurComplement()`, `MatFactorSchurStatus`, `MATSOLVERMUMPS`, `MATSOLVERMKL_PARDISO`
9830: @*/
9831: PetscErrorCode MatFactorCreateSchurComplement(Mat F, Mat *S, MatFactorSchurStatus *status)
9832: {
9833:   PetscFunctionBegin;
9835:   if (S) PetscAssertPointer(S, 2);
9836:   if (status) PetscAssertPointer(status, 3);
9837:   if (S) {
9838:     PetscErrorCode (*f)(Mat, Mat *);

9840:     PetscCall(PetscObjectQueryFunction((PetscObject)F, "MatFactorCreateSchurComplement_C", &f));
9841:     if (f) {
9842:       PetscCall((*f)(F, S));
9843:     } else {
9844:       PetscCall(MatDuplicate(F->schur, MAT_COPY_VALUES, S));
9845:     }
9846:   }
9847:   if (status) *status = F->schur_status;
9848:   PetscFunctionReturn(PETSC_SUCCESS);
9849: }

9851: /*@
9852:   MatFactorGetSchurComplement - Gets access to a Schur complement matrix using the current Schur data within a factored matrix

9854:   Logically Collective

9856:   Input Parameters:
9857: + F      - the factored matrix obtained by calling `MatGetFactor()`
9858: . S      - location where to return the Schur complement, can be `NULL`
9859: - status - the status of the Schur complement matrix, can be `NULL`

9861:   Level: advanced

9863:   Notes:
9864:   You must call `MatFactorSetSchurIS()` before calling this routine.

9866:   Schur complement mode is currently implemented for sequential matrices with factor type of `MATSOLVERMUMPS`

9868:   The routine returns a the Schur Complement stored within the data structures of the solver.

9870:   If `MatFactorInvertSchurComplement()` has previously been called, the returned matrix is actually the inverse of the Schur complement.

9872:   The returned matrix should not be destroyed; the caller should call `MatFactorRestoreSchurComplement()` when the object is no longer needed.

9874:   Use `MatFactorCreateSchurComplement()` to create a copy of the Schur complement matrix that is within a factored matrix

9876:   See `MatCreateSchurComplement()` or `MatGetSchurComplement()` for ways to create virtual or approximate Schur complements.

9878: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorRestoreSchurComplement()`, `MatFactorCreateSchurComplement()`, `MatFactorSchurStatus`
9879: @*/
9880: PetscErrorCode MatFactorGetSchurComplement(Mat F, Mat *S, MatFactorSchurStatus *status)
9881: {
9882:   PetscFunctionBegin;
9884:   if (S) {
9885:     PetscAssertPointer(S, 2);
9886:     *S = F->schur;
9887:   }
9888:   if (status) {
9889:     PetscAssertPointer(status, 3);
9890:     *status = F->schur_status;
9891:   }
9892:   PetscFunctionReturn(PETSC_SUCCESS);
9893: }

9895: static PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat F)
9896: {
9897:   Mat S = F->schur;

9899:   PetscFunctionBegin;
9900:   switch (F->schur_status) {
9901:   case MAT_FACTOR_SCHUR_UNFACTORED: // fall-through
9902:   case MAT_FACTOR_SCHUR_INVERTED:
9903:     if (S) {
9904:       S->ops->solve             = NULL;
9905:       S->ops->matsolve          = NULL;
9906:       S->ops->solvetranspose    = NULL;
9907:       S->ops->matsolvetranspose = NULL;
9908:       S->ops->solveadd          = NULL;
9909:       S->ops->solvetransposeadd = NULL;
9910:       S->factortype             = MAT_FACTOR_NONE;
9911:       PetscCall(PetscFree(S->solvertype));
9912:     }
9913:   case MAT_FACTOR_SCHUR_FACTORED: // fall-through
9914:     break;
9915:   default:
9916:     SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unhandled MatFactorSchurStatus %d", F->schur_status);
9917:   }
9918:   PetscFunctionReturn(PETSC_SUCCESS);
9919: }

9921: /*@
9922:   MatFactorRestoreSchurComplement - Restore the Schur complement matrix object obtained from a call to `MatFactorGetSchurComplement()`

9924:   Logically Collective

9926:   Input Parameters:
9927: + F      - the factored matrix obtained by calling `MatGetFactor()`
9928: . S      - location where the Schur complement is stored
9929: - status - the status of the Schur complement matrix (see `MatFactorSchurStatus`)

9931:   Level: advanced

9933: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorCreateSchurComplement()`, `MatFactorSchurStatus`
9934: @*/
9935: PetscErrorCode MatFactorRestoreSchurComplement(Mat F, Mat *S, MatFactorSchurStatus status)
9936: {
9937:   PetscFunctionBegin;
9939:   if (S) {
9941:     *S = NULL;
9942:   }
9943:   F->schur_status = status;
9944:   PetscCall(MatFactorUpdateSchurStatus_Private(F));
9945:   PetscFunctionReturn(PETSC_SUCCESS);
9946: }

9948: /*@
9949:   MatFactorSolveSchurComplementTranspose - Solve the transpose of the Schur complement system computed during the factorization step

9951:   Logically Collective

9953:   Input Parameters:
9954: + F   - the factored matrix obtained by calling `MatGetFactor()`
9955: . rhs - location where the right-hand side of the Schur complement system is stored
9956: - sol - location where the solution of the Schur complement system has to be returned

9958:   Level: advanced

9960:   Notes:
9961:   The sizes of the vectors should match the size of the Schur complement

9963:   Must be called after `MatFactorSetSchurIS()`

9965: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorSolveSchurComplement()`
9966: @*/
9967: PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat F, Vec rhs, Vec sol)
9968: {
9969:   PetscFunctionBegin;
9976:   PetscCheckSameComm(F, 1, rhs, 2);
9977:   PetscCheckSameComm(F, 1, sol, 3);
9978:   PetscCall(MatFactorFactorizeSchurComplement(F));
9979:   switch (F->schur_status) {
9980:   case MAT_FACTOR_SCHUR_FACTORED:
9981:     PetscCall(MatSolveTranspose(F->schur, rhs, sol));
9982:     break;
9983:   case MAT_FACTOR_SCHUR_INVERTED:
9984:     PetscCall(MatMultTranspose(F->schur, rhs, sol));
9985:     break;
9986:   default:
9987:     SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unhandled MatFactorSchurStatus %d", F->schur_status);
9988:   }
9989:   PetscFunctionReturn(PETSC_SUCCESS);
9990: }

9992: /*@
9993:   MatFactorSolveSchurComplement - Solve the Schur complement system computed during the factorization step

9995:   Logically Collective

9997:   Input Parameters:
9998: + F   - the factored matrix obtained by calling `MatGetFactor()`
9999: . rhs - location where the right-hand side of the Schur complement system is stored
10000: - sol - location where the solution of the Schur complement system has to be returned

10002:   Level: advanced

10004:   Notes:
10005:   The sizes of the vectors should match the size of the Schur complement

10007:   Must be called after `MatFactorSetSchurIS()`

10009: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorSolveSchurComplementTranspose()`
10010: @*/
10011: PetscErrorCode MatFactorSolveSchurComplement(Mat F, Vec rhs, Vec sol)
10012: {
10013:   PetscFunctionBegin;
10020:   PetscCheckSameComm(F, 1, rhs, 2);
10021:   PetscCheckSameComm(F, 1, sol, 3);
10022:   PetscCall(MatFactorFactorizeSchurComplement(F));
10023:   switch (F->schur_status) {
10024:   case MAT_FACTOR_SCHUR_FACTORED:
10025:     PetscCall(MatSolve(F->schur, rhs, sol));
10026:     break;
10027:   case MAT_FACTOR_SCHUR_INVERTED:
10028:     PetscCall(MatMult(F->schur, rhs, sol));
10029:     break;
10030:   default:
10031:     SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unhandled MatFactorSchurStatus %d", F->schur_status);
10032:   }
10033:   PetscFunctionReturn(PETSC_SUCCESS);
10034: }

10036: PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat);
10037: #if PetscDefined(HAVE_CUDA)
10038: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Internal(Mat);
10039: #endif

10041: /* Schur status updated in the interface */
10042: static PetscErrorCode MatFactorInvertSchurComplement_Private(Mat F)
10043: {
10044:   Mat S = F->schur;

10046:   PetscFunctionBegin;
10047:   if (S) {
10048:     PetscMPIInt size;
10049:     PetscBool   isdense, isdensecuda;

10051:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)S), &size));
10052:     PetscCheck(size <= 1, PetscObjectComm((PetscObject)S), PETSC_ERR_SUP, "Not yet implemented");
10053:     PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSEQDENSE, &isdense));
10054:     PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSEQDENSECUDA, &isdensecuda));
10055:     PetscCheck(isdense || isdensecuda, PetscObjectComm((PetscObject)S), PETSC_ERR_SUP, "Not implemented for type %s", ((PetscObject)S)->type_name);
10056:     PetscCall(PetscLogEventBegin(MAT_FactorInvS, F, 0, 0, 0));
10057:     if (isdense) {
10058:       PetscCall(MatSeqDenseInvertFactors_Private(S));
10059:     } else if (isdensecuda) {
10060: #if defined(PETSC_HAVE_CUDA)
10061:       PetscCall(MatSeqDenseCUDAInvertFactors_Internal(S));
10062: #endif
10063:     }
10064:     // HIP??????????????
10065:     PetscCall(PetscLogEventEnd(MAT_FactorInvS, F, 0, 0, 0));
10066:   }
10067:   PetscFunctionReturn(PETSC_SUCCESS);
10068: }

10070: /*@
10071:   MatFactorInvertSchurComplement - Invert the Schur complement matrix computed during the factorization step

10073:   Logically Collective

10075:   Input Parameter:
10076: . F - the factored matrix obtained by calling `MatGetFactor()`

10078:   Level: advanced

10080:   Notes:
10081:   Must be called after `MatFactorSetSchurIS()`.

10083:   Call `MatFactorGetSchurComplement()` or  `MatFactorCreateSchurComplement()` AFTER this call to actually compute the inverse and get access to it.

10085: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorGetSchurComplement()`, `MatFactorCreateSchurComplement()`
10086: @*/
10087: PetscErrorCode MatFactorInvertSchurComplement(Mat F)
10088: {
10089:   PetscFunctionBegin;
10092:   if (F->schur_status == MAT_FACTOR_SCHUR_INVERTED) PetscFunctionReturn(PETSC_SUCCESS);
10093:   PetscCall(MatFactorFactorizeSchurComplement(F));
10094:   PetscCall(MatFactorInvertSchurComplement_Private(F));
10095:   F->schur_status = MAT_FACTOR_SCHUR_INVERTED;
10096:   PetscFunctionReturn(PETSC_SUCCESS);
10097: }

10099: /*@
10100:   MatFactorFactorizeSchurComplement - Factorize the Schur complement matrix computed during the factorization step

10102:   Logically Collective

10104:   Input Parameter:
10105: . F - the factored matrix obtained by calling `MatGetFactor()`

10107:   Level: advanced

10109:   Note:
10110:   Must be called after `MatFactorSetSchurIS()`

10112: .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorInvertSchurComplement()`
10113: @*/
10114: PetscErrorCode MatFactorFactorizeSchurComplement(Mat F)
10115: {
10116:   MatFactorInfo info;

10118:   PetscFunctionBegin;
10121:   if (F->schur_status == MAT_FACTOR_SCHUR_INVERTED || F->schur_status == MAT_FACTOR_SCHUR_FACTORED) PetscFunctionReturn(PETSC_SUCCESS);
10122:   PetscCall(PetscLogEventBegin(MAT_FactorFactS, F, 0, 0, 0));
10123:   PetscCall(PetscMemzero(&info, sizeof(MatFactorInfo)));
10124:   if (F->factortype == MAT_FACTOR_CHOLESKY) { /* LDL^t regarded as Cholesky */
10125:     PetscCall(MatCholeskyFactor(F->schur, NULL, &info));
10126:   } else {
10127:     PetscCall(MatLUFactor(F->schur, NULL, NULL, &info));
10128:   }
10129:   PetscCall(PetscLogEventEnd(MAT_FactorFactS, F, 0, 0, 0));
10130:   F->schur_status = MAT_FACTOR_SCHUR_FACTORED;
10131:   PetscFunctionReturn(PETSC_SUCCESS);
10132: }

10134: /*@
10135:   MatPtAP - Creates the matrix product $C = P^T * A * P$

10137:   Neighbor-wise Collective

10139:   Input Parameters:
10140: + A     - the matrix
10141: . P     - the projection matrix
10142: . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
10143: - fill  - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)), use `PETSC_DETERMINE` or `PETSC_CURRENT` if you do not have a good estimate
10144:           if the result is a dense matrix this is irrelevant

10146:   Output Parameter:
10147: . C - the product matrix

10149:   Level: intermediate

10151:   Notes:
10152:   C will be created and must be destroyed by the user with `MatDestroy()`.

10154:   An alternative approach to this function is to use `MatProductCreate()` and set the desired options before the computation is done

10156:   The deprecated `PETSC_DEFAULT` in `fill` also means use the current value

10158:   Developer Note:
10159:   For matrix types without special implementation the function fallbacks to `MatMatMult()` followed by `MatTransposeMatMult()`.

10161: .seealso: [](ch_matrices), `Mat`, `MatProductCreate()`, `MatMatMult()`, `MatRARt()`
10162: @*/
10163: PetscErrorCode MatPtAP(Mat A, Mat P, MatReuse scall, PetscReal fill, Mat *C)
10164: {
10165:   PetscFunctionBegin;
10166:   if (scall == MAT_REUSE_MATRIX) MatCheckProduct(*C, 5);
10167:   PetscCheck(scall != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Inplace product not supported");

10169:   if (scall == MAT_INITIAL_MATRIX) {
10170:     PetscCall(MatProductCreate(A, P, NULL, C));
10171:     PetscCall(MatProductSetType(*C, MATPRODUCT_PtAP));
10172:     PetscCall(MatProductSetAlgorithm(*C, "default"));
10173:     PetscCall(MatProductSetFill(*C, fill));

10175:     (*C)->product->api_user = PETSC_TRUE;
10176:     PetscCall(MatProductSetFromOptions(*C));
10177:     PetscCheck((*C)->ops->productsymbolic, PetscObjectComm((PetscObject)*C), PETSC_ERR_SUP, "MatProduct %s not supported for A %s and P %s", MatProductTypes[MATPRODUCT_PtAP], ((PetscObject)A)->type_name, ((PetscObject)P)->type_name);
10178:     PetscCall(MatProductSymbolic(*C));
10179:   } else { /* scall == MAT_REUSE_MATRIX */
10180:     PetscCall(MatProductReplaceMats(A, P, NULL, *C));
10181:   }

10183:   PetscCall(MatProductNumeric(*C));
10184:   (*C)->symmetric = A->symmetric;
10185:   (*C)->spd       = A->spd;
10186:   PetscFunctionReturn(PETSC_SUCCESS);
10187: }

10189: /*@
10190:   MatRARt - Creates the matrix product $C = R * A * R^T$

10192:   Neighbor-wise Collective

10194:   Input Parameters:
10195: + A     - the matrix
10196: . R     - the projection matrix
10197: . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
10198: - fill  - expected fill as ratio of nnz(C)/nnz(A), use `PETSC_DETERMINE` or `PETSC_CURRENT` if you do not have a good estimate
10199:           if the result is a dense matrix this is irrelevant

10201:   Output Parameter:
10202: . C - the product matrix

10204:   Level: intermediate

10206:   Notes:
10207:   `C` will be created and must be destroyed by the user with `MatDestroy()`.

10209:   An alternative approach to this function is to use `MatProductCreate()` and set the desired options before the computation is done

10211:   This routine is currently only implemented for pairs of `MATAIJ` matrices and classes
10212:   which inherit from `MATAIJ`. Due to PETSc sparse matrix block row distribution among processes,
10213:   the parallel `MatRARt()` is implemented computing the explicit transpose of `R`, which can be very expensive.
10214:   We recommend using `MatPtAP()` when possible.

10216:   The deprecated `PETSC_DEFAULT` in `fill` also means use the current value

10218: .seealso: [](ch_matrices), `Mat`, `MatProductCreate()`, `MatMatMult()`, `MatPtAP()`
10219: @*/
10220: PetscErrorCode MatRARt(Mat A, Mat R, MatReuse scall, PetscReal fill, Mat *C)
10221: {
10222:   PetscFunctionBegin;
10223:   if (scall == MAT_REUSE_MATRIX) MatCheckProduct(*C, 5);
10224:   PetscCheck(scall != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Inplace product not supported");

10226:   if (scall == MAT_INITIAL_MATRIX) {
10227:     PetscCall(MatProductCreate(A, R, NULL, C));
10228:     PetscCall(MatProductSetType(*C, MATPRODUCT_RARt));
10229:     PetscCall(MatProductSetAlgorithm(*C, "default"));
10230:     PetscCall(MatProductSetFill(*C, fill));

10232:     (*C)->product->api_user = PETSC_TRUE;
10233:     PetscCall(MatProductSetFromOptions(*C));
10234:     PetscCheck((*C)->ops->productsymbolic, PetscObjectComm((PetscObject)*C), PETSC_ERR_SUP, "MatProduct %s not supported for A %s and R %s", MatProductTypes[MATPRODUCT_RARt], ((PetscObject)A)->type_name, ((PetscObject)R)->type_name);
10235:     PetscCall(MatProductSymbolic(*C));
10236:   } else { /* scall == MAT_REUSE_MATRIX */
10237:     PetscCall(MatProductReplaceMats(A, R, NULL, *C));
10238:   }

10240:   PetscCall(MatProductNumeric(*C));
10241:   if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetOption(*C, MAT_SYMMETRIC, PETSC_TRUE));
10242:   PetscFunctionReturn(PETSC_SUCCESS);
10243: }

10245: static PetscErrorCode MatProduct_Private(Mat A, Mat B, MatReuse scall, PetscReal fill, MatProductType ptype, Mat *C)
10246: {
10247:   PetscBool flg = PETSC_TRUE;

10249:   PetscFunctionBegin;
10250:   PetscCheck(scall != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MAT_INPLACE_MATRIX product not supported");
10251:   if (scall == MAT_INITIAL_MATRIX) {
10252:     PetscCall(PetscInfo(A, "Calling MatProduct API with MAT_INITIAL_MATRIX and product type %s\n", MatProductTypes[ptype]));
10253:     PetscCall(MatProductCreate(A, B, NULL, C));
10254:     PetscCall(MatProductSetAlgorithm(*C, MATPRODUCTALGORITHMDEFAULT));
10255:     PetscCall(MatProductSetFill(*C, fill));
10256:   } else { /* scall == MAT_REUSE_MATRIX */
10257:     Mat_Product *product = (*C)->product;

10259:     PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)*C, &flg, MATSEQDENSE, MATMPIDENSE, ""));
10260:     if (flg && product && product->type != ptype) {
10261:       PetscCall(MatProductClear(*C));
10262:       product = NULL;
10263:     }
10264:     PetscCall(PetscInfo(A, "Calling MatProduct API with MAT_REUSE_MATRIX %s product present and product type %s\n", product ? "with" : "without", MatProductTypes[ptype]));
10265:     if (!product) { /* user provide the dense matrix *C without calling MatProductCreate() or reusing it from previous calls */
10266:       PetscCheck(flg, PetscObjectComm((PetscObject)*C), PETSC_ERR_SUP, "Call MatProductCreate() first");
10267:       PetscCall(MatProductCreate_Private(A, B, NULL, *C));
10268:       product        = (*C)->product;
10269:       product->fill  = fill;
10270:       product->clear = PETSC_TRUE;
10271:     } else { /* user may change input matrices A or B when MAT_REUSE_MATRIX */
10272:       flg = PETSC_FALSE;
10273:       PetscCall(MatProductReplaceMats(A, B, NULL, *C));
10274:     }
10275:   }
10276:   if (flg) {
10277:     (*C)->product->api_user = PETSC_TRUE;
10278:     PetscCall(MatProductSetType(*C, ptype));
10279:     PetscCall(MatProductSetFromOptions(*C));
10280:     PetscCall(MatProductSymbolic(*C));
10281:   }
10282:   PetscCall(MatProductNumeric(*C));
10283:   PetscFunctionReturn(PETSC_SUCCESS);
10284: }

10286: /*@
10287:   MatMatMult - Performs matrix-matrix multiplication C=A*B.

10289:   Neighbor-wise Collective

10291:   Input Parameters:
10292: + A     - the left matrix
10293: . B     - the right matrix
10294: . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
10295: - fill  - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use `PETSC_DETERMINE` or `PETSC_CURRENT` if you do not have a good estimate
10296:           if the result is a dense matrix this is irrelevant

10298:   Output Parameter:
10299: . C - the product matrix

10301:   Notes:
10302:   Unless scall is `MAT_REUSE_MATRIX` C will be created.

10304:   `MAT_REUSE_MATRIX` can only be used if the matrices A and B have the same nonzero pattern as in the previous call and C was obtained from a previous
10305:   call to this function with `MAT_INITIAL_MATRIX`.

10307:   To determine the correct fill value, run with `-info` and search for the string "Fill ratio" to see the value actually needed.

10309:   In the special case where matrix `B` (and hence `C`) are dense you can create the correctly sized matrix `C` yourself and then call this routine with `MAT_REUSE_MATRIX`,
10310:   rather than first having `MatMatMult()` create it for you. You can NEVER do this if the matrix `C` is sparse.

10312:   The deprecated `PETSC_DEFAULT` in `fill` also means use the current value

10314:   Example of Usage:
10315: .vb
10316:      MatProductCreate(A,B,NULL,&C);
10317:      MatProductSetType(C,MATPRODUCT_AB);
10318:      MatProductSymbolic(C);
10319:      MatProductNumeric(C); // compute C=A * B
10320:      MatProductReplaceMats(A1,B1,NULL,C); // compute C=A1 * B1
10321:      MatProductNumeric(C);
10322:      MatProductReplaceMats(A2,NULL,NULL,C); // compute C=A2 * B1
10323:      MatProductNumeric(C);
10324: .ve

10326:   Level: intermediate

10328: .seealso: [](ch_matrices), `Mat`, `MatProductType`, `MATPRODUCT_AB`, `MatTransposeMatMult()`, `MatMatTransposeMult()`, `MatPtAP()`, `MatProductCreate()`, `MatProductSymbolic()`, `MatProductReplaceMats()`, `MatProductNumeric()`
10329: @*/
10330: PetscErrorCode MatMatMult(Mat A, Mat B, MatReuse scall, PetscReal fill, Mat *C)
10331: {
10332:   PetscFunctionBegin;
10333:   PetscCall(MatProduct_Private(A, B, scall, fill, MATPRODUCT_AB, C));
10334:   PetscFunctionReturn(PETSC_SUCCESS);
10335: }

10337: /*@
10338:   MatMatTransposeMult - Performs matrix-matrix multiplication $C = A*B^T$.

10340:   Neighbor-wise Collective

10342:   Input Parameters:
10343: + A     - the left matrix
10344: . B     - the right matrix
10345: . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
10346: - fill  - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use `PETSC_DETERMINE` or `PETSC_CURRENT` if not known

10348:   Output Parameter:
10349: . C - the product matrix

10351:   Options Database Key:
10352: . -matmattransmult_mpidense_mpidense_via {allgatherv,cyclic} - Choose between algorithms for `MATMPIDENSE` matrices: the
10353:               first redundantly copies the transposed `B` matrix on each process and requires O(log P) communication complexity;
10354:               the second never stores more than one portion of the `B` matrix at a time but requires O(P) communication complexity.

10356:   Level: intermediate

10358:   Notes:
10359:   C will be created if `MAT_INITIAL_MATRIX` and must be destroyed by the user with `MatDestroy()`.

10361:   `MAT_REUSE_MATRIX` can only be used if the matrices A and B have the same nonzero pattern as in the previous call

10363:   To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
10364:   actually needed.

10366:   This routine is currently only implemented for pairs of `MATSEQAIJ` matrices, for the `MATSEQDENSE` class,
10367:   and for pairs of `MATMPIDENSE` matrices.

10369:   This routine is shorthand for using `MatProductCreate()` with the `MatProductType` of `MATPRODUCT_ABt`

10371:   The deprecated `PETSC_DEFAULT` in `fill` also means use the current value

10373: .seealso: [](ch_matrices), `Mat`, `MatProductCreate()`, `MATPRODUCT_ABt`, `MatMatMult()`, `MatTransposeMatMult()` `MatPtAP()`, `MatProductAlgorithm`, `MatProductType`
10374: @*/
10375: PetscErrorCode MatMatTransposeMult(Mat A, Mat B, MatReuse scall, PetscReal fill, Mat *C)
10376: {
10377:   PetscFunctionBegin;
10378:   PetscCall(MatProduct_Private(A, B, scall, fill, MATPRODUCT_ABt, C));
10379:   if (A == B) PetscCall(MatSetOption(*C, MAT_SYMMETRIC, PETSC_TRUE));
10380:   PetscFunctionReturn(PETSC_SUCCESS);
10381: }

10383: /*@
10384:   MatTransposeMatMult - Performs matrix-matrix multiplication $C = A^T*B$.

10386:   Neighbor-wise Collective

10388:   Input Parameters:
10389: + A     - the left matrix
10390: . B     - the right matrix
10391: . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
10392: - fill  - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use `PETSC_DETERMINE` or `PETSC_CURRENT` if not known

10394:   Output Parameter:
10395: . C - the product matrix

10397:   Level: intermediate

10399:   Notes:
10400:   `C` will be created if `MAT_INITIAL_MATRIX` and must be destroyed by the user with `MatDestroy()`.

10402:   `MAT_REUSE_MATRIX` can only be used if the matrices A and B have the same nonzero pattern as in the previous call.

10404:   This routine is shorthand for using `MatProductCreate()` with the `MatProductType` of `MATPRODUCT_AtB`

10406:   To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value
10407:   actually needed.

10409:   This routine is currently implemented for pairs of `MATAIJ` matrices and pairs of `MATSEQDENSE` matrices and classes
10410:   which inherit from `MATSEQAIJ`.  `C` will be of the same type as the input matrices.

10412:   The deprecated `PETSC_DEFAULT` in `fill` also means use the current value

10414: .seealso: [](ch_matrices), `Mat`, `MatProductCreate()`, `MATPRODUCT_AtB`, `MatMatMult()`, `MatMatTransposeMult()`, `MatPtAP()`
10415: @*/
10416: PetscErrorCode MatTransposeMatMult(Mat A, Mat B, MatReuse scall, PetscReal fill, Mat *C)
10417: {
10418:   PetscFunctionBegin;
10419:   PetscCall(MatProduct_Private(A, B, scall, fill, MATPRODUCT_AtB, C));
10420:   PetscFunctionReturn(PETSC_SUCCESS);
10421: }

10423: /*@
10424:   MatMatMatMult - Performs matrix-matrix-matrix multiplication D=A*B*C.

10426:   Neighbor-wise Collective

10428:   Input Parameters:
10429: + A     - the left matrix
10430: . B     - the middle matrix
10431: . C     - the right matrix
10432: . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
10433: - fill  - expected fill as ratio of nnz(D)/(nnz(A) + nnz(B)+nnz(C)), use `PETSC_DETERMINE` or `PETSC_CURRENT` if you do not have a good estimate
10434:           if the result is a dense matrix this is irrelevant

10436:   Output Parameter:
10437: . D - the product matrix

10439:   Level: intermediate

10441:   Notes:
10442:   Unless `scall` is `MAT_REUSE_MATRIX` `D` will be created.

10444:   `MAT_REUSE_MATRIX` can only be used if the matrices `A`, `B`, and `C` have the same nonzero pattern as in the previous call

10446:   This routine is shorthand for using `MatProductCreate()` with the `MatProductType` of `MATPRODUCT_ABC`

10448:   To determine the correct fill value, run with `-info` and search for the string "Fill ratio" to see the value
10449:   actually needed.

10451:   If you have many matrices with the same non-zero structure to multiply, you
10452:   should use `MAT_REUSE_MATRIX` in all calls but the first

10454:   The deprecated `PETSC_DEFAULT` in `fill` also means use the current value

10456: .seealso: [](ch_matrices), `Mat`, `MatProductCreate()`, `MATPRODUCT_ABC`, `MatMatMult`, `MatPtAP()`, `MatMatTransposeMult()`, `MatTransposeMatMult()`
10457: @*/
10458: PetscErrorCode MatMatMatMult(Mat A, Mat B, Mat C, MatReuse scall, PetscReal fill, Mat *D)
10459: {
10460:   PetscFunctionBegin;
10461:   if (scall == MAT_REUSE_MATRIX) MatCheckProduct(*D, 6);
10462:   PetscCheck(scall != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Inplace product not supported");

10464:   if (scall == MAT_INITIAL_MATRIX) {
10465:     PetscCall(MatProductCreate(A, B, C, D));
10466:     PetscCall(MatProductSetType(*D, MATPRODUCT_ABC));
10467:     PetscCall(MatProductSetAlgorithm(*D, "default"));
10468:     PetscCall(MatProductSetFill(*D, fill));

10470:     (*D)->product->api_user = PETSC_TRUE;
10471:     PetscCall(MatProductSetFromOptions(*D));
10472:     PetscCheck((*D)->ops->productsymbolic, PetscObjectComm((PetscObject)*D), PETSC_ERR_SUP, "MatProduct %s not supported for A %s, B %s and C %s", MatProductTypes[MATPRODUCT_ABC], ((PetscObject)A)->type_name, ((PetscObject)B)->type_name,
10473:                ((PetscObject)C)->type_name);
10474:     PetscCall(MatProductSymbolic(*D));
10475:   } else { /* user may change input matrices when REUSE */
10476:     PetscCall(MatProductReplaceMats(A, B, C, *D));
10477:   }
10478:   PetscCall(MatProductNumeric(*D));
10479:   PetscFunctionReturn(PETSC_SUCCESS);
10480: }

10482: /*@
10483:   MatCreateRedundantMatrix - Create redundant matrices and put them into processors of subcommunicators.

10485:   Collective

10487:   Input Parameters:
10488: + mat      - the matrix
10489: . nsubcomm - the number of subcommunicators (= number of redundant parallel or sequential matrices)
10490: . subcomm  - MPI communicator split from the communicator where mat resides in (or `MPI_COMM_NULL` if nsubcomm is used)
10491: - reuse    - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`

10493:   Output Parameter:
10494: . matredundant - redundant matrix

10496:   Level: advanced

10498:   Notes:
10499:   `MAT_REUSE_MATRIX` can only be used when the nonzero structure of the
10500:   original matrix has not changed from that last call to `MatCreateRedundantMatrix()`.

10502:   This routine creates the duplicated matrices in the subcommunicators; you should NOT create them before
10503:   calling it.

10505:   `PetscSubcommCreate()` can be used to manage the creation of the subcomm but need not be.

10507: .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `PetscSubcommCreate()`, `PetscSubcomm`
10508: @*/
10509: PetscErrorCode MatCreateRedundantMatrix(Mat mat, PetscInt nsubcomm, MPI_Comm subcomm, MatReuse reuse, Mat *matredundant)
10510: {
10511:   MPI_Comm       comm;
10512:   PetscMPIInt    size;
10513:   PetscInt       mloc_sub, nloc_sub, rstart, rend, M = mat->rmap->N, N = mat->cmap->N, bs = mat->rmap->bs;
10514:   Mat_Redundant *redund     = NULL;
10515:   PetscSubcomm   psubcomm   = NULL;
10516:   MPI_Comm       subcomm_in = subcomm;
10517:   Mat           *matseq;
10518:   IS             isrow, iscol;
10519:   PetscBool      newsubcomm = PETSC_FALSE;

10521:   PetscFunctionBegin;
10523:   if (nsubcomm && reuse == MAT_REUSE_MATRIX) {
10524:     PetscAssertPointer(*matredundant, 5);
10526:   }

10528:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
10529:   if (size == 1 || nsubcomm == 1) {
10530:     if (reuse == MAT_INITIAL_MATRIX) {
10531:       PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, matredundant));
10532:     } else {
10533:       PetscCheck(*matredundant != mat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix");
10534:       PetscCall(MatCopy(mat, *matredundant, SAME_NONZERO_PATTERN));
10535:     }
10536:     PetscFunctionReturn(PETSC_SUCCESS);
10537:   }

10539:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
10540:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
10541:   MatCheckPreallocated(mat, 1);

10543:   PetscCall(PetscLogEventBegin(MAT_RedundantMat, mat, 0, 0, 0));
10544:   if (subcomm_in == MPI_COMM_NULL && reuse == MAT_INITIAL_MATRIX) { /* get subcomm if user does not provide subcomm */
10545:     /* create psubcomm, then get subcomm */
10546:     PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
10547:     PetscCallMPI(MPI_Comm_size(comm, &size));
10548:     PetscCheck(nsubcomm >= 1 && nsubcomm <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "nsubcomm must between 1 and %d", size);

10550:     PetscCall(PetscSubcommCreate(comm, &psubcomm));
10551:     PetscCall(PetscSubcommSetNumber(psubcomm, nsubcomm));
10552:     PetscCall(PetscSubcommSetType(psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
10553:     PetscCall(PetscSubcommSetFromOptions(psubcomm));
10554:     PetscCall(PetscCommDuplicate(PetscSubcommChild(psubcomm), &subcomm, NULL));
10555:     newsubcomm = PETSC_TRUE;
10556:     PetscCall(PetscSubcommDestroy(&psubcomm));
10557:   }

10559:   /* get isrow, iscol and a local sequential matrix matseq[0] */
10560:   if (reuse == MAT_INITIAL_MATRIX) {
10561:     mloc_sub = PETSC_DECIDE;
10562:     nloc_sub = PETSC_DECIDE;
10563:     if (bs < 1) {
10564:       PetscCall(PetscSplitOwnership(subcomm, &mloc_sub, &M));
10565:       PetscCall(PetscSplitOwnership(subcomm, &nloc_sub, &N));
10566:     } else {
10567:       PetscCall(PetscSplitOwnershipBlock(subcomm, bs, &mloc_sub, &M));
10568:       PetscCall(PetscSplitOwnershipBlock(subcomm, bs, &nloc_sub, &N));
10569:     }
10570:     PetscCallMPI(MPI_Scan(&mloc_sub, &rend, 1, MPIU_INT, MPI_SUM, subcomm));
10571:     rstart = rend - mloc_sub;
10572:     PetscCall(ISCreateStride(PETSC_COMM_SELF, mloc_sub, rstart, 1, &isrow));
10573:     PetscCall(ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &iscol));
10574:     PetscCall(ISSetIdentity(iscol));
10575:   } else { /* reuse == MAT_REUSE_MATRIX */
10576:     PetscCheck(*matredundant != mat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix");
10577:     /* retrieve subcomm */
10578:     PetscCall(PetscObjectGetComm((PetscObject)*matredundant, &subcomm));
10579:     redund = (*matredundant)->redundant;
10580:     isrow  = redund->isrow;
10581:     iscol  = redund->iscol;
10582:     matseq = redund->matseq;
10583:   }
10584:   PetscCall(MatCreateSubMatrices(mat, 1, &isrow, &iscol, reuse, &matseq));

10586:   /* get matredundant over subcomm */
10587:   if (reuse == MAT_INITIAL_MATRIX) {
10588:     PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm, matseq[0], nloc_sub, reuse, matredundant));

10590:     /* create a supporting struct and attach it to C for reuse */
10591:     PetscCall(PetscNew(&redund));
10592:     (*matredundant)->redundant = redund;
10593:     redund->isrow              = isrow;
10594:     redund->iscol              = iscol;
10595:     redund->matseq             = matseq;
10596:     if (newsubcomm) {
10597:       redund->subcomm = subcomm;
10598:     } else {
10599:       redund->subcomm = MPI_COMM_NULL;
10600:     }
10601:   } else {
10602:     PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm, matseq[0], PETSC_DECIDE, reuse, matredundant));
10603:   }
10604: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
10605:   if (matseq[0]->boundtocpu && matseq[0]->bindingpropagates) {
10606:     PetscCall(MatBindToCPU(*matredundant, PETSC_TRUE));
10607:     PetscCall(MatSetBindingPropagates(*matredundant, PETSC_TRUE));
10608:   }
10609: #endif
10610:   PetscCall(PetscLogEventEnd(MAT_RedundantMat, mat, 0, 0, 0));
10611:   PetscFunctionReturn(PETSC_SUCCESS);
10612: }

10614: /*@C
10615:   MatGetMultiProcBlock - Create multiple 'parallel submatrices' from
10616:   a given `Mat`. Each submatrix can span multiple procs.

10618:   Collective

10620:   Input Parameters:
10621: + mat     - the matrix
10622: . subComm - the sub communicator obtained as if by `MPI_Comm_split(PetscObjectComm((PetscObject)mat))`
10623: - scall   - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`

10625:   Output Parameter:
10626: . subMat - parallel sub-matrices each spanning a given `subcomm`

10628:   Level: advanced

10630:   Notes:
10631:   The submatrix partition across processors is dictated by `subComm` a
10632:   communicator obtained by `MPI_comm_split()` or via `PetscSubcommCreate()`. The `subComm`
10633:   is not restricted to be grouped with consecutive original MPI processes.

10635:   Due the `MPI_Comm_split()` usage, the parallel layout of the submatrices
10636:   map directly to the layout of the original matrix [wrt the local
10637:   row,col partitioning]. So the original 'DiagonalMat' naturally maps
10638:   into the 'DiagonalMat' of the `subMat`, hence it is used directly from
10639:   the `subMat`. However the offDiagMat looses some columns - and this is
10640:   reconstructed with `MatSetValues()`

10642:   This is used by `PCBJACOBI` when a single block spans multiple MPI processes.

10644: .seealso: [](ch_matrices), `Mat`, `MatCreateRedundantMatrix()`, `MatCreateSubMatrices()`, `PCBJACOBI`
10645: @*/
10646: PetscErrorCode MatGetMultiProcBlock(Mat mat, MPI_Comm subComm, MatReuse scall, Mat *subMat)
10647: {
10648:   PetscMPIInt commsize, subCommSize;

10650:   PetscFunctionBegin;
10651:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &commsize));
10652:   PetscCallMPI(MPI_Comm_size(subComm, &subCommSize));
10653:   PetscCheck(subCommSize <= commsize, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "CommSize %d < SubCommZize %d", commsize, subCommSize);

10655:   PetscCheck(scall != MAT_REUSE_MATRIX || *subMat != mat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix");
10656:   PetscCall(PetscLogEventBegin(MAT_GetMultiProcBlock, mat, 0, 0, 0));
10657:   PetscUseTypeMethod(mat, getmultiprocblock, subComm, scall, subMat);
10658:   PetscCall(PetscLogEventEnd(MAT_GetMultiProcBlock, mat, 0, 0, 0));
10659:   PetscFunctionReturn(PETSC_SUCCESS);
10660: }

10662: /*@
10663:   MatGetLocalSubMatrix - Gets a reference to a submatrix specified in local numbering

10665:   Not Collective

10667:   Input Parameters:
10668: + mat   - matrix to extract local submatrix from
10669: . isrow - local row indices for submatrix
10670: - iscol - local column indices for submatrix

10672:   Output Parameter:
10673: . submat - the submatrix

10675:   Level: intermediate

10677:   Notes:
10678:   `submat` should be disposed of with `MatRestoreLocalSubMatrix()`.

10680:   Depending on the format of `mat`, the returned `submat` may not implement `MatMult()`.  Its communicator may be
10681:   the same as `mat`, it may be `PETSC_COMM_SELF`, or some other sub-communictor of `mat`'s.

10683:   `submat` always implements `MatSetValuesLocal()`.  If `isrow` and `iscol` have the same block size, then
10684:   `MatSetValuesBlockedLocal()` will also be implemented.

10686:   `mat` must have had a `ISLocalToGlobalMapping` provided to it with `MatSetLocalToGlobalMapping()`.
10687:   Matrices obtained with `DMCreateMatrix()` generally already have the local to global mapping provided.

10689: .seealso: [](ch_matrices), `Mat`, `MatRestoreLocalSubMatrix()`, `MatCreateLocalRef()`, `MatSetLocalToGlobalMapping()`
10690: @*/
10691: PetscErrorCode MatGetLocalSubMatrix(Mat mat, IS isrow, IS iscol, Mat *submat)
10692: {
10693:   PetscFunctionBegin;
10697:   PetscCheckSameComm(isrow, 2, iscol, 3);
10698:   PetscAssertPointer(submat, 4);
10699:   PetscCheck(mat->rmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must have local to global mapping provided before this call");

10701:   if (mat->ops->getlocalsubmatrix) {
10702:     PetscUseTypeMethod(mat, getlocalsubmatrix, isrow, iscol, submat);
10703:   } else {
10704:     PetscCall(MatCreateLocalRef(mat, isrow, iscol, submat));
10705:   }
10706:   PetscFunctionReturn(PETSC_SUCCESS);
10707: }

10709: /*@
10710:   MatRestoreLocalSubMatrix - Restores a reference to a submatrix specified in local numbering obtained with `MatGetLocalSubMatrix()`

10712:   Not Collective

10714:   Input Parameters:
10715: + mat    - matrix to extract local submatrix from
10716: . isrow  - local row indices for submatrix
10717: . iscol  - local column indices for submatrix
10718: - submat - the submatrix

10720:   Level: intermediate

10722: .seealso: [](ch_matrices), `Mat`, `MatGetLocalSubMatrix()`
10723: @*/
10724: PetscErrorCode MatRestoreLocalSubMatrix(Mat mat, IS isrow, IS iscol, Mat *submat)
10725: {
10726:   PetscFunctionBegin;
10730:   PetscCheckSameComm(isrow, 2, iscol, 3);
10731:   PetscAssertPointer(submat, 4);

10734:   if (mat->ops->restorelocalsubmatrix) {
10735:     PetscUseTypeMethod(mat, restorelocalsubmatrix, isrow, iscol, submat);
10736:   } else {
10737:     PetscCall(MatDestroy(submat));
10738:   }
10739:   *submat = NULL;
10740:   PetscFunctionReturn(PETSC_SUCCESS);
10741: }

10743: /*@
10744:   MatFindZeroDiagonals - Finds all the rows of a matrix that have zero or no diagonal entry in the matrix

10746:   Collective

10748:   Input Parameter:
10749: . mat - the matrix

10751:   Output Parameter:
10752: . is - if any rows have zero diagonals this contains the list of them

10754:   Level: developer

10756: .seealso: [](ch_matrices), `Mat`, `MatMultTranspose()`, `MatMultAdd()`, `MatMultTransposeAdd()`
10757: @*/
10758: PetscErrorCode MatFindZeroDiagonals(Mat mat, IS *is)
10759: {
10760:   PetscFunctionBegin;
10763:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
10764:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");

10766:   if (!mat->ops->findzerodiagonals) {
10767:     Vec                diag;
10768:     const PetscScalar *a;
10769:     PetscInt          *rows;
10770:     PetscInt           rStart, rEnd, r, nrow = 0;

10772:     PetscCall(MatCreateVecs(mat, &diag, NULL));
10773:     PetscCall(MatGetDiagonal(mat, diag));
10774:     PetscCall(MatGetOwnershipRange(mat, &rStart, &rEnd));
10775:     PetscCall(VecGetArrayRead(diag, &a));
10776:     for (r = 0; r < rEnd - rStart; ++r)
10777:       if (a[r] == 0.0) ++nrow;
10778:     PetscCall(PetscMalloc1(nrow, &rows));
10779:     nrow = 0;
10780:     for (r = 0; r < rEnd - rStart; ++r)
10781:       if (a[r] == 0.0) rows[nrow++] = r + rStart;
10782:     PetscCall(VecRestoreArrayRead(diag, &a));
10783:     PetscCall(VecDestroy(&diag));
10784:     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)mat), nrow, rows, PETSC_OWN_POINTER, is));
10785:   } else {
10786:     PetscUseTypeMethod(mat, findzerodiagonals, is);
10787:   }
10788:   PetscFunctionReturn(PETSC_SUCCESS);
10789: }

10791: /*@
10792:   MatFindOffBlockDiagonalEntries - Finds all the rows of a matrix that have entries outside of the main diagonal block (defined by the matrix block size)

10794:   Collective

10796:   Input Parameter:
10797: . mat - the matrix

10799:   Output Parameter:
10800: . is - contains the list of rows with off block diagonal entries

10802:   Level: developer

10804: .seealso: [](ch_matrices), `Mat`, `MatMultTranspose()`, `MatMultAdd()`, `MatMultTransposeAdd()`
10805: @*/
10806: PetscErrorCode MatFindOffBlockDiagonalEntries(Mat mat, IS *is)
10807: {
10808:   PetscFunctionBegin;
10811:   PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
10812:   PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");

10814:   PetscUseTypeMethod(mat, findoffblockdiagonalentries, is);
10815:   PetscFunctionReturn(PETSC_SUCCESS);
10816: }

10818: /*@C
10819:   MatInvertBlockDiagonal - Inverts the block diagonal entries.

10821:   Collective; No Fortran Support

10823:   Input Parameter:
10824: . mat - the matrix

10826:   Output Parameter:
10827: . values - the block inverses in column major order (FORTRAN-like)

10829:   Level: advanced

10831:   Notes:
10832:   The size of the blocks is determined by the block size of the matrix.

10834:   The blocks never overlap between two MPI processes, use `MatInvertVariableBlockEnvelope()` for that case

10836:   The blocks all have the same size, use `MatInvertVariableBlockDiagonal()` for variable block size

10838: .seealso: [](ch_matrices), `Mat`, `MatInvertVariableBlockEnvelope()`, `MatInvertBlockDiagonalMat()`
10839: @*/
10840: PetscErrorCode MatInvertBlockDiagonal(Mat mat, const PetscScalar *values[])
10841: {
10842:   PetscFunctionBegin;
10844:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
10845:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
10846:   PetscUseTypeMethod(mat, invertblockdiagonal, values);
10847:   PetscFunctionReturn(PETSC_SUCCESS);
10848: }

10850: /*@
10851:   MatInvertVariableBlockDiagonal - Inverts the point block diagonal entries.

10853:   Collective; No Fortran Support

10855:   Input Parameters:
10856: + mat     - the matrix
10857: . nblocks - the number of blocks on the process, set with `MatSetVariableBlockSizes()`
10858: - bsizes  - the size of each block on the process, set with `MatSetVariableBlockSizes()`

10860:   Output Parameter:
10861: . values - the block inverses in column major order (FORTRAN-like)

10863:   Level: advanced

10865:   Notes:
10866:   Use `MatInvertBlockDiagonal()` if all blocks have the same size

10868:   The blocks never overlap between two MPI processes, use `MatInvertVariableBlockEnvelope()` for that case

10870: .seealso: [](ch_matrices), `Mat`, `MatInvertBlockDiagonal()`, `MatSetVariableBlockSizes()`, `MatInvertVariableBlockEnvelope()`
10871: @*/
10872: PetscErrorCode MatInvertVariableBlockDiagonal(Mat mat, PetscInt nblocks, const PetscInt bsizes[], PetscScalar values[])
10873: {
10874:   PetscFunctionBegin;
10876:   PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
10877:   PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix");
10878:   PetscUseTypeMethod(mat, invertvariableblockdiagonal, nblocks, bsizes, values);
10879:   PetscFunctionReturn(PETSC_SUCCESS);
10880: }

10882: /*@
10883:   MatInvertBlockDiagonalMat - set the values of matrix C to be the inverted block diagonal of matrix A

10885:   Collective

10887:   Input Parameters:
10888: + A - the matrix
10889: - C - matrix with inverted block diagonal of `A`.  This matrix should be created and may have its type set.

10891:   Level: advanced

10893:   Note:
10894:   The blocksize of the matrix is used to determine the blocks on the diagonal of `C`

10896: .seealso: [](ch_matrices), `Mat`, `MatInvertBlockDiagonal()`
10897: @*/
10898: PetscErrorCode MatInvertBlockDiagonalMat(Mat A, Mat C)
10899: {
10900:   const PetscScalar *vals;
10901:   PetscInt          *dnnz;
10902:   PetscInt           m, rstart, rend, bs, i, j;

10904:   PetscFunctionBegin;
10905:   PetscCall(MatInvertBlockDiagonal(A, &vals));
10906:   PetscCall(MatGetBlockSize(A, &bs));
10907:   PetscCall(MatGetLocalSize(A, &m, NULL));
10908:   PetscCall(MatSetLayouts(C, A->rmap, A->cmap));
10909:   PetscCall(PetscMalloc1(m / bs, &dnnz));
10910:   for (j = 0; j < m / bs; j++) dnnz[j] = 1;
10911:   PetscCall(MatXAIJSetPreallocation(C, bs, dnnz, NULL, NULL, NULL));
10912:   PetscCall(PetscFree(dnnz));
10913:   PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
10914:   PetscCall(MatSetOption(C, MAT_ROW_ORIENTED, PETSC_FALSE));
10915:   for (i = rstart / bs; i < rend / bs; i++) PetscCall(MatSetValuesBlocked(C, 1, &i, 1, &i, &vals[(i - rstart / bs) * bs * bs], INSERT_VALUES));
10916:   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
10917:   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
10918:   PetscCall(MatSetOption(C, MAT_ROW_ORIENTED, PETSC_TRUE));
10919:   PetscFunctionReturn(PETSC_SUCCESS);
10920: }

10922: /*@
10923:   MatTransposeColoringDestroy - Destroys a coloring context for matrix product $C = A*B^T$ that was created
10924:   via `MatTransposeColoringCreate()`.

10926:   Collective

10928:   Input Parameter:
10929: . c - coloring context

10931:   Level: intermediate

10933: .seealso: [](ch_matrices), `Mat`, `MatTransposeColoringCreate()`
10934: @*/
10935: PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring *c)
10936: {
10937:   MatTransposeColoring matcolor = *c;

10939:   PetscFunctionBegin;
10940:   if (!matcolor) PetscFunctionReturn(PETSC_SUCCESS);
10941:   if (--((PetscObject)matcolor)->refct > 0) {
10942:     matcolor = NULL;
10943:     PetscFunctionReturn(PETSC_SUCCESS);
10944:   }

10946:   PetscCall(PetscFree3(matcolor->ncolumns, matcolor->nrows, matcolor->colorforrow));
10947:   PetscCall(PetscFree(matcolor->rows));
10948:   PetscCall(PetscFree(matcolor->den2sp));
10949:   PetscCall(PetscFree(matcolor->colorforcol));
10950:   PetscCall(PetscFree(matcolor->columns));
10951:   if (matcolor->brows > 0) PetscCall(PetscFree(matcolor->lstart));
10952:   PetscCall(PetscHeaderDestroy(c));
10953:   PetscFunctionReturn(PETSC_SUCCESS);
10954: }

10956: /*@
10957:   MatTransColoringApplySpToDen - Given a symbolic matrix product $C = A*B^T$ for which
10958:   a `MatTransposeColoring` context has been created, computes a dense $B^T$ by applying
10959:   `MatTransposeColoring` to sparse `B`.

10961:   Collective

10963:   Input Parameters:
10964: + coloring - coloring context created with `MatTransposeColoringCreate()`
10965: - B        - sparse matrix

10967:   Output Parameter:
10968: . Btdense - dense matrix $B^T$

10970:   Level: developer

10972:   Note:
10973:   These are used internally for some implementations of `MatRARt()`

10975: .seealso: [](ch_matrices), `Mat`, `MatTransposeColoringCreate()`, `MatTransposeColoringDestroy()`, `MatTransColoringApplyDenToSp()`
10976: @*/
10977: PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring coloring, Mat B, Mat Btdense)
10978: {
10979:   PetscFunctionBegin;

10984:   PetscCall((*B->ops->transcoloringapplysptoden)(coloring, B, Btdense));
10985:   PetscFunctionReturn(PETSC_SUCCESS);
10986: }

10988: /*@
10989:   MatTransColoringApplyDenToSp - Given a symbolic matrix product $C_{sp} = A*B^T$ for which
10990:   a `MatTransposeColoring` context has been created and a dense matrix $C_{den} = A*B^T_{dense}$
10991:   in which `B^T_{dens}` is obtained from `MatTransColoringApplySpToDen()`, recover sparse matrix
10992:   $C_{sp}$ from $C_{den}$.

10994:   Collective

10996:   Input Parameters:
10997: + matcoloring - coloring context created with `MatTransposeColoringCreate()`
10998: - Cden        - matrix product of a sparse matrix and a dense matrix Btdense

11000:   Output Parameter:
11001: . Csp - sparse matrix

11003:   Level: developer

11005:   Note:
11006:   These are used internally for some implementations of `MatRARt()`

11008: .seealso: [](ch_matrices), `Mat`, `MatTransposeColoringCreate()`, `MatTransposeColoringDestroy()`, `MatTransColoringApplySpToDen()`
11009: @*/
11010: PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring matcoloring, Mat Cden, Mat Csp)
11011: {
11012:   PetscFunctionBegin;

11017:   PetscCall((*Csp->ops->transcoloringapplydentosp)(matcoloring, Cden, Csp));
11018:   PetscCall(MatAssemblyBegin(Csp, MAT_FINAL_ASSEMBLY));
11019:   PetscCall(MatAssemblyEnd(Csp, MAT_FINAL_ASSEMBLY));
11020:   PetscFunctionReturn(PETSC_SUCCESS);
11021: }

11023: /*@
11024:   MatTransposeColoringCreate - Creates a matrix coloring context for the matrix product $C = A*B^T$.

11026:   Collective

11028:   Input Parameters:
11029: + mat        - the matrix product C
11030: - iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()`

11032:   Output Parameter:
11033: . color - the new coloring context

11035:   Level: intermediate

11037: .seealso: [](ch_matrices), `Mat`, `MatTransposeColoringDestroy()`, `MatTransColoringApplySpToDen()`,
11038:           `MatTransColoringApplyDenToSp()`
11039: @*/
11040: PetscErrorCode MatTransposeColoringCreate(Mat mat, ISColoring iscoloring, MatTransposeColoring *color)
11041: {
11042:   MatTransposeColoring c;
11043:   MPI_Comm             comm;

11045:   PetscFunctionBegin;
11046:   PetscAssertPointer(color, 3);

11048:   PetscCall(PetscLogEventBegin(MAT_TransposeColoringCreate, mat, 0, 0, 0));
11049:   PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
11050:   PetscCall(PetscHeaderCreate(c, MAT_TRANSPOSECOLORING_CLASSID, "MatTransposeColoring", "Matrix product C=A*B^T via coloring", "Mat", comm, MatTransposeColoringDestroy, NULL));
11051:   c->ctype = iscoloring->ctype;
11052:   PetscUseTypeMethod(mat, transposecoloringcreate, iscoloring, c);
11053:   *color = c;
11054:   PetscCall(PetscLogEventEnd(MAT_TransposeColoringCreate, mat, 0, 0, 0));
11055:   PetscFunctionReturn(PETSC_SUCCESS);
11056: }

11058: /*@
11059:   MatGetNonzeroState - Returns a 64-bit integer representing the current state of nonzeros in the matrix. If the
11060:   matrix has had new nonzero locations added to (or removed from) the matrix since the previous call, the value will be larger.

11062:   Not Collective

11064:   Input Parameter:
11065: . mat - the matrix

11067:   Output Parameter:
11068: . state - the current state

11070:   Level: intermediate

11072:   Notes:
11073:   You can only compare states from two different calls to the SAME matrix, you cannot compare calls between
11074:   different matrices

11076:   Use `PetscObjectStateGet()` to check for changes to the numerical values in a matrix

11078:   Use the result of `PetscObjectGetId()` to compare if a previously checked matrix is the same as the current matrix, do not compare object pointers.

11080: .seealso: [](ch_matrices), `Mat`, `PetscObjectStateGet()`, `PetscObjectGetId()`
11081: @*/
11082: PetscErrorCode MatGetNonzeroState(Mat mat, PetscObjectState *state)
11083: {
11084:   PetscFunctionBegin;
11086:   *state = mat->nonzerostate;
11087:   PetscFunctionReturn(PETSC_SUCCESS);
11088: }

11090: /*@
11091:   MatCreateMPIMatConcatenateSeqMat - Creates a single large PETSc matrix by concatenating sequential
11092:   matrices from each processor

11094:   Collective

11096:   Input Parameters:
11097: + comm   - the communicators the parallel matrix will live on
11098: . seqmat - the input sequential matrices
11099: . n      - number of local columns (or `PETSC_DECIDE`)
11100: - reuse  - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`

11102:   Output Parameter:
11103: . mpimat - the parallel matrix generated

11105:   Level: developer

11107:   Note:
11108:   The number of columns of the matrix in EACH processor MUST be the same.

11110: .seealso: [](ch_matrices), `Mat`
11111: @*/
11112: PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm comm, Mat seqmat, PetscInt n, MatReuse reuse, Mat *mpimat)
11113: {
11114:   PetscMPIInt size;

11116:   PetscFunctionBegin;
11117:   PetscCallMPI(MPI_Comm_size(comm, &size));
11118:   if (size == 1) {
11119:     if (reuse == MAT_INITIAL_MATRIX) {
11120:       PetscCall(MatDuplicate(seqmat, MAT_COPY_VALUES, mpimat));
11121:     } else {
11122:       PetscCall(MatCopy(seqmat, *mpimat, SAME_NONZERO_PATTERN));
11123:     }
11124:     PetscFunctionReturn(PETSC_SUCCESS);
11125:   }

11127:   PetscCheck(reuse != MAT_REUSE_MATRIX || seqmat != *mpimat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix");

11129:   PetscCall(PetscLogEventBegin(MAT_Merge, seqmat, 0, 0, 0));
11130:   PetscCall((*seqmat->ops->creatempimatconcatenateseqmat)(comm, seqmat, n, reuse, mpimat));
11131:   PetscCall(PetscLogEventEnd(MAT_Merge, seqmat, 0, 0, 0));
11132:   PetscFunctionReturn(PETSC_SUCCESS);
11133: }

11135: /*@
11136:   MatSubdomainsCreateCoalesce - Creates index subdomains by coalescing adjacent MPI processes' ownership ranges.

11138:   Collective

11140:   Input Parameters:
11141: + A - the matrix to create subdomains from
11142: - N - requested number of subdomains

11144:   Output Parameters:
11145: + n   - number of subdomains resulting on this MPI process
11146: - iss - `IS` list with indices of subdomains on this MPI process

11148:   Level: advanced

11150:   Note:
11151:   The number of subdomains must be smaller than the communicator size

11153: .seealso: [](ch_matrices), `Mat`, `IS`
11154: @*/
11155: PetscErrorCode MatSubdomainsCreateCoalesce(Mat A, PetscInt N, PetscInt *n, IS *iss[])
11156: {
11157:   MPI_Comm    comm, subcomm;
11158:   PetscMPIInt size, rank, color;
11159:   PetscInt    rstart, rend, k;

11161:   PetscFunctionBegin;
11162:   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
11163:   PetscCallMPI(MPI_Comm_size(comm, &size));
11164:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
11165:   PetscCheck(N >= 1 && N < (PetscInt)size, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of subdomains must be > 0 and < %d, got N = %" PetscInt_FMT, size, N);
11166:   *n    = 1;
11167:   k     = ((PetscInt)size) / N + ((PetscInt)size % N > 0); /* There are up to k ranks to a color */
11168:   color = rank / k;
11169:   PetscCallMPI(MPI_Comm_split(comm, color, rank, &subcomm));
11170:   PetscCall(PetscMalloc1(1, iss));
11171:   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
11172:   PetscCall(ISCreateStride(subcomm, rend - rstart, rstart, 1, iss[0]));
11173:   PetscCallMPI(MPI_Comm_free(&subcomm));
11174:   PetscFunctionReturn(PETSC_SUCCESS);
11175: }

11177: /*@
11178:   MatGalerkin - Constructs the coarse grid problem matrix via Galerkin projection.

11180:   If the interpolation and restriction operators are the same, uses `MatPtAP()`.
11181:   If they are not the same, uses `MatMatMatMult()`.

11183:   Once the coarse grid problem is constructed, correct for interpolation operators
11184:   that are not of full rank, which can legitimately happen in the case of non-nested
11185:   geometric multigrid.

11187:   Input Parameters:
11188: + restrct     - restriction operator
11189: . dA          - fine grid matrix
11190: . interpolate - interpolation operator
11191: . reuse       - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX`
11192: - fill        - expected fill, use `PETSC_DETERMINE` or `PETSC_DETERMINE` if you do not have a good estimate

11194:   Output Parameter:
11195: . A - the Galerkin coarse matrix

11197:   Options Database Key:
11198: . -pc_mg_galerkin <both,pmat,mat,none> - for what matrices the Galerkin process should be used

11200:   Level: developer

11202:   Note:
11203:   The deprecated `PETSC_DEFAULT` in `fill` also means use the current value

11205: .seealso: [](ch_matrices), `Mat`, `MatPtAP()`, `MatMatMatMult()`
11206: @*/
11207: PetscErrorCode MatGalerkin(Mat restrct, Mat dA, Mat interpolate, MatReuse reuse, PetscReal fill, Mat *A)
11208: {
11209:   IS  zerorows;
11210:   Vec diag;

11212:   PetscFunctionBegin;
11213:   PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Inplace product not supported");
11214:   /* Construct the coarse grid matrix */
11215:   if (interpolate == restrct) {
11216:     PetscCall(MatPtAP(dA, interpolate, reuse, fill, A));
11217:   } else {
11218:     PetscCall(MatMatMatMult(restrct, dA, interpolate, reuse, fill, A));
11219:   }

11221:   /* If the interpolation matrix is not of full rank, A will have zero rows.
11222:      This can legitimately happen in the case of non-nested geometric multigrid.
11223:      In that event, we set the rows of the matrix to the rows of the identity,
11224:      ignoring the equations (as the RHS will also be zero). */

11226:   PetscCall(MatFindZeroRows(*A, &zerorows));

11228:   if (zerorows != NULL) { /* if there are any zero rows */
11229:     PetscCall(MatCreateVecs(*A, &diag, NULL));
11230:     PetscCall(MatGetDiagonal(*A, diag));
11231:     PetscCall(VecISSet(diag, zerorows, 1.0));
11232:     PetscCall(MatDiagonalSet(*A, diag, INSERT_VALUES));
11233:     PetscCall(VecDestroy(&diag));
11234:     PetscCall(ISDestroy(&zerorows));
11235:   }
11236:   PetscFunctionReturn(PETSC_SUCCESS);
11237: }

11239: /*@C
11240:   MatSetOperation - Allows user to set a matrix operation for any matrix type

11242:   Logically Collective

11244:   Input Parameters:
11245: + mat - the matrix
11246: . op  - the name of the operation
11247: - f   - the function that provides the operation

11249:   Level: developer

11251:   Example Usage:
11252: .vb
11253:   extern PetscErrorCode usermult(Mat, Vec, Vec);

11255:   PetscCall(MatCreateXXX(comm, ..., &A));
11256:   PetscCall(MatSetOperation(A, MATOP_MULT, (PetscVoidFn *)usermult));
11257: .ve

11259:   Notes:
11260:   See the file `include/petscmat.h` for a complete list of matrix
11261:   operations, which all have the form MATOP_<OPERATION>, where
11262:   <OPERATION> is the name (in all capital letters) of the
11263:   user interface routine (e.g., `MatMult()` -> `MATOP_MULT`).

11265:   All user-provided functions (except for `MATOP_DESTROY`) should have the same calling
11266:   sequence as the usual matrix interface routines, since they
11267:   are intended to be accessed via the usual matrix interface
11268:   routines, e.g.,
11269: .vb
11270:   MatMult(Mat, Vec, Vec) -> usermult(Mat, Vec, Vec)
11271: .ve

11273:   In particular each function MUST return `PETSC_SUCCESS` on success and
11274:   nonzero on failure.

11276:   This routine is distinct from `MatShellSetOperation()` in that it can be called on any matrix type.

11278: .seealso: [](ch_matrices), `Mat`, `MatGetOperation()`, `MatCreateShell()`, `MatShellSetContext()`, `MatShellSetOperation()`
11279: @*/
11280: PetscErrorCode MatSetOperation(Mat mat, MatOperation op, void (*f)(void))
11281: {
11282:   PetscFunctionBegin;
11284:   if (op == MATOP_VIEW && !mat->ops->viewnative && f != (void (*)(void))mat->ops->view) mat->ops->viewnative = mat->ops->view;
11285:   (((void (**)(void))mat->ops)[op]) = f;
11286:   PetscFunctionReturn(PETSC_SUCCESS);
11287: }

11289: /*@C
11290:   MatGetOperation - Gets a matrix operation for any matrix type.

11292:   Not Collective

11294:   Input Parameters:
11295: + mat - the matrix
11296: - op  - the name of the operation

11298:   Output Parameter:
11299: . f - the function that provides the operation

11301:   Level: developer

11303:   Example Usage:
11304: .vb
11305:   PetscErrorCode (*usermult)(Mat, Vec, Vec);

11307:   MatGetOperation(A, MATOP_MULT, (void (**)(void))&usermult);
11308: .ve

11310:   Notes:
11311:   See the file include/petscmat.h for a complete list of matrix
11312:   operations, which all have the form MATOP_<OPERATION>, where
11313:   <OPERATION> is the name (in all capital letters) of the
11314:   user interface routine (e.g., `MatMult()` -> `MATOP_MULT`).

11316:   This routine is distinct from `MatShellGetOperation()` in that it can be called on any matrix type.

11318: .seealso: [](ch_matrices), `Mat`, `MatSetOperation()`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()`
11319: @*/
11320: PetscErrorCode MatGetOperation(Mat mat, MatOperation op, void (**f)(void))
11321: {
11322:   PetscFunctionBegin;
11324:   *f = (((void (**)(void))mat->ops)[op]);
11325:   PetscFunctionReturn(PETSC_SUCCESS);
11326: }

11328: /*@
11329:   MatHasOperation - Determines whether the given matrix supports the particular operation.

11331:   Not Collective

11333:   Input Parameters:
11334: + mat - the matrix
11335: - op  - the operation, for example, `MATOP_GET_DIAGONAL`

11337:   Output Parameter:
11338: . has - either `PETSC_TRUE` or `PETSC_FALSE`

11340:   Level: advanced

11342:   Note:
11343:   See `MatSetOperation()` for additional discussion on naming convention and usage of `op`.

11345: .seealso: [](ch_matrices), `Mat`, `MatCreateShell()`, `MatGetOperation()`, `MatSetOperation()`
11346: @*/
11347: PetscErrorCode MatHasOperation(Mat mat, MatOperation op, PetscBool *has)
11348: {
11349:   PetscFunctionBegin;
11351:   PetscAssertPointer(has, 3);
11352:   if (mat->ops->hasoperation) {
11353:     PetscUseTypeMethod(mat, hasoperation, op, has);
11354:   } else {
11355:     if (((void **)mat->ops)[op]) *has = PETSC_TRUE;
11356:     else {
11357:       *has = PETSC_FALSE;
11358:       if (op == MATOP_CREATE_SUBMATRIX) {
11359:         PetscMPIInt size;

11361:         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
11362:         if (size == 1) PetscCall(MatHasOperation(mat, MATOP_CREATE_SUBMATRICES, has));
11363:       }
11364:     }
11365:   }
11366:   PetscFunctionReturn(PETSC_SUCCESS);
11367: }

11369: /*@
11370:   MatHasCongruentLayouts - Determines whether the rows and columns layouts of the matrix are congruent

11372:   Collective

11374:   Input Parameter:
11375: . mat - the matrix

11377:   Output Parameter:
11378: . cong - either `PETSC_TRUE` or `PETSC_FALSE`

11380:   Level: beginner

11382: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatSetSizes()`, `PetscLayout`
11383: @*/
11384: PetscErrorCode MatHasCongruentLayouts(Mat mat, PetscBool *cong)
11385: {
11386:   PetscFunctionBegin;
11389:   PetscAssertPointer(cong, 2);
11390:   if (!mat->rmap || !mat->cmap) {
11391:     *cong = mat->rmap == mat->cmap ? PETSC_TRUE : PETSC_FALSE;
11392:     PetscFunctionReturn(PETSC_SUCCESS);
11393:   }
11394:   if (mat->congruentlayouts == PETSC_DECIDE) { /* first time we compare rows and cols layouts */
11395:     PetscCall(PetscLayoutSetUp(mat->rmap));
11396:     PetscCall(PetscLayoutSetUp(mat->cmap));
11397:     PetscCall(PetscLayoutCompare(mat->rmap, mat->cmap, cong));
11398:     if (*cong) mat->congruentlayouts = 1;
11399:     else mat->congruentlayouts = 0;
11400:   } else *cong = mat->congruentlayouts ? PETSC_TRUE : PETSC_FALSE;
11401:   PetscFunctionReturn(PETSC_SUCCESS);
11402: }

11404: PetscErrorCode MatSetInf(Mat A)
11405: {
11406:   PetscFunctionBegin;
11407:   PetscUseTypeMethod(A, setinf);
11408:   PetscFunctionReturn(PETSC_SUCCESS);
11409: }

11411: /*@
11412:   MatCreateGraph - create a scalar matrix (that is a matrix with one vertex for each block vertex in the original matrix), for use in graph algorithms
11413:   and possibly removes small values from the graph structure.

11415:   Collective

11417:   Input Parameters:
11418: + A       - the matrix
11419: . sym     - `PETSC_TRUE` indicates that the graph should be symmetrized
11420: . scale   - `PETSC_TRUE` indicates that the graph edge weights should be symmetrically scaled with the diagonal entry
11421: . filter  - filter value - < 0: does nothing; == 0: removes only 0.0 entries; otherwise: removes entries with abs(entries) <= value
11422: . num_idx - size of 'index' array
11423: - index   - array of block indices to use for graph strength of connection weight

11425:   Output Parameter:
11426: . graph - the resulting graph

11428:   Level: advanced

11430: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `PCGAMG`
11431: @*/
11432: PetscErrorCode MatCreateGraph(Mat A, PetscBool sym, PetscBool scale, PetscReal filter, PetscInt num_idx, PetscInt index[], Mat *graph)
11433: {
11434:   PetscFunctionBegin;
11438:   PetscAssertPointer(graph, 7);
11439:   PetscCall(PetscLogEventBegin(MAT_CreateGraph, A, 0, 0, 0));
11440:   PetscUseTypeMethod(A, creategraph, sym, scale, filter, num_idx, index, graph);
11441:   PetscCall(PetscLogEventEnd(MAT_CreateGraph, A, 0, 0, 0));
11442:   PetscFunctionReturn(PETSC_SUCCESS);
11443: }

11445: /*@
11446:   MatEliminateZeros - eliminate the nondiagonal zero entries in place from the nonzero structure of a sparse `Mat` in place,
11447:   meaning the same memory is used for the matrix, and no new memory is allocated.

11449:   Collective

11451:   Input Parameters:
11452: + A    - the matrix
11453: - keep - if for a given row of `A`, the diagonal coefficient is zero, indicates whether it should be left in the structure or eliminated as well

11455:   Level: intermediate

11457:   Developer Note:
11458:   The entries in the sparse matrix data structure are shifted to fill in the unneeded locations in the data. Thus the end
11459:   of the arrays in the data structure are unneeded.

11461: .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateGraph()`, `MatFilter()`
11462: @*/
11463: PetscErrorCode MatEliminateZeros(Mat A, PetscBool keep)
11464: {
11465:   PetscFunctionBegin;
11467:   PetscUseTypeMethod(A, eliminatezeros, keep);
11468:   PetscFunctionReturn(PETSC_SUCCESS);
11469: }