Actual source code: ex3.c
petsc-3.7.1 2016-05-15
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";
7: /* Addendum: piggy-backing on this example to test KSPChebyshev methods */
9: #include <petscksp.h>
13: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
14: {
16: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
17: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
18: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
19: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
20: return(0);
21: }
24: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
25: {
27: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
28: return(0);
29: }
33: int main(int argc,char **args)
34: {
35: Mat C;
36: PetscMPIInt rank,size;
37: PetscInt i,m = 5,N,start,end,M,its;
38: PetscScalar val,Ke[16],r[4];
39: PetscReal x,y,h,norm;
41: PetscInt idx[4],count,*rows;
42: Vec u,ustar,b;
43: KSP ksp;
44: PetscBool viewkspest = PETSC_FALSE;
46: PetscInitialize(&argc,&args,(char*)0,help);
47: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
48: PetscOptionsGetBool(NULL,NULL,"-ksp_est_view",&viewkspest,NULL);
49: N = (m+1)*(m+1); /* dimension of matrix */
50: M = m*m; /* number of elements */
51: h = 1.0/m; /* mesh width */
52: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
53: MPI_Comm_size(PETSC_COMM_WORLD,&size);
55: /* Create stiffness matrix */
56: MatCreate(PETSC_COMM_WORLD,&C);
57: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
58: MatSetFromOptions(C);
59: MatSetUp(C);
60: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
61: end = start + M/size + ((M%size) > rank);
63: /* Assemble matrix */
64: FormElementStiffness(h*h,Ke); /* element stiffness for Laplacian */
65: for (i=start; i<end; i++) {
66: /* location of lower left corner of element */
67: x = h*(i % m); y = h*(i/m);
68: /* node numbers for the four corners of element */
69: idx[0] = (m+1)*(i/m) + (i % m);
70: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
71: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
72: }
73: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
74: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
76: /* Create right-hand-side and solution vectors */
77: VecCreate(PETSC_COMM_WORLD,&u);
78: VecSetSizes(u,PETSC_DECIDE,N);
79: VecSetFromOptions(u);
80: PetscObjectSetName((PetscObject)u,"Approx. Solution");
81: VecDuplicate(u,&b);
82: PetscObjectSetName((PetscObject)b,"Right hand side");
83: VecDuplicate(b,&ustar);
84: VecSet(u,0.0);
85: VecSet(b,0.0);
87: /* Assemble right-hand-side vector */
88: for (i=start; i<end; i++) {
89: /* location of lower left corner of element */
90: x = h*(i % m); y = h*(i/m);
91: /* node numbers for the four corners of element */
92: idx[0] = (m+1)*(i/m) + (i % m);
93: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
94: FormElementRhs(x,y,h*h,r);
95: VecSetValues(b,4,idx,r,ADD_VALUES);
96: }
97: VecAssemblyBegin(b);
98: VecAssemblyEnd(b);
100: /* Modify matrix and right-hand-side for Dirichlet boundary conditions */
101: PetscMalloc1(4*m,&rows);
102: for (i=0; i<m+1; i++) {
103: rows[i] = i; /* bottom */
104: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
105: }
106: count = m+1; /* left side */
107: for (i=m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
109: count = 2*m; /* left side */
110: for (i=2*m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
111: for (i=0; i<4*m; i++) {
112: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
113: val = y;
114: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
115: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
116: }
117: MatZeroRows(C,4*m,rows,1.0,0,0);
119: PetscFree(rows);
120: VecAssemblyBegin(u);
121: VecAssemblyEnd(u);
122: VecAssemblyBegin(b);
123: VecAssemblyEnd(b);
125: { Mat A;
126: MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);
127: MatDestroy(&C);
128: MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&C);
129: MatDestroy(&A);
130: }
132: /* Solve linear system */
133: KSPCreate(PETSC_COMM_WORLD,&ksp);
134: KSPSetOperators(ksp,C,C);
135: KSPSetFromOptions(ksp);
136: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
137: KSPSolve(ksp,b,u);
139: if (viewkspest) {
140: KSP kspest;
142: KSPChebyshevEstEigGetKSP(ksp,&kspest);
143: if (kspest) {KSPView(kspest,PETSC_VIEWER_STDOUT_WORLD);}
144: }
146: /* Check error */
147: VecGetOwnershipRange(ustar,&start,&end);
148: for (i=start; i<end; i++) {
149: x = h*(i % (m+1)); y = h*(i/(m+1));
150: val = y;
151: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
152: }
153: VecAssemblyBegin(ustar);
154: VecAssemblyEnd(ustar);
155: VecAXPY(u,-1.0,ustar);
156: VecNorm(u,NORM_2,&norm);
157: KSPGetIterationNumber(ksp,&its);
158: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %D\n",(double)(norm*h),its);
160: /* Free work space */
161: KSPDestroy(&ksp);
162: VecDestroy(&ustar);
163: VecDestroy(&u);
164: VecDestroy(&b);
165: MatDestroy(&C);
166: PetscFinalize();
167: return 0;
168: }