Actual source code: precon.c

  1: #define PETSCKSP_DLL

  3: /*
  4:     The PC (preconditioner) interface routines, callable by users.
  5: */
 6:  #include private/pcimpl.h

  8: /* Logging support */
  9: PetscCookie PC_COOKIE = 0;
 10: PetscEvent  PC_SetUp = 0, PC_SetUpOnBlocks = 0, PC_Apply = 0, PC_ApplyCoarse = 0, PC_ApplyMultiple = 0, PC_ApplySymmetricLeft = 0;
 11: PetscEvent  PC_ApplySymmetricRight = 0, PC_ModifySubMatrices = 0;

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

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

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

 61:    Collective on PC

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

 66:    Level: developer

 68: .keywords: PC, destroy

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

 78:   if (--pc->refct > 0) return(0);

 80:   /* if memory was published with AMS then destroy it */
 81:   PetscObjectDepublish(pc);

 83:   if (pc->ops->destroy)       { (*pc->ops->destroy)(pc);}
 84:   if (pc->diagonalscaleright) {VecDestroy(pc->diagonalscaleright);}
 85:   if (pc->diagonalscaleleft)  {VecDestroy(pc->diagonalscaleleft);}

 87:   if (pc->pmat) {MatDestroy(pc->pmat);}
 88:   if (pc->mat) {MatDestroy(pc->mat);}

 90:   PetscHeaderDestroy(pc);
 91:   return(0);
 92: }

 96: /*@C
 97:    PCDiagonalScale - Indicates if the preconditioner applies an additional left and right
 98:       scaling as needed by certain time-stepping codes.

100:    Collective on PC

102:    Input Parameter:
103: .  pc - the preconditioner context

105:    Output Parameter:
106: .  flag - PETSC_TRUE if it applies the scaling

108:    Level: developer

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

114: .keywords: PC

116: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScaleSet()
117: @*/
118: PetscErrorCode  PCDiagonalScale(PC pc,PetscTruth *flag)
119: {
123:   *flag = pc->diagonalscale;
124:   return(0);
125: }

129: /*@
130:    PCDiagonalScaleSet - Indicates the left scaling to use to apply an additional left and right
131:       scaling as needed by certain time-stepping codes.

133:    Collective on PC

135:    Input Parameters:
136: +  pc - the preconditioner context
137: -  s - scaling vector

139:    Level: intermediate

141:    Notes: 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:    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.

147: .keywords: PC

149: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScale()
150: @*/
151: PetscErrorCode  PCDiagonalScaleSet(PC pc,Vec s)
152: {

158:   pc->diagonalscale     = PETSC_TRUE;
159:   PetscObjectReference((PetscObject)s);
160:   if (pc->diagonalscaleleft) {
161:     VecDestroy(pc->diagonalscaleleft);
162:   }
163:   pc->diagonalscaleleft = s;
164:   if (!pc->diagonalscaleright) {
165:     VecDuplicate(s,&pc->diagonalscaleright);
166:   }
167:   VecCopy(s,pc->diagonalscaleright);
168:   VecReciprocal(pc->diagonalscaleright);
169:   return(0);
170: }

174: /*@C
175:    PCDiagonalScaleLeft - Indicates the left scaling to use to apply an additional left and right
176:       scaling as needed by certain time-stepping codes.

178:    Collective on PC

180:    Input Parameters:
181: +  pc - the preconditioner context
182: .  in - input vector
183: +  out - scaled vector (maybe the same as in)

185:    Level: intermediate

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

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

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

195: .keywords: PC

197: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
198: @*/
199: PetscErrorCode  PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
200: {

207:   if (pc->diagonalscale) {
208:     VecPointwiseMult(out,pc->diagonalscaleleft,in);
209:   } else if (in != out) {
210:     VecCopy(in,out);
211:   }
212:   return(0);
213: }

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

220:    Collective on PC

222:    Input Parameters:
223: +  pc - the preconditioner context
224: .  in - input vector
225: +  out - scaled vector (maybe the same as in)

227:    Level: intermediate

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

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

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

237: .keywords: PC

239: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
240: @*/
241: PetscErrorCode  PCDiagonalScaleRight(PC pc,Vec in,Vec out)
242: {

249:   if (pc->diagonalscale) {
250:     VecPointwiseMult(out,pc->diagonalscaleright,in);
251:   } else if (in != out) {
252:     VecCopy(in,out);
253:   }
254:   return(0);
255: }

259: static PetscErrorCode PCPublish_Petsc(PetscObject obj)
260: {
262:   return(0);
263: }

