Actual source code: precon.c
1: #define PETSCKSP_DLL
3: /*
4: The PC (preconditioner) interface routines, callable by users.
5: */
6: #include private/pcimpl.h
8: /* Logging support */
9: PetscCookie PC_COOKIE = 0;
10: PetscEvent PC_SetUp = 0, PC_SetUpOnBlocks = 0, PC_Apply = 0, PC_ApplyCoarse = 0, PC_ApplyMultiple = 0, PC_ApplySymmetricLeft = 0;
11: PetscEvent PC_ApplySymmetricRight = 0, PC_ModifySubMatrices = 0, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks;
15: PetscErrorCode PCGetDefaultType_Private(PC pc,const char* type[])
16: {
18: PetscMPIInt size;
19: PetscTruth flg1,flg2,set,flg3;
22: MPI_Comm_size(((PetscObject)pc)->comm,&size);
23: if (pc->pmat) {
24: PetscErrorCode (*f)(Mat,PetscTruth*,MatReuse,Mat*);
25: PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",(void (**)(void))&f);
26: if (size == 1) {
27: MatHasOperation(pc->pmat,MATOP_ICCFACTOR_SYMBOLIC,&flg1);
28: MatHasOperation(pc->pmat,MATOP_ILUFACTOR_SYMBOLIC,&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: 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: PetscErrorCode PCDestroy(PC pc)
73: {
78: if (--((PetscObject)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->diagonalscaleright) {VecDestroy(pc->diagonalscaleright);}
85: if (pc->diagonalscaleleft) {VecDestroy(pc->diagonalscaleleft);}
87: if (pc->pmat) {MatDestroy(pc->pmat);}
88: if (pc->mat) {MatDestroy(pc->mat);}
90: PetscHeaderDestroy(pc);
91: return(0);
92: }
96: /*@C
97: PCDiagonalScale - Indicates if the preconditioner applies an additional left and right
98: scaling as needed by certain time-stepping codes.
100: Collective on PC
102: Input Parameter:
103: . pc - the preconditioner context
105: Output Parameter:
106: . flag - PETSC_TRUE if it applies the scaling
108: Level: developer
110: Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
111: $ D M A D^{-1} y = D M b for left preconditioning or
112: $ D A M D^{-1} z = D b for right preconditioning
114: .keywords: PC
116: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScaleSet()
117: @*/
118: PetscErrorCode PCDiagonalScale(PC pc,PetscTruth *flag)
119: {
123: *flag = pc->diagonalscale;
124: return(0);
125: }
129: /*@
130: PCDiagonalScaleSet - Indicates the left scaling to use to apply an additional left and right
131: scaling as needed by certain time-stepping codes.
133: Collective on PC
135: Input Parameters:
136: + pc - the preconditioner context
137: - s - scaling vector
139: Level: intermediate
141: Notes: The system solved via the Krylov method is
142: $ D M A D^{-1} y = D M b for left preconditioning or
143: $ D A M D^{-1} z = D b for right preconditioning
145: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
147: .keywords: PC
149: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScale()
150: @*/
151: PetscErrorCode PCDiagonalScaleSet(PC pc,Vec s)
152: {
158: pc->diagonalscale = PETSC_TRUE;
159: PetscObjectReference((PetscObject)s);
160: if (pc->diagonalscaleleft) {
161: VecDestroy(pc->diagonalscaleleft);
162: }
163: pc->diagonalscaleleft = s;
164: if (!pc->diagonalscaleright) {
165: VecDuplicate(s,&pc->diagonalscaleright);
166: }
167: VecCopy(s,pc->diagonalscaleright);
168: VecReciprocal(pc->diagonalscaleright);
169: return(0);
170: }
174: /*@
175: PCDiagonalScaleLeft - Indicates the left scaling to use to apply an additional left and right
176: scaling as needed by certain time-stepping codes.
178: Collective on PC
180: Input Parameters:
181: + pc - the preconditioner context
182: . in - input vector
183: + out - scaled vector (maybe the same as in)
185: Level: intermediate
187: Notes: The system solved via the Krylov method is
188: $ D M A D^{-1} y = D M b for left preconditioning or
189: $ D A M D^{-1} z = D b for right preconditioning
191: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
193: If diagonal scaling is turned off and in is not out then in is copied to out
195: .keywords: PC
197: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
198: @*/
199: PetscErrorCode PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
200: {
207: if (pc->diagonalscale) {
208: VecPointwiseMult(out,pc->diagonalscaleleft,in);
209: } else if (in != out) {
210: VecCopy(in,out);
211: }
212: return(0);
213: }
217: /*@
218: PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
220: Collective on PC
222: Input Parameters:
223: + pc - the preconditioner context
224: . in - input vector
225: + out - scaled vector (maybe the same as in)
227: Level: intermediate
229: Notes: The system solved via the Krylov method is
230: $ D M A D^{-1} y = D M b for left preconditioning or
231: $ D A M D^{-1} z = D b for right preconditioning
233: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
235: If diagonal scaling is turned off and in is not out then in is copied to out
237: .keywords: PC
239: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
240: @*/
241: PetscErrorCode PCDiagonalScaleRight(PC pc,Vec in,Vec out)
242: {
249: if (pc->diagonalscale) {
250: VecPointwiseMult(out,pc->diagonalscaleright,in);
251: } else if (in != out) {
252: VecCopy(in,out);
253: }
254: return(0);
255: }
257: #if 0
260: static PetscErrorCode PCPublish_Petsc(PetscObject obj)
261: {
263: return(0);
264: }
265: #endif
269: /*@
270: PCCreate - Creates a preconditioner context.
272: Collective on MPI_Comm
274: Input Parameter:
275: . comm - MPI communicator
277: Output Parameter:
278: . pc - location to put the preconditioner context
280: Notes:
281: The default preconditioner on one processor is PCILU with 0 fill on more
282: then one it is PCBJACOBI with ILU() on each processor.
284: Level: developer
286: .keywords: PC, create, context
288: .seealso: PCSetUp(), PCApply(), PCDestroy()
289: @*/
290: PetscErrorCode PCCreate(MPI_Comm comm,PC *newpc)
291: {
292: PC pc;
297: *newpc = 0;
298: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
299: PCInitializePackage(PETSC_NULL);
300: #endif
302: PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_COOKIE,-1,"PC",comm,PCDestroy,PCView);
304: pc->mat = 0;
305: pc->pmat = 0;
306: pc->setupcalled = 0;
307: pc->setfromoptionscalled = 0;
308: pc->data = 0;
309: pc->diagonalscale = PETSC_FALSE;
310: pc->diagonalscaleleft = 0;
311: pc->diagonalscaleright = 0;
313: pc->modifysubmatrices = 0;
314: pc->modifysubmatricesP = 0;
315: PetscPublishAll(pc);
316: *newpc = pc;
317: return(0);
319: }
321: /* -------------------------------------------------------------------------------*/
325: /*@
326: PCApply - Applies the preconditioner to a vector.
328: Collective on PC and Vec
330: Input Parameters:
331: + pc - the preconditioner context
332: - x - input vector
334: Output Parameter:
335: . y - output vector
337: Level: developer
339: .keywords: PC, apply
341: .seealso: PCApplyTranspose(), PCApplyBAorAB()
342: @*/
343: PetscErrorCode PCApply(PC pc,Vec x,Vec y)
344: {
351: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
352: if (pc->setupcalled < 2) {
353: PCSetUp(pc);
354: }
355: if (!pc->ops->apply) SETERRQ(PETSC_ERR_SUP,"PC does not have apply");
357: (*pc->ops->apply)(pc,x,y);
359: return(0);
360: }
364: /*@
365: PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
367: Collective on PC and Vec
369: Input Parameters:
370: + pc - the preconditioner context
371: - x - input vector
373: Output Parameter:
374: . y - output vector
376: Notes:
377: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
379: Level: developer
381: .keywords: PC, apply, symmetric, left
383: .seealso: PCApply(), PCApplySymmetricRight()
384: @*/
385: PetscErrorCode PCApplySymmetricLeft(PC pc,Vec x,Vec y)
386: {
393: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
394: if (pc->setupcalled < 2) {
395: PCSetUp(pc);
396: }
397: if (!pc->ops->applysymmetricleft) SETERRQ(PETSC_ERR_SUP,"PC does not have left symmetric apply");
399: (*pc->ops->applysymmetricleft)(pc,x,y);
401: return(0);
402: }
406: /*@
407: PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
409: Collective on PC and Vec
411: Input Parameters:
412: + pc - the preconditioner context
413: - x - input vector
415: Output Parameter:
416: . y - output vector
418: Level: developer
420: Notes:
421: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
423: .keywords: PC, apply, symmetric, right
425: .seealso: PCApply(), PCApplySymmetricLeft()
426: @*/
427: PetscErrorCode PCApplySymmetricRight(PC pc,Vec x,Vec y)
428: {
435: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
436: if (pc->setupcalled < 2) {
437: PCSetUp(pc);
438: }
439: if (!pc->ops->applysymmetricright) SETERRQ(PETSC_ERR_SUP,"PC does not have left symmetric apply");
441: (*pc->ops->applysymmetricright)(pc,x,y);
443: return(0);
444: }
448: /*@
449: PCApplyTranspose - Applies the transpose of preconditioner to a vector.
451: Collective on PC and Vec
453: Input Parameters:
454: + pc - the preconditioner context
455: - x - input vector
457: Output Parameter:
458: . y - output vector
460: Level: developer
462: .keywords: PC, apply, transpose
464: .seealso: PCApply(), PCApplyBAorAB(), PCApplyBAorABTranspose(), PCApplyTransposeExists()
465: @*/
466: PetscErrorCode PCApplyTranspose(PC pc,Vec x,Vec y)
467: {
474: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
475: if (pc->setupcalled < 2) {
476: PCSetUp(pc);
477: }
478: if (!pc->ops->applytranspose) SETERRQ(PETSC_ERR_SUP,"PC does not have apply transpose");
480: (*pc->ops->applytranspose)(pc,x,y);
482: return(0);
483: }
487: /*@
488: PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation
490: Collective on PC and Vec
492: Input Parameters:
493: . pc - the preconditioner context
495: Output Parameter:
496: . flg - PETSC_TRUE if a transpose operation is defined
498: Level: developer
500: .keywords: PC, apply, transpose
502: .seealso: PCApplyTranspose()
503: @*/
504: PetscErrorCode PCApplyTransposeExists(PC pc,PetscTruth *flg)
505: {
509: if (pc->ops->applytranspose) *flg = PETSC_TRUE;
510: else *flg = PETSC_FALSE;
511: return(0);
512: }
516: /*@
517: PCApplyBAorAB - Applies the preconditioner and operator to a vector.
519: Collective on PC and Vec
521: Input Parameters:
522: + pc - the preconditioner context
523: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
524: . x - input vector
525: - work - work vector
527: Output Parameter:
528: . y - output vector
530: Level: developer
532: .keywords: PC, apply, operator
534: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
535: @*/
536: PetscErrorCode PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
537: {
545: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
546: if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) {
547: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric");
548: }
549: if (pc->diagonalscale && side == PC_SYMMETRIC) {
550: SETERRQ(PETSC_ERR_SUP,"Cannot include diagonal scaling with symmetric preconditioner application");
551: }
553: if (pc->setupcalled < 2) {
554: PCSetUp(pc);
555: }
557: if (pc->diagonalscale) {
558: if (pc->ops->applyBA) {
559: Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
560: VecDuplicate(x,&work2);
561: PCDiagonalScaleRight(pc,x,work2);
562: (*pc->ops->applyBA)(pc,side,work2,y,work);
563: PCDiagonalScaleLeft(pc,y,y);
564: VecDestroy(work2);
565: } else if (side == PC_RIGHT) {
566: PCDiagonalScaleRight(pc,x,y);
567: PCApply(pc,y,work);
568: MatMult(pc->mat,work,y);
569: PCDiagonalScaleLeft(pc,y,y);
570: } else if (side == PC_LEFT) {
571: PCDiagonalScaleRight(pc,x,y);
572: MatMult(pc->mat,y,work);
573: PCApply(pc,work,y);
574: PCDiagonalScaleLeft(pc,y,y);
575: } else if (side == PC_SYMMETRIC) {
576: SETERRQ(PETSC_ERR_SUP,"Cannot provide diagonal scaling with symmetric application of preconditioner");
577: }
578: } else {
579: if (pc->ops->applyBA) {
580: (*pc->ops->applyBA)(pc,side,x,y,work);
581: } else if (side == PC_RIGHT) {
582: PCApply(pc,x,work);
583: MatMult(pc->mat,work,y);
584: } else if (side == PC_LEFT) {
585: MatMult(pc->mat,x,work);
586: PCApply(pc,work,y);
587: } else if (side == PC_SYMMETRIC) {
588: /* There's an extra copy here; maybe should provide 2 work vectors instead? */
589: PCApplySymmetricRight(pc,x,work);
590: MatMult(pc->mat,work,y);
591: VecCopy(y,work);
592: PCApplySymmetricLeft(pc,work,y);
593: }
594: }
595: return(0);
596: }
600: /*@
601: PCApplyBAorABTranspose - Applies the transpose of the preconditioner
602: and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
603: not tr(B*A) = tr(A)*tr(B).
605: Collective on PC and Vec
607: Input Parameters:
608: + pc - the preconditioner context
609: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
610: . x - input vector
611: - work - work vector
613: Output Parameter:
614: . y - output vector
616: Level: developer
618: .keywords: PC, apply, operator, transpose
620: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
621: @*/
622: PetscErrorCode PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
623: {
631: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
632: if (pc->ops->applyBAtranspose) {
633: (*pc->ops->applyBAtranspose)(pc,side,x,y,work);
634: return(0);
635: }
636: if (side != PC_LEFT && side != PC_RIGHT) {
637: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
638: }
640: if (pc->setupcalled < 2) {
641: PCSetUp(pc);
642: }
644: if (side == PC_RIGHT) {
645: PCApplyTranspose(pc,x,work);
646: MatMultTranspose(pc->mat,work,y);
647: } else if (side == PC_LEFT) {
648: MatMultTranspose(pc->mat,x,work);
649: PCApplyTranspose(pc,work,y);
650: }
651: /* add support for PC_SYMMETRIC */
652: return(0); /* actually will never get here */
653: }
655: /* -------------------------------------------------------------------------------*/
659: /*@
660: PCApplyRichardsonExists - Determines whether a particular preconditioner has a
661: built-in fast application of Richardson's method.
663: Not Collective
665: Input Parameter:
666: . pc - the preconditioner
668: Output Parameter:
669: . exists - PETSC_TRUE or PETSC_FALSE
671: Level: developer
673: .keywords: PC, apply, Richardson, exists
675: .seealso: PCApplyRichardson()
676: @*/
677: PetscErrorCode PCApplyRichardsonExists(PC pc,PetscTruth *exists)
678: {
682: if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
683: else *exists = PETSC_FALSE;
684: return(0);
685: }
689: /*@
690: PCApplyRichardson - Applies several steps of Richardson iteration with
691: the particular preconditioner. This routine is usually used by the
692: Krylov solvers and not the application code directly.
694: Collective on PC
696: Input Parameters:
697: + pc - the preconditioner context
698: . x - the initial guess
699: . w - one work vector
700: . rtol - relative decrease in residual norm convergence criteria
701: . abstol - absolute residual norm convergence criteria
702: . dtol - divergence residual norm increase criteria
703: - its - the number of iterations to apply.
705: Output Parameter:
706: . y - the solution
708: Notes:
709: Most preconditioners do not support this function. Use the command
710: PCApplyRichardsonExists() to determine if one does.
712: Except for the multigrid PC this routine ignores the convergence tolerances
713: and always runs for the number of iterations
714:
715: Level: developer
717: .keywords: PC, apply, Richardson
719: .seealso: PCApplyRichardsonExists()
720: @*/
721: PetscErrorCode PCApplyRichardson(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
722: {
730: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
731: if (pc->setupcalled < 2) {
732: PCSetUp(pc);
733: }
734: if (!pc->ops->applyrichardson) SETERRQ(PETSC_ERR_SUP,"PC does not have apply richardson");
735: (*pc->ops->applyrichardson)(pc,x,y,w,rtol,abstol,dtol,its);
736: return(0);
737: }
739: /*
740: a setupcall of 0 indicates never setup,
741: 1 needs to be resetup,
742: 2 does not need any changes.
743: */
746: /*@
747: PCSetUp - Prepares for the use of a preconditioner.
749: Collective on PC
751: Input Parameter:
752: . pc - the preconditioner context
754: Level: developer
756: .keywords: PC, setup
758: .seealso: PCCreate(), PCApply(), PCDestroy()
759: @*/
760: PetscErrorCode PCSetUp(PC pc)
761: {
763: const char *def;
768: if (!pc->mat) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
770: if (pc->setupcalled > 1) {
771: PetscInfo(pc,"Setting PC with identical preconditioner\n");
772: return(0);
773: } else if (!pc->setupcalled) {
774: PetscInfo(pc,"Setting up new PC\n");
775: } else if (pc->flag == SAME_NONZERO_PATTERN) {
776: PetscInfo(pc,"Setting up PC with same nonzero pattern\n");
777: } else {
778: PetscInfo(pc,"Setting up PC with different nonzero pattern\n");
779: }
781: if (!((PetscObject)pc)->type_name) {
782: PCGetDefaultType_Private(pc,&def);
783: PCSetType(pc,def);
784: }
787: if (pc->ops->setup) {
788: (*pc->ops->setup)(pc);
789: }
790: pc->setupcalled = 2;
792: return(0);
793: }
797: /*@
798: PCSetUpOnBlocks - Sets up the preconditioner for each block in
799: the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
800: methods.
802: Collective on PC
804: Input Parameters:
805: . pc - the preconditioner context
807: Level: developer
809: .keywords: PC, setup, blocks
811: .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
812: @*/
813: PetscErrorCode PCSetUpOnBlocks(PC pc)
814: {
819: if (!pc->ops->setuponblocks) return(0);
821: (*pc->ops->setuponblocks)(pc);
823: return(0);
824: }
828: /*@C
829: PCSetModifySubMatrices - Sets a user-defined routine for modifying the
830: submatrices that arise within certain subdomain-based preconditioners.
831: The basic submatrices are extracted from the preconditioner matrix as
832: usual; the user can then alter these (for example, to set different boundary
833: conditions for each submatrix) before they are used for the local solves.
835: Collective on PC
837: Input Parameters:
838: + pc - the preconditioner context
839: . func - routine for modifying the submatrices
840: - ctx - optional user-defined context (may be null)
842: Calling sequence of func:
843: $ func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
845: . row - an array of index sets that contain the global row numbers
846: that comprise each local submatrix
847: . col - an array of index sets that contain the global column numbers
848: that comprise each local submatrix
849: . submat - array of local submatrices
850: - ctx - optional user-defined context for private data for the
851: user-defined func routine (may be null)
853: Notes:
854: PCSetModifySubMatrices() MUST be called before KSPSetUp() and
855: KSPSolve().
857: A routine set by PCSetModifySubMatrices() is currently called within
858: the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
859: preconditioners. All other preconditioners ignore this routine.
861: Level: advanced
863: .keywords: PC, set, modify, submatrices
865: .seealso: PCModifySubMatrices()
866: @*/
867: PetscErrorCode PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
868: {
871: pc->modifysubmatrices = func;
872: pc->modifysubmatricesP = ctx;
873: return(0);
874: }
878: /*@C
879: PCModifySubMatrices - Calls an optional user-defined routine within
880: certain preconditioners if one has been set with PCSetModifySubMarices().
882: Collective on PC
884: Input Parameters:
885: + pc - the preconditioner context
886: . nsub - the number of local submatrices
887: . row - an array of index sets that contain the global row numbers
888: that comprise each local submatrix
889: . col - an array of index sets that contain the global column numbers
890: that comprise each local submatrix
891: . submat - array of local submatrices
892: - ctx - optional user-defined context for private data for the
893: user-defined routine (may be null)
895: Output Parameter:
896: . submat - array of local submatrices (the entries of which may
897: have been modified)
899: Notes:
900: The user should NOT generally call this routine, as it will
901: automatically be called within certain preconditioners (currently
902: block Jacobi, additive Schwarz) if set.
904: The basic submatrices are extracted from the preconditioner matrix
905: as usual; the user can then alter these (for example, to set different
906: boundary conditions for each submatrix) before they are used for the
907: local solves.
909: Level: developer
911: .keywords: PC, modify, submatrices
913: .seealso: PCSetModifySubMatrices()
914: @*/
915: PetscErrorCode PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
916: {
921: if (!pc->modifysubmatrices) return(0);
923: (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);
925: return(0);
926: }
930: /*@
931: PCSetOperators - Sets the matrix associated with the linear system and
932: a (possibly) different one associated with the preconditioner.
934: Collective on PC and Mat
936: Input Parameters:
937: + pc - the preconditioner context
938: . Amat - the matrix associated with the linear system
939: . Pmat - the matrix to be used in constructing the preconditioner, usually
940: the same as Amat.
941: - flag - flag indicating information about the preconditioner matrix structure
942: during successive linear solves. This flag is ignored the first time a
943: linear system is solved, and thus is irrelevant when solving just one linear
944: system.
946: Notes:
947: The flag can be used to eliminate unnecessary work in the preconditioner
948: during the repeated solution of linear systems of the same size. The
949: available options are
950: + SAME_PRECONDITIONER -
951: Pmat is identical during successive linear solves.
952: This option is intended for folks who are using
953: different Amat and Pmat matrices and wish to reuse the
954: same preconditioner matrix. For example, this option
955: saves work by not recomputing incomplete factorization
956: for ILU/ICC preconditioners.
957: . SAME_NONZERO_PATTERN -
958: Pmat has the same nonzero structure during
959: successive linear solves.
960: - DIFFERENT_NONZERO_PATTERN -
961: Pmat does not have the same nonzero structure.
963: Caution:
964: If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
965: and does not check the structure of the matrix. If you erroneously
966: claim that the structure is the same when it actually is not, the new
967: preconditioner will not function correctly. Thus, use this optimization
968: feature carefully!
970: If in doubt about whether your preconditioner matrix has changed
971: structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
973: More Notes about Repeated Solution of Linear Systems:
974: PETSc does NOT reset the matrix entries of either Amat or Pmat
975: to zero after a linear solve; the user is completely responsible for
976: matrix assembly. See the routine MatZeroEntries() if desiring to
977: zero all elements of a matrix.
979: Level: intermediate
981: .keywords: PC, set, operators, matrix, linear system
983: .seealso: PCGetOperators(), MatZeroEntries()
984: @*/
985: PetscErrorCode PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
986: {
988: PetscTruth isbjacobi,isrowbs;
997: /*
998: BlockSolve95 cannot use default BJacobi preconditioning
999: */
1000: if (Amat) {
1001: PetscTypeCompare((PetscObject)Amat,MATMPIROWBS,&isrowbs);
1002: if (isrowbs) {
1003: PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);
1004: if (isbjacobi) {
1005: PCSetType(pc,PCILU);
1006: PetscInfo(pc,"Switching default PC to PCILU since BS95 doesn't support PCBJACOBI\n");
1007: }
1008: }
1009: }
1011: /* reference first in case the matrices are the same */
1012: if (Amat) {PetscObjectReference((PetscObject)Amat);}
1013: if (pc->mat) {MatDestroy(pc->mat);}
1014: if (Pmat) {PetscObjectReference((PetscObject)Pmat);}
1015: if (pc->pmat) {MatDestroy(pc->pmat);}
1016: pc->mat = Amat;
1017: 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: Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1047: are created in PC and returned to the user. In this case, if both operators
1048: mat and pmat are requested, two DIFFERENT operators will be returned. If
1049: only one is requested both operators in the PC will be the same (i.e. as
1050: if one had called KSP/PCSetOperators() with the same argument for both Mats).
1051: The user must set the sizes of the returned matrices and their type etc just
1052: as if the user created them with MatCreate(). For example,
1054: $ KSP/PCGetOperators(ksp/pc,&mat,PETSC_NULL,PETSC_NULL); is equivalent to
1055: $ set size, type, etc of mat
1057: $ MatCreate(comm,&mat);
1058: $ KSP/PCSetOperators(ksp/pc,mat,mat,SAME_NONZERO_PATTERN);
1059: $ PetscObjectDereference((PetscObject)mat);
1060: $ set size, type, etc of mat
1062: and
1064: $ KSP/PCGetOperators(ksp/pc,&mat,&pmat,PETSC_NULL); is equivalent to
1065: $ set size, type, etc of mat and pmat
1067: $ MatCreate(comm,&mat);
1068: $ MatCreate(comm,&pmat);
1069: $ KSP/PCSetOperators(ksp/pc,mat,pmat,SAME_NONZERO_PATTERN);
1070: $ PetscObjectDereference((PetscObject)mat);
1071: $ PetscObjectDereference((PetscObject)pmat);
1072: $ set size, type, etc of mat and pmat
1074: The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1075: of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1076: managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1077: at this is when you create a SNES you do not NEED to create a KSP and attach it to
1078: the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1079: you do not need to attach a PC to it (the KSP object manages the PC object for you).
1080: Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1081: it can be created for you?
1082:
1084: .keywords: PC, get, operators, matrix, linear system
1086: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1087: @*/
1088: PetscErrorCode PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1089: {
1094: if (mat) {
1095: if (!pc->mat) {
1096: MatCreate(((PetscObject)pc)->comm,&pc->mat);
1097: if (!pc->pmat && !pmat) { /* user did NOT request pmat, so make same as mat */
1098: pc->pmat = pc->mat;
1099: PetscObjectReference((PetscObject)pc->pmat);
1100: }
1101: }
1102: *mat = pc->mat;
1103: }
1104: if (pmat) {
1105: if (!pc->pmat) {
1106: MatCreate(((PetscObject)pc)->comm,&pc->mat);
1107: if (!pc->mat && !mat) { /* user did NOT request mat, so make same as pmat */
1108: pc->mat = pc->pmat;
1109: PetscObjectReference((PetscObject)pc->mat);
1110: }
1111: }
1112: *pmat = pc->pmat;
1113: }
1114: if (flag) *flag = pc->flag;
1115: return(0);
1116: }
1120: /*@C
1121: PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1122: possibly a different one associated with the preconditioner have been set in the PC.
1124: Not collective, though the results on all processes should be the same
1126: Input Parameter:
1127: . pc - the preconditioner context
1129: Output Parameters:
1130: + mat - the matrix associated with the linear system was set
1131: - pmat - matrix associated with the preconditioner was set, usually the same
1133: Level: intermediate
1135: .keywords: PC, get, operators, matrix, linear system
1137: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1138: @*/
1139: PetscErrorCode PCGetOperatorsSet(PC pc,PetscTruth *mat,PetscTruth *pmat)
1140: {
1143: if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1144: if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1145: return(0);
1146: }
1150: /*@
1151: PCFactorGetMatrix - Gets the factored matrix from the
1152: preconditioner context. This routine is valid only for the LU,
1153: incomplete LU, Cholesky, and incomplete Cholesky methods.
1155: Not Collective on PC though Mat is parallel if PC is parallel
1157: Input Parameters:
1158: . pc - the preconditioner context
1160: Output parameters:
1161: . mat - the factored matrix
1163: Level: advanced
1165: Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1167: .keywords: PC, get, factored, matrix
1168: @*/
1169: PetscErrorCode PCFactorGetMatrix(PC pc,Mat *mat)
1170: {
1176: if (pc->ops->getfactoredmatrix) {
1177: (*pc->ops->getfactoredmatrix)(pc,mat);
1178: }
1179: return(0);
1180: }
1184: /*@C
1185: PCSetOptionsPrefix - Sets the prefix used for searching for all
1186: PC options in the database.
1188: Collective on PC
1190: Input Parameters:
1191: + pc - the preconditioner context
1192: - prefix - the prefix string to prepend to all PC option requests
1194: Notes:
1195: A hyphen (-) must NOT be given at the beginning of the prefix name.
1196: The first character of all runtime options is AUTOMATICALLY the
1197: hyphen.
1199: Level: advanced
1201: .keywords: PC, set, options, prefix, database
1203: .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1204: @*/
1205: PetscErrorCode PCSetOptionsPrefix(PC pc,const char prefix[])
1206: {
1211: PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);
1212: return(0);
1213: }
1217: /*@C
1218: PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1219: PC options in the database.
1221: Collective on PC
1223: Input Parameters:
1224: + pc - the preconditioner context
1225: - prefix - the prefix string to prepend to all PC option requests
1227: Notes:
1228: A hyphen (-) must NOT be given at the beginning of the prefix name.
1229: The first character of all runtime options is AUTOMATICALLY the
1230: hyphen.
1232: Level: advanced
1234: .keywords: PC, append, options, prefix, database
1236: .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1237: @*/
1238: PetscErrorCode PCAppendOptionsPrefix(PC pc,const char prefix[])
1239: {
1244: PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);
1245: return(0);
1246: }
1250: /*@C
1251: PCGetOptionsPrefix - Gets the prefix used for searching for all
1252: PC options in the database.
1254: Not Collective
1256: Input Parameters:
1257: . pc - the preconditioner context
1259: Output Parameters:
1260: . prefix - pointer to the prefix string used, is returned
1262: Notes: On the fortran side, the user should pass in a string 'prifix' of
1263: sufficient length to hold the prefix.
1265: Level: advanced
1267: .keywords: PC, get, options, prefix, database
1269: .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1270: @*/
1271: PetscErrorCode PCGetOptionsPrefix(PC pc,const char *prefix[])
1272: {
1278: PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);
1279: return(0);
1280: }
1284: /*@
1285: PCPreSolve - Optional pre-solve phase, intended for any
1286: preconditioner-specific actions that must be performed before
1287: the iterative solve itself.
1289: Collective on PC
1291: Input Parameters:
1292: + pc - the preconditioner context
1293: - ksp - the Krylov subspace context
1295: Level: developer
1297: Sample of Usage:
1298: .vb
1299: PCPreSolve(pc,ksp);
1300: KSPSolve(ksp,b,x);
1301: PCPostSolve(pc,ksp);
1302: .ve
1304: Notes:
1305: The pre-solve phase is distinct from the PCSetUp() phase.
1307: KSPSolve() calls this directly, so is rarely called by the user.
1309: .keywords: PC, pre-solve
1311: .seealso: PCPostSolve()
1312: @*/
1313: PetscErrorCode PCPreSolve(PC pc,KSP ksp)
1314: {
1316: Vec x,rhs;
1317: Mat A,B;
1322: KSPGetSolution(ksp,&x);
1323: KSPGetRhs(ksp,&rhs);
1324: /*
1325: Scale the system and have the matrices use the scaled form
1326: only if the two matrices are actually the same (and hence
1327: have the same scaling
1328: */
1329: PCGetOperators(pc,&A,&B,PETSC_NULL);
1330: if (A == B) {
1331: MatScaleSystem(pc->mat,rhs,x);
1332: MatUseScaledForm(pc->mat,PETSC_TRUE);
1333: }
1335: if (pc->ops->presolve) {
1336: (*pc->ops->presolve)(pc,ksp,rhs,x);
1337: }
1338: return(0);
1339: }
1343: /*@
1344: PCPostSolve - Optional post-solve phase, intended for any
1345: preconditioner-specific actions that must be performed after
1346: the iterative solve itself.
1348: Collective on PC
1350: Input Parameters:
1351: + pc - the preconditioner context
1352: - ksp - the Krylov subspace context
1354: Sample of Usage:
1355: .vb
1356: PCPreSolve(pc,ksp);
1357: KSPSolve(ksp,b,x);
1358: PCPostSolve(pc,ksp);
1359: .ve
1361: Note:
1362: KSPSolve() calls this routine directly, so it is rarely called by the user.
1364: Level: developer
1366: .keywords: PC, post-solve
1368: .seealso: PCPreSolve(), KSPSolve()
1369: @*/
1370: PetscErrorCode PCPostSolve(PC pc,KSP ksp)
1371: {
1373: Vec x,rhs;
1374: Mat A,B;
1379: KSPGetSolution(ksp,&x);
1380: KSPGetRhs(ksp,&rhs);
1381: if (pc->ops->postsolve) {
1382: (*pc->ops->postsolve)(pc,ksp,rhs,x);
1383: }
1384: /*
1385: Scale the system and have the matrices use the scaled form
1386: only if the two matrices are actually the same (and hence
1387: have the same scaling
1388: */
1389: PCGetOperators(pc,&A,&B,PETSC_NULL);
1390: if (A == B) {
1391: MatUnScaleSystem(pc->mat,rhs,x);
1392: MatUseScaledForm(pc->mat,PETSC_FALSE);
1393: }
1394: return(0);
1395: }
1399: /*@C
1400: PCView - Prints the PC data structure.
1402: Collective on PC
1404: Input Parameters:
1405: + PC - the PC context
1406: - viewer - optional visualization context
1408: Note:
1409: The available visualization contexts include
1410: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
1411: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1412: output where only the first processor opens
1413: the file. All other processors send their
1414: data to the first processor to print.
1416: The user can open an alternative visualization contexts with
1417: PetscViewerASCIIOpen() (output to a specified file).
1419: Level: developer
1421: .keywords: PC, view
1423: .seealso: KSPView(), PetscViewerASCIIOpen()
1424: @*/
1425: PetscErrorCode PCView(PC pc,PetscViewer viewer)
1426: {
1427: PCType cstr;
1428: PetscErrorCode ierr;
1429: PetscTruth mat_exists,iascii,isstring;
1430: PetscViewerFormat format;
1434: if (!viewer) {
1435: PetscViewerASCIIGetStdout(((PetscObject)pc)->comm,&viewer);
1436: }
1440: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
1441: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
1442: if (iascii) {
1443: PetscViewerGetFormat(viewer,&format);
1444: if (((PetscObject)pc)->prefix) {
1445: PetscViewerASCIIPrintf(viewer,"PC Object:(%s)\n",((PetscObject)pc)->prefix);
1446: } else {
1447: PetscViewerASCIIPrintf(viewer,"PC Object:\n");
1448: }
1449: PCGetType(pc,&cstr);
1450: if (cstr) {
1451: PetscViewerASCIIPrintf(viewer," type: %s\n",cstr);
1452: } else {
1453: PetscViewerASCIIPrintf(viewer," type: not yet set\n");
1454: }
1455: if (pc->ops->view) {
1456: PetscViewerASCIIPushTab(viewer);
1457: (*pc->ops->view)(pc,viewer);
1458: PetscViewerASCIIPopTab(viewer);
1459: }
1460: PetscObjectExists((PetscObject)pc->mat,&mat_exists);
1461: if (mat_exists) {
1462: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1463: if (pc->pmat == pc->mat) {
1464: PetscViewerASCIIPrintf(viewer," linear system matrix = precond matrix:\n");
1465: PetscViewerASCIIPushTab(viewer);
1466: MatView(pc->mat,viewer);
1467: PetscViewerASCIIPopTab(viewer);
1468: } else {
1469: PetscObjectExists((PetscObject)pc->pmat,&mat_exists);
1470: if (mat_exists) {
1471: PetscViewerASCIIPrintf(viewer," linear system matrix followed by preconditioner matrix:\n");
1472: } else {
1473: PetscViewerASCIIPrintf(viewer," linear system matrix:\n");
1474: }
1475: PetscViewerASCIIPushTab(viewer);
1476: MatView(pc->mat,viewer);
1477: if (mat_exists) {MatView(pc->pmat,viewer);}
1478: PetscViewerASCIIPopTab(viewer);
1479: }
1480: PetscViewerPopFormat(viewer);
1481: }
1482: } else if (isstring) {
1483: PCGetType(pc,&cstr);
1484: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
1485: if (pc->ops->view) {(*pc->ops->view)(pc,viewer);}
1486: } else {
1487: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1488: }
1489: return(0);
1490: }
1495: /*@
1496: PCSetInitialGuessNonzero - Tells the iterative solver that the
1497: initial guess is nonzero; otherwise PC assumes the initial guess
1498: is to be zero (and thus zeros it out before solving).
1500: Collective on PC
1502: Input Parameters:
1503: + pc - iterative context obtained from PCCreate()
1504: - flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1506: Level: Developer
1508: Notes:
1509: This is a weird function. Since PC's are linear operators on the right hand side they
1510: CANNOT use an initial guess. This function is for the "pass-through" preconditioners
1511: PCKSP, PCREDUNDANT and PCOPENMP and causes the inner KSP object to use the nonzero
1512: initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP.
1515: .keywords: PC, set, initial guess, nonzero
1517: .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
1518: @*/
1519: PetscErrorCode PCSetInitialGuessNonzero(PC pc,PetscTruth flg)
1520: {
1522: pc->nonzero_guess = flg;
1523: return(0);
1524: }
1528: /*@C
1529: PCRegister - See PCRegisterDynamic()
1531: Level: advanced
1532: @*/
1533: PetscErrorCode PCRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PC))
1534: {
1536: char fullname[PETSC_MAX_PATH_LEN];
1540: PetscFListConcat(path,name,fullname);
1541: PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);
1542: return(0);
1543: }
1547: /*@
1548: PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1550: Collective on PC
1552: Input Parameter:
1553: . pc - the preconditioner object
1555: Output Parameter:
1556: . mat - the explict preconditioned operator
1558: Notes:
1559: This computation is done by applying the operators to columns of the
1560: identity matrix.
1562: Currently, this routine uses a dense matrix format when 1 processor
1563: is used and a sparse format otherwise. This routine is costly in general,
1564: and is recommended for use only with relatively small systems.
1566: Level: advanced
1567:
1568: .keywords: PC, compute, explicit, operator
1570: @*/
1571: PetscErrorCode PCComputeExplicitOperator(PC pc,Mat *mat)
1572: {
1573: Vec in,out;
1575: PetscInt i,M,m,*rows,start,end;
1576: PetscMPIInt size;
1577: MPI_Comm comm;
1578: PetscScalar *array,one = 1.0;
1579:
1584: comm = ((PetscObject)pc)->comm;
1585: MPI_Comm_size(comm,&size);
1587: if (!pc->pmat) SETERRQ(PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1588: MatGetVecs(pc->pmat,&in,0);
1589: VecDuplicate(in,&out);
1590: VecGetOwnershipRange(in,&start,&end);
1591: VecGetSize(in,&M);
1592: VecGetLocalSize(in,&m);
1593: PetscMalloc((m+1)*sizeof(PetscInt),&rows);
1594: for (i=0; i<m; i++) {rows[i] = start + i;}
1596: MatCreate(comm,mat);
1597: MatSetSizes(*mat,m,m,M,M);
1598: if (size == 1) {
1599: MatSetType(*mat,MATSEQDENSE);
1600: MatSeqDenseSetPreallocation(*mat,PETSC_NULL);
1601: } else {
1602: MatSetType(*mat,MATMPIAIJ);
1603: MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);
1604: }
1606: for (i=0; i<M; i++) {
1608: VecSet(in,0.0);
1609: VecSetValues(in,1,&i,&one,INSERT_VALUES);
1610: VecAssemblyBegin(in);
1611: VecAssemblyEnd(in);
1613: /* should fix, allowing user to choose side */
1614: PCApply(pc,in,out);
1615:
1616: VecGetArray(out,&array);
1617: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
1618: VecRestoreArray(out,&array);
1620: }
1621: PetscFree(rows);
1622: VecDestroy(out);
1623: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
1624: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
1625: return(0);
1626: }