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 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:   PetscReal     norm;      /* norm of solution error */
 52:   SampleShellPC *shell;    /* user-defined preconditioner context */
 53:   PetscScalar   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() VecSetSizes() and VecSetFromOptions(),
109:         we specify only the vector's global
110:         dimension; the parallel partitioning is determined at runtime. 
111:       - Note: We form 1 vector from scratch and then duplicate as needed.
112:   */
113:   VecCreate(PETSC_COMM_WORLD,&u);
114:   VecSetSizes(u,PETSC_DECIDE,m*n);
115:   VecSetFromOptions(u);
116:   VecDuplicate(u,&b);
117:   VecDuplicate(b,&x);

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

125:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
126:                 Create the linear solver and set various options
127:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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

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

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

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

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

172:   } else {
173:     PCSetType(pc,PCJACOBI);
174:   }

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

185:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
186:                       Solve the linear system
187:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

191:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
192:                       Check solution and clean up
193:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

210:   if (user_defined_pc) {
211:     SampleShellPCDestroy(shell);
212:   }

214:   PetscFinalize();
215:   return 0;

217: }

219: /***********************************************************************/
220: /*          Routines for a user-defined shell preconditioner           */
221: /***********************************************************************/

223: /*
224:    SampleShellPCCreate - This routine creates a user-defined
225:    preconditioner context.

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

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

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

250:    Output Parameter:
251: .  shell - fully set up user-defined preconditioner context

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

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

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

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

280:    Output Parameter:
281: .  y - preconditioned vector

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

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

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

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

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

313:   VecDestroy(shell->diag);
314:   PetscFree(shell);

316:   return 0;
317: }