xref: /petsc/src/ksp/pc/interface/precon.c (revision 0ce6c5a270a8e214277c9ee66b6a35b5056bb115)
1 #define PETSCKSP_DLL
2 
3 /*
4     The PC (preconditioner) interface routines, callable by users.
5 */
6 #include "private/pcimpl.h"            /*I "petscksp.h" I*/
7 
8 /* Logging support */
9 PetscClassId   PC_CLASSID;
10 PetscLogEvent  PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
11 PetscLogEvent  PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks;
12 
13 #undef __FUNCT__
14 #define __FUNCT__ "PCGetDefaultType_Private"
15 PetscErrorCode PCGetDefaultType_Private(PC pc,const char* type[])
16 {
17   PetscErrorCode ierr;
18   PetscMPIInt    size;
19   PetscBool      flg1,flg2,set,flg3;
20 
21   PetscFunctionBegin;
22   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
23   if (pc->pmat) {
24     PetscErrorCode (*f)(Mat,PetscBool *,MatReuse,Mat*);
25     ierr = PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
26     if (size == 1) {
27       ierr = MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ICC,&flg1);CHKERRQ(ierr);
28       ierr = MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ILU,&flg2);CHKERRQ(ierr);
29       ierr = MatIsSymmetricKnown(pc->pmat,&set,&flg3);CHKERRQ(ierr);
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   PetscFunctionReturn(0);
54 }
55 
56 #undef __FUNCT__
57 #define __FUNCT__ "PCReset"
58 /*@
59    PCReset - Resets a PC context to the pcsetupcalled = 0 state and removes any allocated Vecs and Mats
60 
61    Collective on PC
62 
63    Input Parameter:
64 .  pc - the preconditioner context
65 
66    Level: developer
67 
68    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 
70 .keywords: PC, destroy
71 
72 .seealso: PCCreate(), PCSetUp()
73 @*/
74 PetscErrorCode  PCReset(PC pc)
75 {
76   PetscErrorCode ierr;
77 
78   PetscFunctionBegin;
79   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
80   if (pc->dm) {ierr = DMDestroy(pc->dm);CHKERRQ(ierr); pc->dm = 0;}
81 
82   if (pc->diagonalscaleright) {ierr = VecDestroy(pc->diagonalscaleright);CHKERRQ(ierr);pc->diagonalscaleright = 0;}
83   if (pc->diagonalscaleleft)  {ierr = VecDestroy(pc->diagonalscaleleft);CHKERRQ(ierr);pc->diagonalscaleleft = 0;}
84 
85   if (pc->pmat) {ierr = MatDestroy(pc->pmat);CHKERRQ(ierr);pc->pmat = 0;}
86   if (pc->mat) {ierr = MatDestroy(pc->mat);CHKERRQ(ierr);pc->mat = 0;}
87   if (pc->ops->reset) {
88     ierr = (*pc->ops->reset)(pc);
89   }
90   pc->setupcalled = 0;
91   PetscFunctionReturn(0);
92 }
93 
94 #undef __FUNCT__
95 #define __FUNCT__ "PCDestroy"
96 /*@
97    PCDestroy - Destroys PC context that was created with PCCreate().
98 
99    Collective on PC
100 
101    Input Parameter:
102 .  pc - the preconditioner context
103 
104    Level: developer
105 
106 .keywords: PC, destroy
107 
108 .seealso: PCCreate(), PCSetUp()
109 @*/
110 PetscErrorCode  PCDestroy(PC pc)
111 {
112   PetscErrorCode ierr;
113 
114   PetscFunctionBegin;
115   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
116   if (--((PetscObject)pc)->refct > 0) PetscFunctionReturn(0);
117 
118   /* if memory was published with AMS then destroy it */
119   ierr = PetscObjectDepublish(pc);CHKERRQ(ierr);
120   if (pc->dm) {ierr = DMDestroy(pc->dm);CHKERRQ(ierr);}
121 
122   if (pc->ops->destroy)       {ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr);}
123   if (pc->diagonalscaleright) {ierr = VecDestroy(pc->diagonalscaleright);CHKERRQ(ierr);}
124   if (pc->diagonalscaleleft)  {ierr = VecDestroy(pc->diagonalscaleleft);CHKERRQ(ierr);}
125 
126   if (pc->pmat) {ierr = MatDestroy(pc->pmat);CHKERRQ(ierr);}
127   if (pc->mat) {ierr = MatDestroy(pc->mat);CHKERRQ(ierr);}
128 
129   ierr = PetscHeaderDestroy(pc);CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 #undef __FUNCT__
134 #define __FUNCT__ "PCGetDiagonalScale"
135 /*@C
136    PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
137       scaling as needed by certain time-stepping codes.
138 
139    Logically Collective on PC
140 
141    Input Parameter:
142 .  pc - the preconditioner context
143 
144    Output Parameter:
145 .  flag - PETSC_TRUE if it applies the scaling
146 
147    Level: developer
148 
149    Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
150 $           D M A D^{-1} y = D M b  for left preconditioning or
151 $           D A M D^{-1} z = D b for right preconditioning
152 
153 .keywords: PC
154 
155 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCSetDiagonalScale()
156 @*/
157 PetscErrorCode  PCGetDiagonalScale(PC pc,PetscBool  *flag)
158 {
159   PetscFunctionBegin;
160   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
161   PetscValidPointer(flag,2);
162   *flag = pc->diagonalscale;
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "PCSetDiagonalScale"
168 /*@
169    PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
170       scaling as needed by certain time-stepping codes.
171 
172    Logically Collective on PC
173 
174    Input Parameters:
175 +  pc - the preconditioner context
176 -  s - scaling vector
177 
178    Level: intermediate
179 
180    Notes: The system solved via the Krylov method is
181 $           D M A D^{-1} y = D M b  for left preconditioning or
182 $           D A M D^{-1} z = D b for right preconditioning
183 
184    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
185 
186 .keywords: PC
187 
188 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCGetDiagonalScale()
189 @*/
190 PetscErrorCode  PCSetDiagonalScale(PC pc,Vec s)
191 {
192   PetscErrorCode ierr;
193 
194   PetscFunctionBegin;
195   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
196   PetscValidHeaderSpecific(s,VEC_CLASSID,2);
197   pc->diagonalscale     = PETSC_TRUE;
198   ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr);
199   if (pc->diagonalscaleleft) {
200     ierr = VecDestroy(pc->diagonalscaleleft);CHKERRQ(ierr);
201   }
202   pc->diagonalscaleleft = s;
203   if (!pc->diagonalscaleright) {
204     ierr = VecDuplicate(s,&pc->diagonalscaleright);CHKERRQ(ierr);
205   }
206   ierr = VecCopy(s,pc->diagonalscaleright);CHKERRQ(ierr);
207   ierr = VecReciprocal(pc->diagonalscaleright);CHKERRQ(ierr);
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "PCDiagonalScaleLeft"
213 /*@
214    PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
215 
216    Logically Collective on PC
217 
218    Input Parameters:
219 +  pc - the preconditioner context
220 .  in - input vector
221 +  out - scaled vector (maybe the same as in)
222 
223    Level: intermediate
224 
225    Notes: The system solved via the Krylov method is
226 $           D M A D^{-1} y = D M b  for left preconditioning or
227 $           D A M D^{-1} z = D b for right preconditioning
228 
229    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
230 
231    If diagonal scaling is turned off and in is not out then in is copied to out
232 
233 .keywords: PC
234 
235 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
236 @*/
237 PetscErrorCode  PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
238 {
239   PetscErrorCode ierr;
240 
241   PetscFunctionBegin;
242   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
243   PetscValidHeaderSpecific(in,VEC_CLASSID,2);
244   PetscValidHeaderSpecific(out,VEC_CLASSID,3);
245   if (pc->diagonalscale) {
246     ierr = VecPointwiseMult(out,pc->diagonalscaleleft,in);CHKERRQ(ierr);
247   } else if (in != out) {
248     ierr = VecCopy(in,out);CHKERRQ(ierr);
249   }
250   PetscFunctionReturn(0);
251 }
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "PCDiagonalScaleRight"
255 /*@
256    PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
257 
258    Logically Collective on PC
259 
260    Input Parameters:
261 +  pc - the preconditioner context
262 .  in - input vector
263 +  out - scaled vector (maybe the same as in)
264 
265    Level: intermediate
266 
267    Notes: The system solved via the Krylov method is
268 $           D M A D^{-1} y = D M b  for left preconditioning or
269 $           D A M D^{-1} z = D b for right preconditioning
270 
271    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
272 
273    If diagonal scaling is turned off and in is not out then in is copied to out
274 
275 .keywords: PC
276 
277 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
278 @*/
279 PetscErrorCode  PCDiagonalScaleRight(PC pc,Vec in,Vec out)
280 {
281   PetscErrorCode ierr;
282 
283   PetscFunctionBegin;
284   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
285   PetscValidHeaderSpecific(in,VEC_CLASSID,2);
286   PetscValidHeaderSpecific(out,VEC_CLASSID,3);
287   if (pc->diagonalscale) {
288     ierr = VecPointwiseMult(out,pc->diagonalscaleright,in);CHKERRQ(ierr);
289   } else if (in != out) {
290     ierr = VecCopy(in,out);CHKERRQ(ierr);
291   }
292   PetscFunctionReturn(0);
293 }
294 
295 #if 0
296 #undef __FUNCT__
297 #define __FUNCT__ "PCPublish_Petsc"
298 static PetscErrorCode PCPublish_Petsc(PetscObject obj)
299 {
300   PetscFunctionBegin;
301   PetscFunctionReturn(0);
302 }
303 #endif
304 
305 #undef __FUNCT__
306 #define __FUNCT__ "PCCreate"
307 /*@
308    PCCreate - Creates a preconditioner context.
309 
310    Collective on MPI_Comm
311 
312    Input Parameter:
313 .  comm - MPI communicator
314 
315    Output Parameter:
316 .  pc - location to put the preconditioner context
317 
318    Notes:
319    The default preconditioner for sparse matrices is PCILU or PCICC with 0 fill on one process and block Jacobi with PCILU or ICC
320    in parallel. For dense matrices it is always PCNONE.
321 
322    Level: developer
323 
324 .keywords: PC, create, context
325 
326 .seealso: PCSetUp(), PCApply(), PCDestroy()
327 @*/
328 PetscErrorCode  PCCreate(MPI_Comm comm,PC *newpc)
329 {
330   PC             pc;
331   PetscErrorCode ierr;
332 
333   PetscFunctionBegin;
334   PetscValidPointer(newpc,1);
335   *newpc = 0;
336 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
337   ierr = PCInitializePackage(PETSC_NULL);CHKERRQ(ierr);
338 #endif
339 
340   ierr = PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_CLASSID,-1,"PC",comm,PCDestroy,PCView);CHKERRQ(ierr);
341 
342   pc->mat                  = 0;
343   pc->pmat                 = 0;
344   pc->setupcalled          = 0;
345   pc->setfromoptionscalled = 0;
346   pc->data                 = 0;
347   pc->diagonalscale        = PETSC_FALSE;
348   pc->diagonalscaleleft    = 0;
349   pc->diagonalscaleright   = 0;
350 
351   pc->modifysubmatrices   = 0;
352   pc->modifysubmatricesP  = 0;
353   *newpc = pc;
354   PetscFunctionReturn(0);
355 
356 }
357 
358 /* -------------------------------------------------------------------------------*/
359 
360 #undef __FUNCT__
361 #define __FUNCT__ "PCApply"
362 /*@
363    PCApply - Applies the preconditioner to a vector.
364 
365    Collective on PC and Vec
366 
367    Input Parameters:
368 +  pc - the preconditioner context
369 -  x - input vector
370 
371    Output Parameter:
372 .  y - output vector
373 
374    Level: developer
375 
376 .keywords: PC, apply
377 
378 .seealso: PCApplyTranspose(), PCApplyBAorAB()
379 @*/
380 PetscErrorCode  PCApply(PC pc,Vec x,Vec y)
381 {
382   PetscErrorCode ierr;
383 
384   PetscFunctionBegin;
385   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
386   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
387   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
388   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
389   if (pc->setupcalled < 2) {
390     ierr = PCSetUp(pc);CHKERRQ(ierr);
391   }
392   if (!pc->ops->apply) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have apply");
393   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
394   ierr = (*pc->ops->apply)(pc,x,y);CHKERRQ(ierr);
395   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
396   PetscFunctionReturn(0);
397 }
398 
399 #undef __FUNCT__
400 #define __FUNCT__ "PCApplySymmetricLeft"
401 /*@
402    PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
403 
404    Collective on PC and Vec
405 
406    Input Parameters:
407 +  pc - the preconditioner context
408 -  x - input vector
409 
410    Output Parameter:
411 .  y - output vector
412 
413    Notes:
414    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
415 
416    Level: developer
417 
418 .keywords: PC, apply, symmetric, left
419 
420 .seealso: PCApply(), PCApplySymmetricRight()
421 @*/
422 PetscErrorCode  PCApplySymmetricLeft(PC pc,Vec x,Vec y)
423 {
424   PetscErrorCode ierr;
425 
426   PetscFunctionBegin;
427   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
428   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
429   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
430   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
431   if (pc->setupcalled < 2) {
432     ierr = PCSetUp(pc);CHKERRQ(ierr);
433   }
434   if (!pc->ops->applysymmetricleft) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have left symmetric apply");
435   ierr = PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
436   ierr = (*pc->ops->applysymmetricleft)(pc,x,y);CHKERRQ(ierr);
437   ierr = PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
438   PetscFunctionReturn(0);
439 }
440 
441 #undef __FUNCT__
442 #define __FUNCT__ "PCApplySymmetricRight"
443 /*@
444    PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
445 
446    Collective on PC and Vec
447 
448    Input Parameters:
449 +  pc - the preconditioner context
450 -  x - input vector
451 
452    Output Parameter:
453 .  y - output vector
454 
455    Level: developer
456 
457    Notes:
458    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
459 
460 .keywords: PC, apply, symmetric, right
461 
462 .seealso: PCApply(), PCApplySymmetricLeft()
463 @*/
464 PetscErrorCode  PCApplySymmetricRight(PC pc,Vec x,Vec y)
465 {
466   PetscErrorCode ierr;
467 
468   PetscFunctionBegin;
469   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
470   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
471   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
472   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
473   if (pc->setupcalled < 2) {
474     ierr = PCSetUp(pc);CHKERRQ(ierr);
475   }
476   if (!pc->ops->applysymmetricright) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have left symmetric apply");
477   ierr = PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
478   ierr = (*pc->ops->applysymmetricright)(pc,x,y);CHKERRQ(ierr);
479   ierr = PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
480   PetscFunctionReturn(0);
481 }
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "PCApplyTranspose"
485 /*@
486    PCApplyTranspose - Applies the transpose of preconditioner to a vector.
487 
488    Collective on PC and Vec
489 
490    Input Parameters:
491 +  pc - the preconditioner context
492 -  x - input vector
493 
494    Output Parameter:
495 .  y - output vector
496 
497    Level: developer
498 
499 .keywords: PC, apply, transpose
500 
501 .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCApplyTransposeExists()
502 @*/
503 PetscErrorCode  PCApplyTranspose(PC pc,Vec x,Vec y)
504 {
505   PetscErrorCode ierr;
506 
507   PetscFunctionBegin;
508   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
509   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
510   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
511   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
512   if (pc->setupcalled < 2) {
513     ierr = PCSetUp(pc);CHKERRQ(ierr);
514   }
515   if (!pc->ops->applytranspose) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have apply transpose");
516   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
517   ierr = (*pc->ops->applytranspose)(pc,x,y);CHKERRQ(ierr);
518   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
519   PetscFunctionReturn(0);
520 }
521 
522 #undef __FUNCT__
523 #define __FUNCT__ "PCApplyTransposeExists"
524 /*@
525    PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
526 
527    Collective on PC and Vec
528 
529    Input Parameters:
530 .  pc - the preconditioner context
531 
532    Output Parameter:
533 .  flg - PETSC_TRUE if a transpose operation is defined
534 
535    Level: developer
536 
537 .keywords: PC, apply, transpose
538 
539 .seealso: PCApplyTranspose()
540 @*/
541 PetscErrorCode  PCApplyTransposeExists(PC pc,PetscBool  *flg)
542 {
543   PetscFunctionBegin;
544   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
545   PetscValidPointer(flg,2);
546   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
547   else                         *flg = PETSC_FALSE;
548   PetscFunctionReturn(0);
549 }
550 
551 #undef __FUNCT__
552 #define __FUNCT__ "PCApplyBAorAB"
553 /*@
554    PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x.
555 
556    Collective on PC and Vec
557 
558    Input Parameters:
559 +  pc - the preconditioner context
560 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
561 .  x - input vector
562 -  work - work vector
563 
564    Output Parameter:
565 .  y - output vector
566 
567    Level: developer
568 
569 .keywords: PC, apply, operator
570 
571 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
572 @*/
573 PetscErrorCode  PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
574 {
575   PetscErrorCode ierr;
576 
577   PetscFunctionBegin;
578   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
579   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
580   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
581   PetscValidHeaderSpecific(work,VEC_CLASSID,5);
582   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
583   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");
584   if (pc->diagonalscale && side == PC_SYMMETRIC) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");
585 
586   if (pc->setupcalled < 2) {
587     ierr = PCSetUp(pc);CHKERRQ(ierr);
588   }
589 
590   if (pc->diagonalscale) {
591     if (pc->ops->applyBA) {
592       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
593       ierr = VecDuplicate(x,&work2);CHKERRQ(ierr);
594       ierr = PCDiagonalScaleRight(pc,x,work2);CHKERRQ(ierr);
595       ierr = (*pc->ops->applyBA)(pc,side,work2,y,work);CHKERRQ(ierr);
596       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
597       ierr = VecDestroy(work2);CHKERRQ(ierr);
598     } else if (side == PC_RIGHT) {
599       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
600       ierr = PCApply(pc,y,work);CHKERRQ(ierr);
601       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
602       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
603     } else if (side == PC_LEFT) {
604       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
605       ierr = MatMult(pc->mat,y,work);CHKERRQ(ierr);
606       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
607       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
608     } else if (side == PC_SYMMETRIC) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
609   } else {
610     if (pc->ops->applyBA) {
611       ierr = (*pc->ops->applyBA)(pc,side,x,y,work);CHKERRQ(ierr);
612     } else if (side == PC_RIGHT) {
613       ierr = PCApply(pc,x,work);CHKERRQ(ierr);
614       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
615     } else if (side == PC_LEFT) {
616       ierr = MatMult(pc->mat,x,work);CHKERRQ(ierr);
617       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
618     } else if (side == PC_SYMMETRIC) {
619       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
620       ierr = PCApplySymmetricRight(pc,x,work);CHKERRQ(ierr);
621       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
622       ierr = VecCopy(y,work);CHKERRQ(ierr);
623       ierr = PCApplySymmetricLeft(pc,work,y);CHKERRQ(ierr);
624     }
625   }
626   PetscFunctionReturn(0);
627 }
628 
629 #undef __FUNCT__
630 #define __FUNCT__ "PCApplyBAorABTranspose"
631 /*@
632    PCApplyBAorABTranspose - Applies the transpose of the preconditioner
633    and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
634    NOT tr(B*A) = tr(A)*tr(B).
635 
636    Collective on PC and Vec
637 
638    Input Parameters:
639 +  pc - the preconditioner context
640 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
641 .  x - input vector
642 -  work - work vector
643 
644    Output Parameter:
645 .  y - output vector
646 
647 
648    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
649       defined by B'. This is why this has the funny form that it computes tr(B) * tr(A)
650 
651     Level: developer
652 
653 .keywords: PC, apply, operator, transpose
654 
655 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
656 @*/
657 PetscErrorCode  PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
658 {
659   PetscErrorCode ierr;
660 
661   PetscFunctionBegin;
662   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
663   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
664   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
665   PetscValidHeaderSpecific(work,VEC_CLASSID,5);
666   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
667   if (pc->ops->applyBAtranspose) {
668     ierr = (*pc->ops->applyBAtranspose)(pc,side,x,y,work);CHKERRQ(ierr);
669     PetscFunctionReturn(0);
670   }
671   if (side != PC_LEFT && side != PC_RIGHT) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
672 
673   if (pc->setupcalled < 2) {
674     ierr = PCSetUp(pc);CHKERRQ(ierr);
675   }
676 
677   if (side == PC_RIGHT) {
678     ierr = PCApplyTranspose(pc,x,work);CHKERRQ(ierr);
679     ierr = MatMultTranspose(pc->mat,work,y);CHKERRQ(ierr);
680   } else if (side == PC_LEFT) {
681     ierr = MatMultTranspose(pc->mat,x,work);CHKERRQ(ierr);
682     ierr = PCApplyTranspose(pc,work,y);CHKERRQ(ierr);
683   }
684   /* add support for PC_SYMMETRIC */
685   PetscFunctionReturn(0); /* actually will never get here */
686 }
687 
688 /* -------------------------------------------------------------------------------*/
689 
690 #undef __FUNCT__
691 #define __FUNCT__ "PCApplyRichardsonExists"
692 /*@
693    PCApplyRichardsonExists - Determines whether a particular preconditioner has a
694    built-in fast application of Richardson's method.
695 
696    Not Collective
697 
698    Input Parameter:
699 .  pc - the preconditioner
700 
701    Output Parameter:
702 .  exists - PETSC_TRUE or PETSC_FALSE
703 
704    Level: developer
705 
706 .keywords: PC, apply, Richardson, exists
707 
708 .seealso: PCApplyRichardson()
709 @*/
710 PetscErrorCode  PCApplyRichardsonExists(PC pc,PetscBool  *exists)
711 {
712   PetscFunctionBegin;
713   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
714   PetscValidIntPointer(exists,2);
715   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
716   else                          *exists = PETSC_FALSE;
717   PetscFunctionReturn(0);
718 }
719 
720 #undef __FUNCT__
721 #define __FUNCT__ "PCApplyRichardson"
722 /*@
723    PCApplyRichardson - Applies several steps of Richardson iteration with
724    the particular preconditioner. This routine is usually used by the
725    Krylov solvers and not the application code directly.
726 
727    Collective on PC
728 
729    Input Parameters:
730 +  pc  - the preconditioner context
731 .  b   - the right hand side
732 .  w   - one work vector
733 .  rtol - relative decrease in residual norm convergence criteria
734 .  abstol - absolute residual norm convergence criteria
735 .  dtol - divergence residual norm increase criteria
736 .  its - the number of iterations to apply.
737 -  guesszero - if the input x contains nonzero initial guess
738 
739    Output Parameter:
740 +  outits - number of iterations actually used (for SOR this always equals its)
741 .  reason - the reason the apply terminated
742 -  y - the solution (also contains initial guess if guesszero is PETSC_FALSE
743 
744    Notes:
745    Most preconditioners do not support this function. Use the command
746    PCApplyRichardsonExists() to determine if one does.
747 
748    Except for the multigrid PC this routine ignores the convergence tolerances
749    and always runs for the number of iterations
750 
751    Level: developer
752 
753 .keywords: PC, apply, Richardson
754 
755 .seealso: PCApplyRichardsonExists()
756 @*/
757 PetscErrorCode  PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool  guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
758 {
759   PetscErrorCode ierr;
760 
761   PetscFunctionBegin;
762   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
763   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
764   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
765   PetscValidHeaderSpecific(w,VEC_CLASSID,4);
766   if (b == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"b and y must be different vectors");
767   if (pc->setupcalled < 2) {
768     ierr = PCSetUp(pc);CHKERRQ(ierr);
769   }
770   if (!pc->ops->applyrichardson) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have apply richardson");
771   ierr = (*pc->ops->applyrichardson)(pc,b,y,w,rtol,abstol,dtol,its,guesszero,outits,reason);CHKERRQ(ierr);
772   PetscFunctionReturn(0);
773 }
774 
775 /*
776       a setupcall of 0 indicates never setup,
777                      1 needs to be resetup,
778                      2 does not need any changes.
779 */
780 #undef __FUNCT__
781 #define __FUNCT__ "PCSetUp"
782 /*@
783    PCSetUp - Prepares for the use of a preconditioner.
784 
785    Collective on PC
786 
787    Input Parameter:
788 .  pc - the preconditioner context
789 
790    Level: developer
791 
792 .keywords: PC, setup
793 
794 .seealso: PCCreate(), PCApply(), PCDestroy()
795 @*/
796 PetscErrorCode  PCSetUp(PC pc)
797 {
798   PetscErrorCode ierr;
799   const char     *def;
800 
801   PetscFunctionBegin;
802   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
803   if (!pc->mat) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");
804 
805   if (pc->setupcalled > 1) {
806     ierr = PetscInfo(pc,"Setting PC with identical preconditioner\n");CHKERRQ(ierr);
807     PetscFunctionReturn(0);
808   } else if (!pc->setupcalled) {
809     ierr = PetscInfo(pc,"Setting up new PC\n");CHKERRQ(ierr);
810   } else if (pc->flag == SAME_NONZERO_PATTERN) {
811     ierr = PetscInfo(pc,"Setting up PC with same nonzero pattern\n");CHKERRQ(ierr);
812   } else {
813     ierr = PetscInfo(pc,"Setting up PC with different nonzero pattern\n");CHKERRQ(ierr);
814   }
815 
816   if (!((PetscObject)pc)->type_name) {
817     ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr);
818     ierr = PCSetType(pc,def);CHKERRQ(ierr);
819   }
820 
821   ierr = PetscLogEventBegin(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
822   if (pc->ops->setup) {
823     ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr);
824   }
825   pc->setupcalled = 2;
826   ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
827   PetscFunctionReturn(0);
828 }
829 
830 #undef __FUNCT__
831 #define __FUNCT__ "PCSetUpOnBlocks"
832 /*@
833    PCSetUpOnBlocks - Sets up the preconditioner for each block in
834    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
835    methods.
836 
837    Collective on PC
838 
839    Input Parameters:
840 .  pc - the preconditioner context
841 
842    Level: developer
843 
844 .keywords: PC, setup, blocks
845 
846 .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
847 @*/
848 PetscErrorCode  PCSetUpOnBlocks(PC pc)
849 {
850   PetscErrorCode ierr;
851 
852   PetscFunctionBegin;
853   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
854   if (!pc->ops->setuponblocks) PetscFunctionReturn(0);
855   ierr = PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
856   ierr = (*pc->ops->setuponblocks)(pc);CHKERRQ(ierr);
857   ierr = PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
858   PetscFunctionReturn(0);
859 }
860 
861 #undef __FUNCT__
862 #define __FUNCT__ "PCSetModifySubMatrices"
863 /*@C
864    PCSetModifySubMatrices - Sets a user-defined routine for modifying the
865    submatrices that arise within certain subdomain-based preconditioners.
866    The basic submatrices are extracted from the preconditioner matrix as
867    usual; the user can then alter these (for example, to set different boundary
868    conditions for each submatrix) before they are used for the local solves.
869 
870    Logically Collective on PC
871 
872    Input Parameters:
873 +  pc - the preconditioner context
874 .  func - routine for modifying the submatrices
875 -  ctx - optional user-defined context (may be null)
876 
877    Calling sequence of func:
878 $     func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
879 
880 .  row - an array of index sets that contain the global row numbers
881          that comprise each local submatrix
882 .  col - an array of index sets that contain the global column numbers
883          that comprise each local submatrix
884 .  submat - array of local submatrices
885 -  ctx - optional user-defined context for private data for the
886          user-defined func routine (may be null)
887 
888    Notes:
889    PCSetModifySubMatrices() MUST be called before KSPSetUp() and
890    KSPSolve().
891 
892    A routine set by PCSetModifySubMatrices() is currently called within
893    the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
894    preconditioners.  All other preconditioners ignore this routine.
895 
896    Level: advanced
897 
898 .keywords: PC, set, modify, submatrices
899 
900 .seealso: PCModifySubMatrices()
901 @*/
902 PetscErrorCode  PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
903 {
904   PetscFunctionBegin;
905   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
906   pc->modifysubmatrices  = func;
907   pc->modifysubmatricesP = ctx;
908   PetscFunctionReturn(0);
909 }
910 
911 #undef __FUNCT__
912 #define __FUNCT__ "PCModifySubMatrices"
913 /*@C
914    PCModifySubMatrices - Calls an optional user-defined routine within
915    certain preconditioners if one has been set with PCSetModifySubMarices().
916 
917    Collective on PC
918 
919    Input Parameters:
920 +  pc - the preconditioner context
921 .  nsub - the number of local submatrices
922 .  row - an array of index sets that contain the global row numbers
923          that comprise each local submatrix
924 .  col - an array of index sets that contain the global column numbers
925          that comprise each local submatrix
926 .  submat - array of local submatrices
927 -  ctx - optional user-defined context for private data for the
928          user-defined routine (may be null)
929 
930    Output Parameter:
931 .  submat - array of local submatrices (the entries of which may
932             have been modified)
933 
934    Notes:
935    The user should NOT generally call this routine, as it will
936    automatically be called within certain preconditioners (currently
937    block Jacobi, additive Schwarz) if set.
938 
939    The basic submatrices are extracted from the preconditioner matrix
940    as usual; the user can then alter these (for example, to set different
941    boundary conditions for each submatrix) before they are used for the
942    local solves.
943 
944    Level: developer
945 
946 .keywords: PC, modify, submatrices
947 
948 .seealso: PCSetModifySubMatrices()
949 @*/
950 PetscErrorCode  PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
951 {
952   PetscErrorCode ierr;
953 
954   PetscFunctionBegin;
955   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
956   if (!pc->modifysubmatrices) PetscFunctionReturn(0);
957   ierr = PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
958   ierr = (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);CHKERRQ(ierr);
959   ierr = PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
960   PetscFunctionReturn(0);
961 }
962 
963 #undef __FUNCT__
964 #define __FUNCT__ "PCSetOperators"
965 /*@
966    PCSetOperators - Sets the matrix associated with the linear system and
967    a (possibly) different one associated with the preconditioner.
968 
969    Logically Collective on PC and Mat
970 
971    Input Parameters:
972 +  pc - the preconditioner context
973 .  Amat - the matrix associated with the linear system
974 .  Pmat - the matrix to be used in constructing the preconditioner, usually
975           the same as Amat.
976 -  flag - flag indicating information about the preconditioner matrix structure
977    during successive linear solves.  This flag is ignored the first time a
978    linear system is solved, and thus is irrelevant when solving just one linear
979    system.
980 
981    Notes:
982    The flag can be used to eliminate unnecessary work in the preconditioner
983    during the repeated solution of linear systems of the same size.  The
984    available options are
985 +    SAME_PRECONDITIONER -
986        Pmat is identical during successive linear solves.
987        This option is intended for folks who are using
988        different Amat and Pmat matrices and wish to reuse the
989        same preconditioner matrix.  For example, this option
990        saves work by not recomputing incomplete factorization
991        for ILU/ICC preconditioners.
992 .     SAME_NONZERO_PATTERN -
993        Pmat has the same nonzero structure during
994        successive linear solves.
995 -     DIFFERENT_NONZERO_PATTERN -
996        Pmat does not have the same nonzero structure.
997 
998     Passing a PETSC_NULL for Amat or Pmat removes the matrix that is currently used.
999 
1000     If you wish to replace either Amat or Pmat but leave the other one untouched then
1001     first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
1002     on it and then pass it back in in your call to KSPSetOperators().
1003 
1004    Caution:
1005    If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
1006    and does not check the structure of the matrix.  If you erroneously
1007    claim that the structure is the same when it actually is not, the new
1008    preconditioner will not function correctly.  Thus, use this optimization
1009    feature carefully!
1010 
1011    If in doubt about whether your preconditioner matrix has changed
1012    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
1013 
1014    More Notes about Repeated Solution of Linear Systems:
1015    PETSc does NOT reset the matrix entries of either Amat or Pmat
1016    to zero after a linear solve; the user is completely responsible for
1017    matrix assembly.  See the routine MatZeroEntries() if desiring to
1018    zero all elements of a matrix.
1019 
1020    Level: intermediate
1021 
1022 .keywords: PC, set, operators, matrix, linear system
1023 
1024 .seealso: PCGetOperators(), MatZeroEntries()
1025  @*/
1026 PetscErrorCode  PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
1027 {
1028   PetscErrorCode ierr;
1029   PetscInt       m1,n1,m2,n2;
1030 
1031   PetscFunctionBegin;
1032   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1033   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1034   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3);
1035   if (Amat) PetscCheckSameComm(pc,1,Amat,2);
1036   if (Pmat) PetscCheckSameComm(pc,1,Pmat,3);
1037   if (pc->setupcalled && Amat && Pmat) {
1038     ierr = MatGetLocalSize(Amat,&m1,&n1);CHKERRQ(ierr);
1039     ierr = MatGetLocalSize(pc->mat,&m2,&n2);CHKERRQ(ierr);
1040     if (m1 != m2 || n1 != n2) {
1041       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);
1042     }
1043     ierr = MatGetLocalSize(Pmat,&m1,&n1);CHKERRQ(ierr);
1044     ierr = MatGetLocalSize(pc->pmat,&m2,&n2);CHKERRQ(ierr);
1045     if (m1 != m2 || n1 != n2) {
1046       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);
1047     }
1048   }
1049 
1050   /* reference first in case the matrices are the same */
1051   if (Amat) {ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);}
1052   if (pc->mat) {ierr = MatDestroy(pc->mat);CHKERRQ(ierr);}
1053   if (Pmat) {ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);}
1054   if (pc->pmat) {ierr = MatDestroy(pc->pmat);CHKERRQ(ierr);}
1055   pc->mat  = Amat;
1056   pc->pmat = Pmat;
1057 
1058   if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1059     pc->setupcalled = 1;
1060   }
1061   pc->flag = flag;
1062   PetscFunctionReturn(0);
1063 }
1064 
1065 #undef __FUNCT__
1066 #define __FUNCT__ "PCGetOperators"
1067 /*@C
1068    PCGetOperators - Gets the matrix associated with the linear system and
1069    possibly a different one associated with the preconditioner.
1070 
1071    Not collective, though parallel Mats are returned if the PC is parallel
1072 
1073    Input Parameter:
1074 .  pc - the preconditioner context
1075 
1076    Output Parameters:
1077 +  mat - the matrix associated with the linear system
1078 .  pmat - matrix associated with the preconditioner, usually the same
1079           as mat.
1080 -  flag - flag indicating information about the preconditioner
1081           matrix structure.  See PCSetOperators() for details.
1082 
1083    Level: intermediate
1084 
1085    Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1086       are created in PC and returned to the user. In this case, if both operators
1087       mat and pmat are requested, two DIFFERENT operators will be returned. If
1088       only one is requested both operators in the PC will be the same (i.e. as
1089       if one had called KSP/PCSetOperators() with the same argument for both Mats).
1090       The user must set the sizes of the returned matrices and their type etc just
1091       as if the user created them with MatCreate(). For example,
1092 
1093 $         KSP/PCGetOperators(ksp/pc,&mat,PETSC_NULL,PETSC_NULL); is equivalent to
1094 $           set size, type, etc of mat
1095 
1096 $         MatCreate(comm,&mat);
1097 $         KSP/PCSetOperators(ksp/pc,mat,mat,SAME_NONZERO_PATTERN);
1098 $         PetscObjectDereference((PetscObject)mat);
1099 $           set size, type, etc of mat
1100 
1101      and
1102 
1103 $         KSP/PCGetOperators(ksp/pc,&mat,&pmat,PETSC_NULL); is equivalent to
1104 $           set size, type, etc of mat and pmat
1105 
1106 $         MatCreate(comm,&mat);
1107 $         MatCreate(comm,&pmat);
1108 $         KSP/PCSetOperators(ksp/pc,mat,pmat,SAME_NONZERO_PATTERN);
1109 $         PetscObjectDereference((PetscObject)mat);
1110 $         PetscObjectDereference((PetscObject)pmat);
1111 $           set size, type, etc of mat and pmat
1112 
1113     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1114     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1115     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1116     at this is when you create a SNES you do not NEED to create a KSP and attach it to
1117     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1118     you do not need to attach a PC to it (the KSP object manages the PC object for you).
1119     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1120     it can be created for you?
1121 
1122 
1123 .keywords: PC, get, operators, matrix, linear system
1124 
1125 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1126 @*/
1127 PetscErrorCode  PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1128 {
1129   PetscErrorCode ierr;
1130 
1131   PetscFunctionBegin;
1132   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1133   if (mat) {
1134     if (!pc->mat) {
1135       ierr = MatCreate(((PetscObject)pc)->comm,&pc->mat);CHKERRQ(ierr);
1136       if (!pc->pmat && !pmat) { /* user did NOT request pmat, so make same as mat */
1137         pc->pmat = pc->mat;
1138         ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1139       }
1140     }
1141     *mat  = pc->mat;
1142   }
1143   if (pmat) {
1144     if (!pc->pmat) {
1145       ierr = MatCreate(((PetscObject)pc)->comm,&pc->mat);CHKERRQ(ierr);
1146       if (!pc->mat && !mat) { /* user did NOT request mat, so make same as pmat */
1147         pc->mat = pc->pmat;
1148         ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1149       }
1150     }
1151     *pmat = pc->pmat;
1152   }
1153   if (flag) *flag = pc->flag;
1154   PetscFunctionReturn(0);
1155 }
1156 
1157 #undef __FUNCT__
1158 #define __FUNCT__ "PCGetOperatorsSet"
1159 /*@C
1160    PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1161    possibly a different one associated with the preconditioner have been set in the PC.
1162 
1163    Not collective, though the results on all processes should be the same
1164 
1165    Input Parameter:
1166 .  pc - the preconditioner context
1167 
1168    Output Parameters:
1169 +  mat - the matrix associated with the linear system was set
1170 -  pmat - matrix associated with the preconditioner was set, usually the same
1171 
1172    Level: intermediate
1173 
1174 .keywords: PC, get, operators, matrix, linear system
1175 
1176 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1177 @*/
1178 PetscErrorCode  PCGetOperatorsSet(PC pc,PetscBool  *mat,PetscBool  *pmat)
1179 {
1180   PetscFunctionBegin;
1181   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1182   if (mat)  *mat  = (pc->mat)  ? PETSC_TRUE : PETSC_FALSE;
1183   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1184   PetscFunctionReturn(0);
1185 }
1186 
1187 #undef __FUNCT__
1188 #define __FUNCT__ "PCFactorGetMatrix"
1189 /*@
1190    PCFactorGetMatrix - Gets the factored matrix from the
1191    preconditioner context.  This routine is valid only for the LU,
1192    incomplete LU, Cholesky, and incomplete Cholesky methods.
1193 
1194    Not Collective on PC though Mat is parallel if PC is parallel
1195 
1196    Input Parameters:
1197 .  pc - the preconditioner context
1198 
1199    Output parameters:
1200 .  mat - the factored matrix
1201 
1202    Level: advanced
1203 
1204    Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1205 
1206 .keywords: PC, get, factored, matrix
1207 @*/
1208 PetscErrorCode  PCFactorGetMatrix(PC pc,Mat *mat)
1209 {
1210   PetscErrorCode ierr;
1211 
1212   PetscFunctionBegin;
1213   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1214   PetscValidPointer(mat,2);
1215   if (pc->ops->getfactoredmatrix) {
1216     ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr);
1217   } else {
1218     SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC type does not support getting factor matrix");
1219   }
1220   PetscFunctionReturn(0);
1221 }
1222 
1223 #undef __FUNCT__
1224 #define __FUNCT__ "PCSetOptionsPrefix"
1225 /*@C
1226    PCSetOptionsPrefix - Sets the prefix used for searching for all
1227    PC options in the database.
1228 
1229    Logically Collective on PC
1230 
1231    Input Parameters:
1232 +  pc - the preconditioner context
1233 -  prefix - the prefix string to prepend to all PC option requests
1234 
1235    Notes:
1236    A hyphen (-) must NOT be given at the beginning of the prefix name.
1237    The first character of all runtime options is AUTOMATICALLY the
1238    hyphen.
1239 
1240    Level: advanced
1241 
1242 .keywords: PC, set, options, prefix, database
1243 
1244 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1245 @*/
1246 PetscErrorCode  PCSetOptionsPrefix(PC pc,const char prefix[])
1247 {
1248   PetscErrorCode ierr;
1249 
1250   PetscFunctionBegin;
1251   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1252   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1253   PetscFunctionReturn(0);
1254 }
1255 
1256 #undef __FUNCT__
1257 #define __FUNCT__ "PCAppendOptionsPrefix"
1258 /*@C
1259    PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1260    PC options in the database.
1261 
1262    Logically Collective on PC
1263 
1264    Input Parameters:
1265 +  pc - the preconditioner context
1266 -  prefix - the prefix string to prepend to all PC option requests
1267 
1268    Notes:
1269    A hyphen (-) must NOT be given at the beginning of the prefix name.
1270    The first character of all runtime options is AUTOMATICALLY the
1271    hyphen.
1272 
1273    Level: advanced
1274 
1275 .keywords: PC, append, options, prefix, database
1276 
1277 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1278 @*/
1279 PetscErrorCode  PCAppendOptionsPrefix(PC pc,const char prefix[])
1280 {
1281   PetscErrorCode ierr;
1282 
1283   PetscFunctionBegin;
1284   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1285   ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1286   PetscFunctionReturn(0);
1287 }
1288 
1289 #undef __FUNCT__
1290 #define __FUNCT__ "PCGetOptionsPrefix"
1291 /*@C
1292    PCGetOptionsPrefix - Gets the prefix used for searching for all
1293    PC options in the database.
1294 
1295    Not Collective
1296 
1297    Input Parameters:
1298 .  pc - the preconditioner context
1299 
1300    Output Parameters:
1301 .  prefix - pointer to the prefix string used, is returned
1302 
1303    Notes: On the fortran side, the user should pass in a string 'prifix' of
1304    sufficient length to hold the prefix.
1305 
1306    Level: advanced
1307 
1308 .keywords: PC, get, options, prefix, database
1309 
1310 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1311 @*/
1312 PetscErrorCode  PCGetOptionsPrefix(PC pc,const char *prefix[])
1313 {
1314   PetscErrorCode ierr;
1315 
1316   PetscFunctionBegin;
1317   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1318   PetscValidPointer(prefix,2);
1319   ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1320   PetscFunctionReturn(0);
1321 }
1322 
1323 #undef __FUNCT__
1324 #define __FUNCT__ "PCPreSolve"
1325 /*@
1326    PCPreSolve - Optional pre-solve phase, intended for any
1327    preconditioner-specific actions that must be performed before
1328    the iterative solve itself.
1329 
1330    Collective on PC
1331 
1332    Input Parameters:
1333 +  pc - the preconditioner context
1334 -  ksp - the Krylov subspace context
1335 
1336    Level: developer
1337 
1338    Sample of Usage:
1339 .vb
1340     PCPreSolve(pc,ksp);
1341     KSPSolve(ksp,b,x);
1342     PCPostSolve(pc,ksp);
1343 .ve
1344 
1345    Notes:
1346    The pre-solve phase is distinct from the PCSetUp() phase.
1347 
1348    KSPSolve() calls this directly, so is rarely called by the user.
1349 
1350 .keywords: PC, pre-solve
1351 
1352 .seealso: PCPostSolve()
1353 @*/
1354 PetscErrorCode  PCPreSolve(PC pc,KSP ksp)
1355 {
1356   PetscErrorCode ierr;
1357   Vec            x,rhs;
1358   Mat            A,B;
1359 
1360   PetscFunctionBegin;
1361   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1362   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1363   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1364   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1365   /*
1366       Scale the system and have the matrices use the scaled form
1367     only if the two matrices are actually the same (and hence
1368     have the same scaling
1369   */
1370   ierr = PCGetOperators(pc,&A,&B,PETSC_NULL);CHKERRQ(ierr);
1371   if (A == B) {
1372     ierr = MatScaleSystem(pc->mat,rhs,x);CHKERRQ(ierr);
1373     ierr = MatUseScaledForm(pc->mat,PETSC_TRUE);CHKERRQ(ierr);
1374   }
1375 
1376   if (pc->ops->presolve) {
1377     ierr = (*pc->ops->presolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1378   }
1379   PetscFunctionReturn(0);
1380 }
1381 
1382 #undef __FUNCT__
1383 #define __FUNCT__ "PCPostSolve"
1384 /*@
1385    PCPostSolve - Optional post-solve phase, intended for any
1386    preconditioner-specific actions that must be performed after
1387    the iterative solve itself.
1388 
1389    Collective on PC
1390 
1391    Input Parameters:
1392 +  pc - the preconditioner context
1393 -  ksp - the Krylov subspace context
1394 
1395    Sample of Usage:
1396 .vb
1397     PCPreSolve(pc,ksp);
1398     KSPSolve(ksp,b,x);
1399     PCPostSolve(pc,ksp);
1400 .ve
1401 
1402    Note:
1403    KSPSolve() calls this routine directly, so it is rarely called by the user.
1404 
1405    Level: developer
1406 
1407 .keywords: PC, post-solve
1408 
1409 .seealso: PCPreSolve(), KSPSolve()
1410 @*/
1411 PetscErrorCode  PCPostSolve(PC pc,KSP ksp)
1412 {
1413   PetscErrorCode ierr;
1414   Vec            x,rhs;
1415   Mat            A,B;
1416 
1417   PetscFunctionBegin;
1418   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1419   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1420   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1421   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1422   if (pc->ops->postsolve) {
1423     ierr =  (*pc->ops->postsolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1424   }
1425   /*
1426       Scale the system and have the matrices use the scaled form
1427     only if the two matrices are actually the same (and hence
1428     have the same scaling
1429   */
1430   ierr = PCGetOperators(pc,&A,&B,PETSC_NULL);CHKERRQ(ierr);
1431   if (A == B) {
1432     ierr = MatUnScaleSystem(pc->mat,rhs,x);CHKERRQ(ierr);
1433     ierr = MatUseScaledForm(pc->mat,PETSC_FALSE);CHKERRQ(ierr);
1434   }
1435   PetscFunctionReturn(0);
1436 }
1437 
1438 #undef __FUNCT__
1439 #define __FUNCT__ "PCView"
1440 /*@C
1441    PCView - Prints the PC data structure.
1442 
1443    Collective on PC
1444 
1445    Input Parameters:
1446 +  PC - the PC context
1447 -  viewer - optional visualization context
1448 
1449    Note:
1450    The available visualization contexts include
1451 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1452 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1453          output where only the first processor opens
1454          the file.  All other processors send their
1455          data to the first processor to print.
1456 
1457    The user can open an alternative visualization contexts with
1458    PetscViewerASCIIOpen() (output to a specified file).
1459 
1460    Level: developer
1461 
1462 .keywords: PC, view
1463 
1464 .seealso: KSPView(), PetscViewerASCIIOpen()
1465 @*/
1466 PetscErrorCode  PCView(PC pc,PetscViewer viewer)
1467 {
1468   const PCType      cstr;
1469   PetscErrorCode    ierr;
1470   PetscBool         iascii,isstring;
1471   PetscViewerFormat format;
1472 
1473   PetscFunctionBegin;
1474   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1475   if (!viewer) {
1476     ierr = PetscViewerASCIIGetStdout(((PetscObject)pc)->comm,&viewer);CHKERRQ(ierr);
1477   }
1478   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1479   PetscCheckSameComm(pc,1,viewer,2);
1480 
1481   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1482   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
1483   if (iascii) {
1484     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1485     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)pc,viewer,"PC Object");CHKERRQ(ierr);
1486     if (pc->ops->view) {
1487       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1488       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1489       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1490     }
1491     if (pc->mat) {
1492       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
1493       if (pc->pmat == pc->mat) {
1494         ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix = precond matrix:\n");CHKERRQ(ierr);
1495         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1496         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1497         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1498       } else {
1499         if (pc->pmat) {
1500           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr);
1501         } else {
1502           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix:\n");CHKERRQ(ierr);
1503         }
1504         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1505         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1506         if (pc->pmat) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1507         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1508       }
1509       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1510     }
1511   } else if (isstring) {
1512     ierr = PCGetType(pc,&cstr);CHKERRQ(ierr);
1513     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",cstr);CHKERRQ(ierr);
1514     if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);}
1515   } else {
1516     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1517   }
1518   PetscFunctionReturn(0);
1519 }
1520 
1521 
1522 #undef __FUNCT__
1523 #define __FUNCT__ "PCSetInitialGuessNonzero"
1524 /*@
1525    PCSetInitialGuessNonzero - Tells the iterative solver that the
1526    initial guess is nonzero; otherwise PC assumes the initial guess
1527    is to be zero (and thus zeros it out before solving).
1528 
1529    Logically Collective on PC
1530 
1531    Input Parameters:
1532 +  pc - iterative context obtained from PCCreate()
1533 -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1534 
1535    Level: Developer
1536 
1537    Notes:
1538     This is a weird function. Since PC's are linear operators on the right hand side they
1539     CANNOT use an initial guess. This function is for the "pass-through" preconditioners
1540     PCKSP, PCREDUNDANT and PCOPENMP and causes the inner KSP object to use the nonzero
1541     initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP.
1542 
1543 
1544 .keywords: PC, set, initial guess, nonzero
1545 
1546 .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
1547 @*/
1548 PetscErrorCode  PCSetInitialGuessNonzero(PC pc,PetscBool  flg)
1549 {
1550   PetscFunctionBegin;
1551   PetscValidLogicalCollectiveBool(pc,flg,2);
1552   pc->nonzero_guess   = flg;
1553   PetscFunctionReturn(0);
1554 }
1555 
1556 #undef __FUNCT__
1557 #define __FUNCT__ "PCRegister"
1558 /*@C
1559   PCRegister - See PCRegisterDynamic()
1560 
1561   Level: advanced
1562 @*/
1563 PetscErrorCode  PCRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PC))
1564 {
1565   PetscErrorCode ierr;
1566   char           fullname[PETSC_MAX_PATH_LEN];
1567 
1568   PetscFunctionBegin;
1569 
1570   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1571   ierr = PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1572   PetscFunctionReturn(0);
1573 }
1574 
1575 #undef __FUNCT__
1576 #define __FUNCT__ "PCComputeExplicitOperator"
1577 /*@
1578     PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1579 
1580     Collective on PC
1581 
1582     Input Parameter:
1583 .   pc - the preconditioner object
1584 
1585     Output Parameter:
1586 .   mat - the explict preconditioned operator
1587 
1588     Notes:
1589     This computation is done by applying the operators to columns of the
1590     identity matrix.
1591 
1592     Currently, this routine uses a dense matrix format when 1 processor
1593     is used and a sparse format otherwise.  This routine is costly in general,
1594     and is recommended for use only with relatively small systems.
1595 
1596     Level: advanced
1597 
1598 .keywords: PC, compute, explicit, operator
1599 
1600 .seealso: KSPComputeExplicitOperator()
1601 
1602 @*/
1603 PetscErrorCode  PCComputeExplicitOperator(PC pc,Mat *mat)
1604 {
1605   Vec            in,out;
1606   PetscErrorCode ierr;
1607   PetscInt       i,M,m,*rows,start,end;
1608   PetscMPIInt    size;
1609   MPI_Comm       comm;
1610   PetscScalar    *array,one = 1.0;
1611 
1612   PetscFunctionBegin;
1613   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1614   PetscValidPointer(mat,2);
1615 
1616   comm = ((PetscObject)pc)->comm;
1617   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1618 
1619   if (!pc->pmat) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1620   ierr = MatGetVecs(pc->pmat,&in,0);CHKERRQ(ierr);
1621   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
1622   ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr);
1623   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
1624   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
1625   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&rows);CHKERRQ(ierr);
1626   for (i=0; i<m; i++) {rows[i] = start + i;}
1627 
1628   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
1629   ierr = MatSetSizes(*mat,m,m,M,M);CHKERRQ(ierr);
1630   if (size == 1) {
1631     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
1632     ierr = MatSeqDenseSetPreallocation(*mat,PETSC_NULL);CHKERRQ(ierr);
1633   } else {
1634     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
1635     ierr = MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1636   }
1637 
1638   for (i=0; i<M; i++) {
1639 
1640     ierr = VecSet(in,0.0);CHKERRQ(ierr);
1641     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
1642     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
1643     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
1644 
1645     /* should fix, allowing user to choose side */
1646     ierr = PCApply(pc,in,out);CHKERRQ(ierr);
1647 
1648     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
1649     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
1650     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
1651 
1652   }
1653   ierr = PetscFree(rows);CHKERRQ(ierr);
1654   ierr = VecDestroy(out);CHKERRQ(ierr);
1655   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1656   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1657   PetscFunctionReturn(0);
1658 }
1659 
1660