Actual source code: composite.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: /*
  3:       Defines a preconditioner that can consist of a collection of PCs
  4: */
  5: #include <petsc/private/pcimpl.h>
  6: #include <petscksp.h>            /*I "petscksp.h" I*/

  8: typedef struct _PC_CompositeLink *PC_CompositeLink;
  9: struct _PC_CompositeLink {
 10:   PC               pc;
 11:   PC_CompositeLink next;
 12:   PC_CompositeLink previous;
 13: };

 15: typedef struct {
 16:   PC_CompositeLink head;
 17:   PCCompositeType  type;
 18:   Vec              work1;
 19:   Vec              work2;
 20:   PetscScalar      alpha;
 21: } PC_Composite;

 25: static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y)
 26: {
 27:   PetscErrorCode   ierr;
 28:   PC_Composite     *jac = (PC_Composite*)pc->data;
 29:   PC_CompositeLink next = jac->head;
 30:   Mat              mat  = pc->pmat;


 34:   if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");

 36:   /* Set the reuse flag on children PCs */
 37:   while (next) {
 38:     PCSetReusePreconditioner(next->pc,pc->reusepreconditioner);
 39:     next = next->next;
 40:   }
 41:   next = jac->head;

 43:   if (next->next && !jac->work2) { /* allocate second work vector */
 44:     VecDuplicate(jac->work1,&jac->work2);
 45:   }
 46:   if (pc->useAmat) mat = pc->mat;
 47:   PCApply(next->pc,x,y);                      /* y <- B x */
 48:   while (next->next) {
 49:     next = next->next;
 50:     MatMult(mat,y,jac->work1);                /* work1 <- A y */
 51:     VecWAXPY(jac->work2,-1.0,jac->work1,x);   /* work2 <- x - work1 */
 52:     VecSet(jac->work1,0.0);  /* zero since some PC's may not set all entries in the result */
 53:     PCApply(next->pc,jac->work2,jac->work1);  /* work1 <- C work2 */
 54:     VecAXPY(y,1.0,jac->work1);                /* y <- y + work1 = B x + C (x - A B x) = (B + C (1 - A B)) x */
 55:   }
 56:   if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
 57:     while (next->previous) {
 58:       next = next->previous;
 59:       MatMult(mat,y,jac->work1);
 60:       VecWAXPY(jac->work2,-1.0,jac->work1,x);
 61:       VecSet(jac->work1,0.0);  /* zero since some PC's may not set all entries in the result */
 62:       PCApply(next->pc,jac->work2,jac->work1);
 63:       VecAXPY(y,1.0,jac->work1);
 64:     }
 65:   }
 66:   return(0);
 67: }

 71: static PetscErrorCode PCApplyTranspose_Composite_Multiplicative(PC pc,Vec x,Vec y)
 72: {
 73:   PetscErrorCode   ierr;
 74:   PC_Composite     *jac = (PC_Composite*)pc->data;
 75:   PC_CompositeLink next = jac->head;
 76:   Mat              mat  = pc->pmat;

 79:   if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
 80:   if (next->next && !jac->work2) { /* allocate second work vector */
 81:     VecDuplicate(jac->work1,&jac->work2);
 82:   }
 83:   if (pc->useAmat) mat = pc->mat;
 84:   /* locate last PC */
 85:   while (next->next) {
 86:     next = next->next;
 87:   }
 88:   PCApplyTranspose(next->pc,x,y);
 89:   while (next->previous) {
 90:     next = next->previous;
 91:     MatMultTranspose(mat,y,jac->work1);
 92:     VecWAXPY(jac->work2,-1.0,jac->work1,x);
 93:     VecSet(jac->work1,0.0);  /* zero since some PC's may not set all entries in the result */
 94:     PCApplyTranspose(next->pc,jac->work2,jac->work1);
 95:     VecAXPY(y,1.0,jac->work1);
 96:   }
 97:   if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
 98:     next = jac->head;
 99:     while (next->next) {
100:       next = next->next;
101:       MatMultTranspose(mat,y,jac->work1);
102:       VecWAXPY(jac->work2,-1.0,jac->work1,x);
103:       VecSet(jac->work1,0.0);  /* zero since some PC's may not set all entries in the result */
104:       PCApplyTranspose(next->pc,jac->work2,jac->work1);
105:       VecAXPY(y,1.0,jac->work1);
106:     }
107:   }
108:   return(0);
109: }

