Actual source code: ex15.c

  1: /*$Id: ex15.c,v 1.24 2001/08/07 21:30:54 bsmith Exp $*/

  3: static char help[] = "Solves a linear system in parallel with KSP.  Also\n\
  4: illustrates setting a user-defined shell preconditioner and using the\n\
  5: macro __FUNCT__ to define routine names for use in error handling.\n\
  6: Input parameters include:\n\
  7:   -user_defined_pc : Activate a user-defined preconditioner\n\n";

  9: /*T
 10:    Concepts: KSP^basic parallel example
 11:    Concepts: PC^setting a user-defined shell preconditioner
 12:    Concepts: error handling^Using the macro __FUNCT__ to define routine names;
 13:    Processors: n
 14: T*/

 16: /* 
 17:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 18:   automatically includes:
 19:      petsc.h       - base PETSc routines   petscvec.h - vectors
 20:      petscsys.h    - system routines       petscmat.h - matrices
 21:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 22:      petscviewer.h - viewers               petscpc.h  - preconditioners
 23: */
 24:  #include petscksp.h

 26: /* Define context for user-provided preconditioner */
 27: typedef struct {
 28:   Vec diag;
 29: } SampleShellPC;

 31: /* Declare routines for user-provided preconditioner */
 32: extern int SampleShellPCCreate(SampleShellPC**);
 33: extern int SampleShellPCSetUp(SampleShellPC*,Mat,Vec);
 34: extern int SampleShellPCApply(void*,Vec x,Vec y);
 35: extern int SampleShellPCDestroy(SampleShellPC*);

 37: /* 
 38:    User-defined routines.  Note that immediately before each routine below,
 39:    we define the macro __FUNCT__ to be a string containing the routine name.
 40:    If defined, this macro is used in the PETSc error handlers to provide a
 41:    complete traceback of routine names.  All PETSc library routines use this
 42:    macro, and users can optionally employ it as well in their application
 43:    codes.  Note that users can get a traceback of PETSc errors regardless of
 44:    whether they define __FUNCT__ in application codes; this macro merely
 45:    provides the added traceback detail of the application routine names.
 46: */

 50: int main(int argc,char **args)
 51: {
 52:   Vec           x,b,u;   /* approx solution, RHS, exact solution */
 53:   Mat           A;         /* linear system matrix */
 54:   KSP           ksp;      /* linear solver context */
 55:   PC            pc;        /* preconditioner context */
 56:   PetscReal     norm;      /* norm of solution error */
 57:   SampleShellPC *shell;    /* user-defined preconditioner context */
 58:   PetscScalar   v,one = 1.0,none = -1.0;
 59:   int           i,j,I,J,Istart,Iend,ierr,m = 8,n = 7,its;
 60:   PetscTruth    user_defined_pc;

 62:   PetscInitialize(&argc,&args,(char *)0,help);
 63:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 64:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 67:          Compute the matrix and right-hand-side vector that define
 68:          the linear system, Ax = b.
 69:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 70:   /* 
 71:      Create parallel matrix, specifying only its global dimensions.
 72:      When using MatCreate(), the matrix format can be specified at
 73:      runtime. Also, the parallel partioning of the matrix is
 74:      determined by PETSc at runtime.
 75:   */
 76:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&A);
 77:   MatSetFromOptions(A);

 79:   /* 
 80:      Currently, all PETSc parallel matrix formats are partitioned by
 81:      contiguous chunks of rows across the processors.  Determine which
 82:      rows of the matrix are locally owned. 
 83:   */
 84:   MatGetOwnershipRange(A,&Istart,&Iend);

 86:   /* 
 87:      Set matrix elements for the 2-D, five-point stencil in parallel.
 88:       - Each processor needs to insert only elements that it owns
 89:         locally (but any non-local elements will be sent to the
 90:         appropriate processor during matrix assembly). 
 91:       - Always specify global rows and columns of matrix entries.
 92:    */
 93:   for (I=Istart; I<Iend; I++) {
 94:     v = -1.0; i = I/n; j = I - i*n;
 95:     if (i>0)   {J = I - n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 96:     if (i<m-1) {J = I + n; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 97:     if (j>0)   {J = I - 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 98:     if (j<n-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&v,INSERT_VALUES);}
 99:     v = 4.0; MatSetValues(A,1,&I,1,&I,&v,INSERT_VALUES);
100:   }

102:   /* 
103:      Assemble matrix, using the 2-step process:
104:        MatAssemblyBegin(), MatAssemblyEnd()
105:      Computations can be done while messages are in transition
106:      by placing code between these two statements.
107:   */
108:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
109:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

111:   /* 
112:      Create parallel vectors.
113:       - When using VecCreate() VecSetSizes() and VecSetFromOptions(),
114:         we specify only the vector's global
115:         dimension; the parallel partitioning is determined at runtime. 
116:       - Note: We form 1 vector from scratch and then duplicate as needed.
117:   */
118:   VecCreate(PETSC_COMM_WORLD,&u);
119:   VecSetSizes(u,PETSC_DECIDE,m*n);
120:   VecSetFromOptions(u);
121:   VecDuplicate(u,&b);
122:   VecDuplicate(b,&x);

124:   /* 
125:      Set exact solution; then compute right-hand-side vector.
126:   */
127:   VecSet(&one,u);
128:   MatMult(A,u,b);

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
131:                 Create the linear solver and set various options
132:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

134:   /* 
135:      Create linear solver context
136:   */
137:   KSPCreate(PETSC_COMM_WORLD,&ksp);

139:   /* 
140:      Set operators. Here the matrix that defines the linear system
141:      also serves as the preconditioning matrix.
142:   */
143:   KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);

145:   /* 
146:      Set linear solver defaults for this problem (optional).
147:      - By extracting the KSP and PC contexts from the KSP context,
148:        we can then directly call any KSP and PC routines
149:        to set various options.
150:   */
151:   KSPGetPC(ksp,&pc);
152:   KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,
153:          PETSC_DEFAULT);

