Actual source code: redundant.c

  1: #define PETSCKSP_DLL

  3: /*
  4:   This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner.
  5: */
 6:  #include private/pcimpl.h
 7:  #include petscksp.h

  9: typedef struct {
 10:   KSP          ksp;
 11:   PC           pc;                   /* actual preconditioner used on each processor */
 12:   Vec          xsub,ysub;            /* vectors of a subcommunicator to hold parallel vectors of ((PetscObject)pc)->comm */
 13:   Vec          xdup,ydup;            /* parallel vector that congregates xsub or ysub facilitating vector scattering */
 14:   Mat          pmats;                /* matrix and optional preconditioner matrix belong to a subcommunicator */
 15:   VecScatter   scatterin,scatterout; /* scatter used to move all values to each processor group (subcommunicator) */
 16:   PetscTruth   useparallelmat;
 17:   PetscSubcomm psubcomm;
 18:   PetscInt     nsubcomm;           /* num of data structure PetscSubcomm */
 19: } PC_Redundant;

 23: static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer)
 24: {
 25:   PC_Redundant   *red = (PC_Redundant*)pc->data;
 27:   PetscMPIInt    rank;
 28:   PetscTruth     iascii,isstring;
 29:   PetscViewer    subviewer;

 32:   MPI_Comm_rank(((PetscObject)pc)->comm,&rank);
 33:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
 34:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
 35:   if (iascii) {
 36:     if (!red->psubcomm) {
 37:       PetscViewerASCIIPrintf(viewer,"  Redundant preconditioner: Not yet setup\n");
 38:     } else {
 39:       PetscViewerASCIIPrintf(viewer,"  Redundant preconditioner: First (color=0) of %D PCs follows\n",red->nsubcomm);
 40:       PetscViewerGetSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);
 41:       if (red->psubcomm->color) { /* only view first redundant pc */
 42:         PetscViewerASCIIPushTab(viewer);
 43:         KSPView(red->ksp,subviewer);
 44:         PetscViewerASCIIPopTab(viewer);
 45:       }
 46:       PetscViewerRestoreSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);
 47:     }
 48:   } else if (isstring) {
 49:     PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");
 50:   } else {
 51:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name);
 52:   }
 53:   return(0);
 54: }

 56:  #include private/matimpl.h
 59: static PetscErrorCode PCSetUp_Redundant(PC pc)
 60: {
 61:   PC_Redundant   *red = (PC_Redundant*)pc->data;
 63:   PetscInt       mstart,mend,mlocal,m,mlocal_sub,rstart_sub,rend_sub,mloc_sub;
 64:   PetscMPIInt    size;
 65:   MatReuse       reuse = MAT_INITIAL_MATRIX;
 66:   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
 67:   MPI_Comm       comm = ((PetscObject)pc)->comm,subcomm;
 68:   Vec            vec;
 69:   PetscMPIInt    subsize,subrank;
 70:   const char     *prefix;

 73:   MatGetVecs(pc->pmat,&vec,0);
 74:   VecGetSize(vec,&m);

 76:   if (!pc->setupcalled) {
 77:     if (!red->psubcomm) {
 78:       PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);
 79:       PetscLogObjectMemory(pc,sizeof(PetscSubcomm));

 81:       /* create a new PC that processors in each subcomm have copy of */
 82:       subcomm = red->psubcomm->comm;
 83:       KSPCreate(subcomm,&red->ksp);
 84:       PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);
 85:       PetscLogObjectParent(pc,red->ksp);
 86:       KSPSetType(red->ksp,KSPPREONLY);
 87:       KSPGetPC(red->ksp,&red->pc);
 88:       PCSetType(red->pc,PCLU);

 90:       PCGetOptionsPrefix(pc,&prefix);
 91:       KSPSetOptionsPrefix(red->ksp,prefix);
 92:       KSPAppendOptionsPrefix(red->ksp,"redundant_");
 93:     } else {
 94:        subcomm = red->psubcomm->comm;
 95:     }

 97:     /* create working vectors xsub/ysub and xdup/ydup */
 98:     VecGetLocalSize(vec,&mlocal);
 99:     VecGetOwnershipRange(vec,&mstart,&mend);

