Actual source code: gcreate.c

  1: #include <petsc/private/matimpl.h>

  3: #include <../src/mat/impls/aij/seq/aij.h>
  4: #include <../src/mat/impls/aij/mpi/mpiaij.h>

  6: PetscErrorCode MatSetBlockSizes_Default(Mat mat, PetscInt rbs, PetscInt cbs)
  7: {
  8:   PetscFunctionBegin;
  9:   if (!mat->preallocated) PetscFunctionReturn(PETSC_SUCCESS);
 10:   PetscCheck(mat->rmap->bs <= 0 || mat->rmap->bs == rbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot change row block size %" PetscInt_FMT " to %" PetscInt_FMT, mat->rmap->bs, rbs);
 11:   PetscCheck(mat->cmap->bs <= 0 || mat->cmap->bs == cbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot change column block size %" PetscInt_FMT " to %" PetscInt_FMT, mat->cmap->bs, cbs);
 12:   PetscFunctionReturn(PETSC_SUCCESS);
 13: }

 15: PetscErrorCode MatShift_Basic(Mat Y, PetscScalar a)
 16: {
 17:   PetscInt    i, start, end, oldValA = 0, oldValB = 0;
 18:   PetscScalar alpha = a;
 19:   PetscBool   prevoption;
 20:   PetscBool   isSeqAIJDerived, isMPIAIJDerived; // all classes sharing SEQAIJHEADER or MPIAIJHEADER
 21:   Mat         A = NULL, B = NULL;

 23:   PetscFunctionBegin;
 24:   PetscCall(MatGetOption(Y, MAT_NO_OFF_PROC_ENTRIES, &prevoption));
 25:   PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
 26:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isSeqAIJDerived, MATSEQAIJ, MATSEQBAIJ, MATSEQSBAIJ, ""));
 27:   PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isMPIAIJDerived, MATMPIAIJ, MATMPIBAIJ, MATMPISBAIJ, ""));

 29:   if (isSeqAIJDerived) A = Y;
 30:   else if (isMPIAIJDerived) {
 31:     Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)Y->data;
 32:     A                  = mpiaij->A;
 33:     B                  = mpiaij->B;
 34:   }

 36:   if (A) {
 37:     oldValA                        = ((Mat_SeqAIJ *)A->data)->nonew;
 38:     ((Mat_SeqAIJ *)A->data)->nonew = 0; // so that new nonzero locations are allowed
 39:   }
 40:   if (B) {
 41:     oldValB                        = ((Mat_SeqAIJ *)B->data)->nonew;
 42:     ((Mat_SeqAIJ *)B->data)->nonew = 0;
 43:   }

 45:   PetscCall(MatGetOwnershipRange(Y, &start, &end));
 46:   for (i = start; i < end; i++) {
 47:     if (i < Y->cmap->N) PetscCall(MatSetValues(Y, 1, &i, 1, &i, &alpha, ADD_VALUES));
 48:   }
 49:   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
 50:   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
 51:   PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, prevoption));
 52:   if (A) ((Mat_SeqAIJ *)A->data)->nonew = oldValA;
 53:   if (B) ((Mat_SeqAIJ *)B->data)->nonew = oldValB;
 54:   PetscFunctionReturn(PETSC_SUCCESS);
 55: }

 57: /*@
 58:   MatCreate - Creates a matrix where the type is determined
 59:   from either a call to `MatSetType()` or from the options database
 60:   with a call to `MatSetFromOptions()`.

 62:   Collective

 64:   Input Parameter:
 65: . comm - MPI communicator

 67:   Output Parameter:
 68: . A - the matrix

 70:   Options Database Keys:
 71: + -mat_type seqaij   - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()`
 72: . -mat_type mpiaij   - `MATMPIAIJ` type, uses `MatCreateAIJ()`
 73: . -mat_type seqdense - `MATSEQDENSE`, uses `MatCreateSeqDense()`
 74: . -mat_type mpidense - `MATMPIDENSE` type, uses `MatCreateDense()`
 75: . -mat_type seqbaij  - `MATSEQBAIJ` type, uses `MatCreateSeqBAIJ()`
 76: - -mat_type mpibaij  - `MATMPIBAIJ` type, uses `MatCreateBAIJ()`

 78:    See the manpages for particular formats (e.g., `MATSEQAIJ`)
 79:    for additional format-specific options.

 81:   Level: beginner

 83:   Notes:
 84:   The default matrix type is `MATAIJ`, using the routines `MatCreateSeqAIJ()` or
 85:   `MatCreateAIJ()` if you do not set a type in the options database. If you never call
 86:   `MatSetType()` or `MatSetFromOptions()` it will generate an error when you try to use the
 87:   matrix.

 89: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
 90:           `MatCreateSeqDense()`, `MatCreateDense()`,
 91:           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
 92:           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
 93:           `MatConvert()`
 94: @*/
 95: PetscErrorCode MatCreate(MPI_Comm comm, Mat *A)
 96: {
 97:   Mat B;

 99:   PetscFunctionBegin;
100:   PetscAssertPointer(A, 2);
101:   PetscCall(MatInitializePackage());

103:   PetscCall(PetscHeaderCreate(B, MAT_CLASSID, "Mat", "Matrix", "Mat", comm, MatDestroy, MatView));
104:   PetscCall(PetscLayoutCreate(comm, &B->rmap));
105:   PetscCall(PetscLayoutCreate(comm, &B->cmap));
106:   PetscCall(PetscStrallocpy(VECSTANDARD, &B->defaultvectype));
107:   PetscCall(PetscStrallocpy(PETSCRANDER48, &B->defaultrandtype));

109:   B->symmetric                   = PETSC_BOOL3_UNKNOWN;
110:   B->hermitian                   = PETSC_BOOL3_UNKNOWN;
111:   B->structurally_symmetric      = PETSC_BOOL3_UNKNOWN;
112:   B->spd                         = PETSC_BOOL3_UNKNOWN;
113:   B->symmetry_eternal            = PETSC_FALSE;
114:   B->structural_symmetry_eternal = PETSC_FALSE;

116:   B->congruentlayouts = PETSC_DECIDE;
117:   B->preallocated     = PETSC_FALSE;
118: #if defined(PETSC_HAVE_DEVICE)
119:   B->boundtocpu = PETSC_TRUE;
120: #endif
121:   *A = B;
122:   PetscFunctionReturn(PETSC_SUCCESS);
123: }

