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