111: /*
112:     This is very special for a matrix of the form alpha I + R + S
113: where first preconditioner is built from alpha I + S and second from
114: alpha I + R
115: */
118: static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y)
119: {
120:   PetscErrorCode   ierr;
121:   PC_Composite     *jac = (PC_Composite*)pc->data;
122:   PC_CompositeLink next = jac->head;

125:   if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
126:   if (!next->next || next->next->next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs");

128:   /* Set the reuse flag on children PCs */
129:   PCSetReusePreconditioner(next->pc,pc->reusepreconditioner);
130:   PCSetReusePreconditioner(next->next->pc,pc->reusepreconditioner);

132:   PCApply(next->pc,x,jac->work1);
133:   PCApply(next->next->pc,jac->work1,y);
134:   return(0);
135: }

139: static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y)
140: {
141:   PetscErrorCode   ierr;
142:   PC_Composite     *jac = (PC_Composite*)pc->data;
143:   PC_CompositeLink next = jac->head;

146:   if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");

148:   /* Set the reuse flag on children PCs */
149:   while (next) {
150:     PCSetReusePreconditioner(next->pc,pc->reusepreconditioner);
151:     next = next->next;
152:   }
153:   next = jac->head;

155:   PCApply(next->pc,x,y);
156:   while (next->next) {
157:     next = next->next;
158:     VecSet(jac->work1,0.0);  /* zero since some PC's may not set all entries in the result */
159:     PCApply(next->pc,x,jac->work1);
160:     VecAXPY(y,1.0,jac->work1);
161:   }
162:   return(0);
163: }

167: static PetscErrorCode PCApplyTranspose_Composite_Additive(PC pc,Vec x,Vec y)
168: {
169:   PetscErrorCode   ierr;
170:   PC_Composite     *jac = (PC_Composite*)pc->data;
171:   PC_CompositeLink next = jac->head;

174:   if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
175:   PCApplyTranspose(next->pc,x,y);
176:   while (next->next) {
177:     next = next->next;
178:     VecSet(jac->work1,0.0);  /* zero since some PC's may not set all entries in the result */
179:     PCApplyTranspose(next->pc,x,jac->work1);
180:     VecAXPY(y,1.0,jac->work1);
181:   }
182:   return(0);
183: }

187: static PetscErrorCode PCSetUp_Composite(PC pc)
188: {
189:   PetscErrorCode   ierr;
190:   PC_Composite     *jac = (PC_Composite*)pc->data;
191:   PC_CompositeLink next = jac->head;

194:   if (!jac->work1) {
195:     MatCreateVecs(pc->pmat,&jac->work1,0);
196:   }
197:   while (next) {
198:     PCSetOperators(next->pc,pc->mat,pc->pmat);
199:     next = next->next;
200:   }
201:   return(0);
202: }

206: static PetscErrorCode PCReset_Composite(PC pc)
207: {
208:   PC_Composite     *jac = (PC_Composite*)pc->data;
209:   PetscErrorCode   ierr;
210:   PC_CompositeLink next = jac->head;

213:   while (next) {
214:     PCReset(next->pc);
215:     next = next->next;
216:   }
217:   VecDestroy(&jac->work1);
218:   VecDestroy(&jac->work2);
219:   return(0);
220: }

224: static PetscErrorCode PCDestroy_Composite(PC pc)
225: {
226:   PC_Composite     *jac = (PC_Composite*)pc->data;
227:   PetscErrorCode   ierr;
228:   PC_CompositeLink next = jac->head,next_tmp;

231:   PCReset_Composite(pc);
232:   while (next) {
233:     PCDestroy(&next->pc);
234:     next_tmp = next;
235:     next     = next->next;
236:     PetscFree(next_tmp);
237:   }
238:   PetscFree(pc->data);
239:   return(0);
240: }

244: static PetscErrorCode PCSetFromOptions_Composite(PetscOptions *PetscOptionsObject,PC pc)
245: {
246:   PC_Composite     *jac = (PC_Composite*)pc->data;
247:   PetscErrorCode   ierr;
248:   PetscInt         nmax = 8,i;
249:   PC_CompositeLink next;
250:   char             *pcs[8];
251:   PetscBool        flg;

254:   PetscOptionsHead(PetscOptionsObject,"Composite preconditioner options");
255:   PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);
256:   if (flg) {
257:     PCCompositeSetType(pc,jac->type);
258:   }
259:   PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);
260:   if (flg) {
261:     for (i=0; i<nmax; i++) {
262:       PCCompositeAddPC(pc,pcs[i]);
263:       PetscFree(pcs[i]);   /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */
264:     }
265:   }
266:   PetscOptionsTail();

268:   next = jac->head;
269:   while (next) {
270:     PCSetFromOptions(next->pc);
271:     next = next->next;
272:   }
273:   return(0);
274: }

