1 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 2 3 #include <../src/mat/impls/aij/seq/aij.h> 4 #include <../src/mat/impls/aij/mpi/mpiaij.h> 5 6 PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat, PetscInt rbs, PetscInt cbs) 7 { 8 PetscFunctionBegin; 9 if (!mat->preallocated) PetscFunctionReturn(PETSC_SUCCESS); 10 PetscCheck(mat->rmap->bs <= 0 || mat->rmap->bs == rbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot change row block size %" PetscInt_FMT " to %" PetscInt_FMT, mat->rmap->bs, rbs); 11 PetscCheck(mat->cmap->bs <= 0 || mat->cmap->bs == cbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Cannot change column block size %" PetscInt_FMT " to %" PetscInt_FMT, mat->cmap->bs, cbs); 12 PetscFunctionReturn(PETSC_SUCCESS); 13 } 14 15 PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y, PetscScalar a) 16 { 17 PetscInt i, start, end, oldValA = 0, oldValB = 0; 18 PetscScalar alpha = a; 19 PetscBool prevoption; 20 PetscBool isSeqAIJDerived, isMPIAIJDerived; // all classes sharing SEQAIJHEADER or MPIAIJHEADER 21 Mat A = NULL, B = NULL; 22 23 PetscFunctionBegin; 24 PetscCall(MatGetOption(Y, MAT_NO_OFF_PROC_ENTRIES, &prevoption)); 25 PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 26 PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isSeqAIJDerived, MATSEQAIJ, MATSEQBAIJ, MATSEQSBAIJ, "")); 27 PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &isMPIAIJDerived, MATMPIAIJ, MATMPIBAIJ, MATMPISBAIJ, "")); 28 29 if (isSeqAIJDerived) A = Y; 30 else if (isMPIAIJDerived) { 31 Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)Y->data; 32 A = mpiaij->A; 33 B = mpiaij->B; 34 } 35 36 if (A) { 37 oldValA = ((Mat_SeqAIJ *)(A->data))->nonew; 38 ((Mat_SeqAIJ *)(A->data))->nonew = 0; // so that new nonzero locations are allowed 39 } 40 if (B) { 41 oldValB = ((Mat_SeqAIJ *)(B->data))->nonew; 42 ((Mat_SeqAIJ *)(B->data))->nonew = 0; 43 } 44 45 PetscCall(MatGetOwnershipRange(Y, &start, &end)); 46 for (i = start; i < end; i++) { 47 if (i < Y->cmap->N) PetscCall(MatSetValues(Y, 1, &i, 1, &i, &alpha, ADD_VALUES)); 48 } 49 PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY)); 50 PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY)); 51 PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, prevoption)); 52 if (A) ((Mat_SeqAIJ *)(A->data))->nonew = oldValA; 53 if (B) ((Mat_SeqAIJ *)(B->data))->nonew = oldValB; 54 PetscFunctionReturn(PETSC_SUCCESS); 55 } 56 57 /*@ 58 MatCreate - Creates a matrix where the type is determined 59 from either a call to `MatSetType()` or from the options database 60 with a call to `MatSetFromOptions()`. 61 62 Collective 63 64 Input Parameter: 65 . comm - MPI communicator 66 67 Output Parameter: 68 . A - the matrix 69 70 Options Database Keys: 71 + -mat_type seqaij - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()` 72 . -mat_type mpiaij - `MATMPIAIJ` type, uses `MatCreateAIJ()` 73 . -mat_type seqdense - `MATSEQDENSE`, uses `MatCreateSeqDense()` 74 . -mat_type mpidense - `MATMPIDENSE` type, uses `MatCreateDense()` 75 . -mat_type seqbaij - `MATSEQBAIJ` type, uses `MatCreateSeqBAIJ()` 76 - -mat_type mpibaij - `MATMPIBAIJ` type, uses `MatCreateBAIJ()` 77 78 See the manpages for particular formats (e.g., `MATSEQAIJ`) 79 for additional format-specific options. 80 81 Level: beginner 82 83 Notes: 84 The default matrix type is `MATAIJ`, using the routines `MatCreateSeqAIJ()` or 85 `MatCreateAIJ()` if you do not set a type in the options database. If you never call 86 `MatSetType()` or `MatSetFromOptions()` it will generate an error when you try to use the 87 matrix. 88 89 .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`, 90 `MatCreateSeqDense()`, `MatCreateDense()`, 91 `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, 92 `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`, 93 `MatConvert()` 94 @*/ 95 PetscErrorCode MatCreate(MPI_Comm comm, Mat *A) 96 { 97 Mat B; 98 99 PetscFunctionBegin; 100 PetscAssertPointer(A, 2); 101 102 *A = NULL; 103 PetscCall(MatInitializePackage()); 104 105 PetscCall(PetscHeaderCreate(B, MAT_CLASSID, "Mat", "Matrix", "Mat", comm, MatDestroy, MatView)); 106 PetscCall(PetscLayoutCreate(comm, &B->rmap)); 107 PetscCall(PetscLayoutCreate(comm, &B->cmap)); 108 PetscCall(PetscStrallocpy(VECSTANDARD, &B->defaultvectype)); 109 PetscCall(PetscStrallocpy(PETSCRANDER48, &B->defaultrandtype)); 110 111 B->symmetric = PETSC_BOOL3_UNKNOWN; 112 B->hermitian = PETSC_BOOL3_UNKNOWN; 113 B->structurally_symmetric = PETSC_BOOL3_UNKNOWN; 114 B->spd = PETSC_BOOL3_UNKNOWN; 115 B->symmetry_eternal = PETSC_FALSE; 116 B->structural_symmetry_eternal = PETSC_FALSE; 117 118 B->congruentlayouts = PETSC_DECIDE; 119 B->preallocated = PETSC_FALSE; 120 #if defined(PETSC_HAVE_DEVICE) 121 B->boundtocpu = PETSC_TRUE; 122 #endif 123 *A = B; 124 PetscFunctionReturn(PETSC_SUCCESS); 125 } 126 127 /*@C 128 MatCreateFromOptions - Creates a matrix whose type is set from the options database 129 130 Collective 131 132 Input Parameters: 133 + comm - MPI communicator 134 . prefix - [optional] prefix for the options database 135 . bs - the blocksize (commonly 1) 136 . m - the local number of rows (or `PETSC_DECIDE`) 137 . n - the local number of columns (or `PETSC_DECIDE` or `PETSC_DETERMINE`) 138 . M - the global number of rows (or `PETSC_DETERMINE`) 139 - N - the global number of columns (or `PETSC_DETERMINE`) 140 141 Output Parameter: 142 . A - the matrix 143 144 Options Database Key: 145 . -mat_type - see `MatType`, for example `aij`, `aijcusparse`, `baij`, `sbaij`, dense, defaults to `aij` 146 147 Level: beginner 148 149 .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`, 150 `MatCreateSeqDense()`, `MatCreateDense()`, 151 `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, 152 `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`, 153 `MatConvert()`, `MatCreate()` 154 @*/ 155 PetscErrorCode MatCreateFromOptions(MPI_Comm comm, const char *prefix, PetscInt bs, PetscInt m, PetscInt n, PetscInt M, PetscInt N, Mat *A) 156 { 157 PetscFunctionBegin; 158 PetscAssertPointer(A, 8); 159 PetscCall(MatCreate(comm, A)); 160 if (prefix) PetscCall(MatSetOptionsPrefix(*A, prefix)); 161 PetscCall(MatSetBlockSize(*A, bs)); 162 PetscCall(MatSetSizes(*A, m, n, M, N)); 163 PetscCall(MatSetFromOptions(*A)); 164 PetscFunctionReturn(PETSC_SUCCESS); 165 } 166 167 /*@ 168 MatSetErrorIfFailure - Causes `Mat` to generate an immediate error, for example a zero pivot, is detected. 169 170 Logically Collective 171 172 Input Parameters: 173 + mat - matrix obtained from `MatCreate()` 174 - flg - `PETSC_TRUE` indicates you want the error generated 175 176 Level: advanced 177 178 Note: 179 If this flag is not set then the matrix operation will note the error and continue. The error may cause a later `PC` or `KSP` error 180 or result in a `KSPConvergedReason` indicating the method did not converge. 181 182 .seealso: [](ch_matrices), `Mat`, `PCSetErrorIfFailure()`, `KSPConvergedReason`, `SNESConvergedReason` 183 @*/ 184 PetscErrorCode MatSetErrorIfFailure(Mat mat, PetscBool flg) 185 { 186 PetscFunctionBegin; 187 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 188 PetscValidLogicalCollectiveBool(mat, flg, 2); 189 mat->erroriffailure = flg; 190 PetscFunctionReturn(PETSC_SUCCESS); 191 } 192 193 /*@ 194 MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 195 196 Collective 197 198 Input Parameters: 199 + A - the matrix 200 . m - number of local rows (or `PETSC_DECIDE`) 201 . n - number of local columns (or `PETSC_DECIDE`) 202 . M - number of global rows (or `PETSC_DETERMINE`) 203 - N - number of global columns (or `PETSC_DETERMINE`) 204 205 Level: beginner 206 207 Notes: 208 `m` (`n`) and `M` (`N`) cannot be both `PETSC_DECIDE` 209 If one processor calls this with `M` (`N`) of `PETSC_DECIDE` then all processors must, otherwise the program will hang. 210 211 If `PETSC_DECIDE` is not used for the arguments 'm' and 'n', then the 212 user must ensure that they are chosen to be compatible with the 213 vectors. To do this, one first considers the matrix-vector product 214 'y = A x'. The `m` that is used in the above routine must match the 215 local size used in the vector creation routine `VecCreateMPI()` for 'y'. 216 Likewise, the `n` used must match that used as the local size in 217 `VecCreateMPI()` for 'x'. 218 219 You cannot change the sizes once they have been set. 220 221 The sizes must be set before `MatSetUp()` or MatXXXSetPreallocation() is called. 222 223 .seealso: [](ch_matrices), `Mat`, `MatGetSize()`, `PetscSplitOwnership()` 224 @*/ 225 PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) 226 { 227 PetscFunctionBegin; 228 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 229 PetscValidLogicalCollectiveInt(A, M, 4); 230 PetscValidLogicalCollectiveInt(A, N, 5); 231 PetscCheck(M <= 0 || m <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local row size %" PetscInt_FMT " cannot be larger than global row size %" PetscInt_FMT, m, M); 232 PetscCheck(N <= 0 || n <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local column size %" PetscInt_FMT " cannot be larger than global column size %" PetscInt_FMT, n, N); 233 PetscCheck((A->rmap->n < 0 || A->rmap->N < 0) || (A->rmap->n == m && (M <= 0 || A->rmap->N == M)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot change/reset row sizes to %" PetscInt_FMT " local %" PetscInt_FMT " global after previously setting them to %" PetscInt_FMT " local %" PetscInt_FMT " global", m, M, 234 A->rmap->n, A->rmap->N); 235 PetscCheck((A->cmap->n < 0 || A->cmap->N < 0) || (A->cmap->n == n && (N <= 0 || A->cmap->N == N)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot change/reset column sizes to %" PetscInt_FMT " local %" PetscInt_FMT " global after previously setting them to %" PetscInt_FMT " local %" PetscInt_FMT " global", n, N, 236 A->cmap->n, A->cmap->N); 237 A->rmap->n = m; 238 A->cmap->n = n; 239 A->rmap->N = M > -1 ? M : A->rmap->N; 240 A->cmap->N = N > -1 ? N : A->cmap->N; 241 PetscFunctionReturn(PETSC_SUCCESS); 242 } 243 244 /*@ 245 MatSetFromOptions - Creates a matrix where the type is determined 246 from the options database. 247 248 Collective 249 250 Input Parameter: 251 . B - the matrix 252 253 Options Database Keys: 254 + -mat_type seqaij - `MATSEQAIJ` type, uses `MatCreateSeqAIJ()` 255 . -mat_type mpiaij - `MATMPIAIJ` type, uses `MatCreateAIJ()` 256 . -mat_type seqdense - `MATSEQDENSE` type, uses `MatCreateSeqDense()` 257 . -mat_type mpidense - `MATMPIDENSE`, uses `MatCreateDense()` 258 . -mat_type seqbaij - `MATSEQBAIJ`, uses `MatCreateSeqBAIJ()` 259 - -mat_type mpibaij - `MATMPIBAIJ`, uses `MatCreateBAIJ()` 260 261 See the manpages for particular formats (e.g., `MATSEQAIJ`) 262 for additional format-specific options. 263 264 Level: beginner 265 266 Notes: 267 Generates a parallel MPI matrix if the communicator has more than one processor. The default 268 matrix type is `MATAIJ`, using the routines `MatCreateSeqAIJ()` and `MatCreateAIJ()` if you 269 do not select a type in the options database. 270 271 .seealso: [](ch_matrices), `Mat`, `MatCreateSeqAIJ()`, `MatCreateAIJ()`, 272 `MatCreateSeqDense()`, `MatCreateDense()`, 273 `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, 274 `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`, 275 `MatConvert()` 276 @*/ 277 PetscErrorCode MatSetFromOptions(Mat B) 278 { 279 const char *deft = MATAIJ; 280 char type[256]; 281 PetscBool flg, set; 282 PetscInt bind_below = 0; 283 284 PetscFunctionBegin; 285 PetscValidHeaderSpecific(B, MAT_CLASSID, 1); 286 287 PetscObjectOptionsBegin((PetscObject)B); 288 289 if (B->rmap->bs < 0) { 290 PetscInt newbs = -1; 291 PetscCall(PetscOptionsInt("-mat_block_size", "Set the blocksize used to store the matrix", "MatSetBlockSize", newbs, &newbs, &flg)); 292 if (flg) { 293 PetscCall(PetscLayoutSetBlockSize(B->rmap, newbs)); 294 PetscCall(PetscLayoutSetBlockSize(B->cmap, newbs)); 295 } 296 } 297 298 PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, deft, type, 256, &flg)); 299 if (flg) { 300 PetscCall(MatSetType(B, type)); 301 } else if (!((PetscObject)B)->type_name) { 302 PetscCall(MatSetType(B, deft)); 303 } 304 305 PetscCall(PetscOptionsName("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", &B->checksymmetryonassembly)); 306 PetscCall(PetscOptionsReal("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", B->checksymmetrytol, &B->checksymmetrytol, NULL)); 307 PetscCall(PetscOptionsBool("-mat_null_space_test", "Checks if provided null space is correct in MatAssemblyEnd()", "MatSetNullSpaceTest", B->checknullspaceonassembly, &B->checknullspaceonassembly, NULL)); 308 PetscCall(PetscOptionsBool("-mat_error_if_failure", "Generate an error if an error occurs when factoring the matrix", "MatSetErrorIfFailure", B->erroriffailure, &B->erroriffailure, NULL)); 309 310 PetscTryTypeMethod(B, setfromoptions, PetscOptionsObject); 311 312 flg = PETSC_FALSE; 313 PetscCall(PetscOptionsBool("-mat_new_nonzero_location_err", "Generate an error if new nonzeros are created in the matrix structure (useful to test preallocation)", "MatSetOption", flg, &flg, &set)); 314 if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, flg)); 315 flg = PETSC_FALSE; 316 PetscCall(PetscOptionsBool("-mat_new_nonzero_allocation_err", "Generate an error if new nonzeros are allocated in the matrix structure (useful to test preallocation)", "MatSetOption", flg, &flg, &set)); 317 if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, flg)); 318 flg = PETSC_FALSE; 319 PetscCall(PetscOptionsBool("-mat_ignore_zero_entries", "For AIJ/IS matrices this will stop zero values from creating a zero location in the matrix", "MatSetOption", flg, &flg, &set)); 320 if (set) PetscCall(MatSetOption(B, MAT_IGNORE_ZERO_ENTRIES, flg)); 321 322 flg = PETSC_FALSE; 323 PetscCall(PetscOptionsBool("-mat_form_explicit_transpose", "Hint to form an explicit transpose for operations like MatMultTranspose", "MatSetOption", flg, &flg, &set)); 324 if (set) PetscCall(MatSetOption(B, MAT_FORM_EXPLICIT_TRANSPOSE, flg)); 325 326 /* Bind to CPU if below a user-specified size threshold. 327 * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types, 328 * and putting it here makes is more maintainable than duplicating this for all. */ 329 PetscCall(PetscOptionsInt("-mat_bind_below", "Set the size threshold (in local rows) below which the Mat is bound to the CPU", "MatBindToCPU", bind_below, &bind_below, &flg)); 330 if (flg && B->rmap->n < bind_below) PetscCall(MatBindToCPU(B, PETSC_TRUE)); 331 332 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 333 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)B, PetscOptionsObject)); 334 PetscOptionsEnd(); 335 PetscFunctionReturn(PETSC_SUCCESS); 336 } 337 338 /*@C 339 MatXAIJSetPreallocation - set preallocation for serial and parallel `MATAIJ`, `MATBAIJ`, and `MATSBAIJ` matrices and their unassembled versions. 340 341 Collective 342 343 Input Parameters: 344 + A - matrix being preallocated 345 . bs - block size 346 . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 347 . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 348 . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 349 - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 350 351 Level: beginner 352 353 .seealso: [](ch_matrices), `Mat`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, 354 `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, 355 `PetscSplitOwnership()` 356 @*/ 357 PetscErrorCode MatXAIJSetPreallocation(Mat A, PetscInt bs, const PetscInt dnnz[], const PetscInt onnz[], const PetscInt dnnzu[], const PetscInt onnzu[]) 358 { 359 PetscInt cbs; 360 void (*aij)(void); 361 void (*is)(void); 362 void (*hyp)(void) = NULL; 363 364 PetscFunctionBegin; 365 if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 366 PetscCall(MatSetBlockSize(A, bs)); 367 } 368 PetscCall(PetscLayoutSetUp(A->rmap)); 369 PetscCall(PetscLayoutSetUp(A->cmap)); 370 PetscCall(MatGetBlockSizes(A, &bs, &cbs)); 371 /* these routines assumes bs == cbs, this should be checked somehow */ 372 PetscCall(MatSeqBAIJSetPreallocation(A, bs, 0, dnnz)); 373 PetscCall(MatMPIBAIJSetPreallocation(A, bs, 0, dnnz, 0, onnz)); 374 PetscCall(MatSeqSBAIJSetPreallocation(A, bs, 0, dnnzu)); 375 PetscCall(MatMPISBAIJSetPreallocation(A, bs, 0, dnnzu, 0, onnzu)); 376 /* 377 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 378 good before going on with it. 379 */ 380 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatMPIAIJSetPreallocation_C", &aij)); 381 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatISSetPreallocation_C", &is)); 382 #if defined(PETSC_HAVE_HYPRE) 383 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatHYPRESetPreallocation_C", &hyp)); 384 #endif 385 if (!aij && !is && !hyp) PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", &aij)); 386 if (aij || is || hyp) { 387 if (bs == cbs && bs == 1) { 388 PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz)); 389 PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz, 0, onnz)); 390 PetscCall(MatISSetPreallocation(A, 0, dnnz, 0, onnz)); 391 #if defined(PETSC_HAVE_HYPRE) 392 PetscCall(MatHYPRESetPreallocation(A, 0, dnnz, 0, onnz)); 393 #endif 394 } else { /* Convert block-row precallocation to scalar-row */ 395 PetscInt i, m, *sdnnz, *sonnz; 396 PetscCall(MatGetLocalSize(A, &m, NULL)); 397 PetscCall(PetscMalloc2((!!dnnz) * m, &sdnnz, (!!onnz) * m, &sonnz)); 398 for (i = 0; i < m; i++) { 399 if (dnnz) sdnnz[i] = dnnz[i / bs] * cbs; 400 if (onnz) sonnz[i] = onnz[i / bs] * cbs; 401 } 402 PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL)); 403 PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 404 PetscCall(MatISSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 405 #if defined(PETSC_HAVE_HYPRE) 406 PetscCall(MatHYPRESetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 407 #endif 408 PetscCall(PetscFree2(sdnnz, sonnz)); 409 } 410 } 411 PetscFunctionReturn(PETSC_SUCCESS); 412 } 413 414 /* 415 Merges some information from Cs header to A; the C object is then destroyed 416 417 This is somewhat different from MatHeaderReplace() it would be nice to merge the code 418 */ 419 PetscErrorCode MatHeaderMerge(Mat A, Mat *C) 420 { 421 PetscInt refct; 422 PetscOps Abops; 423 struct _MatOps Aops; 424 char *mtype, *mname, *mprefix; 425 Mat_Product *product; 426 Mat_Redundant *redundant; 427 PetscObjectState state; 428 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 431 PetscValidHeaderSpecific(*C, MAT_CLASSID, 2); 432 if (A == *C) PetscFunctionReturn(PETSC_SUCCESS); 433 PetscCheckSameComm(A, 1, *C, 2); 434 /* save the parts of A we need */ 435 Abops = ((PetscObject)A)->bops[0]; 436 Aops = A->ops[0]; 437 refct = ((PetscObject)A)->refct; 438 mtype = ((PetscObject)A)->type_name; 439 mname = ((PetscObject)A)->name; 440 state = ((PetscObject)A)->state; 441 mprefix = ((PetscObject)A)->prefix; 442 product = A->product; 443 redundant = A->redundant; 444 445 /* zero these so the destroy below does not free them */ 446 ((PetscObject)A)->type_name = NULL; 447 ((PetscObject)A)->name = NULL; 448 449 /* 450 free all the interior data structures from mat 451 cannot use PetscUseTypeMethod(A,destroy); because compiler 452 thinks it may print NULL type_name and name 453 */ 454 PetscTryTypeMethod(A, destroy); 455 456 PetscCall(PetscFree(A->defaultvectype)); 457 PetscCall(PetscFree(A->defaultrandtype)); 458 PetscCall(PetscLayoutDestroy(&A->rmap)); 459 PetscCall(PetscLayoutDestroy(&A->cmap)); 460 PetscCall(PetscFunctionListDestroy(&((PetscObject)A)->qlist)); 461 PetscCall(PetscObjectListDestroy(&((PetscObject)A)->olist)); 462 PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A)); 463 464 /* copy C over to A */ 465 PetscCall(PetscFree(A->factorprefix)); 466 PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat))); 467 468 /* return the parts of A we saved */ 469 ((PetscObject)A)->bops[0] = Abops; 470 A->ops[0] = Aops; 471 ((PetscObject)A)->refct = refct; 472 ((PetscObject)A)->type_name = mtype; 473 ((PetscObject)A)->name = mname; 474 ((PetscObject)A)->prefix = mprefix; 475 ((PetscObject)A)->state = state + 1; 476 A->product = product; 477 A->redundant = redundant; 478 479 /* since these two are copied into A we do not want them destroyed in C */ 480 ((PetscObject)*C)->qlist = NULL; 481 ((PetscObject)*C)->olist = NULL; 482 483 PetscCall(PetscHeaderDestroy(C)); 484 PetscFunctionReturn(PETSC_SUCCESS); 485 } 486 /* 487 Replace A's header with that of C; the C object is then destroyed 488 489 This is essentially code moved from MatDestroy() 490 491 This is somewhat different from MatHeaderMerge() it would be nice to merge the code 492 493 Used in DM hence is declared PETSC_EXTERN 494 */ 495 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A, Mat *C) 496 { 497 PetscInt refct; 498 PetscObjectState state; 499 struct _p_Mat buffer; 500 MatStencilInfo stencil; 501 502 PetscFunctionBegin; 503 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 504 PetscValidHeaderSpecific(*C, MAT_CLASSID, 2); 505 if (A == *C) PetscFunctionReturn(PETSC_SUCCESS); 506 PetscCheckSameComm(A, 1, *C, 2); 507 PetscCheck(((PetscObject)*C)->refct == 1, PetscObjectComm((PetscObject)C), PETSC_ERR_ARG_WRONGSTATE, "Object C has refct %" PetscInt_FMT " > 1, would leave hanging reference", ((PetscObject)*C)->refct); 508 509 /* swap C and A */ 510 refct = ((PetscObject)A)->refct; 511 state = ((PetscObject)A)->state; 512 stencil = A->stencil; 513 PetscCall(PetscMemcpy(&buffer, A, sizeof(struct _p_Mat))); 514 PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat))); 515 PetscCall(PetscMemcpy(*C, &buffer, sizeof(struct _p_Mat))); 516 ((PetscObject)A)->refct = refct; 517 ((PetscObject)A)->state = state + 1; 518 A->stencil = stencil; 519 520 ((PetscObject)*C)->refct = 1; 521 PetscCall(MatShellSetOperation(*C, MATOP_DESTROY, (void (*)(void))NULL)); 522 PetscCall(MatDestroy(C)); 523 PetscFunctionReturn(PETSC_SUCCESS); 524 } 525 526 /*@ 527 MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 528 529 Logically Collective 530 531 Input Parameters: 532 + A - the matrix 533 - flg - bind to the CPU if value of `PETSC_TRUE` 534 535 Level: intermediate 536 537 .seealso: [](ch_matrices), `Mat`, `MatBoundToCPU()` 538 @*/ 539 PetscErrorCode MatBindToCPU(Mat A, PetscBool flg) 540 { 541 PetscFunctionBegin; 542 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 543 PetscValidLogicalCollectiveBool(A, flg, 2); 544 #if defined(PETSC_HAVE_DEVICE) 545 if (A->boundtocpu == flg) PetscFunctionReturn(PETSC_SUCCESS); 546 A->boundtocpu = flg; 547 PetscTryTypeMethod(A, bindtocpu, flg); 548 #endif 549 PetscFunctionReturn(PETSC_SUCCESS); 550 } 551 552 /*@ 553 MatBoundToCPU - query if a matrix is bound to the CPU 554 555 Input Parameter: 556 . A - the matrix 557 558 Output Parameter: 559 . flg - the logical flag 560 561 Level: intermediate 562 563 .seealso: [](ch_matrices), `Mat`, `MatBindToCPU()` 564 @*/ 565 PetscErrorCode MatBoundToCPU(Mat A, PetscBool *flg) 566 { 567 PetscFunctionBegin; 568 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 569 PetscAssertPointer(flg, 2); 570 #if defined(PETSC_HAVE_DEVICE) 571 *flg = A->boundtocpu; 572 #else 573 *flg = PETSC_TRUE; 574 #endif 575 PetscFunctionReturn(PETSC_SUCCESS); 576 } 577 578 PetscErrorCode MatSetValuesCOO_Basic(Mat A, const PetscScalar coo_v[], InsertMode imode) 579 { 580 IS is_coo_i, is_coo_j; 581 const PetscInt *coo_i, *coo_j; 582 PetscInt n, n_i, n_j; 583 PetscScalar zero = 0.; 584 585 PetscFunctionBegin; 586 PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_i", (PetscObject *)&is_coo_i)); 587 PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_j", (PetscObject *)&is_coo_j)); 588 PetscCheck(is_coo_i, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_i IS"); 589 PetscCheck(is_coo_j, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_j IS"); 590 PetscCall(ISGetLocalSize(is_coo_i, &n_i)); 591 PetscCall(ISGetLocalSize(is_coo_j, &n_j)); 592 PetscCheck(n_i == n_j, PETSC_COMM_SELF, PETSC_ERR_COR, "Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT, n_i, n_j); 593 PetscCall(ISGetIndices(is_coo_i, &coo_i)); 594 PetscCall(ISGetIndices(is_coo_j, &coo_j)); 595 if (imode != ADD_VALUES) PetscCall(MatZeroEntries(A)); 596 for (n = 0; n < n_i; n++) PetscCall(MatSetValue(A, coo_i[n], coo_j[n], coo_v ? coo_v[n] : zero, ADD_VALUES)); 597 PetscCall(ISRestoreIndices(is_coo_i, &coo_i)); 598 PetscCall(ISRestoreIndices(is_coo_j, &coo_j)); 599 PetscFunctionReturn(PETSC_SUCCESS); 600 } 601 602 PetscErrorCode MatSetPreallocationCOO_Basic(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) 603 { 604 Mat preallocator; 605 IS is_coo_i, is_coo_j; 606 PetscScalar zero = 0.0; 607 608 PetscFunctionBegin; 609 PetscCall(PetscLayoutSetUp(A->rmap)); 610 PetscCall(PetscLayoutSetUp(A->cmap)); 611 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator)); 612 PetscCall(MatSetType(preallocator, MATPREALLOCATOR)); 613 PetscCall(MatSetSizes(preallocator, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 614 PetscCall(MatSetLayouts(preallocator, A->rmap, A->cmap)); 615 PetscCall(MatSetUp(preallocator)); 616 for (PetscCount n = 0; n < ncoo; n++) PetscCall(MatSetValue(preallocator, coo_i[n], coo_j[n], zero, INSERT_VALUES)); 617 PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY)); 618 PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY)); 619 PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A)); 620 PetscCall(MatDestroy(&preallocator)); 621 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); 622 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo, coo_i, PETSC_COPY_VALUES, &is_coo_i)); 623 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo, coo_j, PETSC_COPY_VALUES, &is_coo_j)); 624 PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_i", (PetscObject)is_coo_i)); 625 PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_j", (PetscObject)is_coo_j)); 626 PetscCall(ISDestroy(&is_coo_i)); 627 PetscCall(ISDestroy(&is_coo_j)); 628 PetscFunctionReturn(PETSC_SUCCESS); 629 } 630 631 /*@C 632 MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices 633 634 Collective 635 636 Input Parameters: 637 + A - matrix being preallocated 638 . ncoo - number of entries 639 . coo_i - row indices 640 - coo_j - column indices 641 642 Level: beginner 643 644 Notes: 645 The indices `coo_i` and `coo_j` may be modified within this function. The caller should not rely on them 646 having any specific value after this function returns. The arrays can be freed or reused immediately 647 after this function returns. 648 649 Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed 650 but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries 651 are allowed and will be properly added or inserted to the matrix, unless the matrix option `MAT_IGNORE_OFF_PROC_ENTRIES` 652 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. 653 654 If you just want to create a sequential AIJ matrix (`MATSEQAIJ`), and your matrix entries in COO format are unique, you can also use 655 `MatCreateSeqAIJFromTriple()`. But that is not recommended for iterative applications. 656 657 .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, 658 `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`, 659 `DMSetMatrixPreallocateSkip()`, `MatCreateSeqAIJFromTriple()` 660 @*/ 661 PetscErrorCode MatSetPreallocationCOO(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) 662 { 663 PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL; 664 665 PetscFunctionBegin; 666 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 667 PetscValidType(A, 1); 668 if (ncoo) PetscAssertPointer(coo_i, 3); 669 if (ncoo) PetscAssertPointer(coo_j, 4); 670 PetscCall(PetscLayoutSetUp(A->rmap)); 671 PetscCall(PetscLayoutSetUp(A->cmap)); 672 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOO_C", &f)); 673 674 PetscCall(PetscLogEventBegin(MAT_PreallCOO, A, 0, 0, 0)); 675 if (f) { 676 PetscCall((*f)(A, ncoo, coo_i, coo_j)); 677 } else { /* allow fallback, very slow */ 678 PetscCall(MatSetPreallocationCOO_Basic(A, ncoo, coo_i, coo_j)); 679 } 680 PetscCall(PetscLogEventEnd(MAT_PreallCOO, A, 0, 0, 0)); 681 A->preallocated = PETSC_TRUE; 682 A->nonzerostate++; 683 PetscFunctionReturn(PETSC_SUCCESS); 684 } 685 686 /*@C 687 MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices 688 689 Collective 690 691 Input Parameters: 692 + A - matrix being preallocated 693 . ncoo - number of entries 694 . coo_i - row indices (local numbering; may be modified) 695 - coo_j - column indices (local numbering; may be modified) 696 697 Level: beginner 698 699 Notes: 700 The local indices are translated using the local to global mapping, thus `MatSetLocalToGlobalMapping()` must have been 701 called prior to this function. For matrices created with `DMCreateMatrix()` the local to global mapping is often already provided. 702 703 The indices `coo_i` and `coo_j` may be modified within this function. They might be translated to corresponding global 704 indices, but the caller should not rely on them having any specific value after this function returns. The arrays 705 can be freed or reused immediately after this function returns. 706 707 Entries can be repeated, see `MatSetValuesCOO()`. Entries with negative row or column indices are allowed 708 but will be ignored. The corresponding entries in `MatSetValuesCOO()` will be ignored too. Remote entries 709 are allowed and will be properly added or inserted to the matrix. 710 711 .seealso: [](ch_matrices), `Mat`, `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, 712 `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`, 713 `DMSetMatrixPreallocateSkip()` 714 @*/ 715 PetscErrorCode MatSetPreallocationCOOLocal(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) 716 { 717 PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL; 718 719 PetscFunctionBegin; 720 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 721 PetscValidType(A, 1); 722 if (ncoo) PetscAssertPointer(coo_i, 3); 723 if (ncoo) PetscAssertPointer(coo_j, 4); 724 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); 725 PetscCall(PetscLayoutSetUp(A->rmap)); 726 PetscCall(PetscLayoutSetUp(A->cmap)); 727 728 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", &f)); 729 if (f) { 730 PetscCall((*f)(A, ncoo, coo_i, coo_j)); 731 A->nonzerostate++; 732 } else { 733 ISLocalToGlobalMapping ltog_row, ltog_col; 734 PetscCall(MatGetLocalToGlobalMapping(A, <og_row, <og_col)); 735 if (ltog_row) PetscCall(ISLocalToGlobalMappingApply(ltog_row, ncoo, coo_i, coo_i)); 736 if (ltog_col) PetscCall(ISLocalToGlobalMappingApply(ltog_col, ncoo, coo_j, coo_j)); 737 PetscCall(MatSetPreallocationCOO(A, ncoo, coo_i, coo_j)); 738 } 739 A->preallocated = PETSC_TRUE; 740 PetscFunctionReturn(PETSC_SUCCESS); 741 } 742 743 /*@ 744 MatSetValuesCOO - set values at once in a matrix preallocated using `MatSetPreallocationCOO()` 745 746 Collective 747 748 Input Parameters: 749 + A - matrix being preallocated 750 . coo_v - the matrix values (can be `NULL`) 751 - imode - the insert mode 752 753 Level: beginner 754 755 Notes: 756 The values must follow the order of the indices prescribed with `MatSetPreallocationCOO()` or `MatSetPreallocationCOOLocal()`. 757 758 When repeated entries are specified in the COO indices the `coo_v` values are first properly summed, regardless of the value of imode. 759 The imode flag indicates if coo_v must be added to the current values of the matrix (`ADD_VALUES`) or overwritten (`INSERT_VALUES`). 760 761 `MatAssemblyBegin()` and `MatAssemblyEnd()` do not need to be called after this routine. It automatically handles the assembly process. 762 763 .seealso: [](ch_matrices), `Mat`, `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES` 764 @*/ 765 PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode) 766 { 767 PetscErrorCode (*f)(Mat, const PetscScalar[], InsertMode) = NULL; 768 PetscBool oldFlg; 769 770 PetscFunctionBegin; 771 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 772 PetscValidType(A, 1); 773 MatCheckPreallocated(A, 1); 774 PetscValidLogicalCollectiveEnum(A, imode, 3); 775 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetValuesCOO_C", &f)); 776 PetscCall(PetscLogEventBegin(MAT_SetVCOO, A, 0, 0, 0)); 777 if (f) { 778 PetscCall((*f)(A, coo_v, imode)); // all known COO implementations do not use MatStash. They do their own off-proc communication 779 PetscCall(MatGetOption(A, MAT_NO_OFF_PROC_ENTRIES, &oldFlg)); 780 PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); // set A->nooffprocentries to avoid costly MatStash scatter in MatAssembly 781 } else { 782 PetscCall(MatSetValuesCOO_Basic(A, coo_v, imode)); // fall back to MatSetValues, which might use MatStash 783 } 784 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 785 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 786 if (f) PetscCall(MatSetOption(A, MAT_NO_OFF_PROC_ENTRIES, oldFlg)); 787 PetscCall(PetscLogEventEnd(MAT_SetVCOO, A, 0, 0, 0)); 788 PetscFunctionReturn(PETSC_SUCCESS); 789 } 790 791 /*@ 792 MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects 793 794 Input Parameters: 795 + A - the matrix 796 - flg - flag indicating whether the boundtocpu flag should be propagated 797 798 Level: developer 799 800 Notes: 801 If the value of flg is set to true, the following will occur 802 + `MatCreateSubMatrices()` and `MatCreateRedundantMatrix()` - bind created matrices to CPU if the input matrix is bound to the CPU. 803 - `MatCreateVecs()` - bind created vectors to CPU if the input matrix is bound to the CPU. 804 805 The bindingpropagates flag itself is also propagated by the above routines. 806 807 Developer Notes: 808 If the fine-scale `DMDA` has the `-dm_bind_below` option set to true, then `DMCreateInterpolationScale()` calls `MatSetBindingPropagates()` 809 on the restriction/interpolation operator to set the bindingpropagates flag to true. 810 811 .seealso: [](ch_matrices), `Mat`, `VecSetBindingPropagates()`, `MatGetBindingPropagates()` 812 @*/ 813 PetscErrorCode MatSetBindingPropagates(Mat A, PetscBool flg) 814 { 815 PetscFunctionBegin; 816 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 817 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 818 A->bindingpropagates = flg; 819 #endif 820 PetscFunctionReturn(PETSC_SUCCESS); 821 } 822 823 /*@ 824 MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects 825 826 Input Parameter: 827 . A - the matrix 828 829 Output Parameter: 830 . flg - flag indicating whether the boundtocpu flag will be propagated 831 832 Level: developer 833 834 .seealso: [](ch_matrices), `Mat`, `MatSetBindingPropagates()` 835 @*/ 836 PetscErrorCode MatGetBindingPropagates(Mat A, PetscBool *flg) 837 { 838 PetscFunctionBegin; 839 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 840 PetscAssertPointer(flg, 2); 841 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 842 *flg = A->bindingpropagates; 843 #else 844 *flg = PETSC_FALSE; 845 #endif 846 PetscFunctionReturn(PETSC_SUCCESS); 847 } 848