1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 27807a1faSBarry Smith 346533700Sstefano_zampini PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat,PetscInt rbs, PetscInt cbs) 446533700Sstefano_zampini { 546533700Sstefano_zampini PetscFunctionBegin; 65c577a9aSstefano_zampini if (!mat->preallocated) PetscFunctionReturn(0); 7aed4548fSBarry 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); 8aed4548fSBarry 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); 946533700Sstefano_zampini PetscFunctionReturn(0); 1046533700Sstefano_zampini } 1146533700Sstefano_zampini 12db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y,PetscScalar a) 137d68702bSBarry Smith { 147d68702bSBarry Smith PetscInt i,start,end; 157d68702bSBarry Smith PetscScalar alpha = a; 167d68702bSBarry Smith PetscBool prevoption; 177d68702bSBarry Smith 187d68702bSBarry Smith PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(MatGetOption(Y,MAT_NO_OFF_PROC_ENTRIES,&prevoption)); 209566063dSJacob Faibussowitsch PetscCall(MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 219566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(Y,&start,&end)); 227d68702bSBarry Smith for (i=start; i<end; i++) { 23ab6153dcSStefano Zampini if (i < Y->cmap->N) { 249566063dSJacob Faibussowitsch PetscCall(MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES)); 257d68702bSBarry Smith } 26ab6153dcSStefano Zampini } 279566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY)); 289566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY)); 299566063dSJacob Faibussowitsch PetscCall(MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption)); 307d68702bSBarry Smith PetscFunctionReturn(0); 317d68702bSBarry Smith } 327d68702bSBarry Smith 3305869f15SSatish Balay /*@ 3469dd0797SLois Curfman McInnes MatCreate - Creates a matrix where the type is determined 357e5f4302SBarry Smith from either a call to MatSetType() or from the options database 367e5f4302SBarry Smith with a call to MatSetFromOptions(). The default matrix type is 3769b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ() 387e5f4302SBarry Smith if you do not set a type in the options database. If you never 397e5f4302SBarry Smith call MatSetType() or MatSetFromOptions() it will generate an 40f8ab6608SSatish Balay error when you try to use the matrix. 4183e1b59cSLois Curfman McInnes 42d083f849SBarry Smith Collective 43cb13003dSBarry Smith 44f69a0ea3SMatthew Knepley Input Parameter: 45f69a0ea3SMatthew Knepley . comm - MPI communicator 467807a1faSBarry Smith 477807a1faSBarry Smith Output Parameter: 48dc401e71SLois Curfman McInnes . A - the matrix 49e0b365e2SLois Curfman McInnes 50273d9f13SBarry Smith Options Database Keys: 51273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 5269b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 53273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 5469b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 55273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 5669b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 57e0b365e2SLois Curfman McInnes 5883e1b59cSLois Curfman McInnes Even More Options Database Keys: 5983e1b59cSLois Curfman McInnes See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 6083e1b59cSLois Curfman McInnes for additional format-specific options. 61e0b365e2SLois Curfman McInnes 62273d9f13SBarry Smith Level: beginner 63273d9f13SBarry Smith 64db781477SPatrick Sanan .seealso: `MatCreateSeqAIJ()`, `MatCreateAIJ()`, 65db781477SPatrick Sanan `MatCreateSeqDense()`, `MatCreateDense()`, 66db781477SPatrick Sanan `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, 67db781477SPatrick Sanan `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`, 68db781477SPatrick Sanan `MatConvert()` 69273d9f13SBarry Smith @*/ 707087cfbeSBarry Smith PetscErrorCode MatCreate(MPI_Comm comm,Mat *A) 71273d9f13SBarry Smith { 72273d9f13SBarry Smith Mat B; 73273d9f13SBarry Smith 74273d9f13SBarry Smith PetscFunctionBegin; 75f69a0ea3SMatthew Knepley PetscValidPointer(A,2); 7697f1f81fSBarry Smith 770298fd71SBarry Smith *A = NULL; 789566063dSJacob Faibussowitsch PetscCall(MatInitializePackage()); 798ba1e511SMatthew Knepley 809566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView)); 819566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm,&B->rmap)); 829566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm,&B->cmap)); 839566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(VECSTANDARD,&B->defaultvectype)); 8426fbe8dcSKarl Rupp 8594342113SStefano Zampini B->congruentlayouts = PETSC_DECIDE; 86273d9f13SBarry Smith B->preallocated = PETSC_FALSE; 876f3d89d0SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 886f3d89d0SStefano Zampini B->boundtocpu = PETSC_TRUE; 896f3d89d0SStefano Zampini #endif 90273d9f13SBarry Smith *A = B; 91273d9f13SBarry Smith PetscFunctionReturn(0); 92273d9f13SBarry Smith } 93273d9f13SBarry Smith 94422a814eSBarry Smith /*@ 9584d44b13SHong Zhang MatSetErrorIfFailure - Causes Mat to generate an error, for example a zero pivot, is detected. 96422a814eSBarry Smith 97422a814eSBarry Smith Logically Collective on Mat 98422a814eSBarry Smith 99422a814eSBarry Smith Input Parameters: 100422a814eSBarry Smith + mat - matrix obtained from MatCreate() 101422a814eSBarry Smith - flg - PETSC_TRUE indicates you want the error generated 102422a814eSBarry Smith 103422a814eSBarry Smith Level: advanced 104422a814eSBarry Smith 105db781477SPatrick Sanan .seealso: `PCSetErrorIfFailure()` 106422a814eSBarry Smith @*/ 10784d44b13SHong Zhang PetscErrorCode MatSetErrorIfFailure(Mat mat,PetscBool flg) 108422a814eSBarry Smith { 109422a814eSBarry Smith PetscFunctionBegin; 110422a814eSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 111422a814eSBarry Smith PetscValidLogicalCollectiveBool(mat,flg,2); 11284d44b13SHong Zhang mat->erroriffailure = flg; 113422a814eSBarry Smith PetscFunctionReturn(0); 114422a814eSBarry Smith } 115422a814eSBarry Smith 116f69a0ea3SMatthew Knepley /*@ 117f69a0ea3SMatthew Knepley MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 118f69a0ea3SMatthew Knepley 119f69a0ea3SMatthew Knepley Collective on Mat 120f69a0ea3SMatthew Knepley 121f69a0ea3SMatthew Knepley Input Parameters: 122f69a0ea3SMatthew Knepley + A - the matrix 123f69a0ea3SMatthew Knepley . m - number of local rows (or PETSC_DECIDE) 124f69a0ea3SMatthew Knepley . n - number of local columns (or PETSC_DECIDE) 125f69a0ea3SMatthew Knepley . M - number of global rows (or PETSC_DETERMINE) 126f69a0ea3SMatthew Knepley - N - number of global columns (or PETSC_DETERMINE) 127f69a0ea3SMatthew Knepley 128f69a0ea3SMatthew Knepley Notes: 129f69a0ea3SMatthew Knepley m (n) and M (N) cannot be both PETSC_DECIDE 130f69a0ea3SMatthew Knepley If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang. 131f69a0ea3SMatthew Knepley 132f69a0ea3SMatthew Knepley If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 133f69a0ea3SMatthew Knepley user must ensure that they are chosen to be compatible with the 134f69a0ea3SMatthew Knepley vectors. To do this, one first considers the matrix-vector product 135f69a0ea3SMatthew Knepley 'y = A x'. The 'm' that is used in the above routine must match the 136f69a0ea3SMatthew Knepley local size used in the vector creation routine VecCreateMPI() for 'y'. 137f69a0ea3SMatthew Knepley Likewise, the 'n' used must match that used as the local size in 138f69a0ea3SMatthew Knepley VecCreateMPI() for 'x'. 139f69a0ea3SMatthew Knepley 140f73d5cc4SBarry Smith You cannot change the sizes once they have been set. 141f73d5cc4SBarry Smith 142f73d5cc4SBarry Smith The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called. 143f73d5cc4SBarry Smith 144f69a0ea3SMatthew Knepley Level: beginner 145f69a0ea3SMatthew Knepley 146db781477SPatrick Sanan .seealso: `MatGetSize()`, `PetscSplitOwnership()` 147f69a0ea3SMatthew Knepley @*/ 1487087cfbeSBarry Smith PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) 149f69a0ea3SMatthew Knepley { 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); 156aed4548fSBarry Smith 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,A->rmap->n,A->rmap->N); 157aed4548fSBarry Smith 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,A->cmap->n,A->cmap->N); 158d0f46423SBarry Smith A->rmap->n = m; 159d0f46423SBarry Smith A->cmap->n = n; 16059cb773eSBarry Smith A->rmap->N = M > -1 ? M : A->rmap->N; 16159cb773eSBarry Smith A->cmap->N = N > -1 ? N : A->cmap->N; 162f69a0ea3SMatthew Knepley PetscFunctionReturn(0); 163f69a0ea3SMatthew Knepley } 164f69a0ea3SMatthew Knepley 16505869f15SSatish Balay /*@ 166273d9f13SBarry Smith MatSetFromOptions - Creates a matrix where the type is determined 167273d9f13SBarry Smith from the options database. Generates a parallel MPI matrix if the 168273d9f13SBarry Smith communicator has more than one processor. The default matrix type is 16969b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if 1707e5f4302SBarry Smith you do not select a type in the options database. 171273d9f13SBarry Smith 172273d9f13SBarry Smith Collective on Mat 173273d9f13SBarry Smith 174273d9f13SBarry Smith Input Parameter: 175273d9f13SBarry Smith . A - the matrix 176273d9f13SBarry Smith 177273d9f13SBarry Smith Options Database Keys: 178273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 17969b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 180273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 18169b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 182273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 18369b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 184273d9f13SBarry Smith 185273d9f13SBarry Smith Even More Options Database Keys: 186273d9f13SBarry Smith See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 187273d9f13SBarry Smith for additional format-specific options. 188bd9ce289SLois Curfman McInnes 1891d69843bSLois Curfman McInnes Level: beginner 1901d69843bSLois Curfman McInnes 191db781477SPatrick Sanan .seealso: `MatCreateSeqAIJ(()`, `MatCreateAIJ()`, 192db781477SPatrick Sanan `MatCreateSeqDense()`, `MatCreateDense()`, 193db781477SPatrick Sanan `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, 194db781477SPatrick Sanan `MatCreateSeqSBAIJ()`, `MatCreateSBAIJ()`, 195db781477SPatrick Sanan `MatConvert()` 1967807a1faSBarry Smith @*/ 1977087cfbeSBarry Smith PetscErrorCode MatSetFromOptions(Mat B) 1987807a1faSBarry Smith { 199f3be49caSLisandro Dalcin const char *deft = MATAIJ; 200f3be49caSLisandro Dalcin char type[256]; 20169df5c0cSJed Brown PetscBool flg,set; 20216e04d98SRichard Tran Mills PetscInt bind_below = 0; 203dbb450caSBarry Smith 2043a40ed3dSBarry Smith PetscFunctionBegin; 2050700a824SBarry Smith PetscValidHeaderSpecific(B,MAT_CLASSID,1); 206f3be49caSLisandro Dalcin 207d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)B); 208535b19f3SBarry Smith 209535b19f3SBarry Smith if (B->rmap->bs < 0) { 210535b19f3SBarry Smith PetscInt newbs = -1; 2119566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg)); 212535b19f3SBarry Smith if (flg) { 2139566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(B->rmap,newbs)); 2149566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(B->cmap,newbs)); 215535b19f3SBarry Smith } 216535b19f3SBarry Smith } 217535b19f3SBarry Smith 2189566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg)); 219273d9f13SBarry Smith if (flg) { 2209566063dSJacob Faibussowitsch PetscCall(MatSetType(B,type)); 221f3be49caSLisandro Dalcin } else if (!((PetscObject)B)->type_name) { 2229566063dSJacob Faibussowitsch PetscCall(MatSetType(B,deft)); 223273d9f13SBarry Smith } 224f3be49caSLisandro Dalcin 2259566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly)); 2269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL)); 2279566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL)); 2289566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_error_if_failure","Generate an error if an error occurs when factoring the matrix","MatSetErrorIfFailure",B->erroriffailure,&B->erroriffailure,NULL)); 229840d65ccSBarry Smith 230f3be49caSLisandro Dalcin if (B->ops->setfromoptions) { 2319566063dSJacob Faibussowitsch PetscCall((*B->ops->setfromoptions)(PetscOptionsObject,B)); 232593c1905SSatish Balay } 233f3be49caSLisandro Dalcin 23469df5c0cSJed Brown flg = PETSC_FALSE; 2359566063dSJacob 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)); 2369566063dSJacob Faibussowitsch if (set) PetscCall(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg)); 23769df5c0cSJed Brown flg = PETSC_FALSE; 2389566063dSJacob 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)); 2399566063dSJacob Faibussowitsch if (set) PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg)); 240478db826SMatthew G. Knepley flg = PETSC_FALSE; 2419566063dSJacob 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)); 2429566063dSJacob Faibussowitsch if (set) PetscCall(MatSetOption(B,MAT_IGNORE_ZERO_ENTRIES,flg)); 24369df5c0cSJed Brown 2441a2c6b5cSJunchao Zhang flg = PETSC_FALSE; 2459566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_form_explicit_transpose","Hint to form an explicit transpose for operations like MatMultTranspose","MatSetOption",flg,&flg,&set)); 2469566063dSJacob Faibussowitsch if (set) PetscCall(MatSetOption(B,MAT_FORM_EXPLICIT_TRANSPOSE,flg)); 2471a2c6b5cSJunchao Zhang 24816e04d98SRichard Tran Mills /* Bind to CPU if below a user-specified size threshold. 24916e04d98SRichard Tran Mills * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types, 25016e04d98SRichard Tran Mills * and putting it here makes is more maintainable than duplicating this for all. */ 2519566063dSJacob 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)); 25216e04d98SRichard Tran Mills if (flg && B->rmap->n < bind_below) { 2539566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(B,PETSC_TRUE)); 25416e04d98SRichard Tran Mills } 25516e04d98SRichard Tran Mills 2565d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B)); 258d0609cedSBarry Smith PetscOptionsEnd(); 2593a40ed3dSBarry Smith PetscFunctionReturn(0); 2607807a1faSBarry Smith } 2617807a1faSBarry Smith 262987010e7SBarry Smith /*@C 263e8bd9bafSStefano Zampini MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 26463562e91SJed Brown 26563562e91SJed Brown Collective on Mat 26663562e91SJed Brown 2674165533cSJose E. Roman Input Parameters: 26863562e91SJed Brown + A - matrix being preallocated 26963562e91SJed Brown . bs - block size 27041319c1dSStefano Zampini . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 27141319c1dSStefano Zampini . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 27241319c1dSStefano Zampini . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 27341319c1dSStefano Zampini - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 27463562e91SJed Brown 27563562e91SJed Brown Level: beginner 27663562e91SJed Brown 277db781477SPatrick Sanan .seealso: `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, 278db781477SPatrick Sanan `PetscSplitOwnership()` 27963562e91SJed Brown @*/ 2804cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[]) 28163562e91SJed Brown { 28241319c1dSStefano Zampini PetscInt cbs; 28363562e91SJed Brown void (*aij)(void); 284e8bd9bafSStefano Zampini void (*is)(void); 285990279feSStefano Zampini void (*hyp)(void) = NULL; 28663562e91SJed Brown 28763562e91SJed Brown PetscFunctionBegin; 28841319c1dSStefano Zampini if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 2899566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A,bs)); 29041319c1dSStefano Zampini } 2919566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 2929566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 2939566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A,&bs,&cbs)); 29441319c1dSStefano Zampini /* these routines assumes bs == cbs, this should be checked somehow */ 2959566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A,bs,0,dnnz)); 2969566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz)); 2979566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu)); 2989566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu)); 29963562e91SJed Brown /* 300e8bd9bafSStefano 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 30163562e91SJed Brown good before going on with it. 30263562e91SJed Brown */ 3039566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij)); 3049566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is)); 305990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 3069566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp)); 307990279feSStefano Zampini #endif 308990279feSStefano Zampini if (!aij && !is && !hyp) { 3099566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij)); 31063562e91SJed Brown } 311990279feSStefano Zampini if (aij || is || hyp) { 31241319c1dSStefano Zampini if (bs == cbs && bs == 1) { 3139566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,0,dnnz)); 3149566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz)); 3159566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(A,0,dnnz,0,onnz)); 316990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 3179566063dSJacob Faibussowitsch PetscCall(MatHYPRESetPreallocation(A,0,dnnz,0,onnz)); 318990279feSStefano Zampini #endif 3193e5f4774SJed Brown } else { /* Convert block-row precallocation to scalar-row */ 32063562e91SJed Brown PetscInt i,m,*sdnnz,*sonnz; 3219566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 3229566063dSJacob Faibussowitsch PetscCall(PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz)); 323dec54756SJed Brown for (i=0; i<m; i++) { 32441319c1dSStefano Zampini if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs; 32541319c1dSStefano Zampini if (onnz) sonnz[i] = onnz[i/bs] * cbs; 32663562e91SJed Brown } 3279566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL)); 3289566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL)); 3299566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL)); 330990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 3319566063dSJacob Faibussowitsch PetscCall(MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL)); 332990279feSStefano Zampini #endif 3339566063dSJacob Faibussowitsch PetscCall(PetscFree2(sdnnz,sonnz)); 33463562e91SJed Brown } 33563562e91SJed Brown } 33663562e91SJed Brown PetscFunctionReturn(0); 33763562e91SJed Brown } 33863562e91SJed Brown 339273d9f13SBarry Smith /* 340eb6b5d47SBarry Smith Merges some information from Cs header to A; the C object is then destroyed 341d0f46423SBarry Smith 342d0f46423SBarry Smith This is somewhat different from MatHeaderReplace() it would be nice to merge the code 343273d9f13SBarry Smith */ 34428be2f97SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat *C) 345273d9f13SBarry Smith { 346d44834fbSBarry Smith PetscInt refct; 34773107ff1SLisandro Dalcin PetscOps Abops; 34873107ff1SLisandro Dalcin struct _MatOps Aops; 3494768301cSVaclav Hapla char *mtype,*mname,*mprefix; 3504222ddf1SHong Zhang Mat_Product *product; 35133e6eea4SJose E. Roman Mat_Redundant *redundant; 352d4a972cbSStefano Zampini PetscObjectState state; 353273d9f13SBarry Smith 354273d9f13SBarry Smith PetscFunctionBegin; 3551dc04de0SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3561dc04de0SStefano Zampini PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 3571dc04de0SStefano Zampini if (A == *C) PetscFunctionReturn(0); 3581dc04de0SStefano Zampini PetscCheckSameComm(A,1,*C,2); 359273d9f13SBarry Smith /* save the parts of A we need */ 36073107ff1SLisandro Dalcin Abops = ((PetscObject)A)->bops[0]; 36173107ff1SLisandro Dalcin Aops = A->ops[0]; 3627adad957SLisandro Dalcin refct = ((PetscObject)A)->refct; 3635c9eb25fSBarry Smith mtype = ((PetscObject)A)->type_name; 3645c9eb25fSBarry Smith mname = ((PetscObject)A)->name; 365d4a972cbSStefano Zampini state = ((PetscObject)A)->state; 3664768301cSVaclav Hapla mprefix = ((PetscObject)A)->prefix; 3674222ddf1SHong Zhang product = A->product; 36833e6eea4SJose E. Roman redundant = A->redundant; 36930735b05SKris Buschelman 3705c9eb25fSBarry Smith /* zero these so the destroy below does not free them */ 371f4259b30SLisandro Dalcin ((PetscObject)A)->type_name = NULL; 372f4259b30SLisandro Dalcin ((PetscObject)A)->name = NULL; 3735c9eb25fSBarry Smith 3747c99f97cSSatish Balay /* free all the interior data structures from mat */ 3759566063dSJacob Faibussowitsch PetscCall((*A->ops->destroy)(A)); 3767c99f97cSSatish Balay 3779566063dSJacob Faibussowitsch PetscCall(PetscFree(A->defaultvectype)); 3789566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&A->rmap)); 3799566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&A->cmap)); 3809566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&((PetscObject)A)->qlist)); 3819566063dSJacob Faibussowitsch PetscCall(PetscObjectListDestroy(&((PetscObject)A)->olist)); 3829566063dSJacob Faibussowitsch PetscCall(PetscComposedQuantitiesDestroy((PetscObject)A)); 383273d9f13SBarry Smith 384273d9f13SBarry Smith /* copy C over to A */ 385*26cc229bSBarry Smith PetscCall(PetscFree(A->factorprefix)); 3869566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A,*C,sizeof(struct _p_Mat))); 387273d9f13SBarry Smith 388273d9f13SBarry Smith /* return the parts of A we saved */ 38973107ff1SLisandro Dalcin ((PetscObject)A)->bops[0] = Abops; 39073107ff1SLisandro Dalcin A->ops[0] = Aops; 3917adad957SLisandro Dalcin ((PetscObject)A)->refct = refct; 3927adad957SLisandro Dalcin ((PetscObject)A)->type_name = mtype; 3937adad957SLisandro Dalcin ((PetscObject)A)->name = mname; 3944768301cSVaclav Hapla ((PetscObject)A)->prefix = mprefix; 395d4a972cbSStefano Zampini ((PetscObject)A)->state = state + 1; 3964222ddf1SHong Zhang A->product = product; 39733e6eea4SJose E. Roman A->redundant = redundant; 398273d9f13SBarry Smith 3995c9eb25fSBarry Smith /* since these two are copied into A we do not want them destroyed in C */ 400f4259b30SLisandro Dalcin ((PetscObject)*C)->qlist = NULL; 401f4259b30SLisandro Dalcin ((PetscObject)*C)->olist = NULL; 40226fbe8dcSKarl Rupp 4039566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(C)); 404273d9f13SBarry Smith PetscFunctionReturn(0); 405273d9f13SBarry Smith } 4068ab5b326SKris Buschelman /* 407eb6b5d47SBarry Smith Replace A's header with that of C; the C object is then destroyed 408d0f46423SBarry Smith 409eb6b5d47SBarry Smith This is essentially code moved from MatDestroy() 410eb6b5d47SBarry Smith 411eb6b5d47SBarry Smith This is somewhat different from MatHeaderMerge() it would be nice to merge the code 412b30237c6SBarry Smith 413b30237c6SBarry Smith Used in DM hence is declared PETSC_EXTERN 4148ab5b326SKris Buschelman */ 41528be2f97SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C) 4168ab5b326SKris Buschelman { 41727b31e29SJed Brown PetscInt refct; 418fefd9316SJose E. Roman PetscObjectState state; 41928be2f97SBarry Smith struct _p_Mat buffer; 42081fa06acSBarry Smith MatStencilInfo stencil; 4218ab5b326SKris Buschelman 4228ab5b326SKris Buschelman PetscFunctionBegin; 42327b31e29SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 42428be2f97SBarry Smith PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 42528be2f97SBarry Smith if (A == *C) PetscFunctionReturn(0); 42628be2f97SBarry Smith PetscCheckSameComm(A,1,*C,2); 427aed4548fSBarry 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); 4286d7c1e57SBarry Smith 42928be2f97SBarry Smith /* swap C and A */ 43027b31e29SJed Brown refct = ((PetscObject)A)->refct; 431fefd9316SJose E. Roman state = ((PetscObject)A)->state; 43281fa06acSBarry Smith stencil = A->stencil; 4339566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(&buffer,A,sizeof(struct _p_Mat))); 4349566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(A,*C,sizeof(struct _p_Mat))); 4359566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat))); 43627b31e29SJed Brown ((PetscObject)A)->refct = refct; 437fefd9316SJose E. Roman ((PetscObject)A)->state = state + 1; 43881fa06acSBarry Smith A->stencil = stencil; 43926fbe8dcSKarl Rupp 440c32d4117SBarry Smith ((PetscObject)*C)->refct = 1; 4419566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL)); 4429566063dSJacob Faibussowitsch PetscCall(MatDestroy(C)); 4438ab5b326SKris Buschelman PetscFunctionReturn(0); 4448ab5b326SKris Buschelman } 445e7e92044SBarry Smith 446e7e92044SBarry Smith /*@ 447b470e4b4SRichard Tran Mills MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 448e7e92044SBarry Smith 4492216c58aSStefano Zampini Logically collective on Mat 4502216c58aSStefano Zampini 451e7e92044SBarry Smith Input Parameters: 452e7e92044SBarry Smith + A - the matrix 453b470e4b4SRichard Tran Mills - flg - bind to the CPU if value of PETSC_TRUE 454e7e92044SBarry Smith 45590ea27d8SSatish Balay Level: intermediate 4562216c58aSStefano Zampini 457db781477SPatrick Sanan .seealso: `MatBoundToCPU()` 458e7e92044SBarry Smith @*/ 459b470e4b4SRichard Tran Mills PetscErrorCode MatBindToCPU(Mat A,PetscBool flg) 460e7e92044SBarry Smith { 4617d871021SStefano Zampini PetscFunctionBegin; 4622ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 4632ffa8ee7SStefano Zampini PetscValidLogicalCollectiveBool(A,flg,2); 4642216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE) 465b470e4b4SRichard Tran Mills if (A->boundtocpu == flg) PetscFunctionReturn(0); 466b470e4b4SRichard Tran Mills A->boundtocpu = flg; 46792f9df4aSRichard Tran Mills if (A->ops->bindtocpu) { 4689566063dSJacob Faibussowitsch PetscCall((*A->ops->bindtocpu)(A,flg)); 469e7e92044SBarry Smith } 4702216c58aSStefano Zampini #endif 471e7e92044SBarry Smith PetscFunctionReturn(0); 4722216c58aSStefano Zampini } 4732216c58aSStefano Zampini 4742216c58aSStefano Zampini /*@ 4752216c58aSStefano Zampini MatBoundToCPU - query if a matrix is bound to the CPU 4762216c58aSStefano Zampini 4772216c58aSStefano Zampini Input Parameter: 4782216c58aSStefano Zampini . A - the matrix 4792216c58aSStefano Zampini 4802216c58aSStefano Zampini Output Parameter: 4812216c58aSStefano Zampini . flg - the logical flag 4822216c58aSStefano Zampini 4832216c58aSStefano Zampini Level: intermediate 4842216c58aSStefano Zampini 485db781477SPatrick Sanan .seealso: `MatBindToCPU()` 4862216c58aSStefano Zampini @*/ 4872216c58aSStefano Zampini PetscErrorCode MatBoundToCPU(Mat A,PetscBool *flg) 4882216c58aSStefano Zampini { 4892ffa8ee7SStefano Zampini PetscFunctionBegin; 4902ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 491dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg,2); 4922216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE) 4932216c58aSStefano Zampini *flg = A->boundtocpu; 4942216c58aSStefano Zampini #else 4952216c58aSStefano Zampini *flg = PETSC_TRUE; 4967d871021SStefano Zampini #endif 4972216c58aSStefano Zampini PetscFunctionReturn(0); 498e7e92044SBarry Smith } 4997e8381f9SStefano Zampini 5007e8381f9SStefano Zampini PetscErrorCode MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode) 5017e8381f9SStefano Zampini { 5027e8381f9SStefano Zampini IS is_coo_i,is_coo_j; 5037e8381f9SStefano Zampini const PetscInt *coo_i,*coo_j; 5047e8381f9SStefano Zampini PetscInt n,n_i,n_j; 5057e8381f9SStefano Zampini PetscScalar zero = 0.; 5067e8381f9SStefano Zampini 5077e8381f9SStefano Zampini PetscFunctionBegin; 5089566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A,"__PETSc_coo_i",(PetscObject*)&is_coo_i)); 5099566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)A,"__PETSc_coo_j",(PetscObject*)&is_coo_j)); 51028b400f6SJacob Faibussowitsch PetscCheck(is_coo_i,PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_i IS"); 51128b400f6SJacob Faibussowitsch PetscCheck(is_coo_j,PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_j IS"); 5129566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is_coo_i,&n_i)); 5139566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is_coo_j,&n_j)); 51408401ef6SPierre Jolivet PetscCheck(n_i == n_j,PETSC_COMM_SELF,PETSC_ERR_COR,"Wrong local size %" PetscInt_FMT " != %" PetscInt_FMT,n_i,n_j); 5159566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is_coo_i,&coo_i)); 5169566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is_coo_j,&coo_j)); 517e61fc153SStefano Zampini if (imode != ADD_VALUES) { 5189566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(A)); 519e61fc153SStefano Zampini } 5207e8381f9SStefano Zampini for (n = 0; n < n_i; n++) { 5219566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,coo_i[n],coo_j[n],coo_v ? coo_v[n] : zero,ADD_VALUES)); 5227e8381f9SStefano Zampini } 5239566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is_coo_i,&coo_i)); 5249566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is_coo_j,&coo_j)); 5257e8381f9SStefano Zampini PetscFunctionReturn(0); 5267e8381f9SStefano Zampini } 5277e8381f9SStefano Zampini 52882a78a4eSJed Brown PetscErrorCode MatSetPreallocationCOO_Basic(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 5297e8381f9SStefano Zampini { 5307e8381f9SStefano Zampini Mat preallocator; 5317e8381f9SStefano Zampini IS is_coo_i,is_coo_j; 5327e8381f9SStefano Zampini PetscScalar zero = 0.0; 5337e8381f9SStefano Zampini 5347e8381f9SStefano Zampini PetscFunctionBegin; 5359566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 5369566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 5379566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&preallocator)); 5389566063dSJacob Faibussowitsch PetscCall(MatSetType(preallocator,MATPREALLOCATOR)); 5399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(preallocator,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 5409566063dSJacob Faibussowitsch PetscCall(MatSetLayouts(preallocator,A->rmap,A->cmap)); 5419566063dSJacob Faibussowitsch PetscCall(MatSetUp(preallocator)); 54282a78a4eSJed Brown for (PetscCount n = 0; n < ncoo; n++) { 5439566063dSJacob Faibussowitsch PetscCall(MatSetValue(preallocator,coo_i[n],coo_j[n],zero,INSERT_VALUES)); 5447e8381f9SStefano Zampini } 5459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY)); 5469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY)); 5479566063dSJacob Faibussowitsch PetscCall(MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A)); 5489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&preallocator)); 5492c71b3e2SJacob 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); 5509566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_i,PETSC_COPY_VALUES,&is_coo_i)); 5519566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_j,PETSC_COPY_VALUES,&is_coo_j)); 5529566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A,"__PETSc_coo_i",(PetscObject)is_coo_i)); 5539566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)A,"__PETSc_coo_j",(PetscObject)is_coo_j)); 5549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_coo_i)); 5559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is_coo_j)); 5567e8381f9SStefano Zampini PetscFunctionReturn(0); 5577e8381f9SStefano Zampini } 5587e8381f9SStefano Zampini 55956856777SBarry Smith /*@C 560c3dd2894SJed Brown MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries with global indices 5617e8381f9SStefano Zampini 5627e8381f9SStefano Zampini Collective on Mat 5637e8381f9SStefano Zampini 5644165533cSJose E. Roman Input Parameters: 5657e8381f9SStefano Zampini + A - matrix being preallocated 56642550becSJunchao Zhang . ncoo - number of entries 5677e8381f9SStefano Zampini . coo_i - row indices 5687e8381f9SStefano Zampini - coo_j - column indices 5697e8381f9SStefano Zampini 5707e8381f9SStefano Zampini Level: beginner 5717e8381f9SStefano Zampini 572394ed5ebSJunchao Zhang Notes: 573394ed5ebSJunchao Zhang Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed 574394ed5ebSJunchao Zhang but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries 575394ed5ebSJunchao Zhang are allowed and will be properly added or inserted to the matrix, unless the matrix option MAT_IGNORE_OFF_PROC_ENTRIES 576394ed5ebSJunchao 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. 5777e8381f9SStefano Zampini 578735d7f90SBarry Smith The arrays coo_i and coo_j may be freed immediately after calling this function. 579735d7f90SBarry Smith 580db781477SPatrick Sanan .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOOLocal()`, `DMSetMatrixPreallocateSkip()` 5817e8381f9SStefano Zampini @*/ 58282a78a4eSJed Brown PetscErrorCode MatSetPreallocationCOO(Mat A,PetscCount ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 5837e8381f9SStefano Zampini { 58482a78a4eSJed Brown PetscErrorCode (*f)(Mat,PetscCount,const PetscInt[],const PetscInt[]) = NULL; 5857e8381f9SStefano Zampini 5867e8381f9SStefano Zampini PetscFunctionBegin; 5877e8381f9SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 5887e8381f9SStefano Zampini PetscValidType(A,1); 5897e8381f9SStefano Zampini if (ncoo) PetscValidIntPointer(coo_i,3); 5907e8381f9SStefano Zampini if (ncoo) PetscValidIntPointer(coo_j,4); 5919566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 5929566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 5939566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOO_C",&f)); 594cbc6b225SStefano Zampini 5959566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_PreallCOO,A,0,0,0)); 5967e8381f9SStefano Zampini if (f) { 5979566063dSJacob Faibussowitsch PetscCall((*f)(A,ncoo,coo_i,coo_j)); 5987e8381f9SStefano Zampini } else { /* allow fallback, very slow */ 5999566063dSJacob Faibussowitsch PetscCall(MatSetPreallocationCOO_Basic(A,ncoo,coo_i,coo_j)); 6007e8381f9SStefano Zampini } 6019566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_PreallCOO,A,0,0,0)); 6026834774dSStefano Zampini A->preallocated = PETSC_TRUE; 603cbc6b225SStefano Zampini A->nonzerostate++; 6047e8381f9SStefano Zampini PetscFunctionReturn(0); 6057e8381f9SStefano Zampini } 6067e8381f9SStefano Zampini 60756856777SBarry Smith /*@C 608c3dd2894SJed Brown MatSetPreallocationCOOLocal - set preallocation for matrices using a coordinate format of the entries with local indices 609c3dd2894SJed Brown 610c3dd2894SJed Brown Collective on Mat 611c3dd2894SJed Brown 612c3dd2894SJed Brown Input Parameters: 613c3dd2894SJed Brown + A - matrix being preallocated 614c3dd2894SJed Brown . ncoo - number of entries 615c3dd2894SJed Brown . coo_i - row indices (local numbering; may be modified) 616c3dd2894SJed Brown - coo_j - column indices (local numbering; may be modified) 617c3dd2894SJed Brown 618c3dd2894SJed Brown Level: beginner 619c3dd2894SJed Brown 620c3dd2894SJed Brown Notes: 621c3dd2894SJed Brown The local indices are translated using the local to global mapping, thus MatSetLocalToGlobalMapping() must have been 622c3dd2894SJed Brown called prior to this function. 623c3dd2894SJed Brown 624c3dd2894SJed Brown The indices coo_i and coo_j may be modified within this function. They might be translated to corresponding global 625735d7f90SBarry Smith indices, but the caller should not rely on them having any specific value after this function returns. The arrays 626735d7f90SBarry Smith can be freed or reused immediately after this function returns. 627c3dd2894SJed Brown 628394ed5ebSJunchao Zhang Entries can be repeated, see MatSetValuesCOO(). Entries with negative row or column indices are allowed 629394ed5ebSJunchao Zhang but will be ignored. The corresponding entries in MatSetValuesCOO() will be ignored too. Remote entries 630394ed5ebSJunchao Zhang are allowed and will be properly added or inserted to the matrix. 631c3dd2894SJed Brown 632db781477SPatrick Sanan .seealso: `MatSetValuesCOO()`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatSeqBAIJSetPreallocation()`, `MatMPIBAIJSetPreallocation()`, `MatSeqSBAIJSetPreallocation()`, `MatMPISBAIJSetPreallocation()`, `MatSetPreallocationCOO()`, `DMSetMatrixPreallocateSkip()` 633c3dd2894SJed Brown @*/ 63482a78a4eSJed Brown PetscErrorCode MatSetPreallocationCOOLocal(Mat A,PetscCount ncoo,PetscInt coo_i[],PetscInt coo_j[]) 635c3dd2894SJed Brown { 6366834774dSStefano Zampini PetscErrorCode (*f)(Mat,PetscCount,PetscInt[],PetscInt[]) = NULL; 637c3dd2894SJed Brown 638c3dd2894SJed Brown PetscFunctionBegin; 639c3dd2894SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 640c3dd2894SJed Brown PetscValidType(A,1); 641c3dd2894SJed Brown if (ncoo) PetscValidIntPointer(coo_i,3); 642c3dd2894SJed Brown if (ncoo) PetscValidIntPointer(coo_j,4); 6436834774dSStefano 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); 6449566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->rmap)); 6459566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(A->cmap)); 646cbc6b225SStefano Zampini 6479566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOOLocal_C",&f)); 6486834774dSStefano Zampini if (f) { 6499566063dSJacob Faibussowitsch PetscCall((*f)(A,ncoo,coo_i,coo_j)); 650cbc6b225SStefano Zampini A->nonzerostate++; 6516834774dSStefano Zampini } else { 652cbc6b225SStefano Zampini ISLocalToGlobalMapping ltog_row,ltog_col; 6539566063dSJacob Faibussowitsch PetscCall(MatGetLocalToGlobalMapping(A,<og_row,<og_col)); 6549566063dSJacob Faibussowitsch if (ltog_row) PetscCall(ISLocalToGlobalMappingApply(ltog_row,ncoo,coo_i,coo_i)); 6559566063dSJacob Faibussowitsch if (ltog_col) PetscCall(ISLocalToGlobalMappingApply(ltog_col,ncoo,coo_j,coo_j)); 6569566063dSJacob Faibussowitsch PetscCall(MatSetPreallocationCOO(A,ncoo,coo_i,coo_j)); 6576834774dSStefano Zampini } 6586834774dSStefano Zampini A->preallocated = PETSC_TRUE; 659c3dd2894SJed Brown PetscFunctionReturn(0); 660c3dd2894SJed Brown } 661c3dd2894SJed Brown 662c3dd2894SJed Brown /*@ 663bfcc3627SStefano Zampini MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO() 6647e8381f9SStefano Zampini 6657e8381f9SStefano Zampini Collective on Mat 6667e8381f9SStefano Zampini 6674165533cSJose E. Roman Input Parameters: 6687e8381f9SStefano Zampini + A - matrix being preallocated 669bfcc3627SStefano Zampini . coo_v - the matrix values (can be NULL) 6707e8381f9SStefano Zampini - imode - the insert mode 6717e8381f9SStefano Zampini 6727e8381f9SStefano Zampini Level: beginner 6737e8381f9SStefano Zampini 674c3dd2894SJed Brown Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO() or MatSetPreallocationCOOLocal(). 675735d7f90SBarry Smith When repeated entries are specified in the COO indices the coo_v values are first properly summed, regardless of the value of imode. 676bfcc3627SStefano Zampini The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES). 677735d7f90SBarry Smith MatAssemblyBegin() and MatAssemblyEnd() do not need to be called after this routine. It automatically handles the assembly process. 6787e8381f9SStefano Zampini 679db781477SPatrick Sanan .seealso: `MatSetPreallocationCOO()`, `MatSetPreallocationCOOLocal()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES` 6807e8381f9SStefano Zampini @*/ 6817e8381f9SStefano Zampini PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode) 6827e8381f9SStefano Zampini { 6837e8381f9SStefano Zampini PetscErrorCode (*f)(Mat,const PetscScalar[],InsertMode) = NULL; 6847e8381f9SStefano Zampini 6857e8381f9SStefano Zampini PetscFunctionBegin; 6867e8381f9SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 6877e8381f9SStefano Zampini PetscValidType(A,1); 6887e8381f9SStefano Zampini MatCheckPreallocated(A,1); 689bfcc3627SStefano Zampini PetscValidLogicalCollectiveEnum(A,imode,3); 6909566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)A,"MatSetValuesCOO_C",&f)); 6919566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_SetVCOO,A,0,0,0)); 6927e8381f9SStefano Zampini if (f) { 6939566063dSJacob Faibussowitsch PetscCall((*f)(A,coo_v,imode)); 6947e8381f9SStefano Zampini } else { /* allow fallback */ 6959566063dSJacob Faibussowitsch PetscCall(MatSetValuesCOO_Basic(A,coo_v,imode)); 6967e8381f9SStefano Zampini } 6979566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_SetVCOO,A,0,0,0)); 6989566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 6999566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 7007e8381f9SStefano Zampini PetscFunctionReturn(0); 7017e8381f9SStefano Zampini } 70265a9ecf2SRichard Tran Mills 70365a9ecf2SRichard Tran Mills /*@ 70465a9ecf2SRichard 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 70565a9ecf2SRichard Tran Mills 70665a9ecf2SRichard Tran Mills Input Parameters: 70765a9ecf2SRichard Tran Mills + A - the matrix 70865a9ecf2SRichard Tran Mills - flg - flag indicating whether the boundtocpu flag should be propagated 70965a9ecf2SRichard Tran Mills 71065a9ecf2SRichard Tran Mills Level: developer 71165a9ecf2SRichard Tran Mills 71265a9ecf2SRichard Tran Mills Notes: 71365a9ecf2SRichard Tran Mills If the value of flg is set to true, the following will occur: 71465a9ecf2SRichard Tran Mills 71565a9ecf2SRichard Tran Mills MatCreateSubMatrices() and MatCreateRedundantMatrix() will bind created matrices to CPU if the input matrix is bound to the CPU. 71665a9ecf2SRichard Tran Mills MatCreateVecs() will bind created vectors to CPU if the input matrix is bound to the CPU. 71765a9ecf2SRichard Tran Mills The bindingpropagates flag itself is also propagated by the above routines. 71865a9ecf2SRichard Tran Mills 71965a9ecf2SRichard Tran Mills Developer Notes: 72065a9ecf2SRichard Tran Mills If the fine-scale DMDA has the -dm_bind_below option set to true, then DMCreateInterpolationScale() calls MatSetBindingPropagates() 72165a9ecf2SRichard Tran Mills on the restriction/interpolation operator to set the bindingpropagates flag to true. 72265a9ecf2SRichard Tran Mills 723db781477SPatrick Sanan .seealso: `VecSetBindingPropagates()`, `MatGetBindingPropagates()` 72465a9ecf2SRichard Tran Mills @*/ 72565a9ecf2SRichard Tran Mills PetscErrorCode MatSetBindingPropagates(Mat A,PetscBool flg) 72665a9ecf2SRichard Tran Mills { 72765a9ecf2SRichard Tran Mills PetscFunctionBegin; 72865a9ecf2SRichard Tran Mills PetscValidHeaderSpecific(A,MAT_CLASSID,1); 72965a9ecf2SRichard Tran Mills #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 73065a9ecf2SRichard Tran Mills A->bindingpropagates = flg; 73165a9ecf2SRichard Tran Mills #endif 73265a9ecf2SRichard Tran Mills PetscFunctionReturn(0); 73365a9ecf2SRichard Tran Mills } 734e9c74fd6SRichard Tran Mills 735e9c74fd6SRichard Tran Mills /*@ 736e9c74fd6SRichard 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 737e9c74fd6SRichard Tran Mills 738e9c74fd6SRichard Tran Mills Input Parameter: 739e9c74fd6SRichard Tran Mills . A - the matrix 740e9c74fd6SRichard Tran Mills 741e9c74fd6SRichard Tran Mills Output Parameter: 742e9c74fd6SRichard Tran Mills . flg - flag indicating whether the boundtocpu flag will be propagated 743e9c74fd6SRichard Tran Mills 744e9c74fd6SRichard Tran Mills Level: developer 745e9c74fd6SRichard Tran Mills 746db781477SPatrick Sanan .seealso: `MatSetBindingPropagates()` 747e9c74fd6SRichard Tran Mills @*/ 748e9c74fd6SRichard Tran Mills PetscErrorCode MatGetBindingPropagates(Mat A,PetscBool *flg) 749e9c74fd6SRichard Tran Mills { 750e9c74fd6SRichard Tran Mills PetscFunctionBegin; 751e9c74fd6SRichard Tran Mills PetscValidHeaderSpecific(A,MAT_CLASSID,1); 752e9c74fd6SRichard Tran Mills PetscValidBoolPointer(flg,2); 753e9c74fd6SRichard Tran Mills #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 754e9c74fd6SRichard Tran Mills *flg = A->bindingpropagates; 755e9c74fd6SRichard Tran Mills #else 756e9c74fd6SRichard Tran Mills *flg = PETSC_FALSE; 757e9c74fd6SRichard Tran Mills #endif 758e9c74fd6SRichard Tran Mills PetscFunctionReturn(0); 759e9c74fd6SRichard Tran Mills } 760