267: /*@
268:    PCCreate - Creates a preconditioner context.

270:    Collective on MPI_Comm

272:    Input Parameter:
273: .  comm - MPI communicator 

275:    Output Parameter:
276: .  pc - location to put the preconditioner context

278:    Notes:
279:    The default preconditioner on one processor is PCILU with 0 fill on more 
280:    then one it is PCBJACOBI with ILU() on each processor.

282:    Level: developer

284: .keywords: PC, create, context

286: .seealso: PCSetUp(), PCApply(), PCDestroy()
287: @*/
288: PetscErrorCode  PCCreate(MPI_Comm comm,PC *newpc)
289: {
290:   PC             pc;

295:   *newpc = 0;
296: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
297:   PCInitializePackage(PETSC_NULL);
298: #endif

300:   PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_COOKIE,-1,"PC",comm,PCDestroy,PCView);
301:   pc->bops->publish        = PCPublish_Petsc;

303:   pc->mat                  = 0;
304:   pc->pmat                 = 0;
305:   pc->setupcalled          = 0;
306:   pc->setfromoptionscalled = 0;
307:   pc->data                 = 0;
308:   pc->diagonalscale        = PETSC_FALSE;
309:   pc->diagonalscaleleft    = 0;
310:   pc->diagonalscaleright   = 0;

312:   pc->modifysubmatrices   = 0;
313:   pc->modifysubmatricesP  = 0;
314:   PetscPublishAll(pc);
315:   *newpc = pc;
316:   return(0);

318: }

320: /* -------------------------------------------------------------------------------*/

324: /*@
325:    PCApply - Applies the preconditioner to a vector.

327:    Collective on PC and Vec

329:    Input Parameters:
330: +  pc - the preconditioner context
331: -  x - input vector

333:    Output Parameter:
334: .  y - output vector

336:    Level: developer

338: .keywords: PC, apply

340: .seealso: PCApplyTranspose(), PCApplyBAorAB()
341: @*/
342: PetscErrorCode  PCApply(PC pc,Vec x,Vec y)
343: {

350:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
351:   if (pc->setupcalled < 2) {
352:     PCSetUp(pc);
353:   }
354:   if (!pc->ops->apply) SETERRQ(PETSC_ERR_SUP,"PC does not have apply");
356:   (*pc->ops->apply)(pc,x,y);
358:   return(0);
359: }

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

366:    Collective on PC and Vec

368:    Input Parameters:
369: +  pc - the preconditioner context
370: -  x - input vector

372:    Output Parameter:
373: .  y - output vector

375:    Notes:
376:    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.

378:    Level: developer

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

382: .seealso: PCApply(), PCApplySymmetricRight()
383: @*/
384: PetscErrorCode  PCApplySymmetricLeft(PC pc,Vec x,Vec y)
385: {

392:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
393:   if (pc->setupcalled < 2) {
394:     PCSetUp(pc);
395:   }
396:   if (!pc->ops->applysymmetricleft) SETERRQ(PETSC_ERR_SUP,"PC does not have left symmetric apply");
398:   (*pc->ops->applysymmetricleft)(pc,x,y);
400:   return(0);
401: }

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

408:    Collective on PC and Vec

410:    Input Parameters:
411: +  pc - the preconditioner context
412: -  x - input vector

414:    Output Parameter:
415: .  y - output vector

417:    Level: developer

419:    Notes:
420:    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.

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

424: .seealso: PCApply(), PCApplySymmetricLeft()
425: @*/
426: PetscErrorCode  PCApplySymmetricRight(PC pc,Vec x,Vec y)
427: {

434:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
435:   if (pc->setupcalled < 2) {
436:     PCSetUp(pc);
437:   }
438:   if (!pc->ops->applysymmetricright) SETERRQ(PETSC_ERR_SUP,"PC does not have left symmetric apply");
440:   (*pc->ops->applysymmetricright)(pc,x,y);
442:   return(0);
443: }

447: /*@
448:    PCApplyTranspose - Applies the transpose of preconditioner to a vector.

450:    Collective on PC and Vec

452:    Input Parameters:
453: +  pc - the preconditioner context
454: -  x - input vector

456:    Output Parameter:
457: .  y - output vector

459:    Level: developer

461: .keywords: PC, apply, transpose

463: .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCHasApplyTranspose()
464: @*/
465: PetscErrorCode  PCApplyTranspose(PC pc,Vec x,Vec y)
466: {

473:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
474:   if (pc->setupcalled < 2) {
475:     PCSetUp(pc);
476:   }
477:   if (!pc->ops->applytranspose) SETERRQ(PETSC_ERR_SUP,"PC does not have apply transpose");
479:   (*pc->ops->applytranspose)(pc,x,y);
481:   return(0);
482: }

