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