1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 27807a1faSBarry Smith 39371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat, PetscInt rbs, PetscInt cbs) { 446533700Sstefano_zampini PetscFunctionBegin; 55c577a9aSstefano_zampini if (!mat->preallocated) PetscFunctionReturn(0); 6aed4548fSBarry Smith 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); 7aed4548fSBarry Smith 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); 846533700Sstefano_zampini PetscFunctionReturn(0); 946533700Sstefano_zampini } 1046533700Sstefano_zampini 119371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y, PetscScalar a) { 127d68702bSBarry Smith PetscInt i, start, end; 137d68702bSBarry Smith PetscScalar alpha = a; 147d68702bSBarry Smith PetscBool prevoption; 157d68702bSBarry Smith 167d68702bSBarry Smith PetscFunctionBegin; 179566063dSJacob Faibussowitsch PetscCall(MatGetOption(Y, MAT_NO_OFF_PROC_ENTRIES, &prevoption)); 189566063dSJacob Faibussowitsch PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE)); 199566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Y, &start, &end)); 207d68702bSBarry Smith for (i = start; i < end; i++) { 21*48a46eb9SPierre Jolivet if (i < Y->cmap->N) PetscCall(MatSetValues(Y, 1, &i, 1, &i, &alpha, ADD_VALUES)); 22ab6153dcSStefano Zampini } 239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY)); 249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY)); 259566063dSJacob Faibussowitsch PetscCall(MatSetOption(Y, MAT_NO_OFF_PROC_ENTRIES, prevoption)); 267d68702bSBarry Smith PetscFunctionReturn(0); 277d68702bSBarry Smith } 287d68702bSBarry Smith 2905869f15SSatish Balay /*@ 3069dd0797SLois Curfman McInnes MatCreate - Creates a matrix where the type is determined 317e5f4302SBarry Smith from either a call to MatSetType() or from the options database 327e5f4302SBarry Smith with a call to MatSetFromOptions(). The default matrix type is 3369b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ() 347e5f4302SBarry Smith if you do not set a type in the options database. If you never 357e5f4302SBarry Smith call MatSetType() or MatSetFromOptions() it will generate an 36f8ab6608SSatish Balay error when you try to use the matrix. 3783e1b59cSLois Curfman McInnes 38d083f849SBarry Smith Collective 39cb13003dSBarry Smith 40f69a0ea3SMatthew Knepley Input Parameter: 41f69a0ea3SMatthew Knepley . comm - MPI communicator 427807a1faSBarry Smith 437807a1faSBarry Smith Output Parameter: 44dc401e71SLois Curfman McInnes . A - the matrix 45e0b365e2SLois Curfman McInnes 46273d9f13SBarry Smith Options Database Keys: 47273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 4869b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 49273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 5069b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 51273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 5269b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 53e0b365e2SLois Curfman McInnes 5483e1b59cSLois Curfman McInnes Even More Options Database Keys: 5583e1b59cSLois Curfman McInnes See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 5683e1b59cSLois Curfman McInnes for additional format-specific options. 57e0b365e2SLois Curfman McInnes 58273d9f13SBarry Smith Level: beginner 59273d9f13SBarry Smith 60db781477SPatrick Sanan .seealso: `MatCreateSeqAIJ()`, `MatCreateAIJ()`, 61db781477SPatrick Sanan `MatCreateSeqDense()`, `MatCreateDense()`, 62db781477SPatrick Sanan `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, 63db781477SPatrick Sanan `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`, 64db781477SPatrick Sanan `MatConvert()` 65273d9f13SBarry Smith @*/ 669371c9d4SSatish Balay PetscErrorCode MatCreate(MPI_Comm comm, Mat *A) { 67273d9f13SBarry Smith Mat B; 68273d9f13SBarry Smith 69273d9f13SBarry Smith PetscFunctionBegin; 70f69a0ea3SMatthew Knepley PetscValidPointer(A, 2); 7197f1f81fSBarry Smith 720298fd71SBarry Smith *A = NULL; 739566063dSJacob Faibussowitsch PetscCall(MatInitializePackage()); 748ba1e511SMatthew Knepley 759566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(B, MAT_CLASSID, "Mat", "Matrix", "Mat", comm, MatDestroy, MatView)); 769566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &B->rmap)); 779566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &B->cmap)); 789566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(VECSTANDARD, &B->defaultvectype)); 7926fbe8dcSKarl Rupp 80b94d7dedSBarry Smith B->symmetric = PETSC_BOOL3_UNKNOWN; 81b94d7dedSBarry Smith B->hermitian = PETSC_BOOL3_UNKNOWN; 82b94d7dedSBarry Smith B->structurally_symmetric = PETSC_BOOL3_UNKNOWN; 83b94d7dedSBarry Smith B->spd = PETSC_BOOL3_UNKNOWN; 84b94d7dedSBarry Smith B->symmetry_eternal = PETSC_FALSE; 85b94d7dedSBarry Smith B->structural_symmetry_eternal = PETSC_FALSE; 86b94d7dedSBarry Smith 8794342113SStefano Zampini B->congruentlayouts = PETSC_DECIDE; 88273d9f13SBarry Smith B->preallocated = PETSC_FALSE; 896f3d89d0SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 906f3d89d0SStefano Zampini B->boundtocpu = PETSC_TRUE; 916f3d89d0SStefano Zampini #endif 92273d9f13SBarry Smith *A = B; 93273d9f13SBarry Smith PetscFunctionReturn(0); 94273d9f13SBarry Smith } 95273d9f13SBarry Smith 96422a814eSBarry Smith /*@ 9784d44b13SHong Zhang MatSetErrorIfFailure - Causes Mat to generate an error, for example a zero pivot, is detected. 98422a814eSBarry Smith 99422a814eSBarry Smith Logically Collective on Mat 100422a814eSBarry Smith 101422a814eSBarry Smith Input Parameters: 102422a814eSBarry Smith + mat - matrix obtained from MatCreate() 103422a814eSBarry Smith - flg - PETSC_TRUE indicates you want the error generated 104422a814eSBarry Smith 105422a814eSBarry Smith Level: advanced 106422a814eSBarry Smith 107db781477SPatrick Sanan .seealso: `PCSetErrorIfFailure()` 108422a814eSBarry Smith @*/ 1099371c9d4SSatish Balay PetscErrorCode MatSetErrorIfFailure(Mat mat, PetscBool flg) { 110422a814eSBarry Smith PetscFunctionBegin; 111422a814eSBarry Smith PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 112422a814eSBarry Smith PetscValidLogicalCollectiveBool(mat, flg, 2); 11384d44b13SHong Zhang mat->erroriffailure = flg; 114422a814eSBarry Smith PetscFunctionReturn(0); 115422a814eSBarry Smith } 116422a814eSBarry Smith 117f69a0ea3SMatthew Knepley /*@ 118f69a0ea3SMatthew Knepley MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 119f69a0ea3SMatthew Knepley 120f69a0ea3SMatthew Knepley Collective on Mat 121f69a0ea3SMatthew Knepley 122f69a0ea3SMatthew Knepley Input Parameters: 123f69a0ea3SMatthew Knepley + A - the matrix 124f69a0ea3SMatthew Knepley . m - number of local rows (or PETSC_DECIDE) 125f69a0ea3SMatthew Knepley . n - number of local columns (or PETSC_DECIDE) 126f69a0ea3SMatthew Knepley . M - number of global rows (or PETSC_DETERMINE) 127f69a0ea3SMatthew Knepley - N - number of global columns (or PETSC_DETERMINE) 128f69a0ea3SMatthew Knepley 129f69a0ea3SMatthew Knepley Notes: 130f69a0ea3SMatthew Knepley m (n) and M (N) cannot be both PETSC_DECIDE 131f69a0ea3SMatthew Knepley If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang. 132f69a0ea3SMatthew Knepley 133f69a0ea3SMatthew Knepley If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 134f69a0ea3SMatthew Knepley user must ensure that they are chosen to be compatible with the 135f69a0ea3SMatthew Knepley vectors. To do this, one first considers the matrix-vector product 136f69a0ea3SMatthew Knepley 'y = A x'. The 'm' that is used in the above routine must match the 137f69a0ea3SMatthew Knepley local size used in the vector creation routine VecCreateMPI() for 'y'. 138f69a0ea3SMatthew Knepley Likewise, the 'n' used must match that used as the local size in 139f69a0ea3SMatthew Knepley VecCreateMPI() for 'x'. 140f69a0ea3SMatthew Knepley 141f73d5cc4SBarry Smith You cannot change the sizes once they have been set. 142f73d5cc4SBarry Smith 143f73d5cc4SBarry Smith The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called. 144f73d5cc4SBarry Smith 145f69a0ea3SMatthew Knepley Level: beginner 146f69a0ea3SMatthew Knepley 147db781477SPatrick Sanan .seealso: `MatGetSize()`, `PetscSplitOwnership()` 148f69a0ea3SMatthew Knepley @*/ 1499371c9d4SSatish Balay PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) { 150f69a0ea3SMatthew Knepley PetscFunctionBegin; 1510700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 152a69c7061SStefano Zampini PetscValidLogicalCollectiveInt(A, M, 4); 153a69c7061SStefano Zampini PetscValidLogicalCollectiveInt(A, N, 5); 154aed4548fSBarry Smith 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); 155aed4548fSBarry Smith 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); 1569371c9d4SSatish Balay 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, 1579371c9d4SSatish Balay A->rmap->n, A->rmap->N); 1589371c9d4SSatish Balay 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, 1599371c9d4SSatish Balay A->cmap->n, A->cmap->N); 160d0f46423SBarry Smith A->rmap->n = m; 161d0f46423SBarry Smith A->cmap->n = n; 16259cb773eSBarry Smith A->rmap->N = M > -1 ? M : A->rmap->N; 16359cb773eSBarry Smith A->cmap->N = N > -1 ? N : A->cmap->N; 164f69a0ea3SMatthew Knepley PetscFunctionReturn(0); 165f69a0ea3SMatthew Knepley } 166f69a0ea3SMatthew Knepley 16705869f15SSatish Balay /*@ 168273d9f13SBarry Smith MatSetFromOptions - Creates a matrix where the type is determined 169273d9f13SBarry Smith from the options database. Generates a parallel MPI matrix if the 170273d9f13SBarry Smith communicator has more than one processor. The default matrix type is 17169b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if 1727e5f4302SBarry Smith you do not select a type in the options database. 173273d9f13SBarry Smith 174273d9f13SBarry Smith Collective on Mat 175273d9f13SBarry Smith 176273d9f13SBarry Smith Input Parameter: 177273d9f13SBarry Smith . A - the matrix 178273d9f13SBarry Smith 179273d9f13SBarry Smith Options Database Keys: 180273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 18169b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 182273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 18369b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 184273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 18569b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 186273d9f13SBarry Smith 187273d9f13SBarry Smith Even More Options Database Keys: 188273d9f13SBarry Smith See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 189273d9f13SBarry Smith for additional format-specific options. 190bd9ce289SLois Curfman McInnes 1911d69843bSLois Curfman McInnes Level: beginner 1921d69843bSLois Curfman McInnes 193db781477SPatrick Sanan .seealso: `MatCreateSeqAIJ(()`, `MatCreateAIJ()`, 194db781477SPatrick Sanan `MatCreateSeqDense()`, `MatCreateDense()`, 195db781477SPatrick Sanan `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, 196db781477SPatrick Sanan `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`, 197db781477SPatrick Sanan `MatConvert()` 1987807a1faSBarry Smith @*/ 1999371c9d4SSatish Balay PetscErrorCode MatSetFromOptions(Mat B) { 200f3be49caSLisandro Dalcin const char *deft = MATAIJ; 201f3be49caSLisandro Dalcin char type[256]; 20269df5c0cSJed Brown PetscBool flg, set; 20316e04d98SRichard Tran Mills PetscInt bind_below = 0; 204dbb450caSBarry Smith 2053a40ed3dSBarry Smith PetscFunctionBegin; 2060700a824SBarry Smith PetscValidHeaderSpecific(B, MAT_CLASSID, 1); 207f3be49caSLisandro Dalcin 208d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)B); 209535b19f3SBarry Smith 210535b19f3SBarry Smith if (B->rmap->bs < 0) { 211535b19f3SBarry Smith PetscInt newbs = -1; 2129566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_block_size", "Set the blocksize used to store the matrix", "MatSetBlockSize", newbs, &newbs, &flg)); 213535b19f3SBarry Smith if (flg) { 2149566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(B->rmap, newbs)); 2159566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(B->cmap, newbs)); 216535b19f3SBarry Smith } 217535b19f3SBarry Smith } 218535b19f3SBarry Smith 2199566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-mat_type", "Matrix type", "MatSetType", MatList, deft, type, 256, &flg)); 220273d9f13SBarry Smith if (flg) { 2219566063dSJacob Faibussowitsch PetscCall(MatSetType(B, type)); 222f3be49caSLisandro Dalcin } else if (!((PetscObject)B)->type_name) { 2239566063dSJacob Faibussowitsch PetscCall(MatSetType(B, deft)); 224273d9f13SBarry Smith } 225f3be49caSLisandro Dalcin 2269566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", &B->checksymmetryonassembly)); 2279566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_is_symmetric", "Checks if mat is symmetric on MatAssemblyEnd()", "MatIsSymmetric", B->checksymmetrytol, &B->checksymmetrytol, NULL)); 2289566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_null_space_test", "Checks if provided null space is correct in MatAssemblyEnd()", "MatSetNullSpaceTest", B->checknullspaceonassembly, &B->checknullspaceonassembly, NULL)); 2299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_error_if_failure", "Generate an error if an error occurs when factoring the matrix", "MatSetErrorIfFailure", B->erroriffailure, &B->erroriffailure, NULL)); 230840d65ccSBarry Smith 231dbbe0bcdSBarry Smith PetscTryTypeMethod(B, setfromoptions, PetscOptionsObject); 232f3be49caSLisandro Dalcin 23369df5c0cSJed Brown flg = PETSC_FALSE; 2349566063dSJacob Faibussowitsch 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)); 2359566063dSJacob Faibussowitsch if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, flg)); 23669df5c0cSJed Brown flg = PETSC_FALSE; 2379566063dSJacob Faibussowitsch 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)); 2389566063dSJacob Faibussowitsch if (set) PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, flg)); 239478db826SMatthew G. Knepley flg = PETSC_FALSE; 2409566063dSJacob Faibussowitsch 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)); 2419566063dSJacob Faibussowitsch if (set) PetscCall(MatSetOption(B, MAT_IGNORE_ZERO_ENTRIES, flg)); 24269df5c0cSJed Brown 2431a2c6b5cSJunchao Zhang flg = PETSC_FALSE; 2449566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_form_explicit_transpose", "Hint to form an explicit transpose for operations like MatMultTranspose", "MatSetOption", flg, &flg, &set)); 2459566063dSJacob Faibussowitsch if (set) PetscCall(MatSetOption(B, MAT_FORM_EXPLICIT_TRANSPOSE, flg)); 2461a2c6b5cSJunchao Zhang 24716e04d98SRichard Tran Mills /* Bind to CPU if below a user-specified size threshold. 24816e04d98SRichard Tran Mills * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types, 24916e04d98SRichard Tran Mills * and putting it here makes is more maintainable than duplicating this for all. */ 2509566063dSJacob Faibussowitsch 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)); 251*48a46eb9SPierre Jolivet if (flg && B->rmap->n < bind_below) PetscCall(MatBindToCPU(B, PETSC_TRUE)); 25216e04d98SRichard Tran Mills 2535d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 254dbbe0bcdSBarry Smith PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)B, PetscOptionsObject)); 255d0609cedSBarry Smith PetscOptionsEnd(); 2563a40ed3dSBarry Smith PetscFunctionReturn(0); 2577807a1faSBarry Smith } 2587807a1faSBarry Smith 259987010e7SBarry Smith /*@C 260e8bd9bafSStefano Zampini MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 26163562e91SJed Brown 26263562e91SJed Brown Collective on Mat 26363562e91SJed Brown 2644165533cSJose E. Roman Input Parameters: 26563562e91SJed Brown + A - matrix being preallocated 26663562e91SJed Brown . bs - block size 26741319c1dSStefano Zampini . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 26841319c1dSStefano Zampini . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 26941319c1dSStefano Zampini . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 27041319c1dSStefano Zampini - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 27163562e91SJed Brown 27263562e91SJed Brown Level: beginner 27363562e91SJed Brown 274db781477SPatrick Sanan .seealso: `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, 275db781477SPatrick Sanan `PetscSplitOwnership()` 27663562e91SJed Brown @*/ 2779371c9d4SSatish Balay PetscErrorCode MatXAIJSetPreallocation(Mat A, PetscInt bs, const PetscInt dnnz[], const PetscInt onnz[], const PetscInt dnnzu[], const PetscInt onnzu[]) { 27841319c1dSStefano Zampini PetscInt cbs; 27963562e91SJed Brown void (*aij)(void); 280e8bd9bafSStefano Zampini void (*is)(void); 281990279feSStefano Zampini void (*hyp)(void) = NULL; 28263562e91SJed Brown 28363562e91SJed Brown PetscFunctionBegin; 28441319c1dSStefano Zampini if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 2859566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A, bs)); 28641319c1dSStefano Zampini } 2879566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 2889566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 2899566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A, &bs, &cbs)); 29041319c1dSStefano Zampini /* these routines assumes bs == cbs, this should be checked somehow */ 2919566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A, bs, 0, dnnz)); 2929566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A, bs, 0, dnnz, 0, onnz)); 2939566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(A, bs, 0, dnnzu)); 2949566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A, bs, 0, dnnzu, 0, onnzu)); 29563562e91SJed Brown /* 296e8bd9bafSStefano Zampini 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 29763562e91SJed Brown good before going on with it. 29863562e91SJed Brown */ 2999566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatMPIAIJSetPreallocation_C", &aij)); 3009566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatISSetPreallocation_C", &is)); 301990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 3029566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatHYPRESetPreallocation_C", &hyp)); 303990279feSStefano Zampini #endif 304*48a46eb9SPierre Jolivet if (!aij && !is && !hyp) PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSeqAIJSetPreallocation_C", &aij)); 305990279feSStefano Zampini if (aij || is || hyp) { 30641319c1dSStefano Zampini if (bs == cbs && bs == 1) { 3079566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz)); 3089566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz, 0, onnz)); 3099566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(A, 0, dnnz, 0, onnz)); 310990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 3119566063dSJacob Faibussowitsch PetscCall(MatHYPRESetPreallocation(A, 0, dnnz, 0, onnz)); 312990279feSStefano Zampini #endif 3133e5f4774SJed Brown } else { /* Convert block-row precallocation to scalar-row */ 31463562e91SJed Brown PetscInt i, m, *sdnnz, *sonnz; 3159566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 3169566063dSJacob Faibussowitsch PetscCall(PetscMalloc2((!!dnnz) * m, &sdnnz, (!!onnz) * m, &sonnz)); 317dec54756SJed Brown for (i = 0; i < m; i++) { 31841319c1dSStefano Zampini if (dnnz) sdnnz[i] = dnnz[i / bs] * cbs; 31941319c1dSStefano Zampini if (onnz) sonnz[i] = onnz[i / bs] * cbs; 32063562e91SJed Brown } 3219566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL)); 3229566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 3239566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 324990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 3259566063dSJacob Faibussowitsch PetscCall(MatHYPRESetPreallocation(A, 0, dnnz ? sdnnz : NULL, 0, onnz ? sonnz : NULL)); 326990279feSStefano Zampini #endif 3279566063dSJacob Faibussowitsch PetscCall(PetscFree2(sdnnz, sonnz)); 32863562e91SJed Brown } 32963562e91SJed Brown } 33063562e91SJed Brown PetscFunctionReturn(0); 33163562e91SJed Brown } 33263562e91SJed Brown 333273d9f13SBarry Smith /* 334eb6b5d47SBarry Smith Merges some information from Cs header to A; the C object is then destroyed 335d0f46423SBarry Smith 336d0f46423SBarry Smith This is somewhat different from MatHeaderReplace() it would be nice to merge the code 337273d9f13SBarry Smith */ 3389371c9d4SSatish Balay PetscErrorCode MatHeaderMerge(Mat A, Mat *C) { 339d44834fbSBarry Smith PetscInt refct; 34073107ff1SLisandro Dalcin PetscOps Abops; 34173107ff1SLisandro Dalcin struct _MatOps Aops; 3424768301cSVaclav Hapla char *mtype, *mname, *mprefix; 3434222ddf1SHong Zhang Mat_Product *product; 34433e6eea4SJose E. Roman Mat_Redundant *redundant; 345d4a972cbSStefano Zampini PetscObjectState state; 346273d9f13SBarry Smith 347273d9f13SBarry Smith PetscFunctionBegin; 3481dc04de0SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 3491dc04de0SStefano Zampini PetscValidHeaderSpecific(*C, MAT_CLASSID, 2); 3501dc04de0SStefano Zampini if (A == *C) PetscFunctionReturn(0); 3511dc04de0SStefano Zampini PetscCheckSameComm(A, 1, *C, 2); 352273d9f13SBarry Smith /* save the parts of A we need */ 35373107ff1SLisandro Dalcin Abops = ((PetscObject)A)->bops[0]; 35473107ff1SLisandro Dalcin Aops = A->ops[0]; 3557adad957SLisandro Dalcin refct = ((PetscObject)A)->refct; 3565c9eb25fSBarry Smith mtype = ((PetscObject)A)->type_name; 3575c9eb25fSBarry Smith mname = ((PetscObject)A)->name; 358d4a972cbSStefano Zampini state = ((PetscObject)A)->state; 3594768301cSVaclav Hapla mprefix = ((PetscObject)A)->prefix; 3604222ddf1SHong Zhang product = A->product; 36133e6eea4SJose E. Roman redundant = A->redundant; 36230735b05SKris Buschelman 3635c9eb25fSBarry Smith /* zero these so the destroy below does not free them */ 364f4259b30SLisandro Dalcin ((PetscObject)A)->type_name = NULL; 365f4259b30SLisandro Dalcin ((PetscObject)A)->name = NULL; 3665c9eb25fSBarry Smith 367dbbe0bcdSBarry Smith /* 368dbbe0bcdSBarry Smith free all the interior data structures from mat 369dbbe0bcdSBarry Smith cannot use PetscUseTypeMethod(A,destroy); because compiler 370dbbe0bcdSBarry Smith thinks it may print NULL type_name and name 371dbbe0bcdSBarry Smith */ 372dbbe0bcdSBarry Smith PetscTryTypeMethod(A, destroy); 3737c99f97cSSatish Balay 3749566063dSJacob Faibussowitsch PetscCall(PetscFree(A->defaultvectype)); 3759566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&A->rmap)); 3769566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&A->cmap)); 3779566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&((PetscObject)A)->qlist)); 3789566063dSJacob Faibussowitsch PetscCall(PetscObjectListDestroy(&((PetscObject)A)->olist)); 3799566063dSJacob Faibussowitsch PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A)); 380273d9f13SBarry Smith 381273d9f13SBarry Smith /* copy C over to A */ 38226cc229bSBarry Smith PetscCall(PetscFree(A->factorprefix)); 3839566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat))); 384273d9f13SBarry Smith 385273d9f13SBarry Smith /* return the parts of A we saved */ 38673107ff1SLisandro Dalcin ((PetscObject)A)->bops[0] = Abops; 38773107ff1SLisandro Dalcin A->ops[0] = Aops; 3887adad957SLisandro Dalcin ((PetscObject)A)->refct = refct; 3897adad957SLisandro Dalcin ((PetscObject)A)->type_name = mtype; 3907adad957SLisandro Dalcin ((PetscObject)A)->name = mname; 3914768301cSVaclav Hapla ((PetscObject)A)->prefix = mprefix; 392d4a972cbSStefano Zampini ((PetscObject)A)->state = state + 1; 3934222ddf1SHong Zhang A->product = product; 39433e6eea4SJose E. Roman A->redundant = redundant; 395273d9f13SBarry Smith 3965c9eb25fSBarry Smith /* since these two are copied into A we do not want them destroyed in C */ 397f4259b30SLisandro Dalcin ((PetscObject)*C)->qlist = NULL; 398f4259b30SLisandro Dalcin ((PetscObject)*C)->olist = NULL; 39926fbe8dcSKarl Rupp 4009566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(C)); 401273d9f13SBarry Smith PetscFunctionReturn(0); 402273d9f13SBarry Smith } 4038ab5b326SKris Buschelman /* 404eb6b5d47SBarry Smith Replace A's header with that of C; the C object is then destroyed 405d0f46423SBarry Smith 406eb6b5d47SBarry Smith This is essentially code moved from MatDestroy() 407eb6b5d47SBarry Smith 408eb6b5d47SBarry Smith This is somewhat different from MatHeaderMerge() it would be nice to merge the code 409b30237c6SBarry Smith 410b30237c6SBarry Smith Used in DM hence is declared PETSC_EXTERN 4118ab5b326SKris Buschelman */ 4129371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A, Mat *C) { 41327b31e29SJed Brown PetscInt refct; 414fefd9316SJose E. Roman PetscObjectState state; 41528be2f97SBarry Smith struct _p_Mat buffer; 41681fa06acSBarry Smith MatStencilInfo stencil; 4178ab5b326SKris Buschelman 4188ab5b326SKris Buschelman PetscFunctionBegin; 41927b31e29SJed Brown PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 42028be2f97SBarry Smith PetscValidHeaderSpecific(*C, MAT_CLASSID, 2); 42128be2f97SBarry Smith if (A == *C) PetscFunctionReturn(0); 42228be2f97SBarry Smith PetscCheckSameComm(A, 1, *C, 2); 423aed4548fSBarry Smith 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); 4246d7c1e57SBarry Smith 42528be2f97SBarry Smith /* swap C and A */ 42627b31e29SJed Brown refct = ((PetscObject)A)->refct; 427fefd9316SJose E. Roman state = ((PetscObject)A)->state; 42881fa06acSBarry Smith stencil = A->stencil; 4299566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(&buffer, A, sizeof(struct _p_Mat))); 4309566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A, *C, sizeof(struct _p_Mat))); 4319566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(*C, &buffer, sizeof(struct _p_Mat))); 43227b31e29SJed Brown ((PetscObject)A)->refct = refct; 433fefd9316SJose E. Roman ((PetscObject)A)->state = state + 1; 43481fa06acSBarry Smith A->stencil = stencil; 43526fbe8dcSKarl Rupp 436c32d4117SBarry Smith ((PetscObject)*C)->refct = 1; 4379566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*C, MATOP_DESTROY, (void (*)(void))NULL)); 4389566063dSJacob Faibussowitsch PetscCall(MatDestroy(C)); 4398ab5b326SKris Buschelman PetscFunctionReturn(0); 4408ab5b326SKris Buschelman } 441e7e92044SBarry Smith 442e7e92044SBarry Smith /*@ 443b470e4b4SRichard Tran Mills MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 444e7e92044SBarry Smith 4452216c58aSStefano Zampini Logically collective on Mat 4462216c58aSStefano Zampini 447e7e92044SBarry Smith Input Parameters: 448e7e92044SBarry Smith + A - the matrix 449b470e4b4SRichard Tran Mills - flg - bind to the CPU if value of PETSC_TRUE 450e7e92044SBarry Smith 45190ea27d8SSatish Balay Level: intermediate 4522216c58aSStefano Zampini 453db781477SPatrick Sanan .seealso: `MatBoundToCPU()` 454e7e92044SBarry Smith @*/ 4559371c9d4SSatish Balay PetscErrorCode MatBindToCPU(Mat A, PetscBool flg) { 4567d871021SStefano Zampini PetscFunctionBegin; 4572ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 4582ffa8ee7SStefano Zampini PetscValidLogicalCollectiveBool(A, flg, 2); 4592216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE) 460b470e4b4SRichard Tran Mills if (A->boundtocpu == flg) PetscFunctionReturn(0); 461b470e4b4SRichard Tran Mills A->boundtocpu = flg; 462dbbe0bcdSBarry Smith PetscTryTypeMethod(A, bindtocpu, flg); 4632216c58aSStefano Zampini #endif 464e7e92044SBarry Smith PetscFunctionReturn(0); 4652216c58aSStefano Zampini } 4662216c58aSStefano Zampini 4672216c58aSStefano Zampini /*@ 4682216c58aSStefano Zampini MatBoundToCPU - query if a matrix is bound to the CPU 4692216c58aSStefano Zampini 4702216c58aSStefano Zampini Input Parameter: 4712216c58aSStefano Zampini . A - the matrix 4722216c58aSStefano Zampini 4732216c58aSStefano Zampini Output Parameter: 4742216c58aSStefano Zampini . flg - the logical flag 4752216c58aSStefano Zampini 4762216c58aSStefano Zampini Level: intermediate 4772216c58aSStefano Zampini 478db781477SPatrick Sanan .seealso: `MatBindToCPU()` 4792216c58aSStefano Zampini @*/ 4809371c9d4SSatish Balay PetscErrorCode MatBoundToCPU(Mat A, PetscBool *flg) { 4812ffa8ee7SStefano Zampini PetscFunctionBegin; 4822ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 483dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg, 2); 4842216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE) 4852216c58aSStefano Zampini *flg = A->boundtocpu; 4862216c58aSStefano Zampini #else 4872216c58aSStefano Zampini *flg = PETSC_TRUE; 4887d871021SStefano Zampini #endif 4892216c58aSStefano Zampini PetscFunctionReturn(0); 490e7e92044SBarry Smith } 4917e8381f9SStefano Zampini 4929371c9d4SSatish Balay PetscErrorCode MatSetValuesCOO_Basic(Mat A, const PetscScalar coo_v[], InsertMode imode) { 4937e8381f9SStefano Zampini IS is_coo_i, is_coo_j; 4947e8381f9SStefano Zampini const PetscInt *coo_i, *coo_j; 4957e8381f9SStefano Zampini PetscInt n, n_i, n_j; 4967e8381f9SStefano Zampini PetscScalar zero = 0.; 4977e8381f9SStefano Zampini 4987e8381f9SStefano Zampini PetscFunctionBegin; 4999566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_i", (PetscObject *)&is_coo_i)); 5009566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A, "__PETSc_coo_j", (PetscObject *)&is_coo_j)); 50128b400f6SJacob Faibussowitsch PetscCheck(is_coo_i, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_i IS"); 50228b400f6SJacob Faibussowitsch PetscCheck(is_coo_j, PetscObjectComm((PetscObject)A), PETSC_ERR_COR, "Missing coo_j IS"); 5039566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is_coo_i, &n_i)); 5049566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is_coo_j, &n_j)); 50508401ef6SPierre Jolivet PetscCheck(n_i == n_j, PETSC_COMM_SELF, PETSC_ERR_COR, "Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT, n_i, n_j); 5069566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is_coo_i, &coo_i)); 5079566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is_coo_j, &coo_j)); 508*48a46eb9SPierre Jolivet if (imode != ADD_VALUES) PetscCall(MatZeroEntries(A)); 509*48a46eb9SPierre Jolivet for (n = 0; n < n_i; n++) PetscCall(MatSetValue(A, coo_i[n], coo_j[n], coo_v ? coo_v[n] : zero, ADD_VALUES)); 5109566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is_coo_i, &coo_i)); 5119566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is_coo_j, &coo_j)); 5127e8381f9SStefano Zampini PetscFunctionReturn(0); 5137e8381f9SStefano Zampini } 5147e8381f9SStefano Zampini 5159371c9d4SSatish Balay PetscErrorCode MatSetPreallocationCOO_Basic(Mat A, PetscCount ncoo, const PetscInt coo_i[], const PetscInt coo_j[]) { 5167e8381f9SStefano Zampini Mat preallocator; 5177e8381f9SStefano Zampini IS is_coo_i, is_coo_j; 5187e8381f9SStefano Zampini PetscScalar zero = 0.0; 5197e8381f9SStefano Zampini 5207e8381f9SStefano Zampini PetscFunctionBegin; 5219566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 5229566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 5239566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator)); 5249566063dSJacob Faibussowitsch PetscCall(MatSetType(preallocator, MATPREALLOCATOR)); 5259566063dSJacob Faibussowitsch PetscCall(MatSetSizes(preallocator, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 5269566063dSJacob Faibussowitsch PetscCall(MatSetLayouts(preallocator, A->rmap, A->cmap)); 5279566063dSJacob Faibussowitsch PetscCall(MatSetUp(preallocator)); 528*48a46eb9SPierre Jolivet for (PetscCount n = 0; n < ncoo; n++) PetscCall(MatSetValue(preallocator, coo_i[n], coo_j[n], zero, INSERT_VALUES)); 5299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY)); 5309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY)); 5319566063dSJacob Faibussowitsch PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A)); 5329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&preallocator)); 5332c71b3e2SJacob Faibussowitsch 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); 5349566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo, coo_i, PETSC_COPY_VALUES, &is_coo_i)); 5359566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncoo, coo_j, PETSC_COPY_VALUES, &is_coo_j)); 5369566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_i", (PetscObject)is_coo_i)); 5379566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A, "__PETSc_coo_j", (PetscObject)is_coo_j)); 5389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_coo_i)); 5399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_coo_j)); 5407e8381f9SStefano Zampini PetscFunctionReturn(0); 5417e8381f9SStefano Zampini } 5427e8381f9SStefano Zampini 54356856777SBarry Smith /*@C 544c3dd2894SJed Brown MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices 5457e8381f9SStefano Zampini 5467e8381f9SStefano Zampini Collective on Mat 5477e8381f9SStefano Zampini 5484165533cSJose E. Roman Input Parameters: 5497e8381f9SStefano Zampini + A - matrix being preallocated 55042550becSJunchao Zhang . ncoo - number of entries 5517e8381f9SStefano Zampini . coo_i - row indices 5527e8381f9SStefano Zampini - coo_j - column indices 5537e8381f9SStefano Zampini 5547e8381f9SStefano Zampini Level: beginner 5557e8381f9SStefano Zampini 556394ed5ebSJunchao Zhang Notes: 557e8729f6fSJunchao Zhang The indices coo_i and coo_j may be modified within this function. The caller should not rely on them 558e8729f6fSJunchao Zhang having any specific value after this function returns. The arrays can be freed or reused immediately 559e8729f6fSJunchao Zhang after this function returns. 560e8729f6fSJunchao Zhang 561394ed5ebSJunchao Zhang Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed 562394ed5ebSJunchao Zhang but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries 563394ed5ebSJunchao Zhang are allowed and will be properly added or inserted to the matrix, unless the matrix option MAT_IGNORE_OFF_PROC_ENTRIES 564394ed5ebSJunchao Zhang 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. 5657e8381f9SStefano Zampini 566db781477SPatrick Sanan .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`, `DMSetMatrixPreallocateSkip()` 5677e8381f9SStefano Zampini @*/ 5689371c9d4SSatish Balay PetscErrorCode MatSetPreallocationCOO(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) { 56982a78a4eSJed Brown PetscErrorCode (*f)(Mat, PetscCount, const PetscInt[], const PetscInt[]) = NULL; 5707e8381f9SStefano Zampini 5717e8381f9SStefano Zampini PetscFunctionBegin; 5727e8381f9SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5737e8381f9SStefano Zampini PetscValidType(A, 1); 5747e8381f9SStefano Zampini if (ncoo) PetscValidIntPointer(coo_i, 3); 5757e8381f9SStefano Zampini if (ncoo) PetscValidIntPointer(coo_j, 4); 5769566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 5779566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 5789566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOO_C", &f)); 579cbc6b225SStefano Zampini 5809566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_PreallCOO, A, 0, 0, 0)); 5817e8381f9SStefano Zampini if (f) { 5829566063dSJacob Faibussowitsch PetscCall((*f)(A, ncoo, coo_i, coo_j)); 5837e8381f9SStefano Zampini } else { /* allow fallback, very slow */ 5849566063dSJacob Faibussowitsch PetscCall(MatSetPreallocationCOO_Basic(A, ncoo, coo_i, coo_j)); 5857e8381f9SStefano Zampini } 5869566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_PreallCOO, A, 0, 0, 0)); 5876834774dSStefano Zampini A->preallocated = PETSC_TRUE; 588cbc6b225SStefano Zampini A->nonzerostate++; 5897e8381f9SStefano Zampini PetscFunctionReturn(0); 5907e8381f9SStefano Zampini } 5917e8381f9SStefano Zampini 59256856777SBarry Smith /*@C 593c3dd2894SJed Brown MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices 594c3dd2894SJed Brown 595c3dd2894SJed Brown Collective on Mat 596c3dd2894SJed Brown 597c3dd2894SJed Brown Input Parameters: 598c3dd2894SJed Brown + A - matrix being preallocated 599c3dd2894SJed Brown . ncoo - number of entries 600c3dd2894SJed Brown . coo_i - row indices (local numbering; may be modified) 601c3dd2894SJed Brown - coo_j - column indices (local numbering; may be modified) 602c3dd2894SJed Brown 603c3dd2894SJed Brown Level: beginner 604c3dd2894SJed Brown 605c3dd2894SJed Brown Notes: 606c3dd2894SJed Brown The local indices are translated using the local to global mapping, thus MatSetLocalToGlobalMapping() must have been 607c3dd2894SJed Brown called prior to this function. 608c3dd2894SJed Brown 609c3dd2894SJed Brown The indices coo_i and coo_j may be modified within this function. They might be translated to corresponding global 610735d7f90SBarry Smith indices, but the caller should not rely on them having any specific value after this function returns. The arrays 611735d7f90SBarry Smith can be freed or reused immediately after this function returns. 612c3dd2894SJed Brown 613394ed5ebSJunchao Zhang Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed 614394ed5ebSJunchao Zhang but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries 615394ed5ebSJunchao Zhang are allowed and will be properly added or inserted to the matrix. 616c3dd2894SJed Brown 617db781477SPatrick Sanan .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`, `DMSetMatrixPreallocateSkip()` 618c3dd2894SJed Brown @*/ 6199371c9d4SSatish Balay PetscErrorCode MatSetPreallocationCOOLocal(Mat A, PetscCount ncoo, PetscInt coo_i[], PetscInt coo_j[]) { 6206834774dSStefano Zampini PetscErrorCode (*f)(Mat, PetscCount, PetscInt[], PetscInt[]) = NULL; 621c3dd2894SJed Brown 622c3dd2894SJed Brown PetscFunctionBegin; 623c3dd2894SJed Brown PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 624c3dd2894SJed Brown PetscValidType(A, 1); 625c3dd2894SJed Brown if (ncoo) PetscValidIntPointer(coo_i, 3); 626c3dd2894SJed Brown if (ncoo) PetscValidIntPointer(coo_j, 4); 6276834774dSStefano Zampini 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); 6289566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 6299566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 630cbc6b225SStefano Zampini 6319566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetPreallocationCOOLocal_C", &f)); 6326834774dSStefano Zampini if (f) { 6339566063dSJacob Faibussowitsch PetscCall((*f)(A, ncoo, coo_i, coo_j)); 634cbc6b225SStefano Zampini A->nonzerostate++; 6356834774dSStefano Zampini } else { 636cbc6b225SStefano Zampini ISLocalToGlobalMapping ltog_row, ltog_col; 6379566063dSJacob Faibussowitsch PetscCall(MatGetLocalToGlobalMapping(A, <og_row, <og_col)); 6389566063dSJacob Faibussowitsch if (ltog_row) PetscCall(ISLocalToGlobalMappingApply(ltog_row, ncoo, coo_i, coo_i)); 6399566063dSJacob Faibussowitsch if (ltog_col) PetscCall(ISLocalToGlobalMappingApply(ltog_col, ncoo, coo_j, coo_j)); 6409566063dSJacob Faibussowitsch PetscCall(MatSetPreallocationCOO(A, ncoo, coo_i, coo_j)); 6416834774dSStefano Zampini } 6426834774dSStefano Zampini A->preallocated = PETSC_TRUE; 643c3dd2894SJed Brown PetscFunctionReturn(0); 644c3dd2894SJed Brown } 645c3dd2894SJed Brown 646c3dd2894SJed Brown /*@ 647bfcc3627SStefano Zampini MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO() 6487e8381f9SStefano Zampini 6497e8381f9SStefano Zampini Collective on Mat 6507e8381f9SStefano Zampini 6514165533cSJose E. Roman Input Parameters: 6527e8381f9SStefano Zampini + A - matrix being preallocated 653bfcc3627SStefano Zampini . coo_v - the matrix values (can be NULL) 6547e8381f9SStefano Zampini - imode - the insert mode 6557e8381f9SStefano Zampini 6567e8381f9SStefano Zampini Level: beginner 6577e8381f9SStefano Zampini 658c3dd2894SJed Brown Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO() or MatSetPreallocationCOOLocal(). 659735d7f90SBarry Smith When repeated entries are specified in the COO indices the coo_v values are first properly summed, regardless of the value of imode. 660bfcc3627SStefano Zampini The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES). 661735d7f90SBarry Smith MatAssemblyBegin() and MatAssemblyEnd() do not need to be called after this routine. It automatically handles the assembly process. 6627e8381f9SStefano Zampini 663db781477SPatrick Sanan .seealso: `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES` 6647e8381f9SStefano Zampini @*/ 6659371c9d4SSatish Balay PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode) { 6667e8381f9SStefano Zampini PetscErrorCode (*f)(Mat, const PetscScalar[], InsertMode) = NULL; 6677e8381f9SStefano Zampini 6687e8381f9SStefano Zampini PetscFunctionBegin; 6697e8381f9SStefano Zampini PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 6707e8381f9SStefano Zampini PetscValidType(A, 1); 6717e8381f9SStefano Zampini MatCheckPreallocated(A, 1); 672bfcc3627SStefano Zampini PetscValidLogicalCollectiveEnum(A, imode, 3); 6739566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatSetValuesCOO_C", &f)); 6749566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_SetVCOO, A, 0, 0, 0)); 6757e8381f9SStefano Zampini if (f) { 6769566063dSJacob Faibussowitsch PetscCall((*f)(A, coo_v, imode)); 6777e8381f9SStefano Zampini } else { /* allow fallback */ 6789566063dSJacob Faibussowitsch PetscCall(MatSetValuesCOO_Basic(A, coo_v, imode)); 6797e8381f9SStefano Zampini } 6809566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_SetVCOO, A, 0, 0, 0)); 6819566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 6829566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 6837e8381f9SStefano Zampini PetscFunctionReturn(0); 6847e8381f9SStefano Zampini } 68565a9ecf2SRichard Tran Mills 68665a9ecf2SRichard Tran Mills /*@ 68765a9ecf2SRichard Tran Mills MatSetBindingPropagates - Sets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects 68865a9ecf2SRichard Tran Mills 68965a9ecf2SRichard Tran Mills Input Parameters: 69065a9ecf2SRichard Tran Mills + A - the matrix 69165a9ecf2SRichard Tran Mills - flg - flag indicating whether the boundtocpu flag should be propagated 69265a9ecf2SRichard Tran Mills 69365a9ecf2SRichard Tran Mills Level: developer 69465a9ecf2SRichard Tran Mills 69565a9ecf2SRichard Tran Mills Notes: 69665a9ecf2SRichard Tran Mills If the value of flg is set to true, the following will occur: 69765a9ecf2SRichard Tran Mills 69865a9ecf2SRichard Tran Mills MatCreateSubMatrices() and MatCreateRedundantMatrix() will bind created matrices to CPU if the input matrix is bound to the CPU. 69965a9ecf2SRichard Tran Mills MatCreateVecs() will bind created vectors to CPU if the input matrix is bound to the CPU. 70065a9ecf2SRichard Tran Mills The bindingpropagates flag itself is also propagated by the above routines. 70165a9ecf2SRichard Tran Mills 70265a9ecf2SRichard Tran Mills Developer Notes: 70365a9ecf2SRichard Tran Mills If the fine-scale DMDA has the -dm_bind_below option set to true, then DMCreateInterpolationScale() calls MatSetBindingPropagates() 70465a9ecf2SRichard Tran Mills on the restriction/interpolation operator to set the bindingpropagates flag to true. 70565a9ecf2SRichard Tran Mills 706db781477SPatrick Sanan .seealso: `VecSetBindingPropagates()`, `MatGetBindingPropagates()` 70765a9ecf2SRichard Tran Mills @*/ 7089371c9d4SSatish Balay PetscErrorCode MatSetBindingPropagates(Mat A, PetscBool flg) { 70965a9ecf2SRichard Tran Mills PetscFunctionBegin; 71065a9ecf2SRichard Tran Mills PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 71165a9ecf2SRichard Tran Mills #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 71265a9ecf2SRichard Tran Mills A->bindingpropagates = flg; 71365a9ecf2SRichard Tran Mills #endif 71465a9ecf2SRichard Tran Mills PetscFunctionReturn(0); 71565a9ecf2SRichard Tran Mills } 716e9c74fd6SRichard Tran Mills 717e9c74fd6SRichard Tran Mills /*@ 718e9c74fd6SRichard Tran Mills MatGetBindingPropagates - Gets whether the state of being bound to the CPU for a GPU matrix type propagates to child and some other associated objects 719e9c74fd6SRichard Tran Mills 720e9c74fd6SRichard Tran Mills Input Parameter: 721e9c74fd6SRichard Tran Mills . A - the matrix 722e9c74fd6SRichard Tran Mills 723e9c74fd6SRichard Tran Mills Output Parameter: 724e9c74fd6SRichard Tran Mills . flg - flag indicating whether the boundtocpu flag will be propagated 725e9c74fd6SRichard Tran Mills 726e9c74fd6SRichard Tran Mills Level: developer 727e9c74fd6SRichard Tran Mills 728db781477SPatrick Sanan .seealso: `MatSetBindingPropagates()` 729e9c74fd6SRichard Tran Mills @*/ 7309371c9d4SSatish Balay PetscErrorCode MatGetBindingPropagates(Mat A, PetscBool *flg) { 731e9c74fd6SRichard Tran Mills PetscFunctionBegin; 732e9c74fd6SRichard Tran Mills PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 733e9c74fd6SRichard Tran Mills PetscValidBoolPointer(flg, 2); 734e9c74fd6SRichard Tran Mills #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 735e9c74fd6SRichard Tran Mills *flg = A->bindingpropagates; 736e9c74fd6SRichard Tran Mills #else 737e9c74fd6SRichard Tran Mills *flg = PETSC_FALSE; 738e9c74fd6SRichard Tran Mills #endif 739e9c74fd6SRichard Tran Mills PetscFunctionReturn(0); 740e9c74fd6SRichard Tran Mills } 741