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 PetscFunctionBeginUser; 18 PetscCall(PetscInitialize(&argc,&argv,0,help)); 19 comm = MPI_COMM_WORLD; 20 PetscCall(PetscMalloc2(m,&dnnz,m,&onnz)); 21 for (i=0; i<m; i++) { 22 dnnz[i] = 1; 23 onnz[i] = 1; 24 } 25 PetscCall(MatCreateAIJ(comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_DECIDE,dnnz,PETSC_DECIDE,onnz,&A)); 26 PetscCall(MatSetFromOptions(A)); 27 PetscCall(MatSetUp(A)); 28 PetscCall(PetscFree2(dnnz,onnz)); 29 30 /* This assembly shrinks memory because we do not insert enough number of values */ 31 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 32 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 33 34 /* MatResetPreallocation restores the memory required by users */ 35 PetscCall(MatResetPreallocation(A)); 36 PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 37 PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 38 PetscCall(MatGetSize(A,&M,&N)); 39 for (i=rstart; i<rend; i++) { 40 PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES)); 41 if (rend<N) { 42 PetscCall(MatSetValue(A,i,rend,1.0,INSERT_VALUES)); 43 } 44 } 45 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 46 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 47 PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 48 PetscCall(MatDestroy(&A)); 49 PetscCall(PetscFinalize()); 50 return 0; 51 } 52 53 /*TEST 54 55 test: 56 suffix: 1 57 58 test: 59 suffix: 2 60 nsize: 2 61 62 TEST*/ 63