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 flg = PETSC_FALSE; 248 ierr = PetscOptionsBool("-mat_form_explicit_transpose","Hint to form an explicit transpose for operations like MatMultTranspose","MatSetOption",flg,&flg,&set);CHKERRQ(ierr); 249 if (set) {ierr = MatSetOption(B,MAT_FORM_EXPLICIT_TRANSPOSE,flg);CHKERRQ(ierr);} 250 251 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 252 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B);CHKERRQ(ierr); 253 ierr = PetscOptionsEnd();CHKERRQ(ierr); 254 PetscFunctionReturn(0); 255 } 256 257 /*@C 258 MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 259 260 Collective on Mat 261 262 Input Arguments: 263 + A - matrix being preallocated 264 . bs - block size 265 . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 266 . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 267 . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 268 - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 269 270 Level: beginner 271 272 .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), 273 PetscSplitOwnership() 274 @*/ 275 PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[]) 276 { 277 PetscErrorCode ierr; 278 PetscInt cbs; 279 void (*aij)(void); 280 void (*is)(void); 281 void (*hyp)(void) = NULL; 282 283 PetscFunctionBegin; 284 if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 285 ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 286 } 287 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 288 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 289 ierr = MatGetBlockSizes(A,&bs,&cbs);CHKERRQ(ierr); 290 /* these routines assumes bs == cbs, this should be checked somehow */ 291 ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr); 292 ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr); 293 ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr); 294 ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr); 295 /* 296 In general, we have to do extra work to preallocate for scalar (AIJ) or unassembled (IS) matrices so we check whether it will do any 297 good before going on with it. 298 */ 299 ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 300 ierr = PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is);CHKERRQ(ierr); 301 #if defined(PETSC_HAVE_HYPRE) 302 ierr = PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp);CHKERRQ(ierr); 303 #endif 304 if (!aij && !is && !hyp) { 305 ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 306 } 307 if (aij || is || hyp) { 308 if (bs == cbs && bs == 1) { 309 ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr); 310 ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 311 ierr = MatISSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 312 #if defined(PETSC_HAVE_HYPRE) 313 ierr = MatHYPRESetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 314 #endif 315 } else { /* Convert block-row precallocation to scalar-row */ 316 PetscInt i,m,*sdnnz,*sonnz; 317 ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 318 ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr); 319 for (i=0; i<m; i++) { 320 if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs; 321 if (onnz) sonnz[i] = onnz[i/bs] * cbs; 322 } 323 ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr); 324 ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 325 ierr = MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 326 #if defined(PETSC_HAVE_HYPRE) 327 ierr = MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 328 #endif 329 ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr); 330 } 331 } 332 PetscFunctionReturn(0); 333 } 334 335 /* 336 Merges some information from Cs header to A; the C object is then destroyed 337 338 This is somewhat different from MatHeaderReplace() it would be nice to merge the code 339 */ 340 PetscErrorCode MatHeaderMerge(Mat A,Mat *C) 341 { 342 PetscErrorCode ierr; 343 PetscInt refct; 344 PetscOps Abops; 345 struct _MatOps Aops; 346 char *mtype,*mname,*mprefix; 347 Mat_Product *product; 348 349 PetscFunctionBegin; 350 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 351 PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 352 if (A == *C) PetscFunctionReturn(0); 353 PetscCheckSameComm(A,1,*C,2); 354 /* save the parts of A we need */ 355 Abops = ((PetscObject)A)->bops[0]; 356 Aops = A->ops[0]; 357 refct = ((PetscObject)A)->refct; 358 mtype = ((PetscObject)A)->type_name; 359 mname = ((PetscObject)A)->name; 360 mprefix = ((PetscObject)A)->prefix; 361 product = A->product; 362 363 /* zero these so the destroy below does not free them */ 364 ((PetscObject)A)->type_name = NULL; 365 ((PetscObject)A)->name = NULL; 366 367 /* free all the interior data structures from mat */ 368 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 369 370 ierr = PetscFree(A->defaultvectype);CHKERRQ(ierr); 371 ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr); 372 ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr); 373 ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr); 374 ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr); 375 376 /* copy C over to A */ 377 ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 378 379 /* return the parts of A we saved */ 380 ((PetscObject)A)->bops[0] = Abops; 381 A->ops[0] = Aops; 382 ((PetscObject)A)->refct = refct; 383 ((PetscObject)A)->type_name = mtype; 384 ((PetscObject)A)->name = mname; 385 ((PetscObject)A)->prefix = mprefix; 386 A->product = product; 387 388 /* since these two are copied into A we do not want them destroyed in C */ 389 ((PetscObject)*C)->qlist = NULL; 390 ((PetscObject)*C)->olist = NULL; 391 392 ierr = PetscHeaderDestroy(C);CHKERRQ(ierr); 393 PetscFunctionReturn(0); 394 } 395 /* 396 Replace A's header with that of C; the C object is then destroyed 397 398 This is essentially code moved from MatDestroy() 399 400 This is somewhat different from MatHeaderMerge() it would be nice to merge the code 401 402 Used in DM hence is declared PETSC_EXTERN 403 */ 404 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C) 405 { 406 PetscErrorCode ierr; 407 PetscInt refct; 408 PetscObjectState state; 409 struct _p_Mat buffer; 410 MatStencilInfo stencil; 411 412 PetscFunctionBegin; 413 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 414 PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 415 if (A == *C) PetscFunctionReturn(0); 416 PetscCheckSameComm(A,1,*C,2); 417 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); 418 419 /* swap C and A */ 420 refct = ((PetscObject)A)->refct; 421 state = ((PetscObject)A)->state; 422 stencil = A->stencil; 423 ierr = PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));CHKERRQ(ierr); 424 ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 425 ierr = PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));CHKERRQ(ierr); 426 ((PetscObject)A)->refct = refct; 427 ((PetscObject)A)->state = state + 1; 428 A->stencil = stencil; 429 430 ((PetscObject)*C)->refct = 1; 431 ierr = MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL);CHKERRQ(ierr); 432 ierr = MatDestroy(C);CHKERRQ(ierr); 433 PetscFunctionReturn(0); 434 } 435 436 /*@ 437 MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 438 439 Input Parameters: 440 + A - the matrix 441 - flg - bind to the CPU if value of PETSC_TRUE 442 443 Level: intermediate 444 @*/ 445 PetscErrorCode MatBindToCPU(Mat A,PetscBool flg) 446 { 447 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 448 PetscErrorCode ierr; 449 450 PetscFunctionBegin; 451 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 452 PetscValidLogicalCollectiveBool(A,flg,2); 453 if (A->boundtocpu == flg) PetscFunctionReturn(0); 454 A->boundtocpu = flg; 455 if (A->ops->bindtocpu) { 456 ierr = (*A->ops->bindtocpu)(A,flg);CHKERRQ(ierr); 457 } 458 PetscFunctionReturn(0); 459 #else 460 PetscFunctionBegin; 461 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 462 PetscValidLogicalCollectiveBool(A,flg,2); 463 PetscFunctionReturn(0); 464 #endif 465 } 466 467 PetscErrorCode MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode) 468 { 469 IS is_coo_i,is_coo_j; 470 const PetscInt *coo_i,*coo_j; 471 PetscInt n,n_i,n_j; 472 PetscScalar zero = 0.; 473 PetscErrorCode ierr; 474 475 PetscFunctionBegin; 476 ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_i",(PetscObject*)&is_coo_i);CHKERRQ(ierr); 477 ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_j",(PetscObject*)&is_coo_j);CHKERRQ(ierr); 478 if (!is_coo_i) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_i IS"); 479 if (!is_coo_j) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_j IS"); 480 ierr = ISGetLocalSize(is_coo_i,&n_i);CHKERRQ(ierr); 481 ierr = ISGetLocalSize(is_coo_j,&n_j);CHKERRQ(ierr); 482 if (n_i != n_j) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_COR,"Wrong local size %D != %D",n_i,n_j); 483 ierr = ISGetIndices(is_coo_i,&coo_i);CHKERRQ(ierr); 484 ierr = ISGetIndices(is_coo_j,&coo_j);CHKERRQ(ierr); 485 if (imode != ADD_VALUES) { 486 ierr = MatZeroEntries(A);CHKERRQ(ierr); 487 } 488 for (n = 0; n < n_i; n++) { 489 ierr = MatSetValue(A,coo_i[n],coo_j[n],coo_v ? coo_v[n] : zero,ADD_VALUES);CHKERRQ(ierr); 490 } 491 ierr = ISRestoreIndices(is_coo_i,&coo_i);CHKERRQ(ierr); 492 ierr = ISRestoreIndices(is_coo_j,&coo_j);CHKERRQ(ierr); 493 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 494 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 495 PetscFunctionReturn(0); 496 } 497 498 PetscErrorCode MatSetPreallocationCOO_Basic(Mat A,PetscInt ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 499 { 500 Mat preallocator; 501 IS is_coo_i,is_coo_j; 502 PetscScalar zero = 0.0; 503 PetscInt n; 504 PetscErrorCode ierr; 505 506 PetscFunctionBegin; 507 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 508 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 509 ierr = MatCreate(PetscObjectComm((PetscObject)A),&preallocator);CHKERRQ(ierr); 510 ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 511 ierr = MatSetSizes(preallocator,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 512 ierr = MatSetLayouts(preallocator,A->rmap,A->cmap);CHKERRQ(ierr); 513 ierr = MatSetUp(preallocator);CHKERRQ(ierr); 514 for (n = 0; n < ncoo; n++) { 515 ierr = MatSetValue(preallocator,coo_i[n],coo_j[n],zero,INSERT_VALUES);CHKERRQ(ierr); 516 } 517 ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 518 ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 519 ierr = MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A);CHKERRQ(ierr); 520 ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 521 ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_i,PETSC_COPY_VALUES,&is_coo_i);CHKERRQ(ierr); 522 ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_j,PETSC_COPY_VALUES,&is_coo_j);CHKERRQ(ierr); 523 ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_i",(PetscObject)is_coo_i);CHKERRQ(ierr); 524 ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_j",(PetscObject)is_coo_j);CHKERRQ(ierr); 525 ierr = ISDestroy(&is_coo_i);CHKERRQ(ierr); 526 ierr = ISDestroy(&is_coo_j);CHKERRQ(ierr); 527 PetscFunctionReturn(0); 528 } 529 530 /*@C 531 MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries 532 533 Collective on Mat 534 535 Input Arguments: 536 + A - matrix being preallocated 537 . ncoo - number of entries in the locally owned part of the parallel matrix 538 . coo_i - row indices 539 - coo_j - column indices 540 541 Level: beginner 542 543 Notes: Entries can be repeated, see MatSetValuesCOO(). Currently optimized for cuSPARSE matrices only. 544 545 .seealso: MatSetValuesCOO(), MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation() 546 @*/ 547 PetscErrorCode MatSetPreallocationCOO(Mat A,PetscInt ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 548 { 549 PetscErrorCode (*f)(Mat,PetscInt,const PetscInt[],const PetscInt[]) = NULL; 550 PetscErrorCode ierr; 551 552 PetscFunctionBegin; 553 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 554 PetscValidType(A,1); 555 if (ncoo) PetscValidIntPointer(coo_i,3); 556 if (ncoo) PetscValidIntPointer(coo_j,4); 557 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 558 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 559 if (PetscDefined(USE_DEBUG)) { 560 PetscInt i; 561 for (i = 0; i < ncoo; i++) { 562 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); 563 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); 564 } 565 } 566 ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOO_C",&f);CHKERRQ(ierr); 567 ierr = PetscLogEventBegin(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr); 568 if (f) { 569 ierr = (*f)(A,ncoo,coo_i,coo_j);CHKERRQ(ierr); 570 } else { /* allow fallback, very slow */ 571 ierr = MatSetPreallocationCOO_Basic(A,ncoo,coo_i,coo_j);CHKERRQ(ierr); 572 } 573 ierr = PetscLogEventEnd(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr); 574 PetscFunctionReturn(0); 575 } 576 577 /*@C 578 MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO() 579 580 Collective on Mat 581 582 Input Arguments: 583 + A - matrix being preallocated 584 . coo_v - the matrix values (can be NULL) 585 - imode - the insert mode 586 587 Level: beginner 588 589 Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO(). 590 When repeated entries are specified in the COO indices the coo_v values are first properly summed. 591 The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES). 592 Currently optimized for cuSPARSE matrices only. 593 Passing coo_v == NULL is equivalent to passing an array of zeros. 594 595 .seealso: MatSetPreallocationCOO(), InsertMode, INSERT_VALUES, ADD_VALUES 596 @*/ 597 PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode) 598 { 599 PetscErrorCode (*f)(Mat,const PetscScalar[],InsertMode) = NULL; 600 PetscErrorCode ierr; 601 602 PetscFunctionBegin; 603 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 604 PetscValidType(A,1); 605 MatCheckPreallocated(A,1); 606 PetscValidLogicalCollectiveEnum(A,imode,3); 607 ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetValuesCOO_C",&f);CHKERRQ(ierr); 608 ierr = PetscLogEventBegin(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr); 609 if (f) { 610 ierr = (*f)(A,coo_v,imode);CHKERRQ(ierr); 611 } else { /* allow fallback */ 612 ierr = MatSetValuesCOO_Basic(A,coo_v,imode);CHKERRQ(ierr); 613 } 614 ierr = PetscLogEventEnd(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr); 615 ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr); 616 PetscFunctionReturn(0); 617 } 618