125: /*@
126:   MatCreateFromOptions - Creates a matrix whose type is set from the options database

128:   Collective

130:   Input Parameters:
131: + comm   - MPI communicator
132: . prefix - [optional] prefix for the options database
133: . bs     - the blocksize (commonly 1)
134: . m      - the local number of rows (or `PETSC_DECIDE`)
135: . n      - the local number of columns (or `PETSC_DECIDE` or `PETSC_DETERMINE`)
136: . M      - the global number of rows (or `PETSC_DETERMINE`)
137: - N      - the global number of columns (or `PETSC_DETERMINE`)

139:   Output Parameter:
140: . A - the matrix

142:   Options Database Key:
143: . -mat_type - see `MatType`, for example `aij`, `aijcusparse`, `baij`, `sbaij`, `dense`, defaults to `aij`

145:   Level: beginner

147: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
148:           `MatCreateSeqDense()`, `MatCreateDense()`,
149:           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
150:           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
151:           `MatConvert()`, `MatCreate()`
152: @*/
153: PetscErrorCode MatCreateFromOptions(MPI_Comm comm, const char *prefix, PetscInt bs, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *A)
154: {
155:   PetscFunctionBegin;
156:   PetscAssertPointer(A, 8);
157:   PetscCall(MatCreate(comm, A));
158:   if (prefix) PetscCall(MatSetOptionsPrefix(*A, prefix));
159:   PetscCall(MatSetBlockSize(*A, bs));
160:   PetscCall(MatSetSizes(*A, m, n, M, N));
161:   PetscCall(MatSetFromOptions(*A));
162:   PetscFunctionReturn(PETSC_SUCCESS);
163: }

165: /*@
166:   MatSetErrorIfFailure - Causes `Mat` to generate an immediate error, for example a zero pivot, is detected.

168:   Logically Collective

170:   Input Parameters:
171: + mat - matrix obtained from `MatCreate()`
172: - flg - `PETSC_TRUE` indicates you want the error generated

174:   Level: advanced

176:   Note:
177:   If this flag is not set then the matrix operation will note the error and continue. The error may cause a later `PC` or `KSP` error
178:   or result in a `KSPConvergedReason` indicating the method did not converge.

180: .seealso: [](ch_matrices), `Mat`, `PCSetErrorIfFailure()`, `KSPConvergedReason`, `SNESConvergedReason`
181: @*/
182: PetscErrorCode MatSetErrorIfFailure(Mat mat, PetscBool flg)
183: {
184:   PetscFunctionBegin;
187:   mat->erroriffailure = flg;
188:   PetscFunctionReturn(PETSC_SUCCESS);
189: }