486: /*@
487:    PCHasApplyTranspose - Test whether the preconditioner has a transpose apply operation

489:    Collective on PC and Vec

491:    Input Parameters:
492: .  pc - the preconditioner context

494:    Output Parameter:
495: .  flg - PETSC_TRUE if a transpose operation is defined

497:    Level: developer

499: .keywords: PC, apply, transpose

501: .seealso: PCApplyTranspose()
502: @*/
503: PetscErrorCode  PCHasApplyTranspose(PC pc,PetscTruth *flg)
504: {
508:   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
509:   else                         *flg = PETSC_FALSE;
510:   return(0);
511: }

515: /*@
516:    PCApplyBAorAB - Applies the preconditioner and operator to a vector. 

518:    Collective on PC and Vec

520:    Input Parameters:
521: +  pc - the preconditioner context
522: .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
523: .  x - input vector
524: -  work - work vector

526:    Output Parameter:
527: .  y - output vector

529:    Level: developer

531: .keywords: PC, apply, operator

533: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
534: @*/
535: PetscErrorCode  PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
536: {

544:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
545:   if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) {
546:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric");
547:   }
548:   if (pc->diagonalscale && side == PC_SYMMETRIC) {
549:     SETERRQ(PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");
550:   }

552:   if (pc->setupcalled < 2) {
553:     PCSetUp(pc);
554:   }

556:   if (pc->diagonalscale) {
557:     if (pc->ops->applyBA) {
558:       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
559:       VecDuplicate(x,&work2);
560:       PCDiagonalScaleRight(pc,x,work2);
561:       (*pc->ops->applyBA)(pc,side,work2,y,work);
562:       PCDiagonalScaleLeft(pc,y,y);
563:       VecDestroy(work2);
564:     } else if (side == PC_RIGHT) {
565:       PCDiagonalScaleRight(pc,x,y);
566:       PCApply(pc,y,work);
567:       MatMult(pc->mat,work,y);
568:       PCDiagonalScaleLeft(pc,y,y);
569:     } else if (side == PC_LEFT) {
570:       PCDiagonalScaleRight(pc,x,y);
571:       MatMult(pc->mat,y,work);
572:       PCApply(pc,work,y);
573:       PCDiagonalScaleLeft(pc,y,y);
574:     } else if (side == PC_SYMMETRIC) {
575:       SETERRQ(PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
576:     }
577:   } else {
578:     if (pc->ops->applyBA) {
579:       (*pc->ops->applyBA)(pc,side,x,y,work);
580:     } else if (side == PC_RIGHT) {
581:       PCApply(pc,x,work);
582:       MatMult(pc->mat,work,y);
583:     } else if (side == PC_LEFT) {
584:       MatMult(pc->mat,x,work);
585:       PCApply(pc,work,y);
586:     } else if (side == PC_SYMMETRIC) {
587:       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
588:       PCApplySymmetricRight(pc,x,work);
589:       MatMult(pc->mat,work,y);
590:       VecCopy(y,work);
591:       PCApplySymmetricLeft(pc,work,y);
592:     }
593:   }
594:   return(0);
595: }

599: /*@ 
600:    PCApplyBAorABTranspose - Applies the transpose of the preconditioner
601:    and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
602:    not tr(B*A) = tr(A)*tr(B).

604:    Collective on PC and Vec

606:    Input Parameters:
607: +  pc - the preconditioner context
608: .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
609: .  x - input vector
610: -  work - work vector

612:    Output Parameter:
613: .  y - output vector

615:    Level: developer

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

619: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
620: @*/
621: PetscErrorCode  PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
622: {

630:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
631:   if (pc->ops->applyBAtranspose) {
632:     (*pc->ops->applyBAtranspose)(pc,side,x,y,work);
633:     return(0);
634:   }
635:   if (side != PC_LEFT && side != PC_RIGHT) {
636:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
637:   }

639:   if (pc->setupcalled < 2) {
640:     PCSetUp(pc);
641:   }

643:   if (side == PC_RIGHT) {
644:     PCApplyTranspose(pc,x,work);
645:     MatMultTranspose(pc->mat,work,y);
646:   } else if (side == PC_LEFT) {
647:     MatMultTranspose(pc->mat,x,work);
648:     PCApplyTranspose(pc,work,y);
649:   }
650:   /* add support for PC_SYMMETRIC */
651:   return(0); /* actually will never get here */
652: }

654: /* -------------------------------------------------------------------------------*/

658: /*@
659:    PCApplyRichardsonExists - Determines whether a particular preconditioner has a 
660:    built-in fast application of Richardson's method.

662:    Not Collective

664:    Input Parameter:
665: .  pc - the preconditioner

667:    Output Parameter:
668: .  exists - PETSC_TRUE or PETSC_FALSE

670:    Level: developer

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

674: .seealso: PCApplyRichardson()
675: @*/
676: PetscErrorCode  PCApplyRichardsonExists(PC pc,PetscTruth *exists)
677: {
681:   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
682:   else                          *exists = PETSC_FALSE;
683:   return(0);
684: }

688: /*@
689:    PCApplyRichardson - Applies several steps of Richardson iteration with 
690:    the particular preconditioner. This routine is usually used by the 
691:    Krylov solvers and not the application code directly.

693:    Collective on PC

695:    Input Parameters:
696: +  pc  - the preconditioner context
697: .  x   - the initial guess 
698: .  w   - one work vector
699: .  rtol - relative decrease in residual norm convergence criteria
700: .  abstol - absolute residual norm convergence criteria
701: .  dtol - divergence residual norm increase criteria
702: -  its - the number of iterations to apply.

704:    Output Parameter:
705: .  y - the solution

707:    Notes: 
708:    Most preconditioners do not support this function. Use the command
709:    PCApplyRichardsonExists() to determine if one does.

711:    Except for the multigrid PC this routine ignores the convergence tolerances
712:    and always runs for the number of iterations
713:  
714:    Level: developer

716: .keywords: PC, apply, Richardson

718: .seealso: PCApplyRichardsonExists()
719: @*/
720: PetscErrorCode  PCApplyRichardson(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
721: {

729:   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
730:   if (pc->setupcalled < 2) {
731:     PCSetUp(pc);
732:   }
733:   if (!pc->ops->applyrichardson) SETERRQ(PETSC_ERR_SUP,"PC does not have apply richardson");
734:   (*pc->ops->applyrichardson)(pc,x,y,w,rtol,abstol,dtol,its);
735:   return(0);
736: }

738: /* 
739:       a setupcall of 0 indicates never setup, 
740:                      1 needs to be resetup,
741:                      2 does not need any changes.
742: */
745: /*@
746:    PCSetUp - Prepares for the use of a preconditioner.

748:    Collective on PC

750:    Input Parameter:
751: .  pc - the preconditioner context

753:    Level: developer

755: .keywords: PC, setup

757: .seealso: PCCreate(), PCApply(), PCDestroy()
758: @*/
759: PetscErrorCode  PCSetUp(PC pc)
760: {
762:   const char     *def;


767:   if (!pc->mat) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}

769:   if (pc->setupcalled > 1) {
770:     PetscInfo(pc,"Setting PC with identical preconditioner\n");
771:     return(0);
772:   } else if (!pc->setupcalled) {
773:     PetscInfo(pc,"Setting up new PC\n");
774:   } else if (pc->flag == SAME_NONZERO_PATTERN) {
775:     PetscInfo(pc,"Setting up PC with same nonzero pattern\n");
776:   } else {
777:     PetscInfo(pc,"Setting up PC with different nonzero pattern\n");
778:   }

780:   if (!pc->type_name) {
781:     PCGetDefaultType_Private(pc,&def);
782:     PCSetType(pc,def);
783:   }

786:   if (pc->ops->setup) {
787:     (*pc->ops->setup)(pc);
788:   }
789:   pc->setupcalled = 2;
791:   return(0);
792: }

796: /*@
797:    PCSetUpOnBlocks - Sets up the preconditioner for each block in
798:    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz 
799:    methods.

801:    Collective on PC

803:    Input Parameters:
804: .  pc - the preconditioner context

806:    Level: developer

808: .keywords: PC, setup, blocks

810: .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
811: @*/
812: PetscErrorCode  PCSetUpOnBlocks(PC pc)
813: {

818:   if (!pc->ops->setuponblocks) return(0);
820:   (*pc->ops->setuponblocks)(pc);
822:   return(0);
823: }

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

834:    Collective on PC

836:    Input Parameters:
837: +  pc - the preconditioner context
838: .  func - routine for modifying the submatrices
839: -  ctx - optional user-defined context (may be null)

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

844: .  row - an array of index sets that contain the global row numbers
845:          that comprise each local submatrix
846: .  col - an array of index sets that contain the global column numbers
847:          that comprise each local submatrix
848: .  submat - array of local submatrices
849: -  ctx - optional user-defined context for private data for the 
850:          user-defined func routine (may be null)

852:    Notes:
853:    PCSetModifySubMatrices() MUST be called before KSPSetUp() and
854:    KSPSolve().

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

860:    Level: advanced

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

864: .seealso: PCModifySubMatrices()
865: @*/
866: PetscErrorCode  PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
867: {
870:   pc->modifysubmatrices  = func;
871:   pc->modifysubmatricesP = ctx;
872:   return(0);
873: }

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

881:    Collective on PC

883:    Input Parameters:
884: +  pc - the preconditioner context
885: .  nsub - the number of local submatrices
886: .  row - an array of index sets that contain the global row numbers
887:          that comprise each local submatrix
888: .  col - an array of index sets that contain the global column numbers
889:          that comprise each local submatrix
890: .  submat - array of local submatrices
891: -  ctx - optional user-defined context for private data for the 
892:          user-defined routine (may be null)

894:    Output Parameter:
895: .  submat - array of local submatrices (the entries of which may
896:             have been modified)

898:    Notes:
899:    The user should NOT generally call this routine, as it will
900:    automatically be called within certain preconditioners (currently
901:    block Jacobi, additive Schwarz) if set.

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

908:    Level: developer

910: .keywords: PC, modify, submatrices

912: .seealso: PCSetModifySubMatrices()
913: @*/
914: PetscErrorCode  PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
915: {

920:   if (!pc->modifysubmatrices) return(0);
922:   (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);
924:   return(0);
925: }

929: /*@
930:    PCSetOperators - Sets the matrix associated with the linear system and 
931:    a (possibly) different one associated with the preconditioner.

933:    Collective on PC and Mat

935:    Input Parameters:
936: +  pc - the preconditioner context
937: .  Amat - the matrix associated with the linear system
938: .  Pmat - the matrix to be used in constructing the preconditioner, usually
939:           the same as Amat. 
940: -  flag - flag indicating information about the preconditioner matrix structure
941:    during successive linear solves.  This flag is ignored the first time a
942:    linear system is solved, and thus is irrelevant when solving just one linear
943:    system.

945:    Notes: 
946:    The flag can be used to eliminate unnecessary work in the preconditioner 
947:    during the repeated solution of linear systems of the same size.  The 
948:    available options are
949: +    SAME_PRECONDITIONER -
950:        Pmat is identical during successive linear solves.
951:        This option is intended for folks who are using
952:        different Amat and Pmat matrices and wish to reuse the
953:        same preconditioner matrix.  For example, this option
954:        saves work by not recomputing incomplete factorization
955:        for ILU/ICC preconditioners.
956: .     SAME_NONZERO_PATTERN -
957:        Pmat has the same nonzero structure during
958:        successive linear solves. 
959: -     DIFFERENT_NONZERO_PATTERN -
960:        Pmat does not have the same nonzero structure.

962:    Caution:
963:    If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
964:    and does not check the structure of the matrix.  If you erroneously
965:    claim that the structure is the same when it actually is not, the new
966:    preconditioner will not function correctly.  Thus, use this optimization
967:    feature carefully!

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

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

978:    Level: intermediate

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

982: .seealso: PCGetOperators(), MatZeroEntries()
983:  @*/
984: PetscErrorCode  PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
985: {
987:   PetscTruth     isbjacobi,isrowbs;


996:   /*
997:       BlockSolve95 cannot use default BJacobi preconditioning
998:   */
999:   if (Amat) {
1000:     PetscTypeCompare((PetscObject)Amat,MATMPIROWBS,&isrowbs);
1001:     if (isrowbs) {
1002:       PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);
1003:       if (isbjacobi) {
1004:         PCSetType(pc,PCILU);
1005:         PetscInfo(pc,"Switching default PC to PCILU since BS95 doesn't support PCBJACOBI\n");
1006:       }
1007:     }
1008:   }

1010:   /* reference first in case the matrices are the same */
1011:   if (Amat) {PetscObjectReference((PetscObject)Amat);}
1012:   if (pc->mat) {MatDestroy(pc->mat);}
1013:   if (Pmat) {PetscObjectReference((PetscObject)Pmat);}
1014:   if (pc->pmat) {MatDestroy(pc->pmat);}
1015:   pc->mat  = Amat;
1016:   pc->pmat = Pmat;

1018:   if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1019:     pc->setupcalled = 1;
1020:   }
1021:   pc->flag = flag;
1022:   return(0);
1023: }

