Actual source code: ex27.c
petsc-3.8.4 2018-03-24
2: static char help[] = "Tests repeated use of assembly for matrices.\n\
3: does nasty case where matrix must be rebuilt.\n\n";
5: #include <petscmat.h>
7: int main(int argc,char **args)
8: {
9: Mat C;
10: PetscInt i,j,m = 5,n = 2,Ii,J;
12: PetscMPIInt rank,size;
13: PetscScalar v;
14: Vec x,y;
16: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
17: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
18: MPI_Comm_size(PETSC_COMM_WORLD,&size);
19: n = 2*size;
21: /* Create the matrix for the five point stencil, YET AGAIN */
22: MatCreate(PETSC_COMM_WORLD,&C);
23: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
24: MatSetFromOptions(C);
25: for (i=0; i<m; i++) {
26: for (j=2*rank; j<2*rank+2; j++) {
27: v = -1.0; Ii = j + n*i;
28: if (i>0) {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
29: if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
30: if (j>0) {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
31: if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
32: v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
33: }
34: }
35: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
36: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
37: for (i=0; i<m; i++) {
38: for (j=2*rank; j<2*rank+2; j++) {
39: v = 1.0; Ii = j + n*i;
40: if (i>0) {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
41: if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
42: if (j>0) {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
43: if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
44: v = -4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
45: }
46: }
47: /* Introduce new nonzero that requires new construction for
48: matrix-vector product */
49: if (rank) {
50: Ii = rank-1; J = m*n-1;
51: MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);
52: }
53: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
54: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
56: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
58: /* Form a couple of vectors to test matrix-vector product */
59: VecCreate(PETSC_COMM_WORLD,&x);
60: VecSetSizes(x,PETSC_DECIDE,m*n);
61: VecSetFromOptions(x);
62: VecDuplicate(x,&y);
63: v = 1.0; VecSet(x,v);
64: MatMult(C,x,y);
66: MatDestroy(&C);
67: VecDestroy(&x);
68: VecDestroy(&y);
69: PetscFinalize();
70: return ierr;
71: }