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