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 PetscCall(MatInitializePackage()); 102 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)); 108 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; 115 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 } 124 125 /*@ 126 MatCreateFromOptions - Creates a matrix whose type is set from the options database 127 128 Collective 129 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`) 138 139 Output Parameter: 140 . A - the matrix 141 142 Options Database Key: 143 . -mat_type - see `MatType`, for example `aij`, `aijcusparse`, `baij`, `sbaij`, dense, defaults to `aij` 144 145 Level: beginner 146 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 } 164 165 /*@ 166 MatSetErrorIfFailure - Causes `Mat` to generate an immediate error, for example a zero pivot, is detected. 167 168 Logically Collective 169 170 Input Parameters: 171 + mat - matrix obtained from `MatCreate()` 172 - flg - `PETSC_TRUE` indicates you want the error generated 173 174 Level: advanced 175 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. 179 180 .seealso: [](ch_matrices), `Mat`, `PCSetErrorIfFailure()`, `KSPConvergedReason`, `SNESConvergedReason` 181 @*/ 182 PetscErrorCode MatSetErrorIfFailure(Mat mat, PetscBool flg) 183 { 184 PetscFunctionBegin; 185 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 186 PetscValidLogicalCollectiveBool(mat, flg, 2); 187 mat->erroriffailure = flg; 188 PetscFunctionReturn(PETSC_SUCCESS); 189 } 190 191 /*@ 192 MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 193 194 Collective 195 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`) 202 203 Level: beginner 204 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. 208 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'. 214 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()`. 217 218 You cannot change the sizes once they have been set. 219 220 The sizes must be set before `MatSetUp()` or MatXXXSetPreallocation() is called. 221 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; 228 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 229 PetscValidLogicalCollectiveInt(A, M, 4); 230 PetscValidLogicalCollectiveInt(A, N, 5); 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 } 243 244 /*@ 245 MatSetFromOptions - Creates a matrix where the type is determined 246 from the options database. 247 248 Collective 249 250 Input Parameter: 251 . B - the matrix 252 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()` 260 261 See the manpages for particular formats (e.g., `MATSEQAIJ`) 262 for additional format-specific options. 263 264 Level: beginner 265 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. 270 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; 283 284 PetscFunctionBegin; 285 PetscValidHeaderSpecific(B, MAT_CLASSID, 1); 286 287 PetscObjectOptionsBegin((PetscObject)B); 288 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 } 294 295 PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, deft, type, 256, &flg)); 296 if (flg) { 297 PetscCall(MatSetType(B, type)); 298 } else if (!((PetscObject)B)->type_name) { 299 PetscCall(MatSetType(B, deft)); 300 } 301 302 PetscCall(PetscOptionsName("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", &B->checksymmetryonassembly)); 303 PetscCall(PetscOptionsReal("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", B->checksymmetrytol, &B->checksymmetrytol, NULL)); 304 PetscCall(PetscOptionsBool("-mat_null_space_test", "Checks if provided null space is correct in MatAssemblyEnd()", "MatSetNullSpaceTest", B->checknullspaceonassembly, &B->checknullspaceonassembly, NULL)); 305 PetscCall(PetscOptionsBool("-mat_error_if_failure", "Generate an error if an error occurs when factoring the matrix", "MatSetErrorIfFailure", B->erroriffailure, &B->erroriffailure, NULL)); 306 307 PetscTryTypeMethod(B, setfromoptions, PetscOptionsObject); 308 309 flg = PETSC_FALSE; 310 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)); 311 if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, flg)); 312 flg = PETSC_FALSE; 313 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)); 314 if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, flg)); 315 flg = PETSC_FALSE; 316 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)); 317 if (set) PetscCall(MatSetOption(B, MAT_IGNORE_ZERO_ENTRIES, flg)); 318 319 flg = PETSC_FALSE; 320 PetscCall(PetscOptionsBool("-mat_form_explicit_transpose", "Hint to form an explicit transpose for operations like MatMultTranspose", "MatSetOption", flg, &flg, &set)); 321 if (set) PetscCall(MatSetOption(B, MAT_FORM_EXPLICIT_TRANSPOSE, flg)); 322 323 /* Bind to CPU if below a user-specified size threshold. 324 * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types, 325 * and putting it here makes is more maintainable than duplicating this for all. */ 326 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)); 327 if (flg && B->rmap->n < bind_below) PetscCall(MatBindToCPU(B, PETSC_TRUE)); 328 329 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 330 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)B, PetscOptionsObject)); 331 PetscOptionsEnd(); 332 PetscFunctionReturn(PETSC_SUCCESS); 333 } 334 335 /*@ 336 MatXAIJSetPreallocation - set preallocation for serial and parallel `MATAIJ`, `MATBAIJ`, and `MATSBAIJ` matrices and their unassembled versions. 337 338 Collective 339 340 Input Parameters: 341 + A - matrix being preallocated 342 . bs - block size 343 . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 344 . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 345 . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 346 - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 347 348 Level: beginner 349 350 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, 351 `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, 352 `PetscSplitOwnership()` 353 @*/ 354 PetscErrorCode MatXAIJSetPreallocation(Mat A, PetscInt bs, const PetscInt dnnz[], const PetscInt onnz[], const PetscInt dnnzu[], const PetscInt onnzu[]) 355 { 356 PetscInt cbs; 357 PetscBool aij, is, hyp; 358 359 PetscFunctionBegin; 360 if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 361 PetscCall(MatSetBlockSize(A, bs)); 362 } 363 PetscCall(PetscLayoutSetUp(A->rmap)); 364 PetscCall(PetscLayoutSetUp(A->cmap)); 365 PetscCall(MatGetBlockSizes(A, &bs, &cbs)); 366 /* these routines assumes bs == cbs, this should be checked somehow */ 367 PetscCall(MatSeqBAIJSetPreallocation(A, bs, 0, dnnz)); 368 PetscCall(MatMPIBAIJSetPreallocation(A, bs, 0, dnnz, 0, onnz)); 369 PetscCall(MatSeqSBAIJSetPreallocation(A, bs, 0, dnnzu)); 370 PetscCall(MatMPISBAIJSetPreallocation(A, bs, 0, dnnzu, 0, onnzu)); 371 /* 372 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 373 good before going on with it. 374 */ 375 PetscCall(PetscObjectHasFunction((PetscObject)A, "MatMPIAIJSetPreallocation_C", &aij)); 376 PetscCall(PetscObjectHasFunction((PetscObject)A, "MatISSetPreallocation_C", &is)); 377 PetscCall(PetscObjectHasFunction((PetscObject)A, "MatHYPRESetPreallocation_C", &hyp)); 378 if (!aij && !is && !hyp) PetscCall(PetscObjectHasFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", &aij)); 379 if (aij || is || hyp) { 380 if (bs == cbs && bs == 1) { 381 PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz)); 382 PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz, 0, onnz)); 383 PetscCall(MatISSetPreallocation(A, 0, dnnz, 0, onnz)); 384 #if defined(PETSC_HAVE_HYPRE) 385 PetscCall(MatHYPRESetPreallocation(A, 0, dnnz, 0, onnz)); 386 #endif 387 } else { /* Convert block-row precallocation to scalar-row */ 388 PetscInt i, m, *sdnnz, *sonnz; 389 PetscCall(MatGetLocalSize(A, &m, NULL)); 390 PetscCall(PetscMalloc2((!!dnnz) * m, &sdnnz, (!!onnz) * m, &sonnz)); 391 for (i = 0; i < m; i++) { 392 if (dnnz) sdnnz[i] = dnnz[i / bs] * cbs; 393 if (onnz) sonnz[i] = onnz[i / bs] * cbs; 394 } 395 PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL)); 396 PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 397 PetscCall(MatISSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 398 #if defined(PETSC_HAVE_HYPRE) 399 PetscCall(MatHYPRESetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 400 #endif 401 PetscCall(PetscFree2(sdnnz, sonnz)); 402 } 403 } 404 PetscFunctionReturn(PETSC_SUCCESS); 405 } 406 407 /*@C 408 MatHeaderMerge - Merges some information from the header of `C` to `A`; the `C` object is then destroyed 409 410 Collective, No Fortran Support 411 412 Input Parameters: 413 + A - a `Mat` being merged into 414 - C - the `Mat` providing the merge information 415 416 Level: developer 417 418 Notes: 419 `A` and `C` must be of the same type. 420 The object list and query function list in `A` are retained, as well as the object name, and prefix. 421 The object state of `A` is increased by 1. 422 423 Developer Note: 424 This is somewhat different from `MatHeaderReplace()`, it would be nice to merge the code 425 426 .seealso: `Mat`, `MatHeaderReplace()` 427 @*/ 428 PetscErrorCode MatHeaderMerge(Mat A, Mat *C) 429 { 430 PetscInt refct; 431 PetscOps Abops; 432 struct _MatOps Aops; 433 char *mtype, *mname, *mprefix; 434 Mat_Product *product; 435 Mat_Redundant *redundant; 436 PetscObjectState state; 437 PetscObjectList olist; 438 PetscFunctionList qlist; 439 440 PetscFunctionBegin; 441 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 442 PetscValidHeaderSpecific(*C, MAT_CLASSID, 2); 443 if (A == *C) PetscFunctionReturn(PETSC_SUCCESS); 444 PetscCheckSameTypeAndComm(A, 1, *C, 2); 445 /* save the parts of A we need */ 446 Abops = ((PetscObject)A)->bops[0]; 447 Aops = A->ops[0]; 448 refct = ((PetscObject)A)->refct; 449 mtype = ((PetscObject)A)->type_name; 450 mname = ((PetscObject)A)->name; 451 state = ((PetscObject)A)->state; 452 mprefix = ((PetscObject)A)->prefix; 453 product = A->product; 454 redundant = A->redundant; 455 qlist = ((PetscObject)A)->qlist; 456 olist = ((PetscObject)A)->olist; 457 458 /* zero these so the destroy below does not free them */ 459 ((PetscObject)A)->type_name = NULL; 460 ((PetscObject)A)->name = NULL; 461 ((PetscObject)A)->qlist = NULL; 462 ((PetscObject)A)->olist = NULL; 463 464 /* 465 free all the interior data structures from mat 466 cannot use PetscUseTypeMethod(A,destroy); because compiler 467 thinks it may print NULL type_name and name 468 */ 469 PetscTryTypeMethod(A, destroy); 470 471 PetscCall(PetscFree(A->defaultvectype)); 472 PetscCall(PetscFree(A->defaultrandtype)); 473 PetscCall(PetscLayoutDestroy(&A->rmap)); 474 PetscCall(PetscLayoutDestroy(&A->cmap)); 475 PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A)); 476 477 /* copy C over to A */ 478 PetscCall(PetscFree(A->factorprefix)); 479 PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat))); 480 481 /* return the parts of A we saved */ 482 ((PetscObject)A)->bops[0] = Abops; 483 A->ops[0] = Aops; 484 ((PetscObject)A)->refct = refct; 485 ((PetscObject)A)->type_name = mtype; 486 ((PetscObject)A)->name = mname; 487 ((PetscObject)A)->prefix = mprefix; 488 ((PetscObject)A)->state = state + 1; 489 A->product = product; 490 A->redundant = redundant; 491 492 /* Append the saved lists */ 493 PetscCall(PetscFunctionListDuplicate(qlist, &((PetscObject)A)->qlist)); 494 PetscCall(PetscObjectListDuplicate(olist, &((PetscObject)A)->olist)); 495 PetscCall(PetscFunctionListDestroy(&qlist)); 496 PetscCall(PetscObjectListDestroy(&olist)); 497 498 /* since these two are copied into A we do not want them destroyed in C */ 499 ((PetscObject)*C)->qlist = NULL; 500 ((PetscObject)*C)->olist = NULL; 501 PetscCall(PetscHeaderDestroy(C)); 502 PetscFunctionReturn(PETSC_SUCCESS); 503 } 504 505 /*@ 506 MatHeaderReplace - Replaces the internal data of matrix `A` by the internal data of matrix `C` while deleting the outer wrapper of `C` 507 508 Input Parameters: 509 + A - a `Mat` whose internal data is to be replaced 510 - C - the `Mat` providing new internal data for `A` 511 512 Level: advanced 513 514 Example Usage\: 515 .vb 516 Mat C; 517 MatCreateSeqAIJWithArrays(..., &C); 518 MatHeaderReplace(A, &C); 519 // C has been destroyed and A contains the matrix entries of C 520 .ve 521 522 Note: 523 This can be used inside a function provided to `SNESSetJacobian()`, `TSSetRHSJacobian()`, or `TSSetIJacobian()` in cases where the user code 524 computes an entirely new sparse matrix (generally with a different matrix nonzero structure/pattern) for each Newton update. 525 It is usually better to reuse the matrix nonzero structure of `A` instead of constructing an entirely new one. 526 527 Developer Note: 528 This is somewhat different from `MatHeaderMerge()` it would be nice to merge the code 529 530 .seealso: `Mat`, `MatHeaderMerge()` 531 @*/ 532 PetscErrorCode MatHeaderReplace(Mat A, Mat *C) 533 { 534 PetscInt refct; 535 PetscObjectState state; 536 struct _p_Mat buffer; 537 MatStencilInfo stencil; 538 539 PetscFunctionBegin; 540 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 541 PetscValidHeaderSpecific(*C, MAT_CLASSID, 2); 542 if (A == *C) PetscFunctionReturn(PETSC_SUCCESS); 543 PetscCheckSameComm(A, 1, *C, 2); 544 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 546 /* swap C and A */ 547 refct = ((PetscObject)A)->refct; 548 state = ((PetscObject)A)->state; 549 stencil = A->stencil; 550 PetscCall(PetscMemcpy(&buffer, A, sizeof(struct _p_Mat))); 551 PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat))); 552 PetscCall(PetscMemcpy(*C, &buffer, sizeof(struct _p_Mat))); 553 ((PetscObject)A)->refct = refct; 554 ((PetscObject)A)->state = state + 1; 555 A->stencil = stencil; 556 557 ((PetscObject)*C)->refct = 1; 558 PetscCall(MatDestroy(C)); 559 PetscFunctionReturn(PETSC_SUCCESS); 560 } 561 562 /*@ 563 MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 564 565 Logically Collective 566 567 Input Parameters: 568 + A - the matrix 569 - flg - bind to the CPU if value of `PETSC_TRUE` 570 571 Level: intermediate 572 573 .seealso: [](ch_matrices), `Mat`, `MatBoundToCPU()` 574 @*/ 575 PetscErrorCode MatBindToCPU(Mat A, PetscBool flg) 576 { 577 PetscFunctionBegin; 578 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 579 PetscValidLogicalCollectiveBool(A, flg, 2); 580 #if defined(PETSC_HAVE_DEVICE) 581 if (A->boundtocpu == flg) PetscFunctionReturn(PETSC_SUCCESS); 582 A->boundtocpu = flg; 583 PetscTryTypeMethod(A, bindtocpu, flg); 584 #endif 585 PetscFunctionReturn(PETSC_SUCCESS); 586 } 587 588 /*@ 589 MatBoundToCPU - query if a matrix is bound to the CPU 590 591 Input Parameter: 592 . A - the matrix 593 594 Output Parameter: 595 . flg - the logical flag 596 597 Level: intermediate 598 599 .seealso: [](ch_matrices), `Mat`, `MatBindToCPU()` 600 @*/ 601 PetscErrorCode MatBoundToCPU(Mat A, PetscBool *flg) 602 { 603 PetscFunctionBegin; 604 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 605 PetscAssertPointer(flg, 2); 606 #if defined(PETSC_HAVE_DEVICE) 607 *flg = A->boundtocpu; 608 #else 609 *flg = PETSC_TRUE; 610 #endif 611 PetscFunctionReturn(PETSC_SUCCESS); 612 } 613 614 PetscErrorCode MatSetValuesCOO_Basic(Mat A, const PetscScalar coo_v[], InsertMode imode) 615 { 616 IS is_coo_i, is_coo_j; 617 const PetscInt *coo_i, *coo_j; 618 PetscInt n, n_i, n_j; 619 PetscScalar zero = 0.; 620 621 PetscFunctionBegin; 622 PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_i", (PetscObject *)&is_coo_i)); 623 PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_j", (PetscObject *)&is_coo_j)); 624 PetscCheck(is_coo_i, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_i IS"); 625 PetscCheck(is_coo_j, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_j IS"); 626 PetscCall(ISGetLocalSize(is_coo_i, &n_i)); 627 PetscCall(ISGetLocalSize(is_coo_j, &n_j)); 628 PetscCheck(n_i == n_j, PETSC_COMM_SELF, PETSC_ERR_COR, "Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT, n_i, n_j); 629 PetscCall(ISGetIndices(is_coo_i, &coo_i)); 630 PetscCall(ISGetIndices(is_coo_j, &coo_j)); 631 if (imode != ADD_VALUES) PetscCall(MatZeroEntries(A)); 632 for (n = 0; n < n_i; n++) PetscCall(MatSetValue(A, coo_i[n], coo_j[n], coo_v ? coo_v[n] : zero, ADD_VALUES)); 633 PetscCall(ISRestoreIndices(is_coo_i, &coo_i)); 634 PetscCall(ISRestoreIndices(is_coo_j, &coo_j)); 635 PetscFunctionReturn(PETSC_SUCCESS); 636 } 637 638 PetscErrorCode MatSetPreallocationCOO_Basic(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) 639 { 640 Mat preallocator; 641 IS is_coo_i, is_coo_j; 642 PetscInt ncoo_i; 643 PetscScalar zero = 0.0; 644 645 PetscFunctionBegin; 646 PetscCall(PetscIntCast(ncoo, &ncoo_i)); 647 PetscCall(PetscLayoutSetUp(A->rmap)); 648 PetscCall(PetscLayoutSetUp(A->cmap)); 649 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator)); 650 PetscCall(MatSetType(preallocator, MATPREALLOCATOR)); 651 PetscCall(MatSetSizes(preallocator, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 652 PetscCall(MatSetLayouts(preallocator, A->rmap, A->cmap)); 653 PetscCall(MatSetUp(preallocator)); 654 for (PetscCount n = 0; n < ncoo; n++) PetscCall(MatSetValue(preallocator, coo_i[n], coo_j[n], zero, INSERT_VALUES)); 655 PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY)); 656 PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY)); 657 PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A)); 658 PetscCall(MatDestroy(&preallocator)); 659 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo_i, coo_i, PETSC_COPY_VALUES, &is_coo_i)); 660 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo_i, coo_j, PETSC_COPY_VALUES, &is_coo_j)); 661 PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_i", (PetscObject)is_coo_i)); 662 PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_j", (PetscObject)is_coo_j)); 663 PetscCall(ISDestroy(&is_coo_i)); 664 PetscCall(ISDestroy(&is_coo_j)); 665 PetscFunctionReturn(PETSC_SUCCESS); 666 } 667 668 /*@C 669 MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices 670 671 Collective 672 673 Input Parameters: 674 + A - matrix being preallocated 675 . ncoo - number of entries 676 . coo_i - row indices 677 - coo_j - column indices 678 679 Level: beginner 680 681 Notes: 682 The indices within `coo_i` and `coo_j` may be modified within this function. The caller should not rely on them 683 having any specific value after this function returns. The arrays can be freed or reused immediately 684 after this function returns. 685 686 Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed 687 but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries 688 are allowed and will be properly added or inserted to the matrix, unless the matrix option `MAT_IGNORE_OFF_PROC_ENTRIES` 689 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 691 If you just want to create a sequential AIJ matrix (`MATSEQAIJ`), and your matrix entries in COO format are unique, you can also use 692 `MatCreateSeqAIJFromTriple()`. But that is not recommended for iterative applications. 693 694 .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, 695 `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`, 696 `DMSetMatrixPreallocateSkip()`, `MatCreateSeqAIJFromTriple()` 697 @*/ 698 PetscErrorCode MatSetPreallocationCOO(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) 699 { 700 PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL; 701 702 PetscFunctionBegin; 703 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 704 PetscValidType(A, 1); 705 if (ncoo) PetscAssertPointer(coo_i, 3); 706 if (ncoo) PetscAssertPointer(coo_j, 4); 707 PetscCall(PetscLayoutSetUp(A->rmap)); 708 PetscCall(PetscLayoutSetUp(A->cmap)); 709 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOO_C", &f)); 710 711 PetscCall(PetscLogEventBegin(MAT_PreallCOO, A, 0, 0, 0)); 712 if (f) { 713 PetscCall((*f)(A, ncoo, coo_i, coo_j)); 714 } else { /* allow fallback, very slow */ 715 PetscCall(MatSetPreallocationCOO_Basic(A, ncoo, coo_i, coo_j)); 716 } 717 PetscCall(PetscLogEventEnd(MAT_PreallCOO, A, 0, 0, 0)); 718 A->preallocated = PETSC_TRUE; 719 A->nonzerostate++; 720 PetscFunctionReturn(PETSC_SUCCESS); 721 } 722 723 /*@C 724 MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices 725 726 Collective 727 728 Input Parameters: 729 + A - matrix being preallocated 730 . ncoo - number of entries 731 . coo_i - row indices (local numbering; may be modified) 732 - coo_j - column indices (local numbering; may be modified) 733 734 Level: beginner 735 736 Notes: 737 The local indices are translated using the local to global mapping, thus `MatSetLocalToGlobalMapping()` must have been 738 called prior to this function. For matrices created with `DMCreateMatrix()` the local to global mapping is often already provided. 739 740 The indices `coo_i` and `coo_j` may be modified within this function. They might be translated to corresponding global 741 indices, but the caller should not rely on them having any specific value after this function returns. The arrays 742 can be freed or reused immediately after this function returns. 743 744 Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed 745 but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries 746 are allowed and will be properly added or inserted to the matrix. 747 748 .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, 749 `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`, 750 `DMSetMatrixPreallocateSkip()` 751 @*/ 752 PetscErrorCode MatSetPreallocationCOOLocal(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) 753 { 754 PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL; 755 756 PetscFunctionBegin; 757 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 758 PetscValidType(A, 1); 759 if (ncoo) PetscAssertPointer(coo_i, 3); 760 if (ncoo) PetscAssertPointer(coo_j, 4); 761 PetscCall(PetscLayoutSetUp(A->rmap)); 762 PetscCall(PetscLayoutSetUp(A->cmap)); 763 764 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", &f)); 765 if (f) { 766 PetscCall((*f)(A, ncoo, coo_i, coo_j)); 767 A->nonzerostate++; 768 } else { 769 PetscInt ncoo_i; 770 ISLocalToGlobalMapping ltog_row, ltog_col; 771 772 PetscCall(MatGetLocalToGlobalMapping(A, <og_row, <og_col)); 773 if (ltog_row) { 774 PetscCall(PetscIntCast(ncoo, &ncoo_i)); 775 PetscCall(ISLocalToGlobalMappingApply(ltog_row, ncoo_i, coo_i, coo_i)); 776 } 777 if (ltog_col) { 778 PetscCall(PetscIntCast(ncoo, &ncoo_i)); 779 PetscCall(ISLocalToGlobalMappingApply(ltog_col, ncoo_i, coo_j, coo_j)); 780 } 781 PetscCall(MatSetPreallocationCOO(A, ncoo, coo_i, coo_j)); 782 } 783 A->preallocated = PETSC_TRUE; 784 PetscFunctionReturn(PETSC_SUCCESS); 785 } 786 787 /*@ 788 MatSetValuesCOO - set values at once in a matrix preallocated using `MatSetPreallocationCOO()` 789 790 Collective 791 792 Input Parameters: 793 + A - matrix being preallocated 794 . coo_v - the matrix values (can be `NULL`) 795 - imode - the insert mode 796 797 Level: beginner 798 799 Notes: 800 The values must follow the order of the indices prescribed with `MatSetPreallocationCOO()` or `MatSetPreallocationCOOLocal()`. 801 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`). 804 805 `MatAssemblyBegin()` and `MatAssemblyEnd()` do not need to be called after this routine. It automatically handles the assembly process. 806 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; 813 814 PetscFunctionBegin; 815 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 816 PetscValidType(A, 1); 817 MatCheckPreallocated(A, 1); 818 PetscValidLogicalCollectiveEnum(A, imode, 3); 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 } 834 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 837 838 Input Parameters: 839 + A - the matrix 840 - flg - flag indicating whether the boundtocpu flag should be propagated 841 842 Level: developer 843 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. 848 849 The bindingpropagates flag itself is also propagated by the above routines. 850 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. 854 855 .seealso: [](ch_matrices), `Mat`, `VecSetBindingPropagates()`, `MatGetBindingPropagates()` 856 @*/ 857 PetscErrorCode MatSetBindingPropagates(Mat A, PetscBool flg) 858 { 859 PetscFunctionBegin; 860 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 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 } 866 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 869 870 Input Parameter: 871 . A - the matrix 872 873 Output Parameter: 874 . flg - flag indicating whether the boundtocpu flag will be propagated 875 876 Level: developer 877 878 .seealso: [](ch_matrices), `Mat`, `MatSetBindingPropagates()` 879 @*/ 880 PetscErrorCode MatGetBindingPropagates(Mat A, PetscBool *flg) 881 { 882 PetscFunctionBegin; 883 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 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 } 892