1027: /*@C
1028:    PCGetOperators - Gets the matrix associated with the linear system and
1029:    possibly a different one associated with the preconditioner.

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

1033:    Input Parameter:
1034: .  pc - the preconditioner context

1036:    Output Parameters:
1037: +  mat - the matrix associated with the linear system
1038: .  pmat - matrix associated with the preconditioner, usually the same
1039:           as mat. 
1040: -  flag - flag indicating information about the preconditioner
1041:           matrix structure.  See PCSetOperators() for details.

1043:    Level: intermediate

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

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

1056: $         MatCreate(comm,&mat);
1057: $         KSP/PCSetOperators(ksp/pc,mat,mat,SAME_NONZERO_PATTERN);
1058: $         PetscObjectDereference((PetscObject)mat);
1059: $           set size, type, etc of mat

1061:      and

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

1066: $         MatCreate(comm,&mat);
1067: $         MatCreate(comm,&pmat);
1068: $         KSP/PCSetOperators(ksp/pc,mat,pmat,SAME_NONZERO_PATTERN);
1069: $         PetscObjectDereference((PetscObject)mat);
1070: $         PetscObjectDereference((PetscObject)pmat);
1071: $           set size, type, etc of mat and pmat

1073:     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1074:     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely 
1075:     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1076:     at this is when you create a SNES you do not NEED to create a KSP and attach it to 
1077:     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1078:     you do not need to attach a PC to it (the KSP object manages the PC object for you).
1079:     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1080:     it can be created for you?
1081:      

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

1085: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1086: @*/
1087: PetscErrorCode  PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1088: {

1093:   if (mat) {
1094:     if (!pc->mat) {
1095:       MatCreate(pc->comm,&pc->mat);
1096:       if (!pc->pmat && !pmat) { /* user did NOT request pmat, so make same as mat */
1097:         pc->pmat = pc->mat;
1098:         PetscObjectReference((PetscObject)pc->pmat);
1099:       }
1100:     }
1101:     *mat  = pc->mat;
1102:   }
1103:   if (pmat) {
1104:     if (!pc->pmat) {
1105:       MatCreate(pc->comm,&pc->mat);
1106:       if (!pc->mat && !mat) { /* user did NOT request mat, so make same as pmat */
1107:         pc->mat = pc->pmat;
1108:         PetscObjectReference((PetscObject)pc->mat);
1109:       }
1110:     }
1111:     *pmat = pc->pmat;
1112:   }
1113:   if (flag) *flag = pc->flag;
1114:   return(0);
1115: }

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

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

