1 /* 2 * 3 * Created on: Sep 25, 2017 4 * Author: Fande Kong 5 */ 6 7 static char help[] = "Illustrate the use of MatResetPreallocation.\n"; 8 9 #include "petscmat.h" 10 11 int main(int argc,char **argv) 12 { 13 Mat A; 14 MPI_Comm comm; 15 PetscInt n=5,m=5,*dnnz,*onnz,i,rstart,rend,M,N; 16 17 PetscCall(PetscInitialize(&argc,&argv,0,help)); 18 comm = MPI_COMM_WORLD; 19 PetscCall(PetscMalloc2(m,&dnnz,m,&onnz)); 20 for (i=0; i<m; i++) { 21 dnnz[i] = 1; 22 onnz[i] = 1; 23 } 24 PetscCall(MatCreateAIJ(comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_DECIDE,dnnz,PETSC_DECIDE,onnz,&A)); 25 PetscCall(MatSetFromOptions(A)); 26 PetscCall(MatSetUp(A)); 27 PetscCall(PetscFree2(dnnz,onnz)); 28 29 /* This assembly shrinks memory because we do not insert enough number of values */ 30 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 31 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 32 33 /* MatResetPreallocation restores the memory required by users */ 34 PetscCall(MatResetPreallocation(A)); 35 PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 36 PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 37 PetscCall(MatGetSize(A,&M,&N)); 38 for (i=rstart; i<rend; i++) { 39 PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES)); 40 if (rend<N) { 41 PetscCall(MatSetValue(A,i,rend,1.0,INSERT_VALUES)); 42 } 43 } 44 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 45 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 46 PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 47 PetscCall(MatDestroy(&A)); 48 PetscCall(PetscFinalize()); 49 return 0; 50 } 51 52 /*TEST 53 54 test: 55 suffix: 1 56 57 test: 58 suffix: 2 59 nsize: 2 60 61 TEST*/ 62