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