1
2 static char help[] = "Bilinear elements on the unit square for Laplacian. To test the parallel\n\
3 matrix assembly, the matrix is intentionally laid out across processors\n\
4 differently from the way it is assembled. Input arguments are:\n\
5 -m <size> : problem size\n\n";
6
7 /* Addendum: piggy-backing on this example to test KSPChebyshev methods */
8
9 #include <petscksp.h>
10
FormElementStiffness(PetscReal H,PetscScalar * Ke)11 int FormElementStiffness(PetscReal H,PetscScalar *Ke)
12 {
13 PetscFunctionBeginUser;
14 Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
15 Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
16 Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
17 Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
18 PetscFunctionReturn(0);
19 }
FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar * r)20 int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
21 {
22 PetscFunctionBeginUser;
23 r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
24 PetscFunctionReturn(0);
25 }
26
main(int argc,char ** args)27 int main(int argc,char **args)
28 {
29 Mat C;
30 PetscMPIInt rank,size;
31 PetscInt i,m = 5,N,start,end,M,its;
32 PetscScalar val,Ke[16],r[4];
33 PetscReal x,y,h,norm;
34 PetscErrorCode ierr;
35 PetscInt idx[4],count,*rows;
36 Vec u,ustar,b;
37 KSP ksp;
38 PetscBool viewkspest = PETSC_FALSE;
39
40 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
41 ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
42 ierr = PetscOptionsGetBool(NULL,NULL,"-ksp_est_view",&viewkspest,NULL);CHKERRQ(ierr);
43 N = (m+1)*(m+1); /* dimension of matrix */
44 M = m*m; /* number of elements */
45 h = 1.0/m; /* mesh width */
46 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
47 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
48
49 /* Create stiffness matrix */
50 ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
51 ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
52 ierr = MatSetFromOptions(C);CHKERRQ(ierr);
53 ierr = MatSetUp(C);CHKERRQ(ierr);
54 start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
55 end = start + M/size + ((M%size) > rank);
56
57 /* Assemble matrix */
58 ierr = FormElementStiffness(h*h,Ke);CHKERRQ(ierr); /* element stiffness for Laplacian */
59 for (i=start; i<end; i++) {
60 /* node numbers for the four corners of element */
61 idx[0] = (m+1)*(i/m) + (i % m);
62 idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
63 ierr = MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);CHKERRQ(ierr);
64 }
65 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
66 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67
68 /* Create right-hand-side and solution vectors */
69 ierr = VecCreate(PETSC_COMM_WORLD,&u);CHKERRQ(ierr);
70 ierr = VecSetSizes(u,PETSC_DECIDE,N);CHKERRQ(ierr);
71 ierr = VecSetFromOptions(u);CHKERRQ(ierr);
72 ierr = PetscObjectSetName((PetscObject)u,"Approx. Solution");CHKERRQ(ierr);
73 ierr = VecDuplicate(u,&b);CHKERRQ(ierr);
74 ierr = PetscObjectSetName((PetscObject)b,"Right hand side");CHKERRQ(ierr);
75 ierr = VecDuplicate(b,&ustar);CHKERRQ(ierr);
76 ierr = VecSet(u,0.0);CHKERRQ(ierr);
77 ierr = VecSet(b,0.0);CHKERRQ(ierr);
78
79 /* Assemble right-hand-side vector */
80 for (i=start; i<end; i++) {
81 /* location of lower left corner of element */
82 x = h*(i % m); y = h*(i/m);
83 /* node numbers for the four corners of element */
84 idx[0] = (m+1)*(i/m) + (i % m);
85 idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
86 ierr = FormElementRhs(x,y,h*h,r);CHKERRQ(ierr);
87 ierr = VecSetValues(b,4,idx,r,ADD_VALUES);CHKERRQ(ierr);
88 }
89 ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
90 ierr = VecAssemblyEnd(b);CHKERRQ(ierr);
91
92 /* Modify matrix and right-hand-side for Dirichlet boundary conditions */
93 ierr = PetscMalloc1(4*m,&rows);CHKERRQ(ierr);
94 for (i=0; i<m+1; i++) {
95 rows[i] = i; /* bottom */
96 rows[3*m - 1 +i] = m*(m+1) + i; /* top */
97 }
98 count = m+1; /* left side */
99 for (i=m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
100
101 count = 2*m; /* left side */
102 for (i=2*m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
103 for (i=0; i<4*m; i++) {
104 val = h*(rows[i]/(m+1));
105 ierr = VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);CHKERRQ(ierr);
106 ierr = VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);CHKERRQ(ierr);
107 }
108 ierr = MatZeroRows(C,4*m,rows,1.0,0,0);CHKERRQ(ierr);
109
110 ierr = PetscFree(rows);CHKERRQ(ierr);
111 ierr = VecAssemblyBegin(u);CHKERRQ(ierr);
112 ierr = VecAssemblyEnd(u);CHKERRQ(ierr);
113 ierr = VecAssemblyBegin(b);CHKERRQ(ierr);
114 ierr = VecAssemblyEnd(b);CHKERRQ(ierr);
115
116 { Mat A;
117 ierr = MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
118 ierr = MatDestroy(&C);CHKERRQ(ierr);
119 ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr);
120 ierr = MatDestroy(&A);CHKERRQ(ierr);
121 }
122
123 /* Solve linear system */
124 ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);
125 ierr = KSPSetOperators(ksp,C,C);CHKERRQ(ierr);
126 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
127 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
128 ierr = KSPSolve(ksp,b,u);CHKERRQ(ierr);
129
130 if (viewkspest) {
131 KSP kspest;
132
133 ierr = KSPChebyshevEstEigGetKSP(ksp,&kspest);CHKERRQ(ierr);
134 if (kspest) {ierr = KSPView(kspest,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
135 }
136
137 /* Check error */
138 ierr = VecGetOwnershipRange(ustar,&start,&end);CHKERRQ(ierr);
139 for (i=start; i<end; i++) {
140 val = h*(i/(m+1));
141 ierr = VecSetValues(ustar,1,&i,&val,INSERT_VALUES);CHKERRQ(ierr);
142 }
143 ierr = VecAssemblyBegin(ustar);CHKERRQ(ierr);
144 ierr = VecAssemblyEnd(ustar);CHKERRQ(ierr);
145 ierr = VecAXPY(u,-1.0,ustar);CHKERRQ(ierr);
146 ierr = VecNorm(u,NORM_2,&norm);CHKERRQ(ierr);
147 ierr = KSPGetIterationNumber(ksp,&its);CHKERRQ(ierr);
148 ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)(norm*h),its);CHKERRQ(ierr);
149
150 /* Free work space */
151 ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
152 ierr = VecDestroy(&ustar);CHKERRQ(ierr);
153 ierr = VecDestroy(&u);CHKERRQ(ierr);
154 ierr = VecDestroy(&b);CHKERRQ(ierr);
155 ierr = MatDestroy(&C);CHKERRQ(ierr);
156 ierr = PetscFinalize();
157 return ierr;
158 }
159
160 /*TEST
161
162 test:
163 args: -pc_type jacobi -ksp_monitor_short -m 5 -ksp_gmres_cgs_refinement_type refine_always
164
165 test:
166 suffix: 2
167 nsize: 2
168 args: -pc_type jacobi -ksp_monitor_short -m 5 -ksp_gmres_cgs_refinement_type refine_always
169
170 test:
171 suffix: 2_kokkos
172 nsize: 2
173 args: -pc_type jacobi -ksp_monitor_short -m 5 -ksp_gmres_cgs_refinement_type refine_always -mat_type aijkokkos -vec_type kokkos
174 output_file: output/ex3_2.out
175 requires: kokkos_kernels
176
177 test:
178 suffix: nocheby
179 args: -ksp_est_view
180
181 test:
182 suffix: chebynoest
183 args: -ksp_est_view -ksp_type chebyshev -ksp_chebyshev_eigenvalues 0.1,1.0
184
185 test:
186 suffix: chebyest
187 args: -ksp_est_view -ksp_type chebyshev -ksp_chebyshev_esteig
188 filter: sed -e "s/Iterations 19/Iterations 20/g"
189
190 TEST*/
191
192