Actual source code: precon.c

  2: /*
  3:     The PC (preconditioner) interface routines, callable by users.
  4: */
  5: #include <private/pcimpl.h>            /*I "petscksp.h" I*/

  7: /* Logging support */
  8: PetscClassId   PC_CLASSID;
  9: PetscLogEvent  PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
 10: PetscLogEvent  PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks;

 14: PetscErrorCode PCGetDefaultType_Private(PC pc,const char* type[])
 15: {
 17:   PetscMPIInt    size;
 18:   PetscBool      flg1,flg2,set,flg3;

 21:   MPI_Comm_size(((PetscObject)pc)->comm,&size);
 22:   if (pc->pmat) {
 23:     PetscErrorCode (*f)(Mat,MatReuse,Mat*);
 24:     PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",(void (**)(void))&f);
 25:     if (size == 1) {
 26:       MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ICC,&flg1);
 27:       MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ILU,&flg2);
 28:       MatIsSymmetricKnown(pc->pmat,&set,&flg3);
 29:       if (flg1 && (!flg2 || (set && flg3))) {
 30:         *type = PCICC;
 31:       } else if (flg2) {
 32:         *type = PCILU;
 33:       } else if (f) { /* likely is a parallel matrix run on one processor */
 34:         *type = PCBJACOBI;
 35:       } else {
 36:         *type = PCNONE;
 37:       }
 38:     } else {
 39:        if (f) {
 40:         *type = PCBJACOBI;
 41:       } else {
 42:         *type = PCNONE;
 43:       }
 44:     }
 45:   } else {
 46:     if (size == 1) {
 47:       *type = PCILU;
 48:     } else {
 49:       *type = PCBJACOBI;
 50:     }
 51:   }
 52:   return(0);
 53: }

 57: /*@
 58:    PCReset - Resets a PC context to the pcsetupcalled = 0 state and removes any allocated Vecs and Mats

 60:    Collective on PC

 62:    Input Parameter:
 63: .  pc - the preconditioner context

 65:    Level: developer

 67:    Notes: This allows a PC to be reused for a different sized linear system but using the same options that have been previously set in the PC

 69: .keywords: PC, destroy

 71: .seealso: PCCreate(), PCSetUp()
 72: @*/
 73: PetscErrorCode  PCReset(PC pc)
 74: {

 79:   if (pc->ops->reset) {
 80:     (*pc->ops->reset)(pc);
 81:   }
 82:   VecDestroy(&pc->diagonalscaleright);
 83:   VecDestroy(&pc->diagonalscaleleft);
 84:   MatDestroy(&pc->pmat);
 85:   MatDestroy(&pc->mat);
 86:   pc->setupcalled = 0;
 87:   return(0);
 88: }

 92: /*@
 93:    PCDestroy - Destroys PC context that was created with PCCreate().

 95:    Collective on PC

 97:    Input Parameter:
 98: .  pc - the preconditioner context

100:    Level: developer

102: .keywords: PC, destroy

104: .seealso: PCCreate(), PCSetUp()
105: @*/
106: PetscErrorCode  PCDestroy(PC *pc)
107: {

111:   if (!*pc) return(0);
113:   if (--((PetscObject)(*pc))->refct > 0) {*pc = 0; return(0);}

115:   PCReset(*pc);

117:   /* if memory was published with AMS then destroy it */
118:   PetscObjectDepublish((*pc));
119:   if ((*pc)->ops->destroy) {(*(*pc)->ops->destroy)((*pc));}
120:   DMDestroy(&(*pc)->dm);
121:   PetscHeaderDestroy(pc);
122:   return(0);
123: }

127: /*@C
128:    PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
129:       scaling as needed by certain time-stepping codes.

131:    Logically Collective on PC

133:    Input Parameter:
134: .  pc - the preconditioner context

136:    Output Parameter:
137: .  flag - PETSC_TRUE if it applies the scaling

139:    Level: developer

141:    Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
142: $           D M A D^{-1} y = D M b  for left preconditioning or
143: $           D A M D^{-1} z = D b for right preconditioning

145: .keywords: PC

147: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCSetDiagonalScale()
148: @*/
149: PetscErrorCode  PCGetDiagonalScale(PC pc,PetscBool  *flag)
150: {
154:   *flag = pc->diagonalscale;
155:   return(0);
156: }

160: /*@
161:    PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
162:       scaling as needed by certain time-stepping codes.

164:    Logically Collective on PC

166:    Input Parameters:
167: +  pc - the preconditioner context
168: -  s - scaling vector

170:    Level: intermediate

172:    Notes: The system solved via the Krylov method is
173: $           D M A D^{-1} y = D M b  for left preconditioning or
174: $           D A M D^{-1} z = D b for right preconditioning

176:    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.

178: .keywords: PC

180: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCGetDiagonalScale()
181: @*/
182: PetscErrorCode  PCSetDiagonalScale(PC pc,Vec s)
183: {

189:   pc->diagonalscale     = PETSC_TRUE;
190:   PetscObjectReference((PetscObject)s);
191:   VecDestroy(&pc->diagonalscaleleft);
192:   pc->diagonalscaleleft = s;
193:   VecDuplicate(s,&pc->diagonalscaleright);
194:   VecCopy(s,pc->diagonalscaleright);
195:   VecReciprocal(pc->diagonalscaleright);
196:   return(0);
197: }

201: /*@
202:    PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.

204:    Logically Collective on PC

206:    Input Parameters:
207: +  pc - the preconditioner context
208: .  in - input vector
209: +  out - scaled vector (maybe the same as in)

211:    Level: intermediate

213:    Notes: The system solved via the Krylov method is
214: $           D M A D^{-1} y = D M b  for left preconditioning or
215: $           D A M D^{-1} z = D b for right preconditioning

217:    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.

219:    If diagonal scaling is turned off and in is not out then in is copied to out

221: .keywords: PC

223: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
224: @*/
225: PetscErrorCode  PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
226: {

233:   if (pc->diagonalscale) {
234:     VecPointwiseMult(out,pc->diagonalscaleleft,in);
235:   } else if (in != out) {
236:     VecCopy(in,out);
237:   }
238:   return(0);
239: }