191: /*@
192:   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility

194:   Collective

196:   Input Parameters:
197: + A - the matrix
198: . m - number of local rows (or `PETSC_DECIDE`)
199: . n - number of local columns (or `PETSC_DECIDE`)
200: . M - number of global rows (or `PETSC_DETERMINE`)
201: - N - number of global columns (or `PETSC_DETERMINE`)

203:   Level: beginner

205:   Notes:
206:   `m` (`n`) and `M` (`N`) cannot be both `PETSC_DECIDE`
207:   If one processor calls this with `M` (`N`) of `PETSC_DECIDE` then all processors must, otherwise the program will hang.

209:   If `PETSC_DECIDE` is not used for the arguments 'm' and 'n', then the
210:   user must ensure that they are chosen to be compatible with the
211:   vectors. To do this, one first considers the matrix-vector product
212:   'y = A x'. The `m` that is used in the above routine must match the
213:   local size of 'y'. Likewise, the `n` used must match the local size of 'x'.

215:   If `m` and `n` are not `PETSC_DECIDE`, then the values determine the `PetscLayout` of the matrix and the ranges returned by
216:   `MatGetOwnershipRange()`,  `MatGetOwnershipRanges()`, `MatGetOwnershipRangeColumn()`, and `MatGetOwnershipRangesColumn()`.

218:   You cannot change the sizes once they have been set.

220:   The sizes must be set before `MatSetUp()` or MatXXXSetPreallocation() is called.

222: .seealso: [](ch_matrices), `Mat`, `MatGetSize()`, `PetscSplitOwnership()`, `MatGetOwnershipRange()`, `MatGetOwnershipRanges()`,
223:           `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRangesColumn()`, `PetscLayout`, `VecSetSizes()`
224: @*/
225: PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
226: {
227:   PetscFunctionBegin;
231:   PetscCheck(M <= 0 || m <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local row size %" PetscInt_FMT " cannot be larger than global row size %" PetscInt_FMT, m, M);
232:   PetscCheck(N <= 0 || n <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local column size %" PetscInt_FMT " cannot be larger than global column size %" PetscInt_FMT, n, N);
233:   PetscCheck((A->rmap->n < 0 || A->rmap->N < 0) || (A->rmap->n == m && (M <= 0 || A->rmap->N == M)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot change/reset row sizes to %" PetscInt_FMT " local %" PetscInt_FMT " global after previously setting them to %" PetscInt_FMT " local %" PetscInt_FMT " global", m, M,
234:              A->rmap->n, A->rmap->N);
235:   PetscCheck((A->cmap->n < 0 || A->cmap->N < 0) || (A->cmap->n == n && (N <= 0 || A->cmap->N == N)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot change/reset column sizes to %" PetscInt_FMT " local %" PetscInt_FMT " global after previously setting them to %" PetscInt_FMT " local %" PetscInt_FMT " global", n, N,
236:              A->cmap->n, A->cmap->N);
237:   A->rmap->n = m;
238:   A->cmap->n = n;
239:   A->rmap->N = M > -1 ? M : A->rmap->N;
240:   A->cmap->N = N > -1 ? N : A->cmap->N;
241:   PetscFunctionReturn(PETSC_SUCCESS);
242: }

244: /*@
245:   MatSetFromOptions - Creates a matrix where the type is determined
246:   from the options database.

248:   Collective

250:   Input Parameter:
251: . B - the matrix

253:   Options Database Keys:
254: + -mat_type seqaij   - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()`
255: . -mat_type mpiaij   - `MATMPIAIJ` type, uses `MatCreateAIJ()`
256: . -mat_type seqdense - `MATSEQDENSE` type, uses `MatCreateSeqDense()`
257: . -mat_type mpidense - `MATMPIDENSE`, uses `MatCreateDense()`
258: . -mat_type seqbaij  - `MATSEQBAIJ`, uses `MatCreateSeqBAIJ()`
259: - -mat_type mpibaij  - `MATMPIBAIJ`, uses `MatCreateBAIJ()`

261:    See the manpages for particular formats (e.g., `MATSEQAIJ`)
262:    for additional format-specific options.

264:   Level: beginner

266:   Notes:
267:   Generates a parallel MPI matrix if the communicator has more than one processor.  The default
268:   matrix type is `MATAIJ`, using the routines `MatCreateSeqAIJ()` and `MatCreateAIJ()` if you
269:   do not select a type in the options database.

271: .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
272:           `MatCreateSeqDense()`, `MatCreateDense()`,
273:           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
274:           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
275:           `MatConvert()`
276: @*/
277: PetscErrorCode MatSetFromOptions(Mat B)
278: {
279:   const char *deft = MATAIJ;
280:   char        type[256];
281:   PetscBool   flg, set;
282:   PetscInt    bind_below = 0, newbs = -1;

284:   PetscFunctionBegin;

287:   PetscObjectOptionsBegin((PetscObject)B);

289:   PetscCall(PetscOptionsInt("-mat_block_size", "Set the blocksize used to store the matrix", "MatSetBlockSize", newbs, &newbs, &flg));
290:   if (flg) {
291:     PetscCall(PetscLayoutSetBlockSize(B->rmap, newbs));
292:     PetscCall(PetscLayoutSetBlockSize(B->cmap, newbs));
293:   }

295:   PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, deft, type, PETSC_STATIC_ARRAY_LENGTH(type), &flg));
296:   if (flg) PetscCall(MatSetType(B, type));
297:   else if (!((PetscObject)B)->type_name) PetscCall(MatSetType(B, deft));

299:   if (newbs > 0) PetscTryTypeMethod(B, setblocksizes, newbs, newbs);

301:   PetscCall(PetscOptionsName("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", &B->checksymmetryonassembly));
302:   PetscCall(PetscOptionsReal("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", B->checksymmetrytol, &B->checksymmetrytol, NULL));
303:   PetscCall(PetscOptionsBool("-mat_null_space_test", "Checks if provided null space is correct in MatAssemblyEnd()", "MatSetNullSpaceTest", B->checknullspaceonassembly, &B->checknullspaceonassembly, NULL));
304:   PetscCall(PetscOptionsBool("-mat_error_if_failure", "Generate an error if an error occurs when factoring the matrix", "MatSetErrorIfFailure", B->erroriffailure, &B->erroriffailure, NULL));

306:   PetscTryTypeMethod(B, setfromoptions, PetscOptionsObject);

308:   flg = PETSC_FALSE;
309:   PetscCall(PetscOptionsBool("-mat_new_nonzero_location_err", "Generate an error if new nonzeros are created in the matrix nonzero structure (useful to test preallocation)", "MatSetOption", flg, &flg, &set));
310:   if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, flg));
311:   flg = PETSC_FALSE;
312:   PetscCall(PetscOptionsBool("-mat_new_nonzero_allocation_err", "Generate an error if new nonzeros are allocated in the matrix nonzero structure (useful to test preallocation)", "MatSetOption", flg, &flg, &set));
313:   if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, flg));
314:   flg = PETSC_FALSE;
315:   PetscCall(PetscOptionsBool("-mat_ignore_zero_entries", "For AIJ/IS matrices this will stop zero values from creating a zero location in the matrix", "MatSetOption", flg, &flg, &set));
316:   if (set) PetscCall(MatSetOption(B, MAT_IGNORE_ZERO_ENTRIES, flg));

