xref: /petsc/src/ksp/pc/interface/precon.c (revision 85e2c93f7d5111e4f7e392c0d6d8b2a65248dd97)
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->diagonalscaleright) {ierr = VecDestroy(pc->diagonalscaleright);CHKERRQ(ierr);}
81   if (pc->diagonalscaleleft)  {ierr = VecDestroy(pc->diagonalscaleleft);CHKERRQ(ierr);}
82   if (pc->pmat) {ierr = MatDestroy(pc->pmat);CHKERRQ(ierr);}
83   if (pc->mat) {ierr = MatDestroy(pc->mat);CHKERRQ(ierr);}
84   if (pc->ops->reset) {
85     ierr = (*pc->ops->reset)(pc);
86   }
87   pc->setupcalled = 0;
88   PetscFunctionReturn(0);
89 }
90 
91 #undef __FUNCT__
92 #define __FUNCT__ "PCDestroy"
93 /*@
94    PCDestroy - Destroys PC context that was created with PCCreate().
95 
96    Collective on PC
97 
98    Input Parameter:
99 .  pc - the preconditioner context
100 
101    Level: developer
102 
103 .keywords: PC, destroy
104 
105 .seealso: PCCreate(), PCSetUp()
106 @*/
107 PetscErrorCode  PCDestroy(PC pc)
108 {
109   PetscErrorCode ierr;
110 
111   PetscFunctionBegin;
112   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
113   if (--((PetscObject)pc)->refct > 0) PetscFunctionReturn(0);
114 
115   ierr = PetscObjectDepublish(pc);CHKERRQ(ierr);
116   ierr = PCReset(pc);CHKERRQ(ierr);
117   if (pc->ops->destroy) {ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr);}
118   if (pc->dm) {ierr = DMDestroy(pc->dm);CHKERRQ(ierr);}
119   ierr = PetscHeaderDestroy(pc);CHKERRQ(ierr);
120   PetscFunctionReturn(0);
121 }
122 
123 #undef __FUNCT__
124 #define __FUNCT__ "PCGetDiagonalScale"
125 /*@C
126    PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
127       scaling as needed by certain time-stepping codes.
128 
129    Logically Collective on PC
130 
131    Input Parameter:
132 .  pc - the preconditioner context
133 
134    Output Parameter:
135 .  flag - PETSC_TRUE if it applies the scaling
136 
137    Level: developer
138 
139    Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
140 $           D M A D^{-1} y = D M b  for left preconditioning or
141 $           D A M D^{-1} z = D b for right preconditioning
142 
143 .keywords: PC
144 
145 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCSetDiagonalScale()
146 @*/
147 PetscErrorCode  PCGetDiagonalScale(PC pc,PetscBool  *flag)
148 {
149   PetscFunctionBegin;
150   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
151   PetscValidPointer(flag,2);
152   *flag = pc->diagonalscale;
153   PetscFunctionReturn(0);
154 }
155 
156 #undef __FUNCT__
157 #define __FUNCT__ "PCSetDiagonalScale"
158 /*@
159    PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
160       scaling as needed by certain time-stepping codes.
161 
162    Logically Collective on PC
163 
164    Input Parameters:
165 +  pc - the preconditioner context
166 -  s - scaling vector
167 
168    Level: intermediate
169 
170    Notes: The system solved via the Krylov method is
171 $           D M A D^{-1} y = D M b  for left preconditioning or
172 $           D A M D^{-1} z = D b for right preconditioning
173 
174    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
175 
176 .keywords: PC
177 
178 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCGetDiagonalScale()
179 @*/
180 PetscErrorCode  PCSetDiagonalScale(PC pc,Vec s)
181 {
182   PetscErrorCode ierr;
183 
184   PetscFunctionBegin;
185   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
186   PetscValidHeaderSpecific(s,VEC_CLASSID,2);
187   pc->diagonalscale     = PETSC_TRUE;
188   ierr = PetscObjectReference((PetscObject)s);CHKERRQ(ierr);
189   if (pc->diagonalscaleleft) {
190     ierr = VecDestroy(pc->diagonalscaleleft);CHKERRQ(ierr);
191   }
192   pc->diagonalscaleleft = s;
193   if (!pc->diagonalscaleright) {
194     ierr = VecDuplicate(s,&pc->diagonalscaleright);CHKERRQ(ierr);
195   }
196   ierr = VecCopy(s,pc->diagonalscaleright);CHKERRQ(ierr);
197   ierr = VecReciprocal(pc->diagonalscaleright);CHKERRQ(ierr);
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "PCDiagonalScaleLeft"
203 /*@
204    PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
205 
206    Logically Collective on PC
207 
208    Input Parameters:
209 +  pc - the preconditioner context
210 .  in - input vector
211 +  out - scaled vector (maybe the same as in)
212 
213    Level: intermediate
214 
215    Notes: The system solved via the Krylov method is
216 $           D M A D^{-1} y = D M b  for left preconditioning or
217 $           D A M D^{-1} z = D b for right preconditioning
218 
219    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
220 
221    If diagonal scaling is turned off and in is not out then in is copied to out
222 
223 .keywords: PC
224 
225 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
226 @*/
227 PetscErrorCode  PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
228 {
229   PetscErrorCode ierr;
230 
231   PetscFunctionBegin;
232   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
233   PetscValidHeaderSpecific(in,VEC_CLASSID,2);
234   PetscValidHeaderSpecific(out,VEC_CLASSID,3);
235   if (pc->diagonalscale) {
236     ierr = VecPointwiseMult(out,pc->diagonalscaleleft,in);CHKERRQ(ierr);
237   } else if (in != out) {
238     ierr = VecCopy(in,out);CHKERRQ(ierr);
239   }
240   PetscFunctionReturn(0);
241 }
242 
243 #undef __FUNCT__
244 #define __FUNCT__ "PCDiagonalScaleRight"
245 /*@
246    PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
247 
248    Logically Collective on PC
249 
250    Input Parameters:
251 +  pc - the preconditioner context
252 .  in - input vector
253 +  out - scaled vector (maybe the same as in)
254 
255    Level: intermediate
256 
257    Notes: The system solved via the Krylov method is
258 $           D M A D^{-1} y = D M b  for left preconditioning or
259 $           D A M D^{-1} z = D b for right preconditioning
260 
261    PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
262 
263    If diagonal scaling is turned off and in is not out then in is copied to out
264 
265 .keywords: PC
266 
267 .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
268 @*/
269 PetscErrorCode  PCDiagonalScaleRight(PC pc,Vec in,Vec out)
270 {
271   PetscErrorCode ierr;
272 
273   PetscFunctionBegin;
274   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
275   PetscValidHeaderSpecific(in,VEC_CLASSID,2);
276   PetscValidHeaderSpecific(out,VEC_CLASSID,3);
277   if (pc->diagonalscale) {
278     ierr = VecPointwiseMult(out,pc->diagonalscaleright,in);CHKERRQ(ierr);
279   } else if (in != out) {
280     ierr = VecCopy(in,out);CHKERRQ(ierr);
281   }
282   PetscFunctionReturn(0);
283 }
284 
285 #if 0
286 #undef __FUNCT__
287 #define __FUNCT__ "PCPublish_Petsc"
288 static PetscErrorCode PCPublish_Petsc(PetscObject obj)
289 {
290   PetscFunctionBegin;
291   PetscFunctionReturn(0);
292 }
293 #endif
294 
295 #undef __FUNCT__
296 #define __FUNCT__ "PCCreate"
297 /*@
298    PCCreate - Creates a preconditioner context.
299 
300    Collective on MPI_Comm
301 
302    Input Parameter:
303 .  comm - MPI communicator
304 
305    Output Parameter:
306 .  pc - location to put the preconditioner context
307 
308    Notes:
309    The default preconditioner for sparse matrices is PCILU or PCICC with 0 fill on one process and block Jacobi with PCILU or ICC
310    in parallel. For dense matrices it is always PCNONE.
311 
312    Level: developer
313 
314 .keywords: PC, create, context
315 
316 .seealso: PCSetUp(), PCApply(), PCDestroy()
317 @*/
318 PetscErrorCode  PCCreate(MPI_Comm comm,PC *newpc)
319 {
320   PC             pc;
321   PetscErrorCode ierr;
322 
323   PetscFunctionBegin;
324   PetscValidPointer(newpc,1);
325   *newpc = 0;
326 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
327   ierr = PCInitializePackage(PETSC_NULL);CHKERRQ(ierr);
328 #endif
329 
330   ierr = PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_CLASSID,-1,"PC",comm,PCDestroy,PCView);CHKERRQ(ierr);
331 
332   pc->mat                  = 0;
333   pc->pmat                 = 0;
334   pc->setupcalled          = 0;
335   pc->setfromoptionscalled = 0;
336   pc->data                 = 0;
337   pc->diagonalscale        = PETSC_FALSE;
338   pc->diagonalscaleleft    = 0;
339   pc->diagonalscaleright   = 0;
340 
341   pc->modifysubmatrices   = 0;
342   pc->modifysubmatricesP  = 0;
343   *newpc = pc;
344   PetscFunctionReturn(0);
345 
346 }
347 
348 /* -------------------------------------------------------------------------------*/
349 
350 #undef __FUNCT__
351 #define __FUNCT__ "PCApply"
352 /*@
353    PCApply - Applies the preconditioner to a vector.
354 
355    Collective on PC and Vec
356 
357    Input Parameters:
358 +  pc - the preconditioner context
359 -  x - input vector
360 
361    Output Parameter:
362 .  y - output vector
363 
364    Level: developer
365 
366 .keywords: PC, apply
367 
368 .seealso: PCApplyTranspose(), PCApplyBAorAB()
369 @*/
370 PetscErrorCode  PCApply(PC pc,Vec x,Vec y)
371 {
372   PetscErrorCode ierr;
373 
374   PetscFunctionBegin;
375   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
376   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
377   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
378   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
379   if (pc->setupcalled < 2) {
380     ierr = PCSetUp(pc);CHKERRQ(ierr);
381   }
382   if (!pc->ops->apply) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have apply");
383   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
384   ierr = (*pc->ops->apply)(pc,x,y);CHKERRQ(ierr);
385   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
386   PetscFunctionReturn(0);
387 }
388 
389 #undef __FUNCT__
390 #define __FUNCT__ "PCApplySymmetricLeft"
391 /*@
392    PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
393 
394    Collective on PC and Vec
395 
396    Input Parameters:
397 +  pc - the preconditioner context
398 -  x - input vector
399 
400    Output Parameter:
401 .  y - output vector
402 
403    Notes:
404    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
405 
406    Level: developer
407 
408 .keywords: PC, apply, symmetric, left
409 
410 .seealso: PCApply(), PCApplySymmetricRight()
411 @*/
412 PetscErrorCode  PCApplySymmetricLeft(PC pc,Vec x,Vec y)
413 {
414   PetscErrorCode ierr;
415 
416   PetscFunctionBegin;
417   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
418   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
419   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
420   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
421   if (pc->setupcalled < 2) {
422     ierr = PCSetUp(pc);CHKERRQ(ierr);
423   }
424   if (!pc->ops->applysymmetricleft) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have left symmetric apply");
425   ierr = PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
426   ierr = (*pc->ops->applysymmetricleft)(pc,x,y);CHKERRQ(ierr);
427   ierr = PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);CHKERRQ(ierr);
428   PetscFunctionReturn(0);
429 }
430 
431 #undef __FUNCT__
432 #define __FUNCT__ "PCApplySymmetricRight"
433 /*@
434    PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
435 
436    Collective on PC and Vec
437 
438    Input Parameters:
439 +  pc - the preconditioner context
440 -  x - input vector
441 
442    Output Parameter:
443 .  y - output vector
444 
445    Level: developer
446 
447    Notes:
448    Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
449 
450 .keywords: PC, apply, symmetric, right
451 
452 .seealso: PCApply(), PCApplySymmetricLeft()
453 @*/
454 PetscErrorCode  PCApplySymmetricRight(PC pc,Vec x,Vec y)
455 {
456   PetscErrorCode ierr;
457 
458   PetscFunctionBegin;
459   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
460   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
461   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
462   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
463   if (pc->setupcalled < 2) {
464     ierr = PCSetUp(pc);CHKERRQ(ierr);
465   }
466   if (!pc->ops->applysymmetricright) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have left symmetric apply");
467   ierr = PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
468   ierr = (*pc->ops->applysymmetricright)(pc,x,y);CHKERRQ(ierr);
469   ierr = PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);CHKERRQ(ierr);
470   PetscFunctionReturn(0);
471 }
472 
473 #undef __FUNCT__
474 #define __FUNCT__ "PCApplyTranspose"
475 /*@
476    PCApplyTranspose - Applies the transpose of preconditioner to a vector.
477 
478    Collective on PC and Vec
479 
480    Input Parameters:
481 +  pc - the preconditioner context
482 -  x - input vector
483 
484    Output Parameter:
485 .  y - output vector
486 
487    Level: developer
488 
489 .keywords: PC, apply, transpose
490 
491 .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCApplyTransposeExists()
492 @*/
493 PetscErrorCode  PCApplyTranspose(PC pc,Vec x,Vec y)
494 {
495   PetscErrorCode ierr;
496 
497   PetscFunctionBegin;
498   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
499   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
500   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
501   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
502   if (pc->setupcalled < 2) {
503     ierr = PCSetUp(pc);CHKERRQ(ierr);
504   }
505   if (!pc->ops->applytranspose) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have apply transpose");
506   ierr = PetscLogEventBegin(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
507   ierr = (*pc->ops->applytranspose)(pc,x,y);CHKERRQ(ierr);
508   ierr = PetscLogEventEnd(PC_Apply,pc,x,y,0);CHKERRQ(ierr);
509   PetscFunctionReturn(0);
510 }
511 
512 #undef __FUNCT__
513 #define __FUNCT__ "PCApplyTransposeExists"
514 /*@
515    PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
516 
517    Collective on PC and Vec
518 
519    Input Parameters:
520 .  pc - the preconditioner context
521 
522    Output Parameter:
523 .  flg - PETSC_TRUE if a transpose operation is defined
524 
525    Level: developer
526 
527 .keywords: PC, apply, transpose
528 
529 .seealso: PCApplyTranspose()
530 @*/
531 PetscErrorCode  PCApplyTransposeExists(PC pc,PetscBool  *flg)
532 {
533   PetscFunctionBegin;
534   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
535   PetscValidPointer(flg,2);
536   if (pc->ops->applytranspose) *flg = PETSC_TRUE;
537   else                         *flg = PETSC_FALSE;
538   PetscFunctionReturn(0);
539 }
540 
541 #undef __FUNCT__
542 #define __FUNCT__ "PCApplyBAorAB"
543 /*@
544    PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x.
545 
546    Collective on PC and Vec
547 
548    Input Parameters:
549 +  pc - the preconditioner context
550 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
551 .  x - input vector
552 -  work - work vector
553 
554    Output Parameter:
555 .  y - output vector
556 
557    Level: developer
558 
559 .keywords: PC, apply, operator
560 
561 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
562 @*/
563 PetscErrorCode  PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
564 {
565   PetscErrorCode ierr;
566 
567   PetscFunctionBegin;
568   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
569   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
570   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
571   PetscValidHeaderSpecific(work,VEC_CLASSID,5);
572   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
573   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");
574   if (pc->diagonalscale && side == PC_SYMMETRIC) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");
575 
576   if (pc->setupcalled < 2) {
577     ierr = PCSetUp(pc);CHKERRQ(ierr);
578   }
579 
580   if (pc->diagonalscale) {
581     if (pc->ops->applyBA) {
582       Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
583       ierr = VecDuplicate(x,&work2);CHKERRQ(ierr);
584       ierr = PCDiagonalScaleRight(pc,x,work2);CHKERRQ(ierr);
585       ierr = (*pc->ops->applyBA)(pc,side,work2,y,work);CHKERRQ(ierr);
586       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
587       ierr = VecDestroy(work2);CHKERRQ(ierr);
588     } else if (side == PC_RIGHT) {
589       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
590       ierr = PCApply(pc,y,work);CHKERRQ(ierr);
591       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
592       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
593     } else if (side == PC_LEFT) {
594       ierr = PCDiagonalScaleRight(pc,x,y);CHKERRQ(ierr);
595       ierr = MatMult(pc->mat,y,work);CHKERRQ(ierr);
596       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
597       ierr = PCDiagonalScaleLeft(pc,y,y);CHKERRQ(ierr);
598     } else if (side == PC_SYMMETRIC) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
599   } else {
600     if (pc->ops->applyBA) {
601       ierr = (*pc->ops->applyBA)(pc,side,x,y,work);CHKERRQ(ierr);
602     } else if (side == PC_RIGHT) {
603       ierr = PCApply(pc,x,work);CHKERRQ(ierr);
604       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
605     } else if (side == PC_LEFT) {
606       ierr = MatMult(pc->mat,x,work);CHKERRQ(ierr);
607       ierr = PCApply(pc,work,y);CHKERRQ(ierr);
608     } else if (side == PC_SYMMETRIC) {
609       /* There's an extra copy here; maybe should provide 2 work vectors instead? */
610       ierr = PCApplySymmetricRight(pc,x,work);CHKERRQ(ierr);
611       ierr = MatMult(pc->mat,work,y);CHKERRQ(ierr);
612       ierr = VecCopy(y,work);CHKERRQ(ierr);
613       ierr = PCApplySymmetricLeft(pc,work,y);CHKERRQ(ierr);
614     }
615   }
616   PetscFunctionReturn(0);
617 }
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "PCApplyBAorABTranspose"
621 /*@
622    PCApplyBAorABTranspose - Applies the transpose of the preconditioner
623    and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
624    NOT tr(B*A) = tr(A)*tr(B).
625 
626    Collective on PC and Vec
627 
628    Input Parameters:
629 +  pc - the preconditioner context
630 .  side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
631 .  x - input vector
632 -  work - work vector
633 
634    Output Parameter:
635 .  y - output vector
636 
637 
638    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
639       defined by B'. This is why this has the funny form that it computes tr(B) * tr(A)
640 
641     Level: developer
642 
643 .keywords: PC, apply, operator, transpose
644 
645 .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
646 @*/
647 PetscErrorCode  PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
648 {
649   PetscErrorCode ierr;
650 
651   PetscFunctionBegin;
652   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
653   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
654   PetscValidHeaderSpecific(y,VEC_CLASSID,4);
655   PetscValidHeaderSpecific(work,VEC_CLASSID,5);
656   if (x == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"x and y must be different vectors");
657   if (pc->ops->applyBAtranspose) {
658     ierr = (*pc->ops->applyBAtranspose)(pc,side,x,y,work);CHKERRQ(ierr);
659     PetscFunctionReturn(0);
660   }
661   if (side != PC_LEFT && side != PC_RIGHT) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
662 
663   if (pc->setupcalled < 2) {
664     ierr = PCSetUp(pc);CHKERRQ(ierr);
665   }
666 
667   if (side == PC_RIGHT) {
668     ierr = PCApplyTranspose(pc,x,work);CHKERRQ(ierr);
669     ierr = MatMultTranspose(pc->mat,work,y);CHKERRQ(ierr);
670   } else if (side == PC_LEFT) {
671     ierr = MatMultTranspose(pc->mat,x,work);CHKERRQ(ierr);
672     ierr = PCApplyTranspose(pc,work,y);CHKERRQ(ierr);
673   }
674   /* add support for PC_SYMMETRIC */
675   PetscFunctionReturn(0); /* actually will never get here */
676 }
677 
678 /* -------------------------------------------------------------------------------*/
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "PCApplyRichardsonExists"
682 /*@
683    PCApplyRichardsonExists - Determines whether a particular preconditioner has a
684    built-in fast application of Richardson's method.
685 
686    Not Collective
687 
688    Input Parameter:
689 .  pc - the preconditioner
690 
691    Output Parameter:
692 .  exists - PETSC_TRUE or PETSC_FALSE
693 
694    Level: developer
695 
696 .keywords: PC, apply, Richardson, exists
697 
698 .seealso: PCApplyRichardson()
699 @*/
700 PetscErrorCode  PCApplyRichardsonExists(PC pc,PetscBool  *exists)
701 {
702   PetscFunctionBegin;
703   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
704   PetscValidIntPointer(exists,2);
705   if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
706   else                          *exists = PETSC_FALSE;
707   PetscFunctionReturn(0);
708 }
709 
710 #undef __FUNCT__
711 #define __FUNCT__ "PCApplyRichardson"
712 /*@
713    PCApplyRichardson - Applies several steps of Richardson iteration with
714    the particular preconditioner. This routine is usually used by the
715    Krylov solvers and not the application code directly.
716 
717    Collective on PC
718 
719    Input Parameters:
720 +  pc  - the preconditioner context
721 .  b   - the right hand side
722 .  w   - one work vector
723 .  rtol - relative decrease in residual norm convergence criteria
724 .  abstol - absolute residual norm convergence criteria
725 .  dtol - divergence residual norm increase criteria
726 .  its - the number of iterations to apply.
727 -  guesszero - if the input x contains nonzero initial guess
728 
729    Output Parameter:
730 +  outits - number of iterations actually used (for SOR this always equals its)
731 .  reason - the reason the apply terminated
732 -  y - the solution (also contains initial guess if guesszero is PETSC_FALSE
733 
734    Notes:
735    Most preconditioners do not support this function. Use the command
736    PCApplyRichardsonExists() to determine if one does.
737 
738    Except for the multigrid PC this routine ignores the convergence tolerances
739    and always runs for the number of iterations
740 
741    Level: developer
742 
743 .keywords: PC, apply, Richardson
744 
745 .seealso: PCApplyRichardsonExists()
746 @*/
747 PetscErrorCode  PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool  guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
748 {
749   PetscErrorCode ierr;
750 
751   PetscFunctionBegin;
752   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
753   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
754   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
755   PetscValidHeaderSpecific(w,VEC_CLASSID,4);
756   if (b == y) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_IDN,"b and y must be different vectors");
757   if (pc->setupcalled < 2) {
758     ierr = PCSetUp(pc);CHKERRQ(ierr);
759   }
760   if (!pc->ops->applyrichardson) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"PC does not have apply richardson");
761   ierr = (*pc->ops->applyrichardson)(pc,b,y,w,rtol,abstol,dtol,its,guesszero,outits,reason);CHKERRQ(ierr);
762   PetscFunctionReturn(0);
763 }
764 
765 /*
766       a setupcall of 0 indicates never setup,
767                      1 needs to be resetup,
768                      2 does not need any changes.
769 */
770 #undef __FUNCT__
771 #define __FUNCT__ "PCSetUp"
772 /*@
773    PCSetUp - Prepares for the use of a preconditioner.
774 
775    Collective on PC
776 
777    Input Parameter:
778 .  pc - the preconditioner context
779 
780    Level: developer
781 
782 .keywords: PC, setup
783 
784 .seealso: PCCreate(), PCApply(), PCDestroy()
785 @*/
786 PetscErrorCode  PCSetUp(PC pc)
787 {
788   PetscErrorCode ierr;
789   const char     *def;
790 
791   PetscFunctionBegin;
792   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
793   if (!pc->mat) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");
794 
795   if (pc->setupcalled > 1) {
796     ierr = PetscInfo(pc,"Setting PC with identical preconditioner\n");CHKERRQ(ierr);
797     PetscFunctionReturn(0);
798   } else if (!pc->setupcalled) {
799     ierr = PetscInfo(pc,"Setting up new PC\n");CHKERRQ(ierr);
800   } else if (pc->flag == SAME_NONZERO_PATTERN) {
801     ierr = PetscInfo(pc,"Setting up PC with same nonzero pattern\n");CHKERRQ(ierr);
802   } else {
803     ierr = PetscInfo(pc,"Setting up PC with different nonzero pattern\n");CHKERRQ(ierr);
804   }
805 
806   if (!((PetscObject)pc)->type_name) {
807     ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr);
808     ierr = PCSetType(pc,def);CHKERRQ(ierr);
809   }
810 
811   ierr = PetscLogEventBegin(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
812   if (pc->ops->setup) {
813     ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr);
814   }
815   pc->setupcalled = 2;
816   ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
817   PetscFunctionReturn(0);
818 }
819 
820 #undef __FUNCT__
821 #define __FUNCT__ "PCSetUpOnBlocks"
822 /*@
823    PCSetUpOnBlocks - Sets up the preconditioner for each block in
824    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
825    methods.
826 
827    Collective on PC
828 
829    Input Parameters:
830 .  pc - the preconditioner context
831 
832    Level: developer
833 
834 .keywords: PC, setup, blocks
835 
836 .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
837 @*/
838 PetscErrorCode  PCSetUpOnBlocks(PC pc)
839 {
840   PetscErrorCode ierr;
841 
842   PetscFunctionBegin;
843   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
844   if (!pc->ops->setuponblocks) PetscFunctionReturn(0);
845   ierr = PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
846   ierr = (*pc->ops->setuponblocks)(pc);CHKERRQ(ierr);
847   ierr = PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
848   PetscFunctionReturn(0);
849 }
850 
851 #undef __FUNCT__
852 #define __FUNCT__ "PCSetModifySubMatrices"
853 /*@C
854    PCSetModifySubMatrices - Sets a user-defined routine for modifying the
855    submatrices that arise within certain subdomain-based preconditioners.
856    The basic submatrices are extracted from the preconditioner matrix as
857    usual; the user can then alter these (for example, to set different boundary
858    conditions for each submatrix) before they are used for the local solves.
859 
860    Logically Collective on PC
861 
862    Input Parameters:
863 +  pc - the preconditioner context
864 .  func - routine for modifying the submatrices
865 -  ctx - optional user-defined context (may be null)
866 
867    Calling sequence of func:
868 $     func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
869 
870 .  row - an array of index sets that contain the global row numbers
871          that comprise each local submatrix
872 .  col - an array of index sets that contain the global column numbers
873          that comprise each local submatrix
874 .  submat - array of local submatrices
875 -  ctx - optional user-defined context for private data for the
876          user-defined func routine (may be null)
877 
878    Notes:
879    PCSetModifySubMatrices() MUST be called before KSPSetUp() and
880    KSPSolve().
881 
882    A routine set by PCSetModifySubMatrices() is currently called within
883    the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
884    preconditioners.  All other preconditioners ignore this routine.
885 
886    Level: advanced
887 
888 .keywords: PC, set, modify, submatrices
889 
890 .seealso: PCModifySubMatrices()
891 @*/
892 PetscErrorCode  PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
893 {
894   PetscFunctionBegin;
895   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
896   pc->modifysubmatrices  = func;
897   pc->modifysubmatricesP = ctx;
898   PetscFunctionReturn(0);
899 }
900 
901 #undef __FUNCT__
902 #define __FUNCT__ "PCModifySubMatrices"
903 /*@C
904    PCModifySubMatrices - Calls an optional user-defined routine within
905    certain preconditioners if one has been set with PCSetModifySubMarices().
906 
907    Collective on PC
908 
909    Input Parameters:
910 +  pc - the preconditioner context
911 .  nsub - the number of local submatrices
912 .  row - an array of index sets that contain the global row numbers
913          that comprise each local submatrix
914 .  col - an array of index sets that contain the global column numbers
915          that comprise each local submatrix
916 .  submat - array of local submatrices
917 -  ctx - optional user-defined context for private data for the
918          user-defined routine (may be null)
919 
920    Output Parameter:
921 .  submat - array of local submatrices (the entries of which may
922             have been modified)
923 
924    Notes:
925    The user should NOT generally call this routine, as it will
926    automatically be called within certain preconditioners (currently
927    block Jacobi, additive Schwarz) if set.
928 
929    The basic submatrices are extracted from the preconditioner matrix
930    as usual; the user can then alter these (for example, to set different
931    boundary conditions for each submatrix) before they are used for the
932    local solves.
933 
934    Level: developer
935 
936 .keywords: PC, modify, submatrices
937 
938 .seealso: PCSetModifySubMatrices()
939 @*/
940 PetscErrorCode  PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
941 {
942   PetscErrorCode ierr;
943 
944   PetscFunctionBegin;
945   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
946   if (!pc->modifysubmatrices) PetscFunctionReturn(0);
947   ierr = PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
948   ierr = (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);CHKERRQ(ierr);
949   ierr = PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
950   PetscFunctionReturn(0);
951 }
952 
953 #undef __FUNCT__
954 #define __FUNCT__ "PCSetOperators"
955 /*@
956    PCSetOperators - Sets the matrix associated with the linear system and
957    a (possibly) different one associated with the preconditioner.
958 
959    Logically Collective on PC and Mat
960 
961    Input Parameters:
962 +  pc - the preconditioner context
963 .  Amat - the matrix associated with the linear system
964 .  Pmat - the matrix to be used in constructing the preconditioner, usually
965           the same as Amat.
966 -  flag - flag indicating information about the preconditioner matrix structure
967    during successive linear solves.  This flag is ignored the first time a
968    linear system is solved, and thus is irrelevant when solving just one linear
969    system.
970 
971    Notes:
972    The flag can be used to eliminate unnecessary work in the preconditioner
973    during the repeated solution of linear systems of the same size.  The
974    available options are
975 +    SAME_PRECONDITIONER -
976        Pmat is identical during successive linear solves.
977        This option is intended for folks who are using
978        different Amat and Pmat matrices and wish to reuse the
979        same preconditioner matrix.  For example, this option
980        saves work by not recomputing incomplete factorization
981        for ILU/ICC preconditioners.
982 .     SAME_NONZERO_PATTERN -
983        Pmat has the same nonzero structure during
984        successive linear solves.
985 -     DIFFERENT_NONZERO_PATTERN -
986        Pmat does not have the same nonzero structure.
987 
988     Passing a PETSC_NULL for Amat or Pmat removes the matrix that is currently used.
989 
990     If you wish to replace either Amat or Pmat but leave the other one untouched then
991     first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
992     on it and then pass it back in in your call to KSPSetOperators().
993 
994    Caution:
995    If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
996    and does not check the structure of the matrix.  If you erroneously
997    claim that the structure is the same when it actually is not, the new
998    preconditioner will not function correctly.  Thus, use this optimization
999    feature carefully!
1000 
1001    If in doubt about whether your preconditioner matrix has changed
1002    structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
1003 
1004    More Notes about Repeated Solution of Linear Systems:
1005    PETSc does NOT reset the matrix entries of either Amat or Pmat
1006    to zero after a linear solve; the user is completely responsible for
1007    matrix assembly.  See the routine MatZeroEntries() if desiring to
1008    zero all elements of a matrix.
1009 
1010    Level: intermediate
1011 
1012 .keywords: PC, set, operators, matrix, linear system
1013 
1014 .seealso: PCGetOperators(), MatZeroEntries()
1015  @*/
1016 PetscErrorCode  PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
1017 {
1018   PetscErrorCode ierr;
1019   PetscInt       m1,n1,m2,n2;
1020 
1021   PetscFunctionBegin;
1022   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1023   if (Amat) PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1024   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_CLASSID,3);
1025   if (Amat) PetscCheckSameComm(pc,1,Amat,2);
1026   if (Pmat) PetscCheckSameComm(pc,1,Pmat,3);
1027   if (pc->setupcalled && Amat && Pmat) {
1028     ierr = MatGetLocalSize(Amat,&m1,&n1);CHKERRQ(ierr);
1029     ierr = MatGetLocalSize(pc->mat,&m2,&n2);CHKERRQ(ierr);
1030     if (m1 != m2 || n1 != n2) {
1031       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);
1032     }
1033     ierr = MatGetLocalSize(Pmat,&m1,&n1);CHKERRQ(ierr);
1034     ierr = MatGetLocalSize(pc->pmat,&m2,&n2);CHKERRQ(ierr);
1035     if (m1 != m2 || n1 != n2) {
1036       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);
1037     }
1038   }
1039 
1040   /* reference first in case the matrices are the same */
1041   if (Amat) {ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);}
1042   if (pc->mat) {ierr = MatDestroy(pc->mat);CHKERRQ(ierr);}
1043   if (Pmat) {ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);}
1044   if (pc->pmat) {ierr = MatDestroy(pc->pmat);CHKERRQ(ierr);}
1045   pc->mat  = Amat;
1046   pc->pmat = Pmat;
1047 
1048   if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1049     pc->setupcalled = 1;
1050   }
1051   pc->flag = flag;
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 #undef __FUNCT__
1056 #define __FUNCT__ "PCGetOperators"
1057 /*@C
1058    PCGetOperators - Gets the matrix associated with the linear system and
1059    possibly a different one associated with the preconditioner.
1060 
1061    Not collective, though parallel Mats are returned if the PC is parallel
1062 
1063    Input Parameter:
1064 .  pc - the preconditioner context
1065 
1066    Output Parameters:
1067 +  mat - the matrix associated with the linear system
1068 .  pmat - matrix associated with the preconditioner, usually the same
1069           as mat.
1070 -  flag - flag indicating information about the preconditioner
1071           matrix structure.  See PCSetOperators() for details.
1072 
1073    Level: intermediate
1074 
1075    Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1076       are created in PC and returned to the user. In this case, if both operators
1077       mat and pmat are requested, two DIFFERENT operators will be returned. If
1078       only one is requested both operators in the PC will be the same (i.e. as
1079       if one had called KSP/PCSetOperators() with the same argument for both Mats).
1080       The user must set the sizes of the returned matrices and their type etc just
1081       as if the user created them with MatCreate(). For example,
1082 
1083 $         KSP/PCGetOperators(ksp/pc,&mat,PETSC_NULL,PETSC_NULL); is equivalent to
1084 $           set size, type, etc of mat
1085 
1086 $         MatCreate(comm,&mat);
1087 $         KSP/PCSetOperators(ksp/pc,mat,mat,SAME_NONZERO_PATTERN);
1088 $         PetscObjectDereference((PetscObject)mat);
1089 $           set size, type, etc of mat
1090 
1091      and
1092 
1093 $         KSP/PCGetOperators(ksp/pc,&mat,&pmat,PETSC_NULL); is equivalent to
1094 $           set size, type, etc of mat and pmat
1095 
1096 $         MatCreate(comm,&mat);
1097 $         MatCreate(comm,&pmat);
1098 $         KSP/PCSetOperators(ksp/pc,mat,pmat,SAME_NONZERO_PATTERN);
1099 $         PetscObjectDereference((PetscObject)mat);
1100 $         PetscObjectDereference((PetscObject)pmat);
1101 $           set size, type, etc of mat and pmat
1102 
1103     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1104     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1105     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1106     at this is when you create a SNES you do not NEED to create a KSP and attach it to
1107     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1108     you do not need to attach a PC to it (the KSP object manages the PC object for you).
1109     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1110     it can be created for you?
1111 
1112 
1113 .keywords: PC, get, operators, matrix, linear system
1114 
1115 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1116 @*/
1117 PetscErrorCode  PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1118 {
1119   PetscErrorCode ierr;
1120 
1121   PetscFunctionBegin;
1122   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1123   if (mat) {
1124     if (!pc->mat) {
1125       if (pc->pmat && !pmat) {  /* pmat has been set, but user did not request it, so use for mat */
1126         pc->mat = pc->pmat;
1127         ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1128       } else {                  /* both mat and pmat are empty */
1129         ierr = MatCreate(((PetscObject)pc)->comm,&pc->mat);CHKERRQ(ierr);
1130         if (!pmat) { /* user did NOT request pmat, so make same as mat */
1131           pc->pmat = pc->mat;
1132           ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1133         }
1134       }
1135     }
1136     *mat  = pc->mat;
1137   }
1138   if (pmat) {
1139     if (!pc->pmat) {
1140       if (pc->mat && !mat) {    /* mat has been set but was not requested, so use for pmat */
1141         pc->pmat = pc->mat;
1142         ierr    = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1143       } else {
1144         ierr = MatCreate(((PetscObject)pc)->comm,&pc->pmat);CHKERRQ(ierr);
1145         if (!mat) { /* user did NOT request mat, so make same as pmat */
1146           pc->mat = pc->pmat;
1147           ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1148         }
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