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