101:     /* get local size of xsub/ysub */
102:     MPI_Comm_size(subcomm,&subsize);
103:     MPI_Comm_rank(subcomm,&subrank);
104:     rstart_sub = pc->pmat->rmap->range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */
105:     if (subrank+1 < subsize){
106:       rend_sub = pc->pmat->rmap->range[red->psubcomm->n*(subrank+1)];
107:     } else {
108:       rend_sub = m;
109:     }
110:     mloc_sub = rend_sub - rstart_sub;
111:     VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);
112:     /* create xsub with empty local arrays, because xdup's arrays will be placed into it */
113:     VecCreateMPIWithArray(subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);

115:     /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it. 
116:        Note: we use communicator dupcomm, not ((PetscObject)pc)->comm! */
117:     VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);
118:     VecCreateMPIWithArray(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);
119: 
120:     /* create vec scatters */
121:     if (!red->scatterin){
122:       IS       is1,is2;
123:       PetscInt *idx1,*idx2,i,j,k;

125:       PetscMalloc2(red->psubcomm->n*mlocal,PetscInt,&idx1,red->psubcomm->n*mlocal,PetscInt,&idx2);
126:       j = 0;
127:       for (k=0; k<red->psubcomm->n; k++){
128:         for (i=mstart; i<mend; i++){
129:           idx1[j]   = i;
130:           idx2[j++] = i + m*k;
131:         }
132:       }
133:       ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,&is1);
134:       ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,&is2);
135:       VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);
136:       ISDestroy(is1);
137:       ISDestroy(is2);

139:       ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);
140:       ISCreateStride(comm,mlocal,mstart,1,&is2);
141:       VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);
142:       ISDestroy(is1);
143:       ISDestroy(is2);
144:       PetscFree2(idx1,idx2);
145:     }
146:   }
147:   VecDestroy(vec);

149:   /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
150:   MPI_Comm_size(comm,&size);
151:   if (size == 1) {
152:     red->useparallelmat = PETSC_FALSE;
153:   }

155:   if (red->useparallelmat) {
156:     if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) {
157:       /* destroy old matrices */
158:       if (red->pmats) {
159:         MatDestroy(red->pmats);
160:       }
161:     } else if (pc->setupcalled == 1) {
162:       reuse = MAT_REUSE_MATRIX;
163:       str   = SAME_NONZERO_PATTERN;
164:     }
165: 
166:     /* grab the parallel matrix and put it into processors of a subcomminicator */
167:     /*--------------------------------------------------------------------------*/
168:     VecGetLocalSize(red->ysub,&mlocal_sub);
169:     MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);
170:     /* tell PC of the subcommunicator its operators */
171:     KSPSetOperators(red->ksp,red->pmats,red->pmats,str);
172:   } else {
173:     KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);
174:   }
175:   if (pc->setfromoptionscalled){
176:     KSPSetFromOptions(red->ksp);
177:   }
178:   KSPSetUp(red->ksp);
179:   return(0);
180: }

184: static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
185: {
186:   PC_Redundant   *red = (PC_Redundant*)pc->data;
188:   PetscScalar    *array;

191:   /* scatter x to xdup */
192:   VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);
193:   VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);
194: 
195:   /* place xdup's local array into xsub */
196:   VecGetArray(red->xdup,&array);
197:   VecPlaceArray(red->xsub,(const PetscScalar*)array);

199:   /* apply preconditioner on each processor */
200:   PCApply(red->pc,red->xsub,red->ysub);
201:   VecResetArray(red->xsub);
202:   VecRestoreArray(red->xdup,&array);
203: 
204:   /* place ysub's local array into ydup */
205:   VecGetArray(red->ysub,&array);
206:   VecPlaceArray(red->ydup,(const PetscScalar*)array);

