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