17807a1faSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 37807a1faSBarry Smith 44a2ae208SSatish Balay #undef __FUNCT__ 5*46533700Sstefano_zampini #define __FUNCT__ "MatSetBlockSizes_Default" 6*46533700Sstefano_zampini PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat mat,PetscInt rbs, PetscInt cbs) 7*46533700Sstefano_zampini { 8*46533700Sstefano_zampini PetscFunctionBegin; 9*46533700Sstefano_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); 10*46533700Sstefano_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); 11*46533700Sstefano_zampini PetscFunctionReturn(0); 12*46533700Sstefano_zampini } 13*46533700Sstefano_zampini 14*46533700Sstefano_zampini #undef __FUNCT__ 157d68702bSBarry Smith #define __FUNCT__ "MatShift_Basic" 16db87b0f2SBarry Smith PETSC_INTERN PetscErrorCode MatShift_Basic(Mat Y,PetscScalar a) 177d68702bSBarry Smith { 187d68702bSBarry Smith PetscErrorCode ierr; 197d68702bSBarry Smith PetscInt i,start,end; 207d68702bSBarry Smith PetscScalar alpha = a; 217d68702bSBarry Smith PetscBool prevoption; 227d68702bSBarry Smith 237d68702bSBarry Smith PetscFunctionBegin; 247d68702bSBarry Smith ierr = MatGetOption(Y,MAT_NO_OFF_PROC_ENTRIES,&prevoption);CHKERRQ(ierr); 257d68702bSBarry Smith ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 267d68702bSBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 277d68702bSBarry Smith for (i=start; i<end; i++) { 287d68702bSBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);CHKERRQ(ierr); 297d68702bSBarry Smith } 307d68702bSBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 317d68702bSBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 327d68702bSBarry Smith ierr = MatSetOption(Y,MAT_NO_OFF_PROC_ENTRIES,prevoption);CHKERRQ(ierr); 337d68702bSBarry Smith PetscFunctionReturn(0); 347d68702bSBarry Smith } 357d68702bSBarry Smith 367d68702bSBarry Smith #undef __FUNCT__ 374a2ae208SSatish Balay #define __FUNCT__ "MatCreate" 3805869f15SSatish Balay /*@ 3969dd0797SLois Curfman McInnes MatCreate - Creates a matrix where the type is determined 407e5f4302SBarry Smith from either a call to MatSetType() or from the options database 417e5f4302SBarry Smith with a call to MatSetFromOptions(). The default matrix type is 4269b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() or MatCreateAIJ() 437e5f4302SBarry Smith if you do not set a type in the options database. If you never 447e5f4302SBarry Smith call MatSetType() or MatSetFromOptions() it will generate an 45f8ab6608SSatish Balay error when you try to use the matrix. 4683e1b59cSLois Curfman McInnes 47cb13003dSBarry Smith Collective on MPI_Comm 48cb13003dSBarry Smith 49f69a0ea3SMatthew Knepley Input Parameter: 50f69a0ea3SMatthew Knepley . comm - MPI communicator 517807a1faSBarry Smith 527807a1faSBarry Smith Output Parameter: 53dc401e71SLois Curfman McInnes . A - the matrix 54e0b365e2SLois Curfman McInnes 55273d9f13SBarry Smith Options Database Keys: 56273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 5769b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 58273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 5969b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 60273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 6169b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 62e0b365e2SLois Curfman McInnes 6383e1b59cSLois Curfman McInnes Even More Options Database Keys: 6483e1b59cSLois Curfman McInnes See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 6583e1b59cSLois Curfman McInnes for additional format-specific options. 66e0b365e2SLois Curfman McInnes 67bd9ce289SLois Curfman McInnes Notes: 68ec6e0d80SSatish Balay 69273d9f13SBarry Smith Level: beginner 70273d9f13SBarry Smith 71a2fc510eSBarry Smith User manual sections: 72a2fc510eSBarry Smith + sec_matcreate 73a2fc510eSBarry Smith - chapter_matrices 74a2fc510eSBarry Smith 75273d9f13SBarry Smith .keywords: matrix, create 76273d9f13SBarry Smith 7769b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ(), MatCreateAIJ(), 7869b1f4b7SBarry Smith MatCreateSeqDense(), MatCreateDense(), 7969b1f4b7SBarry Smith MatCreateSeqBAIJ(), MatCreateBAIJ(), 8069b1f4b7SBarry Smith MatCreateSeqSBAIJ(), MatCreateSBAIJ(), 81273d9f13SBarry Smith MatConvert() 82273d9f13SBarry Smith @*/ 837087cfbeSBarry Smith PetscErrorCode MatCreate(MPI_Comm comm,Mat *A) 84273d9f13SBarry Smith { 85273d9f13SBarry Smith Mat B; 86dfbe8321SBarry Smith PetscErrorCode ierr; 87273d9f13SBarry Smith 88273d9f13SBarry Smith PetscFunctionBegin; 89f69a0ea3SMatthew Knepley PetscValidPointer(A,2); 9097f1f81fSBarry Smith 910298fd71SBarry Smith *A = NULL; 92607a6623SBarry Smith ierr = MatInitializePackage();CHKERRQ(ierr); 938ba1e511SMatthew Knepley 9473107ff1SLisandro Dalcin ierr = PetscHeaderCreate(B,MAT_CLASSID,"Mat","Matrix","Mat",comm,MatDestroy,MatView);CHKERRQ(ierr); 9526283091SBarry Smith ierr = PetscLayoutCreate(comm,&B->rmap);CHKERRQ(ierr); 9626283091SBarry Smith ierr = PetscLayoutCreate(comm,&B->cmap);CHKERRQ(ierr); 9726fbe8dcSKarl Rupp 989ae29715SStefano Zampini B->congruentlayouts = -1; 99273d9f13SBarry Smith B->preallocated = PETSC_FALSE; 100273d9f13SBarry Smith *A = B; 101273d9f13SBarry Smith PetscFunctionReturn(0); 102273d9f13SBarry Smith } 103273d9f13SBarry Smith 1044a2ae208SSatish Balay #undef __FUNCT__ 10584d44b13SHong Zhang #define __FUNCT__ "MatSetErrorIfFailure" 106422a814eSBarry Smith /*@ 10784d44b13SHong Zhang MatSetErrorIfFailure - Causes Mat to generate an error, for example a zero pivot, is detected. 108422a814eSBarry Smith 109422a814eSBarry Smith Logically Collective on Mat 110422a814eSBarry Smith 111422a814eSBarry Smith Input Parameters: 112422a814eSBarry Smith + mat - matrix obtained from MatCreate() 113422a814eSBarry Smith - flg - PETSC_TRUE indicates you want the error generated 114422a814eSBarry Smith 115422a814eSBarry Smith Level: advanced 116422a814eSBarry Smith 117422a814eSBarry Smith .keywords: Mat, set, initial guess, nonzero 118422a814eSBarry Smith 119422a814eSBarry Smith .seealso: PCSetErrorIfFailure() 120422a814eSBarry Smith @*/ 12184d44b13SHong Zhang PetscErrorCode MatSetErrorIfFailure(Mat mat,PetscBool flg) 122422a814eSBarry Smith { 123422a814eSBarry Smith PetscFunctionBegin; 124422a814eSBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 125422a814eSBarry Smith PetscValidLogicalCollectiveBool(mat,flg,2); 12684d44b13SHong Zhang mat->erroriffailure = flg; 127422a814eSBarry Smith PetscFunctionReturn(0); 128422a814eSBarry Smith } 129422a814eSBarry Smith 130422a814eSBarry Smith #undef __FUNCT__ 131f69a0ea3SMatthew Knepley #define __FUNCT__ "MatSetSizes" 132f69a0ea3SMatthew Knepley /*@ 133f69a0ea3SMatthew Knepley MatSetSizes - Sets the local and global sizes, and checks to determine compatibility 134f69a0ea3SMatthew Knepley 135f69a0ea3SMatthew Knepley Collective on Mat 136f69a0ea3SMatthew Knepley 137f69a0ea3SMatthew Knepley Input Parameters: 138f69a0ea3SMatthew Knepley + A - the matrix 139f69a0ea3SMatthew Knepley . m - number of local rows (or PETSC_DECIDE) 140f69a0ea3SMatthew Knepley . n - number of local columns (or PETSC_DECIDE) 141f69a0ea3SMatthew Knepley . M - number of global rows (or PETSC_DETERMINE) 142f69a0ea3SMatthew Knepley - N - number of global columns (or PETSC_DETERMINE) 143f69a0ea3SMatthew Knepley 144f69a0ea3SMatthew Knepley Notes: 145f69a0ea3SMatthew Knepley m (n) and M (N) cannot be both PETSC_DECIDE 146f69a0ea3SMatthew Knepley If one processor calls this with M (N) of PETSC_DECIDE then all processors must, otherwise the program will hang. 147f69a0ea3SMatthew Knepley 148f69a0ea3SMatthew Knepley If PETSC_DECIDE is not used for the arguments 'm' and 'n', then the 149f69a0ea3SMatthew Knepley user must ensure that they are chosen to be compatible with the 150f69a0ea3SMatthew Knepley vectors. To do this, one first considers the matrix-vector product 151f69a0ea3SMatthew Knepley 'y = A x'. The 'm' that is used in the above routine must match the 152f69a0ea3SMatthew Knepley local size used in the vector creation routine VecCreateMPI() for 'y'. 153f69a0ea3SMatthew Knepley Likewise, the 'n' used must match that used as the local size in 154f69a0ea3SMatthew Knepley VecCreateMPI() for 'x'. 155f69a0ea3SMatthew Knepley 156f73d5cc4SBarry Smith You cannot change the sizes once they have been set. 157f73d5cc4SBarry Smith 158f73d5cc4SBarry Smith The sizes must be set before MatSetUp() or MatXXXSetPreallocation() is called. 159f73d5cc4SBarry Smith 160f69a0ea3SMatthew Knepley Level: beginner 161f69a0ea3SMatthew Knepley 162f69a0ea3SMatthew Knepley .seealso: MatGetSize(), PetscSplitOwnership() 163f69a0ea3SMatthew Knepley @*/ 1647087cfbeSBarry Smith PetscErrorCode MatSetSizes(Mat A, PetscInt m, PetscInt n, PetscInt M, PetscInt N) 165f69a0ea3SMatthew Knepley { 166f69a0ea3SMatthew Knepley PetscFunctionBegin; 1670700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 168b28e6c9fSJed Brown if (M > 0) PetscValidLogicalCollectiveInt(A,M,4); 169b28e6c9fSJed Brown if (N > 0) PetscValidLogicalCollectiveInt(A,N,5); 170e32f2f54SBarry 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); 171e32f2f54SBarry 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); 172461878b2SBarry Smith if ((A->rmap->n >= 0 && A->rmap->N >= 0) && (A->rmap->n != m || 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); 173461878b2SBarry Smith if ((A->cmap->n >= 0 && A->cmap->N >= 0) && (A->cmap->n != n || 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); 174d0f46423SBarry Smith A->rmap->n = m; 175d0f46423SBarry Smith A->cmap->n = n; 176d0f46423SBarry Smith A->rmap->N = M; 177d0f46423SBarry Smith A->cmap->N = N; 178f69a0ea3SMatthew Knepley PetscFunctionReturn(0); 179f69a0ea3SMatthew Knepley } 180f69a0ea3SMatthew Knepley 181f69a0ea3SMatthew Knepley #undef __FUNCT__ 1824a2ae208SSatish Balay #define __FUNCT__ "MatSetFromOptions" 18305869f15SSatish Balay /*@ 184273d9f13SBarry Smith MatSetFromOptions - Creates a matrix where the type is determined 185273d9f13SBarry Smith from the options database. Generates a parallel MPI matrix if the 186273d9f13SBarry Smith communicator has more than one processor. The default matrix type is 18769b1f4b7SBarry Smith AIJ, using the routines MatCreateSeqAIJ() and MatCreateAIJ() if 1887e5f4302SBarry Smith you do not select a type in the options database. 189273d9f13SBarry Smith 190273d9f13SBarry Smith Collective on Mat 191273d9f13SBarry Smith 192273d9f13SBarry Smith Input Parameter: 193273d9f13SBarry Smith . A - the matrix 194273d9f13SBarry Smith 195273d9f13SBarry Smith Options Database Keys: 196273d9f13SBarry Smith + -mat_type seqaij - AIJ type, uses MatCreateSeqAIJ() 19769b1f4b7SBarry Smith . -mat_type mpiaij - AIJ type, uses MatCreateAIJ() 198273d9f13SBarry Smith . -mat_type seqdense - dense type, uses MatCreateSeqDense() 19969b1f4b7SBarry Smith . -mat_type mpidense - dense type, uses MatCreateDense() 200273d9f13SBarry Smith . -mat_type seqbaij - block AIJ type, uses MatCreateSeqBAIJ() 20169b1f4b7SBarry Smith - -mat_type mpibaij - block AIJ type, uses MatCreateBAIJ() 202273d9f13SBarry Smith 203273d9f13SBarry Smith Even More Options Database Keys: 204273d9f13SBarry Smith See the manpages for particular formats (e.g., MatCreateSeqAIJ()) 205273d9f13SBarry Smith for additional format-specific options. 206bd9ce289SLois Curfman McInnes 2071d69843bSLois Curfman McInnes Level: beginner 2081d69843bSLois Curfman McInnes 209dc401e71SLois Curfman McInnes .keywords: matrix, create 210e0b365e2SLois Curfman McInnes 21169b1f4b7SBarry Smith .seealso: MatCreateSeqAIJ((), MatCreateAIJ(), 21269b1f4b7SBarry Smith MatCreateSeqDense(), MatCreateDense(), 21369b1f4b7SBarry Smith MatCreateSeqBAIJ(), MatCreateBAIJ(), 21469b1f4b7SBarry Smith MatCreateSeqSBAIJ(), MatCreateSBAIJ(), 215273d9f13SBarry Smith MatConvert() 2167807a1faSBarry Smith @*/ 2177087cfbeSBarry Smith PetscErrorCode MatSetFromOptions(Mat B) 2187807a1faSBarry Smith { 219dfbe8321SBarry Smith PetscErrorCode ierr; 220f3be49caSLisandro Dalcin const char *deft = MATAIJ; 221f3be49caSLisandro Dalcin char type[256]; 22269df5c0cSJed Brown PetscBool flg,set; 223dbb450caSBarry Smith 2243a40ed3dSBarry Smith PetscFunctionBegin; 2250700a824SBarry Smith PetscValidHeaderSpecific(B,MAT_CLASSID,1); 226f3be49caSLisandro Dalcin 2273194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)B);CHKERRQ(ierr); 228535b19f3SBarry Smith 229535b19f3SBarry Smith if (B->rmap->bs < 0) { 230535b19f3SBarry Smith PetscInt newbs = -1; 231535b19f3SBarry Smith ierr = PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatSetBlockSize",newbs,&newbs,&flg);CHKERRQ(ierr); 232535b19f3SBarry Smith if (flg) { 233535b19f3SBarry Smith ierr = PetscLayoutSetBlockSize(B->rmap,newbs);CHKERRQ(ierr); 234535b19f3SBarry Smith ierr = PetscLayoutSetBlockSize(B->cmap,newbs);CHKERRQ(ierr); 235535b19f3SBarry Smith } 236535b19f3SBarry Smith } 237535b19f3SBarry Smith 238a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,type,256,&flg);CHKERRQ(ierr); 239273d9f13SBarry Smith if (flg) { 240f3be49caSLisandro Dalcin ierr = MatSetType(B,type);CHKERRQ(ierr); 241f3be49caSLisandro Dalcin } else if (!((PetscObject)B)->type_name) { 242f3be49caSLisandro Dalcin ierr = MatSetType(B,deft);CHKERRQ(ierr); 243273d9f13SBarry Smith } 244f3be49caSLisandro Dalcin 245840d65ccSBarry Smith ierr = PetscOptionsName("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",&B->checksymmetryonassembly);CHKERRQ(ierr); 24694ae4db5SBarry Smith ierr = PetscOptionsReal("-mat_is_symmetric","Checks if mat is symmetric on MatAssemblyEnd()","MatIsSymmetric",B->checksymmetrytol,&B->checksymmetrytol,NULL);CHKERRQ(ierr); 24794ae4db5SBarry Smith ierr = PetscOptionsBool("-mat_null_space_test","Checks if provided null space is correct in MatAssemblyEnd()","MatSetNullSpaceTest",B->checknullspaceonassembly,&B->checknullspaceonassembly,NULL);CHKERRQ(ierr); 248840d65ccSBarry Smith 249f3be49caSLisandro Dalcin if (B->ops->setfromoptions) { 250e55864a3SBarry Smith ierr = (*B->ops->setfromoptions)(PetscOptionsObject,B);CHKERRQ(ierr); 251593c1905SSatish Balay } 252f3be49caSLisandro Dalcin 25369df5c0cSJed Brown flg = PETSC_FALSE; 25469df5c0cSJed 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); 25569df5c0cSJed Brown if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,flg);CHKERRQ(ierr);} 25669df5c0cSJed Brown flg = PETSC_FALSE; 25769df5c0cSJed 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); 25869df5c0cSJed Brown if (set) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,flg);CHKERRQ(ierr);} 25969df5c0cSJed Brown 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 2664a2ae208SSatish Balay #undef __FUNCT__ 26763562e91SJed Brown #define __FUNCT__ "MatXAIJSetPreallocation" 26863562e91SJed Brown /*@ 269e8bd9bafSStefano Zampini MatXAIJSetPreallocation - set preallocation for serial and parallel AIJ, BAIJ, and SBAIJ matrices and their unassembled versions. 27063562e91SJed Brown 27163562e91SJed Brown Collective on Mat 27263562e91SJed Brown 27363562e91SJed Brown Input Arguments: 27463562e91SJed Brown + A - matrix being preallocated 27563562e91SJed Brown . bs - block size 2763e5f4774SJed Brown . dnnz - number of nonzero blocks per block row of diagonal part of parallel matrix 2773e5f4774SJed Brown . onnz - number of nonzero blocks per block row of off-diagonal part of parallel matrix 2783e5f4774SJed Brown . dnnzu - number of nonzero blocks per block row of upper-triangular part of diagonal part of parallel matrix 2793e5f4774SJed Brown - onnzu - number of nonzero blocks per block row of upper-triangular part of off-diagonal part of parallel matrix 28063562e91SJed Brown 28163562e91SJed Brown Level: beginner 28263562e91SJed Brown 283ab978733SBarry Smith .seealso: MatSeqAIJSetPreallocation(), MatMPIAIJSetPreallocation(), MatSeqBAIJSetPreallocation(), MatMPIBAIJSetPreallocation(), MatSeqSBAIJSetPreallocation(), MatMPISBAIJSetPreallocation(), 284ab978733SBarry Smith PetscSplitOwnership() 28563562e91SJed Brown @*/ 2864cce697cSJed Brown PetscErrorCode MatXAIJSetPreallocation(Mat A,PetscInt bs,const PetscInt dnnz[],const PetscInt onnz[],const PetscInt dnnzu[],const PetscInt onnzu[]) 28763562e91SJed Brown { 28863562e91SJed Brown PetscErrorCode ierr; 28963562e91SJed Brown void (*aij)(void); 290e8bd9bafSStefano Zampini void (*is)(void); 29163562e91SJed Brown 29263562e91SJed Brown PetscFunctionBegin; 293dec54756SJed Brown ierr = MatSetBlockSize(A,bs);CHKERRQ(ierr); 29469bbac97SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 295dec54756SJed Brown ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 296dec54756SJed Brown ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 2973e5f4774SJed Brown ierr = MatSeqBAIJSetPreallocation(A,bs,0,dnnz);CHKERRQ(ierr); 2983e5f4774SJed Brown ierr = MatMPIBAIJSetPreallocation(A,bs,0,dnnz,0,onnz);CHKERRQ(ierr); 2993e5f4774SJed Brown ierr = MatSeqSBAIJSetPreallocation(A,bs,0,dnnzu);CHKERRQ(ierr); 3003e5f4774SJed Brown ierr = MatMPISBAIJSetPreallocation(A,bs,0,dnnzu,0,onnzu);CHKERRQ(ierr); 30163562e91SJed Brown /* 302e8bd9bafSStefano 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 30363562e91SJed Brown good before going on with it. 30463562e91SJed Brown */ 30563562e91SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 306e8bd9bafSStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)A,"MatISSetPreallocation_C",&is);CHKERRQ(ierr); 307e8bd9bafSStefano Zampini if (!aij && !is) { 30863562e91SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 30963562e91SJed Brown } 310e8bd9bafSStefano Zampini if (aij || is) { 31163562e91SJed Brown if (bs == 1) { 3123e5f4774SJed Brown ierr = MatSeqAIJSetPreallocation(A,0,dnnz);CHKERRQ(ierr); 3133e5f4774SJed Brown ierr = MatMPIAIJSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 314e8bd9bafSStefano Zampini ierr = MatISSetPreallocation(A,0,dnnz,0,onnz);CHKERRQ(ierr); 3153e5f4774SJed Brown } else { /* Convert block-row precallocation to scalar-row */ 31663562e91SJed Brown PetscInt i,m,*sdnnz,*sonnz; 3170298fd71SBarry Smith ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 318dcca6d9dSJed Brown ierr = PetscMalloc2((!!dnnz)*m,&sdnnz,(!!onnz)*m,&sonnz);CHKERRQ(ierr); 319dec54756SJed Brown for (i=0; i<m; i++) { 32063562e91SJed Brown if (dnnz) sdnnz[i] = dnnz[i/bs] * bs; 32163562e91SJed Brown if (onnz) sonnz[i] = onnz[i/bs] * bs; 32263562e91SJed Brown } 3230298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL);CHKERRQ(ierr); 3240298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 325e8bd9bafSStefano Zampini ierr = MatISSetPreallocation(A,0,dnnz ? sdnnz : NULL,0,onnz ? sonnz : NULL);CHKERRQ(ierr); 32663562e91SJed Brown ierr = PetscFree2(sdnnz,sonnz);CHKERRQ(ierr); 32763562e91SJed Brown } 32863562e91SJed Brown } 32963562e91SJed Brown PetscFunctionReturn(0); 33063562e91SJed Brown } 33163562e91SJed Brown 332273d9f13SBarry Smith /* 333eb6b5d47SBarry Smith Merges some information from Cs header to A; the C object is then destroyed 334d0f46423SBarry Smith 335d0f46423SBarry Smith This is somewhat different from MatHeaderReplace() it would be nice to merge the code 336273d9f13SBarry Smith */ 3374a2ae208SSatish Balay #undef __FUNCT__ 338eb6b5d47SBarry Smith #define __FUNCT__ "MatHeaderMerge" 33928be2f97SBarry Smith PetscErrorCode MatHeaderMerge(Mat A,Mat *C) 340273d9f13SBarry Smith { 341dfbe8321SBarry Smith PetscErrorCode ierr; 342d44834fbSBarry Smith PetscInt refct; 34373107ff1SLisandro Dalcin PetscOps Abops; 34473107ff1SLisandro Dalcin struct _MatOps Aops; 345273d9f13SBarry Smith char *mtype,*mname; 34630735b05SKris Buschelman void *spptr; 347273d9f13SBarry Smith 348273d9f13SBarry Smith PetscFunctionBegin; 349273d9f13SBarry Smith /* save the parts of A we need */ 35073107ff1SLisandro Dalcin Abops = ((PetscObject)A)->bops[0]; 35173107ff1SLisandro Dalcin Aops = A->ops[0]; 3527adad957SLisandro Dalcin refct = ((PetscObject)A)->refct; 3535c9eb25fSBarry Smith mtype = ((PetscObject)A)->type_name; 3545c9eb25fSBarry Smith mname = ((PetscObject)A)->name; 35530735b05SKris Buschelman spptr = A->spptr; 35630735b05SKris Buschelman 3575c9eb25fSBarry Smith /* zero these so the destroy below does not free them */ 3585c9eb25fSBarry Smith ((PetscObject)A)->type_name = 0; 3595c9eb25fSBarry Smith ((PetscObject)A)->name = 0; 3605c9eb25fSBarry Smith 3617c99f97cSSatish Balay /* free all the interior data structures from mat */ 3627c99f97cSSatish Balay ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 3637c99f97cSSatish Balay 36428be2f97SBarry Smith ierr = PetscFree((*C)->spptr);CHKERRQ(ierr); 36505b42c5fSBarry Smith 3666bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->rmap);CHKERRQ(ierr); 3676bf464f9SBarry Smith ierr = PetscLayoutDestroy(&A->cmap);CHKERRQ(ierr); 368140e18c1SBarry Smith ierr = PetscFunctionListDestroy(&((PetscObject)A)->qlist);CHKERRQ(ierr); 369140e18c1SBarry Smith ierr = PetscObjectListDestroy(&((PetscObject)A)->olist);CHKERRQ(ierr); 370273d9f13SBarry Smith 371273d9f13SBarry Smith /* copy C over to A */ 37228be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 373273d9f13SBarry Smith 374273d9f13SBarry Smith /* return the parts of A we saved */ 37573107ff1SLisandro Dalcin ((PetscObject)A)->bops[0] = Abops; 37673107ff1SLisandro Dalcin A->ops[0] = Aops; 3777adad957SLisandro Dalcin ((PetscObject)A)->refct = refct; 3787adad957SLisandro Dalcin ((PetscObject)A)->type_name = mtype; 3797adad957SLisandro Dalcin ((PetscObject)A)->name = mname; 38030735b05SKris Buschelman A->spptr = spptr; 381273d9f13SBarry Smith 3825c9eb25fSBarry Smith /* since these two are copied into A we do not want them destroyed in C */ 38328be2f97SBarry Smith ((PetscObject)*C)->qlist = 0; 38428be2f97SBarry Smith ((PetscObject)*C)->olist = 0; 38526fbe8dcSKarl Rupp 38628be2f97SBarry Smith ierr = PetscHeaderDestroy(C);CHKERRQ(ierr); 387273d9f13SBarry Smith PetscFunctionReturn(0); 388273d9f13SBarry Smith } 3898ab5b326SKris Buschelman /* 390eb6b5d47SBarry Smith Replace A's header with that of C; the C object is then destroyed 391d0f46423SBarry Smith 392eb6b5d47SBarry Smith This is essentially code moved from MatDestroy() 393eb6b5d47SBarry Smith 394eb6b5d47SBarry Smith This is somewhat different from MatHeaderMerge() it would be nice to merge the code 395b30237c6SBarry Smith 396b30237c6SBarry Smith Used in DM hence is declared PETSC_EXTERN 3978ab5b326SKris Buschelman */ 3988ab5b326SKris Buschelman #undef __FUNCT__ 3998ab5b326SKris Buschelman #define __FUNCT__ "MatHeaderReplace" 40028be2f97SBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat A,Mat *C) 4018ab5b326SKris Buschelman { 4028ab5b326SKris Buschelman PetscErrorCode ierr; 40327b31e29SJed Brown PetscInt refct; 404fefd9316SJose E. Roman PetscObjectState state; 40528be2f97SBarry Smith struct _p_Mat buffer; 4068ab5b326SKris Buschelman 4078ab5b326SKris Buschelman PetscFunctionBegin; 40827b31e29SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 40928be2f97SBarry Smith PetscValidHeaderSpecific(*C,MAT_CLASSID,2); 41028be2f97SBarry Smith if (A == *C) PetscFunctionReturn(0); 41128be2f97SBarry Smith PetscCheckSameComm(A,1,*C,2); 41228be2f97SBarry 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); 4136d7c1e57SBarry Smith 41428be2f97SBarry Smith /* swap C and A */ 41527b31e29SJed Brown refct = ((PetscObject)A)->refct; 416fefd9316SJose E. Roman state = ((PetscObject)A)->state; 41728be2f97SBarry Smith ierr = PetscMemcpy(&buffer,A,sizeof(struct _p_Mat));CHKERRQ(ierr); 41828be2f97SBarry Smith ierr = PetscMemcpy(A,*C,sizeof(struct _p_Mat));CHKERRQ(ierr); 41928be2f97SBarry Smith ierr = PetscMemcpy(*C,&buffer,sizeof(struct _p_Mat));CHKERRQ(ierr); 42027b31e29SJed Brown ((PetscObject)A)->refct = refct; 421fefd9316SJose E. Roman ((PetscObject)A)->state = state + 1; 42226fbe8dcSKarl Rupp 423c32d4117SBarry Smith ((PetscObject)*C)->refct = 1; 42428be2f97SBarry Smith ierr = MatDestroy(C);CHKERRQ(ierr); 4258ab5b326SKris Buschelman PetscFunctionReturn(0); 4268ab5b326SKris Buschelman } 427