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