1125:    Input Parameter:
1126: .  pc - the preconditioner context

1128:    Output Parameters:
1129: +  mat - the matrix associated with the linear system was set
1130: -  pmat - matrix associated with the preconditioner was set, usually the same

1132:    Level: intermediate

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

1136: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1137: @*/
1138: PetscErrorCode  PCGetOperatorsSet(PC pc,PetscTruth *mat,PetscTruth *pmat)
1139: {
1142:   if (mat)  *mat  = (pc->mat)  ? PETSC_TRUE : PETSC_FALSE;
1143:   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1144:   return(0);
1145: }

1149: /*@
1150:    PCGetFactoredMatrix - Gets the factored matrix from the
1151:    preconditioner context.  This routine is valid only for the LU, 
1152:    incomplete LU, Cholesky, and incomplete Cholesky methods.

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

1156:    Input Parameters:
1157: .  pc - the preconditioner context

1159:    Output parameters:
1160: .  mat - the factored matrix

1162:    Level: advanced

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

1166: .keywords: PC, get, factored, matrix
1167: @*/
1168: PetscErrorCode  PCGetFactoredMatrix(PC pc,Mat *mat)
1169: {

1175:   if (pc->ops->getfactoredmatrix) {
1176:     (*pc->ops->getfactoredmatrix)(pc,mat);
1177:   }
1178:   return(0);
1179: }

