Actual source code: ex3.c
1: /*$Id: ex3.c,v 1.73 2001/08/07 21:30:50 bsmith Exp $*/
3: static char help[] = "This example solves a linear system in parallel with SLES. The matrixn
4: uses simple bilinear elements on the unit square. To test the paralleln
5: matrix assembly, the matrix is intentionally laid out across processorsn
6: differently from the way it is assembled. Input arguments are:n
7: -m <size> : problem sizenn";
9: #include petscsles.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: }
20: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
21: {
23: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
24: return(0);
25: }
27: int main(int argc,char **args)
28: {
29: Mat C;
30: int i,m = 5,rank,size,N,start,end,M,its;
31: PetscScalar val,zero = 0.0,one = 1.0,none = -1.0,Ke[16],r[4];
32: PetscReal x,y,h,norm;
33: int ierr,idx[4],count,*rows;
34: Vec u,ustar,b;
35: SLES sles;
36: KSP ksp;
37: IS is;
39: PetscInitialize(&argc,&args,(char *)0,help);
40: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
41: N = (m+1)*(m+1); /* dimension of matrix */
42: M = m*m; /* number of elements */
43: h = 1.0/m; /* mesh width */
44: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
45: MPI_Comm_size(PETSC_COMM_WORLD,&size);
47: /* Create stiffness matrix */
48: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&C);
49: MatSetFromOptions(C);
50: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
51: end = start + M/size + ((M%size) > rank);
53: /* Assemble matrix */
54: FormElementStiffness(h*h,Ke); /* element stiffness for Laplacian */
55: for (i=start; i<end; i++) {
56: /* location of lower left corner of element */
57: x = h*(i % m); y = h*(i/m);
58: /* node numbers for the four corners of element */
59: idx[0] = (m+1)*(i/m) + (i % m);
60: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
61: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
62: }
63: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
64: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
66: /* Create right-hand-side and solution vectors */
67: VecCreate(PETSC_COMM_WORLD,&u);
68: VecSetSizes(u,PETSC_DECIDE,N);
69: VecSetFromOptions(u);
70: PetscObjectSetName((PetscObject)u,"Approx. Solution");
71: VecDuplicate(u,&b);
72: PetscObjectSetName((PetscObject)b,"Right hand side");
73: VecDuplicate(b,&ustar);
74: VecSet(&zero,u);
75: VecSet(&zero,b);
77: /* Assemble right-hand-side vector */
78: for (i=start; i<end; i++) {
79: /* location of lower left corner of element */
80: x = h*(i % m); y = h*(i/m);
81: /* node numbers for the four corners of element */
82: idx[0] = (m+1)*(i/m) + (i % m);
83: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
84: FormElementRhs(x,y,h*h,r);
85: VecSetValues(b,4,idx,r,ADD_VALUES);
86: }
87: VecAssemblyBegin(b);
88: VecAssemblyEnd(b);
90: /* Modify matrix and right-hand-side for Dirichlet boundary conditions */
91: PetscMalloc(4*m*sizeof(int),&rows);
92: for (i=0; i<m+1; i++) {
93: rows[i] = i; /* bottom */
94: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
95: }
96: count = m+1; /* left side */
97: for (i=m+1; i<m*(m+1); i+= m+1) {
98: rows[count++] = i;
99: }
100: count = 2*m; /* left side */
101: for (i=2*m+1; i<m*(m+1); i+= m+1) {
102: rows[count++] = i;
103: }
104: ISCreateGeneral(PETSC_COMM_SELF,4*m,rows,&is);
105: for (i=0; i<4*m; i++) {
106: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
107: val = y;
108: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
109: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
110: }
111: PetscFree(rows);
112: VecAssemblyBegin(u);
113: VecAssemblyEnd(u);
114: VecAssemblyBegin(b);
115: VecAssemblyEnd(b);
117: MatZeroRows(C,is,&one);
118: ISDestroy(is);
121: { Mat A;
122: MatConvert(C,MATSAME,&A);
123: MatDestroy(C);
124: MatConvert(A,MATSAME,&C);
125: MatDestroy(A);
126: }
128: /* Solve linear system */
129: SLESCreate(PETSC_COMM_WORLD,&sles);
130: SLESSetOperators(sles,C,C,DIFFERENT_NONZERO_PATTERN);
131: SLESSetFromOptions(sles);
132: SLESGetKSP(sles,&ksp);
133: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
134: SLESSolve(sles,b,u,&its);
136: /* Check error */
137: VecGetOwnershipRange(ustar,&start,&end);
138: for (i=start; i<end; i++) {
139: x = h*(i % (m+1)); y = h*(i/(m+1));
140: val = y;
141: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
142: }
143: VecAssemblyBegin(ustar);
144: VecAssemblyEnd(ustar);
145: VecAXPY(&none,ustar,u);
146: VecNorm(u,NORM_2,&norm);
147: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %dn",norm*h,its);
149: /* Free work space */
150: SLESDestroy(sles);
151: VecDestroy(ustar);
152: VecDestroy(u);
153: VecDestroy(b);
154: MatDestroy(C);
155: PetscFinalize();
156: return 0;
157: }