Actual source code: ex13.c
1: /*$Id: ex13.c,v 1.29 2001/08/07 21:30:54 bsmith Exp $*/
3: static char help[] = "Solves a variable Poisson problem with SLES.nn";
5: /*T
6: Concepts: SLES^basic sequential example
7: Concepts: SLES^Laplacian, 2d
8: Concepts: Laplacian, 2d
9: Processors: 1
10: T*/
12: /*
13: Include "petscsles.h" so that we can use SLES solvers. Note that this file
14: automatically includes:
15: petsc.h - base PETSc routines petscvec.h - vectors
16: petscsys.h - system routines petscmat.h - matrices
17: petscis.h - index sets petscksp.h - Krylov subspace methods
18: petscviewer.h - viewers petscpc.h - preconditioners
19: */
20: #include petscsles.h
22: /*
23: User-defined context that contains all the data structures used
24: in the linear solution process.
25: */
26: typedef struct {
27: Vec x,b; /* solution vector, right-hand-side vector */
28: Mat A; /* sparse matrix */
29: SLES sles; /* linear solver context */
30: int m,n; /* grid dimensions */
31: PetscScalar hx2,hy2; /* 1/(m+1)*(m+1) and 1/(n+1)*(n+1) */
32: } UserCtx;
34: extern int UserInitializeLinearSolver(int,int,UserCtx *);
35: extern int UserFinalizeLinearSolver(UserCtx *);
36: extern int UserDoLinearSolver(PetscScalar *,UserCtx *userctx,PetscScalar *b,PetscScalar *x);
38: int main(int argc,char **args)
39: {
40: UserCtx userctx;
41: int ierr,m = 6,n = 7,t,tmax = 2,i,I,j,N;
42: PetscScalar *userx,*rho,*solution,*userb,hx,hy,x,y;
43: PetscReal enorm;
45: /*
46: Initialize the PETSc libraries
47: */
48: PetscInitialize(&argc,&args,(char *)0,help);
50: /*
51: The next two lines are for testing only; these allow the user to
52: decide the grid size at runtime.
53: */
54: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
55: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
57: /*
58: Create the empty sparse matrix and linear solver data structures
59: */
60: UserInitializeLinearSolver(m,n,&userctx);
61: N = m*n;
63: /*
64: Allocate arrays to hold the solution to the linear system.
65: This is not normally done in PETSc programs, but in this case,
66: since we are calling these routines from a non-PETSc program, we
67: would like to reuse the data structures from another code. So in
68: the context of a larger application these would be provided by
69: other (non-PETSc) parts of the application code.
70: */
71: PetscMalloc(N*sizeof(PetscScalar),&userx);
72: PetscMalloc(N*sizeof(PetscScalar),&userb);
73: PetscMalloc(N*sizeof(PetscScalar),&solution);
75: /*
76: Allocate an array to hold the coefficients in the elliptic operator
77: */
78: PetscMalloc(N*sizeof(PetscScalar),&rho);
80: /*
81: Fill up the array rho[] with the function rho(x,y) = x; fill the
82: right-hand-side b[] and the solution with a known problem for testing.
83: */
84: hx = 1.0/(m+1);
85: hy = 1.0/(n+1);
86: y = hy;
87: I = 0;
88: for (j=0; j<n; j++) {
89: x = hx;
90: for (i=0; i<m; i++) {
91: rho[I] = x;
92: solution[I] = PetscSinScalar(2.*PETSC_PI*x)*PetscSinScalar(2.*PETSC_PI*y);
93: userb[I] = -2*PETSC_PI*PetscCosScalar(2*PETSC_PI*x)*PetscSinScalar(2*PETSC_PI*y) +
94: 8*PETSC_PI*PETSC_PI*x*PetscSinScalar(2*PETSC_PI*x)*PetscSinScalar(2*PETSC_PI*y);
95: x += hx;
96: I++;
97: }
98: y += hy;
99: }
101: /*
102: Loop over a bunch of timesteps, setting up and solver the linear
103: system for each time-step.
105: Note this is somewhat artificial. It is intended to demonstrate how
106: one may reuse the linear solver stuff in each time-step.
107: */
108: for (t=0; t<tmax; t++) {
109: UserDoLinearSolver(rho,&userctx,userb,userx);
111: /*
112: Compute error: Note that this could (and usually should) all be done
113: using the PETSc vector operations. Here we demonstrate using more
114: standard programming practices to show how they may be mixed with
115: PETSc.
116: */
117: enorm = 0.0;
118: for (i=0; i<N; i++) {
119: enorm += PetscRealPart(PetscConj(solution[i]-userx[i])*(solution[i]-userx[i]));
120: }
121: enorm *= PetscRealPart(hx*hy);
122: printf("m %d n %d error norm %gn",m,n,enorm);
123: }
125: /*
126: We are all finished solving linear systems, so we clean up the
127: data structures.
128: */
129: PetscFree(rho);
130: PetscFree(solution);
131: PetscFree(userx);
132: PetscFree(userb);
133: UserFinalizeLinearSolver(&userctx);
134: PetscFinalize();
136: return 0;
137: }
139: /* ------------------------------------------------------------------------*/
140: int UserInitializeLinearSolver(int m,int n,UserCtx *userctx)
141: {
142: int N,ierr;
144: /*
145: Here we assume use of a grid of size m x n, with all points on the
146: interior of the domain, i.e., we do not include the points corresponding
147: to homogeneous Dirichlet boundary conditions. We assume that the domain
148: is [0,1]x[0,1].
149: */
150: userctx->m = m;
151: userctx->n = n;
152: userctx->hx2 = (m+1)*(m+1);
153: userctx->hy2 = (n+1)*(n+1);
154: N = m*n;
156: /*
157: Create the sparse matrix. Preallocate 5 nonzeros per row.
158: */
159: MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,5,0,&userctx->A);
161: /*
162: Create vectors. Here we create vectors with no memory allocated.
163: This way, we can use the data structures already in the program
164: by using VecPlaceArray() subroutine at a later stage.
165: */
166: VecCreateSeqWithArray(PETSC_COMM_SELF,N,PETSC_NULL,&userctx->b);
167: VecDuplicate(userctx->b,&userctx->x);
169: /*
170: Create linear solver context. This will be used repeatedly for all
171: the linear solves needed.
172: */
173: SLESCreate(PETSC_COMM_SELF,&userctx->sles);
175: return 0;
176: }
178: /*
179: Solves -div (rho grad psi) = F using finite differences.
180: rho is a 2-dimensional array of size m by n, stored in Fortran
181: style by columns. userb is a standard one-dimensional array.
182: */
183: /* ------------------------------------------------------------------------*/
184: int UserDoLinearSolver(PetscScalar *rho,UserCtx *userctx,PetscScalar *userb,PetscScalar *userx)
185: {
186: int ierr,i,j,I,J,m = userctx->m,n = userctx->n,its;
187: Mat A = userctx->A;
188: PC pc;
189: PetscScalar v,hx2 = userctx->hx2,hy2 = userctx->hy2;
191: /*
192: This is not the most efficient way of generating the matrix
193: but let's not worry about it. We should have separate code for
194: the four corners, each edge and then the interior. Then we won't
195: have the slow if-tests inside the loop.
197: Computes the operator
198: -div rho grad
199: on an m by n grid with zero Dirichlet boundary conditions. The rho
200: is assumed to be given on the same grid as the finite difference
201: stencil is applied. For a staggered grid, one would have to change
202: things slightly.
203: */
204: I = 0;
205: for (j=0; j<n; j++) {
206: for (i=0; i<m; i++) {
207: if (j>0) {
208: J = I - m;
209: v = -.5*(rho[I] + rho[J])*hy2;
210: MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);
211: }
212: if (j<n-1) {
213: J = I + m;
214: v = -.5*(rho[I] + rho[J])*hy2;
215: MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);
216: }
217: if (i>0) {
218: J = I - 1;
219: v = -.5*(rho[I] + rho[J])*hx2;
220: MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);
221: }
222: if (i<m-1) {
223: J = I + 1;
224: v = -.5*(rho[I] + rho[J])*hx2;
225: MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);
226: }
227: v = 2.0*rho[I]*(hx2+hy2);
228: MatSetValues(A,1,&I,1,&I,&v,INSERT_VALUES);
229: I++;
230: }
231: }
233: /*
234: Assemble matrix
235: */
236: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
237: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
239: /*
240: Set operators. Here the matrix that defines the linear system
241: also serves as the preconditioning matrix. Since all the matrices
242: will have the same nonzero pattern here, we indicate this so the
243: linear solvers can take advantage of this.
244: */
245: SLESSetOperators(userctx->sles,A,A,SAME_NONZERO_PATTERN);
247: /*
248: Set linear solver defaults for this problem (optional).
249: - Here we set it to use direct LU factorization for the solution
250: */
251: SLESGetPC(userctx->sles,&pc);
252: PCSetType(pc,PCLU);
254: /*
255: Set runtime options, e.g.,
256: -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
257: These options will override those specified above as long as
258: SLESSetFromOptions() is called _after_ any other customization
259: routines.
260:
261: Run the program with the option -help to see all the possible
262: linear solver options.
263: */
264: SLESSetFromOptions(userctx->sles);
266: /*
267: This allows the PETSc linear solvers to compute the solution
268: directly in the user's array rather than in the PETSc vector.
269:
270: This is essentially a hack and not highly recommend unless you
271: are quite comfortable with using PETSc. In general, users should
272: write their entire application using PETSc vectors rather than
273: arrays.
274: */
275: VecPlaceArray(userctx->x,userx);
276: VecPlaceArray(userctx->b,userb);
278: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
279: Solve the linear system
280: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
282: SLESSolve(userctx->sles,userctx->b,userctx->x,&its);
284: /*
285: Put back the PETSc array that belongs in the vector xuserctx->x
286: */
288: return 0;
289: }
291: /* ------------------------------------------------------------------------*/
292: int UserFinalizeLinearSolver(UserCtx *userctx)
293: {
295: /*
296: We are all done and don't need to solve any more linear systems, so
297: we free the work space. All PETSc objects should be destroyed when
298: they are no longer needed.
299: */
300: SLESDestroy(userctx->sles);
301: VecDestroy(userctx->x);
302: VecDestroy(userctx->b);
303: MatDestroy(userctx->A);
304: return 0;
305: }