xref: /petsc/src/ksp/pc/interface/precon.c (revision 28adb3f739167aae2a3fdf9bf501dbd0150de1de)
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 .  x   - the initial guess
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 
709    Output Parameter:
710 +  outits - number of iterations actually used (for SOR this always equals its)
711 .  reason - the reason the apply terminated
712 -  y - the solution
713 
714    Notes:
715    Most preconditioners do not support this function. Use the command
716    PCApplyRichardsonExists() to determine if one does.
717 
718    Except for the multigrid PC this routine ignores the convergence tolerances
719    and always runs for the number of iterations
720 
721    Level: developer
722 
723 .keywords: PC, apply, Richardson
724 
725 .seealso: PCApplyRichardsonExists()
726 @*/
727 PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardson(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscInt *outits,PCRichardsonConvergedReason *reason)
728 {
729   PetscErrorCode ierr;
730 
731   PetscFunctionBegin;
732   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
733   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
734   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
735   PetscValidHeaderSpecific(w,VEC_COOKIE,4);
736   if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
737   if (pc->setupcalled < 2) {
738     ierr = PCSetUp(pc);CHKERRQ(ierr);
739   }
740   if (!pc->ops->applyrichardson) SETERRQ(PETSC_ERR_SUP,"PC does not have apply richardson");
741   ierr = (*pc->ops->applyrichardson)(pc,x,y,w,rtol,abstol,dtol,its,outits,reason);CHKERRQ(ierr);
742   PetscFunctionReturn(0);
743 }
744 
745 /*
746       a setupcall of 0 indicates never setup,
747                      1 needs to be resetup,
748                      2 does not need any changes.
749 */
750 #undef __FUNCT__
751 #define __FUNCT__ "PCSetUp"
752 /*@
753    PCSetUp - Prepares for the use of a preconditioner.
754 
755    Collective on PC
756 
757    Input Parameter:
758 .  pc - the preconditioner context
759 
760    Level: developer
761 
762 .keywords: PC, setup
763 
764 .seealso: PCCreate(), PCApply(), PCDestroy()
765 @*/
766 PetscErrorCode PETSCKSP_DLLEXPORT PCSetUp(PC pc)
767 {
768   PetscErrorCode ierr;
769   const char     *def;
770 
771   PetscFunctionBegin;
772   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
773 
774   if (!pc->mat) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
775 
776   if (pc->setupcalled > 1) {
777     ierr = PetscInfo(pc,"Setting PC with identical preconditioner\n");CHKERRQ(ierr);
778     PetscFunctionReturn(0);
779   } else if (!pc->setupcalled) {
780     ierr = PetscInfo(pc,"Setting up new PC\n");CHKERRQ(ierr);
781   } else if (pc->flag == SAME_NONZERO_PATTERN) {
782     ierr = PetscInfo(pc,"Setting up PC with same nonzero pattern\n");CHKERRQ(ierr);
783   } else {
784     ierr = PetscInfo(pc,"Setting up PC with different nonzero pattern\n");CHKERRQ(ierr);
785   }
786 
787   if (!((PetscObject)pc)->type_name) {
788     ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr);
789     ierr = PCSetType(pc,def);CHKERRQ(ierr);
790   }
791 
792   ierr = PetscLogEventBegin(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
793   if (pc->ops->setup) {
794     ierr = (*pc->ops->setup)(pc);CHKERRQ(ierr);
795   }
796   pc->setupcalled = 2;
797   ierr = PetscLogEventEnd(PC_SetUp,pc,0,0,0);CHKERRQ(ierr);
798   PetscFunctionReturn(0);
799 }
800 
801 #undef __FUNCT__
802 #define __FUNCT__ "PCSetUpOnBlocks"
803 /*@
804    PCSetUpOnBlocks - Sets up the preconditioner for each block in
805    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
806    methods.
807 
808    Collective on PC
809 
810    Input Parameters:
811 .  pc - the preconditioner context
812 
813    Level: developer
814 
815 .keywords: PC, setup, blocks
816 
817 .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
818 @*/
819 PetscErrorCode PETSCKSP_DLLEXPORT PCSetUpOnBlocks(PC pc)
820 {
821   PetscErrorCode ierr;
822 
823   PetscFunctionBegin;
824   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
825   if (!pc->ops->setuponblocks) PetscFunctionReturn(0);
826   ierr = PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
827   ierr = (*pc->ops->setuponblocks)(pc);CHKERRQ(ierr);
828   ierr = PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);CHKERRQ(ierr);
829   PetscFunctionReturn(0);
830 }
831 
832 #undef __FUNCT__
833 #define __FUNCT__ "PCSetModifySubMatrices"
834 /*@C
835    PCSetModifySubMatrices - Sets a user-defined routine for modifying the
836    submatrices that arise within certain subdomain-based preconditioners.
837    The basic submatrices are extracted from the preconditioner matrix as
838    usual; the user can then alter these (for example, to set different boundary
839    conditions for each submatrix) before they are used for the local solves.
840 
841    Collective on PC
842 
843    Input Parameters:
844 +  pc - the preconditioner context
845 .  func - routine for modifying the submatrices
846 -  ctx - optional user-defined context (may be null)
847 
848    Calling sequence of func:
849 $     func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
850 
851 .  row - an array of index sets that contain the global row numbers
852          that comprise each local submatrix
853 .  col - an array of index sets that contain the global column numbers
854          that comprise each local submatrix
855 .  submat - array of local submatrices
856 -  ctx - optional user-defined context for private data for the
857          user-defined func routine (may be null)
858 
859    Notes:
860    PCSetModifySubMatrices() MUST be called before KSPSetUp() and
861    KSPSolve().
862 
863    A routine set by PCSetModifySubMatrices() is currently called within
864    the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
865    preconditioners.  All other preconditioners ignore this routine.
866 
867    Level: advanced
868 
869 .keywords: PC, set, modify, submatrices
870 
871 .seealso: PCModifySubMatrices()
872 @*/
873 PetscErrorCode PETSCKSP_DLLEXPORT PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
874 {
875   PetscFunctionBegin;
876   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
877   pc->modifysubmatrices  = func;
878   pc->modifysubmatricesP = ctx;
879   PetscFunctionReturn(0);
880 }
881 
882 #undef __FUNCT__
883 #define __FUNCT__ "PCModifySubMatrices"
884 /*@C
885    PCModifySubMatrices - Calls an optional user-defined routine within
886    certain preconditioners if one has been set with PCSetModifySubMarices().
887 
888    Collective on PC
889 
890    Input Parameters:
891 +  pc - the preconditioner context
892 .  nsub - the number of local submatrices
893 .  row - an array of index sets that contain the global row numbers
894          that comprise each local submatrix
895 .  col - an array of index sets that contain the global column numbers
896          that comprise each local submatrix
897 .  submat - array of local submatrices
898 -  ctx - optional user-defined context for private data for the
899          user-defined routine (may be null)
900 
901    Output Parameter:
902 .  submat - array of local submatrices (the entries of which may
903             have been modified)
904 
905    Notes:
906    The user should NOT generally call this routine, as it will
907    automatically be called within certain preconditioners (currently
908    block Jacobi, additive Schwarz) if set.
909 
910    The basic submatrices are extracted from the preconditioner matrix
911    as usual; the user can then alter these (for example, to set different
912    boundary conditions for each submatrix) before they are used for the
913    local solves.
914 
915    Level: developer
916 
917 .keywords: PC, modify, submatrices
918 
919 .seealso: PCSetModifySubMatrices()
920 @*/
921 PetscErrorCode PETSCKSP_DLLEXPORT PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
922 {
923   PetscErrorCode ierr;
924 
925   PetscFunctionBegin;
926   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
927   if (!pc->modifysubmatrices) PetscFunctionReturn(0);
928   ierr = PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
929   ierr = (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);CHKERRQ(ierr);
930   ierr = PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);CHKERRQ(ierr);
931   PetscFunctionReturn(0);
932 }
933 
934 #undef __FUNCT__
935 #define __FUNCT__ "PCSetOperators"
936 /*@
937    PCSetOperators - Sets the matrix associated with the linear system and
938    a (possibly) different one associated with the preconditioner.
939 
940    Collective on PC and Mat
941 
942    Input Parameters:
943 +  pc - the preconditioner context
944 .  Amat - the matrix associated with the linear system
945 .  Pmat - the matrix to be used in constructing the preconditioner, usually
946           the same as Amat.
947 -  flag - flag indicating information about the preconditioner matrix structure
948    during successive linear solves.  This flag is ignored the first time a
949    linear system is solved, and thus is irrelevant when solving just one linear
950    system.
951 
952    Notes:
953    The flag can be used to eliminate unnecessary work in the preconditioner
954    during the repeated solution of linear systems of the same size.  The
955    available options are
956 +    SAME_PRECONDITIONER -
957        Pmat is identical during successive linear solves.
958        This option is intended for folks who are using
959        different Amat and Pmat matrices and wish to reuse the
960        same preconditioner matrix.  For example, this option
961        saves work by not recomputing incomplete factorization
962        for ILU/ICC preconditioners.
963 .     SAME_NONZERO_PATTERN -
964        Pmat has the same nonzero structure during
965        successive linear solves.
966 -     DIFFERENT_NONZERO_PATTERN -
967        Pmat does not have the same nonzero structure.
968 
969     Passing a PETSC_NULL for Amat or Pmat removes the matrix that is currently used.
970 
971     If you wish to replace either Amat or Pmat but leave the other one untouched then
972     first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
973     on it and then pass it back in in your call to KSPSetOperators().
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 
1001   PetscFunctionBegin;
1002   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1003   if (Amat) PetscValidHeaderSpecific(Amat,MAT_COOKIE,2);
1004   if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_COOKIE,3);
1005   if (Amat) PetscCheckSameComm(pc,1,Amat,2);
1006   if (Pmat) PetscCheckSameComm(pc,1,Pmat,3);
1007 
1008   /* reference first in case the matrices are the same */
1009   if (Amat) {ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);}
1010   if (pc->mat) {ierr = MatDestroy(pc->mat);CHKERRQ(ierr);}
1011   if (Pmat) {ierr = PetscObjectReference((PetscObject)Pmat);CHKERRQ(ierr);}
1012   if (pc->pmat) {ierr = MatDestroy(pc->pmat);CHKERRQ(ierr);}
1013   pc->mat  = Amat;
1014   pc->pmat = Pmat;
1015 
1016   if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1017     pc->setupcalled = 1;
1018   }
1019   pc->flag = flag;
1020   PetscFunctionReturn(0);
1021 }
1022 
1023 #undef __FUNCT__
1024 #define __FUNCT__ "PCGetOperators"
1025 /*@C
1026    PCGetOperators - Gets the matrix associated with the linear system and
1027    possibly a different one associated with the preconditioner.
1028 
1029    Not collective, though parallel Mats are returned if the PC is parallel
1030 
1031    Input Parameter:
1032 .  pc - the preconditioner context
1033 
1034    Output Parameters:
1035 +  mat - the matrix associated with the linear system
1036 .  pmat - matrix associated with the preconditioner, usually the same
1037           as mat.
1038 -  flag - flag indicating information about the preconditioner
1039           matrix structure.  See PCSetOperators() for details.
1040 
1041    Level: intermediate
1042 
1043    Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1044       are created in PC and returned to the user. In this case, if both operators
1045       mat and pmat are requested, two DIFFERENT operators will be returned. If
1046       only one is requested both operators in the PC will be the same (i.e. as
1047       if one had called KSP/PCSetOperators() with the same argument for both Mats).
1048       The user must set the sizes of the returned matrices and their type etc just
1049       as if the user created them with MatCreate(). For example,
1050 
1051 $         KSP/PCGetOperators(ksp/pc,&mat,PETSC_NULL,PETSC_NULL); is equivalent to
1052 $           set size, type, etc of mat
1053 
1054 $         MatCreate(comm,&mat);
1055 $         KSP/PCSetOperators(ksp/pc,mat,mat,SAME_NONZERO_PATTERN);
1056 $         PetscObjectDereference((PetscObject)mat);
1057 $           set size, type, etc of mat
1058 
1059      and
1060 
1061 $         KSP/PCGetOperators(ksp/pc,&mat,&pmat,PETSC_NULL); is equivalent to
1062 $           set size, type, etc of mat and pmat
1063 
1064 $         MatCreate(comm,&mat);
1065 $         MatCreate(comm,&pmat);
1066 $         KSP/PCSetOperators(ksp/pc,mat,pmat,SAME_NONZERO_PATTERN);
1067 $         PetscObjectDereference((PetscObject)mat);
1068 $         PetscObjectDereference((PetscObject)pmat);
1069 $           set size, type, etc of mat and pmat
1070 
1071     The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1072     of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1073     managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1074     at this is when you create a SNES you do not NEED to create a KSP and attach it to
1075     the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1076     you do not need to attach a PC to it (the KSP object manages the PC object for you).
1077     Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1078     it can be created for you?
1079 
1080 
1081 .keywords: PC, get, operators, matrix, linear system
1082 
1083 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1084 @*/
1085 PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1086 {
1087   PetscErrorCode ierr;
1088 
1089   PetscFunctionBegin;
1090   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1091   if (mat) {
1092     if (!pc->mat) {
1093       ierr = MatCreate(((PetscObject)pc)->comm,&pc->mat);CHKERRQ(ierr);
1094       if (!pc->pmat && !pmat) { /* user did NOT request pmat, so make same as mat */
1095         pc->pmat = pc->mat;
1096         ierr     = PetscObjectReference((PetscObject)pc->pmat);CHKERRQ(ierr);
1097       }
1098     }
1099     *mat  = pc->mat;
1100   }
1101   if (pmat) {
1102     if (!pc->pmat) {
1103       ierr = MatCreate(((PetscObject)pc)->comm,&pc->mat);CHKERRQ(ierr);
1104       if (!pc->mat && !mat) { /* user did NOT request mat, so make same as pmat */
1105         pc->mat = pc->pmat;
1106         ierr    = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr);
1107       }
1108     }
1109     *pmat = pc->pmat;
1110   }
1111   if (flag) *flag = pc->flag;
1112   PetscFunctionReturn(0);
1113 }
1114 
1115 #undef __FUNCT__
1116 #define __FUNCT__ "PCGetOperatorsSet"
1117 /*@C
1118    PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1119    possibly a different one associated with the preconditioner have been set in the PC.
1120 
1121    Not collective, though the results on all processes should be the same
1122 
1123    Input Parameter:
1124 .  pc - the preconditioner context
1125 
1126    Output Parameters:
1127 +  mat - the matrix associated with the linear system was set
1128 -  pmat - matrix associated with the preconditioner was set, usually the same
1129 
1130    Level: intermediate
1131 
1132 .keywords: PC, get, operators, matrix, linear system
1133 
1134 .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1135 @*/
1136 PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperatorsSet(PC pc,PetscTruth *mat,PetscTruth *pmat)
1137 {
1138   PetscFunctionBegin;
1139   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1140   if (mat)  *mat  = (pc->mat)  ? PETSC_TRUE : PETSC_FALSE;
1141   if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1142   PetscFunctionReturn(0);
1143 }
1144 
1145 #undef __FUNCT__
1146 #define __FUNCT__ "PCFactorGetMatrix"
1147 /*@
1148    PCFactorGetMatrix - Gets the factored matrix from the
1149    preconditioner context.  This routine is valid only for the LU,
1150    incomplete LU, Cholesky, and incomplete Cholesky methods.
1151 
1152    Not Collective on PC though Mat is parallel if PC is parallel
1153 
1154    Input Parameters:
1155 .  pc - the preconditioner context
1156 
1157    Output parameters:
1158 .  mat - the factored matrix
1159 
1160    Level: advanced
1161 
1162    Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1163 
1164 .keywords: PC, get, factored, matrix
1165 @*/
1166 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix(PC pc,Mat *mat)
1167 {
1168   PetscErrorCode ierr;
1169 
1170   PetscFunctionBegin;
1171   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1172   PetscValidPointer(mat,2);
1173   if (pc->ops->getfactoredmatrix) {
1174     ierr = (*pc->ops->getfactoredmatrix)(pc,mat);CHKERRQ(ierr);
1175   }
1176   PetscFunctionReturn(0);
1177 }
1178 
1179 #undef __FUNCT__
1180 #define __FUNCT__ "PCSetOptionsPrefix"
1181 /*@C
1182    PCSetOptionsPrefix - Sets the prefix used for searching for all
1183    PC options in the database.
1184 
1185    Collective on PC
1186 
1187    Input Parameters:
1188 +  pc - the preconditioner context
1189 -  prefix - the prefix string to prepend to all PC option requests
1190 
1191    Notes:
1192    A hyphen (-) must NOT be given at the beginning of the prefix name.
1193    The first character of all runtime options is AUTOMATICALLY the
1194    hyphen.
1195 
1196    Level: advanced
1197 
1198 .keywords: PC, set, options, prefix, database
1199 
1200 .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1201 @*/
1202 PetscErrorCode PETSCKSP_DLLEXPORT PCSetOptionsPrefix(PC pc,const char prefix[])
1203 {
1204   PetscErrorCode ierr;
1205 
1206   PetscFunctionBegin;
1207   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1208   ierr = PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1209   PetscFunctionReturn(0);
1210 }
1211 
1212 #undef __FUNCT__
1213 #define __FUNCT__ "PCAppendOptionsPrefix"
1214 /*@C
1215    PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1216    PC options in the database.
1217 
1218    Collective on PC
1219 
1220    Input Parameters:
1221 +  pc - the preconditioner context
1222 -  prefix - the prefix string to prepend to all PC option requests
1223 
1224    Notes:
1225    A hyphen (-) must NOT be given at the beginning of the prefix name.
1226    The first character of all runtime options is AUTOMATICALLY the
1227    hyphen.
1228 
1229    Level: advanced
1230 
1231 .keywords: PC, append, options, prefix, database
1232 
1233 .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1234 @*/
1235 PetscErrorCode PETSCKSP_DLLEXPORT PCAppendOptionsPrefix(PC pc,const char prefix[])
1236 {
1237   PetscErrorCode ierr;
1238 
1239   PetscFunctionBegin;
1240   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1241   ierr = PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1242   PetscFunctionReturn(0);
1243 }
1244 
1245 #undef __FUNCT__
1246 #define __FUNCT__ "PCGetOptionsPrefix"
1247 /*@C
1248    PCGetOptionsPrefix - Gets the prefix used for searching for all
1249    PC options in the database.
1250 
1251    Not Collective
1252 
1253    Input Parameters:
1254 .  pc - the preconditioner context
1255 
1256    Output Parameters:
1257 .  prefix - pointer to the prefix string used, is returned
1258 
1259    Notes: On the fortran side, the user should pass in a string 'prifix' of
1260    sufficient length to hold the prefix.
1261 
1262    Level: advanced
1263 
1264 .keywords: PC, get, options, prefix, database
1265 
1266 .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1267 @*/
1268 PetscErrorCode PETSCKSP_DLLEXPORT PCGetOptionsPrefix(PC pc,const char *prefix[])
1269 {
1270   PetscErrorCode ierr;
1271 
1272   PetscFunctionBegin;
1273   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1274   PetscValidPointer(prefix,2);
1275   ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);CHKERRQ(ierr);
1276   PetscFunctionReturn(0);
1277 }
1278 
1279 #undef __FUNCT__
1280 #define __FUNCT__ "PCPreSolve"
1281 /*@
1282    PCPreSolve - Optional pre-solve phase, intended for any
1283    preconditioner-specific actions that must be performed before
1284    the iterative solve itself.
1285 
1286    Collective on PC
1287 
1288    Input Parameters:
1289 +  pc - the preconditioner context
1290 -  ksp - the Krylov subspace context
1291 
1292    Level: developer
1293 
1294    Sample of Usage:
1295 .vb
1296     PCPreSolve(pc,ksp);
1297     KSPSolve(ksp,b,x);
1298     PCPostSolve(pc,ksp);
1299 .ve
1300 
1301    Notes:
1302    The pre-solve phase is distinct from the PCSetUp() phase.
1303 
1304    KSPSolve() calls this directly, so is rarely called by the user.
1305 
1306 .keywords: PC, pre-solve
1307 
1308 .seealso: PCPostSolve()
1309 @*/
1310 PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC pc,KSP ksp)
1311 {
1312   PetscErrorCode ierr;
1313   Vec            x,rhs;
1314   Mat            A,B;
1315 
1316   PetscFunctionBegin;
1317   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1318   PetscValidHeaderSpecific(ksp,KSP_COOKIE,2);
1319   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1320   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1321   /*
1322       Scale the system and have the matrices use the scaled form
1323     only if the two matrices are actually the same (and hence
1324     have the same scaling
1325   */
1326   ierr = PCGetOperators(pc,&A,&B,PETSC_NULL);CHKERRQ(ierr);
1327   if (A == B) {
1328     ierr = MatScaleSystem(pc->mat,rhs,x);CHKERRQ(ierr);
1329     ierr = MatUseScaledForm(pc->mat,PETSC_TRUE);CHKERRQ(ierr);
1330   }
1331 
1332   if (pc->ops->presolve) {
1333     ierr = (*pc->ops->presolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1334   }
1335   PetscFunctionReturn(0);
1336 }
1337 
1338 #undef __FUNCT__
1339 #define __FUNCT__ "PCPostSolve"
1340 /*@
1341    PCPostSolve - Optional post-solve phase, intended for any
1342    preconditioner-specific actions that must be performed after
1343    the iterative solve itself.
1344 
1345    Collective on PC
1346 
1347    Input Parameters:
1348 +  pc - the preconditioner context
1349 -  ksp - the Krylov subspace context
1350 
1351    Sample of Usage:
1352 .vb
1353     PCPreSolve(pc,ksp);
1354     KSPSolve(ksp,b,x);
1355     PCPostSolve(pc,ksp);
1356 .ve
1357 
1358    Note:
1359    KSPSolve() calls this routine directly, so it is rarely called by the user.
1360 
1361    Level: developer
1362 
1363 .keywords: PC, post-solve
1364 
1365 .seealso: PCPreSolve(), KSPSolve()
1366 @*/
1367 PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC pc,KSP ksp)
1368 {
1369   PetscErrorCode ierr;
1370   Vec            x,rhs;
1371   Mat            A,B;
1372 
1373   PetscFunctionBegin;
1374   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1375   PetscValidHeaderSpecific(ksp,KSP_COOKIE,2);
1376   ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
1377   ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
1378   if (pc->ops->postsolve) {
1379     ierr =  (*pc->ops->postsolve)(pc,ksp,rhs,x);CHKERRQ(ierr);
1380   }
1381   /*
1382       Scale the system and have the matrices use the scaled form
1383     only if the two matrices are actually the same (and hence
1384     have the same scaling
1385   */
1386   ierr = PCGetOperators(pc,&A,&B,PETSC_NULL);CHKERRQ(ierr);
1387   if (A == B) {
1388     ierr = MatUnScaleSystem(pc->mat,rhs,x);CHKERRQ(ierr);
1389     ierr = MatUseScaledForm(pc->mat,PETSC_FALSE);CHKERRQ(ierr);
1390   }
1391   PetscFunctionReturn(0);
1392 }
1393 
1394 #undef __FUNCT__
1395 #define __FUNCT__ "PCView"
1396 /*@C
1397    PCView - Prints the PC data structure.
1398 
1399    Collective on PC
1400 
1401    Input Parameters:
1402 +  PC - the PC context
1403 -  viewer - optional visualization context
1404 
1405    Note:
1406    The available visualization contexts include
1407 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1408 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1409          output where only the first processor opens
1410          the file.  All other processors send their
1411          data to the first processor to print.
1412 
1413    The user can open an alternative visualization contexts with
1414    PetscViewerASCIIOpen() (output to a specified file).
1415 
1416    Level: developer
1417 
1418 .keywords: PC, view
1419 
1420 .seealso: KSPView(), PetscViewerASCIIOpen()
1421 @*/
1422 PetscErrorCode PETSCKSP_DLLEXPORT PCView(PC pc,PetscViewer viewer)
1423 {
1424   const PCType      cstr;
1425   PetscErrorCode    ierr;
1426   PetscTruth        mat_exists,iascii,isstring;
1427   PetscViewerFormat format;
1428 
1429   PetscFunctionBegin;
1430   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1431   if (!viewer) {
1432     ierr = PetscViewerASCIIGetStdout(((PetscObject)pc)->comm,&viewer);CHKERRQ(ierr);
1433   }
1434   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
1435   PetscCheckSameComm(pc,1,viewer,2);
1436 
1437   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1438   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
1439   if (iascii) {
1440     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1441     if (((PetscObject)pc)->prefix) {
1442       ierr = PetscViewerASCIIPrintf(viewer,"PC Object:(%s)\n",((PetscObject)pc)->prefix);CHKERRQ(ierr);
1443     } else {
1444       ierr = PetscViewerASCIIPrintf(viewer,"PC Object:\n");CHKERRQ(ierr);
1445     }
1446     ierr = PCGetType(pc,&cstr);CHKERRQ(ierr);
1447     if (cstr) {
1448       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",cstr);CHKERRQ(ierr);
1449     } else {
1450       ierr = PetscViewerASCIIPrintf(viewer,"  type: not yet set\n");CHKERRQ(ierr);
1451     }
1452     if (pc->ops->view) {
1453       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1454       ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);
1455       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1456     }
1457     ierr = PetscObjectExists((PetscObject)pc->mat,&mat_exists);CHKERRQ(ierr);
1458     if (mat_exists) {
1459       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
1460       if (pc->pmat == pc->mat) {
1461         ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix = precond matrix:\n");CHKERRQ(ierr);
1462         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1463         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1464         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1465       } else {
1466         ierr = PetscObjectExists((PetscObject)pc->pmat,&mat_exists);CHKERRQ(ierr);
1467         if (mat_exists) {
1468           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix followed by preconditioner matrix:\n");CHKERRQ(ierr);
1469         } else {
1470           ierr = PetscViewerASCIIPrintf(viewer,"  linear system matrix:\n");CHKERRQ(ierr);
1471         }
1472         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1473         ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1474         if (mat_exists) {ierr = MatView(pc->pmat,viewer);CHKERRQ(ierr);}
1475         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1476       }
1477       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1478     }
1479   } else if (isstring) {
1480     ierr = PCGetType(pc,&cstr);CHKERRQ(ierr);
1481     ierr = PetscViewerStringSPrintf(viewer," %-7.7s",cstr);CHKERRQ(ierr);
1482     if (pc->ops->view) {ierr = (*pc->ops->view)(pc,viewer);CHKERRQ(ierr);}
1483   } else {
1484     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1485   }
1486   PetscFunctionReturn(0);
1487 }
1488 
1489 
1490 #undef __FUNCT__
1491 #define __FUNCT__ "PCSetInitialGuessNonzero"
1492 /*@
1493    PCSetInitialGuessNonzero - Tells the iterative solver that the
1494    initial guess is nonzero; otherwise PC assumes the initial guess
1495    is to be zero (and thus zeros it out before solving).
1496 
1497    Collective on PC
1498 
1499    Input Parameters:
1500 +  pc - iterative context obtained from PCCreate()
1501 -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1502 
1503    Level: Developer
1504 
1505    Notes:
1506     This is a weird function. Since PC's are linear operators on the right hand side they
1507     CANNOT use an initial guess. This function is for the "pass-through" preconditioners
1508     PCKSP, PCREDUNDANT and PCOPENMP and causes the inner KSP object to use the nonzero
1509     initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP.
1510 
1511 
1512 .keywords: PC, set, initial guess, nonzero
1513 
1514 .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
1515 @*/
1516 PetscErrorCode PETSCKSP_DLLEXPORT PCSetInitialGuessNonzero(PC pc,PetscTruth flg)
1517 {
1518   PetscFunctionBegin;
1519   pc->nonzero_guess   = flg;
1520   PetscFunctionReturn(0);
1521 }
1522 
1523 #undef __FUNCT__
1524 #define __FUNCT__ "PCRegister"
1525 /*@C
1526   PCRegister - See PCRegisterDynamic()
1527 
1528   Level: advanced
1529 @*/
1530 PetscErrorCode PETSCKSP_DLLEXPORT PCRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PC))
1531 {
1532   PetscErrorCode ierr;
1533   char           fullname[PETSC_MAX_PATH_LEN];
1534 
1535   PetscFunctionBegin;
1536 
1537   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
1538   ierr = PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
1539   PetscFunctionReturn(0);
1540 }
1541 
1542 #undef __FUNCT__
1543 #define __FUNCT__ "PCComputeExplicitOperator"
1544 /*@
1545     PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1546 
1547     Collective on PC
1548 
1549     Input Parameter:
1550 .   pc - the preconditioner object
1551 
1552     Output Parameter:
1553 .   mat - the explict preconditioned operator
1554 
1555     Notes:
1556     This computation is done by applying the operators to columns of the
1557     identity matrix.
1558 
1559     Currently, this routine uses a dense matrix format when 1 processor
1560     is used and a sparse format otherwise.  This routine is costly in general,
1561     and is recommended for use only with relatively small systems.
1562 
1563     Level: advanced
1564 
1565 .keywords: PC, compute, explicit, operator
1566 
1567 .seealso: KSPComputeExplicitOperator()
1568 
1569 @*/
1570 PetscErrorCode PETSCKSP_DLLEXPORT PCComputeExplicitOperator(PC pc,Mat *mat)
1571 {
1572   Vec            in,out;
1573   PetscErrorCode ierr;
1574   PetscInt       i,M,m,*rows,start,end;
1575   PetscMPIInt    size;
1576   MPI_Comm       comm;
1577   PetscScalar    *array,one = 1.0;
1578 
1579   PetscFunctionBegin;
1580   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1581   PetscValidPointer(mat,2);
1582 
1583   comm = ((PetscObject)pc)->comm;
1584   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1585 
1586   if (!pc->pmat) SETERRQ(PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1587   ierr = MatGetVecs(pc->pmat,&in,0);CHKERRQ(ierr);
1588   ierr = VecDuplicate(in,&out);CHKERRQ(ierr);
1589   ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr);
1590   ierr = VecGetSize(in,&M);CHKERRQ(ierr);
1591   ierr = VecGetLocalSize(in,&m);CHKERRQ(ierr);
1592   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&rows);CHKERRQ(ierr);
1593   for (i=0; i<m; i++) {rows[i] = start + i;}
1594 
1595   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
1596   ierr = MatSetSizes(*mat,m,m,M,M);CHKERRQ(ierr);
1597   if (size == 1) {
1598     ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr);
1599     ierr = MatSeqDenseSetPreallocation(*mat,PETSC_NULL);CHKERRQ(ierr);
1600   } else {
1601     ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr);
1602     ierr = MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1603   }
1604 
1605   for (i=0; i<M; i++) {
1606 
1607     ierr = VecSet(in,0.0);CHKERRQ(ierr);
1608     ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr);
1609     ierr = VecAssemblyBegin(in);CHKERRQ(ierr);
1610     ierr = VecAssemblyEnd(in);CHKERRQ(ierr);
1611 
1612     /* should fix, allowing user to choose side */
1613     ierr = PCApply(pc,in,out);CHKERRQ(ierr);
1614 
1615     ierr = VecGetArray(out,&array);CHKERRQ(ierr);
1616     ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr);
1617     ierr = VecRestoreArray(out,&array);CHKERRQ(ierr);
1618 
1619   }
1620   ierr = PetscFree(rows);CHKERRQ(ierr);
1621   ierr = VecDestroy(out);CHKERRQ(ierr);
1622   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1623   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1624   PetscFunctionReturn(0);
1625 }
1626 
1627