318:   flg = PETSC_FALSE;
319:   PetscCall(PetscOptionsBool("-mat_form_explicit_transpose", "Hint to form an explicit transpose for operations like MatMultTranspose", "MatSetOption", flg, &flg, &set));
320:   if (set) PetscCall(MatSetOption(B, MAT_FORM_EXPLICIT_TRANSPOSE, flg));

322:   /* Bind to CPU if below a user-specified size threshold.
323:    * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types,
324:    * and putting it here makes is more maintainable than duplicating this for all. */
325:   PetscCall(PetscOptionsInt("-mat_bind_below", "Set the size threshold (in local rows) below which the Mat is bound to the CPU", "MatBindToCPU", bind_below, &bind_below, &flg));
326:   if (flg && B->rmap->n < bind_below) PetscCall(MatBindToCPU(B, PETSC_TRUE));

328:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
329:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)B, PetscOptionsObject));
330:   PetscOptionsEnd();
331:   PetscFunctionReturn(PETSC_SUCCESS);
332: }

334: /*@
335:   MatXAIJSetPreallocation - set preallocation for serial and parallel `MATAIJ`, `MATBAIJ`, and `MATSBAIJ` matrices and their unassembled versions.

337:   Collective

339:   Input Parameters:
340: + A     - matrix being preallocated
341: . bs    - block size
342: . dnnz  - number of nonzero column blocks per block row of diagonal part of parallel matrix
343: . onnz  - number of nonzero column blocks per block row of off-diagonal part of parallel matrix
344: . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix
345: - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix

347:   Level: beginner

349: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`,
350:           `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`,
351:           `PetscSplitOwnership()`
352: @*/
353: PetscErrorCode MatXAIJSetPreallocation(Mat A, PetscInt bs, const PetscInt dnnz[], const PetscInt onnz[], const PetscInt dnnzu[], const PetscInt onnzu[])
354: {
355:   PetscInt  cbs;
356:   PetscBool aij, is, hyp;

358:   PetscFunctionBegin;
359:   if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */
360:     PetscCall(MatSetBlockSize(A, bs));
361:   }
362:   PetscCall(PetscLayoutSetUp(A->rmap));
363:   PetscCall(PetscLayoutSetUp(A->cmap));
364:   PetscCall(MatGetBlockSizes(A, &bs, &cbs));
365:   /* these routines assumes bs == cbs, this should be checked somehow */
366:   PetscCall(MatSeqBAIJSetPreallocation(A, bs, 0, dnnz));
367:   PetscCall(MatMPIBAIJSetPreallocation(A, bs, 0, dnnz, 0, onnz));
368:   PetscCall(MatSeqSBAIJSetPreallocation(A, bs, 0, dnnzu));
369:   PetscCall(MatMPISBAIJSetPreallocation(A, bs, 0, dnnzu, 0, onnzu));
370:   /*
371:     In general, we have to do extra work to preallocate for scalar (AIJ) or unassembled (IS) matrices so we check whether it will do any
372:     good before going on with it.
373:   */
374:   PetscCall(PetscObjectHasFunction((PetscObject)A, "MatMPIAIJSetPreallocation_C", &aij));
375:   PetscCall(PetscObjectHasFunction((PetscObject)A, "MatISSetPreallocation_C", &is));
376:   PetscCall(PetscObjectHasFunction((PetscObject)A, "MatHYPRESetPreallocation_C", &hyp));
377:   if (!aij && !is && !hyp) PetscCall(PetscObjectHasFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", &aij));
378:   if (aij || is || hyp) {
379:     if (bs == cbs && bs == 1) {
380:       PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz));
381:       PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz, 0, onnz));
382:       PetscCall(MatISSetPreallocation(A, 0, dnnz, 0, onnz));
383: #if defined(PETSC_HAVE_HYPRE)
384:       PetscCall(MatHYPRESetPreallocation(A, 0, dnnz, 0, onnz));
385: #endif
386:     } else { /* Convert block-row precallocation to scalar-row */
387:       PetscInt i, m, *sdnnz, *sonnz;
388:       PetscCall(MatGetLocalSize(A, &m, NULL));
389:       PetscCall(PetscMalloc2((!!dnnz) * m, &sdnnz, (!!onnz) * m, &sonnz));
390:       for (i = 0; i < m; i++) {
391:         if (dnnz) sdnnz[i] = dnnz[i / bs] * cbs;
392:         if (onnz) sonnz[i] = onnz[i / bs] * cbs;
393:       }
394:       PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL));
395:       PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
396:       PetscCall(MatISSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
397: #if defined(PETSC_HAVE_HYPRE)
398:       PetscCall(MatHYPRESetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
399: #endif
400:       PetscCall(PetscFree2(sdnnz, sonnz));
401:     }
402:   }
403:   PetscFunctionReturn(PETSC_SUCCESS);
404: }