243: /*@
244:    PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.

246:    Logically Collective on PC

248:    Input Parameters:
249: +  pc - the preconditioner context
250: .  in - input vector
251: +  out - scaled vector (maybe the same as in)

253:    Level: intermediate

255:    Notes: The system solved via the Krylov method is
256: $           D M A D^{-1} y = D M b  for left preconditioning or
257: $           D A M D^{-1} z = D b for right preconditioning

259:    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.

261:    If diagonal scaling is turned off and in is not out then in is copied to out

263: .keywords: PC

265: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
266: @*/
267: PetscErrorCode  PCDiagonalScaleRight(PC pc,Vec in,Vec out)
268: {

275:   if (pc->diagonalscale) {
276:     VecPointwiseMult(out,pc->diagonalscaleright,in);
277:   } else if (in != out) {
278:     VecCopy(in,out);
279:   }
280:   return(0);
281: }

283: #if 0
286: static PetscErrorCode PCPublish_Petsc(PetscObject obj)
287: {
289:   return(0);
290: }
291: #endif

295: /*@
296:    PCCreate - Creates a preconditioner context.

298:    Collective on MPI_Comm

300:    Input Parameter:
301: .  comm - MPI communicator 

303:    Output Parameter:
304: .  pc - location to put the preconditioner context

306:    Notes:
307:    The default preconditioner for sparse matrices is PCILU or PCICC with 0 fill on one process and block Jacobi with PCILU or ICC 
308:    in parallel. For dense matrices it is always PCNONE.

310:    Level: developer

312: .keywords: PC, create, context

314: .seealso: PCSetUp(), PCApply(), PCDestroy()
315: @*/
316: PetscErrorCode  PCCreate(MPI_Comm comm,PC *newpc)
317: {
318:   PC             pc;

323:   *newpc = 0;
324: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
325:   PCInitializePackage(PETSC_NULL);
326: #endif

328:   PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_CLASSID,-1,"PC","Preconditioner","PC",comm,PCDestroy,PCView);

330:   pc->mat                  = 0;
331:   pc->pmat                 = 0;
332:   pc->setupcalled          = 0;
333:   pc->setfromoptionscalled = 0;
334:   pc->data                 = 0;
335:   pc->diagonalscale        = PETSC_FALSE;
336:   pc->diagonalscaleleft    = 0;
337:   pc->diagonalscaleright   = 0;
338:   pc->reuse                = 0;

340:   pc->modifysubmatrices   = 0;
341:   pc->modifysubmatricesP  = 0;
342:   *newpc = pc;
343:   return(0);

345: }

347: /* -------------------------------------------------------------------------------*/

351: /*@
352:    PCApply - Applies the preconditioner to a vector.

354:    Collective on PC and Vec

356:    Input Parameters:
357: +  pc - the preconditioner context
358: -  x - input vector

360:    Output Parameter:
361: .  y - output vector

363:    Level: developer

365: .keywords: PC, apply

367: .seealso: PCApplyTranspose(), PCApplyBAorAB()
368: @*/
369: PetscErrorCode  PCApply(PC pc,Vec x,Vec y)
370: {

377:   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
378:   if (pc->setupcalled < 2) {
379:     PCSetUp(pc);
380:   }
381:   if (!pc->ops->apply) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have apply");
382:   PetscLogEventBegin(PC_Apply,pc,x,y,0);
383:   (*pc->ops->apply)(pc,x,y);
384:   PetscLogEventEnd(PC_Apply,pc,x,y,0);
385:   return(0);
386: }

390: /*@
391:    PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.

393:    Collective on PC and Vec

395:    Input Parameters:
396: +  pc - the preconditioner context
397: -  x - input vector

399:    Output Parameter:
400: .  y - output vector

402:    Notes:
403:    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.

405:    Level: developer

407: .keywords: PC, apply, symmetric, left

409: .seealso: PCApply(), PCApplySymmetricRight()
410: @*/
411: PetscErrorCode  PCApplySymmetricLeft(PC pc,Vec x,Vec y)
412: {

419:   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
420:   if (pc->setupcalled < 2) {
421:     PCSetUp(pc);
422:   }
423:   if (!pc->ops->applysymmetricleft) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have left symmetric apply");
424:   PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);
425:   (*pc->ops->applysymmetricleft)(pc,x,y);
426:   PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);
427:   return(0);
428: }

432: /*@
433:    PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.

435:    Collective on PC and Vec

437:    Input Parameters:
438: +  pc - the preconditioner context
439: -  x - input vector

441:    Output Parameter:
442: .  y - output vector

444:    Level: developer

446:    Notes:
447:    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.

449: .keywords: PC, apply, symmetric, right

451: .seealso: PCApply(), PCApplySymmetricLeft()
452: @*/
453: PetscErrorCode  PCApplySymmetricRight(PC pc,Vec x,Vec y)
454: {

461:   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
462:   if (pc->setupcalled < 2) {
463:     PCSetUp(pc);
464:   }
465:   if (!pc->ops->applysymmetricright) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have left symmetric apply");
466:   PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);
467:   (*pc->ops->applysymmetricright)(pc,x,y);
468:   PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);
469:   return(0);
470: }

