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); 152b28e6c9fSJed Brown if (M > 0) PetscValidLogicalCollectiveInt(A,M,4); 153b28e6c9fSJed Brown if (N > 0) 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);} 24069df5c0cSJed Brown 2415d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 2420633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B);CHKERRQ(ierr); 243f3be49caSLisandro Dalcin ierr = PetscOptionsEnd();CHKERRQ(ierr); 2443a40ed3dSBarry Smith PetscFunctionReturn(0); 2457807a1faSBarry Smith } 2467807a1faSBarry Smith 247987010e7SBarry Smith /*@C 248e8bd9bafSStefano Zampini MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 24963562e91SJed Brown 25063562e91SJed Brown Collective on Mat 25163562e91SJed Brown 25263562e91SJed Brown Input Arguments: 25363562e91SJed Brown + A - matrix being preallocated 25463562e91SJed Brown . bs - block size 25541319c1dSStefano Zampini . dnnz - number of nonzero column blocks per block row of diagonal part of parallel matrix 25641319c1dSStefano Zampini . onnz - number of nonzero column blocks per block row of off-diagonal part of parallel matrix 25741319c1dSStefano Zampini . dnnzu - number of nonzero column blocks per block row of upper-triangular part of diagonal part of parallel matrix 25841319c1dSStefano Zampini - onnzu - number of nonzero column blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 25963562e91SJed Brown 26063562e91SJed Brown Level: beginner 26163562e91SJed Brown 262ab978733SBarry Smith .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), 263ab978733SBarry Smith PetscSplitOwnership() 26463562e91SJed Brown @*/ 2654cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[]) 26663562e91SJed Brown { 26763562e91SJed Brown PetscErrorCode ierr; 26841319c1dSStefano Zampini PetscInt cbs; 26963562e91SJed Brown void (*aij)(void); 270e8bd9bafSStefano Zampini void (*is)(void); 271990279feSStefano Zampini void (*hyp)(void) = NULL; 27263562e91SJed Brown 27363562e91SJed Brown PetscFunctionBegin; 27441319c1dSStefano Zampini if (bs != PETSC_DECIDE) { /* don't mess with an already set block size */ 275dec54756SJed Brown ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 27641319c1dSStefano Zampini } 277dec54756SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 278dec54756SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 27941319c1dSStefano Zampini ierr = MatGetBlockSizes(A,&bs,&cbs);CHKERRQ(ierr); 28041319c1dSStefano Zampini /* these routines assumes bs == cbs, this should be checked somehow */ 2813e5f4774SJed Brown ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr); 2823e5f4774SJed Brown ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr); 2833e5f4774SJed Brown ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr); 2843e5f4774SJed Brown ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr); 28563562e91SJed Brown /* 286e8bd9bafSStefano 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 28763562e91SJed Brown good before going on with it. 28863562e91SJed Brown */ 28963562e91SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 290e8bd9bafSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is);CHKERRQ(ierr); 291990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 292990279feSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatHYPRESetPreallocation_C",&hyp);CHKERRQ(ierr); 293990279feSStefano Zampini #endif 294990279feSStefano Zampini if (!aij && !is && !hyp) { 29563562e91SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 29663562e91SJed Brown } 297990279feSStefano Zampini if (aij || is || hyp) { 29841319c1dSStefano Zampini if (bs == cbs && bs == 1) { 2993e5f4774SJed Brown ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr); 3003e5f4774SJed Brown ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 301e8bd9bafSStefano Zampini ierr = MatISSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 302990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 303990279feSStefano Zampini ierr = MatHYPRESetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 304990279feSStefano Zampini #endif 3053e5f4774SJed Brown } else { /* Convert block-row precallocation to scalar-row */ 30663562e91SJed Brown PetscInt i,m,*sdnnz,*sonnz; 3070298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 308dcca6d9dSJed Brown ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr); 309dec54756SJed Brown for (i=0; i<m; i++) { 31041319c1dSStefano Zampini if (dnnz) sdnnz[i] = dnnz[i/bs] * cbs; 31141319c1dSStefano Zampini if (onnz) sonnz[i] = onnz[i/bs] * cbs; 31263562e91SJed Brown } 3130298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr); 3140298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 315e8bd9bafSStefano Zampini ierr = MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 316990279feSStefano Zampini #if defined(PETSC_HAVE_HYPRE) 317990279feSStefano Zampini ierr = MatHYPRESetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 318990279feSStefano Zampini #endif 31963562e91SJed Brown ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr); 32063562e91SJed Brown } 32163562e91SJed Brown } 32263562e91SJed Brown PetscFunctionReturn(0); 32363562e91SJed Brown } 32463562e91SJed Brown 325273d9f13SBarry Smith /* 326eb6b5d47SBarry Smith Merges some information from Cs header to A; the C object is then destroyed 327d0f46423SBarry Smith 328d0f46423SBarry Smith This is somewhat different from MatHeaderReplace() it would be nice to merge the code 329273d9f13SBarry Smith */ 33028be2f97SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat *C) 331273d9f13SBarry Smith { 332dfbe8321SBarry Smith PetscErrorCode ierr; 333d44834fbSBarry Smith PetscInt refct; 33473107ff1SLisandro Dalcin PetscOps Abops; 33573107ff1SLisandro Dalcin struct _MatOps Aops; 3364768301cSVaclav Hapla char *mtype,*mname,*mprefix; 3374222ddf1SHong Zhang Mat_Product *product; 338273d9f13SBarry Smith 339273d9f13SBarry Smith PetscFunctionBegin; 340273d9f13SBarry Smith /* save the parts of A we need */ 34173107ff1SLisandro Dalcin Abops = ((PetscObject)A)->bops[0]; 34273107ff1SLisandro Dalcin Aops = A->ops[0]; 3437adad957SLisandro Dalcin refct = ((PetscObject)A)->refct; 3445c9eb25fSBarry Smith mtype = ((PetscObject)A)->type_name; 3455c9eb25fSBarry Smith mname = ((PetscObject)A)->name; 3464768301cSVaclav Hapla mprefix = ((PetscObject)A)->prefix; 3474222ddf1SHong Zhang product = A->product; 34830735b05SKris Buschelman 3495c9eb25fSBarry Smith /* zero these so the destroy below does not free them */ 3505c9eb25fSBarry Smith ((PetscObject)A)->type_name = 0; 3515c9eb25fSBarry Smith ((PetscObject)A)->name = 0; 3525c9eb25fSBarry Smith 3537c99f97cSSatish Balay /* free all the interior data structures from mat */ 3547c99f97cSSatish Balay ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 3557c99f97cSSatish Balay 35634136279SStefano Zampini ierr = PetscFree(A->defaultvectype);CHKERRQ(ierr); 3576bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr); 3586bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr); 359140e18c1SBarry Smith ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr); 360140e18c1SBarry Smith ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr); 361273d9f13SBarry Smith 362273d9f13SBarry Smith /* copy C over to A */ 36328be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 364273d9f13SBarry Smith 365273d9f13SBarry Smith /* return the parts of A we saved */ 36673107ff1SLisandro Dalcin ((PetscObject)A)->bops[0] = Abops; 36773107ff1SLisandro Dalcin A->ops[0] = Aops; 3687adad957SLisandro Dalcin ((PetscObject)A)->refct = refct; 3697adad957SLisandro Dalcin ((PetscObject)A)->type_name = mtype; 3707adad957SLisandro Dalcin ((PetscObject)A)->name = mname; 3714768301cSVaclav Hapla ((PetscObject)A)->prefix = mprefix; 3724222ddf1SHong Zhang A->product = product; 373273d9f13SBarry Smith 3745c9eb25fSBarry Smith /* since these two are copied into A we do not want them destroyed in C */ 37528be2f97SBarry Smith ((PetscObject)*C)->qlist = 0; 37628be2f97SBarry Smith ((PetscObject)*C)->olist = 0; 37726fbe8dcSKarl Rupp 37828be2f97SBarry Smith ierr = PetscHeaderDestroy(C);CHKERRQ(ierr); 379273d9f13SBarry Smith PetscFunctionReturn(0); 380273d9f13SBarry Smith } 3818ab5b326SKris Buschelman /* 382eb6b5d47SBarry Smith Replace A's header with that of C; the C object is then destroyed 383d0f46423SBarry Smith 384eb6b5d47SBarry Smith This is essentially code moved from MatDestroy() 385eb6b5d47SBarry Smith 386eb6b5d47SBarry Smith This is somewhat different from MatHeaderMerge() it would be nice to merge the code 387b30237c6SBarry Smith 388b30237c6SBarry Smith Used in DM hence is declared PETSC_EXTERN 3898ab5b326SKris Buschelman */ 39028be2f97SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C) 3918ab5b326SKris Buschelman { 3928ab5b326SKris Buschelman PetscErrorCode ierr; 39327b31e29SJed Brown PetscInt refct; 394fefd9316SJose E. Roman PetscObjectState state; 39528be2f97SBarry Smith struct _p_Mat buffer; 3968ab5b326SKris Buschelman 3978ab5b326SKris Buschelman PetscFunctionBegin; 39827b31e29SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 39928be2f97SBarry Smith PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 40028be2f97SBarry Smith if (A == *C) PetscFunctionReturn(0); 40128be2f97SBarry Smith PetscCheckSameComm(A,1,*C,2); 40228be2f97SBarry 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); 4036d7c1e57SBarry Smith 40428be2f97SBarry Smith /* swap C and A */ 40527b31e29SJed Brown refct = ((PetscObject)A)->refct; 406fefd9316SJose E. Roman state = ((PetscObject)A)->state; 40728be2f97SBarry Smith ierr = PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));CHKERRQ(ierr); 40828be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 40928be2f97SBarry Smith ierr = PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));CHKERRQ(ierr); 41027b31e29SJed Brown ((PetscObject)A)->refct = refct; 411fefd9316SJose E. Roman ((PetscObject)A)->state = state + 1; 41226fbe8dcSKarl Rupp 413c32d4117SBarry Smith ((PetscObject)*C)->refct = 1; 414b117bc16SStefano Zampini ierr = MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))NULL);CHKERRQ(ierr); 41528be2f97SBarry Smith ierr = MatDestroy(C);CHKERRQ(ierr); 4168ab5b326SKris Buschelman PetscFunctionReturn(0); 4178ab5b326SKris Buschelman } 418e7e92044SBarry Smith 419e7e92044SBarry Smith /*@ 420b470e4b4SRichard Tran Mills MatBindToCPU - marks a matrix to temporarily stay on the CPU and perform computations on the CPU 421e7e92044SBarry Smith 422e7e92044SBarry Smith Input Parameters: 423e7e92044SBarry Smith + A - the matrix 424b470e4b4SRichard Tran Mills - flg - bind to the CPU if value of PETSC_TRUE 425e7e92044SBarry Smith 42690ea27d8SSatish Balay Level: intermediate 427e7e92044SBarry Smith @*/ 428b470e4b4SRichard Tran Mills PetscErrorCode MatBindToCPU(Mat A,PetscBool flg) 429e7e92044SBarry Smith { 430e7e92044SBarry Smith #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 431e7e92044SBarry Smith PetscErrorCode ierr; 432e7e92044SBarry Smith 4337d871021SStefano Zampini PetscFunctionBegin; 434*2ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 435*2ffa8ee7SStefano Zampini PetscValidLogicalCollectiveBool(A,flg,2); 436b470e4b4SRichard Tran Mills if (A->boundtocpu == flg) PetscFunctionReturn(0); 437b470e4b4SRichard Tran Mills A->boundtocpu = flg; 43892f9df4aSRichard Tran Mills if (A->ops->bindtocpu) { 43992f9df4aSRichard Tran Mills ierr = (*A->ops->bindtocpu)(A,flg);CHKERRQ(ierr); 440e7e92044SBarry Smith } 441e7e92044SBarry Smith PetscFunctionReturn(0); 4427d871021SStefano Zampini #else 443*2ffa8ee7SStefano Zampini PetscFunctionBegin; 444*2ffa8ee7SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 445*2ffa8ee7SStefano Zampini PetscValidLogicalCollectiveBool(A,flg,2); 446*2ffa8ee7SStefano Zampini PetscFunctionReturn(0); 4477d871021SStefano Zampini #endif 448e7e92044SBarry Smith } 449