Actual source code: ex3.c
1: /*$Id: ex3.c,v 1.33 2001/08/07 21:30:54 bsmith 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(PetscReal,PetscScalar*);
27: extern int FormElementRhs(PetscReal,PetscReal,PetscReal,PetscScalar*);
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: PetscScalar Ke[16]; /* element matrix */
41: PetscScalar r[4]; /* element vector */
42: PetscReal h; /* mesh width */
43: PetscReal norm; /* norm of solution error */
44: PetscReal x,y;
45: PetscScalar 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,&u);
88: VecSetSizes(u,PETSC_DECIDE,N);
89: VecSetFromOptions(u);
90: PetscObjectSetName((PetscObject)u,"Approx. Solution");
91: VecDuplicate(u,&b);
92: PetscObjectSetName((PetscObject)b,"Right hand side");
93: VecDuplicate(b,&ustar);
94: VecSet(&zero,u);
95: VecSet(&zero,b);
97: /*
98: Assemble right-hand-side vector
99: */
100: for (i=start; i<end; i++) {
101: /* location of lower left corner of element */
102: x = h*(i % m); y = h*(i/m);
103: /* node numbers for the four corners of element */
104: idx[0] = (m+1)*(i/m) + (i % m);
105: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
106: FormElementRhs(x,y,h*h,r);
107: VecSetValues(b,4,idx,r,ADD_VALUES);
108: }
109: VecAssemblyBegin(b);
110: VecAssemblyEnd(b);
112: /*
113: Modify matrix and right-hand-side for Dirichlet boundary conditions
114: */
115: PetscMalloc(4*m*sizeof(int),&rows);
116: for (i=0; i<m+1; i++) {
117: rows[i] = i; /* bottom */
118: rows[3*m - 1 +i] = m*(m+1) + i; /* top */
119: }
120: count = m+1; /* left side */
121: for (i=m+1; i<m*(m+1); i+= m+1) {
122: rows[count++] = i;
123: }
124: count = 2*m; /* left side */
125: for (i=2*m+1; i<m*(m+1); i+= m+1) {
126: rows[count++] = i;
127: }
128: ISCreateGeneral(PETSC_COMM_SELF,4*m,rows,&is);
129: for (i=0; i<4*m; i++) {
130: x = h*(rows[i] % (m+1)); y = h*(rows[i]/(m+1));
131: val = y;
132: VecSetValues(u,1,&rows[i],&val,INSERT_VALUES);
133: VecSetValues(b,1,&rows[i],&val,INSERT_VALUES);
134: }
135: PetscFree(rows);
136: VecAssemblyBegin(u);
137: VecAssemblyEnd(u);
138: VecAssemblyBegin(b);
139: VecAssemblyEnd(b);
141: MatZeroRows(A,is,&one);
142: ISDestroy(is);
144: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145: Create the linear solver and set various options
146: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148: SLESCreate(PETSC_COMM_WORLD,&sles);
149: SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
150: SLESGetKSP(sles,&ksp);
151: KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
152: SLESSetFromOptions(sles);
154: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155: Solve the linear system
156: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158: SLESSolve(sles,b,u,&its);
160: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: Check solution and clean up
162: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164: /* Check error */
165: VecGetOwnershipRange(ustar,&start,&end);
166: for (i=start; i<end; i++) {
167: x = h*(i % (m+1)); y = h*(i/(m+1));
168: val = y;
169: VecSetValues(ustar,1,&i,&val,INSERT_VALUES);
170: }
171: VecAssemblyBegin(ustar);
172: VecAssemblyEnd(ustar);
173: VecAXPY(&none,ustar,u);
174: VecNorm(u,NORM_2,&norm);
175: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A Iterations %dn",norm*h,its);
177: /*
178: Free work space. All PETSc objects should be destroyed when they
179: are no longer needed.
180: */
181: SLESDestroy(sles); VecDestroy(u);
182: VecDestroy(ustar); VecDestroy(b);
183: MatDestroy(A);
185: /*
186: Always call PetscFinalize() before exiting a program. This routine
187: - finalizes the PETSc libraries as well as MPI
188: - provides summary and diagnostic information if certain runtime
189: options are chosen (e.g., -log_summary).
190: */
191: PetscFinalize();
192: return 0;
193: }
195: /* --------------------------------------------------------------------- */
196: /* element stiffness for Laplacian */
197: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
198: {
200: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
201: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
202: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
203: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
204: return(0);
205: }
206: /* --------------------------------------------------------------------- */
207: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
208: {
210: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
211: return(0);
212: }