474: /*@
475:    PCApplyTranspose - Applies the transpose of preconditioner to a vector.

477:    Collective on PC and Vec

479:    Input Parameters:
480: +  pc - the preconditioner context
481: -  x - input vector

483:    Output Parameter:
484: .  y - output vector

486:    Notes: For complex numbers this applies the non-Hermitian transpose.

488:    Developer Notes: We need to implement a PCApplyHermitianTranspose()

490:    Level: developer

492: .keywords: PC, apply, transpose

494: .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCApplyTransposeExists()
495: @*/
496: PetscErrorCode  PCApplyTranspose(PC pc,Vec x,Vec y)
497: {

504:   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
505:   if (pc->setupcalled < 2) {
506:     PCSetUp(pc);
507:   }
508:   if (!pc->ops->applytranspose) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have apply transpose");
509:   PetscLogEventBegin(PC_Apply,pc,x,y,0);
510:   (*pc->ops->applytranspose)(pc,x,y);
511:   PetscLogEventEnd(PC_Apply,pc,x,y,0);
512:   return(0);
513: }

517: /*@
518:    PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation

520:    Collective on PC and Vec

522:    Input Parameters:
523: .  pc - the preconditioner context

525:    Output Parameter:
526: .  flg - PETSC_TRUE if a transpose operation is defined

528:    Level: developer

530: .keywords: PC, apply, transpose

532: .seealso: PCApplyTranspose()
533: @*/
534: PetscErrorCode  PCApplyTransposeExists(PC pc,PetscBool  *flg)
535: {
539:   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
540:   else                         *flg = PETSC_FALSE;
541:   return(0);
542: }

546: /*@
547:    PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x.

549:    Collective on PC and Vec

551:    Input Parameters:
552: +  pc - the preconditioner context
553: .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
554: .  x - input vector
555: -  work - work vector

557:    Output Parameter:
558: .  y - output vector

560:    Level: developer

562:    Notes: If the PC has had PCSetDiagonalScale() set then D M A D^{-1} for left preconditioning or  D A M D^{-1} is actually applied. Note that the 
563:    specific KSPSolve() method must also be written to handle the post-solve "correction" for the diagonal scaling.

565: .keywords: PC, apply, operator

567: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
568: @*/
569: PetscErrorCode  PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
570: {

578:   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
579:   if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric");
580:   if (pc->diagonalscale && side == PC_SYMMETRIC) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");

582:   if (pc->setupcalled < 2) {
583:     PCSetUp(pc);
584:   }

586:   if (pc->diagonalscale) {
587:     if (pc->ops->applyBA) {
588:       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
589:       VecDuplicate(x,&work2);
590:       PCDiagonalScaleRight(pc,x,work2);
591:       (*pc->ops->applyBA)(pc,side,work2,y,work);
592:       PCDiagonalScaleLeft(pc,y,y);
593:       VecDestroy(&work2);
594:     } else if (side == PC_RIGHT) {
595:       PCDiagonalScaleRight(pc,x,y);
596:       PCApply(pc,y,work);
597:       MatMult(pc->mat,work,y);
598:       PCDiagonalScaleLeft(pc,y,y);
599:     } else if (side == PC_LEFT) {
600:       PCDiagonalScaleRight(pc,x,y);
601:       MatMult(pc->mat,y,work);
602:       PCApply(pc,work,y);
603:       PCDiagonalScaleLeft(pc,y,y);
604:     } else if (side == PC_SYMMETRIC) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
605:   } else {
606:     if (pc->ops->applyBA) {
607:       (*pc->ops->applyBA)(pc,side,x,y,work);
608:     } else if (side == PC_RIGHT) {
609:       PCApply(pc,x,work);
610:       MatMult(pc->mat,work,y);
611:     } else if (side == PC_LEFT) {
612:       MatMult(pc->mat,x,work);
613:       PCApply(pc,work,y);
614:     } else if (side == PC_SYMMETRIC) {
615:       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
616:       PCApplySymmetricRight(pc,x,work);
617:       MatMult(pc->mat,work,y);
618:       VecCopy(y,work);
619:       PCApplySymmetricLeft(pc,work,y);
620:     }
621:   }
622:   return(0);
623: }

627: /*@ 
628:    PCApplyBAorABTranspose - Applies the transpose of the preconditioner
629:    and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
630:    NOT tr(B*A) = tr(A)*tr(B).

632:    Collective on PC and Vec

634:    Input Parameters:
635: +  pc - the preconditioner context
636: .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
637: .  x - input vector
638: -  work - work vector

640:    Output Parameter:
641: .  y - output vector


644:    Notes: this routine is used internally so that the same Krylov code can be used to solve A x = b and A' x = b, with a preconditioner
645:       defined by B'. This is why this has the funny form that it computes tr(B) * tr(A) 
646:           
647:     Level: developer

649: .keywords: PC, apply, operator, transpose

651: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
652: @*/
653: PetscErrorCode  PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
654: {

662:   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
663:   if (pc->ops->applyBAtranspose) {
664:     (*pc->ops->applyBAtranspose)(pc,side,x,y,work);
665:     return(0);
666:   }
667:   if (side != PC_LEFT && side != PC_RIGHT) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");

669:   if (pc->setupcalled < 2) {
670:     PCSetUp(pc);
671:   }

673:   if (side == PC_RIGHT) {
674:     PCApplyTranspose(pc,x,work);
675:     MatMultTranspose(pc->mat,work,y);
676:   } else if (side == PC_LEFT) {
677:     MatMultTranspose(pc->mat,x,work);
678:     PCApplyTranspose(pc,work,y);
679:   }
680:   /* add support for PC_SYMMETRIC */
681:   return(0); /* actually will never get here */
682: }

684: /* -------------------------------------------------------------------------------*/

688: /*@
689:    PCApplyRichardsonExists - Determines whether a particular preconditioner has a 
690:    built-in fast application of Richardson's method.

692:    Not Collective

694:    Input Parameter:
695: .  pc - the preconditioner

697:    Output Parameter:
698: .  exists - PETSC_TRUE or PETSC_FALSE

700:    Level: developer

702: .keywords: PC, apply, Richardson, exists

704: .seealso: PCApplyRichardson()
705: @*/
706: PetscErrorCode  PCApplyRichardsonExists(PC pc,PetscBool  *exists)
707: {
711:   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
712:   else                          *exists = PETSC_FALSE;
713:   return(0);
714: }

