xref: /petsc/src/mat/tests/ex10.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 static char help[] = "Tests repeated use of assembly for matrices.\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            C;
9   PetscInt       i,j,m = 5,n = 2,Ii,J;
10   PetscMPIInt    rank,size;
11   PetscScalar    v;
12 
13   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
14   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
15   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
16   n    = 2*size;
17 
18   /* create the matrix for the five point stencil, YET AGAIN*/
19   PetscCall(MatCreate(PETSC_COMM_WORLD,&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=2*rank; j<2*rank+2; 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   for (i=0; i<m; i++) {
36     for (j=2*rank; j<2*rank+2; j++) {
37       v = 1.0;  Ii = j + n*i;
38       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
39       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
40       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
41       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
42       v = -4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
43     }
44   }
45   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
46   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
47 
48   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
49 
50   PetscCall(MatDestroy(&C));
51   PetscCall(PetscFinalize());
52   return 0;
53 }
54 
55 /*TEST
56 
57    test:
58 
59 TEST*/
60