155:   /*
156:      Set a user-defined "shell" preconditioner if desired
157:   */
158:   PetscOptionsHasName(PETSC_NULL,"-user_defined_pc",&user_defined_pc);
159:   if (user_defined_pc) {
160:     /* (Required) Indicate to PETSc that we're using a "shell" preconditioner */
161:     PCSetType(pc,PCSHELL);

163:     /* (Optional) Create a context for the user-defined preconditioner; this
164:        context can be used to contain any application-specific data. */
165:     SampleShellPCCreate(&shell);

167:     /* (Required) Set the user-defined routine for applying the preconditioner */
168:     PCShellSetApply(pc,SampleShellPCApply,(void*)shell);

170:     /* (Optional) Set a name for the preconditioner, used for PCView() */
171:     PCShellSetName(pc,"MyPreconditioner");

173:     /* (Optional) Do any setup required for the preconditioner */
174:     SampleShellPCSetUp(shell,A,x);

176:   } else {
177:     PCSetType(pc,PCJACOBI);
178:   }

180:   /* 
181:     Set runtime options, e.g.,
182:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
183:     These options will override those specified above as long as
184:     KSPSetFromOptions() is called _after_ any other customization
185:     routines.
186:   */
187:   KSPSetFromOptions(ksp);

189:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
190:                       Solve the linear system
191:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

193:   KSPSetRhs(ksp,b);
194:   KSPSetSolution(ksp,x);
195:   KSPSolve(ksp);

197:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
198:                       Check solution and clean up
199:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

201:   /* 
202:      Check the error
203:   */
204:   VecAXPY(&none,u,x);
205:   VecNorm(x,NORM_2,&norm);
206:   KSPGetIterationNumber(ksp,&its);
207:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %d\n",norm,its);

209:   /* 
210:      Free work space.  All PETSc objects should be destroyed when they
211:      are no longer needed.
212:   */
213:   KSPDestroy(ksp);
214:   VecDestroy(u);  VecDestroy(x);
215:   VecDestroy(b);  MatDestroy(A);

217:   if (user_defined_pc) {
218:     SampleShellPCDestroy(shell);
219:   }

221:   PetscFinalize();
222:   return 0;

224: }

226: /***********************************************************************/
227: /*          Routines for a user-defined shell preconditioner           */
228: /***********************************************************************/

232: /*
233:    SampleShellPCCreate - This routine creates a user-defined
234:    preconditioner context.

236:    Output Parameter:
237: .  shell - user-defined preconditioner context
238: */
239: int SampleShellPCCreate(SampleShellPC **shell)
240: {
241:   SampleShellPC *newctx;
242:   int           ierr;

244:   PetscNew(SampleShellPC,&newctx);
245:   newctx->diag = 0;
246:   *shell       = newctx;
247:   return 0;
248: }
249: /* ------------------------------------------------------------------- */
252: /*
253:    SampleShellPCSetUp - This routine sets up a user-defined
254:    preconditioner context.  

256:    Input Parameters:
257: .  shell - user-defined preconditioner context
258: .  pmat  - preconditioner matrix
259: .  x     - vector

261:    Output Parameter:
262: .  shell - fully set up user-defined preconditioner context

264:    Notes:
265:    In this example, we define the shell preconditioner to be Jacobi's
266:    method.  Thus, here we create a work vector for storing the reciprocal
267:    of the diagonal of the preconditioner matrix; this vector is then
268:    used within the routine SampleShellPCApply().
269: */
270: int SampleShellPCSetUp(SampleShellPC *shell,Mat pmat,Vec x)
271: {
272:   Vec diag;

275:   VecDuplicate(x,&diag);
276:   MatGetDiagonal(pmat,diag);
277:   VecReciprocal(diag);
278:   shell->diag = diag;

280:   return 0;
281: }
282: /* ------------------------------------------------------------------- */
285: /*
286:    SampleShellPCApply - This routine demonstrates the use of a
287:    user-provided preconditioner.

289:    Input Parameters:
290: .  ctx - optional user-defined context, as set by PCShellSetApply()
291: .  x - input vector

293:    Output Parameter:
294: .  y - preconditioned vector

296:    Notes:
297:    Note that the PCSHELL preconditioner passes a void pointer as the
298:    first input argument.  This can be cast to be the whatever the user
299:    has set (via PCSetShellApply()) the application-defined context to be.

301:    This code implements the Jacobi preconditioner, merely as an
302:    example of working with a PCSHELL.  Note that the Jacobi method
303:    is already provided within PETSc.
304: */
305: int SampleShellPCApply(void *ctx,Vec x,Vec y)
306: {
307:   SampleShellPC *shell = (SampleShellPC*)ctx;
308:   int           ierr;

310:   VecPointwiseMult(x,shell->diag,y);

312:   return 0;
313: }
314: /* ------------------------------------------------------------------- */
317: /*
318:    SampleShellPCDestroy - This routine destroys a user-defined
319:    preconditioner context.

321:    Input Parameter:
322: .  shell - user-defined preconditioner context
323: */
324: int SampleShellPCDestroy(SampleShellPC *shell)
325: {

328:   VecDestroy(shell->diag);
329:   PetscFree(shell);

331:   return 0;
332: }