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