xref: /petsc/src/mat/tests/ex10.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscMPIInt    rank,size;
11c4762a1bSJed Brown   PetscScalar    v;
12c4762a1bSJed Brown 
13*327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
169566063dSJacob Faibussowitsch   PetscCallMPI(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*/
209566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
219566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
229566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
239566063dSJacob Faibussowitsch   PetscCall(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;
279566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
289566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
299566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
309566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
319566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
32c4762a1bSJed Brown     }
33c4762a1bSJed Brown   }
349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
359566063dSJacob Faibussowitsch   PetscCall(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;
399566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
409566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
419566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
429566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
439566063dSJacob Faibussowitsch       v = -4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
44c4762a1bSJed Brown     }
45c4762a1bSJed Brown   }
469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
529566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
53b122ec5aSJacob Faibussowitsch   return 0;
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*TEST
57c4762a1bSJed Brown 
58c4762a1bSJed Brown    test:
59c4762a1bSJed Brown 
60c4762a1bSJed Brown TEST*/
61