208:   /* scatter ydup to y */
209:   VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);
210:   VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);
211:   VecResetArray(red->ydup);
212:   VecRestoreArray(red->ysub,&array);
213:   return(0);
214: }

218: static PetscErrorCode PCDestroy_Redundant(PC pc)
219: {
220:   PC_Redundant   *red = (PC_Redundant*)pc->data;

224:   if (red->scatterin)  {VecScatterDestroy(red->scatterin);}
225:   if (red->scatterout) {VecScatterDestroy(red->scatterout);}
226:   if (red->ysub)       {VecDestroy(red->ysub);}
227:   if (red->xsub)       {VecDestroy(red->xsub);}
228:   if (red->xdup)       {VecDestroy(red->xdup);}
229:   if (red->ydup)       {VecDestroy(red->ydup);}
230:   if (red->pmats) {
231:     MatDestroy(red->pmats);
232:   }
233:   if (red->psubcomm) {PetscSubcommDestroy(red->psubcomm);}
234:   if (red->ksp) {KSPDestroy(red->ksp);}
235:   PetscFree(red);
236:   return(0);
237: }

241: static PetscErrorCode PCSetFromOptions_Redundant(PC pc)
242: {
244:   PC_Redundant   *red = (PC_Redundant*)pc->data;

247:   PetscOptionsHead("Redundant options");
248:   PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);
249:   PetscOptionsTail();
250:   return(0);
251: }

256: PetscErrorCode  PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds)
257: {
258:   PC_Redundant *red = (PC_Redundant*)pc->data;

261:   red->nsubcomm = nreds;
262:   return(0);
263: }

268: /*@
269:    PCRedundantSetNumber - Sets the number of redundant preconditioner contexts.

271:    Collective on PC

273:    Input Parameters:
274: +  pc - the preconditioner context
275: -  nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and
276:                               use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes.

278:    Level: advanced

280: .keywords: PC, redundant solve
281: @*/
282: PetscErrorCode  PCRedundantSetNumber(PC pc,PetscInt nredundant)
283: {
284:   PetscErrorCode ierr,(*f)(PC,PetscInt);

288:   if (nredundant <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant);
289:   PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetNumber_C",(void (**)(void))&f);
290:   if (f) {
291:     (*f)(pc,nredundant);
292:   }
293:   return(0);
294: }

299: PetscErrorCode  PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
300: {
301:   PC_Redundant   *red = (PC_Redundant*)pc->data;

305:   PetscObjectReference((PetscObject)in);
306:   if (red->scatterin) { VecScatterDestroy(red->scatterin); }
307:   red->scatterin  = in;
308:   PetscObjectReference((PetscObject)out);
309:   if (red->scatterout) { VecScatterDestroy(red->scatterout); }
310:   red->scatterout = out;
311:   return(0);
312: }

317: /*@
318:    PCRedundantSetScatter - Sets the scatter used to copy values into the
319:      redundant local solve and the scatter to move them back into the global
320:      vector.

322:    Collective on PC

324:    Input Parameters:
325: +  pc - the preconditioner context
326: .  in - the scatter to move the values in
327: -  out - the scatter to move them out

329:    Level: advanced

331: .keywords: PC, redundant solve
332: @*/
333: PetscErrorCode  PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
334: {
335:   PetscErrorCode ierr,(*f)(PC,VecScatter,VecScatter);

341:   PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetScatter_C",(void (**)(void))&f);
342:   if (f) {
343:     (*f)(pc,in,out);
344:   }
345:   return(0);
346: }

