Actual source code: precon.c
1: /*
2: The PC (preconditioner) interface routines, callable by users.
3: */
4: #include src/ksp/pc/pcimpl.h
6: /* Logging support */
7: PetscCookie PC_COOKIE = 0;
8: PetscEvent PC_SetUp = 0, PC_SetUpOnBlocks = 0, PC_Apply = 0, PC_ApplyCoarse = 0, PC_ApplyMultiple = 0, PC_ApplySymmetricLeft = 0;
9: PetscEvent PC_ApplySymmetricRight = 0, PC_ModifySubMatrices = 0;
13: PetscErrorCode PCGetDefaultType_Private(PC pc,const char* type[])
14: {
16: PetscMPIInt size;
17: PetscTruth flg1,flg2,set,flg3;
20: MPI_Comm_size(pc->comm,&size);
21: if (pc->pmat) {
22: PetscErrorCode (*f)(Mat,PetscTruth*,MatReuse,Mat*);
23: PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",(void (**)(void))&f);
24: if (size == 1) {
25: MatHasOperation(pc->pmat,MATOP_ICCFACTOR_SYMBOLIC,&flg1);
26: MatHasOperation(pc->pmat,MATOP_ILUFACTOR_SYMBOLIC,&flg2);
27: MatIsSymmetricKnown(pc->pmat,&set,&flg3);
28: if (flg1 && (!flg2 || (set && flg3))) {
29: *type = PCICC;
30: } else if (flg2) {
31: *type = PCILU;
32: } else if (f) { /* likely is a parallel matrix run on one processor */
33: *type = PCBJACOBI;
34: } else {
35: *type = PCNONE;
36: }
37: } else {
38: if (f) {
39: *type = PCBJACOBI;
40: } else {
41: *type = PCNONE;
42: }
43: }
44: } else {
45: if (size == 1) {
46: *type = PCILU;
47: } else {
48: *type = PCBJACOBI;
49: }
50: }
51: return(0);
52: }
56: /*@C
57: PCDestroy - Destroys PC context that was created with PCCreate().
59: Collective on PC
61: Input Parameter:
62: . pc - the preconditioner context
64: Level: developer
66: .keywords: PC, destroy
68: .seealso: PCCreate(), PCSetUp()
69: @*/
70: PetscErrorCode PCDestroy(PC pc)
71: {
76: if (--pc->refct > 0) return(0);
78: /* if memory was published with AMS then destroy it */
79: PetscObjectDepublish(pc);
81: if (pc->ops->destroy) { (*pc->ops->destroy)(pc);}
82: if (pc->diagonalscaleright) {VecDestroy(pc->diagonalscaleright);}
83: if (pc->diagonalscaleleft) {VecDestroy(pc->diagonalscaleleft);}
85: PetscLogObjectDestroy(pc);
86: PetscHeaderDestroy(pc);
87: return(0);
88: }
92: /*@C
93: PCDiagonalScale - Indicates if the preconditioner applies an additional left and right
94: scaling as needed by certain time-stepping codes.
96: Collective on PC
98: Input Parameter:
99: . pc - the preconditioner context
101: Output Parameter:
102: . flag - PETSC_TRUE if it applies the scaling
104: Level: developer
106: Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
107: $ D M A D^{-1} y = D M b for left preconditioning or
108: $ D A M D^{-1} z = D b for right preconditioning
110: .keywords: PC
112: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScaleSet()
113: @*/
114: PetscErrorCode PCDiagonalScale(PC pc,PetscTruth *flag)
115: {
119: *flag = pc->diagonalscale;
120: return(0);
121: }
125: /*@
126: PCDiagonalScaleSet - Indicates the left scaling to use to apply an additional left and right
127: scaling as needed by certain time-stepping codes.
129: Collective on PC
131: Input Parameters:
132: + pc - the preconditioner context
133: - s - scaling vector
135: Level: intermediate
137: Notes: The system solved via the Krylov method is
138: $ D M A D^{-1} y = D M b for left preconditioning or
139: $ D A M D^{-1} z = D b for right preconditioning
141: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
143: .keywords: PC
145: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScale()
146: @*/
147: PetscErrorCode PCDiagonalScaleSet(PC pc,Vec s)
148: {
154: pc->diagonalscale = PETSC_TRUE;
155: if (pc->diagonalscaleleft) {
156: VecDestroy(pc->diagonalscaleleft);
157: }
158: pc->diagonalscaleleft = s;
159: PetscObjectReference((PetscObject)s);
160: if (!pc->diagonalscaleright) {
161: VecDuplicate(s,&pc->diagonalscaleright);
162: }
163: VecCopy(s,pc->diagonalscaleright);
164: VecReciprocal(pc->diagonalscaleright);
165: return(0);
166: }
170: /*@C
171: PCDiagonalScaleLeft - Indicates the left scaling to use to apply an additional left and right
172: scaling as needed by certain time-stepping codes.
174: Collective on PC
176: Input Parameters:
177: + pc - the preconditioner context
178: . in - input vector
179: + out - scaled vector (maybe the same as in)
181: Level: intermediate
183: Notes: The system solved via the Krylov method is
184: $ D M A D^{-1} y = D M b for left preconditioning or
185: $ D A M D^{-1} z = D b for right preconditioning
187: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
189: If diagonal scaling is turned off and in is not out then in is copied to out
191: .keywords: PC
193: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
194: @*/
195: PetscErrorCode PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
196: {
203: if (pc->diagonalscale) {
204: VecPointwiseMult(pc->diagonalscaleleft,in,out);
205: } else if (in != out) {
206: VecCopy(in,out);
207: }
208: return(0);
209: }
213: /*@C
214: PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
216: Collective on PC
218: Input Parameters:
219: + pc - the preconditioner context
220: . in - input vector
221: + out - scaled vector (maybe the same as in)
223: Level: intermediate
225: Notes: The system solved via the Krylov method is
226: $ D M A D^{-1} y = D M b for left preconditioning or
227: $ D A M D^{-1} z = D b for right preconditioning
229: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
231: If diagonal scaling is turned off and in is not out then in is copied to out
233: .keywords: PC
235: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
236: @*/
237: PetscErrorCode PCDiagonalScaleRight(PC pc,Vec in,Vec out)
238: {
245: if (pc->diagonalscale) {
246: VecPointwiseMult(pc->diagonalscaleright,in,out);
247: } else if (in != out) {
248: VecCopy(in,out);
249: }
250: return(0);
251: }
255: static PetscErrorCode PCPublish_Petsc(PetscObject obj)
256: {
258: return(0);
259: }
263: /*@C
264: PCCreate - Creates a preconditioner context.
266: Collective on MPI_Comm
268: Input Parameter:
269: . comm - MPI communicator
271: Output Parameter:
272: . pc - location to put the preconditioner context
274: Notes:
275: The default preconditioner on one processor is PCILU with 0 fill on more
276: then one it is PCBJACOBI with ILU() on each processor.
278: Level: developer
280: .keywords: PC, create, context
282: .seealso: PCSetUp(), PCApply(), PCDestroy()
283: @*/
284: PetscErrorCode PCCreate(MPI_Comm comm,PC *newpc)
285: {
286: PC pc;
291: *newpc = 0;
292: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
293: PCInitializePackage(PETSC_NULL);
294: #endif
296: PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_COOKIE,-1,"PC",comm,PCDestroy,PCView);
297: PetscLogObjectCreate(pc);
298: pc->bops->publish = PCPublish_Petsc;
299: pc->mat = 0;
300: pc->setupcalled = 0;
301: pc->data = 0;
302: pc->diagonalscale = PETSC_FALSE;
303: pc->diagonalscaleleft = 0;
304: pc->diagonalscaleright = 0;
306: pc->ops->destroy = 0;
307: pc->ops->apply = 0;
308: pc->ops->applytranspose = 0;
309: pc->ops->applyBA = 0;
310: pc->ops->applyBAtranspose = 0;
311: pc->ops->applyrichardson = 0;
312: pc->ops->view = 0;
313: pc->ops->getfactoredmatrix = 0;
314: pc->ops->applysymmetricright = 0;
315: pc->ops->applysymmetricleft = 0;
316: pc->ops->setuponblocks = 0;
318: pc->modifysubmatrices = 0;
319: pc->modifysubmatricesP = 0;
320: *newpc = pc;
321: PetscPublishAll(pc);
322: return(0);
324: }
326: /* -------------------------------------------------------------------------------*/
330: /*@
331: PCApply - Applies the preconditioner to a vector.
333: Collective on PC and Vec
335: Input Parameters:
336: + pc - the preconditioner context
337: - x - input vector
339: Output Parameter:
340: . y - output vector
342: Level: developer
344: .keywords: PC, apply
346: .seealso: PCApplyTranspose(), PCApplyBAorAB()
347: @*/
348: PetscErrorCode PCApply(PC pc,Vec x,Vec y)
349: {
356: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
358: if (pc->setupcalled < 2) {
359: PCSetUp(pc);
360: }
361: PetscLogEventBegin(PC_Apply,pc,x,y,0);
362: (*pc->ops->apply)(pc,x,y);
363: PetscLogEventEnd(PC_Apply,pc,x,y,0);
364: return(0);
365: }
369: /*@
370: PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
372: Collective on PC and Vec
374: Input Parameters:
375: + pc - the preconditioner context
376: - x - input vector
378: Output Parameter:
379: . y - output vector
381: Notes:
382: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
384: Level: developer
386: .keywords: PC, apply, symmetric, left
388: .seealso: PCApply(), PCApplySymmetricRight()
389: @*/
390: PetscErrorCode PCApplySymmetricLeft(PC pc,Vec x,Vec y)
391: {
398: if (!pc->ops->applysymmetricleft) SETERRQ(PETSC_ERR_SUP,"PC does not have left symmetric apply");
400: if (pc->setupcalled < 2) {
401: PCSetUp(pc);
402: }
404: PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);
405: (*pc->ops->applysymmetricleft)(pc,x,y);
406: PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);
407: return(0);
408: }
412: /*@
413: PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
415: Collective on PC and Vec
417: Input Parameters:
418: + pc - the preconditioner context
419: - x - input vector
421: Output Parameter:
422: . y - output vector
424: Level: developer
426: Notes:
427: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
429: .keywords: PC, apply, symmetric, right
431: .seealso: PCApply(), PCApplySymmetricLeft()
432: @*/
433: PetscErrorCode PCApplySymmetricRight(PC pc,Vec x,Vec y)
434: {
441: if (!pc->ops->applysymmetricright) SETERRQ(PETSC_ERR_SUP,"PC does not have left symmetric apply");
443: if (pc->setupcalled < 2) {
444: PCSetUp(pc);
445: }
447: PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);
448: (*pc->ops->applysymmetricright)(pc,x,y);
449: PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);
450: return(0);
451: }
455: /*@
456: PCApplyTranspose - Applies the transpose of preconditioner to a vector.
458: Collective on PC and Vec
460: Input Parameters:
461: + pc - the preconditioner context
462: - x - input vector
464: Output Parameter:
465: . y - output vector
467: Level: developer
469: .keywords: PC, apply, transpose
471: .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCHasApplyTranspose()
472: @*/
473: PetscErrorCode PCApplyTranspose(PC pc,Vec x,Vec y)
474: {
481: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
482: if (!pc->ops->applytranspose) SETERRQ(PETSC_ERR_SUP," ");
484: if (pc->setupcalled < 2) {
485: PCSetUp(pc);
486: }
488: PetscLogEventBegin(PC_Apply,pc,x,y,0);
489: (*pc->ops->applytranspose)(pc,x,y);
490: PetscLogEventEnd(PC_Apply,pc,x,y,0);
491: return(0);
492: }
496: /*@
497: PCHasApplyTranspose - Test whether the preconditioner has a transpose apply operation
499: Collective on PC and Vec
501: Input Parameters:
502: . pc - the preconditioner context
504: Output Parameter:
505: . flg - PETSC_TRUE if a transpose operation is defined
507: Level: developer
509: .keywords: PC, apply, transpose
511: .seealso: PCApplyTranspose()
512: @*/
513: PetscErrorCode PCHasApplyTranspose(PC pc,PetscTruth *flg)
514: {
518: *flg = (PetscTruth) (pc->ops->applytranspose == 0);
519: return(0);
520: }
524: /*@
525: PCApplyBAorAB - Applies the preconditioner and operator to a vector.
527: Collective on PC and Vec
529: Input Parameters:
530: + pc - the preconditioner context
531: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
532: . x - input vector
533: - work - work vector
535: Output Parameter:
536: . y - output vector
538: Level: developer
540: .keywords: PC, apply, operator
542: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
543: @*/
544: PetscErrorCode PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
545: {
553: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
554: if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) {
555: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric");
556: }
557: if (pc->diagonalscale && side == PC_SYMMETRIC) {
558: SETERRQ(PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");
559: }
561: if (pc->setupcalled < 2) {
562: PCSetUp(pc);
563: }
565: if (pc->diagonalscale) {
566: if (pc->ops->applyBA) {
567: Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
568: VecDuplicate(x,&work2);
569: PCDiagonalScaleRight(pc,x,work2);
570: (*pc->ops->applyBA)(pc,side,work2,y,work);
571: PCDiagonalScaleLeft(pc,y,y);
572: VecDestroy(work2);
573: } else if (side == PC_RIGHT) {
574: PCDiagonalScaleRight(pc,x,y);
575: PCApply(pc,y,work);
576: MatMult(pc->mat,work,y);
577: PCDiagonalScaleLeft(pc,y,y);
578: } else if (side == PC_LEFT) {
579: PCDiagonalScaleRight(pc,x,y);
580: MatMult(pc->mat,y,work);
581: PCApply(pc,work,y);
582: PCDiagonalScaleLeft(pc,y,y);
583: } else if (side == PC_SYMMETRIC) {
584: SETERRQ(PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
585: }
586: } else {
587: if (pc->ops->applyBA) {
588: (*pc->ops->applyBA)(pc,side,x,y,work);
589: } else if (side == PC_RIGHT) {
590: PCApply(pc,x,work);
591: MatMult(pc->mat,work,y);
592: } else if (side == PC_LEFT) {
593: MatMult(pc->mat,x,work);
594: PCApply(pc,work,y);
595: } else if (side == PC_SYMMETRIC) {
596: /* There's an extra copy here; maybe should provide 2 work vectors instead? */
597: PCApplySymmetricRight(pc,x,work);
598: MatMult(pc->mat,work,y);
599: VecCopy(y,work);
600: PCApplySymmetricLeft(pc,work,y);
601: }
602: }
603: return(0);
604: }
608: /*@
609: PCApplyBAorABTranspose - Applies the transpose of the preconditioner
610: and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
611: not tr(B*A) = tr(A)*tr(B).
613: Collective on PC and Vec
615: Input Parameters:
616: + pc - the preconditioner context
617: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
618: . x - input vector
619: - work - work vector
621: Output Parameter:
622: . y - output vector
624: Level: developer
626: .keywords: PC, apply, operator, transpose
628: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
629: @*/
630: PetscErrorCode PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
631: {
639: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
640: if (pc->ops->applyBAtranspose) {
641: (*pc->ops->applyBAtranspose)(pc,side,x,y,work);
642: return(0);
643: }
644: if (side != PC_LEFT && side != PC_RIGHT) {
645: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
646: }
648: if (pc->setupcalled < 2) {
649: PCSetUp(pc);
650: }
652: if (side == PC_RIGHT) {
653: PCApplyTranspose(pc,x,work);
654: MatMultTranspose(pc->mat,work,y);
655: } else if (side == PC_LEFT) {
656: MatMultTranspose(pc->mat,x,work);
657: PCApplyTranspose(pc,work,y);
658: }
659: /* add support for PC_SYMMETRIC */
660: return(0); /* actually will never get here */
661: }
663: /* -------------------------------------------------------------------------------*/
667: /*@
668: PCApplyRichardsonExists - Determines whether a particular preconditioner has a
669: built-in fast application of Richardson's method.
671: Not Collective
673: Input Parameter:
674: . pc - the preconditioner
676: Output Parameter:
677: . exists - PETSC_TRUE or PETSC_FALSE
679: Level: developer
681: .keywords: PC, apply, Richardson, exists
683: .seealso: PCApplyRichardson()
684: @*/
685: PetscErrorCode PCApplyRichardsonExists(PC pc,PetscTruth *exists)
686: {
690: if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
691: else *exists = PETSC_FALSE;
692: return(0);
693: }
697: /*@
698: PCApplyRichardson - Applies several steps of Richardson iteration with
699: the particular preconditioner. This routine is usually used by the
700: Krylov solvers and not the application code directly.
702: Collective on PC
704: Input Parameters:
705: + pc - the preconditioner context
706: . x - the initial guess
707: . w - one work vector
708: . rtol - relative decrease in residual norm convergence criteria
709: . abstol - absolute residual norm convergence criteria
710: . dtol - divergence residual norm increase criteria
711: - its - the number of iterations to apply.
713: Output Parameter:
714: . y - the solution
716: Notes:
717: Most preconditioners do not support this function. Use the command
718: PCApplyRichardsonExists() to determine if one does.
720: Except for the multigrid PC this routine ignores the convergence tolerances
721: and always runs for the number of iterations
722:
723: Level: developer
725: .keywords: PC, apply, Richardson
727: .seealso: PCApplyRichardsonExists()
728: @*/
729: PetscErrorCode PCApplyRichardson(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
730: {
738: if (!pc->ops->applyrichardson) SETERRQ(PETSC_ERR_SUP," ");
740: if (pc->setupcalled < 2) {
741: PCSetUp(pc);
742: }
744: (*pc->ops->applyrichardson)(pc,x,y,w,rtol,abstol,dtol,its);
745: return(0);
746: }
748: /*
749: a setupcall of 0 indicates never setup,
750: 1 needs to be resetup,
751: 2 does not need any changes.
752: */
755: /*@
756: PCSetUp - Prepares for the use of a preconditioner.
758: Collective on PC
760: Input Parameter:
761: . pc - the preconditioner context
763: Level: developer
765: .keywords: PC, setup
767: .seealso: PCCreate(), PCApply(), PCDestroy()
768: @*/
769: PetscErrorCode PCSetUp(PC pc)
770: {
772: const char *def;
777: if (pc->setupcalled > 1) {
778: PetscLogInfo(pc,"PCSetUp:Setting PC with identical preconditioner\n");
779: return(0);
780: } else if (!pc->setupcalled) {
781: PetscLogInfo(pc,"PCSetUp:Setting up new PC\n");
782: } else if (pc->flag == SAME_NONZERO_PATTERN) {
783: PetscLogInfo(pc,"PCSetUp:Setting up PC with same nonzero pattern\n");
784: } else {
785: PetscLogInfo(pc,"PCSetUp:Setting up PC with different nonzero pattern\n");
786: }
788: PetscLogEventBegin(PC_SetUp,pc,0,0,0);
789: if (!pc->mat) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
791: if (!pc->type_name) {
792: PCGetDefaultType_Private(pc,&def);
793: PCSetType(pc,def);
794: }
796: if (pc->ops->setup) {
797: (*pc->ops->setup)(pc);
798: }
799: pc->setupcalled = 2;
800: PetscLogEventEnd(PC_SetUp,pc,0,0,0);
801: return(0);
802: }
806: /*@
807: PCSetUpOnBlocks - Sets up the preconditioner for each block in
808: the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
809: methods.
811: Collective on PC
813: Input Parameters:
814: . pc - the preconditioner context
816: Level: developer
818: .keywords: PC, setup, blocks
820: .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
821: @*/
822: PetscErrorCode PCSetUpOnBlocks(PC pc)
823: {
828: if (!pc->ops->setuponblocks) return(0);
829: PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);
830: (*pc->ops->setuponblocks)(pc);
831: PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);
832: return(0);
833: }
837: /*@C
838: PCSetModifySubMatrices - Sets a user-defined routine for modifying the
839: submatrices that arise within certain subdomain-based preconditioners.
840: The basic submatrices are extracted from the preconditioner matrix as
841: usual; the user can then alter these (for example, to set different boundary
842: conditions for each submatrix) before they are used for the local solves.
844: Collective on PC
846: Input Parameters:
847: + pc - the preconditioner context
848: . func - routine for modifying the submatrices
849: - ctx - optional user-defined context (may be null)
851: Calling sequence of func:
852: $ func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
854: . row - an array of index sets that contain the global row numbers
855: that comprise each local submatrix
856: . col - an array of index sets that contain the global column numbers
857: that comprise each local submatrix
858: . submat - array of local submatrices
859: - ctx - optional user-defined context for private data for the
860: user-defined func routine (may be null)
862: Notes:
863: PCSetModifySubMatrices() MUST be called before KSPSetUp() and
864: KSPSolve().
866: A routine set by PCSetModifySubMatrices() is currently called within
867: the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
868: preconditioners. All other preconditioners ignore this routine.
870: Level: advanced
872: .keywords: PC, set, modify, submatrices
874: .seealso: PCModifySubMatrices()
875: @*/
876: PetscErrorCode PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
877: {
880: pc->modifysubmatrices = func;
881: pc->modifysubmatricesP = ctx;
882: return(0);
883: }
887: /*@C
888: PCModifySubMatrices - Calls an optional user-defined routine within
889: certain preconditioners if one has been set with PCSetModifySubMarices().
891: Collective on PC
893: Input Parameters:
894: + pc - the preconditioner context
895: . nsub - the number of local submatrices
896: . row - an array of index sets that contain the global row numbers
897: that comprise each local submatrix
898: . col - an array of index sets that contain the global column numbers
899: that comprise each local submatrix
900: . submat - array of local submatrices
901: - ctx - optional user-defined context for private data for the
902: user-defined routine (may be null)
904: Output Parameter:
905: . submat - array of local submatrices (the entries of which may
906: have been modified)
908: Notes:
909: The user should NOT generally call this routine, as it will
910: automatically be called within certain preconditioners (currently
911: block Jacobi, additive Schwarz) if set.
913: The basic submatrices are extracted from the preconditioner matrix
914: as usual; the user can then alter these (for example, to set different
915: boundary conditions for each submatrix) before they are used for the
916: local solves.
918: Level: developer
920: .keywords: PC, modify, submatrices
922: .seealso: PCSetModifySubMatrices()
923: @*/
924: PetscErrorCode PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
925: {
929: if (!pc->modifysubmatrices) return(0);
930: PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);
931: (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);
932: PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);
933: return(0);
934: }
938: /*@
939: PCSetOperators - Sets the matrix associated with the linear system and
940: a (possibly) different one associated with the preconditioner.
942: Collective on PC and Mat
944: Input Parameters:
945: + pc - the preconditioner context
946: . Amat - the matrix associated with the linear system
947: . Pmat - the matrix to be used in constructing the preconditioner, usually
948: the same as Amat.
949: - flag - flag indicating information about the preconditioner matrix structure
950: during successive linear solves. This flag is ignored the first time a
951: linear system is solved, and thus is irrelevant when solving just one linear
952: system.
954: Notes:
955: The flag can be used to eliminate unnecessary work in the preconditioner
956: during the repeated solution of linear systems of the same size. The
957: available options are
958: + SAME_PRECONDITIONER -
959: Pmat is identical during successive linear solves.
960: This option is intended for folks who are using
961: different Amat and Pmat matrices and wish to reuse the
962: same preconditioner matrix. For example, this option
963: saves work by not recomputing incomplete factorization
964: for ILU/ICC preconditioners.
965: . SAME_NONZERO_PATTERN -
966: Pmat has the same nonzero structure during
967: successive linear solves.
968: - DIFFERENT_NONZERO_PATTERN -
969: Pmat does not have the same nonzero structure.
971: Caution:
972: If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
973: and does not check the structure of the matrix. If you erroneously
974: claim that the structure is the same when it actually is not, the new
975: preconditioner will not function correctly. Thus, use this optimization
976: feature carefully!
978: If in doubt about whether your preconditioner matrix has changed
979: structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
981: More Notes about Repeated Solution of Linear Systems:
982: PETSc does NOT reset the matrix entries of either Amat or Pmat
983: to zero after a linear solve; the user is completely responsible for
984: matrix assembly. See the routine MatZeroEntries() if desiring to
985: zero all elements of a matrix.
987: Level: intermediate
989: .keywords: PC, set, operators, matrix, linear system
991: .seealso: PCGetOperators(), MatZeroEntries()
992: @*/
993: PetscErrorCode PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
994: {
996: PetscTruth isbjacobi,isrowbs;
1003: /*
1004: BlockSolve95 cannot use default BJacobi preconditioning
1005: */
1006: if (Amat) {
1007: PetscTypeCompare((PetscObject)Amat,MATMPIROWBS,&isrowbs);
1008: if (isrowbs) {
1009: PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);
1010: if (isbjacobi) {
1011: PCSetType(pc,PCILU);
1012: PetscLogInfo(pc,"PCSetOperators:Switching default PC to PCILU since BS95 doesn't support PCBJACOBI\n");
1013: }
1014: }
1015: }
1017: pc->mat = Amat;
1018: pc->pmat = Pmat;
1019: if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1020: pc->setupcalled = 1;
1021: }
1022: pc->flag = flag;
1023: return(0);
1024: }
1028: /*@C
1029: PCGetOperators - Gets the matrix associated with the linear system and
1030: possibly a different one associated with the preconditioner.
1032: Not collective, though parallel Mats are returned if the PC is parallel
1034: Input Parameter:
1035: . pc - the preconditioner context
1037: Output Parameters:
1038: + mat - the matrix associated with the linear system
1039: . pmat - matrix associated with the preconditioner, usually the same
1040: as mat.
1041: - flag - flag indicating information about the preconditioner
1042: matrix structure. See PCSetOperators() for details.
1044: Level: intermediate
1046: .keywords: PC, get, operators, matrix, linear system
1048: .seealso: PCSetOperators()
1049: @*/
1050: PetscErrorCode PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1051: {
1054: if (mat) *mat = pc->mat;
1055: if (pmat) *pmat = pc->pmat;
1056: if (flag) *flag = pc->flag;
1057: return(0);
1058: }
1062: /*@C
1063: PCGetFactoredMatrix - Gets the factored matrix from the
1064: preconditioner context. This routine is valid only for the LU,
1065: incomplete LU, Cholesky, and incomplete Cholesky methods.
1067: Not Collective on PC though Mat is parallel if PC is parallel
1069: Input Parameters:
1070: . pc - the preconditioner context
1072: Output parameters:
1073: . mat - the factored matrix
1075: Level: advanced
1077: Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1079: .keywords: PC, get, factored, matrix
1080: @*/
1081: PetscErrorCode PCGetFactoredMatrix(PC pc,Mat *mat)
1082: {
1088: if (pc->ops->getfactoredmatrix) {
1089: (*pc->ops->getfactoredmatrix)(pc,mat);
1090: }
1091: return(0);
1092: }
1096: /*@C
1097: PCSetOptionsPrefix - Sets the prefix used for searching for all
1098: PC options in the database.
1100: Collective on PC
1102: Input Parameters:
1103: + pc - the preconditioner context
1104: - prefix - the prefix string to prepend to all PC option requests
1106: Notes:
1107: A hyphen (-) must NOT be given at the beginning of the prefix name.
1108: The first character of all runtime options is AUTOMATICALLY the
1109: hyphen.
1111: Level: advanced
1113: .keywords: PC, set, options, prefix, database
1115: .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1116: @*/
1117: PetscErrorCode PCSetOptionsPrefix(PC pc,const char prefix[])
1118: {
1123: PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);
1124: return(0);
1125: }
1129: /*@C
1130: PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1131: PC options in the database.
1133: Collective on PC
1135: Input Parameters:
1136: + pc - the preconditioner context
1137: - prefix - the prefix string to prepend to all PC option requests
1139: Notes:
1140: A hyphen (-) must NOT be given at the beginning of the prefix name.
1141: The first character of all runtime options is AUTOMATICALLY the
1142: hyphen.
1144: Level: advanced
1146: .keywords: PC, append, options, prefix, database
1148: .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1149: @*/
1150: PetscErrorCode PCAppendOptionsPrefix(PC pc,const char prefix[])
1151: {
1156: PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);
1157: return(0);
1158: }
1162: /*@C
1163: PCGetOptionsPrefix - Gets the prefix used for searching for all
1164: PC options in the database.
1166: Not Collective
1168: Input Parameters:
1169: . pc - the preconditioner context
1171: Output Parameters:
1172: . prefix - pointer to the prefix string used, is returned
1174: Notes: On the fortran side, the user should pass in a string 'prifix' of
1175: sufficient length to hold the prefix.
1177: Level: advanced
1179: .keywords: PC, get, options, prefix, database
1181: .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1182: @*/
1183: PetscErrorCode PCGetOptionsPrefix(PC pc,char *prefix[])
1184: {
1190: PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);
1191: return(0);
1192: }
1196: /*@
1197: PCPreSolve - Optional pre-solve phase, intended for any
1198: preconditioner-specific actions that must be performed before
1199: the iterative solve itself.
1201: Collective on PC
1203: Input Parameters:
1204: + pc - the preconditioner context
1205: - ksp - the Krylov subspace context
1207: Level: developer
1209: Sample of Usage:
1210: .vb
1211: PCPreSolve(pc,ksp);
1212: KSPSolve(ksp,b,x);
1213: PCPostSolve(pc,ksp);
1214: .ve
1216: Notes:
1217: The pre-solve phase is distinct from the PCSetUp() phase.
1219: KSPSolve() calls this directly, so is rarely called by the user.
1221: .keywords: PC, pre-solve
1223: .seealso: PCPostSolve()
1224: @*/
1225: PetscErrorCode PCPreSolve(PC pc,KSP ksp)
1226: {
1228: Vec x,rhs;
1229: Mat A,B;
1234: KSPGetSolution(ksp,&x);
1235: KSPGetRhs(ksp,&rhs);
1236: /*
1237: Scale the system and have the matrices use the scaled form
1238: only if the two matrices are actually the same (and hence
1239: have the same scaling
1240: */
1241: PCGetOperators(pc,&A,&B,PETSC_NULL);
1242: if (A == B) {
1243: MatScaleSystem(pc->mat,x,rhs);
1244: MatUseScaledForm(pc->mat,PETSC_TRUE);
1245: }
1247: if (pc->ops->presolve) {
1248: (*pc->ops->presolve)(pc,ksp,x,rhs);
1249: }
1250: return(0);
1251: }
1255: /*@
1256: PCPostSolve - Optional post-solve phase, intended for any
1257: preconditioner-specific actions that must be performed after
1258: the iterative solve itself.
1260: Collective on PC
1262: Input Parameters:
1263: + pc - the preconditioner context
1264: - ksp - the Krylov subspace context
1266: Sample of Usage:
1267: .vb
1268: PCPreSolve(pc,ksp);
1269: KSPSolve(ksp,b,x);
1270: PCPostSolve(pc,ksp);
1271: .ve
1273: Note:
1274: KSPSolve() calls this routine directly, so it is rarely called by the user.
1276: Level: developer
1278: .keywords: PC, post-solve
1280: .seealso: PCPreSolve(), KSPSolve()
1281: @*/
1282: PetscErrorCode PCPostSolve(PC pc,KSP ksp)
1283: {
1285: Vec x,rhs;
1286: Mat A,B;
1291: KSPGetSolution(ksp,&x);
1292: KSPGetRhs(ksp,&rhs);
1293: if (pc->ops->postsolve) {
1294: (*pc->ops->postsolve)(pc,ksp,x,rhs);
1295: }
1297: /*
1298: Scale the system and have the matrices use the scaled form
1299: only if the two matrices are actually the same (and hence
1300: have the same scaling
1301: */
1302: PCGetOperators(pc,&A,&B,PETSC_NULL);
1303: if (A == B) {
1304: MatUnScaleSystem(pc->mat,x,rhs);
1305: MatUseScaledForm(pc->mat,PETSC_FALSE);
1306: }
1307: return(0);
1308: }
1312: /*@C
1313: PCView - Prints the PC data structure.
1315: Collective on PC
1317: Input Parameters:
1318: + PC - the PC context
1319: - viewer - optional visualization context
1321: Note:
1322: The available visualization contexts include
1323: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
1324: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1325: output where only the first processor opens
1326: the file. All other processors send their
1327: data to the first processor to print.
1329: The user can open an alternative visualization contexts with
1330: PetscViewerASCIIOpen() (output to a specified file).
1332: Level: developer
1334: .keywords: PC, view
1336: .seealso: KSPView(), PetscViewerASCIIOpen()
1337: @*/
1338: PetscErrorCode PCView(PC pc,PetscViewer viewer)
1339: {
1340: PCType cstr;
1341: PetscErrorCode ierr;
1342: PetscTruth mat_exists,iascii,isstring;
1343: PetscViewerFormat format;
1347: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(pc->comm);
1351: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
1352: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
1353: if (iascii) {
1354: PetscViewerGetFormat(viewer,&format);
1355: if (pc->prefix) {
1356: PetscViewerASCIIPrintf(viewer,"PC Object:(%s)\n",pc->prefix);
1357: } else {
1358: PetscViewerASCIIPrintf(viewer,"PC Object:\n");
1359: }
1360: PCGetType(pc,&cstr);
1361: if (cstr) {
1362: PetscViewerASCIIPrintf(viewer," type: %s\n",cstr);
1363: } else {
1364: PetscViewerASCIIPrintf(viewer," type: not yet set\n");
1365: }
1366: if (pc->ops->view) {
1367: PetscViewerASCIIPushTab(viewer);
1368: (*pc->ops->view)(pc,viewer);
1369: PetscViewerASCIIPopTab(viewer);
1370: }
1371: PetscObjectExists((PetscObject)pc->mat,&mat_exists);
1372: if (mat_exists) {
1373: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1374: if (pc->pmat == pc->mat) {
1375: PetscViewerASCIIPrintf(viewer," linear system matrix = precond matrix:\n");
1376: PetscViewerASCIIPushTab(viewer);
1377: MatView(pc->mat,viewer);
1378: PetscViewerASCIIPopTab(viewer);
1379: } else {
1380: PetscObjectExists((PetscObject)pc->pmat,&mat_exists);
1381: if (mat_exists) {
1382: PetscViewerASCIIPrintf(viewer," linear system matrix followed by preconditioner matrix:\n");
1383: } else {
1384: PetscViewerASCIIPrintf(viewer," linear system matrix:\n");
1385: }
1386: PetscViewerASCIIPushTab(viewer);
1387: MatView(pc->mat,viewer);
1388: if (mat_exists) {MatView(pc->pmat,viewer);}
1389: PetscViewerASCIIPopTab(viewer);
1390: }
1391: PetscViewerPopFormat(viewer);
1392: }
1393: } else if (isstring) {
1394: PCGetType(pc,&cstr);
1395: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
1396: if (pc->ops->view) {(*pc->ops->view)(pc,viewer);}
1397: } else {
1398: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1399: }
1400: return(0);
1401: }
1405: /*@C
1406: PCRegister - See PCRegisterDynamic()
1408: Level: advanced
1409: @*/
1410: PetscErrorCode PCRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PC))
1411: {
1413: char fullname[PETSC_MAX_PATH_LEN];
1417: PetscFListConcat(path,name,fullname);
1418: PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);
1419: return(0);
1420: }
1424: /*@
1425: PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1427: Collective on PC
1429: Input Parameter:
1430: . pc - the preconditioner object
1432: Output Parameter:
1433: . mat - the explict preconditioned operator
1435: Notes:
1436: This computation is done by applying the operators to columns of the
1437: identity matrix.
1439: Currently, this routine uses a dense matrix format when 1 processor
1440: is used and a sparse format otherwise. This routine is costly in general,
1441: and is recommended for use only with relatively small systems.
1443: Level: advanced
1444:
1445: .keywords: PC, compute, explicit, operator
1447: @*/
1448: PetscErrorCode PCComputeExplicitOperator(PC pc,Mat *mat)
1449: {
1450: Vec in,out;
1452: PetscInt i,M,m,*rows,start,end;
1453: PetscMPIInt size;
1454: MPI_Comm comm;
1455: PetscScalar *array,zero = 0.0,one = 1.0;
1461: comm = pc->comm;
1462: MPI_Comm_size(comm,&size);
1464: if (!pc->pmat) SETERRQ(PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1465: MatGetVecs(pc->pmat,&in,0);
1466: VecDuplicate(in,&out);
1467: VecGetOwnershipRange(in,&start,&end);
1468: VecGetSize(in,&M);
1469: VecGetLocalSize(in,&m);
1470: PetscMalloc((m+1)*sizeof(PetscInt),&rows);
1471: for (i=0; i<m; i++) {rows[i] = start + i;}
1473: MatCreate(comm,m,m,M,M,mat);
1474: if (size == 1) {
1475: MatSetType(*mat,MATSEQDENSE);
1476: MatSeqDenseSetPreallocation(*mat,PETSC_NULL);
1477: } else {
1478: MatSetType(*mat,MATMPIAIJ);
1479: MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);
1480: }
1482: for (i=0; i<M; i++) {
1484: VecSet(&zero,in);
1485: VecSetValues(in,1,&i,&one,INSERT_VALUES);
1486: VecAssemblyBegin(in);
1487: VecAssemblyEnd(in);
1489: /* should fix, allowing user to choose side */
1490: PCApply(pc,in,out);
1491:
1492: VecGetArray(out,&array);
1493: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
1494: VecRestoreArray(out,&array);
1496: }
1497: PetscFree(rows);
1498: VecDestroy(out);
1499: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
1500: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
1501: return(0);
1502: }