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