xref: /petsc/src/mat/tutorials/ex4.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 
11c4762a1bSJed Brown int main(int argc,char **argv)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   Mat             A;
14c4762a1bSJed Brown   MPI_Comm        comm;
15c4762a1bSJed Brown   PetscInt        n=5,m=5,*dnnz,*onnz,i,rstart,rend,M,N;
16c4762a1bSJed Brown   PetscErrorCode  ierr;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
19c4762a1bSJed Brown   comm = MPI_COMM_WORLD;
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(m,&dnnz,m,&onnz));
21c4762a1bSJed Brown   for (i=0; i<m; i++) {
22c4762a1bSJed Brown     dnnz[i] = 1;
23c4762a1bSJed Brown     onnz[i] = 1;
24c4762a1bSJed Brown   }
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateAIJ(comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_DECIDE,dnnz,PETSC_DECIDE,onnz,&A));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(dnnz,onnz));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   /* This assembly shrinks memory because we do not insert enough number of values */
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /* MatResetPreallocation restores the memory required by users */
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatResetPreallocation(A));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&M,&N));
39c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
40*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValue(A,i,i,2.0,INSERT_VALUES));
41c4762a1bSJed Brown     if (rend<N) {
42*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValue(A,i,rend,1.0,INSERT_VALUES));
43c4762a1bSJed Brown     }
44c4762a1bSJed Brown   }
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
49c4762a1bSJed Brown   ierr = PetscFinalize();
50c4762a1bSJed Brown   return ierr;
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown /*TEST
54c4762a1bSJed Brown 
55c4762a1bSJed Brown    test:
56c4762a1bSJed Brown       suffix: 1
57c4762a1bSJed Brown 
58c4762a1bSJed Brown    test:
59c4762a1bSJed Brown       suffix: 2
60c4762a1bSJed Brown       nsize: 2
61c4762a1bSJed Brown 
62c4762a1bSJed Brown TEST*/
63