278: static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer)
279: {
280:   PC_Composite     *jac = (PC_Composite*)pc->data;
281:   PetscErrorCode   ierr;
282:   PC_CompositeLink next = jac->head;
283:   PetscBool        iascii;

286:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
287:   if (iascii) {
288:     PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);
289:     PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");
290:     PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
291:   }
292:   if (iascii) {
293:     PetscViewerASCIIPushTab(viewer);
294:   }
295:   while (next) {
296:     PCView(next->pc,viewer);
297:     next = next->next;
298:   }
299:   if (iascii) {
300:     PetscViewerASCIIPopTab(viewer);
301:     PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
302:   }
303:   return(0);
304: }

306: /* ------------------------------------------------------------------------------*/

310: static PetscErrorCode  PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
311: {
312:   PC_Composite *jac = (PC_Composite*)pc->data;

315:   jac->alpha = alpha;
316:   return(0);
317: }

321: static PetscErrorCode  PCCompositeSetType_Composite(PC pc,PCCompositeType type)
322: {
323:   PC_Composite *jac = (PC_Composite*)pc->data;

326:   if (type == PC_COMPOSITE_ADDITIVE) {
327:     pc->ops->apply          = PCApply_Composite_Additive;
328:     pc->ops->applytranspose = PCApplyTranspose_Composite_Additive;
329:   } else if (type ==  PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
330:     pc->ops->apply          = PCApply_Composite_Multiplicative;
331:     pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative;
332:   } else if (type ==  PC_COMPOSITE_SPECIAL) {
333:     pc->ops->apply          = PCApply_Composite_Special;
334:     pc->ops->applytranspose = NULL;
335:   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type");
336:   jac->type = type;
337:   return(0);
338: }

342: static PetscErrorCode  PCCompositeGetType_Composite(PC pc,PCCompositeType *type)
343: {
344:   PC_Composite *jac = (PC_Composite*)pc->data;

347:   *type = jac->type;
348:   return(0);
349: }

353: static PetscErrorCode  PCCompositeAddPC_Composite(PC pc,PCType type)
354: {
355:   PC_Composite     *jac;
356:   PC_CompositeLink next,ilink;
357:   PetscErrorCode   ierr;
358:   PetscInt         cnt = 0;
359:   const char       *prefix;
360:   char             newprefix[8];

363:   PetscNewLog(pc,&ilink);
364:   ilink->next = 0;
365:   PCCreate(PetscObjectComm((PetscObject)pc),&ilink->pc);
366:   PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->pc);

368:   jac  = (PC_Composite*)pc->data;
369:   next = jac->head;
370:   if (!next) {
371:     jac->head       = ilink;
372:     ilink->previous = NULL;
373:   } else {
374:     cnt++;
375:     while (next->next) {
376:       next = next->next;
377:       cnt++;
378:     }
379:     next->next      = ilink;
380:     ilink->previous = next;
381:   }
382:   PCGetOptionsPrefix(pc,&prefix);
383:   PCSetOptionsPrefix(ilink->pc,prefix);
384:   sprintf(newprefix,"sub_%d_",(int)cnt);
385:   PCAppendOptionsPrefix(ilink->pc,newprefix);
386:   /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
387:   PCSetType(ilink->pc,type);
388:   return(0);
389: }

393: static PetscErrorCode  PCCompositeGetNumberPC_Composite(PC pc,PetscInt *n)
394: {
395:   PC_Composite     *jac;
396:   PC_CompositeLink next;

399:   jac  = (PC_Composite*)pc->data;
400:   next = jac->head;
401:   *n = 0;
402:   while (next) {
403:     next = next->next;
404:     (*n) ++;
405:   }
406:   return(0);
407: }

411: static PetscErrorCode  PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
412: {
413:   PC_Composite     *jac;
414:   PC_CompositeLink next;
415:   PetscInt         i;

418:   jac  = (PC_Composite*)pc->data;
419:   next = jac->head;
420:   for (i=0; i<n; i++) {
421:     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
422:     next = next->next;
423:   }
424:   *subpc = next->pc;
425:   return(0);
426: }

428: /* -------------------------------------------------------------------------------- */
431: /*@
432:    PCCompositeSetType - Sets the type of composite preconditioner.

434:    Logically Collective on PC

436:    Input Parameters:
437: +  pc - the preconditioner context
438: -  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL

440:    Options Database Key:
441: .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type

443:    Level: Developer

445: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
446: @*/
447: PetscErrorCode  PCCompositeSetType(PC pc,PCCompositeType type)
448: {

454:   PetscTryMethod(pc,"PCCompositeSetType_C",(PC,PCCompositeType),(pc,type));
455:   return(0);
456: }

460: /*@
461:    PCCompositeGetType - Gets the type of composite preconditioner.

463:    Logically Collective on PC

465:    Input Parameter:
466: .  pc - the preconditioner context

468:    Output Parameter:
469: .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL

471:    Options Database Key:
472: .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type

474:    Level: Developer

476: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
477: @*/
478: PetscErrorCode  PCCompositeGetType(PC pc,PCCompositeType *type)
479: {

484:   PetscUseMethod(pc,"PCCompositeGetType_C",(PC,PCCompositeType*),(pc,type));
485:   return(0);
486: }