406: /*@C
407:   MatHeaderMerge - Merges some information from the header of `C` to `A`; the `C` object is then destroyed

409:   Collective, No Fortran Support

411:   Input Parameters:
412: + A - a `Mat` being merged into
413: - C - the `Mat` providing the merge information

415:   Level: developer

417:   Notes:
418:   `A` and `C` must be of the same type.
419:   The object list and query function list in `A` are retained, as well as the object name, and prefix.
420:   The object state of `A` is increased by 1.

422:   Developer Note:
423:   This is somewhat different from `MatHeaderReplace()`, it would be nice to merge the code

425: .seealso: `Mat`, `MatHeaderReplace()`
426:  @*/
427: PetscErrorCode MatHeaderMerge(Mat A, Mat *C)
428: {
429:   PetscInt          refct;
430:   PetscOps          Abops;
431:   struct _MatOps    Aops;
432:   char             *mtype, *mname, *mprefix;
433:   Mat_Product      *product;
434:   Mat_Redundant    *redundant;
435:   PetscObjectState  state;
436:   PetscObjectList   olist;
437:   PetscFunctionList qlist;

439:   PetscFunctionBegin;
442:   if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
443:   PetscCheckSameTypeAndComm(A, 1, *C, 2);
444:   /* save the parts of A we need */
445:   Abops     = ((PetscObject)A)->bops[0];
446:   Aops      = A->ops[0];
447:   refct     = ((PetscObject)A)->refct;
448:   mtype     = ((PetscObject)A)->type_name;
449:   mname     = ((PetscObject)A)->name;
450:   state     = ((PetscObject)A)->state;
451:   mprefix   = ((PetscObject)A)->prefix;
452:   product   = A->product;
453:   redundant = A->redundant;
454:   qlist     = ((PetscObject)A)->qlist;
455:   olist     = ((PetscObject)A)->olist;

457:   /* zero these so the destroy below does not free them */
458:   ((PetscObject)A)->type_name = NULL;
459:   ((PetscObject)A)->name      = NULL;
460:   ((PetscObject)A)->qlist     = NULL;
461:   ((PetscObject)A)->olist     = NULL;

463:   /*
464:      free all the interior data structures from mat
465:      cannot use PetscUseTypeMethod(A,destroy); because compiler
466:      thinks it may print NULL type_name and name
467:   */
468:   PetscTryTypeMethod(A, destroy);

470:   PetscCall(PetscFree(A->defaultvectype));
471:   PetscCall(PetscFree(A->defaultrandtype));
472:   PetscCall(PetscLayoutDestroy(&A->rmap));
473:   PetscCall(PetscLayoutDestroy(&A->cmap));
474:   PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A));

476:   /* copy C over to A */
477:   PetscCall(PetscFree(A->factorprefix));
478:   PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));

480:   /* return the parts of A we saved */
481:   ((PetscObject)A)->bops[0]   = Abops;
482:   A->ops[0]                   = Aops;
483:   ((PetscObject)A)->refct     = refct;
484:   ((PetscObject)A)->type_name = mtype;
485:   ((PetscObject)A)->name      = mname;
486:   ((PetscObject)A)->prefix    = mprefix;
487:   ((PetscObject)A)->state     = state + 1;
488:   A->product                  = product;
489:   A->redundant                = redundant;

491:   /* Append the saved lists */
492:   PetscCall(PetscFunctionListDuplicate(qlist, &((PetscObject)A)->qlist));
493:   PetscCall(PetscObjectListDuplicate(olist, &((PetscObject)A)->olist));
494:   PetscCall(PetscFunctionListDestroy(&qlist));
495:   PetscCall(PetscObjectListDestroy(&olist));

497:   /* since these two are copied into A we do not want them destroyed in C */
498:   ((PetscObject)*C)->qlist = NULL;
499:   ((PetscObject)*C)->olist = NULL;
500:   PetscCall(PetscHeaderDestroy(C));
501:   PetscFunctionReturn(PETSC_SUCCESS);
502: }

