Actual source code: ex4.c
petsc-3.3-p6 2013-02-11
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,tol=1.e-14;
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);
46: MatSetUp(C);
48: /* forms the element stiffness for the Laplacian */
49: FormElementStiffness(h*h,Ke);
50: for (i=0; i<M; i++) {
51: /* location of lower left corner of element */
52: x = h*(i % m); y = h*(i/m);
53: /* node numbers for the four corners of element */
54: idx[0] = (m+1)*(i/m) + (i % m);
55: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
56: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
57: }
58: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
59: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
61: /* create right hand side and solution */
63: VecCreateSeq(PETSC_COMM_SELF,N,&u);
64: VecDuplicate(u,&b);
65: VecDuplicate(b,&ustar);
66: VecSet(u,0.0);
67: VecSet(b,0.0);
69: for (i=0; i<M; i++) {
70: /* location of lower left corner of element */
71: x = h*(i % m); y = h*(i/m);
72: /* node numbers for the four corners of element */
73: idx[0] = (m+1)*(i/m) + (i % m);
74: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
75: FormElementRhs(x,y,h*h,r);
76: VecSetValues(b,4,idx,r,ADD_VALUES);
77: }
78: VecAssemblyBegin(b);
79: VecAssemblyEnd(b);
81: /* modify matrix and rhs for Dirichlet boundary conditions */
82: PetscMalloc((4*m+1)*sizeof(PetscInt),&rows);
83: for (i=0; i<m+1; i++) {
84: rows[i] = i; /* bottom */
85: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
86: }
87: count = m+1; /* left side */
88: for (i=m+1; i<m*(m+1); i+= m+1) {
89: rows[count++] = i;
90: }
91: count = 2*m; /* left side */
92: for (i=2*m+1; i<m*(m+1); i+= m+1) {
93: rows[count++] = i;
94: }
95: for (i=0; i<4*m; i++) {
96: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
97: val = y;
98: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
99: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
100: }
101: MatZeroRows(C,4*m,rows,1.0,0,0);
103: PetscFree(rows);
104: VecAssemblyBegin(u);
105: VecAssemblyEnd(u);
106: VecAssemblyBegin(b);
107: VecAssemblyEnd(b);
109: /* solve linear system */
110: KSPCreate(PETSC_COMM_WORLD,&ksp);
111: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
112: KSPSetFromOptions(ksp);
113: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
114: KSPSolve(ksp,b,u);
116: /* check error */
117: for (i=0; i<N; i++) {
118: x = h*(i % (m+1)); y = h*(i/(m+1));
119: val = y;
120: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
121: }
122: VecAssemblyBegin(ustar);
123: VecAssemblyEnd(ustar);
125: VecAXPY(u,-1.0,ustar);
126: VecNorm(u,NORM_2,&norm);
127: KSPGetIterationNumber(ksp,&its);
128: if (norm > tol){
129: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %G Iterations %D\n",norm*h,its);
130: }
132: KSPDestroy(&ksp);
133: VecDestroy(&ustar);
134: VecDestroy(&u);
135: VecDestroy(&b);
136: MatDestroy(&C);
137: PetscFinalize();
138: return 0;
139: }