1 2 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3 4 PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat,PetscInt rbs, PetscInt cbs) 5 { 6 PetscFunctionBegin; 7 if (!mat->preallocated) PetscFunctionReturn(0); 8 if (mat->rmap->bs > 0 && mat->rmap->bs != rbs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot change row block size %D to %D\n",mat->rmap->bs,rbs); 9 if (mat->cmap->bs > 0 && mat->cmap->bs != cbs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot change column block size %D to %D\n",mat->cmap->bs,cbs); 10 PetscFunctionReturn(0); 11 } 12 13 PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y,PetscScalar a) 14 { 15 PetscErrorCode ierr; 16 PetscInt i,start,end; 17 PetscScalar alpha = a; 18 PetscBool prevoption; 19 20 PetscFunctionBegin; 21 ierr = MatGetOption(Y,MAT_NO_OFF_PROC_ENTRIES,&prevoption);CHKERRQ(ierr); 22 ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 23 ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 24 for (i=start; i<end; i++) { 25 if (i < Y->cmap->N) { 26 ierr = MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);CHKERRQ(ierr); 27 } 28 } 29 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 30 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 31 ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption);CHKERRQ(ierr); 32 PetscFunctionReturn(0); 33 } 34 35 /*@ 36 MatCreate - Creates a matrix where the type is determined 37 from either a call to MatSetType() or from the options database 38 with a call to MatSetFromOptions(). The default matrix type is 39 AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ() 40 if you do not set a type in the options database. If you never 41 call MatSetType() or MatSetFromOptions() it will generate an 42 error when you try to use the matrix. 43 44 Collective 45 46 Input Parameter: 47 . comm - MPI communicator 48 49 Output Parameter: 50 . A - the matrix 51 52 Options Database Keys: 53 + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 54 . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 55 . -mat_type seqdense - dense type, uses MatCreateSeqDense() 56 . -mat_type mpidense - dense type, uses MatCreateDense() 57 . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 58 - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 59 60 Even More Options Database Keys: 61 See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 62 for additional format-specific options. 63 64 Level: beginner 65 66 .seealso: MatCreateSeqAIJ(), MatCreateAIJ(), 67 MatCreateSeqDense(), MatCreateDense(), 68 MatCreateSeqBAIJ(), MatCreateBAIJ(), 69 MatCreateSeqSBAIJ(), MatCreateSBAIJ(), 70 MatConvert() 71 @*/ 72 PetscErrorCode MatCreate(MPI_Comm comm,Mat *A) 73 { 74 Mat B; 75 PetscErrorCode ierr; 76 77 PetscFunctionBegin; 78 PetscValidPointer(A,2); 79 80 *A = NULL; 81 ierr = MatInitializePackage();CHKERRQ(ierr); 82 83 ierr = PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);CHKERRQ(ierr); 84 ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr); 85 ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr); 86 ierr = PetscStrallocpy(VECSTANDARD,&B->defaultvectype);CHKERRQ(ierr); 87 88 B->congruentlayouts = PETSC_DECIDE; 89 B->preallocated = PETSC_FALSE; 90 #if defined(PETSC_HAVE_DEVICE) 91 B->boundtocpu = PETSC_TRUE; 92 #endif 93 *A = B; 94 PetscFunctionReturn(0); 95 } 96 97 /*@ 98 MatSetErrorIfFailure - Causes Mat to generate an error, for example a zero pivot, is detected. 99 100 Logically Collective on Mat 101 102 Input Parameters: 103 + mat - matrix obtained from MatCreate() 104 - flg - PETSC_TRUE indicates you want the error generated 105 106 Level: advanced 107 108 .seealso: PCSetErrorIfFailure() 109 @*/ 110 PetscErrorCode MatSetErrorIfFailure(Mat mat,PetscBool flg) 111 { 112 PetscFunctionBegin; 113 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 114 PetscValidLogicalCollectiveBool(mat,flg,2); 115 mat->erroriffailure = flg; 116 PetscFunctionReturn(0); 117 } 118 119 /*@ 120 MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 121 122 Collective on Mat 123 124 Input Parameters: 125 + A - the matrix 126 . m - number of local rows (or PETSC_DECIDE) 127 . n - number of local columns (or PETSC_DECIDE) 128 . M - number of global rows (or PETSC_DETERMINE) 129 - N - number of global columns (or PETSC_DETERMINE) 130 131 Notes: 132 m (n) and M (N) cannot be both PETSC_DECIDE 133 If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang. 134 135 If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 136 user must ensure that they are chosen to be compatible with the 137 vectors. To do this, one first considers the matrix-vector product 138 'y = A x'. The 'm' that is used in the above routine must match the 139 local size used in the vector creation routine VecCreateMPI() for 'y'. 140 Likewise, the 'n' used must match that used as the local size in 141 VecCreateMPI() for 'x'. 142 143 You cannot change the sizes once they have been set. 144 145 The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called. 146 147 Level: beginner 148 149 .seealso: MatGetSize(), PetscSplitOwnership() 150 @*/ 151 PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) 152 { 153 PetscFunctionBegin; 154 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 155 PetscValidLogicalCollectiveInt(A,M,4); 156 PetscValidLogicalCollectiveInt(A,N,5); 157 if (M > 0 && m > M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",m,M); 158 if (N > 0 && n > N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",n,N); 159 if ((A->rmap->n >= 0 && A->rmap->N >= 0) && (A->rmap->n != m || (M > 0 && A->rmap->N != M))) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset row sizes to %D local %D global after previously setting them to %D local %D global",m,M,A->rmap->n,A->rmap->N); 160 if ((A->cmap->n >= 0 && A->cmap->N >= 0) && (A->cmap->n != n || (N > 0 && A->cmap->N != N))) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset column sizes to %D local %D global after previously setting them to %D local %D global",n,N,A->cmap->n,A->cmap->N); 161 A->rmap->n = m; 162 A->cmap->n = n; 163 A->rmap->N = M > -1 ? M : A->rmap->N; 164 A->cmap->N = N > -1 ? N : A->cmap->N; 165 PetscFunctionReturn(0); 166 } 167 168 /*@ 169 MatSetFromOptions - Creates a matrix where the type is determined 170 from the options database. Generates a parallel MPI matrix if the 171 communicator has more than one processor. The default matrix type is 172 AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if 173 you do not select a type in the options database. 174 175 Collective on Mat 176 177 Input Parameter: 178 . A - the matrix 179 180 Options Database Keys: 181 + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 182 . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 183 . -mat_type seqdense - dense type, uses MatCreateSeqDense() 184 . -mat_type mpidense - dense type, uses MatCreateDense() 185 . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 186 - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 187 188 Even More Options Database Keys: 189 See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 190 for additional format-specific options. 191 192 Level: beginner 193 194 .seealso: MatCreateSeqAIJ((), MatCreateAIJ(), 195 MatCreateSeqDense(), MatCreateDense(), 196 MatCreateSeqBAIJ(), MatCreateBAIJ(), 197 MatCreateSeqSBAIJ(), MatCreateSBAIJ(), 198 MatConvert() 199 @*/ 200 PetscErrorCode MatSetFromOptions(Mat B) 201 { 202 PetscErrorCode ierr; 203 const char *deft = MATAIJ; 204 char type[256]; 205 PetscBool flg,set; 206 207 PetscFunctionBegin; 208 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 209 210 ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr); 211 212 if (B->rmap->bs < 0) { 213 PetscInt newbs = -1; 214 ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr); 215 if (flg) { 216 ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr); 217 ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr); 218 } 219 } 220 221 ierr = PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr); 222 if (flg) { 223 ierr = MatSetType(B,type);CHKERRQ(ierr); 224 } else if (!((PetscObject)B)->type_name) { 225 ierr = MatSetType(B,deft);CHKERRQ(ierr); 226 } 227 228 ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr); 229 ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL);CHKERRQ(ierr); 230 ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL);CHKERRQ(ierr); 231 ierr = PetscOptionsBool("-mat_error_if_failure","Generate an error if an error occurs when factoring the matrix","MatSetErrorIfFailure",B->erroriffailure,&B->erroriffailure,NULL);CHKERRQ(ierr); 232 233 if (B->ops->setfromoptions) { 234 ierr = (*B->ops->setfromoptions)(PetscOptionsObject,B);CHKERRQ(ierr); 235 } 236 237 flg = PETSC_FALSE; 238 ierr = 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);CHKERRQ(ierr); 239 if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);} 240 flg = PETSC_FALSE; 241 ierr = 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);CHKERRQ(ierr); 242 if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);} 243 flg = PETSC_FALSE; 244 ierr = 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);CHKERRQ(ierr); 245 if (set) {ierr = MatSetOption(B,MAT_IGNORE_ZERO_ENTRIES,flg);CHKERRQ(ierr);} 246 247 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 248 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B);CHKERRQ(ierr); 249 ierr = PetscOptionsEnd();CHKERRQ(ierr); 250 PetscFunctionReturn(0); 251 } 252 253 /*@C 254 MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 255 256 Collective on Mat 257 258 Input Arguments: 259 + A - matrix being preallocated 260 . bs - block size 261 . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 262 . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 263 . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 264 - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 265 266 Level: beginner 267 268 .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), 269 PetscSplitOwnership() 270 @*/ 271 PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[]) 272 { 273 PetscErrorCode ierr; 274 PetscInt cbs; 275 void (*aij)(void); 276 void (*is)(void); 277 void (*hyp)(void) = NULL; 278 279 PetscFunctionBegin; 280 if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 281 ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 282 } 283 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 284 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 285 ierr = MatGetBlockSizes(A,&bs,&cbs);CHKERRQ(ierr); 286 /* these routines assumes bs == cbs, this should be checked somehow */ 287 ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr); 288 ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr); 289 ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr); 290 ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr); 291 /* 292 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 293 good before going on with it. 294 */ 295 ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 296 ierr = PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is);CHKERRQ(ierr); 297 #if defined(PETSC_HAVE_HYPRE) 298 ierr = PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp);CHKERRQ(ierr); 299 #endif 300 if (!aij && !is && !hyp) { 301 ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 302 } 303 if (aij || is || hyp) { 304 if (bs == cbs && bs == 1) { 305 ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr); 306 ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 307 ierr = MatISSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 308 #if defined(PETSC_HAVE_HYPRE) 309 ierr = MatHYPRESetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 310 #endif 311 } else { /* Convert block-row precallocation to scalar-row */ 312 PetscInt i,m,*sdnnz,*sonnz; 313 ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 314 ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr); 315 for (i=0; i<m; i++) { 316 if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs; 317 if (onnz) sonnz[i] = onnz[i/bs] * cbs; 318 } 319 ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr); 320 ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 321 ierr = MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 322 #if defined(PETSC_HAVE_HYPRE) 323 ierr = MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 324 #endif 325 ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr); 326 } 327 } 328 PetscFunctionReturn(0); 329 } 330 331 /* 332 Merges some information from Cs header to A; the C object is then destroyed 333 334 This is somewhat different from MatHeaderReplace() it would be nice to merge the code 335 */ 336 PetscErrorCode MatHeaderMerge(Mat A,Mat *C) 337 { 338 PetscErrorCode ierr; 339 PetscInt refct; 340 PetscOps Abops; 341 struct _MatOps Aops; 342 char *mtype,*mname,*mprefix; 343 Mat_Product *product; 344 345 PetscFunctionBegin; 346 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 347 PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 348 if (A == *C) PetscFunctionReturn(0); 349 PetscCheckSameComm(A,1,*C,2); 350 /* save the parts of A we need */ 351 Abops = ((PetscObject)A)->bops[0]; 352 Aops = A->ops[0]; 353 refct = ((PetscObject)A)->refct; 354 mtype = ((PetscObject)A)->type_name; 355 mname = ((PetscObject)A)->name; 356 mprefix = ((PetscObject)A)->prefix; 357 product = A->product; 358 359 /* zero these so the destroy below does not free them */ 360 ((PetscObject)A)->type_name = NULL; 361 ((PetscObject)A)->name = NULL; 362 363 /* free all the interior data structures from mat */ 364 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 365 366 ierr = PetscFree(A->defaultvectype);CHKERRQ(ierr); 367 ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr); 368 ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr); 369 ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr); 370 ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr); 371 372 /* copy C over to A */ 373 ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 374 375 /* return the parts of A we saved */ 376 ((PetscObject)A)->bops[0] = Abops; 377 A->ops[0] = Aops; 378 ((PetscObject)A)->refct = refct; 379 ((PetscObject)A)->type_name = mtype; 380 ((PetscObject)A)->name = mname; 381 ((PetscObject)A)->prefix = mprefix; 382 A->product = product; 383 384 /* since these two are copied into A we do not want them destroyed in C */ 385 ((PetscObject)*C)->qlist = NULL; 386 ((PetscObject)*C)->olist = NULL; 387 388 ierr = PetscHeaderDestroy(C);CHKERRQ(ierr); 389 PetscFunctionReturn(0); 390 } 391 /* 392 Replace A's header with that of C; the C object is then destroyed 393 394 This is essentially code moved from MatDestroy() 395 396 This is somewhat different from MatHeaderMerge() it would be nice to merge the code 397 398 Used in DM hence is declared PETSC_EXTERN 399 */ 400 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C) 401 { 402 PetscErrorCode ierr; 403 PetscInt refct; 404 PetscObjectState state; 405 struct _p_Mat buffer; 406 MatStencilInfo stencil; 407 408 PetscFunctionBegin; 409 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 410 PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 411 if (A == *C) PetscFunctionReturn(0); 412 PetscCheckSameComm(A,1,*C,2); 413 if (((PetscObject)*C)->refct != 1) SETERRQ1(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Object C has refct %D > 1, would leave hanging reference",((PetscObject)*C)->refct); 414 415 /* swap C and A */ 416 refct = ((PetscObject)A)->refct; 417 state = ((PetscObject)A)->state; 418 stencil = A->stencil; 419 ierr = PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));CHKERRQ(ierr); 420 ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 421 ierr = PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));CHKERRQ(ierr); 422 ((PetscObject)A)->refct = refct; 423 ((PetscObject)A)->state = state + 1; 424 A->stencil = stencil; 425 426 ((PetscObject)*C)->refct = 1; 427 ierr = MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL);CHKERRQ(ierr); 428 ierr = MatDestroy(C);CHKERRQ(ierr); 429 PetscFunctionReturn(0); 430 } 431 432 /*@ 433 MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 434 435 Input Parameters: 436 + A - the matrix 437 - flg - bind to the CPU if value of PETSC_TRUE 438 439 Level: intermediate 440 @*/ 441 PetscErrorCode MatBindToCPU(Mat A,PetscBool flg) 442 { 443 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 444 PetscErrorCode ierr; 445 446 PetscFunctionBegin; 447 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 448 PetscValidLogicalCollectiveBool(A,flg,2); 449 if (A->boundtocpu == flg) PetscFunctionReturn(0); 450 A->boundtocpu = flg; 451 if (A->ops->bindtocpu) { 452 ierr = (*A->ops->bindtocpu)(A,flg);CHKERRQ(ierr); 453 } 454 PetscFunctionReturn(0); 455 #else 456 PetscFunctionBegin; 457 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 458 PetscValidLogicalCollectiveBool(A,flg,2); 459 PetscFunctionReturn(0); 460 #endif 461 } 462 463 PetscErrorCode MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode) 464 { 465 IS is_coo_i,is_coo_j; 466 const PetscInt *coo_i,*coo_j; 467 PetscInt n,n_i,n_j; 468 PetscScalar zero = 0.; 469 PetscErrorCode ierr; 470 471 PetscFunctionBegin; 472 ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_i",(PetscObject*)&is_coo_i);CHKERRQ(ierr); 473 ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_j",(PetscObject*)&is_coo_j);CHKERRQ(ierr); 474 if (!is_coo_i) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_i IS"); 475 if (!is_coo_j) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_j IS"); 476 ierr = ISGetLocalSize(is_coo_i,&n_i);CHKERRQ(ierr); 477 ierr = ISGetLocalSize(is_coo_j,&n_j);CHKERRQ(ierr); 478 if (n_i != n_j) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_COR,"Wrong local size %D != %D",n_i,n_j); 479 ierr = ISGetIndices(is_coo_i,&coo_i);CHKERRQ(ierr); 480 ierr = ISGetIndices(is_coo_j,&coo_j);CHKERRQ(ierr); 481 if (imode != ADD_VALUES) { 482 ierr = MatZeroEntries(A);CHKERRQ(ierr); 483 } 484 for (n = 0; n < n_i; n++) { 485 ierr = MatSetValue(A,coo_i[n],coo_j[n],coo_v ? coo_v[n] : zero,ADD_VALUES);CHKERRQ(ierr); 486 } 487 ierr = ISRestoreIndices(is_coo_i,&coo_i);CHKERRQ(ierr); 488 ierr = ISRestoreIndices(is_coo_j,&coo_j);CHKERRQ(ierr); 489 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 490 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 491 PetscFunctionReturn(0); 492 } 493 494 PetscErrorCode MatSetPreallocationCOO_Basic(Mat A,PetscInt ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 495 { 496 Mat preallocator; 497 IS is_coo_i,is_coo_j; 498 PetscScalar zero = 0.0; 499 PetscInt n; 500 PetscErrorCode ierr; 501 502 PetscFunctionBegin; 503 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 504 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 505 ierr = MatCreate(PetscObjectComm((PetscObject)A),&preallocator);CHKERRQ(ierr); 506 ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 507 ierr = MatSetSizes(preallocator,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 508 ierr = MatSetLayouts(preallocator,A->rmap,A->cmap);CHKERRQ(ierr); 509 ierr = MatSetUp(preallocator);CHKERRQ(ierr); 510 for (n = 0; n < ncoo; n++) { 511 ierr = MatSetValue(preallocator,coo_i[n],coo_j[n],zero,INSERT_VALUES);CHKERRQ(ierr); 512 } 513 ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 514 ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 515 ierr = MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A);CHKERRQ(ierr); 516 ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 517 ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_i,PETSC_COPY_VALUES,&is_coo_i);CHKERRQ(ierr); 518 ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_j,PETSC_COPY_VALUES,&is_coo_j);CHKERRQ(ierr); 519 ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_i",(PetscObject)is_coo_i);CHKERRQ(ierr); 520 ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_j",(PetscObject)is_coo_j);CHKERRQ(ierr); 521 ierr = ISDestroy(&is_coo_i);CHKERRQ(ierr); 522 ierr = ISDestroy(&is_coo_j);CHKERRQ(ierr); 523 PetscFunctionReturn(0); 524 } 525 526 /*@C 527 MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries 528 529 Collective on Mat 530 531 Input Arguments: 532 + A - matrix being preallocated 533 . ncoo - number of entries in the locally owned part of the parallel matrix 534 . coo_i - row indices 535 - coo_j - column indices 536 537 Level: beginner 538 539 Notes: Entries can be repeated, see MatSetValuesCOO(). Currently optimized for cuSPARSE matrices only. 540 541 .seealso: MatSetValuesCOO(), MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation() 542 @*/ 543 PetscErrorCode MatSetPreallocationCOO(Mat A,PetscInt ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 544 { 545 PetscErrorCode (*f)(Mat,PetscInt,const PetscInt[],const PetscInt[]) = NULL; 546 PetscErrorCode ierr; 547 548 PetscFunctionBegin; 549 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 550 PetscValidType(A,1); 551 if (ncoo) PetscValidIntPointer(coo_i,3); 552 if (ncoo) PetscValidIntPointer(coo_j,4); 553 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 554 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 555 if (PetscDefined(USE_DEBUG)) { 556 PetscInt i; 557 for (i = 0; i < ncoo; i++) { 558 if (coo_i[i] < A->rmap->rstart || coo_i[i] >= A->rmap->rend) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Invalid row index %D! Must be in [%D,%D)",coo_i[i],A->rmap->rstart,A->rmap->rend); 559 if (coo_j[i] < 0 || coo_j[i] >= A->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Invalid col index %D! Must be in [0,%D)",coo_j[i],A->cmap->N); 560 } 561 } 562 ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOO_C",&f);CHKERRQ(ierr); 563 ierr = PetscLogEventBegin(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr); 564 if (f) { 565 ierr = (*f)(A,ncoo,coo_i,coo_j);CHKERRQ(ierr); 566 } else { /* allow fallback, very slow */ 567 ierr = MatSetPreallocationCOO_Basic(A,ncoo,coo_i,coo_j);CHKERRQ(ierr); 568 } 569 ierr = PetscLogEventEnd(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr); 570 PetscFunctionReturn(0); 571 } 572 573 /*@C 574 MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO() 575 576 Collective on Mat 577 578 Input Arguments: 579 + A - matrix being preallocated 580 . coo_v - the matrix values (can be NULL) 581 - imode - the insert mode 582 583 Level: beginner 584 585 Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO(). 586 When repeated entries are specified in the COO indices the coo_v values are first properly summed. 587 The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES). 588 Currently optimized for cuSPARSE matrices only. 589 Passing coo_v == NULL is equivalent to passing an array of zeros. 590 591 .seealso: MatSetPreallocationCOO(), InsertMode, INSERT_VALUES, ADD_VALUES 592 @*/ 593 PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode) 594 { 595 PetscErrorCode (*f)(Mat,const PetscScalar[],InsertMode) = NULL; 596 PetscErrorCode ierr; 597 598 PetscFunctionBegin; 599 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 600 PetscValidType(A,1); 601 MatCheckPreallocated(A,1); 602 PetscValidLogicalCollectiveEnum(A,imode,3); 603 ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetValuesCOO_C",&f);CHKERRQ(ierr); 604 ierr = PetscLogEventBegin(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr); 605 if (f) { 606 ierr = (*f)(A,coo_v,imode);CHKERRQ(ierr); 607 } else { /* allow fallback */ 608 ierr = MatSetValuesCOO_Basic(A,coo_v,imode);CHKERRQ(ierr); 609 } 610 ierr = PetscLogEventEnd(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr); 611 ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr); 612 PetscFunctionReturn(0); 613 } 614