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;

284:   PetscFunctionBegin;

287:   PetscObjectOptionsBegin((PetscObject)B);

289:   if (B->rmap->bs < 0) {
290:     PetscInt newbs = -1;
291:     PetscCall(PetscOptionsInt("-mat_block_size", "Set the blocksize used to store the matrix", "MatSetBlockSize", newbs, &newbs, &flg));
292:     if (flg) {
293:       PetscCall(PetscLayoutSetBlockSize(B->rmap, newbs));
294:       PetscCall(PetscLayoutSetBlockSize(B->cmap, newbs));
295:     }
296:   }

298:   PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, deft, type, 256, &flg));
299:   if (flg) {
300:     PetscCall(MatSetType(B, type));
301:   } else if (!((PetscObject)B)->type_name) {
302:     PetscCall(MatSetType(B, deft));
303:   }

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

310:   PetscTryTypeMethod(B, setfromoptions, PetscOptionsObject);

312:   flg = PETSC_FALSE;
313:   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));
314:   if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, flg));
315:   flg = PETSC_FALSE;
316:   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));
317:   if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, flg));
318:   flg = PETSC_FALSE;
319:   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));
320:   if (set) PetscCall(MatSetOption(B, MAT_IGNORE_ZERO_ENTRIES, flg));

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

326:   /* Bind to CPU if below a user-specified size threshold.
327:    * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types,
328:    * and putting it here makes is more maintainable than duplicating this for all. */
329:   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));
330:   if (flg && B->rmap->n < bind_below) PetscCall(MatBindToCPU(B, PETSC_TRUE));

332:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
333:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)B, PetscOptionsObject));
334:   PetscOptionsEnd();
335:   PetscFunctionReturn(PETSC_SUCCESS);
336: }

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

341:   Collective

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

351:   Level: beginner

