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