718: /*@
719:    PCApplyRichardson - Applies several steps of Richardson iteration with 
720:    the particular preconditioner. This routine is usually used by the 
721:    Krylov solvers and not the application code directly.

723:    Collective on PC

725:    Input Parameters:
726: +  pc  - the preconditioner context
727: .  b   - the right hand side
728: .  w   - one work vector
729: .  rtol - relative decrease in residual norm convergence criteria
730: .  abstol - absolute residual norm convergence criteria
731: .  dtol - divergence residual norm increase criteria
732: .  its - the number of iterations to apply.
733: -  guesszero - if the input x contains nonzero initial guess

735:    Output Parameter:
736: +  outits - number of iterations actually used (for SOR this always equals its)
737: .  reason - the reason the apply terminated
738: -  y - the solution (also contains initial guess if guesszero is PETSC_FALSE

740:    Notes: 
741:    Most preconditioners do not support this function. Use the command
742:    PCApplyRichardsonExists() to determine if one does.

744:    Except for the multigrid PC this routine ignores the convergence tolerances
745:    and always runs for the number of iterations
746:  
747:    Level: developer

749: .keywords: PC, apply, Richardson

751: .seealso: PCApplyRichardsonExists()
752: @*/
753: PetscErrorCode  PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool  guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
754: {

762:   if (b == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"b and y must be different vectors");
763:   if (pc->setupcalled < 2) {
764:     PCSetUp(pc);
765:   }
766:   if (!pc->ops->applyrichardson) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have apply richardson");
767:   (*pc->ops->applyrichardson)(pc,b,y,w,rtol,abstol,dtol,its,guesszero,outits,reason);
768:   return(0);
769: }

771: /* 
772:       a setupcall of 0 indicates never setup, 
773:                      1 needs to be resetup,
774:                      2 does not need any changes.
775: */
778: /*@
779:    PCSetUp - Prepares for the use of a preconditioner.

781:    Collective on PC

783:    Input Parameter:
784: .  pc - the preconditioner context

786:    Level: developer

788: .keywords: PC, setup

790: .seealso: PCCreate(), PCApply(), PCDestroy()
791: @*/
792: PetscErrorCode  PCSetUp(PC pc)
793: {
795:   const char     *def;

799:   if (!pc->mat) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");

801:   if (pc->setupcalled > 1) {
802:     PetscInfo(pc,"Setting PC with identical preconditioner\n");
803:     return(0);
804:   } else if (!pc->setupcalled) {
805:     PetscInfo(pc,"Setting up new PC\n");
806:   } else if (pc->flag == SAME_NONZERO_PATTERN) {
807:     PetscInfo(pc,"Setting up PC with same nonzero pattern\n");
808:   } else {
809:     PetscInfo(pc,"Setting up PC with different nonzero pattern\n");
810:   }

812:   if (!((PetscObject)pc)->type_name) {
813:     PCGetDefaultType_Private(pc,&def);
814:     PCSetType(pc,def);
815:   }

817:   PetscLogEventBegin(PC_SetUp,pc,0,0,0);
818:   if (pc->ops->setup) {
819:     (*pc->ops->setup)(pc);
820:   }
821:   pc->setupcalled = 2;
822:   PetscLogEventEnd(PC_SetUp,pc,0,0,0);
823:   return(0);
824: }

828: /*@
829:    PCSetUpOnBlocks - Sets up the preconditioner for each block in
830:    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz 
831:    methods.

833:    Collective on PC

835:    Input Parameters:
836: .  pc - the preconditioner context

838:    Level: developer

840: .keywords: PC, setup, blocks

842: .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
843: @*/
844: PetscErrorCode  PCSetUpOnBlocks(PC pc)
845: {

850:   if (!pc->ops->setuponblocks) return(0);
851:   PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);
852:   (*pc->ops->setuponblocks)(pc);
853:   PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);
854:   return(0);
855: }

859: /*@C
860:    PCSetModifySubMatrices - Sets a user-defined routine for modifying the
861:    submatrices that arise within certain subdomain-based preconditioners.
862:    The basic submatrices are extracted from the preconditioner matrix as
863:    usual; the user can then alter these (for example, to set different boundary
864:    conditions for each submatrix) before they are used for the local solves.

866:    Logically Collective on PC

868:    Input Parameters:
869: +  pc - the preconditioner context
870: .  func - routine for modifying the submatrices
871: -  ctx - optional user-defined context (may be null)

873:    Calling sequence of func:
874: $     func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);

876: .  row - an array of index sets that contain the global row numbers
877:          that comprise each local submatrix
878: .  col - an array of index sets that contain the global column numbers
879:          that comprise each local submatrix
880: .  submat - array of local submatrices
881: -  ctx - optional user-defined context for private data for the 
882:          user-defined func routine (may be null)

884:    Notes:
885:    PCSetModifySubMatrices() MUST be called before KSPSetUp() and
886:    KSPSolve().

888:    A routine set by PCSetModifySubMatrices() is currently called within
889:    the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
890:    preconditioners.  All other preconditioners ignore this routine.

892:    Level: advanced

894: .keywords: PC, set, modify, submatrices

896: .seealso: PCModifySubMatrices(), PCASMGetSubMatrices()
897: @*/
898: PetscErrorCode  PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
899: {
902:   pc->modifysubmatrices  = func;
903:   pc->modifysubmatricesP = ctx;
904:   return(0);
905: }

