Actual source code: ex4.c
petsc-3.13.4 2020-08-01
2: static char help[] = "Bilinear elements on the unit square for the Laplacian. Input arguments are:\n\
3: -m <size> : problem size\n\n";
5: #include <petscksp.h>
7: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
8: {
9: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
10: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
11: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
12: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
13: return 0;
14: }
15: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
16: {
17: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
18: return 0;
19: }
21: int main(int argc,char **args)
22: {
23: Mat C;
25: PetscInt i,m = 2,N,M,its,idx[4],count,*rows;
26: PetscScalar val,Ke[16],r[4];
27: PetscReal x,y,h,norm;
28: Vec u,ustar,b;
29: KSP ksp;
31: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
32: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
33: N = (m+1)*(m+1); /* dimension of matrix */
34: M = m*m; /* number of elements */
35: h = 1.0/m; /* mesh width */
37: /* create stiffness matrix */
38: MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,9,NULL,&C);
39: MatSetUp(C);
41: /* forms the element stiffness for the Laplacian */
42: FormElementStiffness(h*h,Ke);
43: for (i=0; i<M; i++) {
44: /* node numbers for the four corners of element */
45: idx[0] = (m+1)*(i/m) + (i % m);
46: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
47: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
48: }
49: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
50: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
52: /* create right hand side and solution */
54: VecCreateSeq(PETSC_COMM_SELF,N,&u);
55: VecDuplicate(u,&b);
56: VecDuplicate(b,&ustar);
57: VecSet(u,0.0);
58: VecSet(b,0.0);
60: for (i=0; i<M; i++) {
61: /* location of lower left corner of element */
62: x = h*(i % m); y = h*(i/m);
63: /* node numbers for the four corners of element */
64: idx[0] = (m+1)*(i/m) + (i % m);
65: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
66: FormElementRhs(x,y,h*h,r);
67: VecSetValues(b,4,idx,r,ADD_VALUES);
68: }
69: VecAssemblyBegin(b);
70: VecAssemblyEnd(b);
72: /* modify matrix and rhs for Dirichlet boundary conditions */
73: PetscMalloc1(4*m+1,&rows);
74: for (i=0; i<m+1; i++) {
75: rows[i] = i; /* bottom */
76: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
77: }
78: count = m+1; /* left side */
79: for (i=m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
81: count = 2*m; /* left side */
82: for (i=2*m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
84: for (i=0; i<4*m; i++) {
85: val = h*(rows[i]/(m+1));
86: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
87: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
88: }
89: MatZeroRows(C,4*m,rows,1.0,0,0);
91: PetscFree(rows);
92: VecAssemblyBegin(u);
93: VecAssemblyEnd(u);
94: VecAssemblyBegin(b);
95: VecAssemblyEnd(b);
97: /* solve linear system */
98: KSPCreate(PETSC_COMM_WORLD,&ksp);
99: KSPSetOperators(ksp,C,C);
100: KSPSetFromOptions(ksp);
101: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
102: KSPSolve(ksp,b,u);
104: /* check error */
105: for (i=0; i<N; i++) {
106: val = h*(i/(m+1));
107: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
108: }
109: VecAssemblyBegin(ustar);
110: VecAssemblyEnd(ustar);
112: VecAXPY(u,-1.0,ustar);
113: VecNorm(u,NORM_2,&norm);
114: KSPGetIterationNumber(ksp,&its);
115: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)(norm*h),its);
117: KSPDestroy(&ksp);
118: VecDestroy(&ustar);
119: VecDestroy(&u);
120: VecDestroy(&b);
121: MatDestroy(&C);
122: PetscFinalize();
123: return ierr;
124: }
126: /*TEST
128: test:
129: args: -ksp_monitor_short -m 5 -pc_type jacobi -ksp_gmres_cgs_refinement_type refine_always
131: test:
132: suffix: 3
133: args: -pc_type sor -pc_sor_symmetric -ksp_monitor_short -m 5 -ksp_gmres_cgs_refinement_type refine_always
135: test:
136: suffix: 5
137: args: -pc_type eisenstat -ksp_monitor_short -m 5 -ksp_gmres_cgs_refinement_type refine_always
139: TEST*/