17807a1faSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 37807a1faSBarry Smith 446533700Sstefano_zampini PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat,PetscInt rbs, PetscInt cbs) 546533700Sstefano_zampini { 646533700Sstefano_zampini PetscFunctionBegin; 75c577a9aSstefano_zampini if (!mat->preallocated) PetscFunctionReturn(0); 846533700Sstefano_zampini if (mat->rmap->bs > 0 && mat->rmap->bs != rbs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot change row block size %D to %D\n",mat->rmap->bs,rbs); 946533700Sstefano_zampini if (mat->cmap->bs > 0 && mat->cmap->bs != cbs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot change column block size %D to %D\n",mat->cmap->bs,cbs); 1046533700Sstefano_zampini PetscFunctionReturn(0); 1146533700Sstefano_zampini } 1246533700Sstefano_zampini 13db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y,PetscScalar a) 147d68702bSBarry Smith { 157d68702bSBarry Smith PetscErrorCode ierr; 167d68702bSBarry Smith PetscInt i,start,end; 177d68702bSBarry Smith PetscScalar alpha = a; 187d68702bSBarry Smith PetscBool prevoption; 197d68702bSBarry Smith 207d68702bSBarry Smith PetscFunctionBegin; 217d68702bSBarry Smith ierr = MatGetOption(Y,MAT_NO_OFF_PROC_ENTRIES,&prevoption);CHKERRQ(ierr); 227d68702bSBarry Smith ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 237d68702bSBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 247d68702bSBarry Smith for (i=start; i<end; i++) { 25ab6153dcSStefano Zampini if (i < Y->cmap->N) { 267d68702bSBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);CHKERRQ(ierr); 277d68702bSBarry Smith } 28ab6153dcSStefano Zampini } 297d68702bSBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 307d68702bSBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 317d68702bSBarry Smith ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption);CHKERRQ(ierr); 327d68702bSBarry Smith PetscFunctionReturn(0); 337d68702bSBarry Smith } 347d68702bSBarry Smith 3505869f15SSatish Balay /*@ 3669dd0797SLois Curfman McInnes MatCreate - Creates a matrix where the type is determined 377e5f4302SBarry Smith from either a call to MatSetType() or from the options database 387e5f4302SBarry Smith with a call to MatSetFromOptions(). The default matrix type is 3969b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ() 407e5f4302SBarry Smith if you do not set a type in the options database. If you never 417e5f4302SBarry Smith call MatSetType() or MatSetFromOptions() it will generate an 42f8ab6608SSatish Balay error when you try to use the matrix. 4383e1b59cSLois Curfman McInnes 44d083f849SBarry Smith Collective 45cb13003dSBarry Smith 46f69a0ea3SMatthew Knepley Input Parameter: 47f69a0ea3SMatthew Knepley . comm - MPI communicator 487807a1faSBarry Smith 497807a1faSBarry Smith Output Parameter: 50dc401e71SLois Curfman McInnes . A - the matrix 51e0b365e2SLois Curfman McInnes 52273d9f13SBarry Smith Options Database Keys: 53273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 5469b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 55273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 5669b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 57273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 5869b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 59e0b365e2SLois Curfman McInnes 6083e1b59cSLois Curfman McInnes Even More Options Database Keys: 6183e1b59cSLois Curfman McInnes See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 6283e1b59cSLois Curfman McInnes for additional format-specific options. 63e0b365e2SLois Curfman McInnes 64273d9f13SBarry Smith Level: beginner 65273d9f13SBarry Smith 6669b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ(), MatCreateAIJ(), 6769b1f4b7SBarry Smith MatCreateSeqDense(), MatCreateDense(), 6869b1f4b7SBarry Smith MatCreateSeqBAIJ(), MatCreateBAIJ(), 6969b1f4b7SBarry Smith MatCreateSeqSBAIJ(), MatCreateSBAIJ(), 70273d9f13SBarry Smith MatConvert() 71273d9f13SBarry Smith @*/ 727087cfbeSBarry Smith PetscErrorCode MatCreate(MPI_Comm comm,Mat *A) 73273d9f13SBarry Smith { 74273d9f13SBarry Smith Mat B; 75dfbe8321SBarry Smith PetscErrorCode ierr; 76273d9f13SBarry Smith 77273d9f13SBarry Smith PetscFunctionBegin; 78f69a0ea3SMatthew Knepley PetscValidPointer(A,2); 7997f1f81fSBarry Smith 800298fd71SBarry Smith *A = NULL; 81607a6623SBarry Smith ierr = MatInitializePackage();CHKERRQ(ierr); 828ba1e511SMatthew Knepley 8373107ff1SLisandro Dalcin ierr = PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);CHKERRQ(ierr); 8426283091SBarry Smith ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr); 8526283091SBarry Smith ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr); 8634136279SStefano Zampini ierr = PetscStrallocpy(VECSTANDARD,&B->defaultvectype);CHKERRQ(ierr); 8726fbe8dcSKarl Rupp 8894342113SStefano Zampini B->congruentlayouts = PETSC_DECIDE; 89273d9f13SBarry Smith B->preallocated = PETSC_FALSE; 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 105422a814eSBarry Smith .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 146f69a0ea3SMatthew Knepley .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); 154cb9de970SPierre Jolivet if (M > 0 && m > M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",m,M); 155cb9de970SPierre Jolivet if (N > 0 && n > N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",n,N); 15659cb773eSBarry Smith if ((A->rmap->n >= 0 && A->rmap->N >= 0) && (A->rmap->n != m || (M > 0 && A->rmap->N != M))) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset row sizes to %D local %D global after previously setting them to %D local %D global",m,M,A->rmap->n,A->rmap->N); 15759cb773eSBarry Smith if ((A->cmap->n >= 0 && A->cmap->N >= 0) && (A->cmap->n != n || (N > 0 && A->cmap->N != N))) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset column sizes to %D local %D global after previously setting them to %D local %D global",n,N,A->cmap->n,A->cmap->N); 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 19169b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateAIJ(), 19269b1f4b7SBarry Smith MatCreateSeqDense(), MatCreateDense(), 19369b1f4b7SBarry Smith MatCreateSeqBAIJ(), MatCreateBAIJ(), 19469b1f4b7SBarry Smith MatCreateSeqSBAIJ(), MatCreateSBAIJ(), 195273d9f13SBarry Smith MatConvert() 1967807a1faSBarry Smith @*/ 1977087cfbeSBarry Smith PetscErrorCode MatSetFromOptions(Mat B) 1987807a1faSBarry Smith { 199dfbe8321SBarry Smith PetscErrorCode ierr; 200f3be49caSLisandro Dalcin const char *deft = MATAIJ; 201f3be49caSLisandro Dalcin char type[256]; 20269df5c0cSJed Brown PetscBool flg,set; 203dbb450caSBarry Smith 2043a40ed3dSBarry Smith PetscFunctionBegin; 2050700a824SBarry Smith PetscValidHeaderSpecific(B,MAT_CLASSID,1); 206f3be49caSLisandro Dalcin 2073194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr); 208535b19f3SBarry Smith 209535b19f3SBarry Smith if (B->rmap->bs < 0) { 210535b19f3SBarry Smith PetscInt newbs = -1; 211535b19f3SBarry Smith ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr); 212535b19f3SBarry Smith if (flg) { 213535b19f3SBarry Smith ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr); 214535b19f3SBarry Smith ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr); 215535b19f3SBarry Smith } 216535b19f3SBarry Smith } 217535b19f3SBarry Smith 218a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr); 219273d9f13SBarry Smith if (flg) { 220f3be49caSLisandro Dalcin ierr = MatSetType(B,type);CHKERRQ(ierr); 221f3be49caSLisandro Dalcin } else if (!((PetscObject)B)->type_name) { 222f3be49caSLisandro Dalcin ierr = MatSetType(B,deft);CHKERRQ(ierr); 223273d9f13SBarry Smith } 224f3be49caSLisandro Dalcin 225840d65ccSBarry Smith ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr); 22694ae4db5SBarry Smith ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL);CHKERRQ(ierr); 22794ae4db5SBarry Smith ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL);CHKERRQ(ierr); 2285307bdd6SStefano Zampini ierr = PetscOptionsBool("-mat_error_if_failure","Generate an error if an error occurs when factoring the matrix","MatSetErrorIfFailure",B->erroriffailure,&B->erroriffailure,NULL);CHKERRQ(ierr); 229840d65ccSBarry Smith 230f3be49caSLisandro Dalcin if (B->ops->setfromoptions) { 231e55864a3SBarry Smith ierr = (*B->ops->setfromoptions)(PetscOptionsObject,B);CHKERRQ(ierr); 232593c1905SSatish Balay } 233f3be49caSLisandro Dalcin 23469df5c0cSJed Brown flg = PETSC_FALSE; 23569df5c0cSJed Brown ierr = PetscOptionsBool("-mat_new_nonzero_location_err","Generate an error if new nonzeros are created in the matrix structure (useful to test preallocation)","MatSetOption",flg,&flg,&set);CHKERRQ(ierr); 23669df5c0cSJed Brown if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);} 23769df5c0cSJed Brown flg = PETSC_FALSE; 23869df5c0cSJed Brown ierr = PetscOptionsBool("-mat_new_nonzero_allocation_err","Generate an error if new nonzeros are allocated in the matrix structure (useful to test preallocation)","MatSetOption",flg,&flg,&set);CHKERRQ(ierr); 23969df5c0cSJed Brown if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);} 240478db826SMatthew G. Knepley flg = PETSC_FALSE; 241478db826SMatthew G. Knepley ierr = PetscOptionsBool("-mat_ignore_zero_entries","For AIJ/IS matrices this will stop zero values from creating a zero location in the matrix","MatSetOption",flg,&flg,&set);CHKERRQ(ierr); 242478db826SMatthew G. Knepley if (set) {ierr = MatSetOption(B,MAT_IGNORE_ZERO_ENTRIES,flg);CHKERRQ(ierr);} 24369df5c0cSJed Brown 2445d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 2450633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B);CHKERRQ(ierr); 246f3be49caSLisandro Dalcin ierr = PetscOptionsEnd();CHKERRQ(ierr); 2473a40ed3dSBarry Smith PetscFunctionReturn(0); 2487807a1faSBarry Smith } 2497807a1faSBarry Smith 250987010e7SBarry Smith /*@C 251e8bd9bafSStefano Zampini MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 25263562e91SJed Brown 25363562e91SJed Brown Collective on Mat 25463562e91SJed Brown 25563562e91SJed Brown Input Arguments: 25663562e91SJed Brown + A - matrix being preallocated 25763562e91SJed Brown . bs - block size 25841319c1dSStefano Zampini . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 25941319c1dSStefano Zampini . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 26041319c1dSStefano Zampini . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 26141319c1dSStefano Zampini - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 26263562e91SJed Brown 26363562e91SJed Brown Level: beginner 26463562e91SJed Brown 265ab978733SBarry Smith .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), 266ab978733SBarry Smith PetscSplitOwnership() 26763562e91SJed Brown @*/ 2684cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[]) 26963562e91SJed Brown { 27063562e91SJed Brown PetscErrorCode ierr; 27141319c1dSStefano Zampini PetscInt cbs; 27263562e91SJed Brown void (*aij)(void); 273e8bd9bafSStefano Zampini void (*is)(void); 274990279feSStefano Zampini void (*hyp)(void) = NULL; 27563562e91SJed Brown 27663562e91SJed Brown PetscFunctionBegin; 27741319c1dSStefano Zampini if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 278dec54756SJed Brown ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 27941319c1dSStefano Zampini } 280dec54756SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 281dec54756SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 28241319c1dSStefano Zampini ierr = MatGetBlockSizes(A,&bs,&cbs);CHKERRQ(ierr); 28341319c1dSStefano Zampini /* these routines assumes bs == cbs, this should be checked somehow */ 2843e5f4774SJed Brown ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr); 2853e5f4774SJed Brown ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr); 2863e5f4774SJed Brown ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr); 2873e5f4774SJed Brown ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr); 28863562e91SJed Brown /* 289e8bd9bafSStefano 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 29063562e91SJed Brown good before going on with it. 29163562e91SJed Brown */ 29263562e91SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 293e8bd9bafSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is);CHKERRQ(ierr); 294990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 295990279feSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp);CHKERRQ(ierr); 296990279feSStefano Zampini #endif 297990279feSStefano Zampini if (!aij && !is && !hyp) { 29863562e91SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 29963562e91SJed Brown } 300990279feSStefano Zampini if (aij || is || hyp) { 30141319c1dSStefano Zampini if (bs == cbs && bs == 1) { 3023e5f4774SJed Brown ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr); 3033e5f4774SJed Brown ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 304e8bd9bafSStefano Zampini ierr = MatISSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 305990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 306990279feSStefano Zampini ierr = MatHYPRESetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 307990279feSStefano Zampini #endif 3083e5f4774SJed Brown } else { /* Convert block-row precallocation to scalar-row */ 30963562e91SJed Brown PetscInt i,m,*sdnnz,*sonnz; 3100298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 311dcca6d9dSJed Brown ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr); 312dec54756SJed Brown for (i=0; i<m; i++) { 31341319c1dSStefano Zampini if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs; 31441319c1dSStefano Zampini if (onnz) sonnz[i] = onnz[i/bs] * cbs; 31563562e91SJed Brown } 3160298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr); 3170298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 318e8bd9bafSStefano Zampini ierr = MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 319990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 320990279feSStefano Zampini ierr = MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 321990279feSStefano Zampini #endif 32263562e91SJed Brown ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr); 32363562e91SJed Brown } 32463562e91SJed Brown } 32563562e91SJed Brown PetscFunctionReturn(0); 32663562e91SJed Brown } 32763562e91SJed Brown 328273d9f13SBarry Smith /* 329eb6b5d47SBarry Smith Merges some information from Cs header to A; the C object is then destroyed 330d0f46423SBarry Smith 331d0f46423SBarry Smith This is somewhat different from MatHeaderReplace() it would be nice to merge the code 332273d9f13SBarry Smith */ 33328be2f97SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat *C) 334273d9f13SBarry Smith { 335dfbe8321SBarry Smith PetscErrorCode ierr; 336d44834fbSBarry Smith PetscInt refct; 33773107ff1SLisandro Dalcin PetscOps Abops; 33873107ff1SLisandro Dalcin struct _MatOps Aops; 3394768301cSVaclav Hapla char *mtype,*mname,*mprefix; 3404222ddf1SHong Zhang Mat_Product *product; 341273d9f13SBarry Smith 342273d9f13SBarry Smith PetscFunctionBegin; 3431dc04de0SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3441dc04de0SStefano Zampini PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 3451dc04de0SStefano Zampini if (A == *C) PetscFunctionReturn(0); 3461dc04de0SStefano Zampini PetscCheckSameComm(A,1,*C,2); 347273d9f13SBarry Smith /* save the parts of A we need */ 34873107ff1SLisandro Dalcin Abops = ((PetscObject)A)->bops[0]; 34973107ff1SLisandro Dalcin Aops = A->ops[0]; 3507adad957SLisandro Dalcin refct = ((PetscObject)A)->refct; 3515c9eb25fSBarry Smith mtype = ((PetscObject)A)->type_name; 3525c9eb25fSBarry Smith mname = ((PetscObject)A)->name; 3534768301cSVaclav Hapla mprefix = ((PetscObject)A)->prefix; 3544222ddf1SHong Zhang product = A->product; 35530735b05SKris Buschelman 3565c9eb25fSBarry Smith /* zero these so the destroy below does not free them */ 357f4259b30SLisandro Dalcin ((PetscObject)A)->type_name = NULL; 358f4259b30SLisandro Dalcin ((PetscObject)A)->name = NULL; 3595c9eb25fSBarry Smith 3607c99f97cSSatish Balay /* free all the interior data structures from mat */ 3617c99f97cSSatish Balay ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 3627c99f97cSSatish Balay 36334136279SStefano Zampini ierr = PetscFree(A->defaultvectype);CHKERRQ(ierr); 3646bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr); 3656bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr); 366140e18c1SBarry Smith ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr); 367140e18c1SBarry Smith ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr); 368273d9f13SBarry Smith 369273d9f13SBarry Smith /* copy C over to A */ 37028be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 371273d9f13SBarry Smith 372273d9f13SBarry Smith /* return the parts of A we saved */ 37373107ff1SLisandro Dalcin ((PetscObject)A)->bops[0] = Abops; 37473107ff1SLisandro Dalcin A->ops[0] = Aops; 3757adad957SLisandro Dalcin ((PetscObject)A)->refct = refct; 3767adad957SLisandro Dalcin ((PetscObject)A)->type_name = mtype; 3777adad957SLisandro Dalcin ((PetscObject)A)->name = mname; 3784768301cSVaclav Hapla ((PetscObject)A)->prefix = mprefix; 3794222ddf1SHong Zhang A->product = product; 380273d9f13SBarry Smith 3815c9eb25fSBarry Smith /* since these two are copied into A we do not want them destroyed in C */ 382f4259b30SLisandro Dalcin ((PetscObject)*C)->qlist = NULL; 383f4259b30SLisandro Dalcin ((PetscObject)*C)->olist = NULL; 38426fbe8dcSKarl Rupp 38528be2f97SBarry Smith ierr = PetscHeaderDestroy(C);CHKERRQ(ierr); 386273d9f13SBarry Smith PetscFunctionReturn(0); 387273d9f13SBarry Smith } 3888ab5b326SKris Buschelman /* 389eb6b5d47SBarry Smith Replace A's header with that of C; the C object is then destroyed 390d0f46423SBarry Smith 391eb6b5d47SBarry Smith This is essentially code moved from MatDestroy() 392eb6b5d47SBarry Smith 393eb6b5d47SBarry Smith This is somewhat different from MatHeaderMerge() it would be nice to merge the code 394b30237c6SBarry Smith 395b30237c6SBarry Smith Used in DM hence is declared PETSC_EXTERN 3968ab5b326SKris Buschelman */ 39728be2f97SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C) 3988ab5b326SKris Buschelman { 3998ab5b326SKris Buschelman PetscErrorCode ierr; 40027b31e29SJed Brown PetscInt refct; 401fefd9316SJose E. Roman PetscObjectState state; 40228be2f97SBarry Smith struct _p_Mat buffer; 40381fa06acSBarry Smith MatStencilInfo stencil; 4048ab5b326SKris Buschelman 4058ab5b326SKris Buschelman PetscFunctionBegin; 40627b31e29SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 40728be2f97SBarry Smith PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 40828be2f97SBarry Smith if (A == *C) PetscFunctionReturn(0); 40928be2f97SBarry Smith PetscCheckSameComm(A,1,*C,2); 41028be2f97SBarry Smith if (((PetscObject)*C)->refct != 1) SETERRQ1(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Object C has refct %D > 1, would leave hanging reference",((PetscObject)*C)->refct); 4116d7c1e57SBarry Smith 41228be2f97SBarry Smith /* swap C and A */ 41327b31e29SJed Brown refct = ((PetscObject)A)->refct; 414fefd9316SJose E. Roman state = ((PetscObject)A)->state; 41581fa06acSBarry Smith stencil = A->stencil; 41628be2f97SBarry Smith ierr = PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));CHKERRQ(ierr); 41728be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 41828be2f97SBarry Smith ierr = PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));CHKERRQ(ierr); 41927b31e29SJed Brown ((PetscObject)A)->refct = refct; 420fefd9316SJose E. Roman ((PetscObject)A)->state = state + 1; 42181fa06acSBarry Smith A->stencil = stencil; 42226fbe8dcSKarl Rupp 423c32d4117SBarry Smith ((PetscObject)*C)->refct = 1; 424b117bc16SStefano Zampini ierr = MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL);CHKERRQ(ierr); 42528be2f97SBarry Smith ierr = MatDestroy(C);CHKERRQ(ierr); 4268ab5b326SKris Buschelman PetscFunctionReturn(0); 4278ab5b326SKris Buschelman } 428e7e92044SBarry Smith 429e7e92044SBarry Smith /*@ 430b470e4b4SRichard Tran Mills MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 431e7e92044SBarry Smith 432e7e92044SBarry Smith Input Parameters: 433e7e92044SBarry Smith + A - the matrix 434b470e4b4SRichard Tran Mills - flg - bind to the CPU if value of PETSC_TRUE 435e7e92044SBarry Smith 43690ea27d8SSatish Balay Level: intermediate 437e7e92044SBarry Smith @*/ 438b470e4b4SRichard Tran Mills PetscErrorCode MatBindToCPU(Mat A,PetscBool flg) 439e7e92044SBarry Smith { 440e7e92044SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 441e7e92044SBarry Smith PetscErrorCode ierr; 442e7e92044SBarry Smith 4437d871021SStefano Zampini PetscFunctionBegin; 4442ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 4452ffa8ee7SStefano Zampini PetscValidLogicalCollectiveBool(A,flg,2); 446b470e4b4SRichard Tran Mills if (A->boundtocpu == flg) PetscFunctionReturn(0); 447b470e4b4SRichard Tran Mills A->boundtocpu = flg; 44892f9df4aSRichard Tran Mills if (A->ops->bindtocpu) { 44992f9df4aSRichard Tran Mills ierr = (*A->ops->bindtocpu)(A,flg);CHKERRQ(ierr); 450e7e92044SBarry Smith } 451e7e92044SBarry Smith PetscFunctionReturn(0); 4527d871021SStefano Zampini #else 4532ffa8ee7SStefano Zampini PetscFunctionBegin; 4542ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 4552ffa8ee7SStefano Zampini PetscValidLogicalCollectiveBool(A,flg,2); 4562ffa8ee7SStefano Zampini PetscFunctionReturn(0); 4577d871021SStefano Zampini #endif 458e7e92044SBarry Smith } 4597e8381f9SStefano Zampini 4607e8381f9SStefano Zampini PetscErrorCode MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode) 4617e8381f9SStefano Zampini { 4627e8381f9SStefano Zampini IS is_coo_i,is_coo_j; 4637e8381f9SStefano Zampini const PetscInt *coo_i,*coo_j; 4647e8381f9SStefano Zampini PetscInt n,n_i,n_j; 4657e8381f9SStefano Zampini PetscScalar zero = 0.; 4667e8381f9SStefano Zampini PetscErrorCode ierr; 4677e8381f9SStefano Zampini 4687e8381f9SStefano Zampini PetscFunctionBegin; 4697e8381f9SStefano Zampini ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_i",(PetscObject*)&is_coo_i);CHKERRQ(ierr); 4707e8381f9SStefano Zampini ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_j",(PetscObject*)&is_coo_j);CHKERRQ(ierr); 4717e8381f9SStefano Zampini if (!is_coo_i) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_i IS"); 4727e8381f9SStefano Zampini if (!is_coo_j) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_j IS"); 4737e8381f9SStefano Zampini ierr = ISGetLocalSize(is_coo_i,&n_i);CHKERRQ(ierr); 4747e8381f9SStefano Zampini ierr = ISGetLocalSize(is_coo_j,&n_j);CHKERRQ(ierr); 4757e8381f9SStefano Zampini if (n_i != n_j) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_COR,"Wrong local size %D != %D",n_i,n_j); 4767e8381f9SStefano Zampini ierr = ISGetIndices(is_coo_i,&coo_i);CHKERRQ(ierr); 4777e8381f9SStefano Zampini ierr = ISGetIndices(is_coo_j,&coo_j);CHKERRQ(ierr); 478e61fc153SStefano Zampini if (imode != ADD_VALUES) { 479e61fc153SStefano Zampini ierr = MatZeroEntries(A);CHKERRQ(ierr); 480e61fc153SStefano Zampini } 4817e8381f9SStefano Zampini for (n = 0; n < n_i; n++) { 482e61fc153SStefano Zampini ierr = MatSetValue(A,coo_i[n],coo_j[n],coo_v ? coo_v[n] : zero,ADD_VALUES);CHKERRQ(ierr); 4837e8381f9SStefano Zampini } 4847e8381f9SStefano Zampini ierr = ISRestoreIndices(is_coo_i,&coo_i);CHKERRQ(ierr); 4857e8381f9SStefano Zampini ierr = ISRestoreIndices(is_coo_j,&coo_j);CHKERRQ(ierr); 4867e8381f9SStefano Zampini ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4877e8381f9SStefano Zampini ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 4887e8381f9SStefano Zampini PetscFunctionReturn(0); 4897e8381f9SStefano Zampini } 4907e8381f9SStefano Zampini 4917e8381f9SStefano Zampini PetscErrorCode MatSetPreallocationCOO_Basic(Mat A,PetscInt ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 4927e8381f9SStefano Zampini { 4937e8381f9SStefano Zampini Mat preallocator; 4947e8381f9SStefano Zampini IS is_coo_i,is_coo_j; 4957e8381f9SStefano Zampini PetscScalar zero = 0.0; 4967e8381f9SStefano Zampini PetscInt n; 4977e8381f9SStefano Zampini PetscErrorCode ierr; 4987e8381f9SStefano Zampini 4997e8381f9SStefano Zampini PetscFunctionBegin; 5007e8381f9SStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 5017e8381f9SStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 5027e8381f9SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),&preallocator);CHKERRQ(ierr); 5037e8381f9SStefano Zampini ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 5047e8381f9SStefano Zampini ierr = MatSetSizes(preallocator,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 5057e8381f9SStefano Zampini ierr = MatSetLayouts(preallocator,A->rmap,A->cmap);CHKERRQ(ierr); 5067e8381f9SStefano Zampini ierr = MatSetUp(preallocator);CHKERRQ(ierr); 5077e8381f9SStefano Zampini for (n = 0; n < ncoo; n++) { 5087e8381f9SStefano Zampini ierr = MatSetValue(preallocator,coo_i[n],coo_j[n],zero,INSERT_VALUES);CHKERRQ(ierr); 5097e8381f9SStefano Zampini } 5107e8381f9SStefano Zampini ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5117e8381f9SStefano Zampini ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5127e8381f9SStefano Zampini ierr = MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A);CHKERRQ(ierr); 5137e8381f9SStefano Zampini ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 5147e8381f9SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_i,PETSC_COPY_VALUES,&is_coo_i);CHKERRQ(ierr); 5157e8381f9SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_j,PETSC_COPY_VALUES,&is_coo_j);CHKERRQ(ierr); 5167e8381f9SStefano Zampini ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_i",(PetscObject)is_coo_i);CHKERRQ(ierr); 5177e8381f9SStefano Zampini ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_j",(PetscObject)is_coo_j);CHKERRQ(ierr); 5187e8381f9SStefano Zampini ierr = ISDestroy(&is_coo_i);CHKERRQ(ierr); 5197e8381f9SStefano Zampini ierr = ISDestroy(&is_coo_j);CHKERRQ(ierr); 5207e8381f9SStefano Zampini PetscFunctionReturn(0); 5217e8381f9SStefano Zampini } 5227e8381f9SStefano Zampini 5237e8381f9SStefano Zampini /*@C 5247e8381f9SStefano Zampini MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries 5257e8381f9SStefano Zampini 5267e8381f9SStefano Zampini Collective on Mat 5277e8381f9SStefano Zampini 5287e8381f9SStefano Zampini Input Arguments: 5297e8381f9SStefano Zampini + A - matrix being preallocated 5307e8381f9SStefano Zampini . ncoo - number of entries in the locally owned part of the parallel matrix 5317e8381f9SStefano Zampini . coo_i - row indices 5327e8381f9SStefano Zampini - coo_j - column indices 5337e8381f9SStefano Zampini 5347e8381f9SStefano Zampini Level: beginner 5357e8381f9SStefano Zampini 536*bfcc3627SStefano Zampini Notes: Entries can be repeated, see MatSetValuesCOO(). Currently optimized for cuSPARSE matrices only. 5377e8381f9SStefano Zampini 5387e8381f9SStefano Zampini .seealso: MatSetValuesCOO(), MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation() 5397e8381f9SStefano Zampini @*/ 5407e8381f9SStefano Zampini PetscErrorCode MatSetPreallocationCOO(Mat A,PetscInt ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 5417e8381f9SStefano Zampini { 5427e8381f9SStefano Zampini PetscErrorCode (*f)(Mat,PetscInt,const PetscInt[],const PetscInt[]) = NULL; 5437e8381f9SStefano Zampini PetscErrorCode ierr; 5447e8381f9SStefano Zampini 5457e8381f9SStefano Zampini PetscFunctionBegin; 5467e8381f9SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 5477e8381f9SStefano Zampini PetscValidType(A,1); 5487e8381f9SStefano Zampini if (ncoo) PetscValidIntPointer(coo_i,3); 5497e8381f9SStefano Zampini if (ncoo) PetscValidIntPointer(coo_j,4); 5507e8381f9SStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 5517e8381f9SStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 5527e8381f9SStefano Zampini if (PetscDefined(USE_DEBUG)) { 5537e8381f9SStefano Zampini PetscInt i; 5547e8381f9SStefano Zampini for (i = 0; i < ncoo; i++) { 5557e8381f9SStefano Zampini if (coo_i[i] < A->rmap->rstart || coo_i[i] >= A->rmap->rend) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_USER,"Invalid row index %D! Must be in [%D,%D)",coo_i[i],A->rmap->rstart,A->rmap->rend); 5567e8381f9SStefano Zampini if (coo_j[i] < 0 || coo_j[i] >= A->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Invalid col index %D! Must be in [0,%D)",coo_j[i],A->cmap->N); 5577e8381f9SStefano Zampini } 5587e8381f9SStefano Zampini } 5597e8381f9SStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOO_C",&f);CHKERRQ(ierr); 560*bfcc3627SStefano Zampini ierr = PetscLogEventBegin(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr); 5617e8381f9SStefano Zampini if (f) { 5627e8381f9SStefano Zampini ierr = (*f)(A,ncoo,coo_i,coo_j);CHKERRQ(ierr); 5637e8381f9SStefano Zampini } else { /* allow fallback, very slow */ 5647e8381f9SStefano Zampini ierr = MatSetPreallocationCOO_Basic(A,ncoo,coo_i,coo_j);CHKERRQ(ierr); 5657e8381f9SStefano Zampini } 566*bfcc3627SStefano Zampini ierr = PetscLogEventEnd(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr); 5677e8381f9SStefano Zampini PetscFunctionReturn(0); 5687e8381f9SStefano Zampini } 5697e8381f9SStefano Zampini 5707e8381f9SStefano Zampini /*@C 571*bfcc3627SStefano Zampini MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO() 5727e8381f9SStefano Zampini 5737e8381f9SStefano Zampini Collective on Mat 5747e8381f9SStefano Zampini 5757e8381f9SStefano Zampini Input Arguments: 5767e8381f9SStefano Zampini + A - matrix being preallocated 577*bfcc3627SStefano Zampini . coo_v - the matrix values (can be NULL) 5787e8381f9SStefano Zampini - imode - the insert mode 5797e8381f9SStefano Zampini 5807e8381f9SStefano Zampini Level: beginner 5817e8381f9SStefano Zampini 5827e8381f9SStefano Zampini Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO(). 583*bfcc3627SStefano Zampini When repeated entries are specified in the COO indices the coo_v values are first properly summed. 584*bfcc3627SStefano Zampini The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES). 5857e8381f9SStefano Zampini Currently optimized for cuSPARSE matrices only. 586*bfcc3627SStefano Zampini Passing coo_v == NULL is equivalent to passing an array of zeros. 5877e8381f9SStefano Zampini 5887e8381f9SStefano Zampini .seealso: MatSetPreallocationCOO(), InsertMode, INSERT_VALUES, ADD_VALUES 5897e8381f9SStefano Zampini @*/ 5907e8381f9SStefano Zampini PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode) 5917e8381f9SStefano Zampini { 5927e8381f9SStefano Zampini PetscErrorCode (*f)(Mat,const PetscScalar[],InsertMode) = NULL; 5937e8381f9SStefano Zampini PetscErrorCode ierr; 5947e8381f9SStefano Zampini 5957e8381f9SStefano Zampini PetscFunctionBegin; 5967e8381f9SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 5977e8381f9SStefano Zampini PetscValidType(A,1); 5987e8381f9SStefano Zampini MatCheckPreallocated(A,1); 599*bfcc3627SStefano Zampini PetscValidLogicalCollectiveEnum(A,imode,3); 6007e8381f9SStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetValuesCOO_C",&f);CHKERRQ(ierr); 601*bfcc3627SStefano Zampini ierr = PetscLogEventBegin(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr); 6027e8381f9SStefano Zampini if (f) { 6037e8381f9SStefano Zampini ierr = (*f)(A,coo_v,imode);CHKERRQ(ierr); 6047e8381f9SStefano Zampini } else { /* allow fallback */ 6057e8381f9SStefano Zampini ierr = MatSetValuesCOO_Basic(A,coo_v,imode);CHKERRQ(ierr); 6067e8381f9SStefano Zampini } 607*bfcc3627SStefano Zampini ierr = PetscLogEventEnd(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr); 608e61fc153SStefano Zampini ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr); 6097e8381f9SStefano Zampini PetscFunctionReturn(0); 6107e8381f9SStefano Zampini } 611