xref: /petsc/src/mat/tutorials/ex4.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown  *
3c4762a1bSJed Brown  *  Created on: Sep 25, 2017
4c4762a1bSJed Brown  *      Author: Fande Kong
5c4762a1bSJed Brown  */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown static char help[] = "Illustrate the use of MatResetPreallocation.\n";
8c4762a1bSJed Brown 
9c4762a1bSJed Brown #include "petscmat.h"
10c4762a1bSJed Brown 
119371c9d4SSatish Balay int main(int argc, char **argv) {
12c4762a1bSJed Brown   Mat      A;
13c4762a1bSJed Brown   MPI_Comm comm;
14c4762a1bSJed Brown   PetscInt n = 5, m = 5, *dnnz, *onnz, i, rstart, rend, M, N;
15c4762a1bSJed Brown 
16327415f7SBarry Smith   PetscFunctionBeginUser;
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, help));
18c4762a1bSJed Brown   comm = MPI_COMM_WORLD;
199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(m, &dnnz, m, &onnz));
20c4762a1bSJed Brown   for (i = 0; i < m; i++) {
21c4762a1bSJed Brown     dnnz[i] = 1;
22c4762a1bSJed Brown     onnz[i] = 1;
23c4762a1bSJed Brown   }
249566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(comm, m, n, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DECIDE, dnnz, PETSC_DECIDE, onnz, &A));
259566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
269566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
279566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnnz, onnz));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* This assembly shrinks memory because we do not insert enough number of values */
309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* MatResetPreallocation restores the memory required by users */
349566063dSJacob Faibussowitsch   PetscCall(MatResetPreallocation(A));
359566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
369566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
379566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
38c4762a1bSJed Brown   for (i = rstart; i < rend; i++) {
399566063dSJacob Faibussowitsch     PetscCall(MatSetValue(A, i, i, 2.0, INSERT_VALUES));
40*48a46eb9SPierre Jolivet     if (rend < N) PetscCall(MatSetValue(A, i, rend, 1.0, INSERT_VALUES));
41c4762a1bSJed Brown   }
429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
449566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
469566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
47b122ec5aSJacob Faibussowitsch   return 0;
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50c4762a1bSJed Brown /*TEST
51c4762a1bSJed Brown 
52c4762a1bSJed Brown    test:
53c4762a1bSJed Brown       suffix: 1
54c4762a1bSJed Brown 
55c4762a1bSJed Brown    test:
56c4762a1bSJed Brown       suffix: 2
57c4762a1bSJed Brown       nsize: 2
58c4762a1bSJed Brown 
59c4762a1bSJed Brown TEST*/
60