1183: /*@C
1184:    PCSetOptionsPrefix - Sets the prefix used for searching for all 
1185:    PC options in the database.

1187:    Collective on PC

1189:    Input Parameters:
1190: +  pc - the preconditioner context
1191: -  prefix - the prefix string to prepend to all PC option requests

1193:    Notes:
1194:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1195:    The first character of all runtime options is AUTOMATICALLY the
1196:    hyphen.

1198:    Level: advanced

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

1202: .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1203: @*/
1204: PetscErrorCode  PCSetOptionsPrefix(PC pc,const char prefix[])
1205: {

1210:   PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);
1211:   return(0);
1212: }

1216: /*@C
1217:    PCAppendOptionsPrefix - Appends to the prefix used for searching for all 
1218:    PC options in the database.

1220:    Collective on PC

1222:    Input Parameters:
1223: +  pc - the preconditioner context
1224: -  prefix - the prefix string to prepend to all PC option requests

1226:    Notes:
1227:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1228:    The first character of all runtime options is AUTOMATICALLY the
1229:    hyphen.

1231:    Level: advanced

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

1235: .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1236: @*/
1237: PetscErrorCode  PCAppendOptionsPrefix(PC pc,const char prefix[])
1238: {

1243:   PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);
1244:   return(0);
1245: }

1249: /*@C
1250:    PCGetOptionsPrefix - Gets the prefix used for searching for all 
1251:    PC options in the database.

1253:    Not Collective

1255:    Input Parameters:
1256: .  pc - the preconditioner context

1258:    Output Parameters:
1259: .  prefix - pointer to the prefix string used, is returned

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

1264:    Level: advanced

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

1268: .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1269: @*/
1270: PetscErrorCode  PCGetOptionsPrefix(PC pc,const char *prefix[])
1271: {

1277:   PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);
1278:   return(0);
1279: }