909: /*@C
910:    PCModifySubMatrices - Calls an optional user-defined routine within 
911:    certain preconditioners if one has been set with PCSetModifySubMarices().

913:    Collective on PC

915:    Input Parameters:
916: +  pc - the preconditioner context
917: .  nsub - the number of local submatrices
918: .  row - an array of index sets that contain the global row numbers
919:          that comprise each local submatrix
920: .  col - an array of index sets that contain the global column numbers
921:          that comprise each local submatrix
922: .  submat - array of local submatrices
923: -  ctx - optional user-defined context for private data for the 
924:          user-defined routine (may be null)

926:    Output Parameter:
927: .  submat - array of local submatrices (the entries of which may
928:             have been modified)

930:    Notes:
931:    The user should NOT generally call this routine, as it will
932:    automatically be called within certain preconditioners (currently
933:    block Jacobi, additive Schwarz) if set.

935:    The basic submatrices are extracted from the preconditioner matrix
936:    as usual; the user can then alter these (for example, to set different
937:    boundary conditions for each submatrix) before they are used for the
938:    local solves.

940:    Level: developer

942: .keywords: PC, modify, submatrices

944: .seealso: PCSetModifySubMatrices()
945: @*/
946: PetscErrorCode  PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
947: {

952:   if (!pc->modifysubmatrices) return(0);
953:   PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);
954:   (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);
955:   PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);
956:   return(0);
957: }

961: /*@
962:    PCSetOperators - Sets the matrix associated with the linear system and 
963:    a (possibly) different one associated with the preconditioner.

965:    Logically Collective on PC and Mat

967:    Input Parameters:
968: +  pc - the preconditioner context
969: .  Amat - the matrix associated with the linear system
970: .  Pmat - the matrix to be used in constructing the preconditioner, usually
971:           the same as Amat. 
972: -  flag - flag indicating information about the preconditioner matrix structure
973:    during successive linear solves.  This flag is ignored the first time a
974:    linear system is solved, and thus is irrelevant when solving just one linear
975:    system.

977:    Notes: 
978:    The flag can be used to eliminate unnecessary work in the preconditioner 
979:    during the repeated solution of linear systems of the same size.  The 
980:    available options are
981: +    SAME_PRECONDITIONER -
982:        Pmat is identical during successive linear solves.
983:        This option is intended for folks who are using
984:        different Amat and Pmat matrices and wish to reuse the
985:        same preconditioner matrix.  For example, this option
986:        saves work by not recomputing incomplete factorization
987:        for ILU/ICC preconditioners.
988: .     SAME_NONZERO_PATTERN -
989:        Pmat has the same nonzero structure during
990:        successive linear solves. 
991: -     DIFFERENT_NONZERO_PATTERN -
992:        Pmat does not have the same nonzero structure.

994:     Passing a PETSC_NULL for Amat or Pmat removes the matrix that is currently used.

996:     If you wish to replace either Amat or Pmat but leave the other one untouched then
997:     first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
998:     on it and then pass it back in in your call to KSPSetOperators().

1000:    Caution:
1001:    If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
1002:    and does not check the structure of the matrix.  If you erroneously
1003:    claim that the structure is the same when it actually is not, the new
1004:    preconditioner will not function correctly.  Thus, use this optimization
1005:    feature carefully!

1007:    If in doubt about whether your preconditioner matrix has changed
1008:    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.

1010:    More Notes about Repeated Solution of Linear Systems:
1011:    PETSc does NOT reset the matrix entries of either Amat or Pmat
1012:    to zero after a linear solve; the user is completely responsible for
1013:    matrix assembly.  See the routine MatZeroEntries() if desiring to
1014:    zero all elements of a matrix.

1016:    Level: intermediate

1018: .keywords: PC, set, operators, matrix, linear system

1020: .seealso: PCGetOperators(), MatZeroEntries()
1021:  @*/
1022: PetscErrorCode  PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
1023: {
1025:   PetscInt       m1,n1,m2,n2;

1033:   if (pc->setupcalled && Amat && Pmat) {
1034:     MatGetLocalSize(Amat,&m1,&n1);
1035:     MatGetLocalSize(pc->mat,&m2,&n2);
1036:     if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot change local size of Amat after use old sizes %D %D new sizes %D %D",m2,n2,m1,n1);
1037:     MatGetLocalSize(Pmat,&m1,&n1);
1038:     MatGetLocalSize(pc->pmat,&m2,&n2);
1039:     if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot change local size of Pmat after use old sizes %D %D new sizes %D %D",m2,n2,m1,n1);
1040:   }

1042:   /* reference first in case the matrices are the same */
1043:   if (Amat) {PetscObjectReference((PetscObject)Amat);}
1044:   MatDestroy(&pc->mat);
1045:   if (Pmat) {PetscObjectReference((PetscObject)Pmat);}
1046:   MatDestroy(&pc->pmat);
1047:   pc->mat  = Amat;
1048:   pc->pmat = Pmat;

1050:   if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1051:     pc->setupcalled = 1;
1052:   }
1053:   pc->flag = flag;
1054:   return(0);
1055: }

