Actual source code: ex7.c

  1: /*$Id: ex7.c,v 1.56 2001/04/10 19:36:40 bsmith Exp $*/

  3: static char help[] = "Block Jacobi preconditioner for solving a linear system in parallel with SLES.n
  4: The code indicates then
  5: procedures for setting the particular block sizes and for using differentn
  6: linear solvers on the individual blocks.nn";

  8: /*
  9:    Note:  This example focuses on ways to customize the block Jacobi
 10:    preconditioner. See ex1.c and ex2.c for more detailed comments on
 11:    the basic usage of SLES (including working with matrices and vectors).

 13:    Recall: The block Jacobi method is equivalent to the ASM preconditioner
 14:    with zero overlap.
 15:  */

 17: /*T
 18:    Concepts: SLES^customizing the block Jacobi preconditioner
 19:    Processors: n
 20: T*/

 22: /* 
 23:   Include "petscsles.h" so that we can use SLES solvers.  Note that this file
 24:   automatically includes:
 25:      petsc.h       - base PETSc routines   petscvec.h - vectors
 26:      petscsys.h    - system routines       petscmat.h - matrices
 27:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 28:      petscviewer.h - viewers               petscpc.h  - preconditioners
 29: */
 30:  #include petscsles.h

 32: int main(int argc,char **args)
 33: {
 34:   Vec        x,b,u;      /* approx solution, RHS, exact solution */
 35:   Mat        A;            /* linear system matrix */
 36:   SLES       sles;         /* SLES context */
 37:   SLES       *subsles;     /* array of local SLES contexts on this processor */
 38:   PC         pc;           /* PC context */
 39:   PC         subpc;        /* PC context for subdomain */
 40:   KSP        subksp;       /* KSP context for subdomain */
 41:   double     norm;         /* norm of solution error */
 42:   int        i,j,I,J,ierr,*blks,m = 8,n;
 43:   int        rank,size,its,nlocal,first,Istart,Iend;
 44:   Scalar     v,one = 1.0,none = -1.0;
 45:   PetscTruth isbjacobi,flg;

 47:   PetscInitialize(&argc,&args,(char *)0,help);
 48:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 49:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 50:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 51:   n = m+2;

 53:   /* -------------------------------------------------------------------
 54:          Compute the matrix and right-hand-side vector that define
 55:          the linear system, Ax = b.
 56:      ------------------------------------------------------------------- */

 58:   /* 
 59:      Create and assemble parallel matrix
 60:   */
 61:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&A);
 62:   MatSetFromOptions(A);
 63:   MatGetOwnershipRange(A,&Istart,&Iend);
 64:   for (I=Istart; I<Iend; I++) {
 65:     v = -1.0; i = I/n; j = I - i*n;
 66:     if (i>0)   {J = I - n; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
 67:     if (i<m-1) {J = I + n; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
 68:     if (j>0)   {J = I - 1; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
 69:     if (j<n-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&v,ADD_VALUES);}
 70:     v = 4.0; MatSetValues(A,1,&I,1,&I,&v,ADD_VALUES);
 71:   }
 72:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 73:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 75:   /*
 76:      Create parallel vectors
 77:   */
 78:   VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,m*n,&u);
 79:   VecSetFromOptions(u);
 80:   VecDuplicate(u,&b);
 81:   VecDuplicate(b,&x);

 83:   /*
 84:      Set exact solution; then compute right-hand-side vector.
 85:   */
 86:   VecSet(&one,u);
 87:   MatMult(A,u,b);

 89:   /*
 90:      Create linear solver context
 91:   */
 92:   SLESCreate(PETSC_COMM_WORLD,&sles);

 94:   /* 
 95:      Set operators. Here the matrix that defines the linear system
 96:      also serves as the preconditioning matrix.
 97:   */
 98:   SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);

100:   /*
101:      Set default preconditioner for this program to be block Jacobi.
102:      This choice can be overridden at runtime with the option
103:         -pc_type <type>
104:   */
105:   SLESGetPC(sles,&pc);
106:   PCSetType(pc,PCBJACOBI);


109:   /* -------------------------------------------------------------------
110:                    Define the problem decomposition
111:      ------------------------------------------------------------------- */