351: PetscErrorCode  PCRedundantGetPC_Redundant(PC pc,PC *innerpc)
352: {
354:   PC_Redundant   *red = (PC_Redundant*)pc->data;
355:   MPI_Comm       comm,subcomm;
356:   const char     *prefix;

359:   if (!red->psubcomm) {
360:     PetscObjectGetComm((PetscObject)pc,&comm);
361:     PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);
362:     PetscLogObjectMemory(pc,sizeof(PetscSubcomm));

364:     /* create a new PC that processors in each subcomm have copy of */
365:     subcomm = red->psubcomm->comm;
366:     KSPCreate(subcomm,&red->ksp);
367:     PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);
368:     PetscLogObjectParent(pc,red->ksp);
369:     KSPSetType(red->ksp,KSPPREONLY);
370:     KSPGetPC(red->ksp,&red->pc);
371:     PCSetType(red->pc,PCLU);

373:     PCGetOptionsPrefix(pc,&prefix);
374:     KSPSetOptionsPrefix(red->ksp,prefix);
375:     KSPAppendOptionsPrefix(red->ksp,"redundant_");
376:   }

378:   KSPGetPC(red->ksp,innerpc);
379:   return(0);
380: }

385: /*@
386:    PCRedundantGetPC - Gets the sequential PC created by the redundant PC.

388:    Not Collective

390:    Input Parameter:
391: .  pc - the preconditioner context

393:    Output Parameter:
394: .  innerpc - the sequential PC 

396:    Level: advanced

398: .keywords: PC, redundant solve
399: @*/
400: PetscErrorCode  PCRedundantGetPC(PC pc,PC *innerpc)
401: {
402:   PetscErrorCode ierr,(*f)(PC,PC*);

407:   PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetPC_C",(void (**)(void))&f);
408:   if (f) {
409:     (*f)(pc,innerpc);
410:   }
411:   return(0);
412: }

417: PetscErrorCode  PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
418: {
419:   PC_Redundant *red = (PC_Redundant*)pc->data;

422:   if (mat)  *mat  = red->pmats;
423:   if (pmat) *pmat = red->pmats;
424:   return(0);
425: }

430: /*@
431:    PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix

433:    Not Collective

435:    Input Parameter:
436: .  pc - the preconditioner context

438:    Output Parameters:
439: +  mat - the matrix
440: -  pmat - the (possibly different) preconditioner matrix

442:    Level: advanced

444: .keywords: PC, redundant solve
445: @*/
446: PetscErrorCode  PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
447: {
448:   PetscErrorCode ierr,(*f)(PC,Mat*,Mat*);

454:   PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetOperators_C",(void (**)(void))&f);
455:   if (f) {
456:     (*f)(pc,mat,pmat);
457:   }
458:   return(0);
459: }

461: /* -------------------------------------------------------------------------------------*/
462: /*MC
463:      PCREDUNDANT - Runs a preconditioner for the entire problem on subgroups of processors

465:      Options for the redundant preconditioners can be set with -redundant_pc_xxx

467:   Options Database:
468: .  -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and
469:                               use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes.

471:    Level: intermediate

473: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(),
474:            PCRedundantGetPC(), PCRedundantGetOperators(), PCRedundantSetNumber()
475: M*/

480: PetscErrorCode  PCCreate_Redundant(PC pc)
481: {
483:   PC_Redundant   *red;
484:   PetscMPIInt    size;
485: 
487:   PetscNewLog(pc,PC_Redundant,&red);
488:   MPI_Comm_size(((PetscObject)pc)->comm,&size);
489:   red->nsubcomm       = size;
490:   red->useparallelmat = PETSC_TRUE;
491:   pc->data            = (void*)red;

493:   pc->ops->apply           = PCApply_Redundant;
494:   pc->ops->applytranspose  = 0;
495:   pc->ops->setup           = PCSetUp_Redundant;
496:   pc->ops->destroy         = PCDestroy_Redundant;
497:   pc->ops->setfromoptions  = PCSetFromOptions_Redundant;
498:   pc->ops->view            = PCView_Redundant;
499:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant",
500:                     PCRedundantSetScatter_Redundant);
501:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant",
502:                     PCRedundantSetNumber_Redundant);
503:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant",
504:                     PCRedundantGetPC_Redundant);
505:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant",
506:                     PCRedundantGetOperators_Redundant);
507:   return(0);
508: }