353: .seealso: [](ch_matrices), `Mat`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`,
354:           `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`,
355:           `PetscSplitOwnership()`
356: @*/
357: PetscErrorCode MatXAIJSetPreallocation(Mat A, PetscInt bs, const PetscInt dnnz[], const PetscInt onnz[], const PetscInt dnnzu[], const PetscInt onnzu[])
358: {
359:   PetscInt cbs;
360:   void (*aij)(void);
361:   void (*is)(void);
362:   void (*hyp)(void) = NULL;

364:   PetscFunctionBegin;
365:   if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */
366:     PetscCall(MatSetBlockSize(A, bs));
367:   }
368:   PetscCall(PetscLayoutSetUp(A->rmap));
369:   PetscCall(PetscLayoutSetUp(A->cmap));
370:   PetscCall(MatGetBlockSizes(A, &bs, &cbs));
371:   /* these routines assumes bs == cbs, this should be checked somehow */
372:   PetscCall(MatSeqBAIJSetPreallocation(A, bs, 0, dnnz));
373:   PetscCall(MatMPIBAIJSetPreallocation(A, bs, 0, dnnz, 0, onnz));
374:   PetscCall(MatSeqSBAIJSetPreallocation(A, bs, 0, dnnzu));
375:   PetscCall(MatMPISBAIJSetPreallocation(A, bs, 0, dnnzu, 0, onnzu));
376:   /*
377:     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
378:     good before going on with it.
379:   */
380:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatMPIAIJSetPreallocation_C", &aij));
381:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatISSetPreallocation_C", &is));
382: #if defined(PETSC_HAVE_HYPRE)
383:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatHYPRESetPreallocation_C", &hyp));
384: #endif
385:   if (!aij && !is && !hyp) PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", &aij));
386:   if (aij || is || hyp) {
387:     if (bs == cbs && bs == 1) {
388:       PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz));
389:       PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz, 0, onnz));
390:       PetscCall(MatISSetPreallocation(A, 0, dnnz, 0, onnz));
391: #if defined(PETSC_HAVE_HYPRE)
392:       PetscCall(MatHYPRESetPreallocation(A, 0, dnnz, 0, onnz));
393: #endif
394:     } else { /* Convert block-row precallocation to scalar-row */
395:       PetscInt i, m, *sdnnz, *sonnz;
396:       PetscCall(MatGetLocalSize(A, &m, NULL));
397:       PetscCall(PetscMalloc2((!!dnnz) * m, &sdnnz, (!!onnz) * m, &sonnz));
398:       for (i = 0; i < m; i++) {
399:         if (dnnz) sdnnz[i] = dnnz[i / bs] * cbs;
400:         if (onnz) sonnz[i] = onnz[i / bs] * cbs;
401:       }
402:       PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL));
403:       PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
404:       PetscCall(MatISSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
405: #if defined(PETSC_HAVE_HYPRE)
406:       PetscCall(MatHYPRESetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
407: #endif
408:       PetscCall(PetscFree2(sdnnz, sonnz));
409:     }
410:   }
411:   PetscFunctionReturn(PETSC_SUCCESS);
412: }

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

417:   Collective, No Fortran Support

419:   Input Parameters:
420: + A - a `Mat` being merged into
421: - C - the `Mat` providing the merge information

423:   Level: developer

425:   Notes:
426:   `A` and `C` must be of the same type.
427:   The object list and query function list in `A` are retained, as well as the object name, and prefix.
428:   The object state of `A` is increased by 1.

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

433: .seealso: `Mat`, `MatHeaderReplace()`
434:  @*/
435: PetscErrorCode MatHeaderMerge(Mat A, Mat *C)
436: {
437:   PetscInt          refct;
438:   PetscOps          Abops;
439:   struct _MatOps    Aops;
440:   char             *mtype, *mname, *mprefix;
441:   Mat_Product      *product;
442:   Mat_Redundant    *redundant;
443:   PetscObjectState  state;
444:   PetscObjectList   olist;
445:   PetscFunctionList qlist;

447:   PetscFunctionBegin;
450:   if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
451:   PetscCheckSameTypeAndComm(A, 1, *C, 2);
452:   /* save the parts of A we need */
453:   Abops     = ((PetscObject)A)->bops[0];
454:   Aops      = A->ops[0];
455:   refct     = ((PetscObject)A)->refct;
456:   mtype     = ((PetscObject)A)->type_name;
457:   mname     = ((PetscObject)A)->name;
458:   state     = ((PetscObject)A)->state;
459:   mprefix   = ((PetscObject)A)->prefix;
460:   product   = A->product;
461:   redundant = A->redundant;
462:   qlist     = ((PetscObject)A)->qlist;
463:   olist     = ((PetscObject)A)->olist;

465:   /* zero these so the destroy below does not free them */
466:   ((PetscObject)A)->type_name = NULL;
467:   ((PetscObject)A)->name      = NULL;
468:   ((PetscObject)A)->qlist     = NULL;
469:   ((PetscObject)A)->olist     = NULL;

471:   /*
472:      free all the interior data structures from mat
473:      cannot use PetscUseTypeMethod(A,destroy); because compiler
474:      thinks it may print NULL type_name and name
475:   */
476:   PetscTryTypeMethod(A, destroy);

478:   PetscCall(PetscFree(A->defaultvectype));
479:   PetscCall(PetscFree(A->defaultrandtype));
480:   PetscCall(PetscLayoutDestroy(&A->rmap));
481:   PetscCall(PetscLayoutDestroy(&A->cmap));
482:   PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A));

484:   /* copy C over to A */
485:   PetscCall(PetscFree(A->factorprefix));
486:   PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));

488:   /* return the parts of A we saved */
489:   ((PetscObject)A)->bops[0]   = Abops;
490:   A->ops[0]                   = Aops;
491:   ((PetscObject)A)->refct     = refct;
492:   ((PetscObject)A)->type_name = mtype;
493:   ((PetscObject)A)->name      = mname;
494:   ((PetscObject)A)->prefix    = mprefix;
495:   ((PetscObject)A)->state     = state + 1;
496:   A->product                  = product;
497:   A->redundant                = redundant;

499:   /* Append the saved lists */
500:   PetscCall(PetscFunctionListDuplicate(qlist, &((PetscObject)A)->qlist));
501:   PetscCall(PetscObjectListDuplicate(olist, &((PetscObject)A)->olist));
502:   PetscCall(PetscFunctionListDestroy(&qlist));
503:   PetscCall(PetscObjectListDestroy(&olist));

505:   /* since these two are copied into A we do not want them destroyed in C */
506:   ((PetscObject)*C)->qlist = NULL;
507:   ((PetscObject)*C)->olist = NULL;
508:   PetscCall(PetscHeaderDestroy(C));
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }

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

515:   Input Parameters:
516: + A - a `Mat` whose internal data is to be replaced
517: - C - the `Mat` providing new internal data for `A`

519:   Level: advanced

521:   Example Usage\:
522: .vb
523:   Mat C;
524:   MatCreateSeqAIJWithArrays(..., &C);
525:   MatHeaderReplace(A, &C);
526:   // C has been destroyed and A contains the matrix entries of C
527: .ve

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

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

537: .seealso: `Mat`, `MatHeaderMerge()`
538:  @*/
539: PetscErrorCode MatHeaderReplace(Mat A, Mat *C)
540: {
541:   PetscInt         refct;
542:   PetscObjectState state;
543:   struct _p_Mat    buffer;
544:   MatStencilInfo   stencil;

546:   PetscFunctionBegin;
549:   if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
550:   PetscCheckSameComm(A, 1, *C, 2);
551:   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);

553:   /* swap C and A */
554:   refct   = ((PetscObject)A)->refct;
555:   state   = ((PetscObject)A)->state;
556:   stencil = A->stencil;
557:   PetscCall(PetscMemcpy(&buffer, A, sizeof(struct _p_Mat)));
558:   PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));
559:   PetscCall(PetscMemcpy(*C, &buffer, sizeof(struct _p_Mat)));
560:   ((PetscObject)A)->refct = refct;
561:   ((PetscObject)A)->state = state + 1;
562:   A->stencil              = stencil;

564:   ((PetscObject)*C)->refct = 1;
565:   PetscCall(MatDestroy(C));
566:   PetscFunctionReturn(PETSC_SUCCESS);
567: }

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

572:   Logically Collective

574:   Input Parameters:
575: + A   - the matrix
576: - flg - bind to the CPU if value of `PETSC_TRUE`

578:   Level: intermediate

580: .seealso: [](ch_matrices), `Mat`, `MatBoundToCPU()`
581: @*/
582: PetscErrorCode MatBindToCPU(Mat A, PetscBool flg)
583: {
584:   PetscFunctionBegin;
587: #if defined(PETSC_HAVE_DEVICE)
588:   if (A->boundtocpu == flg) PetscFunctionReturn(PETSC_SUCCESS);
589:   A->boundtocpu = flg;
590:   PetscTryTypeMethod(A, bindtocpu, flg);
591: #endif
592:   PetscFunctionReturn(PETSC_SUCCESS);
593: }

595: /*@
596:   MatBoundToCPU - query if a matrix is bound to the CPU

598:   Input Parameter:
599: . A - the matrix

601:   Output Parameter:
602: . flg - the logical flag

604:   Level: intermediate

606: .seealso: [](ch_matrices), `Mat`, `MatBindToCPU()`
607: @*/
608: PetscErrorCode MatBoundToCPU(Mat A, PetscBool *flg)
609: {
610:   PetscFunctionBegin;
612:   PetscAssertPointer(flg, 2);
613: #if defined(PETSC_HAVE_DEVICE)
614:   *flg = A->boundtocpu;
615: #else
616:   *flg = PETSC_TRUE;
617: #endif
618:   PetscFunctionReturn(PETSC_SUCCESS);
619: }

621: PetscErrorCode MatSetValuesCOO_Basic(Mat A, const PetscScalar coo_v[], InsertMode imode)
622: {
623:   IS              is_coo_i, is_coo_j;
624:   const PetscInt *coo_i, *coo_j;
625:   PetscInt        n, n_i, n_j;
626:   PetscScalar     zero = 0.;

628:   PetscFunctionBegin;
629:   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_i", (PetscObject *)&is_coo_i));
630:   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_j", (PetscObject *)&is_coo_j));
631:   PetscCheck(is_coo_i, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_i IS");
632:   PetscCheck(is_coo_j, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_j IS");
633:   PetscCall(ISGetLocalSize(is_coo_i, &n_i));
634:   PetscCall(ISGetLocalSize(is_coo_j, &n_j));
635:   PetscCheck(n_i == n_j, PETSC_COMM_SELF, PETSC_ERR_COR, "Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT, n_i, n_j);
636:   PetscCall(ISGetIndices(is_coo_i, &coo_i));
637:   PetscCall(ISGetIndices(is_coo_j, &coo_j));
638:   if (imode != ADD_VALUES) PetscCall(MatZeroEntries(A));
639:   for (n = 0; n < n_i; n++) PetscCall(MatSetValue(A, coo_i[n], coo_j[n], coo_v ? coo_v[n] : zero, ADD_VALUES));
640:   PetscCall(ISRestoreIndices(is_coo_i, &coo_i));
641:   PetscCall(ISRestoreIndices(is_coo_j, &coo_j));
642:   PetscFunctionReturn(PETSC_SUCCESS);
643: }

645: PetscErrorCode MatSetPreallocationCOO_Basic(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
646: {
647:   Mat         preallocator;
648:   IS          is_coo_i, is_coo_j;
649:   PetscScalar zero = 0.0;

651:   PetscFunctionBegin;
652:   PetscCheck(ncoo <= PETSC_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "ncoo %" PetscCount_FMT " overflowed PetscInt; configure --with-64-bit-indices or request support", ncoo);
653:   PetscCall(PetscLayoutSetUp(A->rmap));
654:   PetscCall(PetscLayoutSetUp(A->cmap));
655:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator));
656:   PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
657:   PetscCall(MatSetSizes(preallocator, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
658:   PetscCall(MatSetLayouts(preallocator, A->rmap, A->cmap));
659:   PetscCall(MatSetUp(preallocator));
660:   for (PetscCount n = 0; n < ncoo; n++) PetscCall(MatSetValue(preallocator, coo_i[n], coo_j[n], zero, INSERT_VALUES));
661:   PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
662:   PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
663:   PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A));
664:   PetscCall(MatDestroy(&preallocator));
665:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, (PetscInt)ncoo, coo_i, PETSC_COPY_VALUES, &is_coo_i));
666:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, (PetscInt)ncoo, coo_j, PETSC_COPY_VALUES, &is_coo_j));
667:   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_i", (PetscObject)is_coo_i));
668:   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_j", (PetscObject)is_coo_j));
669:   PetscCall(ISDestroy(&is_coo_i));
670:   PetscCall(ISDestroy(&is_coo_j));
671:   PetscFunctionReturn(PETSC_SUCCESS);
672: }

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

677:   Collective

679:   Input Parameters:
680: + A     - matrix being preallocated
681: . ncoo  - number of entries
682: . coo_i - row indices
683: - coo_j - column indices

685:   Level: beginner

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

692:   Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed
693:   but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries
694:   are allowed and will be properly added or inserted to the matrix, unless the matrix option `MAT_IGNORE_OFF_PROC_ENTRIES`
695:   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.

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

700: .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
701:           `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`,
702:           `DMSetMatrixPreallocateSkip()`, `MatCreateSeqAIJFromTriple()`
703: @*/
704: PetscErrorCode MatSetPreallocationCOO(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
705: {
706:   PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;

708:   PetscFunctionBegin;
711:   if (ncoo) PetscAssertPointer(coo_i, 3);
712:   if (ncoo) PetscAssertPointer(coo_j, 4);
713:   PetscCall(PetscLayoutSetUp(A->rmap));
714:   PetscCall(PetscLayoutSetUp(A->cmap));
715:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOO_C", &f));

717:   PetscCall(PetscLogEventBegin(MAT_PreallCOO, A, 0, 0, 0));
718:   if (f) {
719:     PetscCall((*f)(A, ncoo, coo_i, coo_j));
720:   } else { /* allow fallback, very slow */
721:     PetscCall(MatSetPreallocationCOO_Basic(A, ncoo, coo_i, coo_j));
722:   }
723:   PetscCall(PetscLogEventEnd(MAT_PreallCOO, A, 0, 0, 0));
724:   A->preallocated = PETSC_TRUE;
725:   A->nonzerostate++;
726:   PetscFunctionReturn(PETSC_SUCCESS);
727: }

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

732:   Collective

734:   Input Parameters:
735: + A     - matrix being preallocated
736: . ncoo  - number of entries
737: . coo_i - row indices (local numbering; may be modified)
738: - coo_j - column indices (local numbering; may be modified)

740:   Level: beginner

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

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

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

754: .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
755:           `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`,
756:           `DMSetMatrixPreallocateSkip()`
757: @*/
758: PetscErrorCode MatSetPreallocationCOOLocal(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
759: {
760:   PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;

762:   PetscFunctionBegin;
765:   if (ncoo) PetscAssertPointer(coo_i, 3);
766:   if (ncoo) PetscAssertPointer(coo_j, 4);
767:   PetscCheck(ncoo <= PETSC_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "ncoo %" PetscCount_FMT " overflowed PetscInt; configure --with-64-bit-indices or request support", ncoo);
768:   PetscCall(PetscLayoutSetUp(A->rmap));
769:   PetscCall(PetscLayoutSetUp(A->cmap));

771:   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", &f));
772:   if (f) {
773:     PetscCall((*f)(A, ncoo, coo_i, coo_j));
774:     A->nonzerostate++;
775:   } else {
776:     ISLocalToGlobalMapping ltog_row, ltog_col;

778:     PetscCall(MatGetLocalToGlobalMapping(A, &ltog_row, &ltog_col));
779:     if (ltog_row) PetscCall(ISLocalToGlobalMappingApply(ltog_row, (PetscInt)ncoo, coo_i, coo_i));
780:     if (ltog_col) PetscCall(ISLocalToGlobalMappingApply(ltog_col, (PetscInt)ncoo, coo_j, coo_j));
781:     PetscCall(MatSetPreallocationCOO(A, ncoo, coo_i, coo_j));
782:   }
783:   A->preallocated = PETSC_TRUE;
784:   PetscFunctionReturn(PETSC_SUCCESS);
785: }

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

790:   Collective

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

797:   Level: beginner

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

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

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

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

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

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

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

842:   Level: developer

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

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

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

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

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

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

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

876:   Level: developer

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