Actual source code: ex4.c

  1: /*
  2:   Contributed by: Fande Kong
  3:  */

  5: static char help[] = "Illustrate the use of MatResetPreallocation.\n";

  7: #include <petscmat.h>

  9: int main(int argc, char **argv)
 10: {
 11:   Mat      A;
 12:   MPI_Comm comm;
 13:   PetscInt n = 5, m = 5, *dnnz, *onnz, i, rstart, rend, M, N;

 15:   PetscFunctionBeginUser;
 16:   PetscCall(PetscInitialize(&argc, &argv, 0, help));
 17:   comm = MPI_COMM_WORLD;
 18:   PetscCall(PetscMalloc2(m, &dnnz, m, &onnz));
 19:   for (i = 0; i < m; i++) {
 20:     dnnz[i] = 1;
 21:     onnz[i] = 1;
 22:   }
 23:   PetscCall(MatCreateAIJ(comm, m, n, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DECIDE, dnnz, PETSC_DECIDE, onnz, &A));
 24:   PetscCall(MatSetFromOptions(A));
 25:   PetscCall(MatSetUp(A));
 26:   PetscCall(PetscFree2(dnnz, onnz));

 28:   /* since the matrix has never been assembled this reset should do nothing */
 29:   PetscCall(MatResetPreallocation(A));

 31:   /* This assembly shrinks memory because we do not insert enough number of values */
 32:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
 33:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));

 35:   /* MatResetPreallocation() restores the memory required by users */
 36:   PetscCall(MatResetPreallocation(A));
 37:   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
 38:   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
 39:   PetscCall(MatGetSize(A, &M, &N));
 40:   for (i = rstart; i < rend; i++) {
 41:     PetscCall(MatSetValue(A, i, i, 2.0, INSERT_VALUES));
 42:     if (rend < N) PetscCall(MatSetValue(A, i, rend, 1.0, INSERT_VALUES));
 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: }

 52: /*TEST

 54:    test:
 55:       suffix: 1

 57:    test:
 58:       suffix: 2
 59:       nsize: 2

 61: TEST*/