1059: /*@C
1060:    PCGetOperators - Gets the matrix associated with the linear system and
1061:    possibly a different one associated with the preconditioner.

1063:    Not collective, though parallel Mats are returned if the PC is parallel

1065:    Input Parameter:
1066: .  pc - the preconditioner context

1068:    Output Parameters:
1069: +  mat - the matrix associated with the linear system
1070: .  pmat - matrix associated with the preconditioner, usually the same
1071:           as mat. 
1072: -  flag - flag indicating information about the preconditioner
1073:           matrix structure.  See PCSetOperators() for details.

1075:    Level: intermediate

1077:    Notes: Does not increase the reference count of the matrices, so you should not destroy them

1079:    Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1080:       are created in PC and returned to the user. In this case, if both operators
1081:       mat and pmat are requested, two DIFFERENT operators will be returned. If
1082:       only one is requested both operators in the PC will be the same (i.e. as
1083:       if one had called KSP/PCSetOperators() with the same argument for both Mats).
1084:       The user must set the sizes of the returned matrices and their type etc just
1085:       as if the user created them with MatCreate(). For example,

1087: $         KSP/PCGetOperators(ksp/pc,&mat,PETSC_NULL,PETSC_NULL); is equivalent to
1088: $           set size, type, etc of mat

1090: $         MatCreate(comm,&mat);
1091: $         KSP/PCSetOperators(ksp/pc,mat,mat,SAME_NONZERO_PATTERN);
1092: $         PetscObjectDereference((PetscObject)mat);
1093: $           set size, type, etc of mat

1095:      and

1097: $         KSP/PCGetOperators(ksp/pc,&mat,&pmat,PETSC_NULL); is equivalent to
1098: $           set size, type, etc of mat and pmat

1100: $         MatCreate(comm,&mat);
1101: $         MatCreate(comm,&pmat);
1102: $         KSP/PCSetOperators(ksp/pc,mat,pmat,SAME_NONZERO_PATTERN);
1103: $         PetscObjectDereference((PetscObject)mat);
1104: $         PetscObjectDereference((PetscObject)pmat);
1105: $           set size, type, etc of mat and pmat

1107:     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1108:     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely 
1109:     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1110:     at this is when you create a SNES you do not NEED to create a KSP and attach it to 
1111:     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1112:     you do not need to attach a PC to it (the KSP object manages the PC object for you).
1113:     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1114:     it can be created for you?
1115:      

1117: .keywords: PC, get, operators, matrix, linear system

1119: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1120: @*/
1121: PetscErrorCode  PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1122: {

1127:   if (mat) {
1128:     if (!pc->mat) {
1129:       if (pc->pmat && !pmat) {  /* pmat has been set, but user did not request it, so use for mat */
1130:         pc->mat = pc->pmat;
1131:         PetscObjectReference((PetscObject)pc->mat);
1132:       } else {                  /* both mat and pmat are empty */
1133:         MatCreate(((PetscObject)pc)->comm,&pc->mat);
1134:         if (!pmat) { /* user did NOT request pmat, so make same as mat */
1135:           pc->pmat = pc->mat;
1136:           PetscObjectReference((PetscObject)pc->pmat);
1137:         }
1138:       }
1139:     }
1140:     *mat  = pc->mat;
1141:   }
1142:   if (pmat) {
1143:     if (!pc->pmat) {
1144:       if (pc->mat && !mat) {    /* mat has been set but was not requested, so use for pmat */
1145:         pc->pmat = pc->mat;
1146:         PetscObjectReference((PetscObject)pc->pmat);
1147:       } else {
1148:         MatCreate(((PetscObject)pc)->comm,&pc->pmat);
1149:         if (!mat) { /* user did NOT request mat, so make same as pmat */
1150:           pc->mat = pc->pmat;
1151:           PetscObjectReference((PetscObject)pc->mat);
1152:         }
1153:       }
1154:     }
1155:     *pmat = pc->pmat;
1156:   }
1157:   if (flag) *flag = pc->flag;
1158:   return(0);
1159: }

1163: /*@C
1164:    PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1165:    possibly a different one associated with the preconditioner have been set in the PC.

1167:    Not collective, though the results on all processes should be the same

1169:    Input Parameter:
1170: .  pc - the preconditioner context

1172:    Output Parameters:
1173: +  mat - the matrix associated with the linear system was set
1174: -  pmat - matrix associated with the preconditioner was set, usually the same

1176:    Level: intermediate

1178: .keywords: PC, get, operators, matrix, linear system

1180: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1181: @*/
1182: PetscErrorCode  PCGetOperatorsSet(PC pc,PetscBool  *mat,PetscBool  *pmat)
1183: {
1186:   if (mat)  *mat  = (pc->mat)  ? PETSC_TRUE : PETSC_FALSE;
1187:   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1188:   return(0);
1189: }

1193: /*@
1194:    PCFactorGetMatrix - Gets the factored matrix from the
1195:    preconditioner context.  This routine is valid only for the LU, 
1196:    incomplete LU, Cholesky, and incomplete Cholesky methods.

1198:    Not Collective on PC though Mat is parallel if PC is parallel

1200:    Input Parameters:
1201: .  pc - the preconditioner context

1203:    Output parameters:
1204: .  mat - the factored matrix

1206:    Level: advanced

1208:    Notes: Does not increase the reference count for the matrix so DO NOT destroy it

1210: .keywords: PC, get, factored, matrix
1211: @*/
1212: PetscErrorCode  PCFactorGetMatrix(PC pc,Mat *mat)
1213: {

1219:   if (pc->ops->getfactoredmatrix) {
1220:     (*pc->ops->getfactoredmatrix)(pc,mat);
1221:   } else {
1222:     SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC type does not support getting factor matrix");
1223:   }
1224:   return(0);
1225: }

1229: /*@C
1230:    PCSetOptionsPrefix - Sets the prefix used for searching for all 
1231:    PC options in the database.

1233:    Logically Collective on PC

1235:    Input Parameters:
1236: +  pc - the preconditioner context
1237: -  prefix - the prefix string to prepend to all PC option requests

1239:    Notes:
1240:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1241:    The first character of all runtime options is AUTOMATICALLY the
1242:    hyphen.

1244:    Level: advanced

1246: .keywords: PC, set, options, prefix, database

1248: .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1249: @*/
1250: PetscErrorCode  PCSetOptionsPrefix(PC pc,const char prefix[])
1251: {

1256:   PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);
1257:   return(0);
1258: }

1262: /*@C
1263:    PCAppendOptionsPrefix - Appends to the prefix used for searching for all 
1264:    PC options in the database.

1266:    Logically Collective on PC

1268:    Input Parameters:
1269: +  pc - the preconditioner context
1270: -  prefix - the prefix string to prepend to all PC option requests

1272:    Notes:
1273:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1274:    The first character of all runtime options is AUTOMATICALLY the
1275:    hyphen.

1277:    Level: advanced

1279: .keywords: PC, append, options, prefix, database

1281: .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1282: @*/
1283: PetscErrorCode  PCAppendOptionsPrefix(PC pc,const char prefix[])
1284: {

1289:   PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);
1290:   return(0);
1291: }