504: /*@
505:   MatHeaderReplace - Replaces the internal data of matrix `A` by the internal data of matrix `C` while deleting the outer wrapper of `C`

507:   Input Parameters:
508: + A - a `Mat` whose internal data is to be replaced
509: - C - the `Mat` providing new internal data for `A`

511:   Level: advanced

513:   Example Usage\:
514: .vb
515:   Mat C;
516:   MatCreateSeqAIJWithArrays(..., &C);
517:   MatHeaderReplace(A, &C);
518:   // C has been destroyed and A contains the matrix entries of C
519: .ve

521:   Note:
522:   This can be used inside a function provided to `SNESSetJacobian()`, `TSSetRHSJacobian()`, or `TSSetIJacobian()` in cases where the user code
523:   computes an entirely new sparse matrix  (generally with a different matrix nonzero structure/pattern) for each Newton update.
524:   It is usually better to reuse the matrix nonzero structure of `A` instead of constructing an entirely new one.

526:   Developer Note:
527:   This is somewhat different from `MatHeaderMerge()` it would be nice to merge the code

529: .seealso: `Mat`, `MatHeaderMerge()`
530:  @*/
531: PetscErrorCode MatHeaderReplace(Mat A, Mat *C)
532: {
533:   PetscInt         refct;
534:   PetscObjectState state;
535:   struct _p_Mat    buffer;
536:   MatStencilInfo   stencil;

538:   PetscFunctionBegin;
541:   if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
542:   PetscCheckSameComm(A, 1, *C, 2);
543:   PetscCheck(((PetscObject)*C)->refct == 1, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "Object C has refct %" PetscInt_FMT " > 1, would leave hanging reference", ((PetscObject)*C)->refct);

545:   /* swap C and A */
546:   refct   = ((PetscObject)A)->refct;
547:   state   = ((PetscObject)A)->state;
548:   stencil = A->stencil;
549:   PetscCall(PetscMemcpy(&buffer, A, sizeof(struct _p_Mat)));
550:   PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));
551:   PetscCall(PetscMemcpy(*C, &buffer, sizeof(struct _p_Mat)));
552:   ((PetscObject)A)->refct = refct;
553:   ((PetscObject)A)->state = state + 1;
554:   A->stencil              = stencil;

556:   ((PetscObject)*C)->refct = 1;
557:   PetscCall(MatDestroy(C));
558:   PetscFunctionReturn(PETSC_SUCCESS);
559: }

561: /*@
562:   MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU

564:   Logically Collective

566:   Input Parameters:
567: + A   - the matrix
568: - flg - bind to the CPU if value of `PETSC_TRUE`

570:   Level: intermediate

572: .seealso: [](ch_matrices), `Mat`, `MatBoundToCPU()`
573: @*/
574: PetscErrorCode MatBindToCPU(Mat A, PetscBool flg)
575: {
576:   PetscFunctionBegin;
579: #if defined(PETSC_HAVE_DEVICE)
580:   if (A->boundtocpu == flg) PetscFunctionReturn(PETSC_SUCCESS);
581:   A->boundtocpu = flg;
582:   PetscTryTypeMethod(A, bindtocpu, flg);
583: #endif
584:   PetscFunctionReturn(PETSC_SUCCESS);
585: }

587: /*@
588:   MatBoundToCPU - query if a matrix is bound to the CPU

590:   Input Parameter:
591: . A - the matrix

593:   Output Parameter:
594: . flg - the logical flag

596:   Level: intermediate

598: .seealso: [](ch_matrices), `Mat`, `MatBindToCPU()`
599: @*/
600: PetscErrorCode MatBoundToCPU(Mat A, PetscBool *flg)
601: {
602:   PetscFunctionBegin;
604:   PetscAssertPointer(flg, 2);
605: #if defined(PETSC_HAVE_DEVICE)
606:   *flg = A->boundtocpu;
607: #else
608:   *flg = PETSC_TRUE;
609: #endif
610:   PetscFunctionReturn(PETSC_SUCCESS);
611: }

613: PetscErrorCode MatSetValuesCOO_Basic(Mat A, const PetscScalar coo_v[], InsertMode imode)
614: {
615:   IS              is_coo_i, is_coo_j;
616:   const PetscInt *coo_i, *coo_j;
617:   PetscInt        n, n_i, n_j;
618:   PetscScalar     zero = 0.;

620:   PetscFunctionBegin;
621:   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_i", (PetscObject *)&is_coo_i));
622:   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_j", (PetscObject *)&is_coo_j));
623:   PetscCheck(is_coo_i, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_i IS");
624:   PetscCheck(is_coo_j, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_j IS");
625:   PetscCall(ISGetLocalSize(is_coo_i, &n_i));
626:   PetscCall(ISGetLocalSize(is_coo_j, &n_j));
627:   PetscCheck(n_i == n_j, PETSC_COMM_SELF, PETSC_ERR_COR, "Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT, n_i, n_j);
628:   PetscCall(ISGetIndices(is_coo_i, &coo_i));
629:   PetscCall(ISGetIndices(is_coo_j, &coo_j));
630:   if (imode != ADD_VALUES) PetscCall(MatZeroEntries(A));
631:   for (n = 0; n < n_i; n++) PetscCall(MatSetValue(A, coo_i[n], coo_j[n], coo_v ? coo_v[n] : zero, ADD_VALUES));
632:   PetscCall(ISRestoreIndices(is_coo_i, &coo_i));
633:   PetscCall(ISRestoreIndices(is_coo_j, &coo_j));
634:   PetscFunctionReturn(PETSC_SUCCESS);
635: }

637: PetscErrorCode MatSetPreallocationCOO_Basic(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
638: {
639:   Mat         preallocator;
640:   IS          is_coo_i, is_coo_j;
641:   PetscInt    ncoo_i;
642:   PetscScalar zero = 0.0;

644:   PetscFunctionBegin;
645:   PetscCall(PetscIntCast(ncoo, &ncoo_i));
646:   PetscCall(PetscLayoutSetUp(A->rmap));
647:   PetscCall(PetscLayoutSetUp(A->cmap));
648:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator));
649:   PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
650:   PetscCall(MatSetSizes(preallocator, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
651:   PetscCall(MatSetLayouts(preallocator, A->rmap, A->cmap));
652:   PetscCall(MatSetUp(preallocator));
653:   for (PetscCount n = 0; n < ncoo; n++) PetscCall(MatSetValue(preallocator, coo_i[n], coo_j[n], zero, INSERT_VALUES));
654:   PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
655:   PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
656:   PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A));
657:   PetscCall(MatDestroy(&preallocator));
658:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo_i, coo_i, PETSC_COPY_VALUES, &is_coo_i));
659:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo_i, coo_j, PETSC_COPY_VALUES, &is_coo_j));
660:   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_i", (PetscObject)is_coo_i));
661:   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_j", (PetscObject)is_coo_j));
662:   PetscCall(ISDestroy(&is_coo_i));
663:   PetscCall(ISDestroy(&is_coo_j));
664:   PetscFunctionReturn(PETSC_SUCCESS);
665: }

