1c4762a1bSJed Brown /* 2c4762a1bSJed Brown * 3c4762a1bSJed Brown * Created on: Sep 25, 2017 4c4762a1bSJed Brown * Author: Fande Kong 5c4762a1bSJed Brown */ 6c4762a1bSJed Brown 7c4762a1bSJed Brown static char help[] = "Illustrate the use of MatResetPreallocation.\n"; 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include "petscmat.h" 10c4762a1bSJed Brown 11c4762a1bSJed Brown int main(int argc,char **argv) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown Mat A; 14c4762a1bSJed Brown MPI_Comm comm; 15c4762a1bSJed Brown PetscInt n=5,m=5,*dnnz,*onnz,i,rstart,rend,M,N; 16c4762a1bSJed Brown 17*327415f7SBarry Smith PetscFunctionBeginUser; 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 19c4762a1bSJed Brown comm = MPI_COMM_WORLD; 209566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(m,&dnnz,m,&onnz)); 21c4762a1bSJed Brown for (i=0; i<m; i++) { 22c4762a1bSJed Brown dnnz[i] = 1; 23c4762a1bSJed Brown onnz[i] = 1; 24c4762a1bSJed Brown } 259566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_DECIDE,dnnz,PETSC_DECIDE,onnz,&A)); 269566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 279566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 289566063dSJacob Faibussowitsch PetscCall(PetscFree2(dnnz,onnz)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* This assembly shrinks memory because we do not insert enough number of values */ 319566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 329566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* MatResetPreallocation restores the memory required by users */ 359566063dSJacob Faibussowitsch PetscCall(MatResetPreallocation(A)); 369566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 379566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 389566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 39c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 409566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES)); 41c4762a1bSJed Brown if (rend<N) { 429566063dSJacob Faibussowitsch PetscCall(MatSetValue(A,i,rend,1.0,INSERT_VALUES)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown } 459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 479566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 499566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 50b122ec5aSJacob Faibussowitsch return 0; 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 53c4762a1bSJed Brown /*TEST 54c4762a1bSJed Brown 55c4762a1bSJed Brown test: 56c4762a1bSJed Brown suffix: 1 57c4762a1bSJed Brown 58c4762a1bSJed Brown test: 59c4762a1bSJed Brown suffix: 2 60c4762a1bSJed Brown nsize: 2 61c4762a1bSJed Brown 62c4762a1bSJed Brown TEST*/ 63