Actual source code: precon.c
petsc-3.6.1 2015-07-22
2: /*
3: The PC (preconditioner) interface routines, callable by users.
4: */
5: #include <petsc/private/pcimpl.h> /*I "petscksp.h" I*/
6: #include <petscdm.h>
8: /* Logging support */
9: PetscClassId PC_CLASSID;
10: PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
11: PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks, PC_ApplyOnMproc;
15: PetscErrorCode PCGetDefaultType_Private(PC pc,const char *type[])
16: {
18: PetscMPIInt size;
19: PetscBool flg1,flg2,set,flg3;
22: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
23: if (pc->pmat) {
24: PetscErrorCode (*f)(Mat,MatReuse,Mat*);
25: PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",&f);
26: if (size == 1) {
27: MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ICC,&flg1);
28: MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ILU,&flg2);
29: MatIsSymmetricKnown(pc->pmat,&set,&flg3);
30: if (flg1 && (!flg2 || (set && flg3))) {
31: *type = PCICC;
32: } else if (flg2) {
33: *type = PCILU;
34: } else if (f) { /* likely is a parallel matrix run on one processor */
35: *type = PCBJACOBI;
36: } else {
37: *type = PCNONE;
38: }
39: } else {
40: if (f) {
41: *type = PCBJACOBI;
42: } else {
43: *type = PCNONE;
44: }
45: }
46: } else {
47: if (size == 1) {
48: *type = PCILU;
49: } else {
50: *type = PCBJACOBI;
51: }
52: }
53: return(0);
54: }
58: /*@
59: PCReset - Resets a PC context to the pcsetupcalled = 0 state and removes any allocated Vecs and Mats
61: Collective on PC
63: Input Parameter:
64: . pc - the preconditioner context
66: Level: developer
68: 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
70: .keywords: PC, destroy
72: .seealso: PCCreate(), PCSetUp()
73: @*/
74: PetscErrorCode PCReset(PC pc)
75: {
80: if (pc->ops->reset) {
81: (*pc->ops->reset)(pc);
82: }
83: VecDestroy(&pc->diagonalscaleright);
84: VecDestroy(&pc->diagonalscaleleft);
85: MatDestroy(&pc->pmat);
86: MatDestroy(&pc->mat);
88: pc->setupcalled = 0;
89: return(0);
90: }
94: /*@
95: PCDestroy - Destroys PC context that was created with PCCreate().
97: Collective on PC
99: Input Parameter:
100: . pc - the preconditioner context
102: Level: developer
104: .keywords: PC, destroy
106: .seealso: PCCreate(), PCSetUp()
107: @*/
108: PetscErrorCode PCDestroy(PC *pc)
109: {
113: if (!*pc) return(0);
115: if (--((PetscObject)(*pc))->refct > 0) {*pc = 0; return(0);}
117: PCReset(*pc);
119: /* if memory was published with SAWs then destroy it */
120: PetscObjectSAWsViewOff((PetscObject)*pc);
121: if ((*pc)->ops->destroy) {(*(*pc)->ops->destroy)((*pc));}
122: DMDestroy(&(*pc)->dm);
123: PetscHeaderDestroy(pc);
124: return(0);
125: }
129: /*@C
130: PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right
131: scaling as needed by certain time-stepping codes.
133: Logically Collective on PC
135: Input Parameter:
136: . pc - the preconditioner context
138: Output Parameter:
139: . flag - PETSC_TRUE if it applies the scaling
141: Level: developer
143: Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
144: $ D M A D^{-1} y = D M b for left preconditioning or
145: $ D A M D^{-1} z = D b for right preconditioning
147: .keywords: PC
149: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCSetDiagonalScale()
150: @*/
151: PetscErrorCode PCGetDiagonalScale(PC pc,PetscBool *flag)
152: {
156: *flag = pc->diagonalscale;
157: return(0);
158: }
162: /*@
163: PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right
164: scaling as needed by certain time-stepping codes.
166: Logically Collective on PC
168: Input Parameters:
169: + pc - the preconditioner context
170: - s - scaling vector
172: Level: intermediate
174: Notes: The system solved via the Krylov method is
175: $ D M A D^{-1} y = D M b for left preconditioning or
176: $ D A M D^{-1} z = D b for right preconditioning
178: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
180: .keywords: PC
182: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCGetDiagonalScale()
183: @*/
184: PetscErrorCode PCSetDiagonalScale(PC pc,Vec s)
185: {
191: pc->diagonalscale = PETSC_TRUE;
193: PetscObjectReference((PetscObject)s);
194: VecDestroy(&pc->diagonalscaleleft);
196: pc->diagonalscaleleft = s;
198: VecDuplicate(s,&pc->diagonalscaleright);
199: VecCopy(s,pc->diagonalscaleright);
200: VecReciprocal(pc->diagonalscaleright);
201: return(0);
202: }
206: /*@
207: PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes.
209: Logically Collective on PC
211: Input Parameters:
212: + pc - the preconditioner context
213: . in - input vector
214: + out - scaled vector (maybe the same as in)
216: Level: intermediate
218: Notes: The system solved via the Krylov method is
219: $ D M A D^{-1} y = D M b for left preconditioning or
220: $ D A M D^{-1} z = D b for right preconditioning
222: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
224: If diagonal scaling is turned off and in is not out then in is copied to out
226: .keywords: PC
228: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
229: @*/
230: PetscErrorCode PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
231: {
238: if (pc->diagonalscale) {
239: VecPointwiseMult(out,pc->diagonalscaleleft,in);
240: } else if (in != out) {
241: VecCopy(in,out);
242: }
243: return(0);
244: }
248: /*@
249: PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
251: Logically Collective on PC
253: Input Parameters:
254: + pc - the preconditioner context
255: . in - input vector
256: + out - scaled vector (maybe the same as in)
258: Level: intermediate
260: Notes: The system solved via the Krylov method is
261: $ D M A D^{-1} y = D M b for left preconditioning or
262: $ D A M D^{-1} z = D b for right preconditioning
264: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
266: If diagonal scaling is turned off and in is not out then in is copied to out
268: .keywords: PC
270: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
271: @*/
272: PetscErrorCode PCDiagonalScaleRight(PC pc,Vec in,Vec out)
273: {
280: if (pc->diagonalscale) {
281: VecPointwiseMult(out,pc->diagonalscaleright,in);
282: } else if (in != out) {
283: VecCopy(in,out);
284: }
285: return(0);
286: }
290: /*@
291: PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the
292: operator during the preconditioning process it applies the Amat provided to TSSetRHSJacobian(),
293: TSSetIJacobian(), SNESSetJacobian(), KSPSetOperator() or PCSetOperator() not the Pmat.
295: Logically Collective on PC
297: Input Parameters:
298: + pc - the preconditioner context
299: - flg - PETSC_TRUE to use the Amat, PETSC_FALSE to use the Pmat (default is false)
301: Options Database Key:
302: . -pc_use_amat <true,false>
304: Notes:
305: For the common case in which the linear system matrix and the matrix used to construct the
306: preconditioner are identical, this routine is does nothing.
308: Level: intermediate
310: .seealso: PCGetUseAmat(), PCBJACOBI, PGMG, PCFIELDSPLIT, PCCOMPOSITE
311: @*/
312: PetscErrorCode PCSetUseAmat(PC pc,PetscBool flg)
313: {
316: pc->useAmat = flg;
317: return(0);
318: }
322: /*@
323: PCSetErrorIfFailure - Causes PC to generate an error if a FPE, for example a zero pivot, is detected.
325: Logically Collective on PC
327: Input Parameters:
328: + pc - iterative context obtained from PCCreate()
329: - flg - PETSC_TRUE indicates you want the error generated
331: Level: advanced
333: Notes:
334: Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call KSPGetConvergedReason() after a KSPSolve()
335: to determine if it has converged or failed. Or use -ksp_error_if_not_converged to cause the program to terminate as soon as lack of convergence is
336: detected.
338: This is propagated into KSPs used by this PC, which then propagate it into PCs used by those KSPs
340: .keywords: PC, set, initial guess, nonzero
342: .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
343: @*/
344: PetscErrorCode PCSetErrorIfFailure(PC pc,PetscBool flg)
345: {
349: pc->erroriffailure = flg;
350: return(0);
351: }
355: /*@
356: PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the
357: operator during the preconditioning process it applies the Amat provided to TSSetRHSJacobian(),
358: TSSetIJacobian(), SNESSetJacobian(), KSPSetOperator() or PCSetOperator() not the Pmat.
360: Logically Collective on PC
362: Input Parameter:
363: . pc - the preconditioner context
365: Output Parameter:
366: . flg - PETSC_TRUE to use the Amat, PETSC_FALSE to use the Pmat (default is false)
368: Notes:
369: For the common case in which the linear system matrix and the matrix used to construct the
370: preconditioner are identical, this routine is does nothing.
372: Level: intermediate
374: .seealso: PCSetUseAmat(), PCBJACOBI, PGMG, PCFIELDSPLIT, PCCOMPOSITE
375: @*/
376: PetscErrorCode PCGetUseAmat(PC pc,PetscBool *flg)
377: {
380: *flg = pc->useAmat;
381: return(0);
382: }
386: /*@
387: PCCreate - Creates a preconditioner context.
389: Collective on MPI_Comm
391: Input Parameter:
392: . comm - MPI communicator
394: Output Parameter:
395: . pc - location to put the preconditioner context
397: Notes:
398: The default preconditioner for sparse matrices is PCILU or PCICC with 0 fill on one process and block Jacobi with PCILU or ICC
399: in parallel. For dense matrices it is always PCNONE.
401: Level: developer
403: .keywords: PC, create, context
405: .seealso: PCSetUp(), PCApply(), PCDestroy()
406: @*/
407: PetscErrorCode PCCreate(MPI_Comm comm,PC *newpc)
408: {
409: PC pc;
414: *newpc = 0;
415: PCInitializePackage();
417: PetscHeaderCreate(pc,PC_CLASSID,"PC","Preconditioner","PC",comm,PCDestroy,PCView);
419: pc->mat = 0;
420: pc->pmat = 0;
421: pc->setupcalled = 0;
422: pc->setfromoptionscalled = 0;
423: pc->data = 0;
424: pc->diagonalscale = PETSC_FALSE;
425: pc->diagonalscaleleft = 0;
426: pc->diagonalscaleright = 0;
428: pc->modifysubmatrices = 0;
429: pc->modifysubmatricesP = 0;
431: *newpc = pc;
432: return(0);
434: }
436: /* -------------------------------------------------------------------------------*/
440: /*@
441: PCApply - Applies the preconditioner to a vector.
443: Collective on PC and Vec
445: Input Parameters:
446: + pc - the preconditioner context
447: - x - input vector
449: Output Parameter:
450: . y - output vector
452: Level: developer
454: .keywords: PC, apply
456: .seealso: PCApplyTranspose(), PCApplyBAorAB()
457: @*/
458: PetscErrorCode PCApply(PC pc,Vec x,Vec y)
459: {
461: PetscInt m,n,mv,nv;
467: if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
468: if (pc->erroriffailure) {VecValidValues(x,2,PETSC_TRUE);}
469: MatGetLocalSize(pc->mat,&m,&n);
470: VecGetLocalSize(x,&nv);
471: VecGetLocalSize(y,&mv);
472: if (mv != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Preconditioner number of local rows %D does not equal resulting vector number of rows %D",m,mv);
473: if (nv != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Preconditioner number of local columns %D does not equal resulting vector number of rows %D",n,nv);
474: VecLocked(y,3);
476: if (pc->setupcalled < 2) {
477: PCSetUp(pc);
478: }
479: if (!pc->ops->apply) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply");
480: VecLockPush(x);
481: PetscLogEventBegin(PC_Apply,pc,x,y,0);
482: (*pc->ops->apply)(pc,x,y);
483: PetscLogEventEnd(PC_Apply,pc,x,y,0);
484: if (pc->erroriffailure) {VecValidValues(y,3,PETSC_FALSE);}
485: VecLockPop(x);
486: return(0);
487: }
491: /*@
492: PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
494: Collective on PC and Vec
496: Input Parameters:
497: + pc - the preconditioner context
498: - x - input vector
500: Output Parameter:
501: . y - output vector
503: Notes:
504: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
506: Level: developer
508: .keywords: PC, apply, symmetric, left
510: .seealso: PCApply(), PCApplySymmetricRight()
511: @*/
512: PetscErrorCode PCApplySymmetricLeft(PC pc,Vec x,Vec y)
513: {
520: if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
521: if (pc->erroriffailure) {VecValidValues(x,2,PETSC_TRUE);}
522: if (pc->setupcalled < 2) {
523: PCSetUp(pc);
524: }
525: if (!pc->ops->applysymmetricleft) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply");
526: VecLockPush(x);
527: PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);
528: (*pc->ops->applysymmetricleft)(pc,x,y);
529: PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);
530: VecLockPop(x);
531: if (pc->erroriffailure) {VecValidValues(y,3,PETSC_FALSE);}
532: return(0);
533: }
537: /*@
538: PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
540: Collective on PC and Vec
542: Input Parameters:
543: + pc - the preconditioner context
544: - x - input vector
546: Output Parameter:
547: . y - output vector
549: Level: developer
551: Notes:
552: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
554: .keywords: PC, apply, symmetric, right
556: .seealso: PCApply(), PCApplySymmetricLeft()
557: @*/
558: PetscErrorCode PCApplySymmetricRight(PC pc,Vec x,Vec y)
559: {
566: if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
567: if (pc->erroriffailure) {VecValidValues(x,2,PETSC_TRUE);}
568: if (pc->setupcalled < 2) {
569: PCSetUp(pc);
570: }
571: if (!pc->ops->applysymmetricright) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have left symmetric apply");
572: VecLockPush(x);
573: PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);
574: (*pc->ops->applysymmetricright)(pc,x,y);
575: PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);
576: VecLockPop(x);
577: if (pc->erroriffailure) {VecValidValues(y,3,PETSC_FALSE);}
578: return(0);
579: }
583: /*@
584: PCApplyTranspose - Applies the transpose of preconditioner to a vector.
586: Collective on PC and Vec
588: Input Parameters:
589: + pc - the preconditioner context
590: - x - input vector
592: Output Parameter:
593: . y - output vector
595: Notes: For complex numbers this applies the non-Hermitian transpose.
597: Developer Notes: We need to implement a PCApplyHermitianTranspose()
599: Level: developer
601: .keywords: PC, apply, transpose
603: .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCApplyTransposeExists()
604: @*/
605: PetscErrorCode PCApplyTranspose(PC pc,Vec x,Vec y)
606: {
613: if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
614: if (pc->erroriffailure) {VecValidValues(x,2,PETSC_TRUE);}
615: if (pc->setupcalled < 2) {
616: PCSetUp(pc);
617: }
618: if (!pc->ops->applytranspose) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply transpose");
619: VecLockPush(x);
620: PetscLogEventBegin(PC_Apply,pc,x,y,0);
621: (*pc->ops->applytranspose)(pc,x,y);
622: PetscLogEventEnd(PC_Apply,pc,x,y,0);
623: VecLockPop(x);
624: if (pc->erroriffailure) {VecValidValues(y,3,PETSC_FALSE);}
625: return(0);
626: }
630: /*@
631: PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
633: Collective on PC and Vec
635: Input Parameters:
636: . pc - the preconditioner context
638: Output Parameter:
639: . flg - PETSC_TRUE if a transpose operation is defined
641: Level: developer
643: .keywords: PC, apply, transpose
645: .seealso: PCApplyTranspose()
646: @*/
647: PetscErrorCode PCApplyTransposeExists(PC pc,PetscBool *flg)
648: {
652: if (pc->ops->applytranspose) *flg = PETSC_TRUE;
653: else *flg = PETSC_FALSE;
654: return(0);
655: }
659: /*@
660: PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x.
662: Collective on PC and Vec
664: Input Parameters:
665: + pc - the preconditioner context
666: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
667: . x - input vector
668: - work - work vector
670: Output Parameter:
671: . y - output vector
673: Level: developer
675: 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
676: specific KSPSolve() method must also be written to handle the post-solve "correction" for the diagonal scaling.
678: .keywords: PC, apply, operator
680: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
681: @*/
682: PetscErrorCode PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
683: {
691: if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
692: if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric");
693: if (pc->diagonalscale && side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");
694: if (pc->erroriffailure) {VecValidValues(x,3,PETSC_TRUE);}
696: if (pc->setupcalled < 2) {
697: PCSetUp(pc);
698: }
700: if (pc->diagonalscale) {
701: if (pc->ops->applyBA) {
702: Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
703: VecDuplicate(x,&work2);
704: PCDiagonalScaleRight(pc,x,work2);
705: (*pc->ops->applyBA)(pc,side,work2,y,work);
706: PCDiagonalScaleLeft(pc,y,y);
707: VecDestroy(&work2);
708: } else if (side == PC_RIGHT) {
709: PCDiagonalScaleRight(pc,x,y);
710: PCApply(pc,y,work);
711: MatMult(pc->mat,work,y);
712: PCDiagonalScaleLeft(pc,y,y);
713: } else if (side == PC_LEFT) {
714: PCDiagonalScaleRight(pc,x,y);
715: MatMult(pc->mat,y,work);
716: PCApply(pc,work,y);
717: PCDiagonalScaleLeft(pc,y,y);
718: } else if (side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
719: } else {
720: if (pc->ops->applyBA) {
721: (*pc->ops->applyBA)(pc,side,x,y,work);
722: } else if (side == PC_RIGHT) {
723: PCApply(pc,x,work);
724: MatMult(pc->mat,work,y);
725: } else if (side == PC_LEFT) {
726: MatMult(pc->mat,x,work);
727: PCApply(pc,work,y);
728: } else if (side == PC_SYMMETRIC) {
729: /* There's an extra copy here; maybe should provide 2 work vectors instead? */
730: PCApplySymmetricRight(pc,x,work);
731: MatMult(pc->mat,work,y);
732: VecCopy(y,work);
733: PCApplySymmetricLeft(pc,work,y);
734: }
735: }
736: if (pc->erroriffailure) {VecValidValues(y,4,PETSC_FALSE);}
737: return(0);
738: }
742: /*@
743: PCApplyBAorABTranspose - Applies the transpose of the preconditioner
744: and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
745: NOT tr(B*A) = tr(A)*tr(B).
747: Collective on PC and Vec
749: Input Parameters:
750: + pc - the preconditioner context
751: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
752: . x - input vector
753: - work - work vector
755: Output Parameter:
756: . y - output vector
759: 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
760: defined by B'. This is why this has the funny form that it computes tr(B) * tr(A)
762: Level: developer
764: .keywords: PC, apply, operator, transpose
766: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
767: @*/
768: PetscErrorCode PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
769: {
777: if (x == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
778: if (pc->erroriffailure) {VecValidValues(x,3,PETSC_TRUE);}
779: if (pc->ops->applyBAtranspose) {
780: (*pc->ops->applyBAtranspose)(pc,side,x,y,work);
781: if (pc->erroriffailure) {VecValidValues(y,4,PETSC_FALSE);}
782: return(0);
783: }
784: if (side != PC_LEFT && side != PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
786: if (pc->setupcalled < 2) {
787: PCSetUp(pc);
788: }
790: if (side == PC_RIGHT) {
791: PCApplyTranspose(pc,x,work);
792: MatMultTranspose(pc->mat,work,y);
793: } else if (side == PC_LEFT) {
794: MatMultTranspose(pc->mat,x,work);
795: PCApplyTranspose(pc,work,y);
796: }
797: /* add support for PC_SYMMETRIC */
798: if (pc->erroriffailure) {VecValidValues(y,4,PETSC_FALSE);}
799: return(0);
800: }
802: /* -------------------------------------------------------------------------------*/
806: /*@
807: PCApplyRichardsonExists - Determines whether a particular preconditioner has a
808: built-in fast application of Richardson's method.
810: Not Collective
812: Input Parameter:
813: . pc - the preconditioner
815: Output Parameter:
816: . exists - PETSC_TRUE or PETSC_FALSE
818: Level: developer
820: .keywords: PC, apply, Richardson, exists
822: .seealso: PCApplyRichardson()
823: @*/
824: PetscErrorCode PCApplyRichardsonExists(PC pc,PetscBool *exists)
825: {
829: if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
830: else *exists = PETSC_FALSE;
831: return(0);
832: }
836: /*@
837: PCApplyRichardson - Applies several steps of Richardson iteration with
838: the particular preconditioner. This routine is usually used by the
839: Krylov solvers and not the application code directly.
841: Collective on PC
843: Input Parameters:
844: + pc - the preconditioner context
845: . b - the right hand side
846: . w - one work vector
847: . rtol - relative decrease in residual norm convergence criteria
848: . abstol - absolute residual norm convergence criteria
849: . dtol - divergence residual norm increase criteria
850: . its - the number of iterations to apply.
851: - guesszero - if the input x contains nonzero initial guess
853: Output Parameter:
854: + outits - number of iterations actually used (for SOR this always equals its)
855: . reason - the reason the apply terminated
856: - y - the solution (also contains initial guess if guesszero is PETSC_FALSE
858: Notes:
859: Most preconditioners do not support this function. Use the command
860: PCApplyRichardsonExists() to determine if one does.
862: Except for the multigrid PC this routine ignores the convergence tolerances
863: and always runs for the number of iterations
865: Level: developer
867: .keywords: PC, apply, Richardson
869: .seealso: PCApplyRichardsonExists()
870: @*/
871: PetscErrorCode PCApplyRichardson(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
872: {
880: if (b == y) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_IDN,"b and y must be different vectors");
881: if (pc->setupcalled < 2) {
882: PCSetUp(pc);
883: }
884: if (!pc->ops->applyrichardson) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC does not have apply richardson");
885: (*pc->ops->applyrichardson)(pc,b,y,w,rtol,abstol,dtol,its,guesszero,outits,reason);
886: return(0);
887: }
891: /*@
892: PCGetSetUpFailedReason - Gets the reason a PCSetUp() failed or 0 if it did not fail
894: Logically Collective on PC
896: Input Parameter:
897: . pc - the preconditioner context
899: Output Parameter:
900: . reason - the reason it failed, currently only -1
902: Level: advanced
904: .keywords: PC, setup
906: .seealso: PCCreate(), PCApply(), PCDestroy()
907: @*/
908: PetscErrorCode PCGetSetUpFailedReason(PC pc,PetscInt *reason)
909: {
911: if (pc->setupcalled < 0) *reason = pc->setupcalled;
912: else *reason = 0;
913: return(0);
914: }
917: /*
918: a setupcall of 0 indicates never setup,
919: 1 indicates has been previously setup
920: -1 indicates a PCSetUp() was attempted and failed
921: */
924: /*@
925: PCSetUp - Prepares for the use of a preconditioner.
927: Collective on PC
929: Input Parameter:
930: . pc - the preconditioner context
932: Level: developer
934: .keywords: PC, setup
936: .seealso: PCCreate(), PCApply(), PCDestroy()
937: @*/
938: PetscErrorCode PCSetUp(PC pc)
939: {
940: PetscErrorCode ierr;
941: const char *def;
942: PetscObjectState matstate, matnonzerostate;
946: if (!pc->mat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");
948: if (pc->setupcalled && pc->reusepreconditioner) {
949: PetscInfo(pc,"Leaving PC with identical preconditioner since reuse preconditioner is set\n");
950: return(0);
951: }
953: PetscObjectStateGet((PetscObject)pc->pmat,&matstate);
954: MatGetNonzeroState(pc->pmat,&matnonzerostate);
955: if (!pc->setupcalled) {
956: PetscInfo(pc,"Setting up PC for first time\n");
957: pc->flag = DIFFERENT_NONZERO_PATTERN;
958: } else if (matstate == pc->matstate) {
959: PetscInfo(pc,"Leaving PC with identical preconditioner since operator is unchanged\n");
960: return(0);
961: } else {
962: if (matnonzerostate > pc->matnonzerostate) {
963: PetscInfo(pc,"Setting up PC with different nonzero pattern\n");
964: pc->flag = DIFFERENT_NONZERO_PATTERN;
965: } else {
966: PetscInfo(pc,"Setting up PC with same nonzero pattern\n");
967: pc->flag = SAME_NONZERO_PATTERN;
968: }
969: }
970: pc->matstate = matstate;
971: pc->matnonzerostate = matnonzerostate;
973: if (!((PetscObject)pc)->type_name) {
974: PCGetDefaultType_Private(pc,&def);
975: PCSetType(pc,def);
976: }
978: MatSetErrorIfFPE(pc->pmat,pc->erroriffailure);
979: MatSetErrorIfFPE(pc->mat,pc->erroriffailure);
980: PetscLogEventBegin(PC_SetUp,pc,0,0,0);
981: if (pc->ops->setup) {
982: (*pc->ops->setup)(pc);
983: }
984: PetscLogEventEnd(PC_SetUp,pc,0,0,0);
985: if (!pc->setupcalled) pc->setupcalled = 1;
986: return(0);
987: }
991: /*@
992: PCSetUpOnBlocks - Sets up the preconditioner for each block in
993: the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
994: methods.
996: Collective on PC
998: Input Parameters:
999: . pc - the preconditioner context
1001: Level: developer
1003: .keywords: PC, setup, blocks
1005: .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
1006: @*/
1007: PetscErrorCode PCSetUpOnBlocks(PC pc)
1008: {
1013: if (!pc->ops->setuponblocks) return(0);
1014: PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);
1015: (*pc->ops->setuponblocks)(pc);
1016: PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);
1017: return(0);
1018: }
1022: /*@C
1023: PCSetModifySubMatrices - Sets a user-defined routine for modifying the
1024: submatrices that arise within certain subdomain-based preconditioners.
1025: The basic submatrices are extracted from the preconditioner matrix as
1026: usual; the user can then alter these (for example, to set different boundary
1027: conditions for each submatrix) before they are used for the local solves.
1029: Logically Collective on PC
1031: Input Parameters:
1032: + pc - the preconditioner context
1033: . func - routine for modifying the submatrices
1034: - ctx - optional user-defined context (may be null)
1036: Calling sequence of func:
1037: $ func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
1039: . row - an array of index sets that contain the global row numbers
1040: that comprise each local submatrix
1041: . col - an array of index sets that contain the global column numbers
1042: that comprise each local submatrix
1043: . submat - array of local submatrices
1044: - ctx - optional user-defined context for private data for the
1045: user-defined func routine (may be null)
1047: Notes:
1048: PCSetModifySubMatrices() MUST be called before KSPSetUp() and
1049: KSPSolve().
1051: A routine set by PCSetModifySubMatrices() is currently called within
1052: the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
1053: preconditioners. All other preconditioners ignore this routine.
1055: Level: advanced
1057: .keywords: PC, set, modify, submatrices
1059: .seealso: PCModifySubMatrices(), PCASMGetSubMatrices()
1060: @*/
1061: PetscErrorCode PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
1062: {
1065: pc->modifysubmatrices = func;
1066: pc->modifysubmatricesP = ctx;
1067: return(0);
1068: }
1072: /*@C
1073: PCModifySubMatrices - Calls an optional user-defined routine within
1074: certain preconditioners if one has been set with PCSetModifySubMarices().
1076: Collective on PC
1078: Input Parameters:
1079: + pc - the preconditioner context
1080: . nsub - the number of local submatrices
1081: . row - an array of index sets that contain the global row numbers
1082: that comprise each local submatrix
1083: . col - an array of index sets that contain the global column numbers
1084: that comprise each local submatrix
1085: . submat - array of local submatrices
1086: - ctx - optional user-defined context for private data for the
1087: user-defined routine (may be null)
1089: Output Parameter:
1090: . submat - array of local submatrices (the entries of which may
1091: have been modified)
1093: Notes:
1094: The user should NOT generally call this routine, as it will
1095: automatically be called within certain preconditioners (currently
1096: block Jacobi, additive Schwarz) if set.
1098: The basic submatrices are extracted from the preconditioner matrix
1099: as usual; the user can then alter these (for example, to set different
1100: boundary conditions for each submatrix) before they are used for the
1101: local solves.
1103: Level: developer
1105: .keywords: PC, modify, submatrices
1107: .seealso: PCSetModifySubMatrices()
1108: @*/
1109: PetscErrorCode PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
1110: {
1115: if (!pc->modifysubmatrices) return(0);
1116: PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);
1117: (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);
1118: PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);
1119: return(0);
1120: }
1124: /*@
1125: PCSetOperators - Sets the matrix associated with the linear system and
1126: a (possibly) different one associated with the preconditioner.
1128: Logically Collective on PC and Mat
1130: Input Parameters:
1131: + pc - the preconditioner context
1132: . Amat - the matrix that defines the linear system
1133: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
1135: Notes:
1136: Passing a NULL for Amat or Pmat removes the matrix that is currently used.
1138: If you wish to replace either Amat or Pmat but leave the other one untouched then
1139: first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
1140: on it and then pass it back in in your call to KSPSetOperators().
1142: More Notes about Repeated Solution of Linear Systems:
1143: PETSc does NOT reset the matrix entries of either Amat or Pmat
1144: to zero after a linear solve; the user is completely responsible for
1145: matrix assembly. See the routine MatZeroEntries() if desiring to
1146: zero all elements of a matrix.
1148: Level: intermediate
1150: .keywords: PC, set, operators, matrix, linear system
1152: .seealso: PCGetOperators(), MatZeroEntries()
1153: @*/
1154: PetscErrorCode PCSetOperators(PC pc,Mat Amat,Mat Pmat)
1155: {
1156: PetscErrorCode ierr;
1157: PetscInt m1,n1,m2,n2;
1165: if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) {
1166: MatGetLocalSize(Amat,&m1,&n1);
1167: MatGetLocalSize(pc->mat,&m2,&n2);
1168: 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);
1169: MatGetLocalSize(Pmat,&m1,&n1);
1170: MatGetLocalSize(pc->pmat,&m2,&n2);
1171: 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);
1172: }
1174: if (Pmat != pc->pmat) {
1175: /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */
1176: pc->matnonzerostate = -1;
1177: pc->matstate = -1;
1178: }
1180: /* reference first in case the matrices are the same */
1181: if (Amat) {PetscObjectReference((PetscObject)Amat);}
1182: MatDestroy(&pc->mat);
1183: if (Pmat) {PetscObjectReference((PetscObject)Pmat);}
1184: MatDestroy(&pc->pmat);
1185: pc->mat = Amat;
1186: pc->pmat = Pmat;
1187: return(0);
1188: }
1192: /*@
1193: PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed.
1195: Logically Collective on PC
1197: Input Parameters:
1198: + pc - the preconditioner context
1199: - flag - PETSC_TRUE do not compute a new preconditioner, PETSC_FALSE do compute a new preconditioner
1201: Level: intermediate
1203: .seealso: PCGetOperators(), MatZeroEntries(), PCGetReusePreconditioner(), KSPSetReusePreconditioner()
1204: @*/
1205: PetscErrorCode PCSetReusePreconditioner(PC pc,PetscBool flag)
1206: {
1209: pc->reusepreconditioner = flag;
1210: return(0);
1211: }
1215: /*@
1216: PCGetReusePreconditioner - Determines if the PC reuses the current preconditioner even if the operator in the preconditioner has changed.
1218: Not Collective
1220: Input Parameter:
1221: . pc - the preconditioner context
1223: Output Parameter:
1224: . flag - PETSC_TRUE do not compute a new preconditioner, PETSC_FALSE do compute a new preconditioner
1226: Level: intermediate
1228: .seealso: PCGetOperators(), MatZeroEntries(), PCSetReusePreconditioner()
1229: @*/
1230: PetscErrorCode PCGetReusePreconditioner(PC pc,PetscBool *flag)
1231: {
1234: *flag = pc->reusepreconditioner;
1235: return(0);
1236: }
1240: /*@C
1241: PCGetOperators - Gets the matrix associated with the linear system and
1242: possibly a different one associated with the preconditioner.
1244: Not collective, though parallel Mats are returned if the PC is parallel
1246: Input Parameter:
1247: . pc - the preconditioner context
1249: Output Parameters:
1250: + Amat - the matrix defining the linear system
1251: - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat.
1253: Level: intermediate
1255: Notes: Does not increase the reference count of the matrices, so you should not destroy them
1257: Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1258: are created in PC and returned to the user. In this case, if both operators
1259: mat and pmat are requested, two DIFFERENT operators will be returned. If
1260: only one is requested both operators in the PC will be the same (i.e. as
1261: if one had called KSP/PCSetOperators() with the same argument for both Mats).
1262: The user must set the sizes of the returned matrices and their type etc just
1263: as if the user created them with MatCreate(). For example,
1265: $ KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to
1266: $ set size, type, etc of Amat
1268: $ MatCreate(comm,&mat);
1269: $ KSP/PCSetOperators(ksp/pc,Amat,Amat);
1270: $ PetscObjectDereference((PetscObject)mat);
1271: $ set size, type, etc of Amat
1273: and
1275: $ KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to
1276: $ set size, type, etc of Amat and Pmat
1278: $ MatCreate(comm,&Amat);
1279: $ MatCreate(comm,&Pmat);
1280: $ KSP/PCSetOperators(ksp/pc,Amat,Pmat);
1281: $ PetscObjectDereference((PetscObject)Amat);
1282: $ PetscObjectDereference((PetscObject)Pmat);
1283: $ set size, type, etc of Amat and Pmat
1285: The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1286: of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1287: managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1288: at this is when you create a SNES you do not NEED to create a KSP and attach it to
1289: the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1290: you do not need to attach a PC to it (the KSP object manages the PC object for you).
1291: Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1292: it can be created for you?
1295: .keywords: PC, get, operators, matrix, linear system
1297: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1298: @*/
1299: PetscErrorCode PCGetOperators(PC pc,Mat *Amat,Mat *Pmat)
1300: {
1305: if (Amat) {
1306: if (!pc->mat) {
1307: if (pc->pmat && !Pmat) { /* Apmat has been set, but user did not request it, so use for Amat */
1308: pc->mat = pc->pmat;
1309: PetscObjectReference((PetscObject)pc->mat);
1310: } else { /* both Amat and Pmat are empty */
1311: MatCreate(PetscObjectComm((PetscObject)pc),&pc->mat);
1312: if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */
1313: pc->pmat = pc->mat;
1314: PetscObjectReference((PetscObject)pc->pmat);
1315: }
1316: }
1317: }
1318: *Amat = pc->mat;
1319: }
1320: if (Pmat) {
1321: if (!pc->pmat) {
1322: if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */
1323: pc->pmat = pc->mat;
1324: PetscObjectReference((PetscObject)pc->pmat);
1325: } else {
1326: MatCreate(PetscObjectComm((PetscObject)pc),&pc->pmat);
1327: if (!Amat) { /* user did NOT request Amat, so make same as Pmat */
1328: pc->mat = pc->pmat;
1329: PetscObjectReference((PetscObject)pc->mat);
1330: }
1331: }
1332: }
1333: *Pmat = pc->pmat;
1334: }
1335: return(0);
1336: }
1340: /*@C
1341: PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1342: possibly a different one associated with the preconditioner have been set in the PC.
1344: Not collective, though the results on all processes should be the same
1346: Input Parameter:
1347: . pc - the preconditioner context
1349: Output Parameters:
1350: + mat - the matrix associated with the linear system was set
1351: - pmat - matrix associated with the preconditioner was set, usually the same
1353: Level: intermediate
1355: .keywords: PC, get, operators, matrix, linear system
1357: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1358: @*/
1359: PetscErrorCode PCGetOperatorsSet(PC pc,PetscBool *mat,PetscBool *pmat)
1360: {
1363: if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1364: if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1365: return(0);
1366: }
1370: /*@
1371: PCFactorGetMatrix - Gets the factored matrix from the
1372: preconditioner context. This routine is valid only for the LU,
1373: incomplete LU, Cholesky, and incomplete Cholesky methods.
1375: Not Collective on PC though Mat is parallel if PC is parallel
1377: Input Parameters:
1378: . pc - the preconditioner context
1380: Output parameters:
1381: . mat - the factored matrix
1383: Level: advanced
1385: Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1387: .keywords: PC, get, factored, matrix
1388: @*/
1389: PetscErrorCode PCFactorGetMatrix(PC pc,Mat *mat)
1390: {
1396: if (pc->ops->getfactoredmatrix) {
1397: (*pc->ops->getfactoredmatrix)(pc,mat);
1398: } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PC type does not support getting factor matrix");
1399: return(0);
1400: }
1404: /*@C
1405: PCSetOptionsPrefix - Sets the prefix used for searching for all
1406: PC options in the database.
1408: Logically Collective on PC
1410: Input Parameters:
1411: + pc - the preconditioner context
1412: - prefix - the prefix string to prepend to all PC option requests
1414: Notes:
1415: A hyphen (-) must NOT be given at the beginning of the prefix name.
1416: The first character of all runtime options is AUTOMATICALLY the
1417: hyphen.
1419: Level: advanced
1421: .keywords: PC, set, options, prefix, database
1423: .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1424: @*/
1425: PetscErrorCode PCSetOptionsPrefix(PC pc,const char prefix[])
1426: {
1431: PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);
1432: return(0);
1433: }
1437: /*@C
1438: PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1439: PC options in the database.
1441: Logically Collective on PC
1443: Input Parameters:
1444: + pc - the preconditioner context
1445: - prefix - the prefix string to prepend to all PC option requests
1447: Notes:
1448: A hyphen (-) must NOT be given at the beginning of the prefix name.
1449: The first character of all runtime options is AUTOMATICALLY the
1450: hyphen.
1452: Level: advanced
1454: .keywords: PC, append, options, prefix, database
1456: .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1457: @*/
1458: PetscErrorCode PCAppendOptionsPrefix(PC pc,const char prefix[])
1459: {
1464: PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);
1465: return(0);
1466: }
1470: /*@C
1471: PCGetOptionsPrefix - Gets the prefix used for searching for all
1472: PC options in the database.
1474: Not Collective
1476: Input Parameters:
1477: . pc - the preconditioner context
1479: Output Parameters:
1480: . prefix - pointer to the prefix string used, is returned
1482: Notes: On the fortran side, the user should pass in a string 'prifix' of
1483: sufficient length to hold the prefix.
1485: Level: advanced
1487: .keywords: PC, get, options, prefix, database
1489: .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1490: @*/
1491: PetscErrorCode PCGetOptionsPrefix(PC pc,const char *prefix[])
1492: {
1498: PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);
1499: return(0);
1500: }
1504: /*@
1505: PCPreSolve - Optional pre-solve phase, intended for any
1506: preconditioner-specific actions that must be performed before
1507: the iterative solve itself.
1509: Collective on PC
1511: Input Parameters:
1512: + pc - the preconditioner context
1513: - ksp - the Krylov subspace context
1515: Level: developer
1517: Sample of Usage:
1518: .vb
1519: PCPreSolve(pc,ksp);
1520: KSPSolve(ksp,b,x);
1521: PCPostSolve(pc,ksp);
1522: .ve
1524: Notes:
1525: The pre-solve phase is distinct from the PCSetUp() phase.
1527: KSPSolve() calls this directly, so is rarely called by the user.
1529: .keywords: PC, pre-solve
1531: .seealso: PCPostSolve()
1532: @*/
1533: PetscErrorCode PCPreSolve(PC pc,KSP ksp)
1534: {
1536: Vec x,rhs;
1541: pc->presolvedone++;
1542: if (pc->presolvedone > 2) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot embed PCPreSolve() more than twice");
1543: KSPGetSolution(ksp,&x);
1544: KSPGetRhs(ksp,&rhs);
1546: if (pc->ops->presolve) {
1547: (*pc->ops->presolve)(pc,ksp,rhs,x);
1548: }
1549: return(0);
1550: }
1554: /*@
1555: PCPostSolve - Optional post-solve phase, intended for any
1556: preconditioner-specific actions that must be performed after
1557: the iterative solve itself.
1559: Collective on PC
1561: Input Parameters:
1562: + pc - the preconditioner context
1563: - ksp - the Krylov subspace context
1565: Sample of Usage:
1566: .vb
1567: PCPreSolve(pc,ksp);
1568: KSPSolve(ksp,b,x);
1569: PCPostSolve(pc,ksp);
1570: .ve
1572: Note:
1573: KSPSolve() calls this routine directly, so it is rarely called by the user.
1575: Level: developer
1577: .keywords: PC, post-solve
1579: .seealso: PCPreSolve(), KSPSolve()
1580: @*/
1581: PetscErrorCode PCPostSolve(PC pc,KSP ksp)
1582: {
1584: Vec x,rhs;
1589: pc->presolvedone--;
1590: KSPGetSolution(ksp,&x);
1591: KSPGetRhs(ksp,&rhs);
1592: if (pc->ops->postsolve) {
1593: (*pc->ops->postsolve)(pc,ksp,rhs,x);
1594: }
1595: return(0);
1596: }
1600: /*@C
1601: PCLoad - Loads a PC that has been stored in binary with PCView().
1603: Collective on PetscViewer
1605: Input Parameters:
1606: + newdm - the newly loaded PC, this needs to have been created with PCCreate() or
1607: some related function before a call to PCLoad().
1608: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
1610: Level: intermediate
1612: Notes:
1613: The type is determined by the data in the file, any type set into the PC before this call is ignored.
1615: Notes for advanced users:
1616: Most users should not need to know the details of the binary storage
1617: format, since PCLoad() and PCView() completely hide these details.
1618: But for anyone who's interested, the standard binary matrix storage
1619: format is
1620: .vb
1621: has not yet been determined
1622: .ve
1624: .seealso: PetscViewerBinaryOpen(), PCView(), MatLoad(), VecLoad()
1625: @*/
1626: PetscErrorCode PCLoad(PC newdm, PetscViewer viewer)
1627: {
1629: PetscBool isbinary;
1630: PetscInt classid;
1631: char type[256];
1636: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1637: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1639: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
1640: if (classid != PC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not PC next in file");
1641: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
1642: PCSetType(newdm, type);
1643: if (newdm->ops->load) {
1644: (*newdm->ops->load)(newdm,viewer);
1645: }
1646: return(0);
1647: }
1649: #include <petscdraw.h>
1650: #if defined(PETSC_HAVE_SAWS)
1651: #include <petscviewersaws.h>
1652: #endif
1655: /*@C
1656: PCView - Prints the PC data structure.
1658: Collective on PC
1660: Input Parameters:
1661: + PC - the PC context
1662: - viewer - optional visualization context
1664: Note:
1665: The available visualization contexts include
1666: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
1667: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1668: output where only the first processor opens
1669: the file. All other processors send their
1670: data to the first processor to print.
1672: The user can open an alternative visualization contexts with
1673: PetscViewerASCIIOpen() (output to a specified file).
1675: Level: developer
1677: .keywords: PC, view
1679: .seealso: KSPView(), PetscViewerASCIIOpen()
1680: @*/
1681: PetscErrorCode PCView(PC pc,PetscViewer viewer)
1682: {
1683: PCType cstr;
1684: PetscErrorCode ierr;
1685: PetscBool iascii,isstring,isbinary,isdraw;
1686: PetscViewerFormat format;
1687: #if defined(PETSC_HAVE_SAWS)
1688: PetscBool issaws;
1689: #endif
1693: if (!viewer) {
1694: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&viewer);
1695: }
1699: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1700: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
1701: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
1702: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
1703: #if defined(PETSC_HAVE_SAWS)
1704: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
1705: #endif
1707: if (iascii) {
1708: PetscViewerGetFormat(viewer,&format);
1709: PetscObjectPrintClassNamePrefixType((PetscObject)pc,viewer);
1710: if (!pc->setupcalled) {
1711: PetscViewerASCIIPrintf(viewer," PC has not been set up so information may be incomplete\n");
1712: }
1713: if (pc->ops->view) {
1714: PetscViewerASCIIPushTab(viewer);
1715: (*pc->ops->view)(pc,viewer);
1716: PetscViewerASCIIPopTab(viewer);
1717: }
1718: if (pc->mat) {
1719: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1720: if (pc->pmat == pc->mat) {
1721: PetscViewerASCIIPrintf(viewer," linear system matrix = precond matrix:\n");
1722: PetscViewerASCIIPushTab(viewer);
1723: MatView(pc->mat,viewer);
1724: PetscViewerASCIIPopTab(viewer);
1725: } else {
1726: if (pc->pmat) {
1727: PetscViewerASCIIPrintf(viewer," linear system matrix followed by preconditioner matrix:\n");
1728: } else {
1729: PetscViewerASCIIPrintf(viewer," linear system matrix:\n");
1730: }
1731: PetscViewerASCIIPushTab(viewer);
1732: MatView(pc->mat,viewer);
1733: if (pc->pmat) {MatView(pc->pmat,viewer);}
1734: PetscViewerASCIIPopTab(viewer);
1735: }
1736: PetscViewerPopFormat(viewer);
1737: }
1738: } else if (isstring) {
1739: PCGetType(pc,&cstr);
1740: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
1741: if (pc->ops->view) {(*pc->ops->view)(pc,viewer);}
1742: } else if (isbinary) {
1743: PetscInt classid = PC_FILE_CLASSID;
1744: MPI_Comm comm;
1745: PetscMPIInt rank;
1746: char type[256];
1748: PetscObjectGetComm((PetscObject)pc,&comm);
1749: MPI_Comm_rank(comm,&rank);
1750: if (!rank) {
1751: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
1752: PetscStrncpy(type,((PetscObject)pc)->type_name,256);
1753: PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);
1754: }
1755: if (pc->ops->view) {
1756: (*pc->ops->view)(pc,viewer);
1757: }
1758: } else if (isdraw) {
1759: PetscDraw draw;
1760: char str[25];
1761: PetscReal x,y,bottom,h;
1762: PetscInt n;
1764: PetscViewerDrawGetDraw(viewer,0,&draw);
1765: PetscDrawGetCurrentPoint(draw,&x,&y);
1766: if (pc->mat) {
1767: MatGetSize(pc->mat,&n,NULL);
1768: PetscSNPrintf(str,25,"PC: %s (%D)",((PetscObject)pc)->type_name,n);
1769: } else {
1770: PetscSNPrintf(str,25,"PC: %s",((PetscObject)pc)->type_name);
1771: }
1772: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_RED,PETSC_DRAW_BLACK,str,NULL,&h);
1773: bottom = y - h;
1774: PetscDrawPushCurrentPoint(draw,x,bottom);
1775: if (pc->ops->view) {
1776: (*pc->ops->view)(pc,viewer);
1777: }
1778: PetscDrawPopCurrentPoint(draw);
1779: #if defined(PETSC_HAVE_SAWS)
1780: } else if (issaws) {
1781: PetscMPIInt rank;
1783: PetscObjectName((PetscObject)pc);
1784: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
1785: if (!((PetscObject)pc)->amsmem && !rank) {
1786: PetscObjectViewSAWs((PetscObject)pc,viewer);
1787: }
1788: if (pc->mat) {MatView(pc->mat,viewer);}
1789: if (pc->pmat && pc->pmat != pc->mat) {MatView(pc->pmat,viewer);}
1790: #endif
1791: }
1792: return(0);
1793: }
1798: /*@
1799: PCSetInitialGuessNonzero - Tells the iterative solver that the
1800: initial guess is nonzero; otherwise PC assumes the initial guess
1801: is to be zero (and thus zeros it out before solving).
1803: Logically Collective on PC
1805: Input Parameters:
1806: + pc - iterative context obtained from PCCreate()
1807: - flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1809: Level: Developer
1811: Notes:
1812: This is a weird function. Since PC's are linear operators on the right hand side they
1813: CANNOT use an initial guess. This function is for the "pass-through" preconditioners
1814: PCKSP and PCREDUNDANT and causes the inner KSP object to use the nonzero
1815: initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP.
1818: .keywords: PC, set, initial guess, nonzero
1820: .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
1821: @*/
1822: PetscErrorCode PCSetInitialGuessNonzero(PC pc,PetscBool flg)
1823: {
1826: pc->nonzero_guess = flg;
1827: return(0);
1828: }
1832: /*@
1833: PCGetInitialGuessNonzero - Determines if the iterative solver assumes that the
1834: initial guess is nonzero; otherwise PC assumes the initial guess
1835: is to be zero (and thus zeros it out before solving).
1837: Logically Collective on PC
1839: Input Parameter:
1840: . pc - iterative context obtained from PCCreate()
1842: Output Parameter:
1843: . flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1845: Level: Developer
1847: .keywords: PC, set, initial guess, nonzero
1849: .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll(), PCSetInitialGuessNonzero()
1850: @*/
1851: PetscErrorCode PCGetInitialGuessNonzero(PC pc,PetscBool *flg)
1852: {
1854: *flg = pc->nonzero_guess;
1855: return(0);
1856: }
1860: /*@C
1861: PCRegister - Adds a method to the preconditioner package.
1863: Not collective
1865: Input Parameters:
1866: + name_solver - name of a new user-defined solver
1867: - routine_create - routine to create method context
1869: Notes:
1870: PCRegister() may be called multiple times to add several user-defined preconditioners.
1872: Sample usage:
1873: .vb
1874: PCRegister("my_solver", MySolverCreate);
1875: .ve
1877: Then, your solver can be chosen with the procedural interface via
1878: $ PCSetType(pc,"my_solver")
1879: or at runtime via the option
1880: $ -pc_type my_solver
1882: Level: advanced
1884: .keywords: PC, register
1886: .seealso: PCRegisterAll(), PCRegisterDestroy()
1887: @*/
1888: PetscErrorCode PCRegister(const char sname[],PetscErrorCode (*function)(PC))
1889: {
1893: PetscFunctionListAdd(&PCList,sname,function);
1894: return(0);
1895: }
1899: /*@
1900: PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1902: Collective on PC
1904: Input Parameter:
1905: . pc - the preconditioner object
1907: Output Parameter:
1908: . mat - the explict preconditioned operator
1910: Notes:
1911: This computation is done by applying the operators to columns of the
1912: identity matrix.
1914: Currently, this routine uses a dense matrix format when 1 processor
1915: is used and a sparse format otherwise. This routine is costly in general,
1916: and is recommended for use only with relatively small systems.
1918: Level: advanced
1920: .keywords: PC, compute, explicit, operator
1922: .seealso: KSPComputeExplicitOperator()
1924: @*/
1925: PetscErrorCode PCComputeExplicitOperator(PC pc,Mat *mat)
1926: {
1927: Vec in,out;
1929: PetscInt i,M,m,*rows,start,end;
1930: PetscMPIInt size;
1931: MPI_Comm comm;
1932: PetscScalar *array,one = 1.0;
1938: PetscObjectGetComm((PetscObject)pc,&comm);
1939: MPI_Comm_size(comm,&size);
1941: if (!pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1942: MatCreateVecs(pc->pmat,&in,0);
1943: VecDuplicate(in,&out);
1944: VecGetOwnershipRange(in,&start,&end);
1945: VecGetSize(in,&M);
1946: VecGetLocalSize(in,&m);
1947: PetscMalloc1(m+1,&rows);
1948: for (i=0; i<m; i++) rows[i] = start + i;
1950: MatCreate(comm,mat);
1951: MatSetSizes(*mat,m,m,M,M);
1952: if (size == 1) {
1953: MatSetType(*mat,MATSEQDENSE);
1954: MatSeqDenseSetPreallocation(*mat,NULL);
1955: } else {
1956: MatSetType(*mat,MATMPIAIJ);
1957: MatMPIAIJSetPreallocation(*mat,0,NULL,0,NULL);
1958: }
1959: MatSetOption(*mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
1961: for (i=0; i<M; i++) {
1963: VecSet(in,0.0);
1964: VecSetValues(in,1,&i,&one,INSERT_VALUES);
1965: VecAssemblyBegin(in);
1966: VecAssemblyEnd(in);
1968: /* should fix, allowing user to choose side */
1969: PCApply(pc,in,out);
1971: VecGetArray(out,&array);
1972: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
1973: VecRestoreArray(out,&array);
1975: }
1976: PetscFree(rows);
1977: VecDestroy(&out);
1978: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
1979: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
1980: return(0);
1981: }
1985: /*@
1986: PCSetCoordinates - sets the coordinates of all the nodes on the local process
1988: Collective on PC
1990: Input Parameters:
1991: + pc - the solver context
1992: . dim - the dimension of the coordinates 1, 2, or 3
1993: - coords - the coordinates
1995: Level: intermediate
1997: Notes: coords is an array of the 3D coordinates for the nodes on
1998: the local processor. So if there are 108 equation on a processor
1999: for a displacement finite element discretization of elasticity (so
2000: that there are 36 = 108/3 nodes) then the array must have 108
2001: double precision values (ie, 3 * 36). These x y z coordinates
2002: should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x,
2003: ... , N-1.z ].
2005: .seealso: MatSetNearNullSpace
2006: @*/
2007: PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
2008: {
2012: PetscTryMethod(pc,"PCSetCoordinates_C",(PC,PetscInt,PetscInt,PetscReal*),(pc,dim,nloc,coords));
2013: return(0);
2014: }