1c4762a1bSJed Brown /*
2*46233b44SBarry Smith Contributed by: Fande Kong
3c4762a1bSJed Brown */
4c4762a1bSJed Brown
5c4762a1bSJed Brown static char help[] = "Illustrate the use of MatResetPreallocation.\n";
6c4762a1bSJed Brown
7*46233b44SBarry Smith #include <petscmat.h>
8c4762a1bSJed Brown
main(int argc,char ** argv)9d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
10d71ae5a4SJacob Faibussowitsch {
11c4762a1bSJed Brown Mat A;
12c4762a1bSJed Brown MPI_Comm comm;
13c4762a1bSJed Brown PetscInt n = 5, m = 5, *dnnz, *onnz, i, rstart, rend, M, N;
14c4762a1bSJed Brown
15327415f7SBarry Smith PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help));
17c4762a1bSJed Brown comm = MPI_COMM_WORLD;
189566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(m, &dnnz, m, &onnz));
19c4762a1bSJed Brown for (i = 0; i < m; i++) {
20c4762a1bSJed Brown dnnz[i] = 1;
21c4762a1bSJed Brown onnz[i] = 1;
22c4762a1bSJed Brown }
239566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(comm, m, n, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DECIDE, dnnz, PETSC_DECIDE, onnz, &A));
249566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
259566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
269566063dSJacob Faibussowitsch PetscCall(PetscFree2(dnnz, onnz));
27c4762a1bSJed Brown
281f14be2bSBarry Smith /* since the matrix has never been assembled this reset should do nothing */
291f14be2bSBarry Smith PetscCall(MatResetPreallocation(A));
301f14be2bSBarry Smith
31c4762a1bSJed Brown /* This assembly shrinks memory because we do not insert enough number of values */
329566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
339566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
34c4762a1bSJed Brown
351f14be2bSBarry Smith /* MatResetPreallocation() restores the memory required by users */
369566063dSJacob Faibussowitsch PetscCall(MatResetPreallocation(A));
379566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
389566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
399566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N));
40c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
419566063dSJacob Faibussowitsch PetscCall(MatSetValue(A, i, i, 2.0, INSERT_VALUES));
4248a46eb9SPierre Jolivet if (rend < N) PetscCall(MatSetValue(A, i, rend, 1.0, INSERT_VALUES));
43c4762a1bSJed Brown }
449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
469566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
489566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
49b122ec5aSJacob Faibussowitsch return 0;
50c4762a1bSJed Brown }
51c4762a1bSJed Brown
52c4762a1bSJed Brown /*TEST
53c4762a1bSJed Brown
54c4762a1bSJed Brown test:
55c4762a1bSJed Brown suffix: 1
56c4762a1bSJed Brown
57c4762a1bSJed Brown test:
58c4762a1bSJed Brown suffix: 2
59c4762a1bSJed Brown nsize: 2
60c4762a1bSJed Brown
61c4762a1bSJed Brown TEST*/
62