1295: /*@C
1296:    PCGetOptionsPrefix - Gets the prefix used for searching for all 
1297:    PC options in the database.

1299:    Not Collective

1301:    Input Parameters:
1302: .  pc - the preconditioner context

1304:    Output Parameters:
1305: .  prefix - pointer to the prefix string used, is returned

1307:    Notes: On the fortran side, the user should pass in a string 'prifix' of
1308:    sufficient length to hold the prefix.

1310:    Level: advanced

1312: .keywords: PC, get, options, prefix, database

1314: .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1315: @*/
1316: PetscErrorCode  PCGetOptionsPrefix(PC pc,const char *prefix[])
1317: {

1323:   PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);
1324:   return(0);
1325: }

1329: /*@
1330:    PCPreSolve - Optional pre-solve phase, intended for any
1331:    preconditioner-specific actions that must be performed before 
1332:    the iterative solve itself.

1334:    Collective on PC

1336:    Input Parameters:
1337: +  pc - the preconditioner context
1338: -  ksp - the Krylov subspace context

1340:    Level: developer

1342:    Sample of Usage:
1343: .vb
1344:     PCPreSolve(pc,ksp);
1345:     KSPSolve(ksp,b,x);
1346:     PCPostSolve(pc,ksp);
1347: .ve

1349:    Notes:
1350:    The pre-solve phase is distinct from the PCSetUp() phase.

1352:    KSPSolve() calls this directly, so is rarely called by the user.

1354: .keywords: PC, pre-solve

1356: .seealso: PCPostSolve()
1357: @*/
1358: PetscErrorCode  PCPreSolve(PC pc,KSP ksp)
1359: {
1361:   Vec            x,rhs;
1362:   Mat            A,B;

1367:   KSPGetSolution(ksp,&x);
1368:   KSPGetRhs(ksp,&rhs);
1369:   /*
1370:       Scale the system and have the matrices use the scaled form
1371:     only if the two matrices are actually the same (and hence
1372:     have the same scaling
1373:   */
1374:   PCGetOperators(pc,&A,&B,PETSC_NULL);
1375:   if (A == B) {
1376:     MatScaleSystem(pc->mat,rhs,x);
1377:     MatUseScaledForm(pc->mat,PETSC_TRUE);
1378:   }

1380:   if (pc->ops->presolve) {
1381:     (*pc->ops->presolve)(pc,ksp,rhs,x);
1382:   }
1383:   return(0);
1384: }

1388: /*@
1389:    PCPostSolve - Optional post-solve phase, intended for any
1390:    preconditioner-specific actions that must be performed after
1391:    the iterative solve itself.

1393:    Collective on PC

1395:    Input Parameters:
1396: +  pc - the preconditioner context
1397: -  ksp - the Krylov subspace context

1399:    Sample of Usage:
1400: .vb
1401:     PCPreSolve(pc,ksp);
1402:     KSPSolve(ksp,b,x);
1403:     PCPostSolve(pc,ksp);
1404: .ve

1406:    Note:
1407:    KSPSolve() calls this routine directly, so it is rarely called by the user.

1409:    Level: developer

1411: .keywords: PC, post-solve

1413: .seealso: PCPreSolve(), KSPSolve()
1414: @*/
1415: PetscErrorCode  PCPostSolve(PC pc,KSP ksp)
1416: {
1418:   Vec            x,rhs;
1419:   Mat            A,B;

1424:   KSPGetSolution(ksp,&x);
1425:   KSPGetRhs(ksp,&rhs);
1426:   if (pc->ops->postsolve) {
1427:      (*pc->ops->postsolve)(pc,ksp,rhs,x);
1428:   }
1429:   /*
1430:       Scale the system and have the matrices use the scaled form
1431:     only if the two matrices are actually the same (and hence
1432:     have the same scaling
1433:   */
1434:   PCGetOperators(pc,&A,&B,PETSC_NULL);
1435:   if (A == B) {
1436:     MatUnScaleSystem(pc->mat,rhs,x);
1437:     MatUseScaledForm(pc->mat,PETSC_FALSE);
1438:   }
1439:   return(0);
1440: }

1444: /*@C
1445:    PCView - Prints the PC data structure.

1447:    Collective on PC

1449:    Input Parameters:
1450: +  PC - the PC context
1451: -  viewer - optional visualization context

1453:    Note:
1454:    The available visualization contexts include
1455: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1456: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1457:          output where only the first processor opens
1458:          the file.  All other processors send their 
1459:          data to the first processor to print. 

1461:    The user can open an alternative visualization contexts with
1462:    PetscViewerASCIIOpen() (output to a specified file).

1464:    Level: developer

1466: .keywords: PC, view

1468: .seealso: KSPView(), PetscViewerASCIIOpen()
1469: @*/
1470: PetscErrorCode  PCView(PC pc,PetscViewer viewer)
1471: {
1472:   const PCType      cstr;
1473:   PetscErrorCode    ierr;
1474:   PetscBool         iascii,isstring;
1475:   PetscViewerFormat format;

1479:   if (!viewer) {
1480:     PetscViewerASCIIGetStdout(((PetscObject)pc)->comm,&viewer);
1481:   }

1485:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1486:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
1487:   if (iascii) {
1488:     PetscViewerGetFormat(viewer,&format);
1489:     PetscObjectPrintClassNamePrefixType((PetscObject)pc,viewer,"PC Object");
1490:     if (pc->ops->view) {
1491:       PetscViewerASCIIPushTab(viewer);
1492:       (*pc->ops->view)(pc,viewer);
1493:       PetscViewerASCIIPopTab(viewer);
1494:     }
1495:     if (pc->mat) {
1496:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1497:       if (pc->pmat == pc->mat) {
1498:         PetscViewerASCIIPrintf(viewer,"  linear system matrix = precond matrix:\n");
1499:         PetscViewerASCIIPushTab(viewer);
1500:         MatView(pc->mat,viewer);
1501:         PetscViewerASCIIPopTab(viewer);
1502:       } else {
1503:         if (pc->pmat) {
1504:           PetscViewerASCIIPrintf(viewer,"  linear system matrix followed by preconditioner matrix:\n");
1505:         } else {
1506:           PetscViewerASCIIPrintf(viewer,"  linear system matrix:\n");
1507:         }
1508:         PetscViewerASCIIPushTab(viewer);
1509:         MatView(pc->mat,viewer);
1510:         if (pc->pmat) {MatView(pc->pmat,viewer);}
1511:         PetscViewerASCIIPopTab(viewer);
1512:       }
1513:       PetscViewerPopFormat(viewer);
1514:     }
1515:   } else if (isstring) {
1516:     PCGetType(pc,&cstr);
1517:     PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
1518:     if (pc->ops->view) {(*pc->ops->view)(pc,viewer);}
1519:   } else {
1520:     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1521:   }
1522:   return(0);
1523: }


