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) PetscCall(MatSetValue(A, i, rend, 1.0, INSERT_VALUES)); 42 } 43 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 44 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 45 PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 46 PetscCall(MatDestroy(&A)); 47 PetscCall(PetscFinalize()); 48 return 0; 49 } 50 51 /*TEST 52 53 test: 54 suffix: 1 55 56 test: 57 suffix: 2 58 nsize: 2 59 60 TEST*/ 61