113:   /* 
114:      Call PCBJacobiSetTotalBlocks() to set individually the size of
115:      each block in the preconditioner.  This could also be done with
116:      the runtime option
117:          -pc_bjacobi_blocks <blocks>
118:      Also, see the command PCBJacobiSetLocalBlocks() to set the
119:      local blocks.

121:       Note: The default decomposition is 1 block per processor.
122:   */
123:   PetscMalloc(m*sizeof(int),&blks);
124:   for (i=0; i<m; i++) blks[i] = n;
125:   PCBJacobiSetTotalBlocks(pc,m,blks);
126:   PetscFree(blks);


129:   /* -------------------------------------------------------------------
130:                Set the linear solvers for the subblocks
131:      ------------------------------------------------------------------- */

133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
134:        Basic method, should be sufficient for the needs of most users.
135:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

137:      By default, the block Jacobi method uses the same solver on each
138:      block of the problem.  To set the same solver options on all blocks, 
139:      use the prefix -sub before the usual PC and KSP options, e.g.,
140:           -sub_pc_type <pc> -sub_ksp_type <ksp> -sub_ksp_rtol 1.e-4
141:   */

143:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
144:         Advanced method, setting different solvers for various blocks.
145:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

147:      Note that each block's SLES context is completely independent of
148:      the others, and the full range of uniprocessor SLES options is
149:      available for each block. The following section of code is intended
150:      to be a simple illustration of setting different linear solvers for
151:      the individual blocks.  These choices are obviously not recommended
152:      for solving this particular problem.
153:   */
154:   PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);
155:   if (isbjacobi) {
156:     /* 
157:        Call SLESSetUp() to set the block Jacobi data structures (including
158:        creation of an internal SLES context for each block).

160:        Note: SLESSetUp() MUST be called before PCBJacobiGetSubSLES().
161:     */
162:     SLESSetUp(sles,x,b);

164:     /*
165:        Extract the array of SLES contexts for the local blocks
166:     */
167:     PCBJacobiGetSubSLES(pc,&nlocal,&first,&subsles);

169:     /*
170:        Loop over the local blocks, setting various SLES options
171:        for each block.  
172:     */
173:     for (i=0; i<nlocal; i++) {
174:       SLESGetPC(subsles[i],&subpc);
175:       SLESGetKSP(subsles[i],&subksp);
176:       if (!rank) {
177:         if (i%2) {
178:           PCSetType(subpc,PCILU);
179:         } else {
180:           PCSetType(subpc,PCNONE);
181:           KSPSetType(subksp,KSPBCGS);
182:           KSPSetTolerances(subksp,1.e-6,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
183:         }
184:       } else {
185:         PCSetType(subpc,PCJACOBI);
186:         KSPSetType(subksp,KSPGMRES);
187:         KSPSetTolerances(subksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
188:       }
189:     }
190:   }

192:   /* -------------------------------------------------------------------
193:                       Solve the linear system
194:      ------------------------------------------------------------------- */

196:   /* 
197:     Set runtime options
198:   */
199:   SLESSetFromOptions(sles);

201:   /*
202:      Solve the linear system
203:   */
204:   SLESSolve(sles,b,x,&its);

206:   /*
207:      View info about the solver
208:   */
209:   PetscOptionsHasName(PETSC_NULL,"-noslesview",&flg);
210:   if (!flg) {
211:     SLESView(sles,PETSC_VIEWER_STDOUT_WORLD);
212:   }

214:   /* -------------------------------------------------------------------
215:                       Check solution and clean up
216:      ------------------------------------------------------------------- */

218:   /*
219:      Check the error
220:   */
221:   VecAXPY(&none,u,x);
222:   VecNorm(x,NORM_2,&norm);
223:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A iterations %dn",norm,its);

225:   /* 
226:      Free work space.  All PETSc objects should be destroyed when they
227:      are no longer needed.
228:   */
229:   SLESDestroy(sles);
230:   VecDestroy(u);  VecDestroy(x);
231:   VecDestroy(b);  MatDestroy(A);
232:   PetscFinalize();
233:   return 0;
234: }