xref: /petsc/src/mat/utils/gcreate.c (revision 647c55fd8ae18e0846d428867ca031dc04cbd941)
1 #include <petsc/private/matimpl.h> /*I "petscmat.h"  I*/
2 
3 #include <../src/mat/impls/aij/seq/aij.h>
4 #include <../src/mat/impls/aij/mpi/mpiaij.h>
5 
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 }
14 
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;
22 
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, ""));
28 
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   }
35 
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   }
44 
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 }
56 
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()`.
61 
62   Collective
63 
64   Input Parameter:
65 . comm - MPI communicator
66 
67   Output Parameter:
68 . A - the matrix
69 
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()`
77 
78    See the manpages for particular formats (e.g., `MATSEQAIJ`)
79    for additional format-specific options.
80 
81   Level: beginner
82 
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.
88 
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;
98 
99   PetscFunctionBegin;
100   PetscAssertPointer(A, 2);
101 
102   *A = NULL;
103   PetscCall(MatInitializePackage());
104 
105   PetscCall(PetscHeaderCreate(B, MAT_CLASSID, "Mat", "Matrix", "Mat", comm, MatDestroy, MatView));
106   PetscCall(PetscLayoutCreate(comm, &B->rmap));
107   PetscCall(PetscLayoutCreate(comm, &B->cmap));
108   PetscCall(PetscStrallocpy(VECSTANDARD, &B->defaultvectype));
109   PetscCall(PetscStrallocpy(PETSCRANDER48, &B->defaultrandtype));
110 
111   B->symmetric                   = PETSC_BOOL3_UNKNOWN;
112   B->hermitian                   = PETSC_BOOL3_UNKNOWN;
113   B->structurally_symmetric      = PETSC_BOOL3_UNKNOWN;
114   B->spd                         = PETSC_BOOL3_UNKNOWN;
115   B->symmetry_eternal            = PETSC_FALSE;
116   B->structural_symmetry_eternal = PETSC_FALSE;
117 
118   B->congruentlayouts = PETSC_DECIDE;
119   B->preallocated     = PETSC_FALSE;
120 #if defined(PETSC_HAVE_DEVICE)
121   B->boundtocpu = PETSC_TRUE;
122 #endif
123   *A = B;
124   PetscFunctionReturn(PETSC_SUCCESS);
125 }
126 
127 /*@C
128   MatCreateFromOptions - Creates a matrix whose type is set from the options database
129 
130   Collective
131 
132   Input Parameters:
133 + comm   - MPI communicator
134 . prefix - [optional] prefix for the options database
135 . bs     - the blocksize (commonly 1)
136 . m      - the local number of rows (or `PETSC_DECIDE`)
137 . n      - the local number of columns (or `PETSC_DECIDE` or `PETSC_DETERMINE`)
138 . M      - the global number of rows (or `PETSC_DETERMINE`)
139 - N      - the global number of columns (or `PETSC_DETERMINE`)
140 
141   Output Parameter:
142 . A - the matrix
143 
144   Options Database Key:
145 . -mat_type - see `MatType`, for example `aij`, `aijcusparse`, `baij`, `sbaij`, dense, defaults to `aij`
146 
147   Level: beginner
148 
149 .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
150           `MatCreateSeqDense()`, `MatCreateDense()`,
151           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
152           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
153           `MatConvert()`, `MatCreate()`
154 @*/
155 PetscErrorCode MatCreateFromOptions(MPI_Comm comm, const char *prefix, PetscInt bs, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *A)
156 {
157   PetscFunctionBegin;
158   PetscAssertPointer(A, 8);
159   PetscCall(MatCreate(comm, A));
160   if (prefix) PetscCall(MatSetOptionsPrefix(*A, prefix));
161   PetscCall(MatSetBlockSize(*A, bs));
162   PetscCall(MatSetSizes(*A, m, n, M, N));
163   PetscCall(MatSetFromOptions(*A));
164   PetscFunctionReturn(PETSC_SUCCESS);
165 }
166 
167 /*@
168   MatSetErrorIfFailure - Causes `Mat` to generate an immediate error, for example a zero pivot, is detected.
169 
170   Logically Collective
171 
172   Input Parameters:
173 + mat - matrix obtained from `MatCreate()`
174 - flg - `PETSC_TRUE` indicates you want the error generated
175 
176   Level: advanced
177 
178   Note:
179   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
180   or result in a `KSPConvergedReason` indicating the method did not converge.
181 
182 .seealso: [](ch_matrices), `Mat`, `PCSetErrorIfFailure()`, `KSPConvergedReason`, `SNESConvergedReason`
183 @*/
184 PetscErrorCode MatSetErrorIfFailure(Mat mat, PetscBool flg)
185 {
186   PetscFunctionBegin;
187   PetscValidHeaderSpecific(mat, MAT_CLASSID, 1);
188   PetscValidLogicalCollectiveBool(mat, flg, 2);
189   mat->erroriffailure = flg;
190   PetscFunctionReturn(PETSC_SUCCESS);
191 }
192 
193 /*@
194   MatSetSizes - Sets the local and global sizes, and checks to determine compatibility
195 
196   Collective
197 
198   Input Parameters:
199 + A - the matrix
200 . m - number of local rows (or `PETSC_DECIDE`)
201 . n - number of local columns (or `PETSC_DECIDE`)
202 . M - number of global rows (or `PETSC_DETERMINE`)
203 - N - number of global columns (or `PETSC_DETERMINE`)
204 
205   Level: beginner
206 
207   Notes:
208   `m` (`n`) and `M` (`N`) cannot be both `PETSC_DECIDE`
209   If one processor calls this with `M` (`N`) of `PETSC_DECIDE` then all processors must, otherwise the program will hang.
210 
211   If `PETSC_DECIDE` is not used for the arguments 'm' and 'n', then the
212   user must ensure that they are chosen to be compatible with the
213   vectors. To do this, one first considers the matrix-vector product
214   'y = A x'. The `m` that is used in the above routine must match the
215   local size used in the vector creation routine `VecCreateMPI()` for 'y'.
216   Likewise, the `n` used must match that used as the local size in
217   `VecCreateMPI()` for 'x'.
218 
219   If `m` and `n` are not `PETSC_DECIDE`, then the values determine the `PetscLayout` of the matrix and the ranges returned by
220   `MatGetOwnershipRange()`,  `MatGetOwnershipRanges()`, `MatGetOwnershipRangeColumn()`, and `MatGetOwnershipRangesColumn()`.
221 
222   You cannot change the sizes once they have been set.
223 
224   The sizes must be set before `MatSetUp()` or MatXXXSetPreallocation() is called.
225 
226 .seealso: [](ch_matrices), `Mat`, `MatGetSize()`, `PetscSplitOwnership()`, `MatGetOwnershipRange()`, `MatGetOwnershipRanges()`,
227           `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRangesColumn()`, `PetscLayout`, `VecSetSizes()`
228 @*/
229 PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N)
230 {
231   PetscFunctionBegin;
232   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
233   PetscValidLogicalCollectiveInt(A, M, 4);
234   PetscValidLogicalCollectiveInt(A, N, 5);
235   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);
236   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);
237   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,
238              A->rmap->n, A->rmap->N);
239   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,
240              A->cmap->n, A->cmap->N);
241   A->rmap->n = m;
242   A->cmap->n = n;
243   A->rmap->N = M > -1 ? M : A->rmap->N;
244   A->cmap->N = N > -1 ? N : A->cmap->N;
245   PetscFunctionReturn(PETSC_SUCCESS);
246 }
247 
248 /*@
249   MatSetFromOptions - Creates a matrix where the type is determined
250   from the options database.
251 
252   Collective
253 
254   Input Parameter:
255 . B - the matrix
256 
257   Options Database Keys:
258 + -mat_type seqaij   - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()`
259 . -mat_type mpiaij   - `MATMPIAIJ` type, uses `MatCreateAIJ()`
260 . -mat_type seqdense - `MATSEQDENSE` type, uses `MatCreateSeqDense()`
261 . -mat_type mpidense - `MATMPIDENSE`, uses `MatCreateDense()`
262 . -mat_type seqbaij  - `MATSEQBAIJ`, uses `MatCreateSeqBAIJ()`
263 - -mat_type mpibaij  - `MATMPIBAIJ`, uses `MatCreateBAIJ()`
264 
265    See the manpages for particular formats (e.g., `MATSEQAIJ`)
266    for additional format-specific options.
267 
268   Level: beginner
269 
270   Notes:
271   Generates a parallel MPI matrix if the communicator has more than one processor.  The default
272   matrix type is `MATAIJ`, using the routines `MatCreateSeqAIJ()` and `MatCreateAIJ()` if you
273   do not select a type in the options database.
274 
275 .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`,
276           `MatCreateSeqDense()`, `MatCreateDense()`,
277           `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`,
278           `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`,
279           `MatConvert()`
280 @*/
281 PetscErrorCode MatSetFromOptions(Mat B)
282 {
283   const char *deft = MATAIJ;
284   char        type[256];
285   PetscBool   flg, set;
286   PetscInt    bind_below = 0;
287 
288   PetscFunctionBegin;
289   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
290 
291   PetscObjectOptionsBegin((PetscObject)B);
292 
293   if (B->rmap->bs < 0) {
294     PetscInt newbs = -1;
295     PetscCall(PetscOptionsInt("-mat_block_size", "Set the blocksize used to store the matrix", "MatSetBlockSize", newbs, &newbs, &flg));
296     if (flg) {
297       PetscCall(PetscLayoutSetBlockSize(B->rmap, newbs));
298       PetscCall(PetscLayoutSetBlockSize(B->cmap, newbs));
299     }
300   }
301 
302   PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, deft, type, 256, &flg));
303   if (flg) {
304     PetscCall(MatSetType(B, type));
305   } else if (!((PetscObject)B)->type_name) {
306     PetscCall(MatSetType(B, deft));
307   }
308 
309   PetscCall(PetscOptionsName("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", &B->checksymmetryonassembly));
310   PetscCall(PetscOptionsReal("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", B->checksymmetrytol, &B->checksymmetrytol, NULL));
311   PetscCall(PetscOptionsBool("-mat_null_space_test", "Checks if provided null space is correct in MatAssemblyEnd()", "MatSetNullSpaceTest", B->checknullspaceonassembly, &B->checknullspaceonassembly, NULL));
312   PetscCall(PetscOptionsBool("-mat_error_if_failure", "Generate an error if an error occurs when factoring the matrix", "MatSetErrorIfFailure", B->erroriffailure, &B->erroriffailure, NULL));
313 
314   PetscTryTypeMethod(B, setfromoptions, PetscOptionsObject);
315 
316   flg = PETSC_FALSE;
317   PetscCall(PetscOptionsBool("-mat_new_nonzero_location_err", "Generate an error if new nonzeros are created in the matrix structure (useful to test preallocation)", "MatSetOption", flg, &flg, &set));
318   if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, flg));
319   flg = PETSC_FALSE;
320   PetscCall(PetscOptionsBool("-mat_new_nonzero_allocation_err", "Generate an error if new nonzeros are allocated in the matrix structure (useful to test preallocation)", "MatSetOption", flg, &flg, &set));
321   if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, flg));
322   flg = PETSC_FALSE;
323   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));
324   if (set) PetscCall(MatSetOption(B, MAT_IGNORE_ZERO_ENTRIES, flg));
325 
326   flg = PETSC_FALSE;
327   PetscCall(PetscOptionsBool("-mat_form_explicit_transpose", "Hint to form an explicit transpose for operations like MatMultTranspose", "MatSetOption", flg, &flg, &set));
328   if (set) PetscCall(MatSetOption(B, MAT_FORM_EXPLICIT_TRANSPOSE, flg));
329 
330   /* Bind to CPU if below a user-specified size threshold.
331    * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types,
332    * and putting it here makes is more maintainable than duplicating this for all. */
333   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));
334   if (flg && B->rmap->n < bind_below) PetscCall(MatBindToCPU(B, PETSC_TRUE));
335 
336   /* process any options handlers added with PetscObjectAddOptionsHandler() */
337   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)B, PetscOptionsObject));
338   PetscOptionsEnd();
339   PetscFunctionReturn(PETSC_SUCCESS);
340 }
341 
342 /*@C
343   MatXAIJSetPreallocation - set preallocation for serial and parallel `MATAIJ`, `MATBAIJ`, and `MATSBAIJ` matrices and their unassembled versions.
344 
345   Collective
346 
347   Input Parameters:
348 + A     - matrix being preallocated
349 . bs    - block size
350 . dnnz  - number of nonzero column blocks per block row of diagonal part of parallel matrix
351 . onnz  - number of nonzero column blocks per block row of off-diagonal part of parallel matrix
352 . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix
353 - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix
354 
355   Level: beginner
356 
357 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`,
358           `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`,
359           `PetscSplitOwnership()`
360 @*/
361 PetscErrorCode MatXAIJSetPreallocation(Mat A, PetscInt bs, const PetscInt dnnz[], const PetscInt onnz[], const PetscInt dnnzu[], const PetscInt onnzu[])
362 {
363   PetscInt cbs;
364   void (*aij)(void);
365   void (*is)(void);
366   void (*hyp)(void) = NULL;
367 
368   PetscFunctionBegin;
369   if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */
370     PetscCall(MatSetBlockSize(A, bs));
371   }
372   PetscCall(PetscLayoutSetUp(A->rmap));
373   PetscCall(PetscLayoutSetUp(A->cmap));
374   PetscCall(MatGetBlockSizes(A, &bs, &cbs));
375   /* these routines assumes bs == cbs, this should be checked somehow */
376   PetscCall(MatSeqBAIJSetPreallocation(A, bs, 0, dnnz));
377   PetscCall(MatMPIBAIJSetPreallocation(A, bs, 0, dnnz, 0, onnz));
378   PetscCall(MatSeqSBAIJSetPreallocation(A, bs, 0, dnnzu));
379   PetscCall(MatMPISBAIJSetPreallocation(A, bs, 0, dnnzu, 0, onnzu));
380   /*
381     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
382     good before going on with it.
383   */
384   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatMPIAIJSetPreallocation_C", &aij));
385   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatISSetPreallocation_C", &is));
386 #if defined(PETSC_HAVE_HYPRE)
387   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatHYPRESetPreallocation_C", &hyp));
388 #endif
389   if (!aij && !is && !hyp) PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", &aij));
390   if (aij || is || hyp) {
391     if (bs == cbs && bs == 1) {
392       PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz));
393       PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz, 0, onnz));
394       PetscCall(MatISSetPreallocation(A, 0, dnnz, 0, onnz));
395 #if defined(PETSC_HAVE_HYPRE)
396       PetscCall(MatHYPRESetPreallocation(A, 0, dnnz, 0, onnz));
397 #endif
398     } else { /* Convert block-row precallocation to scalar-row */
399       PetscInt i, m, *sdnnz, *sonnz;
400       PetscCall(MatGetLocalSize(A, &m, NULL));
401       PetscCall(PetscMalloc2((!!dnnz) * m, &sdnnz, (!!onnz) * m, &sonnz));
402       for (i = 0; i < m; i++) {
403         if (dnnz) sdnnz[i] = dnnz[i / bs] * cbs;
404         if (onnz) sonnz[i] = onnz[i / bs] * cbs;
405       }
406       PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL));
407       PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
408       PetscCall(MatISSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
409 #if defined(PETSC_HAVE_HYPRE)
410       PetscCall(MatHYPRESetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL));
411 #endif
412       PetscCall(PetscFree2(sdnnz, sonnz));
413     }
414   }
415   PetscFunctionReturn(PETSC_SUCCESS);
416 }
417 
418 /*@C
419   MatHeaderMerge - Merges some information from the header of `C` to `A`; the `C` object is then destroyed
420 
421   Collective, No Fortran Support
422 
423   Input Parameters:
424 + A - a `Mat` being merged into
425 - C - the `Mat` providing the merge information
426 
427   Level: developer
428 
429   Notes:
430   `A` and `C` must be of the same type.
431   The object list and query function list in `A` are retained, as well as the object name, and prefix.
432   The object state of `A` is increased by 1.
433 
434   Developer Note:
435   This is somewhat different from `MatHeaderReplace()`, it would be nice to merge the code
436 
437 .seealso: `Mat`, `MatHeaderReplace()`
438  @*/
439 PetscErrorCode MatHeaderMerge(Mat A, Mat *C)
440 {
441   PetscInt          refct;
442   PetscOps          Abops;
443   struct _MatOps    Aops;
444   char             *mtype, *mname, *mprefix;
445   Mat_Product      *product;
446   Mat_Redundant    *redundant;
447   PetscObjectState  state;
448   PetscObjectList   olist;
449   PetscFunctionList qlist;
450 
451   PetscFunctionBegin;
452   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
453   PetscValidHeaderSpecific(*C, MAT_CLASSID, 2);
454   if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
455   PetscCheckSameTypeAndComm(A, 1, *C, 2);
456   /* save the parts of A we need */
457   Abops     = ((PetscObject)A)->bops[0];
458   Aops      = A->ops[0];
459   refct     = ((PetscObject)A)->refct;
460   mtype     = ((PetscObject)A)->type_name;
461   mname     = ((PetscObject)A)->name;
462   state     = ((PetscObject)A)->state;
463   mprefix   = ((PetscObject)A)->prefix;
464   product   = A->product;
465   redundant = A->redundant;
466   qlist     = ((PetscObject)A)->qlist;
467   olist     = ((PetscObject)A)->olist;
468 
469   /* zero these so the destroy below does not free them */
470   ((PetscObject)A)->type_name = NULL;
471   ((PetscObject)A)->name      = NULL;
472   ((PetscObject)A)->qlist     = NULL;
473   ((PetscObject)A)->olist     = NULL;
474 
475   /*
476      free all the interior data structures from mat
477      cannot use PetscUseTypeMethod(A,destroy); because compiler
478      thinks it may print NULL type_name and name
479   */
480   PetscTryTypeMethod(A, destroy);
481 
482   PetscCall(PetscFree(A->defaultvectype));
483   PetscCall(PetscFree(A->defaultrandtype));
484   PetscCall(PetscLayoutDestroy(&A->rmap));
485   PetscCall(PetscLayoutDestroy(&A->cmap));
486   PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A));
487 
488   /* copy C over to A */
489   PetscCall(PetscFree(A->factorprefix));
490   PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));
491 
492   /* return the parts of A we saved */
493   ((PetscObject)A)->bops[0]   = Abops;
494   A->ops[0]                   = Aops;
495   ((PetscObject)A)->refct     = refct;
496   ((PetscObject)A)->type_name = mtype;
497   ((PetscObject)A)->name      = mname;
498   ((PetscObject)A)->prefix    = mprefix;
499   ((PetscObject)A)->state     = state + 1;
500   A->product                  = product;
501   A->redundant                = redundant;
502 
503   /* Append the saved lists */
504   PetscCall(PetscFunctionListDuplicate(qlist, &((PetscObject)A)->qlist));
505   PetscCall(PetscObjectListDuplicate(olist, &((PetscObject)A)->olist));
506   PetscCall(PetscFunctionListDestroy(&qlist));
507   PetscCall(PetscObjectListDestroy(&olist));
508 
509   /* since these two are copied into A we do not want them destroyed in C */
510   ((PetscObject)*C)->qlist = NULL;
511   ((PetscObject)*C)->olist = NULL;
512   PetscCall(PetscHeaderDestroy(C));
513   PetscFunctionReturn(PETSC_SUCCESS);
514 }
515 
516 /*@
517   MatHeaderReplace - Replaces the internal data of matrix `A` by the internal data of matrix `C` while deleting the outer wrapper of `C`
518 
519   Input Parameters:
520 + A - a `Mat` whose internal data is to be replaced
521 - C - the `Mat` providing new internal data for `A`
522 
523   Level: advanced
524 
525   Example Usage\:
526 .vb
527   Mat C;
528   MatCreateSeqAIJWithArrays(..., &C);
529   MatHeaderReplace(A, &C);
530   // C has been destroyed and A contains the matrix entries of C
531 .ve
532 
533   Note:
534   This can be used inside a function provided to `SNESSetJacobian()`, `TSSetRHSJacobian()`, or `TSSetIJacobian()` in cases where the user code computes an entirely new sparse matrix
535   (generally with a different nonzero pattern) for each Newton update. It is usually better to reuse the matrix structure of `A` instead of constructing an entirely new one.
536 
537   Developer Note:
538   This is somewhat different from `MatHeaderMerge()` it would be nice to merge the code
539 
540 .seealso: `Mat`, `MatHeaderMerge()`
541  @*/
542 PetscErrorCode MatHeaderReplace(Mat A, Mat *C)
543 {
544   PetscInt         refct;
545   PetscObjectState state;
546   struct _p_Mat    buffer;
547   MatStencilInfo   stencil;
548 
549   PetscFunctionBegin;
550   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
551   PetscValidHeaderSpecific(*C, MAT_CLASSID, 2);
552   if (A == *C) PetscFunctionReturn(PETSC_SUCCESS);
553   PetscCheckSameComm(A, 1, *C, 2);
554   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);
555 
556   /* swap C and A */
557   refct   = ((PetscObject)A)->refct;
558   state   = ((PetscObject)A)->state;
559   stencil = A->stencil;
560   PetscCall(PetscMemcpy(&buffer, A, sizeof(struct _p_Mat)));
561   PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat)));
562   PetscCall(PetscMemcpy(*C, &buffer, sizeof(struct _p_Mat)));
563   ((PetscObject)A)->refct = refct;
564   ((PetscObject)A)->state = state + 1;
565   A->stencil              = stencil;
566 
567   ((PetscObject)*C)->refct = 1;
568   PetscCall(MatDestroy(C));
569   PetscFunctionReturn(PETSC_SUCCESS);
570 }
571 
572 /*@
573   MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU
574 
575   Logically Collective
576 
577   Input Parameters:
578 + A   - the matrix
579 - flg - bind to the CPU if value of `PETSC_TRUE`
580 
581   Level: intermediate
582 
583 .seealso: [](ch_matrices), `Mat`, `MatBoundToCPU()`
584 @*/
585 PetscErrorCode MatBindToCPU(Mat A, PetscBool flg)
586 {
587   PetscFunctionBegin;
588   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
589   PetscValidLogicalCollectiveBool(A, flg, 2);
590 #if defined(PETSC_HAVE_DEVICE)
591   if (A->boundtocpu == flg) PetscFunctionReturn(PETSC_SUCCESS);
592   A->boundtocpu = flg;
593   PetscTryTypeMethod(A, bindtocpu, flg);
594 #endif
595   PetscFunctionReturn(PETSC_SUCCESS);
596 }
597 
598 /*@
599   MatBoundToCPU - query if a matrix is bound to the CPU
600 
601   Input Parameter:
602 . A - the matrix
603 
604   Output Parameter:
605 . flg - the logical flag
606 
607   Level: intermediate
608 
609 .seealso: [](ch_matrices), `Mat`, `MatBindToCPU()`
610 @*/
611 PetscErrorCode MatBoundToCPU(Mat A, PetscBool *flg)
612 {
613   PetscFunctionBegin;
614   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
615   PetscAssertPointer(flg, 2);
616 #if defined(PETSC_HAVE_DEVICE)
617   *flg = A->boundtocpu;
618 #else
619   *flg = PETSC_TRUE;
620 #endif
621   PetscFunctionReturn(PETSC_SUCCESS);
622 }
623 
624 PetscErrorCode MatSetValuesCOO_Basic(Mat A, const PetscScalar coo_v[], InsertMode imode)
625 {
626   IS              is_coo_i, is_coo_j;
627   const PetscInt *coo_i, *coo_j;
628   PetscInt        n, n_i, n_j;
629   PetscScalar     zero = 0.;
630 
631   PetscFunctionBegin;
632   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_i", (PetscObject *)&is_coo_i));
633   PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_j", (PetscObject *)&is_coo_j));
634   PetscCheck(is_coo_i, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_i IS");
635   PetscCheck(is_coo_j, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_j IS");
636   PetscCall(ISGetLocalSize(is_coo_i, &n_i));
637   PetscCall(ISGetLocalSize(is_coo_j, &n_j));
638   PetscCheck(n_i == n_j, PETSC_COMM_SELF, PETSC_ERR_COR, "Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT, n_i, n_j);
639   PetscCall(ISGetIndices(is_coo_i, &coo_i));
640   PetscCall(ISGetIndices(is_coo_j, &coo_j));
641   if (imode != ADD_VALUES) PetscCall(MatZeroEntries(A));
642   for (n = 0; n < n_i; n++) PetscCall(MatSetValue(A, coo_i[n], coo_j[n], coo_v ? coo_v[n] : zero, ADD_VALUES));
643   PetscCall(ISRestoreIndices(is_coo_i, &coo_i));
644   PetscCall(ISRestoreIndices(is_coo_j, &coo_j));
645   PetscFunctionReturn(PETSC_SUCCESS);
646 }
647 
648 PetscErrorCode MatSetPreallocationCOO_Basic(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
649 {
650   Mat         preallocator;
651   IS          is_coo_i, is_coo_j;
652   PetscScalar zero = 0.0;
653 
654   PetscFunctionBegin;
655   PetscCall(PetscLayoutSetUp(A->rmap));
656   PetscCall(PetscLayoutSetUp(A->cmap));
657   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator));
658   PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
659   PetscCall(MatSetSizes(preallocator, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
660   PetscCall(MatSetLayouts(preallocator, A->rmap, A->cmap));
661   PetscCall(MatSetUp(preallocator));
662   for (PetscCount n = 0; n < ncoo; n++) PetscCall(MatSetValue(preallocator, coo_i[n], coo_j[n], zero, INSERT_VALUES));
663   PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
664   PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
665   PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A));
666   PetscCall(MatDestroy(&preallocator));
667   PetscCheck(ncoo <= PETSC_MAX_INT, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "ncoo %" PetscCount_FMT " overflowed PetscInt; configure --with-64-bit-indices or request support", ncoo);
668   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo, coo_i, PETSC_COPY_VALUES, &is_coo_i));
669   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo, coo_j, PETSC_COPY_VALUES, &is_coo_j));
670   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_i", (PetscObject)is_coo_i));
671   PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_j", (PetscObject)is_coo_j));
672   PetscCall(ISDestroy(&is_coo_i));
673   PetscCall(ISDestroy(&is_coo_j));
674   PetscFunctionReturn(PETSC_SUCCESS);
675 }
676 
677 /*@C
678   MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices
679 
680   Collective
681 
682   Input Parameters:
683 + A     - matrix being preallocated
684 . ncoo  - number of entries
685 . coo_i - row indices
686 - coo_j - column indices
687 
688   Level: beginner
689 
690   Notes:
691   The indices within `coo_i` and `coo_j` may be modified within this function. The caller should not rely on them
692   having any specific value after this function returns. The arrays can be freed or reused immediately
693   after this function returns.
694 
695   Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed
696   but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries
697   are allowed and will be properly added or inserted to the matrix, unless the matrix option `MAT_IGNORE_OFF_PROC_ENTRIES`
698   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.
699 
700   If you just want to create a sequential AIJ matrix (`MATSEQAIJ`), and your matrix entries in COO format are unique, you can also use
701   `MatCreateSeqAIJFromTriple()`. But that is not recommended for iterative applications.
702 
703 .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
704           `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`,
705           `DMSetMatrixPreallocateSkip()`, `MatCreateSeqAIJFromTriple()`
706 @*/
707 PetscErrorCode MatSetPreallocationCOO(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
708 {
709   PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;
710 
711   PetscFunctionBegin;
712   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
713   PetscValidType(A, 1);
714   if (ncoo) PetscAssertPointer(coo_i, 3);
715   if (ncoo) PetscAssertPointer(coo_j, 4);
716   PetscCall(PetscLayoutSetUp(A->rmap));
717   PetscCall(PetscLayoutSetUp(A->cmap));
718   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOO_C", &f));
719 
720   PetscCall(PetscLogEventBegin(MAT_PreallCOO, A, 0, 0, 0));
721   if (f) {
722     PetscCall((*f)(A, ncoo, coo_i, coo_j));
723   } else { /* allow fallback, very slow */
724     PetscCall(MatSetPreallocationCOO_Basic(A, ncoo, coo_i, coo_j));
725   }
726   PetscCall(PetscLogEventEnd(MAT_PreallCOO, A, 0, 0, 0));
727   A->preallocated = PETSC_TRUE;
728   A->nonzerostate++;
729   PetscFunctionReturn(PETSC_SUCCESS);
730 }
731 
732 /*@C
733   MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices
734 
735   Collective
736 
737   Input Parameters:
738 + A     - matrix being preallocated
739 . ncoo  - number of entries
740 . coo_i - row indices (local numbering; may be modified)
741 - coo_j - column indices (local numbering; may be modified)
742 
743   Level: beginner
744 
745   Notes:
746   The local indices are translated using the local to global mapping, thus `MatSetLocalToGlobalMapping()` must have been
747   called prior to this function. For matrices created with `DMCreateMatrix()` the local to global mapping is often already provided.
748 
749   The indices `coo_i` and `coo_j` may be modified within this function. They might be translated to corresponding global
750   indices, but the caller should not rely on them having any specific value after this function returns. The arrays
751   can be freed or reused immediately after this function returns.
752 
753   Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed
754   but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries
755   are allowed and will be properly added or inserted to the matrix.
756 
757 .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`,
758           `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`,
759           `DMSetMatrixPreallocateSkip()`
760 @*/
761 PetscErrorCode MatSetPreallocationCOOLocal(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[])
762 {
763   PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL;
764 
765   PetscFunctionBegin;
766   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
767   PetscValidType(A, 1);
768   if (ncoo) PetscAssertPointer(coo_i, 3);
769   if (ncoo) PetscAssertPointer(coo_j, 4);
770   PetscCheck(ncoo <= PETSC_MAX_INT, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "ncoo %" PetscCount_FMT " overflowed PetscInt; configure --with-64-bit-indices or request support", ncoo);
771   PetscCall(PetscLayoutSetUp(A->rmap));
772   PetscCall(PetscLayoutSetUp(A->cmap));
773 
774   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", &f));
775   if (f) {
776     PetscCall((*f)(A, ncoo, coo_i, coo_j));
777     A->nonzerostate++;
778   } else {
779     ISLocalToGlobalMapping ltog_row, ltog_col;
780     PetscCall(MatGetLocalToGlobalMapping(A, &ltog_row, &ltog_col));
781     if (ltog_row) PetscCall(ISLocalToGlobalMappingApply(ltog_row, ncoo, coo_i, coo_i));
782     if (ltog_col) PetscCall(ISLocalToGlobalMappingApply(ltog_col, ncoo, coo_j, coo_j));
783     PetscCall(MatSetPreallocationCOO(A, ncoo, coo_i, coo_j));
784   }
785   A->preallocated = PETSC_TRUE;
786   PetscFunctionReturn(PETSC_SUCCESS);
787 }
788 
789 /*@
790   MatSetValuesCOO - set values at once in a matrix preallocated using `MatSetPreallocationCOO()`
791 
792   Collective
793 
794   Input Parameters:
795 + A     - matrix being preallocated
796 . coo_v - the matrix values (can be `NULL`)
797 - imode - the insert mode
798 
799   Level: beginner
800 
801   Notes:
802   The values must follow the order of the indices prescribed with `MatSetPreallocationCOO()` or `MatSetPreallocationCOOLocal()`.
803 
804   When repeated entries are specified in the COO indices the `coo_v` values are first properly summed, regardless of the value of imode.
805   The imode flag indicates if coo_v must be added to the current values of the matrix (`ADD_VALUES`) or overwritten (`INSERT_VALUES`).
806 
807   `MatAssemblyBegin()` and `MatAssemblyEnd()` do not need to be called after this routine. It automatically handles the assembly process.
808 
809 .seealso: [](ch_matrices), `Mat`, `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`
810 @*/
811 PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode)
812 {
813   PetscErrorCode (*f)(Mat, const PetscScalar[], InsertMode) = NULL;
814   PetscBool oldFlg;
815 
816   PetscFunctionBegin;
817   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
818   PetscValidType(A, 1);
819   MatCheckPreallocated(A, 1);
820   PetscValidLogicalCollectiveEnum(A, imode, 3);
821   PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetValuesCOO_C", &f));
822   PetscCall(PetscLogEventBegin(MAT_SetVCOO, A, 0, 0, 0));
823   if (f) {
824     PetscCall((*f)(A, coo_v, imode)); // all known COO implementations do not use MatStash. They do their own off-proc communication
825     PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &oldFlg));
826     PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); // set A->nooffprocentries to avoid costly MatStash scatter in MatAssembly
827   } else {
828     PetscCall(MatSetValuesCOO_Basic(A, coo_v, imode)); // fall back to MatSetValues, which might use MatStash
829   }
830   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
831   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
832   if (f) PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, oldFlg));
833   PetscCall(PetscLogEventEnd(MAT_SetVCOO, A, 0, 0, 0));
834   PetscFunctionReturn(PETSC_SUCCESS);
835 }
836 
837 /*@
838   MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects
839 
840   Input Parameters:
841 + A   - the matrix
842 - flg - flag indicating whether the boundtocpu flag should be propagated
843 
844   Level: developer
845 
846   Notes:
847   If the value of flg is set to true, the following will occur
848 +   `MatCreateSubMatrices()` and `MatCreateRedundantMatrix()` - bind created matrices to CPU if the input matrix is bound to the CPU.
849 -   `MatCreateVecs()` - bind created vectors to CPU if the input matrix is bound to the CPU.
850 
851   The bindingpropagates flag itself is also propagated by the above routines.
852 
853   Developer Notes:
854   If the fine-scale `DMDA` has the `-dm_bind_below` option set to true, then `DMCreateInterpolationScale()` calls `MatSetBindingPropagates()`
855   on the restriction/interpolation operator to set the bindingpropagates flag to true.
856 
857 .seealso: [](ch_matrices), `Mat`, `VecSetBindingPropagates()`, `MatGetBindingPropagates()`
858 @*/
859 PetscErrorCode MatSetBindingPropagates(Mat A, PetscBool flg)
860 {
861   PetscFunctionBegin;
862   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
863 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
864   A->bindingpropagates = flg;
865 #endif
866   PetscFunctionReturn(PETSC_SUCCESS);
867 }
868 
869 /*@
870   MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects
871 
872   Input Parameter:
873 . A - the matrix
874 
875   Output Parameter:
876 . flg - flag indicating whether the boundtocpu flag will be propagated
877 
878   Level: developer
879 
880 .seealso: [](ch_matrices), `Mat`, `MatSetBindingPropagates()`
881 @*/
882 PetscErrorCode MatGetBindingPropagates(Mat A, PetscBool *flg)
883 {
884   PetscFunctionBegin;
885   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
886   PetscAssertPointer(flg, 2);
887 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
888   *flg = A->bindingpropagates;
889 #else
890   *flg = PETSC_FALSE;
891 #endif
892   PetscFunctionReturn(PETSC_SUCCESS);
893 }
894