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++) { 257d68702bSBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);CHKERRQ(ierr); 267d68702bSBarry Smith } 277d68702bSBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 287d68702bSBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 297d68702bSBarry Smith ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption);CHKERRQ(ierr); 307d68702bSBarry Smith PetscFunctionReturn(0); 317d68702bSBarry Smith } 327d68702bSBarry Smith 3305869f15SSatish Balay /*@ 3469dd0797SLois Curfman McInnes MatCreate - Creates a matrix where the type is determined 357e5f4302SBarry Smith from either a call to MatSetType() or from the options database 367e5f4302SBarry Smith with a call to MatSetFromOptions(). The default matrix type is 3769b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ() 387e5f4302SBarry Smith if you do not set a type in the options database. If you never 397e5f4302SBarry Smith call MatSetType() or MatSetFromOptions() it will generate an 40f8ab6608SSatish Balay error when you try to use the matrix. 4183e1b59cSLois Curfman McInnes 42cb13003dSBarry Smith Collective on MPI_Comm 43cb13003dSBarry Smith 44f69a0ea3SMatthew Knepley Input Parameter: 45f69a0ea3SMatthew Knepley . comm - MPI communicator 467807a1faSBarry Smith 477807a1faSBarry Smith Output Parameter: 48dc401e71SLois Curfman McInnes . A - the matrix 49e0b365e2SLois Curfman McInnes 50273d9f13SBarry Smith Options Database Keys: 51273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 5269b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 53273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 5469b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 55273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 5669b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 57e0b365e2SLois Curfman McInnes 5883e1b59cSLois Curfman McInnes Even More Options Database Keys: 5983e1b59cSLois Curfman McInnes See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 6083e1b59cSLois Curfman McInnes for additional format-specific options. 61e0b365e2SLois Curfman McInnes 62bd9ce289SLois Curfman McInnes Notes: 63ec6e0d80SSatish Balay 64273d9f13SBarry Smith Level: beginner 65273d9f13SBarry Smith 66a2fc510eSBarry Smith User manual sections: 67a2fc510eSBarry Smith + sec_matcreate 68a2fc510eSBarry Smith - chapter_matrices 69a2fc510eSBarry Smith 70273d9f13SBarry Smith .keywords: matrix, create 71273d9f13SBarry Smith 7269b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ(), MatCreateAIJ(), 7369b1f4b7SBarry Smith MatCreateSeqDense(), MatCreateDense(), 7469b1f4b7SBarry Smith MatCreateSeqBAIJ(), MatCreateBAIJ(), 7569b1f4b7SBarry Smith MatCreateSeqSBAIJ(), MatCreateSBAIJ(), 76273d9f13SBarry Smith MatConvert() 77273d9f13SBarry Smith @*/ 787087cfbeSBarry Smith PetscErrorCode MatCreate(MPI_Comm comm,Mat *A) 79273d9f13SBarry Smith { 80273d9f13SBarry Smith Mat B; 81dfbe8321SBarry Smith PetscErrorCode ierr; 82273d9f13SBarry Smith 83273d9f13SBarry Smith PetscFunctionBegin; 84f69a0ea3SMatthew Knepley PetscValidPointer(A,2); 8597f1f81fSBarry Smith 860298fd71SBarry Smith *A = NULL; 87607a6623SBarry Smith ierr = MatInitializePackage();CHKERRQ(ierr); 888ba1e511SMatthew Knepley 8973107ff1SLisandro Dalcin ierr = PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);CHKERRQ(ierr); 9026283091SBarry Smith ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr); 9126283091SBarry Smith ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr); 9226fbe8dcSKarl Rupp 939ae29715SStefano Zampini B->congruentlayouts = -1; 94273d9f13SBarry Smith B->preallocated = PETSC_FALSE; 95273d9f13SBarry Smith *A = B; 96273d9f13SBarry Smith PetscFunctionReturn(0); 97273d9f13SBarry Smith } 98273d9f13SBarry Smith 99422a814eSBarry Smith /*@ 10084d44b13SHong Zhang MatSetErrorIfFailure - Causes Mat to generate an error, for example a zero pivot, is detected. 101422a814eSBarry Smith 102422a814eSBarry Smith Logically Collective on Mat 103422a814eSBarry Smith 104422a814eSBarry Smith Input Parameters: 105422a814eSBarry Smith + mat - matrix obtained from MatCreate() 106422a814eSBarry Smith - flg - PETSC_TRUE indicates you want the error generated 107422a814eSBarry Smith 108422a814eSBarry Smith Level: advanced 109422a814eSBarry Smith 110422a814eSBarry Smith .keywords: Mat, set, initial guess, nonzero 111422a814eSBarry Smith 112422a814eSBarry Smith .seealso: PCSetErrorIfFailure() 113422a814eSBarry Smith @*/ 11484d44b13SHong Zhang PetscErrorCode MatSetErrorIfFailure(Mat mat,PetscBool flg) 115422a814eSBarry Smith { 116422a814eSBarry Smith PetscFunctionBegin; 117422a814eSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 118422a814eSBarry Smith PetscValidLogicalCollectiveBool(mat,flg,2); 11984d44b13SHong Zhang mat->erroriffailure = flg; 120422a814eSBarry Smith PetscFunctionReturn(0); 121422a814eSBarry Smith } 122422a814eSBarry Smith 123f69a0ea3SMatthew Knepley /*@ 124f69a0ea3SMatthew Knepley MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 125f69a0ea3SMatthew Knepley 126f69a0ea3SMatthew Knepley Collective on Mat 127f69a0ea3SMatthew Knepley 128f69a0ea3SMatthew Knepley Input Parameters: 129f69a0ea3SMatthew Knepley + A - the matrix 130f69a0ea3SMatthew Knepley . m - number of local rows (or PETSC_DECIDE) 131f69a0ea3SMatthew Knepley . n - number of local columns (or PETSC_DECIDE) 132f69a0ea3SMatthew Knepley . M - number of global rows (or PETSC_DETERMINE) 133f69a0ea3SMatthew Knepley - N - number of global columns (or PETSC_DETERMINE) 134f69a0ea3SMatthew Knepley 135f69a0ea3SMatthew Knepley Notes: 136f69a0ea3SMatthew Knepley m (n) and M (N) cannot be both PETSC_DECIDE 137f69a0ea3SMatthew Knepley If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang. 138f69a0ea3SMatthew Knepley 139f69a0ea3SMatthew Knepley If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 140f69a0ea3SMatthew Knepley user must ensure that they are chosen to be compatible with the 141f69a0ea3SMatthew Knepley vectors. To do this, one first considers the matrix-vector product 142f69a0ea3SMatthew Knepley 'y = A x'. The 'm' that is used in the above routine must match the 143f69a0ea3SMatthew Knepley local size used in the vector creation routine VecCreateMPI() for 'y'. 144f69a0ea3SMatthew Knepley Likewise, the 'n' used must match that used as the local size in 145f69a0ea3SMatthew Knepley VecCreateMPI() for 'x'. 146f69a0ea3SMatthew Knepley 147f73d5cc4SBarry Smith You cannot change the sizes once they have been set. 148f73d5cc4SBarry Smith 149f73d5cc4SBarry Smith The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called. 150f73d5cc4SBarry Smith 151f69a0ea3SMatthew Knepley Level: beginner 152f69a0ea3SMatthew Knepley 153f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership() 154f69a0ea3SMatthew Knepley @*/ 1557087cfbeSBarry Smith PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) 156f69a0ea3SMatthew Knepley { 157f69a0ea3SMatthew Knepley PetscFunctionBegin; 1580700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 159b28e6c9fSJed Brown if (M > 0) PetscValidLogicalCollectiveInt(A,M,4); 160b28e6c9fSJed Brown if (N > 0) PetscValidLogicalCollectiveInt(A,N,5); 161e32f2f54SBarry Smith if (M > 0 && m > M) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local column size %D cannot be larger than global column size %D",m,M); 162e32f2f54SBarry Smith if (N > 0 && n > N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local row size %D cannot be larger than global row size %D",n,N); 163*59cb773eSBarry 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); 164*59cb773eSBarry 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); 165d0f46423SBarry Smith A->rmap->n = m; 166d0f46423SBarry Smith A->cmap->n = n; 167*59cb773eSBarry Smith A->rmap->N = M > -1 ? M : A->rmap->N; 168*59cb773eSBarry Smith A->cmap->N = N > -1 ? N : A->cmap->N; 169f69a0ea3SMatthew Knepley PetscFunctionReturn(0); 170f69a0ea3SMatthew Knepley } 171f69a0ea3SMatthew Knepley 17205869f15SSatish Balay /*@ 173273d9f13SBarry Smith MatSetFromOptions - Creates a matrix where the type is determined 174273d9f13SBarry Smith from the options database. Generates a parallel MPI matrix if the 175273d9f13SBarry Smith communicator has more than one processor. The default matrix type is 17669b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if 1777e5f4302SBarry Smith you do not select a type in the options database. 178273d9f13SBarry Smith 179273d9f13SBarry Smith Collective on Mat 180273d9f13SBarry Smith 181273d9f13SBarry Smith Input Parameter: 182273d9f13SBarry Smith . A - the matrix 183273d9f13SBarry Smith 184273d9f13SBarry Smith Options Database Keys: 185273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 18669b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 187273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 18869b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 189273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 19069b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 191273d9f13SBarry Smith 192273d9f13SBarry Smith Even More Options Database Keys: 193273d9f13SBarry Smith See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 194273d9f13SBarry Smith for additional format-specific options. 195bd9ce289SLois Curfman McInnes 1961d69843bSLois Curfman McInnes Level: beginner 1971d69843bSLois Curfman McInnes 198dc401e71SLois Curfman McInnes .keywords: matrix, create 199e0b365e2SLois Curfman McInnes 20069b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateAIJ(), 20169b1f4b7SBarry Smith MatCreateSeqDense(), MatCreateDense(), 20269b1f4b7SBarry Smith MatCreateSeqBAIJ(), MatCreateBAIJ(), 20369b1f4b7SBarry Smith MatCreateSeqSBAIJ(), MatCreateSBAIJ(), 204273d9f13SBarry Smith MatConvert() 2057807a1faSBarry Smith @*/ 2067087cfbeSBarry Smith PetscErrorCode MatSetFromOptions(Mat B) 2077807a1faSBarry Smith { 208dfbe8321SBarry Smith PetscErrorCode ierr; 209f3be49caSLisandro Dalcin const char *deft = MATAIJ; 210f3be49caSLisandro Dalcin char type[256]; 21169df5c0cSJed Brown PetscBool flg,set; 212dbb450caSBarry Smith 2133a40ed3dSBarry Smith PetscFunctionBegin; 2140700a824SBarry Smith PetscValidHeaderSpecific(B,MAT_CLASSID,1); 215f3be49caSLisandro Dalcin 2163194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr); 217535b19f3SBarry Smith 218535b19f3SBarry Smith if (B->rmap->bs < 0) { 219535b19f3SBarry Smith PetscInt newbs = -1; 220535b19f3SBarry Smith ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr); 221535b19f3SBarry Smith if (flg) { 222535b19f3SBarry Smith ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr); 223535b19f3SBarry Smith ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr); 224535b19f3SBarry Smith } 225535b19f3SBarry Smith } 226535b19f3SBarry Smith 227a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr); 228273d9f13SBarry Smith if (flg) { 229f3be49caSLisandro Dalcin ierr = MatSetType(B,type);CHKERRQ(ierr); 230f3be49caSLisandro Dalcin } else if (!((PetscObject)B)->type_name) { 231f3be49caSLisandro Dalcin ierr = MatSetType(B,deft);CHKERRQ(ierr); 232273d9f13SBarry Smith } 233f3be49caSLisandro Dalcin 234840d65ccSBarry Smith ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr); 23594ae4db5SBarry Smith ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL);CHKERRQ(ierr); 23694ae4db5SBarry Smith ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL);CHKERRQ(ierr); 237840d65ccSBarry Smith 238f3be49caSLisandro Dalcin if (B->ops->setfromoptions) { 239e55864a3SBarry Smith ierr = (*B->ops->setfromoptions)(PetscOptionsObject,B);CHKERRQ(ierr); 240593c1905SSatish Balay } 241f3be49caSLisandro Dalcin 24269df5c0cSJed Brown flg = PETSC_FALSE; 24369df5c0cSJed 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); 24469df5c0cSJed Brown if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);} 24569df5c0cSJed Brown flg = PETSC_FALSE; 24669df5c0cSJed 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); 24769df5c0cSJed Brown if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);} 24869df5c0cSJed Brown 2495d973c19SBarry Smith /* process any options handlers added with PetscObjectAddOptionsHandler() */ 2500633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)B);CHKERRQ(ierr); 251f3be49caSLisandro Dalcin ierr = PetscOptionsEnd();CHKERRQ(ierr); 2523a40ed3dSBarry Smith PetscFunctionReturn(0); 2537807a1faSBarry Smith } 2547807a1faSBarry Smith 255987010e7SBarry Smith /*@C 256e8bd9bafSStefano Zampini MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 25763562e91SJed Brown 25863562e91SJed Brown Collective on Mat 25963562e91SJed Brown 26063562e91SJed Brown Input Arguments: 26163562e91SJed Brown + A - matrix being preallocated 26263562e91SJed Brown . bs - block size 2633e5f4774SJed Brown . dnnz - number of nonzero blocks per block row of diagonal part of parallel matrix 2643e5f4774SJed Brown . onnz - number of nonzero blocks per block row of off-diagonal part of parallel matrix 2653e5f4774SJed Brown . dnnzu - number of nonzero blocks per block row of upper-triangular part of diagonal part of parallel matrix 2663e5f4774SJed Brown - onnzu - number of nonzero blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 26763562e91SJed Brown 26863562e91SJed Brown Level: beginner 26963562e91SJed Brown 270ab978733SBarry Smith .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), 271ab978733SBarry Smith PetscSplitOwnership() 27263562e91SJed Brown @*/ 2734cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[]) 27463562e91SJed Brown { 27563562e91SJed Brown PetscErrorCode ierr; 27663562e91SJed Brown void (*aij)(void); 277e8bd9bafSStefano Zampini void (*is)(void); 27863562e91SJed Brown 27963562e91SJed Brown PetscFunctionBegin; 280dec54756SJed Brown ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 28169bbac97SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 282dec54756SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 283dec54756SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 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); 294e8bd9bafSStefano Zampini if (!aij && !is) { 29563562e91SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 29663562e91SJed Brown } 297e8bd9bafSStefano Zampini if (aij || is) { 29863562e91SJed Brown if (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); 3023e5f4774SJed Brown } else { /* Convert block-row precallocation to scalar-row */ 30363562e91SJed Brown PetscInt i,m,*sdnnz,*sonnz; 3040298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 305dcca6d9dSJed Brown ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr); 306dec54756SJed Brown for (i=0; i<m; i++) { 30763562e91SJed Brown if (dnnz) sdnnz[i] = dnnz[i/bs] * bs; 30863562e91SJed Brown if (onnz) sonnz[i] = onnz[i/bs] * bs; 30963562e91SJed Brown } 3100298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr); 3110298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 312e8bd9bafSStefano Zampini ierr = MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 31363562e91SJed Brown ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr); 31463562e91SJed Brown } 31563562e91SJed Brown } 31663562e91SJed Brown PetscFunctionReturn(0); 31763562e91SJed Brown } 31863562e91SJed Brown 319273d9f13SBarry Smith /* 320eb6b5d47SBarry Smith Merges some information from Cs header to A; the C object is then destroyed 321d0f46423SBarry Smith 322d0f46423SBarry Smith This is somewhat different from MatHeaderReplace() it would be nice to merge the code 323273d9f13SBarry Smith */ 32428be2f97SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat *C) 325273d9f13SBarry Smith { 326dfbe8321SBarry Smith PetscErrorCode ierr; 327d44834fbSBarry Smith PetscInt refct; 32873107ff1SLisandro Dalcin PetscOps Abops; 32973107ff1SLisandro Dalcin struct _MatOps Aops; 330273d9f13SBarry Smith char *mtype,*mname; 331273d9f13SBarry Smith 332273d9f13SBarry Smith PetscFunctionBegin; 333273d9f13SBarry Smith /* save the parts of A we need */ 33473107ff1SLisandro Dalcin Abops = ((PetscObject)A)->bops[0]; 33573107ff1SLisandro Dalcin Aops = A->ops[0]; 3367adad957SLisandro Dalcin refct = ((PetscObject)A)->refct; 3375c9eb25fSBarry Smith mtype = ((PetscObject)A)->type_name; 3385c9eb25fSBarry Smith mname = ((PetscObject)A)->name; 33930735b05SKris Buschelman 3405c9eb25fSBarry Smith /* zero these so the destroy below does not free them */ 3415c9eb25fSBarry Smith ((PetscObject)A)->type_name = 0; 3425c9eb25fSBarry Smith ((PetscObject)A)->name = 0; 3435c9eb25fSBarry Smith 3447c99f97cSSatish Balay /* free all the interior data structures from mat */ 3457c99f97cSSatish Balay ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 3467c99f97cSSatish Balay 3476bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr); 3486bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr); 349140e18c1SBarry Smith ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr); 350140e18c1SBarry Smith ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr); 351273d9f13SBarry Smith 352273d9f13SBarry Smith /* copy C over to A */ 35328be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 354273d9f13SBarry Smith 355273d9f13SBarry Smith /* return the parts of A we saved */ 35673107ff1SLisandro Dalcin ((PetscObject)A)->bops[0] = Abops; 35773107ff1SLisandro Dalcin A->ops[0] = Aops; 3587adad957SLisandro Dalcin ((PetscObject)A)->refct = refct; 3597adad957SLisandro Dalcin ((PetscObject)A)->type_name = mtype; 3607adad957SLisandro Dalcin ((PetscObject)A)->name = mname; 361273d9f13SBarry Smith 3625c9eb25fSBarry Smith /* since these two are copied into A we do not want them destroyed in C */ 36328be2f97SBarry Smith ((PetscObject)*C)->qlist = 0; 36428be2f97SBarry Smith ((PetscObject)*C)->olist = 0; 36526fbe8dcSKarl Rupp 36628be2f97SBarry Smith ierr = PetscHeaderDestroy(C);CHKERRQ(ierr); 367273d9f13SBarry Smith PetscFunctionReturn(0); 368273d9f13SBarry Smith } 3698ab5b326SKris Buschelman /* 370eb6b5d47SBarry Smith Replace A's header with that of C; the C object is then destroyed 371d0f46423SBarry Smith 372eb6b5d47SBarry Smith This is essentially code moved from MatDestroy() 373eb6b5d47SBarry Smith 374eb6b5d47SBarry Smith This is somewhat different from MatHeaderMerge() it would be nice to merge the code 375b30237c6SBarry Smith 376b30237c6SBarry Smith Used in DM hence is declared PETSC_EXTERN 3778ab5b326SKris Buschelman */ 37828be2f97SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C) 3798ab5b326SKris Buschelman { 3808ab5b326SKris Buschelman PetscErrorCode ierr; 38127b31e29SJed Brown PetscInt refct; 382fefd9316SJose E. Roman PetscObjectState state; 38328be2f97SBarry Smith struct _p_Mat buffer; 3848ab5b326SKris Buschelman 3858ab5b326SKris Buschelman PetscFunctionBegin; 38627b31e29SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 38728be2f97SBarry Smith PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 38828be2f97SBarry Smith if (A == *C) PetscFunctionReturn(0); 38928be2f97SBarry Smith PetscCheckSameComm(A,1,*C,2); 39028be2f97SBarry 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); 3916d7c1e57SBarry Smith 39228be2f97SBarry Smith /* swap C and A */ 39327b31e29SJed Brown refct = ((PetscObject)A)->refct; 394fefd9316SJose E. Roman state = ((PetscObject)A)->state; 39528be2f97SBarry Smith ierr = PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));CHKERRQ(ierr); 39628be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 39728be2f97SBarry Smith ierr = PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));CHKERRQ(ierr); 39827b31e29SJed Brown ((PetscObject)A)->refct = refct; 399fefd9316SJose E. Roman ((PetscObject)A)->state = state + 1; 40026fbe8dcSKarl Rupp 401c32d4117SBarry Smith ((PetscObject)*C)->refct = 1; 40228be2f97SBarry Smith ierr = MatDestroy(C);CHKERRQ(ierr); 4038ab5b326SKris Buschelman PetscFunctionReturn(0); 4048ab5b326SKris Buschelman } 405