Actual source code: ex3.c
petsc-3.4.2 2013-07-02
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: #include <petscksp.h>
11: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
12: {
14: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
15: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
16: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
17: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
18: return(0);
19: }
22: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
23: {
25: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
26: return(0);
27: }
31: int main(int argc,char **args)
32: {
33: Mat C;
34: PetscMPIInt rank,size;
35: PetscInt i,m = 5,N,start,end,M,its;
36: PetscScalar val,Ke[16],r[4];
37: PetscReal x,y,h,norm,tol=1.e-14;
39: PetscInt idx[4],count,*rows;
40: Vec u,ustar,b;
41: KSP ksp;
43: PetscInitialize(&argc,&args,(char*)0,help);
44: PetscOptionsGetInt(NULL,"-m",&m,NULL);
45: N = (m+1)*(m+1); /* dimension of matrix */
46: M = m*m; /* number of elements */
47: h = 1.0/m; /* mesh width */
48: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
49: MPI_Comm_size(PETSC_COMM_WORLD,&size);
51: /* Create stiffness matrix */
52: MatCreate(PETSC_COMM_WORLD,&C);
53: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
54: MatSetFromOptions(C);
55: MatSetUp(C);
56: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
57: end = start + M/size + ((M%size) > rank);
59: /* Assemble matrix */
60: FormElementStiffness(h*h,Ke); /* element stiffness for Laplacian */
61: for (i=start; i<end; i++) {
62: /* location of lower left corner of element */
63: x = h*(i % m); y = h*(i/m);
64: /* node numbers for the four corners of element */
65: idx[0] = (m+1)*(i/m) + (i % m);
66: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
67: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
68: }
69: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
70: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
72: /* Create right-hand-side and solution vectors */
73: VecCreate(PETSC_COMM_WORLD,&u);
74: VecSetSizes(u,PETSC_DECIDE,N);
75: VecSetFromOptions(u);
76: PetscObjectSetName((PetscObject)u,"Approx. Solution");
77: VecDuplicate(u,&b);
78: PetscObjectSetName((PetscObject)b,"Right hand side");
79: VecDuplicate(b,&ustar);
80: VecSet(u,0.0);
81: VecSet(b,0.0);
83: /* Assemble right-hand-side vector */
84: for (i=start; i<end; i++) {
85: /* location of lower left corner of element */
86: x = h*(i % m); y = h*(i/m);
87: /* node numbers for the four corners of element */
88: idx[0] = (m+1)*(i/m) + (i % m);
89: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
90: FormElementRhs(x,y,h*h,r);
91: VecSetValues(b,4,idx,r,ADD_VALUES);
92: }
93: VecAssemblyBegin(b);
94: VecAssemblyEnd(b);
96: /* Modify matrix and right-hand-side for Dirichlet boundary conditions */
97: PetscMalloc(4*m*sizeof(PetscInt),&rows);
98: for (i=0; i<m+1; i++) {
99: rows[i] = i; /* bottom */
100: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
101: }
102: count = m+1; /* left side */
103: for (i=m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
105: count = 2*m; /* left side */
106: for (i=2*m+1; i<m*(m+1); i+= m+1) rows[count++] = i;
107: for (i=0; i<4*m; i++) {
108: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
109: val = y;
110: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
111: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
112: }
113: MatZeroRows(C,4*m,rows,1.0,0,0);
115: PetscFree(rows);
116: VecAssemblyBegin(u);
117: VecAssemblyEnd(u);
118: VecAssemblyBegin(b);
119: VecAssemblyEnd(b);
121: { Mat A;
122: MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);
123: MatDestroy(&C);
124: MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&C);
125: MatDestroy(&A);
126: }
128: /* Solve linear system */
129: KSPCreate(PETSC_COMM_WORLD,&ksp);
130: KSPSetOperators(ksp,C,C,DIFFERENT_NONZERO_PATTERN);
131: KSPSetFromOptions(ksp);
132: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
133: KSPSolve(ksp,b,u);
135: /* Check error */
136: VecGetOwnershipRange(ustar,&start,&end);
137: for (i=start; i<end; i++) {
138: x = h*(i % (m+1)); y = h*(i/(m+1));
139: val = y;
140: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
141: }
142: VecAssemblyBegin(ustar);
143: VecAssemblyEnd(ustar);
144: VecAXPY(u,-1.0,ustar);
145: VecNorm(u,NORM_2,&norm);
146: KSPGetIterationNumber(ksp,&its);
147: if (norm > tol) {
148: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %G Iterations %D\n",norm*h,its);
149: }
151: /* Free work space */
152: KSPDestroy(&ksp);
153: VecDestroy(&ustar);
154: VecDestroy(&u);
155: VecDestroy(&b);
156: MatDestroy(&C);
157: PetscFinalize();
158: return 0;
159: }