667: /*@C
668:   MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices

670:   Collective

672:   Input Parameters:
673: + A     - matrix being preallocated
674: . ncoo  - number of entries
675: . coo_i - row indices
676: - coo_j - column indices

678:   Level: beginner

680:   Notes:
681:   The indices within `coo_i` and `coo_j` may be modified within this function. The caller should not rely on them
682:   having any specific value after this function returns. The arrays can be freed or reused immediately
683:   after this function returns.

685:   Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed
686:   but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries
687:   are allowed and will be properly added or inserted to the matrix, unless the matrix option `MAT_IGNORE_OFF_PROC_ENTRIES`
688:   is set, in which case remote entries are ignored, or `MAT_NO_OFF_PROC_ENTRIES` is set, in which case an error will be generated.

690:   If you just want to create a sequential AIJ matrix (`MATSEQAIJ`), and your matrix entries in COO format are unique, you can also use
691:   `MatCreateSeqAIJFromTriple()`. But that is not recommended for iterative applications.

693: .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
694:           `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`,
695:           `DMSetMatrixPreallocateSkip()`, `MatCreateSeqAIJFromTriple()`
696: @*/
697: PetscErrorCode MatSetPreallocationCOO(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
698: {
699:   PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;

701:   PetscFunctionBegin;
704:   if (ncoo) PetscAssertPointer(coo_i, 3);
705:   if (ncoo) PetscAssertPointer(coo_j, 4);
706:   PetscCall(PetscLayoutSetUp(A->rmap));
707:   PetscCall(PetscLayoutSetUp(A->cmap));
708:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOO_C", &f));

710:   PetscCall(PetscLogEventBegin(MAT_PreallCOO, A, 0, 0, 0));
711:   if (f) {
712:     PetscCall((*f)(A, ncoo, coo_i, coo_j));
713:   } else { /* allow fallback, very slow */
714:     PetscCall(MatSetPreallocationCOO_Basic(A, ncoo, coo_i, coo_j));
715:   }
716:   PetscCall(PetscLogEventEnd(MAT_PreallCOO, A, 0, 0, 0));
717:   A->preallocated = PETSC_TRUE;
718:   A->nonzerostate++;
719:   PetscFunctionReturn(PETSC_SUCCESS);
720: }

722: /*@C
723:   MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices

725:   Collective

727:   Input Parameters:
728: + A     - matrix being preallocated
729: . ncoo  - number of entries
730: . coo_i - row indices (local numbering; may be modified)
731: - coo_j - column indices (local numbering; may be modified)

733:   Level: beginner

735:   Notes:
736:   The local indices are translated using the local to global mapping, thus `MatSetLocalToGlobalMapping()` must have been
737:   called prior to this function. For matrices created with `DMCreateMatrix()` the local to global mapping is often already provided.

739:   The indices `coo_i` and `coo_j` may be modified within this function. They might be translated to corresponding global
740:   indices, but the caller should not rely on them having any specific value after this function returns. The arrays
741:   can be freed or reused immediately after this function returns.

743:   Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed
744:   but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries
745:   are allowed and will be properly added or inserted to the matrix.

747: .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
748:           `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`,
749:           `DMSetMatrixPreallocateSkip()`
750: @*/
751: PetscErrorCode MatSetPreallocationCOOLocal(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
752: {
753:   PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;

755:   PetscFunctionBegin;
758:   if (ncoo) PetscAssertPointer(coo_i, 3);
759:   if (ncoo) PetscAssertPointer(coo_j, 4);
760:   PetscCall(PetscLayoutSetUp(A->rmap));
761:   PetscCall(PetscLayoutSetUp(A->cmap));

763:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", &f));
764:   if (f) {
765:     PetscCall((*f)(A, ncoo, coo_i, coo_j));
766:     A->nonzerostate++;
767:   } else {
768:     PetscInt               ncoo_i;
769:     ISLocalToGlobalMapping ltog_row, ltog_col;

771:     PetscCall(MatGetLocalToGlobalMapping(A, &ltog_row, &ltog_col));
772:     if (ltog_row) {
773:       PetscCall(PetscIntCast(ncoo, &ncoo_i));
774:       PetscCall(ISLocalToGlobalMappingApply(ltog_row, ncoo_i, coo_i, coo_i));
775:     }
776:     if (ltog_col) {
777:       PetscCall(PetscIntCast(ncoo, &ncoo_i));
778:       PetscCall(ISLocalToGlobalMappingApply(ltog_col, ncoo_i, coo_j, coo_j));
779:     }
780:     PetscCall(MatSetPreallocationCOO(A, ncoo, coo_i, coo_j));
781:   }
782:   A->preallocated = PETSC_TRUE;
783:   PetscFunctionReturn(PETSC_SUCCESS);
784: }

786: /*@
787:   MatSetValuesCOO - set values at once in a matrix preallocated using `MatSetPreallocationCOO()`

789:   Collective

791:   Input Parameters:
792: + A     - matrix being preallocated
793: . coo_v - the matrix values (can be `NULL`)
794: - imode - the insert mode

796:   Level: beginner

798:   Notes:
799:   The values must follow the order of the indices prescribed with `MatSetPreallocationCOO()` or `MatSetPreallocationCOOLocal()`.

801:   When repeated entries are specified in the COO indices the `coo_v` values are first properly summed, regardless of the value of imode.
802:   The imode flag indicates if coo_v must be added to the current values of the matrix (`ADD_VALUES`) or overwritten (`INSERT_VALUES`).

804:   `MatAssemblyBegin()` and `MatAssemblyEnd()` do not need to be called after this routine. It automatically handles the assembly process.

806: .seealso: [](ch_matrices), `Mat`, `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`
807: @*/
808: PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode)
809: {
810:   PetscErrorCode (*f)(Mat, const PetscScalar[], InsertMode) = NULL;
811:   PetscBool oldFlg;

813:   PetscFunctionBegin;
816:   MatCheckPreallocated(A, 1);
818:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetValuesCOO_C", &f));
819:   PetscCall(PetscLogEventBegin(MAT_SetVCOO, A, 0, 0, 0));
820:   if (f) {
821:     PetscCall((*f)(A, coo_v, imode)); // all known COO implementations do not use MatStash. They do their own off-proc communication
822:     PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &oldFlg));
823:     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); // set A->nooffprocentries to avoid costly MatStash scatter in MatAssembly
824:   } else {
825:     PetscCall(MatSetValuesCOO_Basic(A, coo_v, imode)); // fall back to MatSetValues, which might use MatStash
826:   }
827:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
828:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
829:   if (f) PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, oldFlg));
830:   PetscCall(PetscLogEventEnd(MAT_SetVCOO, A, 0, 0, 0));
831:   PetscFunctionReturn(PETSC_SUCCESS);
832: }