1528: /*@
1529:    PCSetInitialGuessNonzero - Tells the iterative solver that the 
1530:    initial guess is nonzero; otherwise PC assumes the initial guess
1531:    is to be zero (and thus zeros it out before solving).

1533:    Logically Collective on PC

1535:    Input Parameters:
1536: +  pc - iterative context obtained from PCCreate()
1537: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

1539:    Level: Developer

1541:    Notes:
1542:     This is a weird function. Since PC's are linear operators on the right hand side they
1543:     CANNOT use an initial guess. This function is for the "pass-through" preconditioners
1544:     PCKSP, PCREDUNDANT and PCHMPI and causes the inner KSP object to use the nonzero
1545:     initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP.


1548: .keywords: PC, set, initial guess, nonzero

1550: .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
1551: @*/
1552: PetscErrorCode  PCSetInitialGuessNonzero(PC pc,PetscBool  flg)
1553: {
1556:   pc->nonzero_guess   = flg;
1557:   return(0);
1558: }

1562: /*@C
1563:   PCRegister - See PCRegisterDynamic()

1565:   Level: advanced
1566: @*/
1567: PetscErrorCode  PCRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PC))
1568: {
1570:   char           fullname[PETSC_MAX_PATH_LEN];


1574:   PetscFListConcat(path,name,fullname);
1575:   PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);
1576:   return(0);
1577: }

1581: /*@
1582:     PCComputeExplicitOperator - Computes the explicit preconditioned operator.  

1584:     Collective on PC

1586:     Input Parameter:
1587: .   pc - the preconditioner object

1589:     Output Parameter:
1590: .   mat - the explict preconditioned operator

1592:     Notes:
1593:     This computation is done by applying the operators to columns of the 
1594:     identity matrix.

1596:     Currently, this routine uses a dense matrix format when 1 processor
1597:     is used and a sparse format otherwise.  This routine is costly in general,
1598:     and is recommended for use only with relatively small systems.

1600:     Level: advanced
1601:    
1602: .keywords: PC, compute, explicit, operator

1604: .seealso: KSPComputeExplicitOperator()

1606: @*/
1607: PetscErrorCode  PCComputeExplicitOperator(PC pc,Mat *mat)
1608: {
1609:   Vec            in,out;
1611:   PetscInt       i,M,m,*rows,start,end;
1612:   PetscMPIInt    size;
1613:   MPI_Comm       comm;
1614:   PetscScalar    *array,one = 1.0;
1615: 

1620:   comm = ((PetscObject)pc)->comm;
1621:   MPI_Comm_size(comm,&size);

1623:   if (!pc->pmat) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1624:   MatGetVecs(pc->pmat,&in,0);
1625:   VecDuplicate(in,&out);
1626:   VecGetOwnershipRange(in,&start,&end);
1627:   VecGetSize(in,&M);
1628:   VecGetLocalSize(in,&m);
1629:   PetscMalloc((m+1)*sizeof(PetscInt),&rows);
1630:   for (i=0; i<m; i++) {rows[i] = start + i;}

1632:   MatCreate(comm,mat);
1633:   MatSetSizes(*mat,m,m,M,M);
1634:   if (size == 1) {
1635:     MatSetType(*mat,MATSEQDENSE);
1636:     MatSeqDenseSetPreallocation(*mat,PETSC_NULL);
1637:   } else {
1638:     MatSetType(*mat,MATMPIAIJ);
1639:     MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);
1640:   }

1642:   for (i=0; i<M; i++) {

1644:     VecSet(in,0.0);
1645:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
1646:     VecAssemblyBegin(in);
1647:     VecAssemblyEnd(in);

1649:     /* should fix, allowing user to choose side */
1650:     PCApply(pc,in,out);
1651: 
1652:     VecGetArray(out,&array);
1653:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
1654:     VecRestoreArray(out,&array);

1656:   }
1657:   PetscFree(rows);
1658:   VecDestroy(&out);
1659:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
1660:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
1661:   return(0);
1662: }

1666: /*@
1667:    PCSetCoordinates - sets the coordinates of all the nodes on the local process

1669:    Collective on PC

1671:    Input Parameters:
1672: +  pc - the solver context
1673: .  dim - the dimension of the coordinates 1, 2, or 3
1674: -  coords - the coordinates

1676:    Level: intermediate

1678:    Notes: coords is an array of the 3D coordinates for the nodes on
1679:    the local processor.  So if there are 108 equation on a processor
1680:    for a displacement finite element discretization of elasticity (so
1681:    that there are 36 = 108/3 nodes) then the array must have 108
1682:    double precision values (ie, 3 * 36).  These x y z coordinates
1683:    should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
1684:    ... , N-1.z ].

1686: .seealso: PCPROMETHEUS
1687: @*/
1688: PetscErrorCode  PCSetCoordinates(PC pc,PetscInt dim,PetscReal *coords)
1689: {

1693:   PetscTryMethod(pc,"PCSetCoordinates_C",(PC,PetscInt,PetscReal*),(pc,dim,coords));
1694:   return(0);
1695: }