Actual source code: ex15.c

  1: /*$Id: ex15.c,v 1.22 2001/03/23 23:23:55 balay Exp $*/

  3: static char help[] = "Solves a linear system in parallel with SLES.  Alson
  4: illustrates setting a user-defined shell preconditioner and using then
  5: Input parameters include:n
  6:   -user_defined_pc : Activate a user-defined preconditionernn";

  8: /*T
  9:    Concepts: SLES^basic parallel example
 10:    Concepts: PC^setting a user-defined shell preconditioner
 11:    Processors: n
 12: T*/

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

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

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

 35: /* 
 36:    User-defined routines.  Note that immediately before each routine below,
 37:    If defined, this macro is used in the PETSc error handlers to provide a
 38:    complete traceback of routine names.  All PETSc library routines use this
 39:    macro, and users can optionally employ it as well in their application
 40:    codes.  Note that users can get a traceback of PETSc errors regardless of
 41:    provides the added traceback detail of the application routine names.
 42: */

 44: int main(int argc,char **args)
 45: {
 46:   Vec           x,b,u;   /* approx solution, RHS, exact solution */
 47:   Mat           A;         /* linear system matrix */
 48:   SLES          sles;      /* linear solver context */
 49:   PC            pc;        /* preconditioner context */
 50:   KSP           ksp;       /* Krylov subspace method context */
 51:   double        norm;      /* norm of solution error */
 52:   SampleShellPC *shell;    /* user-defined preconditioner context */
 53:   Scalar        v,one = 1.0,none = -1.0;
 54:   int           i,j,I,J,Istart,Iend,ierr,m = 8,n = 7,its;
 55:   PetscTruth    user_defined_pc;

 57:   PetscInitialize(&argc,&args,(char *)0,help);
 58:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 59:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);

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

 74:   /* 
 75:      Currently, all PETSc parallel matrix formats are partitioned by
 76:      contiguous chunks of rows across the processors.  Determine which
 77:      rows of the matrix are locally owned. 
 78:   */
 79:   MatGetOwnershipRange(A,&Istart,&Iend);

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

 97:   /* 
 98:      Assemble matrix, using the 2-step process:
 99:        MatAssemblyBegin(), MatAssemblyEnd()
100:      Computations can be done while messages are in transition
101:      by placing code between these two statements.
102:   */
103:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
104:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

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

117:   /* 
118:      Set exact solution; then compute right-hand-side vector.
119:   */
120:   VecSet(&one,u);
121:   MatMult(A,u,b);

123:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
124:                 Create the linear solver and set various options
125:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

127:   /* 
128:      Create linear solver context
129:   */
130:   SLESCreate(PETSC_COMM_WORLD,&sles);

132:   /* 
133:      Set operators. Here the matrix that defines the linear system
134:      also serves as the preconditioning matrix.
135:   */
136:   SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);

138:   /* 
139:      Set linear solver defaults for this problem (optional).
140:      - By extracting the KSP and PC contexts from the SLES context,
141:        we can then directly call any KSP and PC routines
142:        to set various options.
143:   */
144:   SLESGetKSP(sles,&ksp);
145:   SLESGetPC(sles,&pc);
146:   KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,
147:          PETSC_DEFAULT);

149:   /*
150:      Set a user-defined "shell" preconditioner if desired
151:   */
152:   PetscOptionsHasName(PETSC_NULL,"-user_defined_pc",&user_defined_pc);
153:   if (user_defined_pc) {
154:     /* (Required) Indicate to PETSc that we're using a "shell" preconditioner */
155:     PCSetType(pc,PCSHELL);

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

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

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

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

170:   } else {
171:     PCSetType(pc,PCJACOBI);
172:   }

174:   /* 
175:     Set runtime options, e.g.,
176:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
177:     These options will override those specified above as long as
178:     SLESSetFromOptions() is called _after_ any other customization
179:     routines.
180:   */
181:   SLESSetFromOptions(sles);

183:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
184:                       Solve the linear system
185:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

187:   SLESSolve(sles,b,x,&its);

189:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
190:                       Check solution and clean up
191:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

193:   /* 
194:      Check the error
195:   */
196:   VecAXPY(&none,u,x);
197:   VecNorm(x,NORM_2,&norm);
198:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %dn",norm,its);

200:   /* 
201:      Free work space.  All PETSc objects should be destroyed when they
202:      are no longer needed.
203:   */
204:   SLESDestroy(sles);
205:   VecDestroy(u);  VecDestroy(x);
206:   VecDestroy(b);  MatDestroy(A);

208:   if (user_defined_pc) {
209:     SampleShellPCDestroy(shell);
210:   }

212:   PetscFinalize();
213:   return 0;

215: }

217: /***********************************************************************/
218: /*          Routines for a user-defined shell preconditioner           */
219: /***********************************************************************/

221: /*
222:    SampleShellPCCreate - This routine creates a user-defined
223:    preconditioner context.

225:    Output Parameter:
226: .  shell - user-defined preconditioner context
227: */
228: int SampleShellPCCreate(SampleShellPC **shell)
229: {
230:   SampleShellPC *newctx;
231:   int           ierr;

233:   ierr         = PetscNew(SampleShellPC,&newctx);
234:   newctx->diag = 0;
235:   *shell       = newctx;
236:   return 0;
237: }
238: /* ------------------------------------------------------------------- */
239: /*
240:    SampleShellPCSetUp - This routine sets up a user-defined
241:    preconditioner context.  

243:    Input Parameters:
244: .  shell - user-defined preconditioner context
245: .  pmat  - preconditioner matrix
246: .  x     - vector

248:    Output Parameter:
249: .  shell - fully set up user-defined preconditioner context

251:    Notes:
252:    In this example, we define the shell preconditioner to be Jacobi's
253:    method.  Thus, here we create a work vector for storing the reciprocal
254:    of the diagonal of the preconditioner matrix; this vector is then
255:    used within the routine SampleShellPCApply().
256: */
257: int SampleShellPCSetUp(SampleShellPC *shell,Mat pmat,Vec x)
258: {
259:   Vec diag;

262:   VecDuplicate(x,&diag);
263:   MatGetDiagonal(pmat,diag);
264:   VecReciprocal(diag);
265:   shell->diag = diag;

267:   return 0;
268: }
269: /* ------------------------------------------------------------------- */
270: /*
271:    SampleShellPCApply - This routine demonstrates the use of a
272:    user-provided preconditioner.

274:    Input Parameters:
275: .  ctx - optional user-defined context, as set by PCShellSetApply()
276: .  x - input vector

278:    Output Parameter:
279: .  y - preconditioned vector

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

286:    This code implements the Jacobi preconditioner, merely as an
287:    example of working with a PCSHELL.  Note that the Jacobi method
288:    is already provided within PETSc.
289: */
290: int SampleShellPCApply(void *ctx,Vec x,Vec y)
291: {
292:   SampleShellPC *shell = (SampleShellPC*)ctx;
293:   int           ierr;

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

297:   return 0;
298: }
299: /* ------------------------------------------------------------------- */
300: /*
301:    SampleShellPCDestroy - This routine destroys a user-defined
302:    preconditioner context.

304:    Input Parameter:
305: .  shell - user-defined preconditioner context
306: */
307: int SampleShellPCDestroy(SampleShellPC *shell)
308: {

311:   VecDestroy(shell->diag);
312:   PetscFree(shell);

314:   return 0;
315: }