xref: /petsc/src/mat/tests/ex10.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests repeated use of assembly for matrices.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            C;
9c4762a1bSJed Brown   PetscInt       i,j,m = 5,n = 2,Ii,J;
10c4762a1bSJed Brown   PetscErrorCode ierr;
11c4762a1bSJed Brown   PetscMPIInt    rank,size;
12c4762a1bSJed Brown   PetscScalar    v;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
17c4762a1bSJed Brown   n    = 2*size;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
24c4762a1bSJed Brown   for (i=0; i<m; i++) {
25c4762a1bSJed Brown     for (j=2*rank; j<2*rank+2; j++) {
26c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
27*5f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
28*5f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
29*5f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
30*5f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
31*5f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
32c4762a1bSJed Brown     }
33c4762a1bSJed Brown   }
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
36c4762a1bSJed Brown   for (i=0; i<m; i++) {
37c4762a1bSJed Brown     for (j=2*rank; j<2*rank+2; j++) {
38c4762a1bSJed Brown       v = 1.0;  Ii = j + n*i;
39*5f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
40*5f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
41*5f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
42*5f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
43*5f80ce2aSJacob Faibussowitsch       v = -4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
44c4762a1bSJed Brown     }
45c4762a1bSJed Brown   }
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
48c4762a1bSJed Brown 
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
50c4762a1bSJed Brown 
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
52c4762a1bSJed Brown   ierr = PetscFinalize();
53c4762a1bSJed Brown   return ierr;
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*TEST
57c4762a1bSJed Brown 
58c4762a1bSJed Brown    test:
59c4762a1bSJed Brown 
60c4762a1bSJed Brown TEST*/
61