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 Mat A; 13 MPI_Comm comm; 14 PetscInt n = 5, m = 5, *dnnz, *onnz, i, rstart, rend, M, N; 15 16 PetscFunctionBeginUser; 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) PetscCall(MatSetValue(A, i, rend, 1.0, INSERT_VALUES)); 41 } 42 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 43 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 44 PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 45 PetscCall(MatDestroy(&A)); 46 PetscCall(PetscFinalize()); 47 return 0; 48 } 49 50 /*TEST 51 52 test: 53 suffix: 1 54 55 test: 56 suffix: 2 57 nsize: 2 58 59 TEST*/ 60