Actual source code: ex15.c

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

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

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

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

 30: /* Declare routines for user-provided preconditioner */

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

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

 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,&A);
 77:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 78:   MatSetFromOptions(A);

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

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

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

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

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

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

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

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

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

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

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

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

172:     /* (Optional) Set user-defined function to free objects used by custom preconditioner */
173:     PCShellSetDestroy(pc,SampleShellPCDestroy);

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

178:     /* (Optional) Do any setup required for the preconditioner */
179:     /* Note: This function could be set with PCShellSetSetUp and it would be called when necessary */
180:     SampleShellPCSetUp(pc,A,x);

182:   } else {
183:     PCSetType(pc,PCJACOBI);
184:   }

186:   /* 
187:     Set runtime options, e.g.,
188:         -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
189:     These options will override those specified above as long as
190:     KSPSetFromOptions() is called _after_ any other customization
191:     routines.
192:   */
193:   KSPSetFromOptions(ksp);

195:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
196:                       Solve the linear system
197:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

199:   KSPSolve(ksp,b,x);

201:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
202:                       Check solution and clean up
203:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

205:   /* 
206:      Check the error
207:   */
208:   VecAXPY(x,none,u);
209:   VecNorm(x,NORM_2,&norm);
210:   KSPGetIterationNumber(ksp,&its);
211:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %D\n",norm,its);

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

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: PetscErrorCode SampleShellPCCreate(SampleShellPC **shell)
240: {
241:   SampleShellPC  *newctx;

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: .  pc    - preconditioner object
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: PetscErrorCode SampleShellPCSetUp(PC pc,Mat pmat,Vec x)
271: {
272:   SampleShellPC  *shell;
273:   Vec            diag;

276:   PCShellGetContext(pc,(void**)&shell);
277:   VecDuplicate(x,&diag);
278:   MatGetDiagonal(pmat,diag);
279:   VecReciprocal(diag);
280:   shell->diag = diag;

282:   return 0;
283: }
284: /* ------------------------------------------------------------------- */
287: /*
288:    SampleShellPCApply - This routine demonstrates the use of a
289:    user-provided preconditioner.

291:    Input Parameters:
292: +  pc - preconditioner object
293: -  x - input vector

295:    Output Parameter:
296: .  y - preconditioned vector

298:    Notes:
299:    This code implements the Jacobi preconditioner, merely as an
300:    example of working with a PCSHELL.  Note that the Jacobi method
301:    is already provided within PETSc.
302: */
303: PetscErrorCode SampleShellPCApply(PC pc,Vec x,Vec y)
304: {
305:   SampleShellPC   *shell;
306:   PetscErrorCode  ierr;

308:   PCShellGetContext(pc,(void**)&shell);
309:   VecPointwiseMult(y,x,shell->diag);

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

320:    Input Parameter:
321: .  shell - user-defined preconditioner context
322: */
323: PetscErrorCode SampleShellPCDestroy(PC pc)
324: {
325:   SampleShellPC *shell;

328:   PCShellGetContext(pc,(void**)&shell);
329:   VecDestroy(shell->diag);
330:   PetscFree(shell);

332:   return 0;
333: }