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; 906f3d89d0SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 916f3d89d0SStefano Zampini B->boundtocpu = PETSC_TRUE; 926f3d89d0SStefano Zampini #endif 93273d9f13SBarry Smith *A = B; 94273d9f13SBarry Smith PetscFunctionReturn(0); 95273d9f13SBarry Smith } 96273d9f13SBarry Smith 97422a814eSBarry Smith /*@ 9884d44b13SHong Zhang MatSetErrorIfFailure - Causes Mat to generate an error, for example a zero pivot, is detected. 99422a814eSBarry Smith 100422a814eSBarry Smith Logically Collective on Mat 101422a814eSBarry Smith 102422a814eSBarry Smith Input Parameters: 103422a814eSBarry Smith + mat - matrix obtained from MatCreate() 104422a814eSBarry Smith - flg - PETSC_TRUE indicates you want the error generated 105422a814eSBarry Smith 106422a814eSBarry Smith Level: advanced 107422a814eSBarry Smith 108422a814eSBarry Smith .seealso: PCSetErrorIfFailure() 109422a814eSBarry Smith @*/ 11084d44b13SHong Zhang PetscErrorCode MatSetErrorIfFailure(Mat mat,PetscBool flg) 111422a814eSBarry Smith { 112422a814eSBarry Smith PetscFunctionBegin; 113422a814eSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 114422a814eSBarry Smith PetscValidLogicalCollectiveBool(mat,flg,2); 11584d44b13SHong Zhang mat->erroriffailure = flg; 116422a814eSBarry Smith PetscFunctionReturn(0); 117422a814eSBarry Smith } 118422a814eSBarry Smith 119f69a0ea3SMatthew Knepley /*@ 120f69a0ea3SMatthew Knepley MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 121f69a0ea3SMatthew Knepley 122f69a0ea3SMatthew Knepley Collective on Mat 123f69a0ea3SMatthew Knepley 124f69a0ea3SMatthew Knepley Input Parameters: 125f69a0ea3SMatthew Knepley + A - the matrix 126f69a0ea3SMatthew Knepley . m - number of local rows (or PETSC_DECIDE) 127f69a0ea3SMatthew Knepley . n - number of local columns (or PETSC_DECIDE) 128f69a0ea3SMatthew Knepley . M - number of global rows (or PETSC_DETERMINE) 129f69a0ea3SMatthew Knepley - N - number of global columns (or PETSC_DETERMINE) 130f69a0ea3SMatthew Knepley 131f69a0ea3SMatthew Knepley Notes: 132f69a0ea3SMatthew Knepley m (n) and M (N) cannot be both PETSC_DECIDE 133f69a0ea3SMatthew Knepley If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang. 134f69a0ea3SMatthew Knepley 135f69a0ea3SMatthew Knepley If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 136f69a0ea3SMatthew Knepley user must ensure that they are chosen to be compatible with the 137f69a0ea3SMatthew Knepley vectors. To do this, one first considers the matrix-vector product 138f69a0ea3SMatthew Knepley 'y = A x'. The 'm' that is used in the above routine must match the 139f69a0ea3SMatthew Knepley local size used in the vector creation routine VecCreateMPI() for 'y'. 140f69a0ea3SMatthew Knepley Likewise, the 'n' used must match that used as the local size in 141f69a0ea3SMatthew Knepley VecCreateMPI() for 'x'. 142f69a0ea3SMatthew Knepley 143f73d5cc4SBarry Smith You cannot change the sizes once they have been set. 144f73d5cc4SBarry Smith 145f73d5cc4SBarry Smith The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called. 146f73d5cc4SBarry Smith 147f69a0ea3SMatthew Knepley Level: beginner 148f69a0ea3SMatthew Knepley 149f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership() 150f69a0ea3SMatthew Knepley @*/ 1517087cfbeSBarry Smith PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) 152f69a0ea3SMatthew Knepley { 153f69a0ea3SMatthew Knepley PetscFunctionBegin; 1540700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 155a69c7061SStefano Zampini PetscValidLogicalCollectiveInt(A,M,4); 156a69c7061SStefano Zampini PetscValidLogicalCollectiveInt(A,N,5); 157cb9de970SPierre 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); 158cb9de970SPierre 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); 15959cb773eSBarry 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); 16059cb773eSBarry 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); 161d0f46423SBarry Smith A->rmap->n = m; 162d0f46423SBarry Smith A->cmap->n = n; 16359cb773eSBarry Smith A->rmap->N = M > -1 ? M : A->rmap->N; 16459cb773eSBarry Smith A->cmap->N = N > -1 ? N : A->cmap->N; 165f69a0ea3SMatthew Knepley PetscFunctionReturn(0); 166f69a0ea3SMatthew Knepley } 167f69a0ea3SMatthew Knepley 16805869f15SSatish Balay /*@ 169273d9f13SBarry Smith MatSetFromOptions - Creates a matrix where the type is determined 170273d9f13SBarry Smith from the options database. Generates a parallel MPI matrix if the 171273d9f13SBarry Smith communicator has more than one processor. The default matrix type is 17269b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if 1737e5f4302SBarry Smith you do not select a type in the options database. 174273d9f13SBarry Smith 175273d9f13SBarry Smith Collective on Mat 176273d9f13SBarry Smith 177273d9f13SBarry Smith Input Parameter: 178273d9f13SBarry Smith . A - the matrix 179273d9f13SBarry Smith 180273d9f13SBarry Smith Options Database Keys: 181273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 18269b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 183273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 18469b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 185273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 18669b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 187273d9f13SBarry Smith 188273d9f13SBarry Smith Even More Options Database Keys: 189273d9f13SBarry Smith See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 190273d9f13SBarry Smith for additional format-specific options. 191bd9ce289SLois Curfman McInnes 1921d69843bSLois Curfman McInnes Level: beginner 1931d69843bSLois Curfman McInnes 19469b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateAIJ(), 19569b1f4b7SBarry Smith MatCreateSeqDense(), MatCreateDense(), 19669b1f4b7SBarry Smith MatCreateSeqBAIJ(), MatCreateBAIJ(), 19769b1f4b7SBarry Smith MatCreateSeqSBAIJ(), MatCreateSBAIJ(), 198273d9f13SBarry Smith MatConvert() 1997807a1faSBarry Smith @*/ 2007087cfbeSBarry Smith PetscErrorCode MatSetFromOptions(Mat B) 2017807a1faSBarry Smith { 202dfbe8321SBarry Smith PetscErrorCode ierr; 203f3be49caSLisandro Dalcin const char *deft = MATAIJ; 204f3be49caSLisandro Dalcin char type[256]; 20569df5c0cSJed Brown PetscBool flg,set; 20616e04d98SRichard Tran Mills PetscInt bind_below = 0; 207dbb450caSBarry Smith 2083a40ed3dSBarry Smith PetscFunctionBegin; 2090700a824SBarry Smith PetscValidHeaderSpecific(B,MAT_CLASSID,1); 210f3be49caSLisandro Dalcin 2113194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr); 212535b19f3SBarry Smith 213535b19f3SBarry Smith if (B->rmap->bs < 0) { 214535b19f3SBarry Smith PetscInt newbs = -1; 215535b19f3SBarry Smith ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr); 216535b19f3SBarry Smith if (flg) { 217535b19f3SBarry Smith ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr); 218535b19f3SBarry Smith ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr); 219535b19f3SBarry Smith } 220535b19f3SBarry Smith } 221535b19f3SBarry Smith 222a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr); 223273d9f13SBarry Smith if (flg) { 224f3be49caSLisandro Dalcin ierr = MatSetType(B,type);CHKERRQ(ierr); 225f3be49caSLisandro Dalcin } else if (!((PetscObject)B)->type_name) { 226f3be49caSLisandro Dalcin ierr = MatSetType(B,deft);CHKERRQ(ierr); 227273d9f13SBarry Smith } 228f3be49caSLisandro Dalcin 229840d65ccSBarry Smith ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr); 23094ae4db5SBarry Smith ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL);CHKERRQ(ierr); 23194ae4db5SBarry Smith ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL);CHKERRQ(ierr); 2325307bdd6SStefano 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); 233840d65ccSBarry Smith 234f3be49caSLisandro Dalcin if (B->ops->setfromoptions) { 235e55864a3SBarry Smith ierr = (*B->ops->setfromoptions)(PetscOptionsObject,B);CHKERRQ(ierr); 236593c1905SSatish Balay } 237f3be49caSLisandro Dalcin 23869df5c0cSJed Brown flg = PETSC_FALSE; 23969df5c0cSJed 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); 24069df5c0cSJed Brown if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);} 24169df5c0cSJed Brown flg = PETSC_FALSE; 24269df5c0cSJed 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); 24369df5c0cSJed Brown if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);} 244478db826SMatthew G. Knepley flg = PETSC_FALSE; 245478db826SMatthew 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); 246478db826SMatthew G. Knepley if (set) {ierr = MatSetOption(B,MAT_IGNORE_ZERO_ENTRIES,flg);CHKERRQ(ierr);} 24769df5c0cSJed Brown 2481a2c6b5cSJunchao Zhang flg = PETSC_FALSE; 2491a2c6b5cSJunchao Zhang ierr = PetscOptionsBool("-mat_form_explicit_transpose","Hint to form an explicit transpose for operations like MatMultTranspose","MatSetOption",flg,&flg,&set);CHKERRQ(ierr); 2501a2c6b5cSJunchao Zhang if (set) {ierr = MatSetOption(B,MAT_FORM_EXPLICIT_TRANSPOSE,flg);CHKERRQ(ierr);} 2511a2c6b5cSJunchao Zhang 25216e04d98SRichard Tran Mills /* Bind to CPU if below a user-specified size threshold. 25316e04d98SRichard Tran Mills * This perhaps belongs in the options for the GPU Mat types, but MatBindToCPU() does nothing when called on non-GPU types, 25416e04d98SRichard Tran Mills * and putting it here makes is more maintainable than duplicating this for all. */ 25516e04d98SRichard Tran Mills ierr = PetscOptionsInt("-mat_bind_below","Set the size threshold (in local rows) below which the Mat is bound to the CPU","MatBindToCPU",bind_below,&bind_below,&flg);CHKERRQ(ierr); 25616e04d98SRichard Tran Mills if (flg && B->rmap->n < bind_below) { 25716e04d98SRichard Tran Mills ierr = MatBindToCPU(B,PETSC_TRUE);CHKERRQ(ierr); 25816e04d98SRichard Tran Mills } 25916e04d98SRichard Tran Mills 2605d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 2610633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B);CHKERRQ(ierr); 262f3be49caSLisandro Dalcin ierr = PetscOptionsEnd();CHKERRQ(ierr); 2633a40ed3dSBarry Smith PetscFunctionReturn(0); 2647807a1faSBarry Smith } 2657807a1faSBarry Smith 266987010e7SBarry Smith /*@C 267e8bd9bafSStefano Zampini MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 26863562e91SJed Brown 26963562e91SJed Brown Collective on Mat 27063562e91SJed Brown 2714165533cSJose E. Roman Input Parameters: 27263562e91SJed Brown + A - matrix being preallocated 27363562e91SJed Brown . bs - block size 27441319c1dSStefano Zampini . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 27541319c1dSStefano Zampini . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 27641319c1dSStefano Zampini . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 27741319c1dSStefano Zampini - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 27863562e91SJed Brown 27963562e91SJed Brown Level: beginner 28063562e91SJed Brown 281ab978733SBarry Smith .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), 282ab978733SBarry Smith PetscSplitOwnership() 28363562e91SJed Brown @*/ 2844cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[]) 28563562e91SJed Brown { 28663562e91SJed Brown PetscErrorCode ierr; 28741319c1dSStefano Zampini PetscInt cbs; 28863562e91SJed Brown void (*aij)(void); 289e8bd9bafSStefano Zampini void (*is)(void); 290990279feSStefano Zampini void (*hyp)(void) = NULL; 29163562e91SJed Brown 29263562e91SJed Brown PetscFunctionBegin; 29341319c1dSStefano Zampini if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 294dec54756SJed Brown ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 29541319c1dSStefano Zampini } 296dec54756SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 297dec54756SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 29841319c1dSStefano Zampini ierr = MatGetBlockSizes(A,&bs,&cbs);CHKERRQ(ierr); 29941319c1dSStefano Zampini /* these routines assumes bs == cbs, this should be checked somehow */ 3003e5f4774SJed Brown ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr); 3013e5f4774SJed Brown ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr); 3023e5f4774SJed Brown ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr); 3033e5f4774SJed Brown ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr); 30463562e91SJed Brown /* 305e8bd9bafSStefano 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 30663562e91SJed Brown good before going on with it. 30763562e91SJed Brown */ 30863562e91SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 309e8bd9bafSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is);CHKERRQ(ierr); 310990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 311990279feSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp);CHKERRQ(ierr); 312990279feSStefano Zampini #endif 313990279feSStefano Zampini if (!aij && !is && !hyp) { 31463562e91SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 31563562e91SJed Brown } 316990279feSStefano Zampini if (aij || is || hyp) { 31741319c1dSStefano Zampini if (bs == cbs && bs == 1) { 3183e5f4774SJed Brown ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr); 3193e5f4774SJed Brown ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 320e8bd9bafSStefano Zampini ierr = MatISSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 321990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 322990279feSStefano Zampini ierr = MatHYPRESetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 323990279feSStefano Zampini #endif 3243e5f4774SJed Brown } else { /* Convert block-row precallocation to scalar-row */ 32563562e91SJed Brown PetscInt i,m,*sdnnz,*sonnz; 3260298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 327dcca6d9dSJed Brown ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr); 328dec54756SJed Brown for (i=0; i<m; i++) { 32941319c1dSStefano Zampini if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs; 33041319c1dSStefano Zampini if (onnz) sonnz[i] = onnz[i/bs] * cbs; 33163562e91SJed Brown } 3320298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr); 3330298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 334e8bd9bafSStefano Zampini ierr = MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 335990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 336990279feSStefano Zampini ierr = MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 337990279feSStefano Zampini #endif 33863562e91SJed Brown ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr); 33963562e91SJed Brown } 34063562e91SJed Brown } 34163562e91SJed Brown PetscFunctionReturn(0); 34263562e91SJed Brown } 34363562e91SJed Brown 344273d9f13SBarry Smith /* 345eb6b5d47SBarry Smith Merges some information from Cs header to A; the C object is then destroyed 346d0f46423SBarry Smith 347d0f46423SBarry Smith This is somewhat different from MatHeaderReplace() it would be nice to merge the code 348273d9f13SBarry Smith */ 34928be2f97SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat *C) 350273d9f13SBarry Smith { 351dfbe8321SBarry Smith PetscErrorCode ierr; 352d44834fbSBarry Smith PetscInt refct; 35373107ff1SLisandro Dalcin PetscOps Abops; 35473107ff1SLisandro Dalcin struct _MatOps Aops; 3554768301cSVaclav Hapla char *mtype,*mname,*mprefix; 3564222ddf1SHong Zhang Mat_Product *product; 357*d4a972cbSStefano Zampini PetscObjectState state; 358273d9f13SBarry Smith 359273d9f13SBarry Smith PetscFunctionBegin; 3601dc04de0SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3611dc04de0SStefano Zampini PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 3621dc04de0SStefano Zampini if (A == *C) PetscFunctionReturn(0); 3631dc04de0SStefano Zampini PetscCheckSameComm(A,1,*C,2); 364273d9f13SBarry Smith /* save the parts of A we need */ 36573107ff1SLisandro Dalcin Abops = ((PetscObject)A)->bops[0]; 36673107ff1SLisandro Dalcin Aops = A->ops[0]; 3677adad957SLisandro Dalcin refct = ((PetscObject)A)->refct; 3685c9eb25fSBarry Smith mtype = ((PetscObject)A)->type_name; 3695c9eb25fSBarry Smith mname = ((PetscObject)A)->name; 370*d4a972cbSStefano Zampini state = ((PetscObject)A)->state; 3714768301cSVaclav Hapla mprefix = ((PetscObject)A)->prefix; 3724222ddf1SHong Zhang product = A->product; 37330735b05SKris Buschelman 3745c9eb25fSBarry Smith /* zero these so the destroy below does not free them */ 375f4259b30SLisandro Dalcin ((PetscObject)A)->type_name = NULL; 376f4259b30SLisandro Dalcin ((PetscObject)A)->name = NULL; 3775c9eb25fSBarry Smith 3787c99f97cSSatish Balay /* free all the interior data structures from mat */ 3797c99f97cSSatish Balay ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 3807c99f97cSSatish Balay 38134136279SStefano Zampini ierr = PetscFree(A->defaultvectype);CHKERRQ(ierr); 3826bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr); 3836bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr); 384140e18c1SBarry Smith ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr); 385140e18c1SBarry Smith ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr); 386*d4a972cbSStefano Zampini ierr = PetscComposedQuantitiesDestroy((PetscObject)A);CHKERRQ(ierr); 387273d9f13SBarry Smith 388273d9f13SBarry Smith /* copy C over to A */ 38928be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 390273d9f13SBarry Smith 391273d9f13SBarry Smith /* return the parts of A we saved */ 39273107ff1SLisandro Dalcin ((PetscObject)A)->bops[0] = Abops; 39373107ff1SLisandro Dalcin A->ops[0] = Aops; 3947adad957SLisandro Dalcin ((PetscObject)A)->refct = refct; 3957adad957SLisandro Dalcin ((PetscObject)A)->type_name = mtype; 3967adad957SLisandro Dalcin ((PetscObject)A)->name = mname; 3974768301cSVaclav Hapla ((PetscObject)A)->prefix = mprefix; 398*d4a972cbSStefano Zampini ((PetscObject)A)->state = state + 1; 3994222ddf1SHong Zhang A->product = product; 400273d9f13SBarry Smith 4015c9eb25fSBarry Smith /* since these two are copied into A we do not want them destroyed in C */ 402f4259b30SLisandro Dalcin ((PetscObject)*C)->qlist = NULL; 403f4259b30SLisandro Dalcin ((PetscObject)*C)->olist = NULL; 40426fbe8dcSKarl Rupp 40528be2f97SBarry Smith ierr = PetscHeaderDestroy(C);CHKERRQ(ierr); 406273d9f13SBarry Smith PetscFunctionReturn(0); 407273d9f13SBarry Smith } 4088ab5b326SKris Buschelman /* 409eb6b5d47SBarry Smith Replace A's header with that of C; the C object is then destroyed 410d0f46423SBarry Smith 411eb6b5d47SBarry Smith This is essentially code moved from MatDestroy() 412eb6b5d47SBarry Smith 413eb6b5d47SBarry Smith This is somewhat different from MatHeaderMerge() it would be nice to merge the code 414b30237c6SBarry Smith 415b30237c6SBarry Smith Used in DM hence is declared PETSC_EXTERN 4168ab5b326SKris Buschelman */ 41728be2f97SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C) 4188ab5b326SKris Buschelman { 4198ab5b326SKris Buschelman PetscErrorCode ierr; 42027b31e29SJed Brown PetscInt refct; 421fefd9316SJose E. Roman PetscObjectState state; 42228be2f97SBarry Smith struct _p_Mat buffer; 42381fa06acSBarry Smith MatStencilInfo stencil; 4248ab5b326SKris Buschelman 4258ab5b326SKris Buschelman PetscFunctionBegin; 42627b31e29SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 42728be2f97SBarry Smith PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 42828be2f97SBarry Smith if (A == *C) PetscFunctionReturn(0); 42928be2f97SBarry Smith PetscCheckSameComm(A,1,*C,2); 43028be2f97SBarry 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); 4316d7c1e57SBarry Smith 43228be2f97SBarry Smith /* swap C and A */ 43327b31e29SJed Brown refct = ((PetscObject)A)->refct; 434fefd9316SJose E. Roman state = ((PetscObject)A)->state; 43581fa06acSBarry Smith stencil = A->stencil; 43628be2f97SBarry Smith ierr = PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));CHKERRQ(ierr); 43728be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 43828be2f97SBarry Smith ierr = PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));CHKERRQ(ierr); 43927b31e29SJed Brown ((PetscObject)A)->refct = refct; 440fefd9316SJose E. Roman ((PetscObject)A)->state = state + 1; 44181fa06acSBarry Smith A->stencil = stencil; 44226fbe8dcSKarl Rupp 443c32d4117SBarry Smith ((PetscObject)*C)->refct = 1; 444b117bc16SStefano Zampini ierr = MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL);CHKERRQ(ierr); 44528be2f97SBarry Smith ierr = MatDestroy(C);CHKERRQ(ierr); 4468ab5b326SKris Buschelman PetscFunctionReturn(0); 4478ab5b326SKris Buschelman } 448e7e92044SBarry Smith 449e7e92044SBarry Smith /*@ 450b470e4b4SRichard Tran Mills MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 451e7e92044SBarry Smith 4522216c58aSStefano Zampini Logically collective on Mat 4532216c58aSStefano Zampini 454e7e92044SBarry Smith Input Parameters: 455e7e92044SBarry Smith + A - the matrix 456b470e4b4SRichard Tran Mills - flg - bind to the CPU if value of PETSC_TRUE 457e7e92044SBarry Smith 45890ea27d8SSatish Balay Level: intermediate 4592216c58aSStefano Zampini 4602216c58aSStefano Zampini .seealso: MatBoundToCPU() 461e7e92044SBarry Smith @*/ 462b470e4b4SRichard Tran Mills PetscErrorCode MatBindToCPU(Mat A,PetscBool flg) 463e7e92044SBarry Smith { 4647d871021SStefano Zampini PetscFunctionBegin; 4652ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 4662ffa8ee7SStefano Zampini PetscValidLogicalCollectiveBool(A,flg,2); 4672216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE) 468b470e4b4SRichard Tran Mills if (A->boundtocpu == flg) PetscFunctionReturn(0); 469b470e4b4SRichard Tran Mills A->boundtocpu = flg; 47092f9df4aSRichard Tran Mills if (A->ops->bindtocpu) { 4712216c58aSStefano Zampini PetscErrorCode ierr; 47292f9df4aSRichard Tran Mills ierr = (*A->ops->bindtocpu)(A,flg);CHKERRQ(ierr); 473e7e92044SBarry Smith } 4742216c58aSStefano Zampini #endif 475e7e92044SBarry Smith PetscFunctionReturn(0); 4762216c58aSStefano Zampini } 4772216c58aSStefano Zampini 4782216c58aSStefano Zampini /*@ 4792216c58aSStefano Zampini MatBoundToCPU - query if a matrix is bound to the CPU 4802216c58aSStefano Zampini 4812216c58aSStefano Zampini Input Parameter: 4822216c58aSStefano Zampini . A - the matrix 4832216c58aSStefano Zampini 4842216c58aSStefano Zampini Output Parameter: 4852216c58aSStefano Zampini . flg - the logical flag 4862216c58aSStefano Zampini 4872216c58aSStefano Zampini Level: intermediate 4882216c58aSStefano Zampini 4892216c58aSStefano Zampini .seealso: MatBindToCPU() 4902216c58aSStefano Zampini @*/ 4912216c58aSStefano Zampini PetscErrorCode MatBoundToCPU(Mat A,PetscBool *flg) 4922216c58aSStefano Zampini { 4932ffa8ee7SStefano Zampini PetscFunctionBegin; 4942ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 4952216c58aSStefano Zampini PetscValidPointer(flg,2); 4962216c58aSStefano Zampini #if defined(PETSC_HAVE_DEVICE) 4972216c58aSStefano Zampini *flg = A->boundtocpu; 4982216c58aSStefano Zampini #else 4992216c58aSStefano Zampini *flg = PETSC_TRUE; 5007d871021SStefano Zampini #endif 5012216c58aSStefano Zampini PetscFunctionReturn(0); 502e7e92044SBarry Smith } 5037e8381f9SStefano Zampini 5047e8381f9SStefano Zampini PetscErrorCode MatSetValuesCOO_Basic(Mat A,const PetscScalar coo_v[],InsertMode imode) 5057e8381f9SStefano Zampini { 5067e8381f9SStefano Zampini IS is_coo_i,is_coo_j; 5077e8381f9SStefano Zampini const PetscInt *coo_i,*coo_j; 5087e8381f9SStefano Zampini PetscInt n,n_i,n_j; 5097e8381f9SStefano Zampini PetscScalar zero = 0.; 5107e8381f9SStefano Zampini PetscErrorCode ierr; 5117e8381f9SStefano Zampini 5127e8381f9SStefano Zampini PetscFunctionBegin; 5137e8381f9SStefano Zampini ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_i",(PetscObject*)&is_coo_i);CHKERRQ(ierr); 5147e8381f9SStefano Zampini ierr = PetscObjectQuery((PetscObject)A,"__PETSc_coo_j",(PetscObject*)&is_coo_j);CHKERRQ(ierr); 5157e8381f9SStefano Zampini if (!is_coo_i) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_i IS"); 5167e8381f9SStefano Zampini if (!is_coo_j) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_COR,"Missing coo_j IS"); 5177e8381f9SStefano Zampini ierr = ISGetLocalSize(is_coo_i,&n_i);CHKERRQ(ierr); 5187e8381f9SStefano Zampini ierr = ISGetLocalSize(is_coo_j,&n_j);CHKERRQ(ierr); 5197e8381f9SStefano Zampini if (n_i != n_j) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_COR,"Wrong local size %D != %D",n_i,n_j); 5207e8381f9SStefano Zampini ierr = ISGetIndices(is_coo_i,&coo_i);CHKERRQ(ierr); 5217e8381f9SStefano Zampini ierr = ISGetIndices(is_coo_j,&coo_j);CHKERRQ(ierr); 522e61fc153SStefano Zampini if (imode != ADD_VALUES) { 523e61fc153SStefano Zampini ierr = MatZeroEntries(A);CHKERRQ(ierr); 524e61fc153SStefano Zampini } 5257e8381f9SStefano Zampini for (n = 0; n < n_i; n++) { 526e61fc153SStefano Zampini ierr = MatSetValue(A,coo_i[n],coo_j[n],coo_v ? coo_v[n] : zero,ADD_VALUES);CHKERRQ(ierr); 5277e8381f9SStefano Zampini } 5287e8381f9SStefano Zampini ierr = ISRestoreIndices(is_coo_i,&coo_i);CHKERRQ(ierr); 5297e8381f9SStefano Zampini ierr = ISRestoreIndices(is_coo_j,&coo_j);CHKERRQ(ierr); 5307e8381f9SStefano Zampini ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5317e8381f9SStefano Zampini ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5327e8381f9SStefano Zampini PetscFunctionReturn(0); 5337e8381f9SStefano Zampini } 5347e8381f9SStefano Zampini 5357e8381f9SStefano Zampini PetscErrorCode MatSetPreallocationCOO_Basic(Mat A,PetscInt ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 5367e8381f9SStefano Zampini { 5377e8381f9SStefano Zampini Mat preallocator; 5387e8381f9SStefano Zampini IS is_coo_i,is_coo_j; 5397e8381f9SStefano Zampini PetscScalar zero = 0.0; 5407e8381f9SStefano Zampini PetscInt n; 5417e8381f9SStefano Zampini PetscErrorCode ierr; 5427e8381f9SStefano Zampini 5437e8381f9SStefano Zampini PetscFunctionBegin; 5447e8381f9SStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 5457e8381f9SStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 5467e8381f9SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),&preallocator);CHKERRQ(ierr); 5477e8381f9SStefano Zampini ierr = MatSetType(preallocator,MATPREALLOCATOR);CHKERRQ(ierr); 5487e8381f9SStefano Zampini ierr = MatSetSizes(preallocator,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 5497e8381f9SStefano Zampini ierr = MatSetLayouts(preallocator,A->rmap,A->cmap);CHKERRQ(ierr); 5507e8381f9SStefano Zampini ierr = MatSetUp(preallocator);CHKERRQ(ierr); 5517e8381f9SStefano Zampini for (n = 0; n < ncoo; n++) { 5527e8381f9SStefano Zampini ierr = MatSetValue(preallocator,coo_i[n],coo_j[n],zero,INSERT_VALUES);CHKERRQ(ierr); 5537e8381f9SStefano Zampini } 5547e8381f9SStefano Zampini ierr = MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5557e8381f9SStefano Zampini ierr = MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5567e8381f9SStefano Zampini ierr = MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A);CHKERRQ(ierr); 5577e8381f9SStefano Zampini ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 5587e8381f9SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_i,PETSC_COPY_VALUES,&is_coo_i);CHKERRQ(ierr); 5597e8381f9SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,ncoo,coo_j,PETSC_COPY_VALUES,&is_coo_j);CHKERRQ(ierr); 5607e8381f9SStefano Zampini ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_i",(PetscObject)is_coo_i);CHKERRQ(ierr); 5617e8381f9SStefano Zampini ierr = PetscObjectCompose((PetscObject)A,"__PETSc_coo_j",(PetscObject)is_coo_j);CHKERRQ(ierr); 5627e8381f9SStefano Zampini ierr = ISDestroy(&is_coo_i);CHKERRQ(ierr); 5637e8381f9SStefano Zampini ierr = ISDestroy(&is_coo_j);CHKERRQ(ierr); 5647e8381f9SStefano Zampini PetscFunctionReturn(0); 5657e8381f9SStefano Zampini } 5667e8381f9SStefano Zampini 567b5f22f13SPierre Jolivet /*@ 5687e8381f9SStefano Zampini MatSetPreallocationCOO - set preallocation for matrices using a coordinate format of the entries 5697e8381f9SStefano Zampini 5707e8381f9SStefano Zampini Collective on Mat 5717e8381f9SStefano Zampini 5724165533cSJose E. Roman Input Parameters: 5737e8381f9SStefano Zampini + A - matrix being preallocated 5747e8381f9SStefano Zampini . ncoo - number of entries in the locally owned part of the parallel matrix 5757e8381f9SStefano Zampini . coo_i - row indices 5767e8381f9SStefano Zampini - coo_j - column indices 5777e8381f9SStefano Zampini 5787e8381f9SStefano Zampini Level: beginner 5797e8381f9SStefano Zampini 580bfcc3627SStefano Zampini Notes: Entries can be repeated, see MatSetValuesCOO(). Currently optimized for cuSPARSE matrices only. 5817e8381f9SStefano Zampini 5827e8381f9SStefano Zampini .seealso: MatSetValuesCOO(), MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation() 5837e8381f9SStefano Zampini @*/ 5847e8381f9SStefano Zampini PetscErrorCode MatSetPreallocationCOO(Mat A,PetscInt ncoo,const PetscInt coo_i[],const PetscInt coo_j[]) 5857e8381f9SStefano Zampini { 5867e8381f9SStefano Zampini PetscErrorCode (*f)(Mat,PetscInt,const PetscInt[],const PetscInt[]) = NULL; 5877e8381f9SStefano Zampini PetscErrorCode ierr; 5887e8381f9SStefano Zampini 5897e8381f9SStefano Zampini PetscFunctionBegin; 5907e8381f9SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 5917e8381f9SStefano Zampini PetscValidType(A,1); 5927e8381f9SStefano Zampini if (ncoo) PetscValidIntPointer(coo_i,3); 5937e8381f9SStefano Zampini if (ncoo) PetscValidIntPointer(coo_j,4); 5947e8381f9SStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 5957e8381f9SStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 5967e8381f9SStefano Zampini if (PetscDefined(USE_DEBUG)) { 5977e8381f9SStefano Zampini PetscInt i; 5987e8381f9SStefano Zampini for (i = 0; i < ncoo; i++) { 5997e8381f9SStefano 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); 6007e8381f9SStefano 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); 6017e8381f9SStefano Zampini } 6027e8381f9SStefano Zampini } 6037e8381f9SStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetPreallocationCOO_C",&f);CHKERRQ(ierr); 604bfcc3627SStefano Zampini ierr = PetscLogEventBegin(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr); 6057e8381f9SStefano Zampini if (f) { 6067e8381f9SStefano Zampini ierr = (*f)(A,ncoo,coo_i,coo_j);CHKERRQ(ierr); 6077e8381f9SStefano Zampini } else { /* allow fallback, very slow */ 6087e8381f9SStefano Zampini ierr = MatSetPreallocationCOO_Basic(A,ncoo,coo_i,coo_j);CHKERRQ(ierr); 6097e8381f9SStefano Zampini } 610bfcc3627SStefano Zampini ierr = PetscLogEventEnd(MAT_PreallCOO,A,0,0,0);CHKERRQ(ierr); 6117e8381f9SStefano Zampini PetscFunctionReturn(0); 6127e8381f9SStefano Zampini } 6137e8381f9SStefano Zampini 614b5f22f13SPierre Jolivet /*@ 615bfcc3627SStefano Zampini MatSetValuesCOO - set values at once in a matrix preallocated using MatSetPreallocationCOO() 6167e8381f9SStefano Zampini 6177e8381f9SStefano Zampini Collective on Mat 6187e8381f9SStefano Zampini 6194165533cSJose E. Roman Input Parameters: 6207e8381f9SStefano Zampini + A - matrix being preallocated 621bfcc3627SStefano Zampini . coo_v - the matrix values (can be NULL) 6227e8381f9SStefano Zampini - imode - the insert mode 6237e8381f9SStefano Zampini 6247e8381f9SStefano Zampini Level: beginner 6257e8381f9SStefano Zampini 6267e8381f9SStefano Zampini Notes: The values must follow the order of the indices prescribed with MatSetPreallocationCOO(). 627bfcc3627SStefano Zampini When repeated entries are specified in the COO indices the coo_v values are first properly summed. 628bfcc3627SStefano Zampini The imode flag indicates if coo_v must be added to the current values of the matrix (ADD_VALUES) or overwritten (INSERT_VALUES). 6297e8381f9SStefano Zampini Currently optimized for cuSPARSE matrices only. 630bfcc3627SStefano Zampini Passing coo_v == NULL is equivalent to passing an array of zeros. 6317e8381f9SStefano Zampini 6327e8381f9SStefano Zampini .seealso: MatSetPreallocationCOO(), InsertMode, INSERT_VALUES, ADD_VALUES 6337e8381f9SStefano Zampini @*/ 6347e8381f9SStefano Zampini PetscErrorCode MatSetValuesCOO(Mat A, const PetscScalar coo_v[], InsertMode imode) 6357e8381f9SStefano Zampini { 6367e8381f9SStefano Zampini PetscErrorCode (*f)(Mat,const PetscScalar[],InsertMode) = NULL; 6377e8381f9SStefano Zampini PetscErrorCode ierr; 6387e8381f9SStefano Zampini 6397e8381f9SStefano Zampini PetscFunctionBegin; 6407e8381f9SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 6417e8381f9SStefano Zampini PetscValidType(A,1); 6427e8381f9SStefano Zampini MatCheckPreallocated(A,1); 643bfcc3627SStefano Zampini PetscValidLogicalCollectiveEnum(A,imode,3); 6447e8381f9SStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatSetValuesCOO_C",&f);CHKERRQ(ierr); 645bfcc3627SStefano Zampini ierr = PetscLogEventBegin(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr); 6467e8381f9SStefano Zampini if (f) { 6477e8381f9SStefano Zampini ierr = (*f)(A,coo_v,imode);CHKERRQ(ierr); 6487e8381f9SStefano Zampini } else { /* allow fallback */ 6497e8381f9SStefano Zampini ierr = MatSetValuesCOO_Basic(A,coo_v,imode);CHKERRQ(ierr); 6507e8381f9SStefano Zampini } 651bfcc3627SStefano Zampini ierr = PetscLogEventEnd(MAT_SetVCOO,A,0,0,0);CHKERRQ(ierr); 652e61fc153SStefano Zampini ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr); 6537e8381f9SStefano Zampini PetscFunctionReturn(0); 6547e8381f9SStefano Zampini } 65565a9ecf2SRichard Tran Mills 65665a9ecf2SRichard Tran Mills /*@ 65765a9ecf2SRichard 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 65865a9ecf2SRichard Tran Mills 65965a9ecf2SRichard Tran Mills Input Parameters: 66065a9ecf2SRichard Tran Mills + A - the matrix 66165a9ecf2SRichard Tran Mills - flg - flag indicating whether the boundtocpu flag should be propagated 66265a9ecf2SRichard Tran Mills 66365a9ecf2SRichard Tran Mills Level: developer 66465a9ecf2SRichard Tran Mills 66565a9ecf2SRichard Tran Mills Notes: 66665a9ecf2SRichard Tran Mills If the value of flg is set to true, the following will occur: 66765a9ecf2SRichard Tran Mills 66865a9ecf2SRichard Tran Mills MatCreateSubMatrices() and MatCreateRedundantMatrix() will bind created matrices to CPU if the input matrix is bound to the CPU. 66965a9ecf2SRichard Tran Mills MatCreateVecs() will bind created vectors to CPU if the input matrix is bound to the CPU. 67065a9ecf2SRichard Tran Mills The bindingpropagates flag itself is also propagated by the above routines. 67165a9ecf2SRichard Tran Mills 67265a9ecf2SRichard Tran Mills Developer Notes: 67365a9ecf2SRichard Tran Mills If the fine-scale DMDA has the -dm_bind_below option set to true, then DMCreateInterpolationScale() calls MatSetBindingPropagates() 67465a9ecf2SRichard Tran Mills on the restriction/interpolation operator to set the bindingpropagates flag to true. 67565a9ecf2SRichard Tran Mills 676e9c74fd6SRichard Tran Mills .seealso: VecSetBindingPropagates(), MatGetBindingPropagates() 67765a9ecf2SRichard Tran Mills @*/ 67865a9ecf2SRichard Tran Mills PetscErrorCode MatSetBindingPropagates(Mat A,PetscBool flg) 67965a9ecf2SRichard Tran Mills { 68065a9ecf2SRichard Tran Mills PetscFunctionBegin; 68165a9ecf2SRichard Tran Mills PetscValidHeaderSpecific(A,MAT_CLASSID,1); 68265a9ecf2SRichard Tran Mills #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 68365a9ecf2SRichard Tran Mills A->bindingpropagates = flg; 68465a9ecf2SRichard Tran Mills #endif 68565a9ecf2SRichard Tran Mills PetscFunctionReturn(0); 68665a9ecf2SRichard Tran Mills } 687e9c74fd6SRichard Tran Mills 688e9c74fd6SRichard Tran Mills /*@ 689e9c74fd6SRichard 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 690e9c74fd6SRichard Tran Mills 691e9c74fd6SRichard Tran Mills Input Parameter: 692e9c74fd6SRichard Tran Mills . A - the matrix 693e9c74fd6SRichard Tran Mills 694e9c74fd6SRichard Tran Mills Output Parameter: 695e9c74fd6SRichard Tran Mills . flg - flag indicating whether the boundtocpu flag will be propagated 696e9c74fd6SRichard Tran Mills 697e9c74fd6SRichard Tran Mills Level: developer 698e9c74fd6SRichard Tran Mills 699e9c74fd6SRichard Tran Mills .seealso: MatSetBindingPropagates() 700e9c74fd6SRichard Tran Mills @*/ 701e9c74fd6SRichard Tran Mills PetscErrorCode MatGetBindingPropagates(Mat A,PetscBool *flg) 702e9c74fd6SRichard Tran Mills { 703e9c74fd6SRichard Tran Mills PetscFunctionBegin; 704e9c74fd6SRichard Tran Mills PetscValidHeaderSpecific(A,MAT_CLASSID,1); 705e9c74fd6SRichard Tran Mills PetscValidBoolPointer(flg,2); 706e9c74fd6SRichard Tran Mills #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 707e9c74fd6SRichard Tran Mills *flg = A->bindingpropagates; 708e9c74fd6SRichard Tran Mills #else 709e9c74fd6SRichard Tran Mills *flg = PETSC_FALSE; 710e9c74fd6SRichard Tran Mills #endif 711e9c74fd6SRichard Tran Mills PetscFunctionReturn(0); 712e9c74fd6SRichard Tran Mills } 713