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