Actual source code: ex4.c
2: static char help[] = "Solves a linear system with KSP. The matrix uses simple\n\
3: bilinear elements on the unit square. Input arguments are:\n\
4: -m <size> : problem size\n\n";
6: #include petscksp.h
10: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
11: {
12: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
13: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
14: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
15: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
16: return 0;
17: }
20: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
21: {
22: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
23: return 0;
24: }
28: int main(int argc,char **args)
29: {
30: Mat C;
32: PetscInt i,m = 2,N,M,its,idx[4],count,*rows;
33: PetscScalar val,Ke[16],r[4];
34: PetscReal x,y,h,norm;
35: Vec u,ustar,b;
36: KSP ksp;
38: PetscInitialize(&argc,&args,(char *)0,help);
39: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
40: N = (m+1)*(m+1); /* dimension of matrix */
41: M = m*m; /* number of elements */
42: h = 1.0/m; /* mesh width */
44: /* create stiffness matrix */
45: MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,9,PETSC_NULL,&C);
47: /* forms the element stiffness for the Laplacian */
48: FormElementStiffness(h*h,Ke);
49: for (i=0; i<M; i++) {
50: /* location of lower left corner of element */
51: x = h*(i % m); y = h*(i/m);
52: /* node numbers for the four corners of element */
53: idx[0] = (m+1)*(i/m) + (i % m);
54: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
55: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
56: }
57: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
58: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
60: /* create right hand side and solution */
62: VecCreateSeq(PETSC_COMM_SELF,N,&u);
63: VecDuplicate(u,&b);
64: VecDuplicate(b,&ustar);
65: VecSet(u,0.0);
66: VecSet(b,0.0);
68: for (i=0; i<M; i++) {
69: /* location of lower left corner of element */
70: x = h*(i % m); y = h*(i/m);
71: /* node numbers for the four corners of element */
72: idx[0] = (m+1)*(i/m) + (i % m);
73: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
74: FormElementRhs(x,y,h*h,r);
75: VecSetValues(b,4,idx,r,ADD_VALUES);
76: }
77: VecAssemblyBegin(b);
78: VecAssemblyEnd(b);
80: /* modify matrix and rhs for Dirichlet boundary conditions */
81: PetscMalloc((4*m+1)*sizeof(PetscInt),&rows);
82: for (i=0; i<m+1; i++) {
83: rows[i] = i; /* bottom */
84: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
85: }
86: count = m+1; /* left side */
87: for (i=m+1; i<m*(m+1); i+= m+1) {
88: rows[count++] = i;
89: }
90: count = 2*m; /* left side */
91: for (i=2*m+1; i<m*(m+1); i+= m+1) {
92: rows[count++] = i;
93: }
94: for (i=0; i<4*m; i++) {
95: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
96: val = y;
97: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
98: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
99: }
100: MatZeroRows(C,4*m,rows,1.0);
102: PetscFree(rows);
103: VecAssemblyBegin(u);
104: VecAssemblyEnd(u);
105: VecAssemblyBegin(b);
106: VecAssemblyEnd(b);
108: /* solve linear system */
109: KSPCreate(PETSC_COMM_WORLD,&ksp);
110: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
111: KSPSetFromOptions(ksp);
112: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
113: KSPSolve(ksp,b,u);
115: /* check error */
116: for (i=0; i<N; i++) {
117: x = h*(i % (m+1)); y = h*(i/(m+1));
118: val = y;
119: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
120: }
121: VecAssemblyBegin(ustar);
122: VecAssemblyEnd(ustar);
124: VecAXPY(u,-1.0,ustar);
125: VecNorm(u,NORM_2,&norm);
126: KSPGetIterationNumber(ksp,&its);
127: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %D\n",norm*h,its);
129: KSPDestroy(ksp);
130: VecDestroy(ustar);
131: VecDestroy(u);
132: VecDestroy(b);
133: MatDestroy(C);
134: PetscFinalize();
135: return 0;
136: }