Actual source code: ex8.c

  2: static char help[] = "Illustrates use of the preconditioner ASM.\n\
  3: The Additive Schwarz Method for solving a linear system in parallel with KSP.  The\n\
  4: code indicates the procedure for setting user-defined subdomains.  Input\n\
  5: parameters include:\n\
  6:   -user_set_subdomain_solvers:  User explicitly sets subdomain solvers\n\
  7:   -user_set_subdomains:  Activate user-defined subdomains\n\n";

  9: /*
 10:    Note:  This example focuses on setting the subdomains for the ASM 
 11:    preconditioner for a problem on a 2D rectangular grid.  See ex1.c
 12:    and ex2.c for more detailed comments on the basic usage of KSP
 13:    (including working with matrices and vectors).

 15:    The ASM preconditioner is fully parallel, but currently the routine
 16:    PCASMCreateSubdomains2D(), which is used in this example to demonstrate
 17:    user-defined subdomains (activated via -user_set_subdomains), is
 18:    uniprocessor only.

 20:    This matrix in this linear system arises from the discretized Laplacian,
 21:    and thus is not very interesting in terms of experimenting with variants
 22:    of the ASM preconditioner.  
 23: */

 25: /*T
 26:    Concepts: KSP^Additive Schwarz Method (ASM) with user-defined subdomains
 27:    Processors: n
 28: T*/

 30: /* 
 31:   Include "petscksp.h" so that we can use KSP solvers.  Note that this file
 32:   automatically includes:
 33:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 34:      petscmat.h - matrices
 35:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 36:      petscviewer.h - viewers               petscpc.h  - preconditioners
 37: */
 38: #include <petscksp.h>

 42: int main(int argc,char **args)
 43: {
 44:   Vec            x,b,u;                 /* approx solution, RHS, exact solution */
 45:   Mat            A;                       /* linear system matrix */
 46:   KSP            ksp;                    /* linear solver context */
 47:   PC             pc;                      /* PC context */
 48:   IS             *is,*is_local;           /* array of index sets that define the subdomains */
 49:   PetscInt       overlap = 1;             /* width of subdomain overlap */
 50:   PetscInt       Nsub;                    /* number of subdomains */
 51:   PetscInt       m = 15,n = 17;          /* mesh dimensions in x- and y- directions */
 52:   PetscInt       M = 2,N = 1;            /* number of subdomains in x- and y- directions */
 53:   PetscInt       i,j,Ii,J,Istart,Iend;
 55:   PetscMPIInt    size;
 56:   PetscBool      flg;
 57:   PetscBool      user_subdomains = PETSC_FALSE;
 58:   PetscScalar    v, one = 1.0;
 59:   PetscReal      e;

 61:   PetscInitialize(&argc,&args,(char *)0,help);
 62:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 63:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 64:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 65:   PetscOptionsGetInt(PETSC_NULL,"-Mdomains",&M,PETSC_NULL);
 66:   PetscOptionsGetInt(PETSC_NULL,"-Ndomains",&N,PETSC_NULL);
 67:   PetscOptionsGetInt(PETSC_NULL,"-overlap",&overlap,PETSC_NULL);
 68:   PetscOptionsGetBool(PETSC_NULL,"-user_set_subdomains",&user_subdomains,PETSC_NULL);

 70:   /* -------------------------------------------------------------------
 71:          Compute the matrix and right-hand-side vector that define
 72:          the linear system, Ax = b.
 73:      ------------------------------------------------------------------- */

 75:   /* 
 76:      Assemble the matrix for the five point stencil, YET AGAIN 
 77:   */
 78:   MatCreate(PETSC_COMM_WORLD,&A);
 79:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 80:   MatSetFromOptions(A);
 81:   MatGetOwnershipRange(A,&Istart,&Iend);
 82:   for (Ii=Istart; Ii<Iend; Ii++) {
 83:     v = -1.0; i = Ii/n; j = Ii - i*n;
 84:     if (i>0)   {J = Ii - n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 85:     if (i<m-1) {J = Ii + n; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 86:     if (j>0)   {J = Ii - 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 87:     if (j<n-1) {J = Ii + 1; MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);}
 88:     v = 4.0; MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 89:   }
 90:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 91:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 93:   /* 
 94:      Create and set vectors 
 95:   */
 96:   VecCreate(PETSC_COMM_WORLD,&b);
 97:   VecSetSizes(b,PETSC_DECIDE,m*n);
 98:   VecSetFromOptions(b);
 99:   VecDuplicate(b,&u);
100:   VecDuplicate(b,&x);
101:   VecSet(u,one);
102:   MatMult(A,u,b);

104:   /* 
105:      Create linear solver context 
106:   */
107:   KSPCreate(PETSC_COMM_WORLD,&ksp);

109:   /* 
110:      Set operators. Here the matrix that defines the linear system
111:      also serves as the preconditioning matrix.
112:   */
113:   KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);

115:   /* 
116:      Set the default preconditioner for this program to be ASM
117:   */
118:   KSPGetPC(ksp,&pc);
119:   PCSetType(pc,PCASM);

121:   /* -------------------------------------------------------------------
122:                   Define the problem decomposition
123:      ------------------------------------------------------------------- */

125:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
126:        Basic method, should be sufficient for the needs of many users.
127:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

129:      Set the overlap, using the default PETSc decomposition via
130:          PCASMSetOverlap(pc,overlap);
131:      Could instead use the option -pc_asm_overlap <ovl> 

133:      Set the total number of blocks via -pc_asm_blocks <blks>
134:      Note:  The ASM default is to use 1 block per processor.  To
135:      experiment on a single processor with various overlaps, you
136:      must specify use of multiple blocks!
137:   */

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
140:        More advanced method, setting user-defined subdomains
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

143:      Firstly, create index sets that define the subdomains.  The utility
144:      routine PCASMCreateSubdomains2D() is a simple example (that currently
145:      supports 1 processor only!).  More generally, the user should write
146:      a custom routine for a particular problem geometry.

148:      Then call either PCASMSetLocalSubdomains() or PCASMSetTotalSubdomains()
149:      to set the subdomains for the ASM preconditioner.
150:   */

152:   if (!user_subdomains) { /* basic version */
153:     PCASMSetOverlap(pc,overlap);
154:   } else { /* advanced version */
155:     if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"PCASMCreateSubdomains() is currently a uniprocessor routine only!");
156:     PCASMCreateSubdomains2D(m,n,M,N,1,overlap,&Nsub,&is,&is_local);
157:     PCASMSetLocalSubdomains(pc,Nsub,is,is_local);
158:     PetscOptionsGetBool(PETSC_NULL,"-subdomain_view",&flg,PETSC_NULL);
159:     if (flg){
160:       printf("Nmesh points: %d x %d; subdomain partition: %d x %d; overlap: %d; Nsub: %d\n",m,n,M,N,overlap,Nsub);
161:       printf("IS:\n");
162:       for (i=0; i<Nsub; i++){
163:         printf("  IS[%d]\n",i);
164:         ISView(is[i],PETSC_VIEWER_STDOUT_SELF);
165:       }
166:       printf("IS_local:\n");
167:       for (i=0; i<Nsub; i++){
168:         printf("  IS_local[%d]\n",i);
169:         ISView(is_local[i],PETSC_VIEWER_STDOUT_SELF);
170:       }
171:     }
172:   }

174:   /* -------------------------------------------------------------------
175:                 Set the linear solvers for the subblocks
176:      ------------------------------------------------------------------- */

178:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
179:        Basic method, should be sufficient for the needs of most users.
180:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

182:      By default, the ASM preconditioner uses the same solver on each
183:      block of the problem.  To set the same solver options on all blocks,
184:      use the prefix -sub before the usual PC and KSP options, e.g.,
185:           -sub_pc_type <pc> -sub_ksp_type <ksp> -sub_ksp_rtol 1.e-4

187:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
188:         Advanced method, setting different solvers for various blocks.
189:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 

191:      Note that each block's KSP context is completely independent of
192:      the others, and the full range of uniprocessor KSP options is
193:      available for each block.

195:      - Use PCASMGetSubKSP() to extract the array of KSP contexts for
196:        the local blocks.
197:      - See ex7.c for a simple example of setting different linear solvers
198:        for the individual blocks for the block Jacobi method (which is
199:        equivalent to the ASM method with zero overlap).
200:   */

202:   flg  = PETSC_FALSE;
203:   PetscOptionsGetBool(PETSC_NULL,"-user_set_subdomain_solvers",&flg,PETSC_NULL);
204:   if (flg) {
205:     KSP        *subksp;       /* array of KSP contexts for local subblocks */
206:     PetscInt   nlocal,first;  /* number of local subblocks, first local subblock */
207:     PC         subpc;          /* PC context for subblock */
208:     PetscBool  isasm;

210:     PetscPrintf(PETSC_COMM_WORLD,"User explicitly sets subdomain solvers.\n");

212:     /* 
213:        Set runtime options
214:     */
215:     KSPSetFromOptions(ksp);

217:     /* 
218:        Flag an error if PCTYPE is changed from the runtime options
219:      */
220:     PetscTypeCompare((PetscObject)pc,PCASM,&isasm);
221:     if (!isasm) SETERRQ(PETSC_COMM_WORLD,1,"Cannot Change the PCTYPE when manually changing the subdomain solver settings");

223:     /* 
224:        Call KSPSetUp() to set the block Jacobi data structures (including
225:        creation of an internal KSP context for each block).

227:        Note: KSPSetUp() MUST be called before PCASMGetSubKSP().
228:     */
229:     KSPSetUp(ksp);

231:     /*
232:        Extract the array of KSP contexts for the local blocks
233:     */
234:     PCASMGetSubKSP(pc,&nlocal,&first,&subksp);

236:     /*
237:        Loop over the local blocks, setting various KSP options
238:        for each block.  
239:     */
240:     for (i=0; i<nlocal; i++) {
241:       KSPGetPC(subksp[i],&subpc);
242:       PCSetType(subpc,PCILU);
243:       KSPSetType(subksp[i],KSPGMRES);
244:       KSPSetTolerances(subksp[i],1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
245:     }
246:   } else {
247:     /* 
248:        Set runtime options
249:     */
250:     KSPSetFromOptions(ksp);
251:   }

253:   /* -------------------------------------------------------------------
254:                       Solve the linear system
255:      ------------------------------------------------------------------- */

257:   KSPSolve(ksp,b,x);

259:   /* -------------------------------------------------------------------
260:                       Compare result to the exact solution
261:      ------------------------------------------------------------------- */
262:   VecAXPY(x,-1.0,u);
263:   VecNorm(x,NORM_INFINITY, &e);

265:   flg  = PETSC_FALSE;
266:   PetscOptionsGetBool(PETSC_NULL,"-print_error",&flg,PETSC_NULL);
267:   if(flg) {
268:     PetscPrintf(PETSC_COMM_WORLD, "Infinity norm of the error: %G\n", e);
269:   }

271:   /* 
272:      Free work space.  All PETSc objects should be destroyed when they
273:      are no longer needed.
274:   */

276:   if (user_subdomains) {
277:     for (i=0; i<Nsub; i++) {
278:       ISDestroy(&is[i]);
279:       ISDestroy(&is_local[i]);
280:     }
281:     PetscFree(is);
282:     PetscFree(is_local);
283:   }
284:   KSPDestroy(&ksp);
285:   VecDestroy(&u);
286:   VecDestroy(&x);
287:   VecDestroy(&b);
288:   MatDestroy(&A);
289:   PetscFinalize();
290:   return 0;
291: }