xref: /petsc/src/mat/tests/ex8.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 static char help[] = "Tests automatic allocation of matrix storage space.\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            C;
9   PetscInt       i,j,m = 3,n = 3,Ii,J;
10   PetscScalar    v;
11   MatInfo        info;
12 
13   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
14   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
15   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
16 
17   /* create the matrix for the five point stencil, YET AGAIN */
18   PetscCall(MatCreate(PETSC_COMM_SELF,&C));
19   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
20   PetscCall(MatSetFromOptions(C));
21   PetscCall(MatSetUp(C));
22   for (i=0; i<m; i++) {
23     for (j=0; j<n; j++) {
24       v = -1.0;  Ii = j + n*i;
25       if (i>0)   {J=Ii-n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
26       if (i<m-1) {J=Ii+n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
27       if (j>0)   {J=Ii-1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
28       if (j<n-1) {J=Ii+1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
29       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
30     }
31   }
32   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
33   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
34   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF));
35 
36   PetscCall(MatGetInfo(C,MAT_LOCAL,&info));
37   PetscCall(PetscPrintf(PETSC_COMM_SELF,"matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated));
38 
39   PetscCall(MatDestroy(&C));
40   PetscCall(PetscFinalize());
41   return 0;
42 }
43 
44 /*TEST
45 
46    test:
47 
48 TEST*/
49