Actual source code: ex3.c
1: /*$Id: ex3.c,v 1.29 2001/03/23 23:23:55 balay Exp $*/
3: static char help[] = "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: /*T
10: Concepts: SLES^basic parallel example
11: Concepts: Matrices^inserting elements by blocks
12: Processors: n
13: T*/
15: /*
16: Include "petscsles.h" so that we can use SLES solvers. Note that this file
17: automatically includes:
18: petsc.h - base PETSc routines petscvec.h - vectors
19: petscsys.h - system routines petscmat.h - matrices
20: petscis.h - index sets petscksp.h - Krylov subspace methods
21: petscviewer.h - viewers petscpc.h - preconditioners
22: */
23: #include "petscsles.h"
25: /* Declare user-defined routines */
26: extern int FormElementStiffness(double,Scalar*);
27: extern int FormElementRhs(double,double,double,Scalar*);
29: int main(int argc,char **args)
30: {
31: Vec u,b,ustar; /* approx solution, RHS, exact solution */
32: Mat A; /* linear system matrix */
33: SLES sles; /* linear solver context */
34: KSP ksp; /* Krylov subspace method context */
35: IS is; /* index set - used for boundary conditions */
36: int N; /* dimension of system (global) */
37: int M; /* number of elements (global) */
38: int rank; /* processor rank */
39: int size; /* size of communicator */
40: Scalar Ke[16]; /* element matrix */
41: Scalar r[4]; /* element vector */
42: double h; /* mesh width */
43: double norm; /* norm of solution error */
44: double x,y;
45: Scalar val,zero = 0.0,one = 1.0,none = -1.0;
46: int ierr,idx[4],count,*rows,i,m = 5,start,end,its;
48: PetscInitialize(&argc,&args,(char *)0,help);
49: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
50: N = (m+1)*(m+1);
51: M = m*m;
52: h = 1.0/m;
53: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
54: MPI_Comm_size(PETSC_COMM_WORLD,&size);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Compute the matrix and right-hand-side vector that define
58: the linear system, Au = b.
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: /*
62: Create stiffness matrix
63: */
64: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&A);
65: MatSetFromOptions(A);
66: start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank);
67: end = start + M/size + ((M%size) > rank);
69: /*
70: Assemble matrix
71: */
72: FormElementStiffness(h*h,Ke);
73: for (i=start; i<end; i++) {
74: /* location of lower left corner of element */
75: x = h*(i % m); y = h*(i/m);
76: /* node numbers for the four corners of element */
77: idx[0] = (m+1)*(i/m) + (i % m);
78: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
79: MatSetValues(A,4,idx,4,idx,Ke,ADD_VALUES);
80: }
81: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
82: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
84: /*
85: Create right-hand-side and solution vectors
86: */
87: VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,N,&u);
88: VecSetFromOptions(u);
89: PetscObjectSetName((PetscObject)u,"Approx. Solution");
90: VecDuplicate(u,&b);
91: PetscObjectSetName((PetscObject)b,"Right hand side");
92: VecDuplicate(b,&ustar);
93: VecSet(&zero,u);
94: VecSet(&zero,b);
96: /*
97: Assemble right-hand-side vector
98: */
99: for (i=start; i<end; i++) {
100: /* location of lower left corner of element */
101: x = h*(i % m); y = h*(i/m);
102: /* node numbers for the four corners of element */
103: idx[0] = (m+1)*(i/m) + (i % m);
104: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
105: FormElementRhs(x,y,h*h,r);
106: VecSetValues(b,4,idx,r,ADD_VALUES);
107: }
108: VecAssemblyBegin(b);
109: VecAssemblyEnd(b);
111: /*
112: Modify matrix and right-hand-side for Dirichlet boundary conditions
113: */
114: PetscMalloc(4*m*sizeof(int),&rows);
115: for (i=0; i<m+1; i++) {
116: rows[i] = i; /* bottom */
117: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
118: }
119: count = m+1; /* left side */
120: for (i=m+1; i<m*(m+1); i+= m+1) {
121: rows[count++] = i;
122: }
123: count = 2*m; /* left side */
124: for (i=2*m+1; i<m*(m+1); i+= m+1) {
125: rows[count++] = i;
126: }
127: ISCreateGeneral(PETSC_COMM_SELF,4*m,rows,&is);
128: for (i=0; i<4*m; i++) {
129: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
130: val = y;
131: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
132: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
133: }
134: PetscFree(rows);
135: VecAssemblyBegin(u);
136: VecAssemblyEnd(u);
137: VecAssemblyBegin(b);
138: VecAssemblyEnd(b);
140: MatZeroRows(A,is,&one);
141: ISDestroy(is);
143: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144: Create the linear solver and set various options
145: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147: SLESCreate(PETSC_COMM_WORLD,&sles);
148: SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
149: SLESGetKSP(sles,&ksp);
150: KSPSetInitialGuessNonzero(ksp);
151: SLESSetFromOptions(sles);
153: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154: Solve the linear system
155: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157: SLESSolve(sles,b,u,&its);
159: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: Check solution and clean up
161: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: /* Check error */
164: VecGetOwnershipRange(ustar,&start,&end);
165: for (i=start; i<end; i++) {
166: x = h*(i % (m+1)); y = h*(i/(m+1));
167: val = y;
168: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
169: }
170: VecAssemblyBegin(ustar);
171: VecAssemblyEnd(ustar);
172: VecAXPY(&none,ustar,u);
173: VecNorm(u,NORM_2,&norm);
174: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %dn",norm*h,its);
176: /*
177: Free work space. All PETSc objects should be destroyed when they
178: are no longer needed.
179: */
180: SLESDestroy(sles); VecDestroy(u);
181: VecDestroy(ustar); VecDestroy(b);
182: MatDestroy(A);
184: /*
185: Always call PetscFinalize() before exiting a program. This routine
186: - finalizes the PETSc libraries as well as MPI
187: - provides summary and diagnostic information if certain runtime
188: options are chosen (e.g., -log_summary).
189: */
190: PetscFinalize();
191: return 0;
192: }
194: /* --------------------------------------------------------------------- */
195: /* element stiffness for Laplacian */
196: int FormElementStiffness(double H,Scalar *Ke)
197: {
199: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
200: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
201: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
202: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
203: return(0);
204: }
205: /* --------------------------------------------------------------------- */
206: int FormElementRhs(double x,double y,double H,Scalar *r)
207: {
209: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
210: return(0);
211: }