xref: /petsc/src/mat/tests/ex8.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
15   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
16   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
17 
18   /* create the matrix for the five point stencil, YET AGAIN */
19   PetscCall(MatCreate(PETSC_COMM_SELF,&C));
20   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
21   PetscCall(MatSetFromOptions(C));
22   PetscCall(MatSetUp(C));
23   for (i=0; i<m; i++) {
24     for (j=0; j<n; j++) {
25       v = -1.0;  Ii = j + n*i;
26       if (i>0)   {J=Ii-n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
27       if (i<m-1) {J=Ii+n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
28       if (j>0)   {J=Ii-1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
29       if (j<n-1) {J=Ii+1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
30       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
31     }
32   }
33   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
34   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
35   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF));
36 
37   PetscCall(MatGetInfo(C,MAT_LOCAL,&info));
38   PetscCall(PetscPrintf(PETSC_COMM_SELF,"matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated));
39 
40   PetscCall(MatDestroy(&C));
41   PetscCall(PetscFinalize());
42   return 0;
43 }
44 
45 /*TEST
46 
47    test:
48 
49 TEST*/
50