Actual source code: precon.c
1: /*$Id: precon.c,v 1.216 2001/08/21 21:03:13 bsmith Exp $*/
2: /*
3: The PC (preconditioner) interface routines, callable by users.
4: */
5: #include src/sles/pc/pcimpl.h
7: /* Logging support */
8: int PC_COOKIE;
9: int PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
10: int PC_ApplySymmetricRight, PC_ModifySubMatrices;
12: EXTERN int SLESInitializePackage(const char[]);
16: /*@C
17: PCNullSpaceAttach - attaches a null space to a preconditioner object.
18: This null space will be removed from the resulting vector whenever
19: the preconditioner is applied.
21: Collective on PC
23: Input Parameters:
24: + pc - the preconditioner context
25: - nullsp - the null space object
27: Level: developer
29: Notes:
30: Overwrites any previous null space that may have been attached
32: Users manual sections:
33: . Section 4.15 Solving Singular Systems
35: .keywords: PC, destroy, null space
37: .seealso: PCCreate(), PCSetUp(), MatNullSpaceCreate(), MatNullSpace
39: @*/
40: int PCNullSpaceAttach(PC pc,MatNullSpace nullsp)
41: {
48: if (pc->nullsp) {
49: MatNullSpaceDestroy(pc->nullsp);
50: }
51: pc->nullsp = nullsp;
52: PetscObjectReference((PetscObject)nullsp);
53: return(0);
54: }
58: /*@C
59: PCDestroy - Destroys PC context that was created with PCCreate().
61: Collective on PC
63: Input Parameter:
64: . pc - the preconditioner context
66: Level: developer
68: .keywords: PC, destroy
70: .seealso: PCCreate(), PCSetUp()
71: @*/
72: int PCDestroy(PC pc)
73: {
78: if (--pc->refct > 0) return(0);
80: /* if memory was published with AMS then destroy it */
81: PetscObjectDepublish(pc);
83: if (pc->ops->destroy) { (*pc->ops->destroy)(pc);}
84: if (pc->nullsp) {MatNullSpaceDestroy(pc->nullsp);}
85: if (pc->diagonalscaleright) {VecDestroy(pc->diagonalscaleright);}
86: if (pc->diagonalscaleleft) {VecDestroy(pc->diagonalscaleleft);}
88: PetscLogObjectDestroy(pc);
89: PetscHeaderDestroy(pc);
90: return(0);
91: }
95: /*@C
96: PCDiagonalScale - Indicates if the preconditioner applies an additional left and right
97: scaling as needed by certain time-stepping codes.
99: Collective on PC
101: Input Parameter:
102: . pc - the preconditioner context
104: Output Parameter:
105: . flag - PETSC_TRUE if it applies the scaling
107: Level: developer
109: Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
110: $ D M A D^{-1} y = D M b for left preconditioning or
111: $ D A M D^{-1} z = D b for right preconditioning
113: .keywords: PC
115: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScaleSet()
116: @*/
117: int PCDiagonalScale(PC pc,PetscTruth *flag)
118: {
121: *flag = pc->diagonalscale;
122: return(0);
123: }
127: /*@
128: PCDiagonalScaleSet - Indicates the left scaling to use to apply an additional left and right
129: scaling as needed by certain time-stepping codes.
131: Collective on PC
133: Input Parameters:
134: + pc - the preconditioner context
135: - s - scaling vector
137: Level: intermediate
139: Notes: The system solved via the Krylov method is
140: $ D M A D^{-1} y = D M b for left preconditioning or
141: $ D A M D^{-1} z = D b for right preconditioning
143: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
145: .keywords: PC
147: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScale()
148: @*/
149: int PCDiagonalScaleSet(PC pc,Vec s)
150: {
156: pc->diagonalscale = PETSC_TRUE;
157: if (pc->diagonalscaleleft) {
158: VecDestroy(pc->diagonalscaleleft);
159: }
160: pc->diagonalscaleleft = s;
161: PetscObjectReference((PetscObject)s);
162: if (!pc->diagonalscaleright) {
163: VecDuplicate(s,&pc->diagonalscaleright);
164: }
165: VecCopy(s,pc->diagonalscaleright);
166: VecReciprocal(pc->diagonalscaleright);
167: return(0);
168: }
172: /*@C
173: PCDiagonalScaleLeft - Indicates the left scaling to use to apply an additional left and right
174: scaling as needed by certain time-stepping codes.
176: Collective on PC
178: Input Parameters:
179: + pc - the preconditioner context
180: . in - input vector
181: + out - scaled vector (maybe the same as in)
183: Level: intermediate
185: Notes: The system solved via the Krylov method is
186: $ D M A D^{-1} y = D M b for left preconditioning or
187: $ D A M D^{-1} z = D b for right preconditioning
189: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
191: If diagonal scaling is turned off and in is not out then in is copied to out
193: .keywords: PC
195: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
196: @*/
197: int PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
198: {
205: if (pc->diagonalscale) {
206: VecPointwiseMult(pc->diagonalscaleleft,in,out);
207: } else if (in != out) {
208: VecCopy(in,out);
209: }
210: return(0);
211: }
215: /*@C
216: PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
218: Collective on PC
220: Input Parameters:
221: + pc - the preconditioner context
222: . in - input vector
223: + out - scaled vector (maybe the same as in)
225: Level: intermediate
227: Notes: The system solved via the Krylov method is
228: $ D M A D^{-1} y = D M b for left preconditioning or
229: $ D A M D^{-1} z = D b for right preconditioning
231: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
233: If diagonal scaling is turned off and in is not out then in is copied to out
235: .keywords: PC
237: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
238: @*/
239: int PCDiagonalScaleRight(PC pc,Vec in,Vec out)
240: {
247: if (pc->diagonalscale) {
248: VecPointwiseMult(pc->diagonalscaleright,in,out);
249: } else if (in != out) {
250: VecCopy(in,out);
251: }
252: return(0);
253: }
257: static int PCPublish_Petsc(PetscObject obj)
258: {
259: #if defined(PETSC_HAVE_AMS)
260: PC v = (PC) obj;
261: int ierr;
262: #endif
266: #if defined(PETSC_HAVE_AMS)
267: /* if it is already published then return */
268: if (v->amem >=0) return(0);
270: PetscObjectPublishBaseBegin(obj);
271: PetscObjectPublishBaseEnd(obj);
272: #endif
274: return(0);
275: }
279: /*@C
280: PCCreate - Creates a preconditioner context.
282: Collective on MPI_Comm
284: Input Parameter:
285: . comm - MPI communicator
287: Output Parameter:
288: . pc - location to put the preconditioner context
290: Notes:
291: The default preconditioner on one processor is PCILU with 0 fill on more
292: then one it is PCBJACOBI with ILU() on each processor.
294: Level: developer
296: .keywords: PC, create, context
298: .seealso: PCSetUp(), PCApply(), PCDestroy()
299: @*/
300: int PCCreate(MPI_Comm comm,PC *newpc)
301: {
302: PC pc;
307: *newpc = 0;
308: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
309: SLESInitializePackage(PETSC_NULL);
310: #endif
312: PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_COOKIE,-1,"PC",comm,PCDestroy,PCView);
313: PetscLogObjectCreate(pc);
314: pc->bops->publish = PCPublish_Petsc;
315: pc->vec = 0;
316: pc->mat = 0;
317: pc->setupcalled = 0;
318: pc->nullsp = 0;
319: pc->data = 0;
320: pc->diagonalscale = PETSC_FALSE;
321: pc->diagonalscaleleft = 0;
322: pc->diagonalscaleright = 0;
324: pc->ops->destroy = 0;
325: pc->ops->apply = 0;
326: pc->ops->applytranspose = 0;
327: pc->ops->applyBA = 0;
328: pc->ops->applyBAtranspose = 0;
329: pc->ops->applyrichardson = 0;
330: pc->ops->view = 0;
331: pc->ops->getfactoredmatrix = 0;
332: pc->ops->applysymmetricright = 0;
333: pc->ops->applysymmetricleft = 0;
334: pc->ops->setuponblocks = 0;
336: pc->modifysubmatrices = 0;
337: pc->modifysubmatricesP = 0;
338: *newpc = pc;
339: PetscPublishAll(pc);
340: return(0);
342: }
344: /* -------------------------------------------------------------------------------*/
348: /*@
349: PCApply - Applies the preconditioner to a vector.
351: Collective on PC and Vec
353: Input Parameters:
354: + pc - the preconditioner context
355: - x - input vector
357: Output Parameter:
358: . y - output vector
360: Level: developer
362: .keywords: PC, apply
364: .seealso: PCApplyTranspose(), PCApplyBAorAB()
365: @*/
366: int PCApply(PC pc,Vec x,Vec y,PCSide side)
367: {
368: int ierr;
374: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
376: if (pc->setupcalled < 2) {
377: PCSetUp(pc);
378: }
380: /* Remove null space from input vector y */
381: if (side == PC_RIGHT && pc->nullsp) {
382: Vec tmp;
383: MatNullSpaceRemove(pc->nullsp,x,&tmp);
384: x = tmp;
385: SETERRQ(1,"Cannot deflate out null space if using right preconditioner!");
386: }
388: PetscLogEventBegin(PC_Apply,pc,x,y,0);
389: (*pc->ops->apply)(pc,x,y);
390: PetscLogEventEnd(PC_Apply,pc,x,y,0);
392: /* Remove null space from preconditioned vector y */
393: if (side == PC_LEFT && pc->nullsp) {
394: MatNullSpaceRemove(pc->nullsp,y,PETSC_NULL);
395: }
396: return(0);
397: }
401: /*@
402: PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
404: Collective on PC and Vec
406: Input Parameters:
407: + pc - the preconditioner context
408: - x - input vector
410: Output Parameter:
411: . y - output vector
413: Notes:
414: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
416: Level: developer
418: .keywords: PC, apply, symmetric, left
420: .seealso: PCApply(), PCApplySymmetricRight()
421: @*/
422: int PCApplySymmetricLeft(PC pc,Vec x,Vec y)
423: {
430: if (!pc->ops->applysymmetricleft) SETERRQ(1,"PC does not have left symmetric apply");
432: if (pc->setupcalled < 2) {
433: PCSetUp(pc);
434: }
436: PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);
437: (*pc->ops->applysymmetricleft)(pc,x,y);
438: PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);
439: return(0);
440: }
444: /*@
445: PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
447: Collective on PC and Vec
449: Input Parameters:
450: + pc - the preconditioner context
451: - x - input vector
453: Output Parameter:
454: . y - output vector
456: Level: developer
458: Notes:
459: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
461: .keywords: PC, apply, symmetric, right
463: .seealso: PCApply(), PCApplySymmetricLeft()
464: @*/
465: int PCApplySymmetricRight(PC pc,Vec x,Vec y)
466: {
473: if (!pc->ops->applysymmetricright) SETERRQ(1,"PC does not have left symmetric apply");
475: if (pc->setupcalled < 2) {
476: PCSetUp(pc);
477: }
479: PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);
480: (*pc->ops->applysymmetricright)(pc,x,y);
481: PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);
482: return(0);
483: }
487: /*@
488: PCApplyTranspose - Applies the transpose of preconditioner to a vector.
490: Collective on PC and Vec
492: Input Parameters:
493: + pc - the preconditioner context
494: - x - input vector
496: Output Parameter:
497: . y - output vector
499: Level: developer
501: .keywords: PC, apply, transpose
503: .seealso: PCApplyTranspose(), PCApplyBAorAB(), PCApplyBAorABTranspose()
504: @*/
505: int PCApplyTranspose(PC pc,Vec x,Vec y)
506: {
513: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
514: if (!pc->ops->applytranspose) SETERRQ(PETSC_ERR_SUP," ");
516: if (pc->setupcalled < 2) {
517: PCSetUp(pc);
518: }
520: PetscLogEventBegin(PC_Apply,pc,x,y,0);
521: (*pc->ops->applytranspose)(pc,x,y);
522: PetscLogEventEnd(PC_Apply,pc,x,y,0);
523: return(0);
524: }
528: /*@
529: PCApplyBAorAB - Applies the preconditioner and operator to a vector.
531: Collective on PC and Vec
533: Input Parameters:
534: + pc - the preconditioner context
535: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
536: . x - input vector
537: - work - work vector
539: Output Parameter:
540: . y - output vector
542: Level: developer
544: .keywords: PC, apply, operator
546: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
547: @*/
548: int PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
549: {
550: int ierr;
557: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
558: if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) {
559: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric");
560: }
561: if (pc->diagonalscale && side == PC_SYMMETRIC) {
562: SETERRQ(1,"Cannot include diagonal scaling with symmetric preconditioner application");
563: }
565: if (pc->setupcalled < 2) {
566: PCSetUp(pc);
567: }
569: if (pc->diagonalscale) {
570: if (pc->ops->applyBA) {
571: Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
572: VecDuplicate(x,&work2);
573: PCDiagonalScaleRight(pc,x,work2);
574: (*pc->ops->applyBA)(pc,side,work2,y,work);
575: /* Remove null space from preconditioned vector y */
576: if (pc->nullsp) {
577: MatNullSpaceRemove(pc->nullsp,y,PETSC_NULL);
578: }
579: PCDiagonalScaleLeft(pc,y,y);
580: VecDestroy(work2);
581: } else if (side == PC_RIGHT) {
582: PCDiagonalScaleRight(pc,x,y);
583: PCApply(pc,y,work,side);
584: MatMult(pc->mat,work,y);
585: PCDiagonalScaleLeft(pc,y,y);
586: } else if (side == PC_LEFT) {
587: PCDiagonalScaleRight(pc,x,y);
588: MatMult(pc->mat,y,work);
589: PCApply(pc,work,y,side);
590: PCDiagonalScaleLeft(pc,y,y);
591: } else if (side == PC_SYMMETRIC) {
592: SETERRQ(1,"Cannot provide diagonal scaling with symmetric application of preconditioner");
593: }
594: } else {
595: if (pc->ops->applyBA) {
596: (*pc->ops->applyBA)(pc,side,x,y,work);
597: /* Remove null space from preconditioned vector y */
598: if (pc->nullsp) {
599: MatNullSpaceRemove(pc->nullsp,y,PETSC_NULL);
600: }
601: } else if (side == PC_RIGHT) {
602: PCApply(pc,x,work,side);
603: MatMult(pc->mat,work,y);
604: } else if (side == PC_LEFT) {
605: MatMult(pc->mat,x,work);
606: PCApply(pc,work,y,side);
607: } else if (side == PC_SYMMETRIC) {
608: /* There's an extra copy here; maybe should provide 2 work vectors instead? */
609: PCApplySymmetricRight(pc,x,work);
610: MatMult(pc->mat,work,y);
611: VecCopy(y,work);
612: PCApplySymmetricLeft(pc,work,y);
613: }
614: }
615: return(0);
616: }
620: /*@
621: PCApplyBAorABTranspose - Applies the transpose of the preconditioner
622: and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
623: not tr(B*A) = tr(A)*tr(B).
625: Collective on PC and Vec
627: Input Parameters:
628: + pc - the preconditioner context
629: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
630: . x - input vector
631: - work - work vector
633: Output Parameter:
634: . y - output vector
636: Level: developer
638: .keywords: PC, apply, operator, transpose
640: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
641: @*/
642: int PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
643: {
651: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
652: if (pc->ops->applyBAtranspose) {
653: (*pc->ops->applyBAtranspose)(pc,side,x,y,work);
654: return(0);
655: }
656: if (side != PC_LEFT && side != PC_RIGHT) {
657: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
658: }
660: if (pc->setupcalled < 2) {
661: PCSetUp(pc);
662: }
664: if (side == PC_RIGHT) {
665: PCApplyTranspose(pc,x,work);
666: MatMultTranspose(pc->mat,work,y);
667: } else if (side == PC_LEFT) {
668: MatMultTranspose(pc->mat,x,work);
669: PCApplyTranspose(pc,work,y);
670: }
671: /* add support for PC_SYMMETRIC */
672: return(0); /* actually will never get here */
673: }
675: /* -------------------------------------------------------------------------------*/
679: /*@
680: PCApplyRichardsonExists - Determines whether a particular preconditioner has a
681: built-in fast application of Richardson's method.
683: Not Collective
685: Input Parameter:
686: . pc - the preconditioner
688: Output Parameter:
689: . exists - PETSC_TRUE or PETSC_FALSE
691: Level: developer
693: .keywords: PC, apply, Richardson, exists
695: .seealso: PCApplyRichardson()
696: @*/
697: int PCApplyRichardsonExists(PC pc,PetscTruth *exists)
698: {
702: if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
703: else *exists = PETSC_FALSE;
704: return(0);
705: }
709: /*@
710: PCApplyRichardson - Applies several steps of Richardson iteration with
711: the particular preconditioner. This routine is usually used by the
712: Krylov solvers and not the application code directly.
714: Collective on PC
716: Input Parameters:
717: + pc - the preconditioner context
718: . x - the initial guess
719: . w - one work vector
720: . rtol - relative decrease in residual norm convergence criteria
721: . atol - absolute residual norm convergence criteria
722: . dtol - divergence residual norm increase criteria
723: - its - the number of iterations to apply.
725: Output Parameter:
726: . y - the solution
728: Notes:
729: Most preconditioners do not support this function. Use the command
730: PCApplyRichardsonExists() to determine if one does.
732: Except for the multigrid PC this routine ignores the convergence tolerances
733: and always runs for the number of iterations
734:
735: Level: developer
737: .keywords: PC, apply, Richardson
739: .seealso: PCApplyRichardsonExists()
740: @*/
741: int PCApplyRichardson(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal atol, PetscReal dtol,int its)
742: {
750: if (!pc->ops->applyrichardson) SETERRQ(PETSC_ERR_SUP," ");
752: if (pc->setupcalled < 2) {
753: PCSetUp(pc);
754: }
756: (*pc->ops->applyrichardson)(pc,x,y,w,rtol,atol,dtol,its);
757: return(0);
758: }
760: /*
761: a setupcall of 0 indicates never setup,
762: 1 needs to be resetup,
763: 2 does not need any changes.
764: */
767: /*@
768: PCSetUp - Prepares for the use of a preconditioner.
770: Collective on PC
772: Input Parameter:
773: . pc - the preconditioner context
775: Level: developer
777: .keywords: PC, setup
779: .seealso: PCCreate(), PCApply(), PCDestroy()
780: @*/
781: int PCSetUp(PC pc)
782: {
784: PetscTruth flg;
789: if (pc->setupcalled > 1) {
790: PetscLogInfo(pc,"PCSetUp:Setting PC with identical preconditioner\n");
791: return(0);
792: } else if (pc->setupcalled == 0) {
793: PetscLogInfo(pc,"PCSetUp:Setting up new PC\n");
794: } else if (pc->flag == SAME_NONZERO_PATTERN) {
795: PetscLogInfo(pc,"PCSetUp:Setting up PC with same nonzero pattern\n");
796: } else {
797: PetscLogInfo(pc,"PCSetUp:Setting up PC with different nonzero pattern\n");
798: }
800: PetscLogEventBegin(PC_SetUp,pc,0,0,0);
801: if (!pc->vec) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Vector must be set first");}
802: if (!pc->mat) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
804: if (!pc->type_name) {
805: int size;
807: MPI_Comm_size(pc->comm,&size);
808: if (size == 1) {
809: MatHasOperation(pc->pmat,MATOP_ICCFACTOR_SYMBOLIC,&flg);
810: if (flg) { /* for sbaij mat */
811: PCSetType(pc,PCICC);
812: } else {
813: PCSetType(pc,PCILU);
814: }
815: } else {
816: PCSetType(pc,PCBJACOBI);
817: }
818: }
819: if (pc->ops->setup) {
820: (*pc->ops->setup)(pc);
821: }
822: pc->setupcalled = 2;
823: if (pc->nullsp) {
824: PetscTruth test;
825: PetscOptionsHasName(pc->prefix,"-pc_test_null_space",&test);
826: if (test) {
827: MatNullSpaceTest(pc->nullsp,pc->mat);
828: }
829: }
830: PetscLogEventEnd(PC_SetUp,pc,0,0,0);
831: return(0);
832: }
836: /*@
837: PCSetUpOnBlocks - Sets up the preconditioner for each block in
838: the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
839: methods.
841: Collective on PC
843: Input Parameters:
844: . pc - the preconditioner context
846: Level: developer
848: .keywords: PC, setup, blocks
850: .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
851: @*/
852: int PCSetUpOnBlocks(PC pc)
853: {
858: if (!pc->ops->setuponblocks) return(0);
859: PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);
860: (*pc->ops->setuponblocks)(pc);
861: PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);
862: return(0);
863: }
867: /*@C
868: PCSetModifySubMatrices - Sets a user-defined routine for modifying the
869: submatrices that arise within certain subdomain-based preconditioners.
870: The basic submatrices are extracted from the preconditioner matrix as
871: usual; the user can then alter these (for example, to set different boundary
872: conditions for each submatrix) before they are used for the local solves.
874: Collective on PC
876: Input Parameters:
877: + pc - the preconditioner context
878: . func - routine for modifying the submatrices
879: - ctx - optional user-defined context (may be null)
881: Calling sequence of func:
882: $ func (PC pc,int nsub,IS *row,IS *col,Mat *submat,void *ctx);
884: . row - an array of index sets that contain the global row numbers
885: that comprise each local submatrix
886: . col - an array of index sets that contain the global column numbers
887: that comprise each local submatrix
888: . submat - array of local submatrices
889: - ctx - optional user-defined context for private data for the
890: user-defined func routine (may be null)
892: Notes:
893: PCSetModifySubMatrices() MUST be called before SLESSetUp() and
894: SLESSolve().
896: A routine set by PCSetModifySubMatrices() is currently called within
897: the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
898: preconditioners. All other preconditioners ignore this routine.
900: Level: advanced
902: .keywords: PC, set, modify, submatrices
904: .seealso: PCModifySubMatrices()
905: @*/
906: int PCSetModifySubMatrices(PC pc,int(*func)(PC,int,const IS[],const IS[],Mat[],void*),void *ctx)
907: {
910: pc->modifysubmatrices = func;
911: pc->modifysubmatricesP = ctx;
912: return(0);
913: }
917: /*@C
918: PCModifySubMatrices - Calls an optional user-defined routine within
919: certain preconditioners if one has been set with PCSetModifySubMarices().
921: Collective on PC
923: Input Parameters:
924: + pc - the preconditioner context
925: . nsub - the number of local submatrices
926: . row - an array of index sets that contain the global row numbers
927: that comprise each local submatrix
928: . col - an array of index sets that contain the global column numbers
929: that comprise each local submatrix
930: . submat - array of local submatrices
931: - ctx - optional user-defined context for private data for the
932: user-defined routine (may be null)
934: Output Parameter:
935: . submat - array of local submatrices (the entries of which may
936: have been modified)
938: Notes:
939: The user should NOT generally call this routine, as it will
940: automatically be called within certain preconditioners (currently
941: block Jacobi, additive Schwarz) if set.
943: The basic submatrices are extracted from the preconditioner matrix
944: as usual; the user can then alter these (for example, to set different
945: boundary conditions for each submatrix) before they are used for the
946: local solves.
948: Level: developer
950: .keywords: PC, modify, submatrices
952: .seealso: PCSetModifySubMatrices()
953: @*/
954: int PCModifySubMatrices(PC pc,int nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
955: {
959: if (!pc->modifysubmatrices) return(0);
960: PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);
961: (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);
962: PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);
963: return(0);
964: }
968: /*@
969: PCSetOperators - Sets the matrix associated with the linear system and
970: a (possibly) different one associated with the preconditioner.
972: Collective on PC and Mat
974: Input Parameters:
975: + pc - the preconditioner context
976: . Amat - the matrix associated with the linear system
977: . Pmat - the matrix to be used in constructing the preconditioner, usually
978: the same as Amat.
979: - flag - flag indicating information about the preconditioner matrix structure
980: during successive linear solves. This flag is ignored the first time a
981: linear system is solved, and thus is irrelevant when solving just one linear
982: system.
984: Notes:
985: The flag can be used to eliminate unnecessary work in the preconditioner
986: during the repeated solution of linear systems of the same size. The
987: available options are
988: + SAME_PRECONDITIONER -
989: Pmat is identical during successive linear solves.
990: This option is intended for folks who are using
991: different Amat and Pmat matrices and wish to reuse the
992: same preconditioner matrix. For example, this option
993: saves work by not recomputing incomplete factorization
994: for ILU/ICC preconditioners.
995: . SAME_NONZERO_PATTERN -
996: Pmat has the same nonzero structure during
997: successive linear solves.
998: - DIFFERENT_NONZERO_PATTERN -
999: Pmat does not have the same nonzero structure.
1001: Caution:
1002: If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
1003: and does not check the structure of the matrix. If you erroneously
1004: claim that the structure is the same when it actually is not, the new
1005: preconditioner will not function correctly. Thus, use this optimization
1006: feature carefully!
1008: If in doubt about whether your preconditioner matrix has changed
1009: structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
1011: More Notes about Repeated Solution of Linear Systems:
1012: PETSc does NOT reset the matrix entries of either Amat or Pmat
1013: to zero after a linear solve; the user is completely responsible for
1014: matrix assembly. See the routine MatZeroEntries() if desiring to
1015: zero all elements of a matrix.
1017: Level: intermediate
1019: .keywords: PC, set, operators, matrix, linear system
1021: .seealso: PCGetOperators(), MatZeroEntries()
1022: @*/
1023: int PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
1024: {
1025: int ierr;
1026: PetscTruth isbjacobi,isrowbs;
1033: /*
1034: BlockSolve95 cannot use default BJacobi preconditioning
1035: */
1036: PetscTypeCompare((PetscObject)Amat,MATMPIROWBS,&isrowbs);
1037: if (isrowbs) {
1038: PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);
1039: if (isbjacobi) {
1040: PCSetType(pc,PCILU);
1041: PetscLogInfo(pc,"PCSetOperators:Switching default PC to PCILU since BS95 doesn't support PCBJACOBI\n");
1042: }
1043: }
1045: pc->mat = Amat;
1046: pc->pmat = Pmat;
1047: if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1048: pc->setupcalled = 1;
1049: }
1050: pc->flag = flag;
1051: return(0);
1052: }
1056: /*@C
1057: PCGetOperators - Gets the matrix associated with the linear system and
1058: possibly a different one associated with the preconditioner.
1060: Not collective, though parallel Mats are returned if the PC is parallel
1062: Input Parameter:
1063: . pc - the preconditioner context
1065: Output Parameters:
1066: + mat - the matrix associated with the linear system
1067: . pmat - matrix associated with the preconditioner, usually the same
1068: as mat.
1069: - flag - flag indicating information about the preconditioner
1070: matrix structure. See PCSetOperators() for details.
1072: Level: intermediate
1074: .keywords: PC, get, operators, matrix, linear system
1076: .seealso: PCSetOperators()
1077: @*/
1078: int PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1079: {
1082: if (mat) *mat = pc->mat;
1083: if (pmat) *pmat = pc->pmat;
1084: if (flag) *flag = pc->flag;
1085: return(0);
1086: }
1090: /*@
1091: PCSetVector - Sets a vector associated with the preconditioner.
1093: Collective on PC and Vec
1095: Input Parameters:
1096: + pc - the preconditioner context
1097: - vec - the vector
1099: Notes:
1100: The vector must be set so that the preconditioner knows what type
1101: of vector to allocate if necessary.
1103: Level: intermediate
1105: .keywords: PC, set, vector
1107: .seealso: PCGetVector()
1109: @*/
1110: int PCSetVector(PC pc,Vec vec)
1111: {
1116: pc->vec = vec;
1117: return(0);
1118: }
1122: /*@
1123: PCGetVector - Gets a vector associated with the preconditioner; if the
1124: vector was not get set it will return a 0 pointer.
1126: Not collective, but vector is shared by all processors that share the PC
1128: Input Parameter:
1129: . pc - the preconditioner context
1131: Output Parameter:
1132: . vec - the vector
1134: Level: intermediate
1136: .keywords: PC, get, vector
1138: .seealso: PCSetVector()
1140: @*/
1141: int PCGetVector(PC pc,Vec *vec)
1142: {
1145: *vec = pc->vec;
1146: return(0);
1147: }
1151: /*@C
1152: PCGetFactoredMatrix - Gets the factored matrix from the
1153: preconditioner context. This routine is valid only for the LU,
1154: incomplete LU, Cholesky, and incomplete Cholesky methods.
1156: Not Collective on PC though Mat is parallel if PC is parallel
1158: Input Parameters:
1159: . pc - the preconditioner context
1161: Output parameters:
1162: . mat - the factored matrix
1164: Level: advanced
1166: .keywords: PC, get, factored, matrix
1167: @*/
1168: int PCGetFactoredMatrix(PC pc,Mat *mat)
1169: {
1174: if (pc->ops->getfactoredmatrix) {
1175: (*pc->ops->getfactoredmatrix)(pc,mat);
1176: }
1177: return(0);
1178: }
1182: /*@C
1183: PCSetOptionsPrefix - Sets the prefix used for searching for all
1184: PC options in the database.
1186: Collective on PC
1188: Input Parameters:
1189: + pc - the preconditioner context
1190: - prefix - the prefix string to prepend to all PC option requests
1192: Notes:
1193: A hyphen (-) must NOT be given at the beginning of the prefix name.
1194: The first character of all runtime options is AUTOMATICALLY the
1195: hyphen.
1197: Level: advanced
1199: .keywords: PC, set, options, prefix, database
1201: .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1202: @*/
1203: int PCSetOptionsPrefix(PC pc,const char prefix[])
1204: {
1209: PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);
1210: return(0);
1211: }
1215: /*@C
1216: PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1217: PC options in the database.
1219: Collective on PC
1221: Input Parameters:
1222: + pc - the preconditioner context
1223: - prefix - the prefix string to prepend to all PC option requests
1225: Notes:
1226: A hyphen (-) must NOT be given at the beginning of the prefix name.
1227: The first character of all runtime options is AUTOMATICALLY the
1228: hyphen.
1230: Level: advanced
1232: .keywords: PC, append, options, prefix, database
1234: .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1235: @*/
1236: int PCAppendOptionsPrefix(PC pc,const char prefix[])
1237: {
1242: PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);
1243: return(0);
1244: }
1248: /*@C
1249: PCGetOptionsPrefix - Gets the prefix used for searching for all
1250: PC options in the database.
1252: Not Collective
1254: Input Parameters:
1255: . pc - the preconditioner context
1257: Output Parameters:
1258: . prefix - pointer to the prefix string used, is returned
1260: Notes: On the fortran side, the user should pass in a string 'prifix' of
1261: sufficient length to hold the prefix.
1263: Level: advanced
1265: .keywords: PC, get, options, prefix, database
1267: .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1268: @*/
1269: int PCGetOptionsPrefix(PC pc,char *prefix[])
1270: {
1275: PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);
1276: return(0);
1277: }
1281: /*@
1282: PCPreSolve - Optional pre-solve phase, intended for any
1283: preconditioner-specific actions that must be performed before
1284: the iterative solve itself.
1286: Collective on PC
1288: Input Parameters:
1289: + pc - the preconditioner context
1290: - ksp - the Krylov subspace context
1292: Level: developer
1294: Sample of Usage:
1295: .vb
1296: PCPreSolve(pc,ksp);
1297: KSPSolve(ksp,its);
1298: PCPostSolve(pc,ksp);
1299: .ve
1301: Notes:
1302: The pre-solve phase is distinct from the PCSetUp() phase.
1304: SLESSolve() calls this directly, so is rarely called by the user.
1306: .keywords: PC, pre-solve
1308: .seealso: PCPostSolve()
1309: @*/
1310: int PCPreSolve(PC pc,KSP ksp)
1311: {
1313: Vec x,rhs;
1314: Mat A,B;
1319: KSPGetSolution(ksp,&x);
1320: KSPGetRhs(ksp,&rhs);
1321: /*
1322: Scale the system and have the matrices use the scaled form
1323: only if the two matrices are actually the same (and hence
1324: have the same scaling
1325: */
1326: PCGetOperators(pc,&A,&B,PETSC_NULL);
1327: if (A == B) {
1328: MatScaleSystem(pc->mat,x,rhs);
1329: MatUseScaledForm(pc->mat,PETSC_TRUE);
1330: }
1332: if (pc->ops->presolve) {
1333: (*pc->ops->presolve)(pc,ksp,x,rhs);
1334: }
1335: return(0);
1336: }
1340: /*@
1341: PCPostSolve - Optional post-solve phase, intended for any
1342: preconditioner-specific actions that must be performed after
1343: the iterative solve itself.
1345: Collective on PC
1347: Input Parameters:
1348: + pc - the preconditioner context
1349: - ksp - the Krylov subspace context
1351: Sample of Usage:
1352: .vb
1353: PCPreSolve(pc,ksp);
1354: KSPSolve(ksp,its);
1355: PCPostSolve(pc,ksp);
1356: .ve
1358: Note:
1359: SLESSolve() calls this routine directly, so it is rarely called by the user.
1361: Level: developer
1363: .keywords: PC, post-solve
1365: .seealso: PCPreSolve(), SLESSolve()
1366: @*/
1367: int PCPostSolve(PC pc,KSP ksp)
1368: {
1370: Vec x,rhs;
1371: Mat A,B;
1375: KSPGetSolution(ksp,&x);
1376: KSPGetRhs(ksp,&rhs);
1377: if (pc->ops->postsolve) {
1378: (*pc->ops->postsolve)(pc,ksp,x,rhs);
1379: }
1381: /*
1382: Scale the system and have the matrices use the scaled form
1383: only if the two matrices are actually the same (and hence
1384: have the same scaling
1385: */
1386: PCGetOperators(pc,&A,&B,PETSC_NULL);
1387: if (A == B) {
1388: MatUnScaleSystem(pc->mat,x,rhs);
1389: MatUseScaledForm(pc->mat,PETSC_FALSE);
1390: }
1391: return(0);
1392: }
1396: /*@C
1397: PCView - Prints the PC data structure.
1399: Collective on PC
1401: Input Parameters:
1402: + PC - the PC context
1403: - viewer - optional visualization context
1405: Note:
1406: The available visualization contexts include
1407: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
1408: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1409: output where only the first processor opens
1410: the file. All other processors send their
1411: data to the first processor to print.
1413: The user can open an alternative visualization contexts with
1414: PetscViewerASCIIOpen() (output to a specified file).
1416: Level: developer
1418: .keywords: PC, view
1420: .seealso: KSPView(), PetscViewerASCIIOpen()
1421: @*/
1422: int PCView(PC pc,PetscViewer viewer)
1423: {
1424: PCType cstr;
1425: int ierr;
1426: PetscTruth mat_exists,isascii,isstring;
1427: PetscViewerFormat format;
1431: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(pc->comm);
1435: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
1436: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
1437: if (isascii) {
1438: PetscViewerGetFormat(viewer,&format);
1439: if (pc->prefix) {
1440: PetscViewerASCIIPrintf(viewer,"PC Object:(%s)\n",pc->prefix);
1441: } else {
1442: PetscViewerASCIIPrintf(viewer,"PC Object:\n");
1443: }
1444: PCGetType(pc,&cstr);
1445: if (cstr) {
1446: PetscViewerASCIIPrintf(viewer," type: %s\n",cstr);
1447: } else {
1448: PetscViewerASCIIPrintf(viewer," type: not yet set\n");
1449: }
1450: if (pc->ops->view) {
1451: PetscViewerASCIIPushTab(viewer);
1452: (*pc->ops->view)(pc,viewer);
1453: PetscViewerASCIIPopTab(viewer);
1454: }
1455: PetscObjectExists((PetscObject)pc->mat,&mat_exists);
1456: if (mat_exists) {
1457: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1458: if (pc->pmat == pc->mat) {
1459: PetscViewerASCIIPrintf(viewer," linear system matrix = precond matrix:\n");
1460: PetscViewerASCIIPushTab(viewer);
1461: MatView(pc->mat,viewer);
1462: PetscViewerASCIIPopTab(viewer);
1463: } else {
1464: PetscObjectExists((PetscObject)pc->pmat,&mat_exists);
1465: if (mat_exists) {
1466: PetscViewerASCIIPrintf(viewer," linear system matrix followed by preconditioner matrix:\n");
1467: } else {
1468: PetscViewerASCIIPrintf(viewer," linear system matrix:\n");
1469: }
1470: PetscViewerASCIIPushTab(viewer);
1471: MatView(pc->mat,viewer);
1472: if (mat_exists) {MatView(pc->pmat,viewer);}
1473: PetscViewerASCIIPopTab(viewer);
1474: }
1475: PetscViewerPopFormat(viewer);
1476: }
1477: } else if (isstring) {
1478: PCGetType(pc,&cstr);
1479: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
1480: if (pc->ops->view) {(*pc->ops->view)(pc,viewer);}
1481: } else {
1482: SETERRQ1(1,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1483: }
1484: return(0);
1485: }
1489: /*@C
1490: PCRegister - See PCRegisterDynamic()
1492: Level: advanced
1493: @*/
1494: int PCRegister(const char sname[],const char path[],const char name[],int (*function)(PC))
1495: {
1496: int ierr;
1497: char fullname[256];
1501: PetscFListConcat(path,name,fullname);
1502: PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);
1503: return(0);
1504: }
1508: /*@
1509: PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1511: Collective on PC
1513: Input Parameter:
1514: . pc - the preconditioner object
1516: Output Parameter:
1517: . mat - the explict preconditioned operator
1519: Notes:
1520: This computation is done by applying the operators to columns of the
1521: identity matrix.
1523: Currently, this routine uses a dense matrix format when 1 processor
1524: is used and a sparse format otherwise. This routine is costly in general,
1525: and is recommended for use only with relatively small systems.
1527: Level: advanced
1528:
1529: .keywords: PC, compute, explicit, operator
1531: @*/
1532: int PCComputeExplicitOperator(PC pc,Mat *mat)
1533: {
1534: Vec in,out;
1535: int ierr,i,M,m,size,*rows,start,end;
1536: MPI_Comm comm;
1537: PetscScalar *array,zero = 0.0,one = 1.0;
1543: comm = pc->comm;
1544: MPI_Comm_size(comm,&size);
1546: PCGetVector(pc,&in);
1547: VecDuplicate(in,&out);
1548: VecGetOwnershipRange(in,&start,&end);
1549: VecGetSize(in,&M);
1550: VecGetLocalSize(in,&m);
1551: PetscMalloc((m+1)*sizeof(int),&rows);
1552: for (i=0; i<m; i++) {rows[i] = start + i;}
1554: if (size == 1) {
1555: MatCreateSeqDense(comm,M,M,PETSC_NULL,mat);
1556: } else {
1557: MatCreateMPIAIJ(comm,m,m,M,M,0,0,0,0,mat);
1558: }
1560: for (i=0; i<M; i++) {
1562: VecSet(&zero,in);
1563: VecSetValues(in,1,&i,&one,INSERT_VALUES);
1564: VecAssemblyBegin(in);
1565: VecAssemblyEnd(in);
1567: /* should fix, allowing user to choose side */
1568: PCApply(pc,in,out,PC_LEFT);
1569:
1570: VecGetArray(out,&array);
1571: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
1572: VecRestoreArray(out,&array);
1574: }
1575: PetscFree(rows);
1576: VecDestroy(out);
1577: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
1578: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
1579: return(0);
1580: }