834: /*@
835:   MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects

837:   Input Parameters:
838: + A   - the matrix
839: - flg - flag indicating whether the boundtocpu flag should be propagated

841:   Level: developer

843:   Notes:
844:   If the value of flg is set to true, the following will occur
845: +   `MatCreateSubMatrices()` and `MatCreateRedundantMatrix()` - bind created matrices to CPU if the input matrix is bound to the CPU.
846: -   `MatCreateVecs()` - bind created vectors to CPU if the input matrix is bound to the CPU.

848:   The bindingpropagates flag itself is also propagated by the above routines.

850:   Developer Notes:
851:   If the fine-scale `DMDA` has the `-dm_bind_below` option set to true, then `DMCreateInterpolationScale()` calls `MatSetBindingPropagates()`
852:   on the restriction/interpolation operator to set the bindingpropagates flag to true.

854: .seealso: [](ch_matrices), `Mat`, `VecSetBindingPropagates()`, `MatGetBindingPropagates()`
855: @*/
856: PetscErrorCode MatSetBindingPropagates(Mat A, PetscBool flg)
857: {
858:   PetscFunctionBegin;
860: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
861:   A->bindingpropagates = flg;
862: #endif
863:   PetscFunctionReturn(PETSC_SUCCESS);
864: }

866: /*@
867:   MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects

869:   Input Parameter:
870: . A - the matrix

872:   Output Parameter:
873: . flg - flag indicating whether the boundtocpu flag will be propagated

875:   Level: developer

877: .seealso: [](ch_matrices), `Mat`, `MatSetBindingPropagates()`
878: @*/
879: PetscErrorCode MatGetBindingPropagates(Mat A, PetscBool *flg)
880: {
881:   PetscFunctionBegin;
883:   PetscAssertPointer(flg, 2);
884: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
885:   *flg = A->bindingpropagates;
886: #else
887:   *flg = PETSC_FALSE;
888: #endif
889:   PetscFunctionReturn(PETSC_SUCCESS);
890: }