1 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 2 3 PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat, PetscInt rbs, PetscInt cbs) { 4 PetscFunctionBegin; 5 if (!mat->preallocated) PetscFunctionReturn(0); 6 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); 7 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); 8 PetscFunctionReturn(0); 9 } 10 11 PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y, PetscScalar a) { 12 PetscInt i, start, end; 13 PetscScalar alpha = a; 14 PetscBool prevoption; 15 16 PetscFunctionBegin; 17 PetscCall(MatGetOption(Y, MAT_NO_OFF_PROC_ENTRIES, &prevoption)); 18 PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 19 PetscCall(MatGetOwnershipRange(Y, &start, &end)); 20 for (i = start; i < end; i++) { 21 if (i < Y->cmap->N) PetscCall(MatSetValues(Y, 1, &i, 1, &i, &alpha, ADD_VALUES)); 22 } 23 PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY)); 24 PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY)); 25 PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, prevoption)); 26 PetscFunctionReturn(0); 27 } 28 29 /*@ 30 MatCreate - Creates a matrix where the type is determined 31 from either a call to `MatSetType()` or from the options database 32 with a call to `MatSetFromOptions()`. The default matrix type is 33 `MATAIJ`, using the routines `MatCreateSeqAIJ()` or `MatCreateAIJ()` 34 if you do not set a type in the options database. If you never 35 call `MatSetType()` or `MatSetFromOptions()` it will generate an 36 error when you try to use the matrix. 37 38 Collective 39 40 Input Parameter: 41 . comm - MPI communicator 42 43 Output Parameter: 44 . A - the matrix 45 46 Options Database Keys: 47 + -mat_type seqaij - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()` 48 . -mat_type mpiaij - `MATMPIAIJ` type, uses `MatCreateAIJ()` 49 . -mat_type seqdense - `MATSEQDENSE`, uses `MatCreateSeqDense()` 50 . -mat_type mpidense - `MATMPIDENSE` type, uses `MatCreateDense()` 51 . -mat_type seqbaij - `MATSEQBAIJ` type, uses `MatCreateSeqBAIJ()` 52 - -mat_type mpibaij - `MATMPIBAIJ` type, uses `MatCreateBAIJ()` 53 54 Even More Options Database Keys: 55 See the manpages for particular formats (e.g., `MatCreateSeqAIJ()`) 56 for additional format-specific options. 57 58 Level: beginner 59 60 `.seealso: `MatCreateSeqAIJ()`, `MatCreateAIJ()`, 61 `MatCreateSeqDense()`, `MatCreateDense()`, 62 `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, 63 `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`, 64 `MatConvert()` 65 @*/ 66 PetscErrorCode MatCreate(MPI_Comm comm, Mat *A) { 67 Mat B; 68 69 PetscFunctionBegin; 70 PetscValidPointer(A, 2); 71 72 *A = NULL; 73 PetscCall(MatInitializePackage()); 74 75 PetscCall(PetscHeaderCreate(B, MAT_CLASSID, "Mat", "Matrix", "Mat", comm, MatDestroy, MatView)); 76 PetscCall(PetscLayoutCreate(comm, &B->rmap)); 77 PetscCall(PetscLayoutCreate(comm, &B->cmap)); 78 PetscCall(PetscStrallocpy(VECSTANDARD, &B->defaultvectype)); 79 PetscCall(PetscStrallocpy(PETSCRANDER48, &B->defaultrandtype)); 80 81 B->symmetric = PETSC_BOOL3_UNKNOWN; 82 B->hermitian = PETSC_BOOL3_UNKNOWN; 83 B->structurally_symmetric = PETSC_BOOL3_UNKNOWN; 84 B->spd = PETSC_BOOL3_UNKNOWN; 85 B->symmetry_eternal = PETSC_FALSE; 86 B->structural_symmetry_eternal = PETSC_FALSE; 87 88 B->congruentlayouts = PETSC_DECIDE; 89 B->preallocated = PETSC_FALSE; 90 #if defined(PETSC_HAVE_DEVICE) 91 B->boundtocpu = PETSC_TRUE; 92 #endif 93 *A = B; 94 PetscFunctionReturn(0); 95 } 96 97 /*@ 98 MatSetErrorIfFailure - Causes `Mat` to generate an immediate error, for example a zero pivot, is detected. 99 100 Logically Collective on mat 101 102 Input Parameters: 103 + mat - matrix obtained from `MatCreate()` 104 - flg - `PETSC_TRUE` indicates you want the error generated 105 106 Level: advanced 107 108 Note: 109 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 110 or result in a `KSPConvergedReason` indicating the method did not converge. 111 112 .seealso: `PCSetErrorIfFailure()`, `KSPConvergedReason`, `SNESConvergedReason` 113 @*/ 114 PetscErrorCode MatSetErrorIfFailure(Mat mat, PetscBool flg) { 115 PetscFunctionBegin; 116 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 117 PetscValidLogicalCollectiveBool(mat, flg, 2); 118 mat->erroriffailure = flg; 119 PetscFunctionReturn(0); 120 } 121 122 /*@ 123 MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 124 125 Collective on Mat 126 127 Input Parameters: 128 + A - the matrix 129 . m - number of local rows (or `PETSC_DECIDE`) 130 . n - number of local columns (or `PETSC_DECIDE`) 131 . M - number of global rows (or `PETSC_DETERMINE`) 132 - N - number of global columns (or `PETSC_DETERMINE`) 133 134 Notes: 135 m (n) and M (N) cannot be both `PETSC_DECIDE` 136 If one processor calls this with M (N) of `PETSC_DECIDE` then all processors must, otherwise the program will hang. 137 138 If `PETSC_DECIDE` is not used for the arguments 'm' and 'n', then the 139 user must ensure that they are chosen to be compatible with the 140 vectors. To do this, one first considers the matrix-vector product 141 'y = A x'. The 'm' that is used in the above routine must match the 142 local size used in the vector creation routine VecCreateMPI() for 'y'. 143 Likewise, the 'n' used must match that used as the local size in 144 `VecCreateMPI()` for 'x'. 145 146 You cannot change the sizes once they have been set. 147 148 The sizes must be set before `MatSetUp()` or MatXXXSetPreallocation() is called. 149 150 Level: beginner 151 152 .seealso: `MatGetSize()`, `PetscSplitOwnership()` 153 @*/ 154 PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) { 155 PetscFunctionBegin; 156 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 157 PetscValidLogicalCollectiveInt(A, M, 4); 158 PetscValidLogicalCollectiveInt(A, N, 5); 159 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); 160 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); 161 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, 162 A->rmap->n, A->rmap->N); 163 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, 164 A->cmap->n, A->cmap->N); 165 A->rmap->n = m; 166 A->cmap->n = n; 167 A->rmap->N = M > -1 ? M : A->rmap->N; 168 A->cmap->N = N > -1 ? N : A->cmap->N; 169 PetscFunctionReturn(0); 170 } 171 172 /*@ 173 MatSetFromOptions - Creates a matrix where the type is determined 174 from the options database. Generates a parallel MPI matrix if the 175 communicator has more than one processor. The default matrix type is 176 `MATAIJ`, using the routines `MatCreateSeqAIJ()` and `MatCreateAIJ()` if 177 you do not select a type in the options database. 178 179 Collective on B 180 181 Input Parameter: 182 . A - the matrix 183 184 Options Database Keys: 185 + -mat_type seqaij - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()` 186 . -mat_type mpiaij - `MATMPIAIJ` type, uses `MatCreateAIJ()` 187 . -mat_type seqdense - `MATSEQDENSE` type, uses `MatCreateSeqDense()` 188 . -mat_type mpidense - `MATMPIDENSE`, uses `MatCreateDense()` 189 . -mat_type seqbaij - `MATSEQBAIJ`, uses `MatCreateSeqBAIJ()` 190 - -mat_type mpibaij - `MATMPIBAIJ`, uses `MatCreateBAIJ()` 191 192 Even More Options Database Keys: 193 See the manpages for particular formats (e.g., `MatCreateSeqAIJ()`) 194 for additional format-specific options. 195 196 Level: beginner 197 198 .seealso: `MatCreateSeqAIJ(()`, `MatCreateAIJ()`, 199 `MatCreateSeqDense()`, `MatCreateDense()`, 200 `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, 201 `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`, 202 `MatConvert()` 203 @*/ 204 PetscErrorCode MatSetFromOptions(Mat B) { 205 const char *deft = MATAIJ; 206 char type[256]; 207 PetscBool flg, set; 208 PetscInt bind_below = 0; 209 210 PetscFunctionBegin; 211 PetscValidHeaderSpecific(B, MAT_CLASSID, 1); 212 213 PetscObjectOptionsBegin((PetscObject)B); 214 215 if (B->rmap->bs < 0) { 216 PetscInt newbs = -1; 217 PetscCall(PetscOptionsInt("-mat_block_size", "Set the blocksize used to store the matrix", "MatSetBlockSize", newbs, &newbs, &flg)); 218 if (flg) { 219 PetscCall(PetscLayoutSetBlockSize(B->rmap, newbs)); 220 PetscCall(PetscLayoutSetBlockSize(B->cmap, newbs)); 221 } 222 } 223 224 PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, deft, type, 256, &flg)); 225 if (flg) { 226 PetscCall(MatSetType(B, type)); 227 } else if (!((PetscObject)B)->type_name) { 228 PetscCall(MatSetType(B, deft)); 229 } 230 231 PetscCall(PetscOptionsName("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", &B->checksymmetryonassembly)); 232 PetscCall(PetscOptionsReal("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", B->checksymmetrytol, &B->checksymmetrytol, NULL)); 233 PetscCall(PetscOptionsBool("-mat_null_space_test", "Checks if provided null space is correct in MatAssemblyEnd()", "MatSetNullSpaceTest", B->checknullspaceonassembly, &B->checknullspaceonassembly, NULL)); 234 PetscCall(PetscOptionsBool("-mat_error_if_failure", "Generate an error if an error occurs when factoring the matrix", "MatSetErrorIfFailure", B->erroriffailure, &B->erroriffailure, NULL)); 235 236 PetscTryTypeMethod(B, setfromoptions, PetscOptionsObject); 237 238 flg = PETSC_FALSE; 239 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)); 240 if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, flg)); 241 flg = PETSC_FALSE; 242 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)); 243 if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, flg)); 244 flg = PETSC_FALSE; 245 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)); 246 if (set) PetscCall(MatSetOption(B, MAT_IGNORE_ZERO_ENTRIES, flg)); 247 248 flg = PETSC_FALSE; 249 PetscCall(PetscOptionsBool("-mat_form_explicit_transpose", "Hint to form an explicit transpose for operations like MatMultTranspose", "MatSetOption", flg, &flg, &set)); 250 if (set) PetscCall(MatSetOption(B, MAT_FORM_EXPLICIT_TRANSPOSE, flg)); 251 252 /* Bind to CPU if below a user-specified size threshold. 253 * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types, 254 * and putting it here makes is more maintainable than duplicating this for all. */ 255 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)); 256 if (flg && B->rmap->n < bind_below) PetscCall(MatBindToCPU(B, PETSC_TRUE)); 257 258 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 259 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)B, PetscOptionsObject)); 260 PetscOptionsEnd(); 261 PetscFunctionReturn(0); 262 } 263 264 /*@C 265 MatXAIJSetPreallocation - set preallocation for serial and parallel `MATAIJ`, `MATBAIJ`, and `MATSBAIJ` matrices and their unassembled versions. 266 267 Collective on A 268 269 Input Parameters: 270 + A - matrix being preallocated 271 . bs - block size 272 . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 273 . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 274 . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 275 - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 276 277 Level: beginner 278 279 .seealso: `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, 280 `PetscSplitOwnership()` 281 @*/ 282 PetscErrorCode MatXAIJSetPreallocation(Mat A, PetscInt bs, const PetscInt dnnz[], const PetscInt onnz[], const PetscInt dnnzu[], const PetscInt onnzu[]) { 283 PetscInt cbs; 284 void (*aij)(void); 285 void (*is)(void); 286 void (*hyp)(void) = NULL; 287 288 PetscFunctionBegin; 289 if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 290 PetscCall(MatSetBlockSize(A, bs)); 291 } 292 PetscCall(PetscLayoutSetUp(A->rmap)); 293 PetscCall(PetscLayoutSetUp(A->cmap)); 294 PetscCall(MatGetBlockSizes(A, &bs, &cbs)); 295 /* these routines assumes bs == cbs, this should be checked somehow */ 296 PetscCall(MatSeqBAIJSetPreallocation(A, bs, 0, dnnz)); 297 PetscCall(MatMPIBAIJSetPreallocation(A, bs, 0, dnnz, 0, onnz)); 298 PetscCall(MatSeqSBAIJSetPreallocation(A, bs, 0, dnnzu)); 299 PetscCall(MatMPISBAIJSetPreallocation(A, bs, 0, dnnzu, 0, onnzu)); 300 /* 301 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 302 good before going on with it. 303 */ 304 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatMPIAIJSetPreallocation_C", &aij)); 305 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatISSetPreallocation_C", &is)); 306 #if defined(PETSC_HAVE_HYPRE) 307 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatHYPRESetPreallocation_C", &hyp)); 308 #endif 309 if (!aij && !is && !hyp) PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", &aij)); 310 if (aij || is || hyp) { 311 if (bs == cbs && bs == 1) { 312 PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz)); 313 PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz, 0, onnz)); 314 PetscCall(MatISSetPreallocation(A, 0, dnnz, 0, onnz)); 315 #if defined(PETSC_HAVE_HYPRE) 316 PetscCall(MatHYPRESetPreallocation(A, 0, dnnz, 0, onnz)); 317 #endif 318 } else { /* Convert block-row precallocation to scalar-row */ 319 PetscInt i, m, *sdnnz, *sonnz; 320 PetscCall(MatGetLocalSize(A, &m, NULL)); 321 PetscCall(PetscMalloc2((!!dnnz) * m, &sdnnz, (!!onnz) * m, &sonnz)); 322 for (i = 0; i < m; i++) { 323 if (dnnz) sdnnz[i] = dnnz[i / bs] * cbs; 324 if (onnz) sonnz[i] = onnz[i / bs] * cbs; 325 } 326 PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL)); 327 PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 328 PetscCall(MatISSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 329 #if defined(PETSC_HAVE_HYPRE) 330 PetscCall(MatHYPRESetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 331 #endif 332 PetscCall(PetscFree2(sdnnz, sonnz)); 333 } 334 } 335 PetscFunctionReturn(0); 336 } 337 338 /* 339 Merges some information from Cs header to A; the C object is then destroyed 340 341 This is somewhat different from MatHeaderReplace() it would be nice to merge the code 342 */ 343 PetscErrorCode MatHeaderMerge(Mat A, Mat *C) { 344 PetscInt refct; 345 PetscOps Abops; 346 struct _MatOps Aops; 347 char *mtype, *mname, *mprefix; 348 Mat_Product *product; 349 Mat_Redundant *redundant; 350 PetscObjectState state; 351 352 PetscFunctionBegin; 353 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 354 PetscValidHeaderSpecific(*C, MAT_CLASSID, 2); 355 if (A == *C) PetscFunctionReturn(0); 356 PetscCheckSameComm(A, 1, *C, 2); 357 /* save the parts of A we need */ 358 Abops = ((PetscObject)A)->bops[0]; 359 Aops = A->ops[0]; 360 refct = ((PetscObject)A)->refct; 361 mtype = ((PetscObject)A)->type_name; 362 mname = ((PetscObject)A)->name; 363 state = ((PetscObject)A)->state; 364 mprefix = ((PetscObject)A)->prefix; 365 product = A->product; 366 redundant = A->redundant; 367 368 /* zero these so the destroy below does not free them */ 369 ((PetscObject)A)->type_name = NULL; 370 ((PetscObject)A)->name = NULL; 371 372 /* 373 free all the interior data structures from mat 374 cannot use PetscUseTypeMethod(A,destroy); because compiler 375 thinks it may print NULL type_name and name 376 */ 377 PetscTryTypeMethod(A, destroy); 378 379 PetscCall(PetscFree(A->defaultvectype)); 380 PetscCall(PetscFree(A->defaultrandtype)); 381 PetscCall(PetscLayoutDestroy(&A->rmap)); 382 PetscCall(PetscLayoutDestroy(&A->cmap)); 383 PetscCall(PetscFunctionListDestroy(&((PetscObject)A)->qlist)); 384 PetscCall(PetscObjectListDestroy(&((PetscObject)A)->olist)); 385 PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A)); 386 387 /* copy C over to A */ 388 PetscCall(PetscFree(A->factorprefix)); 389 PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat))); 390 391 /* return the parts of A we saved */ 392 ((PetscObject)A)->bops[0] = Abops; 393 A->ops[0] = Aops; 394 ((PetscObject)A)->refct = refct; 395 ((PetscObject)A)->type_name = mtype; 396 ((PetscObject)A)->name = mname; 397 ((PetscObject)A)->prefix = mprefix; 398 ((PetscObject)A)->state = state + 1; 399 A->product = product; 400 A->redundant = redundant; 401 402 /* since these two are copied into A we do not want them destroyed in C */ 403 ((PetscObject)*C)->qlist = NULL; 404 ((PetscObject)*C)->olist = NULL; 405 406 PetscCall(PetscHeaderDestroy(C)); 407 PetscFunctionReturn(0); 408 } 409 /* 410 Replace A's header with that of C; the C object is then destroyed 411 412 This is essentially code moved from MatDestroy() 413 414 This is somewhat different from MatHeaderMerge() it would be nice to merge the code 415 416 Used in DM hence is declared PETSC_EXTERN 417 */ 418 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A, Mat *C) { 419 PetscInt refct; 420 PetscObjectState state; 421 struct _p_Mat buffer; 422 MatStencilInfo stencil; 423 424 PetscFunctionBegin; 425 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 426 PetscValidHeaderSpecific(*C, MAT_CLASSID, 2); 427 if (A == *C) PetscFunctionReturn(0); 428 PetscCheckSameComm(A, 1, *C, 2); 429 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); 430 431 /* swap C and A */ 432 refct = ((PetscObject)A)->refct; 433 state = ((PetscObject)A)->state; 434 stencil = A->stencil; 435 PetscCall(PetscMemcpy(&buffer, A, sizeof(struct _p_Mat))); 436 PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat))); 437 PetscCall(PetscMemcpy(*C, &buffer, sizeof(struct _p_Mat))); 438 ((PetscObject)A)->refct = refct; 439 ((PetscObject)A)->state = state + 1; 440 A->stencil = stencil; 441 442 ((PetscObject)*C)->refct = 1; 443 PetscCall(MatShellSetOperation(*C, MATOP_DESTROY, (void (*)(void))NULL)); 444 PetscCall(MatDestroy(C)); 445 PetscFunctionReturn(0); 446 } 447 448 /*@ 449 MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 450 451 Logically collective on A 452 453 Input Parameters: 454 + A - the matrix 455 - flg - bind to the CPU if value of `PETSC_TRUE` 456 457 Level: intermediate 458 459 .seealso: `MatBoundToCPU()` 460 @*/ 461 PetscErrorCode MatBindToCPU(Mat A, PetscBool flg) { 462 PetscFunctionBegin; 463 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 464 PetscValidLogicalCollectiveBool(A, flg, 2); 465 #if defined(PETSC_HAVE_DEVICE) 466 if (A->boundtocpu == flg) PetscFunctionReturn(0); 467 A->boundtocpu = flg; 468 PetscTryTypeMethod(A, bindtocpu, flg); 469 #endif 470 PetscFunctionReturn(0); 471 } 472 473 /*@ 474 MatBoundToCPU - query if a matrix is bound to the CPU 475 476 Input Parameter: 477 . A - the matrix 478 479 Output Parameter: 480 . flg - the logical flag 481 482 Level: intermediate 483 484 .seealso: `MatBindToCPU()` 485 @*/ 486 PetscErrorCode MatBoundToCPU(Mat A, PetscBool *flg) { 487 PetscFunctionBegin; 488 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 489 PetscValidBoolPointer(flg, 2); 490 #if defined(PETSC_HAVE_DEVICE) 491 *flg = A->boundtocpu; 492 #else 493 *flg = PETSC_TRUE; 494 #endif 495 PetscFunctionReturn(0); 496 } 497 498 PetscErrorCode MatSetValuesCOO_Basic(Mat A, const PetscScalar coo_v[], InsertMode imode) { 499 IS is_coo_i, is_coo_j; 500 const PetscInt *coo_i, *coo_j; 501 PetscInt n, n_i, n_j; 502 PetscScalar zero = 0.; 503 504 PetscFunctionBegin; 505 PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_i", (PetscObject *)&is_coo_i)); 506 PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_j", (PetscObject *)&is_coo_j)); 507 PetscCheck(is_coo_i, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_i IS"); 508 PetscCheck(is_coo_j, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_j IS"); 509 PetscCall(ISGetLocalSize(is_coo_i, &n_i)); 510 PetscCall(ISGetLocalSize(is_coo_j, &n_j)); 511 PetscCheck(n_i == n_j, PETSC_COMM_SELF, PETSC_ERR_COR, "Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT, n_i, n_j); 512 PetscCall(ISGetIndices(is_coo_i, &coo_i)); 513 PetscCall(ISGetIndices(is_coo_j, &coo_j)); 514 if (imode != ADD_VALUES) PetscCall(MatZeroEntries(A)); 515 for (n = 0; n < n_i; n++) PetscCall(MatSetValue(A, coo_i[n], coo_j[n], coo_v ? coo_v[n] : zero, ADD_VALUES)); 516 PetscCall(ISRestoreIndices(is_coo_i, &coo_i)); 517 PetscCall(ISRestoreIndices(is_coo_j, &coo_j)); 518 PetscFunctionReturn(0); 519 } 520 521 PetscErrorCode MatSetPreallocationCOO_Basic(Mat A, PetscCount ncoo, const PetscInt coo_i[], const PetscInt coo_j[]) { 522 Mat preallocator; 523 IS is_coo_i, is_coo_j; 524 PetscScalar zero = 0.0; 525 526 PetscFunctionBegin; 527 PetscCall(PetscLayoutSetUp(A->rmap)); 528 PetscCall(PetscLayoutSetUp(A->cmap)); 529 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator)); 530 PetscCall(MatSetType(preallocator, MATPREALLOCATOR)); 531 PetscCall(MatSetSizes(preallocator, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 532 PetscCall(MatSetLayouts(preallocator, A->rmap, A->cmap)); 533 PetscCall(MatSetUp(preallocator)); 534 for (PetscCount n = 0; n < ncoo; n++) PetscCall(MatSetValue(preallocator, coo_i[n], coo_j[n], zero, INSERT_VALUES)); 535 PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY)); 536 PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY)); 537 PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A)); 538 PetscCall(MatDestroy(&preallocator)); 539 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); 540 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo, coo_i, PETSC_COPY_VALUES, &is_coo_i)); 541 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo, coo_j, PETSC_COPY_VALUES, &is_coo_j)); 542 PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_i", (PetscObject)is_coo_i)); 543 PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_j", (PetscObject)is_coo_j)); 544 PetscCall(ISDestroy(&is_coo_i)); 545 PetscCall(ISDestroy(&is_coo_j)); 546 PetscFunctionReturn(0); 547 } 548 549 /*@C 550 MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices 551 552 Collective on A 553 554 Input Parameters: 555 + A - matrix being preallocated 556 . ncoo - number of entries 557 . coo_i - row indices 558 - coo_j - column indices 559 560 Level: beginner 561 562 Notes: 563 The indices coo_i and coo_j may be modified within this function. The caller should not rely on them 564 having any specific value after this function returns. The arrays can be freed or reused immediately 565 after this function returns. 566 567 Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed 568 but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries 569 are allowed and will be properly added or inserted to the matrix, unless the matrix option MAT_IGNORE_OFF_PROC_ENTRIES 570 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. 571 572 .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`, `DMSetMatrixPreallocateSkip()` 573 @*/ 574 PetscErrorCode MatSetPreallocationCOO(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) { 575 PetscErrorCode (*f)(Mat, PetscCount, const PetscInt[], const PetscInt[]) = NULL; 576 577 PetscFunctionBegin; 578 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 579 PetscValidType(A, 1); 580 if (ncoo) PetscValidIntPointer(coo_i, 3); 581 if (ncoo) PetscValidIntPointer(coo_j, 4); 582 PetscCall(PetscLayoutSetUp(A->rmap)); 583 PetscCall(PetscLayoutSetUp(A->cmap)); 584 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOO_C", &f)); 585 586 PetscCall(PetscLogEventBegin(MAT_PreallCOO, A, 0, 0, 0)); 587 if (f) { 588 PetscCall((*f)(A, ncoo, coo_i, coo_j)); 589 } else { /* allow fallback, very slow */ 590 PetscCall(MatSetPreallocationCOO_Basic(A, ncoo, coo_i, coo_j)); 591 } 592 PetscCall(PetscLogEventEnd(MAT_PreallCOO, A, 0, 0, 0)); 593 A->preallocated = PETSC_TRUE; 594 A->nonzerostate++; 595 PetscFunctionReturn(0); 596 } 597 598 /*@C 599 MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices 600 601 Collective on A 602 603 Input Parameters: 604 + A - matrix being preallocated 605 . ncoo - number of entries 606 . coo_i - row indices (local numbering; may be modified) 607 - coo_j - column indices (local numbering; may be modified) 608 609 Level: beginner 610 611 Notes: 612 The local indices are translated using the local to global mapping, thus `MatSetLocalToGlobalMapping()` must have been 613 called prior to this function. For matrices created with `DMCreateMatrix()` the local to global mapping is often already provided. 614 615 The indices coo_i and coo_j may be modified within this function. They might be translated to corresponding global 616 indices, but the caller should not rely on them having any specific value after this function returns. The arrays 617 can be freed or reused immediately after this function returns. 618 619 Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed 620 but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries 621 are allowed and will be properly added or inserted to the matrix. 622 623 .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`, `DMSetMatrixPreallocateSkip()` 624 @*/ 625 PetscErrorCode MatSetPreallocationCOOLocal(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) { 626 PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL; 627 628 PetscFunctionBegin; 629 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 630 PetscValidType(A, 1); 631 if (ncoo) PetscValidIntPointer(coo_i, 3); 632 if (ncoo) PetscValidIntPointer(coo_j, 4); 633 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); 634 PetscCall(PetscLayoutSetUp(A->rmap)); 635 PetscCall(PetscLayoutSetUp(A->cmap)); 636 637 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", &f)); 638 if (f) { 639 PetscCall((*f)(A, ncoo, coo_i, coo_j)); 640 A->nonzerostate++; 641 } else { 642 ISLocalToGlobalMapping ltog_row, ltog_col; 643 PetscCall(MatGetLocalToGlobalMapping(A, <og_row, <og_col)); 644 if (ltog_row) PetscCall(ISLocalToGlobalMappingApply(ltog_row, ncoo, coo_i, coo_i)); 645 if (ltog_col) PetscCall(ISLocalToGlobalMappingApply(ltog_col, ncoo, coo_j, coo_j)); 646 PetscCall(MatSetPreallocationCOO(A, ncoo, coo_i, coo_j)); 647 } 648 A->preallocated = PETSC_TRUE; 649 PetscFunctionReturn(0); 650 } 651 652 /*@ 653 MatSetValuesCOO - set values at once in a matrix preallocated using `MatSetPreallocationCOO()` 654 655 Collective on A 656 657 Input Parameters: 658 + A - matrix being preallocated 659 . coo_v - the matrix values (can be NULL) 660 - imode - the insert mode 661 662 Level: beginner 663 664 Notes: 665 The values must follow the order of the indices prescribed with `MatSetPreallocationCOO()` or `MatSetPreallocationCOOLocal()`. 666 667 When repeated entries are specified in the COO indices the coo_v values are first properly summed, regardless of the value of imode. 668 The imode flag indicates if coo_v must be added to the current values of the matrix (`ADD_VALUES`) or overwritten (`INSERT_VALUES`). 669 670 `MatAssemblyBegin()` and `MatAssemblyEnd()` do not need to be called after this routine. It automatically handles the assembly process. 671 672 .seealso: `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES` 673 @*/ 674 PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode) { 675 PetscErrorCode (*f)(Mat, const PetscScalar[], InsertMode) = NULL; 676 677 PetscFunctionBegin; 678 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 679 PetscValidType(A, 1); 680 MatCheckPreallocated(A, 1); 681 PetscValidLogicalCollectiveEnum(A, imode, 3); 682 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetValuesCOO_C", &f)); 683 PetscCall(PetscLogEventBegin(MAT_SetVCOO, A, 0, 0, 0)); 684 if (f) { 685 PetscCall((*f)(A, coo_v, imode)); 686 } else { /* allow fallback */ 687 PetscCall(MatSetValuesCOO_Basic(A, coo_v, imode)); 688 } 689 PetscCall(PetscLogEventEnd(MAT_SetVCOO, A, 0, 0, 0)); 690 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 691 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 692 PetscFunctionReturn(0); 693 } 694 695 /*@ 696 MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects 697 698 Input Parameters: 699 + A - the matrix 700 - flg - flag indicating whether the boundtocpu flag should be propagated 701 702 Level: developer 703 704 Notes: 705 If the value of flg is set to true, the following will occur: 706 707 `MatCreateSubMatrices()` and `MatCreateRedundantMatrix()` will bind created matrices to CPU if the input matrix is bound to the CPU. 708 709 `MatCreateVecs()` will bind created vectors to CPU if the input matrix is bound to the CPU. 710 The bindingpropagates flag itself is also propagated by the above routines. 711 712 Developer Note: 713 If the fine-scale `DMDA `has the -dm_bind_below option set to true, then `DMCreateInterpolationScale()` calls `MatSetBindingPropagates()` 714 on the restriction/interpolation operator to set the bindingpropagates flag to true. 715 716 .seealso: `VecSetBindingPropagates()`, `MatGetBindingPropagates()` 717 @*/ 718 PetscErrorCode MatSetBindingPropagates(Mat A, PetscBool flg) { 719 PetscFunctionBegin; 720 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 721 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 722 A->bindingpropagates = flg; 723 #endif 724 PetscFunctionReturn(0); 725 } 726 727 /*@ 728 MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects 729 730 Input Parameter: 731 . A - the matrix 732 733 Output Parameter: 734 . flg - flag indicating whether the boundtocpu flag will be propagated 735 736 Level: developer 737 738 .seealso: `MatSetBindingPropagates()` 739 @*/ 740 PetscErrorCode MatGetBindingPropagates(Mat A, PetscBool *flg) { 741 PetscFunctionBegin; 742 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 743 PetscValidBoolPointer(flg, 2); 744 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 745 *flg = A->bindingpropagates; 746 #else 747 *flg = PETSC_FALSE; 748 #endif 749 PetscFunctionReturn(0); 750 } 751