490: /*@
491:    PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
492:      for alphaI + R + S

494:    Logically Collective on PC

496:    Input Parameter:
497: +  pc - the preconditioner context
498: -  alpha - scale on identity

500:    Level: Developer

502: .keywords: PC, set, type, composite preconditioner, additive, multiplicative
503: @*/
504: PetscErrorCode  PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
505: {

511:   PetscTryMethod(pc,"PCCompositeSpecialSetAlpha_C",(PC,PetscScalar),(pc,alpha));
512:   return(0);
513: }

517: /*@C
518:    PCCompositeAddPC - Adds another PC to the composite PC.

520:    Collective on PC

522:    Input Parameters:
523: +  pc - the preconditioner context
524: -  type - the type of the new preconditioner

526:    Level: Developer

528: .keywords: PC, composite preconditioner, add
529: @*/
530: PetscErrorCode  PCCompositeAddPC(PC pc,PCType type)
531: {

536:   PetscTryMethod(pc,"PCCompositeAddPC_C",(PC,PCType),(pc,type));
537:   return(0);
538: }

542: /*@
543:    PCCompositeGetNumberPC - Gets the number of PC objects in the composite PC.

545:    Not Collective

547:    Input Parameter:
548: .  pc - the preconditioner context

550:    Output Parameter:
551: .  num - the number of sub pcs

553:    Level: Developer

555: .keywords: PC, get, composite preconditioner, sub preconditioner

557: .seealso: PCCompositeGetPC()
558: @*/
559: PetscErrorCode  PCCompositeGetNumberPC(PC pc,PetscInt *num)
560: {

566:   PetscUseMethod(pc,"PCCompositeGetNumberPC_C",(PC,PetscInt*),(pc,num));
567:   return(0);
568: }

572: /*@
573:    PCCompositeGetPC - Gets one of the PC objects in the composite PC.

575:    Not Collective

577:    Input Parameter:
578: +  pc - the preconditioner context
579: -  n - the number of the pc requested

581:    Output Parameters:
582: .  subpc - the PC requested

584:    Level: Developer

586: .keywords: PC, get, composite preconditioner, sub preconditioner

588: .seealso: PCCompositeAddPC(), PCCompositeGetNumberPC()
589: @*/
590: PetscErrorCode  PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
591: {

597:   PetscUseMethod(pc,"PCCompositeGetPC_C",(PC,PetscInt,PC*),(pc,n,subpc));
598:   return(0);
599: }

601: /* -------------------------------------------------------------------------------------------*/

603: /*MC
604:      PCCOMPOSITE - Build a preconditioner by composing together several preconditioners

606:    Options Database Keys:
607: +  -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
608: .  -pc_use_amat - Activates PCSetUseAmat()
609: -  -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose

611:    Level: intermediate

613:    Concepts: composing solvers

615:    Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
616:           inner PCs to be PCKSP.
617:           Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
618:           the incorrect answer) unless you use KSPFGMRES as the outer Krylov method


621: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
622:            PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
623:            PCCompositeGetPC(), PCSetUseAmat()

625: M*/

629: PETSC_EXTERN PetscErrorCode PCCreate_Composite(PC pc)
630: {
632:   PC_Composite   *jac;

635:   PetscNewLog(pc,&jac);

637:   pc->ops->apply           = PCApply_Composite_Additive;
638:   pc->ops->applytranspose  = PCApplyTranspose_Composite_Additive;
639:   pc->ops->setup           = PCSetUp_Composite;
640:   pc->ops->reset           = PCReset_Composite;
641:   pc->ops->destroy         = PCDestroy_Composite;
642:   pc->ops->setfromoptions  = PCSetFromOptions_Composite;
643:   pc->ops->view            = PCView_Composite;
644:   pc->ops->applyrichardson = 0;

646:   pc->data   = (void*)jac;
647:   jac->type  = PC_COMPOSITE_ADDITIVE;
648:   jac->work1 = 0;
649:   jac->work2 = 0;
650:   jac->head  = 0;

652:   PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSetType_C",PCCompositeSetType_Composite);
653:   PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetType_C",PCCompositeGetType_Composite);
654:   PetscObjectComposeFunction((PetscObject)pc,"PCCompositeAddPC_C",PCCompositeAddPC_Composite);
655:   PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetNumberPC_C",PCCompositeGetNumberPC_Composite);
656:   PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetPC_C",PCCompositeGetPC_Composite);
657:   PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",PCCompositeSpecialSetAlpha_Composite);
658:   return(0);
659: }