Actual source code: ex4.c
petsc-3.7.1 2016-05-15
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>
9: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
10: {
11: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
12: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
13: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
14: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
15: return 0;
16: }
19: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
20: {
21: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
22: return 0;
23: }
27: int main(int argc,char **args)
28: {
29: Mat C;
31: PetscInt i,m = 2,N,M,its,idx[4],count,*rows;
32: PetscScalar val,Ke[16],r[4];
33: PetscReal x,y,h,norm;
34: Vec u,ustar,b;
35: KSP ksp;
37: PetscInitialize(&argc,&args,(char*)0,help);
38: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
39: N = (m+1)*(m+1); /* dimension of matrix */
40: M = m*m; /* number of elements */
41: h = 1.0/m; /* mesh width */
43: /* create stiffness matrix */
44: MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,9,NULL,&C);
45: MatSetUp(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: PetscMalloc1(4*m+1,&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) rows[count++] = i;
89: count = 2*m; /* left side */
90: for (i=2*m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
92: for (i=0; i<4*m; i++) {
93: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
94: val = y;
95: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
96: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
97: }
98: MatZeroRows(C,4*m,rows,1.0,0,0);
100: PetscFree(rows);
101: VecAssemblyBegin(u);
102: VecAssemblyEnd(u);
103: VecAssemblyBegin(b);
104: VecAssemblyEnd(b);
106: /* solve linear system */
107: KSPCreate(PETSC_COMM_WORLD,&ksp);
108: KSPSetOperators(ksp,C,C);
109: KSPSetFromOptions(ksp);
110: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
111: KSPSolve(ksp,b,u);
113: /* check error */
114: for (i=0; i<N; i++) {
115: x = h*(i % (m+1)); y = h*(i/(m+1));
116: val = y;
117: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
118: }
119: VecAssemblyBegin(ustar);
120: VecAssemblyEnd(ustar);
122: VecAXPY(u,-1.0,ustar);
123: VecNorm(u,NORM_2,&norm);
124: KSPGetIterationNumber(ksp,&its);
125: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)(norm*h),its);
127: KSPDestroy(&ksp);
128: VecDestroy(&ustar);
129: VecDestroy(&u);
130: VecDestroy(&b);
131: MatDestroy(&C);
132: PetscFinalize();
133: return 0;
134: }