1283: /*@
1284:    PCPreSolve - Optional pre-solve phase, intended for any
1285:    preconditioner-specific actions that must be performed before 
1286:    the iterative solve itself.

1288:    Collective on PC

1290:    Input Parameters:
1291: +  pc - the preconditioner context
1292: -  ksp - the Krylov subspace context

1294:    Level: developer

1296:    Sample of Usage:
1297: .vb
1298:     PCPreSolve(pc,ksp);
1299:     KSPSolve(ksp,b,x);
1300:     PCPostSolve(pc,ksp);
1301: .ve

1303:    Notes:
1304:    The pre-solve phase is distinct from the PCSetUp() phase.

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

1308: .keywords: PC, pre-solve

1310: .seealso: PCPostSolve()
1311: @*/
1312: PetscErrorCode  PCPreSolve(PC pc,KSP ksp)
1313: {
1315:   Vec            x,rhs;
1316:   Mat            A,B;

1321:   KSPGetSolution(ksp,&x);
1322:   KSPGetRhs(ksp,&rhs);
1323:   /*
1324:       Scale the system and have the matrices use the scaled form
1325:     only if the two matrices are actually the same (and hence
1326:     have the same scaling
1327:   */
1328:   PCGetOperators(pc,&A,&B,PETSC_NULL);
1329:   if (A == B) {
1330:     MatScaleSystem(pc->mat,rhs,x);
1331:     MatUseScaledForm(pc->mat,PETSC_TRUE);
1332:   }

1334:   if (pc->ops->presolve) {
1335:     (*pc->ops->presolve)(pc,ksp,rhs,x);
1336:   }
1337:   return(0);
1338: }

1342: /*@
1343:    PCPostSolve - Optional post-solve phase, intended for any
1344:    preconditioner-specific actions that must be performed after
1345:    the iterative solve itself.

1347:    Collective on PC

1349:    Input Parameters:
1350: +  pc - the preconditioner context
1351: -  ksp - the Krylov subspace context

1353:    Sample of Usage:
1354: .vb
1355:     PCPreSolve(pc,ksp);
1356:     KSPSolve(ksp,b,x);
1357:     PCPostSolve(pc,ksp);
1358: .ve

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

1363:    Level: developer

1365: .keywords: PC, post-solve

1367: .seealso: PCPreSolve(), KSPSolve()
1368: @*/
1369: PetscErrorCode  PCPostSolve(PC pc,KSP ksp)
1370: {
1372:   Vec            x,rhs;
1373:   Mat            A,B;

1378:   KSPGetSolution(ksp,&x);
1379:   KSPGetRhs(ksp,&rhs);
1380:   if (pc->ops->postsolve) {
1381:      (*pc->ops->postsolve)(pc,ksp,rhs,x);
1382:   }
1383:   /*
1384:       Scale the system and have the matrices use the scaled form
1385:     only if the two matrices are actually the same (and hence
1386:     have the same scaling
1387:   */
1388:   PCGetOperators(pc,&A,&B,PETSC_NULL);
1389:   if (A == B) {
1390:     MatUnScaleSystem(pc->mat,rhs,x);
1391:     MatUseScaledForm(pc->mat,PETSC_FALSE);
1392:   }
1393:   return(0);
1394: }

