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) PetscCall((*B->ops->setfromoptions)(PetscOptionsObject,B)); 231 232 flg = PETSC_FALSE; 233 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)); 234 if (set) PetscCall(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg)); 235 flg = PETSC_FALSE; 236 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)); 237 if (set) PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg)); 238 flg = PETSC_FALSE; 239 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)); 240 if (set) PetscCall(MatSetOption(B,MAT_IGNORE_ZERO_ENTRIES,flg)); 241 242 flg = PETSC_FALSE; 243 PetscCall(PetscOptionsBool("-mat_form_explicit_transpose","Hint to form an explicit transpose for operations like MatMultTranspose","MatSetOption",flg,&flg,&set)); 244 if (set) PetscCall(MatSetOption(B,MAT_FORM_EXPLICIT_TRANSPOSE,flg)); 245 246 /* Bind to CPU if below a user-specified size threshold. 247 * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types, 248 * and putting it here makes is more maintainable than duplicating this for all. */ 249 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)); 250 if (flg && B->rmap->n < bind_below) { 251 PetscCall(MatBindToCPU(B,PETSC_TRUE)); 252 } 253 254 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 255 PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B)); 256 PetscOptionsEnd(); 257 PetscFunctionReturn(0); 258 } 259 260 /*@C 261 MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 262 263 Collective on Mat 264 265 Input Parameters: 266 + A - matrix being preallocated 267 . bs - block size 268 . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 269 . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 270 . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 271 - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 272 273 Level: beginner 274 275 .seealso: `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, 276 `PetscSplitOwnership()` 277 @*/ 278 PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[]) 279 { 280 PetscInt cbs; 281 void (*aij)(void); 282 void (*is)(void); 283 void (*hyp)(void) = NULL; 284 285 PetscFunctionBegin; 286 if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 287 PetscCall(MatSetBlockSize(A,bs)); 288 } 289 PetscCall(PetscLayoutSetUp(A->rmap)); 290 PetscCall(PetscLayoutSetUp(A->cmap)); 291 PetscCall(MatGetBlockSizes(A,&bs,&cbs)); 292 /* these routines assumes bs == cbs, this should be checked somehow */ 293 PetscCall(MatSeqBAIJSetPreallocation(A,bs,0,dnnz)); 294 PetscCall(MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz)); 295 PetscCall(MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu)); 296 PetscCall(MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu)); 297 /* 298 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 299 good before going on with it. 300 */ 301 PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij)); 302 PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is)); 303 #if defined(PETSC_HAVE_HYPRE) 304 PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp)); 305 #endif 306 if (!aij && !is && !hyp) { 307 PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij)); 308 } 309 if (aij || is || hyp) { 310 if (bs == cbs && bs == 1) { 311 PetscCall(MatSeqAIJSetPreallocation(A,0,dnnz)); 312 PetscCall(MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz)); 313 PetscCall(MatISSetPreallocation(A,0,dnnz,0,onnz)); 314 #if defined(PETSC_HAVE_HYPRE) 315 PetscCall(MatHYPRESetPreallocation(A,0,dnnz,0,onnz)); 316 #endif 317 } else { /* Convert block-row precallocation to scalar-row */ 318 PetscInt i,m,*sdnnz,*sonnz; 319 PetscCall(MatGetLocalSize(A,&m,NULL)); 320 PetscCall(PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz)); 321 for (i=0; i<m; i++) { 322 if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs; 323 if (onnz) sonnz[i] = onnz[i/bs] * cbs; 324 } 325 PetscCall(MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL)); 326 PetscCall(MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL)); 327 PetscCall(MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL)); 328 #if defined(PETSC_HAVE_HYPRE) 329 PetscCall(MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL)); 330 #endif 331 PetscCall(PetscFree2(sdnnz,sonnz)); 332 } 333 } 334 PetscFunctionReturn(0); 335 } 336 337 /* 338 Merges some information from Cs header to A; the C object is then destroyed 339 340 This is somewhat different from MatHeaderReplace() it would be nice to merge the code 341 */ 342 PetscErrorCode MatHeaderMerge(Mat A,Mat *C) 343 { 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 /* free all the interior data structures from mat */ 373 PetscCall((*A->ops->destroy)(A)); 374 375 PetscCall(PetscFree(A->defaultvectype)); 376 PetscCall(PetscLayoutDestroy(&A->rmap)); 377 PetscCall(PetscLayoutDestroy(&A->cmap)); 378 PetscCall(PetscFunctionListDestroy(&((PetscObject)A)->qlist)); 379 PetscCall(PetscObjectListDestroy(&((PetscObject)A)->olist)); 380 PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A)); 381 382 /* copy C over to A */ 383 PetscCall(PetscFree(A->factorprefix)); 384 PetscCall(PetscMemcpy(A,*C,sizeof(struct _p_Mat))); 385 386 /* return the parts of A we saved */ 387 ((PetscObject)A)->bops[0] = Abops; 388 A->ops[0] = Aops; 389 ((PetscObject)A)->refct = refct; 390 ((PetscObject)A)->type_name = mtype; 391 ((PetscObject)A)->name = mname; 392 ((PetscObject)A)->prefix = mprefix; 393 ((PetscObject)A)->state = state + 1; 394 A->product = product; 395 A->redundant = redundant; 396 397 /* since these two are copied into A we do not want them destroyed in C */ 398 ((PetscObject)*C)->qlist = NULL; 399 ((PetscObject)*C)->olist = NULL; 400 401 PetscCall(PetscHeaderDestroy(C)); 402 PetscFunctionReturn(0); 403 } 404 /* 405 Replace A's header with that of C; the C object is then destroyed 406 407 This is essentially code moved from MatDestroy() 408 409 This is somewhat different from MatHeaderMerge() it would be nice to merge the code 410 411 Used in DM hence is declared PETSC_EXTERN 412 */ 413 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C) 414 { 415 PetscInt refct; 416 PetscObjectState state; 417 struct _p_Mat buffer; 418 MatStencilInfo stencil; 419 420 PetscFunctionBegin; 421 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 422 PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 423 if (A == *C) PetscFunctionReturn(0); 424 PetscCheckSameComm(A,1,*C,2); 425 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); 426 427 /* swap C and A */ 428 refct = ((PetscObject)A)->refct; 429 state = ((PetscObject)A)->state; 430 stencil = A->stencil; 431 PetscCall(PetscMemcpy(&buffer,A,sizeof(struct _p_Mat))); 432 PetscCall(PetscMemcpy(A,*C,sizeof(struct _p_Mat))); 433 PetscCall(PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat))); 434 ((PetscObject)A)->refct = refct; 435 ((PetscObject)A)->state = state + 1; 436 A->stencil = stencil; 437 438 ((PetscObject)*C)->refct = 1; 439 PetscCall(MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL)); 440 PetscCall(MatDestroy(C)); 441 PetscFunctionReturn(0); 442 } 443 444 /*@ 445 MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 446 447 Logically collective on Mat 448 449 Input Parameters: 450 + A - the matrix 451 - flg - bind to the CPU if value of PETSC_TRUE 452 453 Level: intermediate 454 455 .seealso: `MatBoundToCPU()` 456 @*/ 457 PetscErrorCode MatBindToCPU(Mat A,PetscBool flg) 458 { 459 PetscFunctionBegin; 460 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 461 PetscValidLogicalCollectiveBool(A,flg,2); 462 #if defined(PETSC_HAVE_DEVICE) 463 if (A->boundtocpu == flg) PetscFunctionReturn(0); 464 A->boundtocpu = flg; 465 if (A->ops->bindtocpu) PetscCall((*A->ops->bindtocpu)(A,flg)); 466 #endif 467 PetscFunctionReturn(0); 468 } 469 470 /*@ 471 MatBoundToCPU - query if a matrix is bound to the CPU 472 473 Input Parameter: 474 . A - the matrix 475 476 Output Parameter: 477 . flg - the logical flag 478 479 Level: intermediate 480 481 .seealso: `MatBindToCPU()` 482 @*/ 483 PetscErrorCode MatBoundToCPU(Mat A,PetscBool *flg) 484 { 485 PetscFunctionBegin; 486 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 487 PetscValidBoolPointer(flg,2); 488 #if defined(PETSC_HAVE_DEVICE) 489 *flg = A->boundtocpu; 490 #else 491 *flg = PETSC_TRUE; 492 #endif 493 PetscFunctionReturn(0); 494 } 495 496 PetscErrorCode MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode) 497 { 498 IS is_coo_i,is_coo_j; 499 const PetscInt *coo_i,*coo_j; 500 PetscInt n,n_i,n_j; 501 PetscScalar zero = 0.; 502 503 PetscFunctionBegin; 504 PetscCall(PetscObjectQuery((PetscObject)A,"__PETSc_coo_i",(PetscObject*)&is_coo_i)); 505 PetscCall(PetscObjectQuery((PetscObject)A,"__PETSc_coo_j",(PetscObject*)&is_coo_j)); 506 PetscCheck(is_coo_i,PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_i IS"); 507 PetscCheck(is_coo_j,PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_j IS"); 508 PetscCall(ISGetLocalSize(is_coo_i,&n_i)); 509 PetscCall(ISGetLocalSize(is_coo_j,&n_j)); 510 PetscCheck(n_i == n_j,PETSC_COMM_SELF,PETSC_ERR_COR,"Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT,n_i,n_j); 511 PetscCall(ISGetIndices(is_coo_i,&coo_i)); 512 PetscCall(ISGetIndices(is_coo_j,&coo_j)); 513 if (imode != ADD_VALUES) { 514 PetscCall(MatZeroEntries(A)); 515 } 516 for (n = 0; n < n_i; n++) { 517 PetscCall(MatSetValue(A,coo_i[n],coo_j[n],coo_v ? coo_v[n] : zero,ADD_VALUES)); 518 } 519 PetscCall(ISRestoreIndices(is_coo_i,&coo_i)); 520 PetscCall(ISRestoreIndices(is_coo_j,&coo_j)); 521 PetscFunctionReturn(0); 522 } 523 524 PetscErrorCode MatSetPreallocationCOO_Basic(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 525 { 526 Mat preallocator; 527 IS is_coo_i,is_coo_j; 528 PetscScalar zero = 0.0; 529 530 PetscFunctionBegin; 531 PetscCall(PetscLayoutSetUp(A->rmap)); 532 PetscCall(PetscLayoutSetUp(A->cmap)); 533 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&preallocator)); 534 PetscCall(MatSetType(preallocator,MATPREALLOCATOR)); 535 PetscCall(MatSetSizes(preallocator,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 536 PetscCall(MatSetLayouts(preallocator,A->rmap,A->cmap)); 537 PetscCall(MatSetUp(preallocator)); 538 for (PetscCount n = 0; n < ncoo; n++) { 539 PetscCall(MatSetValue(preallocator,coo_i[n],coo_j[n],zero,INSERT_VALUES)); 540 } 541 PetscCall(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY)); 542 PetscCall(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY)); 543 PetscCall(MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A)); 544 PetscCall(MatDestroy(&preallocator)); 545 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); 546 PetscCall(ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_i,PETSC_COPY_VALUES,&is_coo_i)); 547 PetscCall(ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_j,PETSC_COPY_VALUES,&is_coo_j)); 548 PetscCall(PetscObjectCompose((PetscObject)A,"__PETSc_coo_i",(PetscObject)is_coo_i)); 549 PetscCall(PetscObjectCompose((PetscObject)A,"__PETSc_coo_j",(PetscObject)is_coo_j)); 550 PetscCall(ISDestroy(&is_coo_i)); 551 PetscCall(ISDestroy(&is_coo_j)); 552 PetscFunctionReturn(0); 553 } 554 555 /*@C 556 MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices 557 558 Collective on Mat 559 560 Input Parameters: 561 + A - matrix being preallocated 562 . ncoo - number of entries 563 . coo_i - row indices 564 - coo_j - column indices 565 566 Level: beginner 567 568 Notes: 569 Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed 570 but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries 571 are allowed and will be properly added or inserted to the matrix, unless the matrix option MAT_IGNORE_OFF_PROC_ENTRIES 572 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. 573 574 The arrays coo_i and coo_j may be freed immediately after calling this function. 575 576 .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`, `DMSetMatrixPreallocateSkip()` 577 @*/ 578 PetscErrorCode MatSetPreallocationCOO(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 579 { 580 PetscErrorCode (*f)(Mat,PetscCount,const PetscInt[],const PetscInt[]) = NULL; 581 582 PetscFunctionBegin; 583 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 584 PetscValidType(A,1); 585 if (ncoo) PetscValidIntPointer(coo_i,3); 586 if (ncoo) PetscValidIntPointer(coo_j,4); 587 PetscCall(PetscLayoutSetUp(A->rmap)); 588 PetscCall(PetscLayoutSetUp(A->cmap)); 589 PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOO_C",&f)); 590 591 PetscCall(PetscLogEventBegin(MAT_PreallCOO,A,0,0,0)); 592 if (f) { 593 PetscCall((*f)(A,ncoo,coo_i,coo_j)); 594 } else { /* allow fallback, very slow */ 595 PetscCall(MatSetPreallocationCOO_Basic(A,ncoo,coo_i,coo_j)); 596 } 597 PetscCall(PetscLogEventEnd(MAT_PreallCOO,A,0,0,0)); 598 A->preallocated = PETSC_TRUE; 599 A->nonzerostate++; 600 PetscFunctionReturn(0); 601 } 602 603 /*@C 604 MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices 605 606 Collective on Mat 607 608 Input Parameters: 609 + A - matrix being preallocated 610 . ncoo - number of entries 611 . coo_i - row indices (local numbering; may be modified) 612 - coo_j - column indices (local numbering; may be modified) 613 614 Level: beginner 615 616 Notes: 617 The local indices are translated using the local to global mapping, thus MatSetLocalToGlobalMapping() must have been 618 called prior to this function. 619 620 The indices coo_i and coo_j may be modified within this function. They might be translated to corresponding global 621 indices, but the caller should not rely on them having any specific value after this function returns. The arrays 622 can be freed or reused immediately after this function returns. 623 624 Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed 625 but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries 626 are allowed and will be properly added or inserted to the matrix. 627 628 .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`, `DMSetMatrixPreallocateSkip()` 629 @*/ 630 PetscErrorCode MatSetPreallocationCOOLocal(Mat A,PetscCount ncoo,PetscInt coo_i[],PetscInt coo_j[]) 631 { 632 PetscErrorCode (*f)(Mat,PetscCount,PetscInt[],PetscInt[]) = NULL; 633 634 PetscFunctionBegin; 635 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 636 PetscValidType(A,1); 637 if (ncoo) PetscValidIntPointer(coo_i,3); 638 if (ncoo) PetscValidIntPointer(coo_j,4); 639 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); 640 PetscCall(PetscLayoutSetUp(A->rmap)); 641 PetscCall(PetscLayoutSetUp(A->cmap)); 642 643 PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOOLocal_C",&f)); 644 if (f) { 645 PetscCall((*f)(A,ncoo,coo_i,coo_j)); 646 A->nonzerostate++; 647 } else { 648 ISLocalToGlobalMapping ltog_row,ltog_col; 649 PetscCall(MatGetLocalToGlobalMapping(A,<og_row,<og_col)); 650 if (ltog_row) PetscCall(ISLocalToGlobalMappingApply(ltog_row,ncoo,coo_i,coo_i)); 651 if (ltog_col) PetscCall(ISLocalToGlobalMappingApply(ltog_col,ncoo,coo_j,coo_j)); 652 PetscCall(MatSetPreallocationCOO(A,ncoo,coo_i,coo_j)); 653 } 654 A->preallocated = PETSC_TRUE; 655 PetscFunctionReturn(0); 656 } 657 658 /*@ 659 MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO() 660 661 Collective on Mat 662 663 Input Parameters: 664 + A - matrix being preallocated 665 . coo_v - the matrix values (can be NULL) 666 - imode - the insert mode 667 668 Level: beginner 669 670 Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO() or MatSetPreallocationCOOLocal(). 671 When repeated entries are specified in the COO indices the coo_v values are first properly summed, regardless of the value of imode. 672 The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES). 673 MatAssemblyBegin() and MatAssemblyEnd() do not need to be called after this routine. It automatically handles the assembly process. 674 675 .seealso: `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES` 676 @*/ 677 PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode) 678 { 679 PetscErrorCode (*f)(Mat,const PetscScalar[],InsertMode) = NULL; 680 681 PetscFunctionBegin; 682 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 683 PetscValidType(A,1); 684 MatCheckPreallocated(A,1); 685 PetscValidLogicalCollectiveEnum(A,imode,3); 686 PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSetValuesCOO_C",&f)); 687 PetscCall(PetscLogEventBegin(MAT_SetVCOO,A,0,0,0)); 688 if (f) { 689 PetscCall((*f)(A,coo_v,imode)); 690 } else { /* allow fallback */ 691 PetscCall(MatSetValuesCOO_Basic(A,coo_v,imode)); 692 } 693 PetscCall(PetscLogEventEnd(MAT_SetVCOO,A,0,0,0)); 694 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 695 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 696 PetscFunctionReturn(0); 697 } 698 699 /*@ 700 MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects 701 702 Input Parameters: 703 + A - the matrix 704 - flg - flag indicating whether the boundtocpu flag should be propagated 705 706 Level: developer 707 708 Notes: 709 If the value of flg is set to true, the following will occur: 710 711 MatCreateSubMatrices() and MatCreateRedundantMatrix() will bind created matrices to CPU if the input matrix is bound to the CPU. 712 MatCreateVecs() will bind created vectors to CPU if the input matrix is bound to the CPU. 713 The bindingpropagates flag itself is also propagated by the above routines. 714 715 Developer Notes: 716 If the fine-scale DMDA has the -dm_bind_below option set to true, then DMCreateInterpolationScale() calls MatSetBindingPropagates() 717 on the restriction/interpolation operator to set the bindingpropagates flag to true. 718 719 .seealso: `VecSetBindingPropagates()`, `MatGetBindingPropagates()` 720 @*/ 721 PetscErrorCode MatSetBindingPropagates(Mat A,PetscBool flg) 722 { 723 PetscFunctionBegin; 724 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 725 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 726 A->bindingpropagates = flg; 727 #endif 728 PetscFunctionReturn(0); 729 } 730 731 /*@ 732 MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects 733 734 Input Parameter: 735 . A - the matrix 736 737 Output Parameter: 738 . flg - flag indicating whether the boundtocpu flag will be propagated 739 740 Level: developer 741 742 .seealso: `MatSetBindingPropagates()` 743 @*/ 744 PetscErrorCode MatGetBindingPropagates(Mat A,PetscBool *flg) 745 { 746 PetscFunctionBegin; 747 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 748 PetscValidBoolPointer(flg,2); 749 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 750 *flg = A->bindingpropagates; 751 #else 752 *flg = PETSC_FALSE; 753 #endif 754 PetscFunctionReturn(0); 755 } 756