1398: /*@C
1399:    PCView - Prints the PC data structure.

1401:    Collective on PC

1403:    Input Parameters:
1404: +  PC - the PC context
1405: -  viewer - optional visualization context

1407:    Note:
1408:    The available visualization contexts include
1409: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1410: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1411:          output where only the first processor opens
1412:          the file.  All other processors send their 
1413:          data to the first processor to print. 

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

1418:    Level: developer

1420: .keywords: PC, view

1422: .seealso: KSPView(), PetscViewerASCIIOpen()
1423: @*/
1424: PetscErrorCode  PCView(PC pc,PetscViewer viewer)
1425: {
1426:   PCType            cstr;
1427:   PetscErrorCode    ierr;
1428:   PetscTruth        mat_exists,iascii,isstring;
1429:   PetscViewerFormat format;

1433:   if (!viewer) {
1434:     PetscViewerASCIIGetStdout(pc->comm,&viewer);
1435:   }

1439:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
1440:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
1441:   if (iascii) {
1442:     PetscViewerGetFormat(viewer,&format);
1443:     if (pc->prefix) {
1444:       PetscViewerASCIIPrintf(viewer,"PC Object:(%s)\n",pc->prefix);
1445:     } else {
1446:       PetscViewerASCIIPrintf(viewer,"PC Object:\n");
1447:     }
1448:     PCGetType(pc,&cstr);
1449:     if (cstr) {
1450:       PetscViewerASCIIPrintf(viewer,"  type: %s\n",cstr);
1451:     } else {
1452:       PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");
1453:     }
1454:     if (pc->ops->view) {
1455:       PetscViewerASCIIPushTab(viewer);
1456:       (*pc->ops->view)(pc,viewer);
1457:       PetscViewerASCIIPopTab(viewer);
1458:     }
1459:     PetscObjectExists((PetscObject)pc->mat,&mat_exists);
1460:     if (mat_exists) {
1461:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1462:       if (pc->pmat == pc->mat) {
1463:         PetscViewerASCIIPrintf(viewer,"  linear system matrix = precond matrix:\n");
1464:         PetscViewerASCIIPushTab(viewer);
1465:         MatView(pc->mat,viewer);
1466:         PetscViewerASCIIPopTab(viewer);
1467:       } else {
1468:         PetscObjectExists((PetscObject)pc->pmat,&mat_exists);
1469:         if (mat_exists) {
1470:           PetscViewerASCIIPrintf(viewer,"  linear system matrix followed by preconditioner matrix:\n");
1471:         } else {
1472:           PetscViewerASCIIPrintf(viewer,"  linear system matrix:\n");
1473:         }
1474:         PetscViewerASCIIPushTab(viewer);
1475:         MatView(pc->mat,viewer);
1476:         if (mat_exists) {MatView(pc->pmat,viewer);}
1477:         PetscViewerASCIIPopTab(viewer);
1478:       }
1479:       PetscViewerPopFormat(viewer);
1480:     }
1481:   } else if (isstring) {
1482:     PCGetType(pc,&cstr);
1483:     PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
1484:     if (pc->ops->view) {(*pc->ops->view)(pc,viewer);}
1485:   } else {
1486:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1487:   }
1488:   return(0);
1489: }

1493: /*@C
1494:   PCRegister - See PCRegisterDynamic()

1496:   Level: advanced
1497: @*/
1498: PetscErrorCode  PCRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PC))
1499: {
1501:   char           fullname[PETSC_MAX_PATH_LEN];


1505:   PetscFListConcat(path,name,fullname);
1506:   PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);
1507:   return(0);
1508: }

1512: /*@
1513:     PCComputeExplicitOperator - Computes the explicit preconditioned operator.  

1515:     Collective on PC

1517:     Input Parameter:
1518: .   pc - the preconditioner object

1520:     Output Parameter:
1521: .   mat - the explict preconditioned operator

1523:     Notes:
1524:     This computation is done by applying the operators to columns of the 
1525:     identity matrix.

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

1531:     Level: advanced
1532:    
1533: .keywords: PC, compute, explicit, operator

1535: @*/
1536: PetscErrorCode  PCComputeExplicitOperator(PC pc,Mat *mat)
1537: {
1538:   Vec            in,out;
1540:   PetscInt       i,M,m,*rows,start,end;
1541:   PetscMPIInt    size;
1542:   MPI_Comm       comm;
1543:   PetscScalar    *array,one = 1.0;
1544: 

1549:   comm = pc->comm;
1550:   MPI_Comm_size(comm,&size);

1552:   if (!pc->pmat) SETERRQ(PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1553:   MatGetVecs(pc->pmat,&in,0);
1554:   VecDuplicate(in,&out);
1555:   VecGetOwnershipRange(in,&start,&end);
1556:   VecGetSize(in,&M);
1557:   VecGetLocalSize(in,&m);
1558:   PetscMalloc((m+1)*sizeof(PetscInt),&rows);
1559:   for (i=0; i<m; i++) {rows[i] = start + i;}

1561:   MatCreate(comm,mat);
1562:   MatSetSizes(*mat,m,m,M,M);
1563:   if (size == 1) {
1564:     MatSetType(*mat,MATSEQDENSE);
1565:     MatSeqDenseSetPreallocation(*mat,PETSC_NULL);
1566:   } else {
1567:     MatSetType(*mat,MATMPIAIJ);
1568:     MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);
1569:   }

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

1573:     VecSet(in,0.0);
1574:     VecSetValues(in,1,&i,&one,INSERT_VALUES);
1575:     VecAssemblyBegin(in);
1576:     VecAssemblyEnd(in);

1578:     /* should fix, allowing user to choose side */
1579:     PCApply(pc,in,out);
1580: 
1581:     VecGetArray(out,&array);
1582:     MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
1583:     VecRestoreArray(out,&array);

1585:   }
1586:   PetscFree(rows);
1587:   VecDestroy(out);
1588:   MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
1589:   MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
1590:   return(0);
1591: }