Actual source code: precon.c
1: /*$Id: precon.c,v 1.213 2001/04/10 19:36:06 bsmith Exp bsmith $*/
2: /*
3: The PC (preconditioner) interface routines, callable by users.
4: */
5: #include src/sles/pc/pcimpl.h
7: /*@C
8: PCNullSpaceAttach - attaches a null space to a preconditioner object.
9: This null space will be removed from the resulting vector whenever
10: the preconditioner is applied.
12: Collective on PC
14: Input Parameters:
15: + pc - the preconditioner context
16: - nullsp - the null space object
18: Level: developer
20: Notes:
21: Overwrites any previous null space that may have been attached
23: .keywords: PC, destroy, null space
25: .seealso: PCCreate(), PCSetUp()
26: @*/
27: int PCNullSpaceAttach(PC pc,MatNullSpace nullsp)
28: {
35: if (pc->nullsp) {
36: MatNullSpaceDestroy(pc->nullsp);
37: }
38: pc->nullsp = nullsp;
39: PetscObjectReference((PetscObject)nullsp);
40: return(0);
41: }
43: /*@C
44: PCDestroy - Destroys PC context that was created with PCCreate().
46: Collective on PC
48: Input Parameter:
49: . pc - the preconditioner context
51: Level: developer
53: .keywords: PC, destroy
55: .seealso: PCCreate(), PCSetUp()
56: @*/
57: int PCDestroy(PC pc)
58: {
63: if (--pc->refct > 0) return(0);
65: /* if memory was published with AMS then destroy it */
66: PetscObjectDepublish(pc);
68: if (pc->ops->destroy) { (*pc->ops->destroy)(pc);}
69: if (pc->nullsp) {MatNullSpaceDestroy(pc->nullsp);}
70: if (pc->diagonalscaleright) {VecDestroy(pc->diagonalscaleright);}
71: if (pc->diagonalscaleleft) {VecDestroy(pc->diagonalscaleleft);}
73: PetscLogObjectDestroy(pc);
74: PetscHeaderDestroy(pc);
75: return(0);
76: }
78: /*@C
79: PCDiagonalScale - Indicates if the preconditioner applies an additional left and right
80: scaling as needed by certain time-stepping codes.
82: Collective on PC
84: Input Parameter:
85: . pc - the preconditioner context
87: Output Parameter:
88: . flag - PETSC_TRUE if it applies the scaling
90: Level: developer
92: Notes: If this returns PETSC_TRUE then the system solved via the Krylov method is
93: $ D M A D^{-1} y = D M b for left preconditioning or
94: $ D A M D^{-1} z = D b for right preconditioning
96: .keywords: PC
98: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScaleSet()
99: @*/
100: int PCDiagonalScale(PC pc,PetscTruth *flag)
101: {
104: *flag = pc->diagonalscale;
105: return(0);
106: }
108: /*@C
109: PCDiagonalScaleSet - Indicates the left scaling to use to apply an additional left and right
110: scaling as needed by certain time-stepping codes.
112: Collective on PC
114: Input Parameters:
115: + pc - the preconditioner context
116: - s - scaling vector
118: Level: intermediate
120: Notes: The system solved via the Krylov method is
121: $ D M A D^{-1} y = D M b for left preconditioning or
122: $ D A M D^{-1} z = D b for right preconditioning
124: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
126: .keywords: PC
128: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleRight(), PCDiagonalScale()
129: @*/
130: int PCDiagonalScaleSet(PC pc,Vec s)
131: {
137: pc->diagonalscale = PETSC_TRUE;
138: if (pc->diagonalscaleleft) {
139: VecDestroy(pc->diagonalscaleleft);
140: }
141: pc->diagonalscaleleft = s;
142: ierr = PetscObjectReference((PetscObject)s);
143: if (!pc->diagonalscaleright) {
144: VecDuplicate(s,&pc->diagonalscaleright);
145: }
146: VecCopy(s,pc->diagonalscaleright);
147: VecReciprocal(pc->diagonalscaleright);
148: return(0);
149: }
151: /*@C
152: PCDiagonalScaleLeft - Indicates the left scaling to use to apply an additional left and right
153: scaling as needed by certain time-stepping codes.
155: Collective on PC
157: Input Parameters:
158: + pc - the preconditioner context
159: . in - input vector
160: + out - scaled vector (maybe the same as in)
162: Level: intermediate
164: Notes: The system solved via the Krylov method is
165: $ D M A D^{-1} y = D M b for left preconditioning or
166: $ D A M D^{-1} z = D b for right preconditioning
168: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
170: If diagonal scaling is turned off and in is not out then in is copied to out
172: .keywords: PC
174: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleSet(), PCDiagonalScaleRight(), PCDiagonalScale()
175: @*/
176: int PCDiagonalScaleLeft(PC pc,Vec in,Vec out)
177: {
184: if (pc->diagonalscale) {
185: VecPointwiseMult(pc->diagonalscaleleft,in,out);
186: } else if (in != out) {
187: VecCopy(in,out);
188: }
189: return(0);
190: }
192: /*@C
193: PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes.
195: Collective on PC
197: Input Parameters:
198: + pc - the preconditioner context
199: . in - input vector
200: + out - scaled vector (maybe the same as in)
202: Level: intermediate
204: Notes: The system solved via the Krylov method is
205: $ D M A D^{-1} y = D M b for left preconditioning or
206: $ D A M D^{-1} z = D b for right preconditioning
208: PCDiagonalScaleLeft() scales a vector by D. PCDiagonalScaleRight() scales a vector by D^{-1}.
210: If diagonal scaling is turned off and in is not out then in is copied to out
212: .keywords: PC
214: .seealso: PCCreate(), PCSetUp(), PCDiagonalScaleLeft(), PCDiagonalScaleSet(), PCDiagonalScale()
215: @*/
216: int PCDiagonalScaleRight(PC pc,Vec in,Vec out)
217: {
224: if (pc->diagonalscale) {
225: VecPointwiseMult(pc->diagonalscaleright,in,out);
226: } else if (in != out) {
227: VecCopy(in,out);
228: }
229: return(0);
230: }
232: static int PCPublish_Petsc(PetscObject obj)
233: {
234: #if defined(PETSC_HAVE_AMS)
235: PC v = (PC) obj;
236: int ierr;
237: #endif
241: #if defined(PETSC_HAVE_AMS)
242: /* if it is already published then return */
243: if (v->amem >=0) return(0);
245: PetscObjectPublishBaseBegin(obj);
246: PetscObjectPublishBaseEnd(obj);
247: #endif
249: return(0);
250: }
252: /*@C
253: PCCreate - Creates a preconditioner context.
255: Collective on MPI_Comm
257: Input Parameter:
258: . comm - MPI communicator
260: Output Parameter:
261: . pc - location to put the preconditioner context
263: Notes:
264: The default preconditioner on one processor is PCILU with 0 fill on more
265: then one it is PCBJACOBI with ILU() on each processor.
267: Level: developer
269: .keywords: PC, create, context
271: .seealso: PCSetUp(), PCApply(), PCDestroy()
272: @*/
273: int PCCreate(MPI_Comm comm,PC *newpc)
274: {
275: PC pc;
279: *newpc = 0;
281: PetscHeaderCreate(pc,_p_PC,struct _PCOps,PC_COOKIE,-1,"PC",comm,PCDestroy,PCView);
282: PetscLogObjectCreate(pc);
283: pc->bops->publish = PCPublish_Petsc;
284: pc->vec = 0;
285: pc->mat = 0;
286: pc->setupcalled = 0;
287: pc->nullsp = 0;
288: pc->data = 0;
289: pc->diagonalscale = PETSC_FALSE;
290: pc->diagonalscaleleft = 0;
291: pc->diagonalscaleright = 0;
293: pc->ops->destroy = 0;
294: pc->ops->apply = 0;
295: pc->ops->applytranspose = 0;
296: pc->ops->applyBA = 0;
297: pc->ops->applyBAtranspose = 0;
298: pc->ops->applyrichardson = 0;
299: pc->ops->view = 0;
300: pc->ops->getfactoredmatrix = 0;
301: pc->ops->applysymmetricright = 0;
302: pc->ops->applysymmetricleft = 0;
303: pc->ops->setuponblocks = 0;
305: pc->modifysubmatrices = 0;
306: pc->modifysubmatricesP = 0;
307: *newpc = pc;
308: PetscPublishAll(pc);
309: return(0);
311: }
313: /* -------------------------------------------------------------------------------*/
315: /*@
316: PCApply - Applies the preconditioner to a vector.
318: Collective on PC and Vec
320: Input Parameters:
321: + pc - the preconditioner context
322: - x - input vector
324: Output Parameter:
325: . y - output vector
327: Level: developer
329: .keywords: PC, apply
331: .seealso: PCApplyTranspose(), PCApplyBAorAB()
332: @*/
333: int PCApply(PC pc,Vec x,Vec y)
334: {
335: int ierr;
341: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
343: if (pc->setupcalled < 2) {
344: PCSetUp(pc);
345: }
347: PetscLogEventBegin(PC_Apply,pc,x,y,0);
348: (*pc->ops->apply)(pc,x,y);
349: PetscLogEventEnd(PC_Apply,pc,x,y,0);
351: /* Remove null space from preconditioned vector y */
352: if (pc->nullsp) {
353: MatNullSpaceRemove(pc->nullsp,y,PETSC_NULL);
354: }
355: return(0);
356: }
358: /*@
359: PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector.
361: Collective on PC and Vec
363: Input Parameters:
364: + pc - the preconditioner context
365: - x - input vector
367: Output Parameter:
368: . y - output vector
370: Notes:
371: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
373: Level: developer
375: .keywords: PC, apply, symmetric, left
377: .seealso: PCApply(), PCApplySymmetricRight()
378: @*/
379: int PCApplySymmetricLeft(PC pc,Vec x,Vec y)
380: {
387: if (!pc->ops->applysymmetricleft) SETERRQ(1,"PC does not have left symmetric apply");
389: if (pc->setupcalled < 2) {
390: PCSetUp(pc);
391: }
393: PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);
394: (*pc->ops->applysymmetricleft)(pc,x,y);
395: PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);
396: return(0);
397: }
399: /*@
400: PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector.
402: Collective on PC and Vec
404: Input Parameters:
405: + pc - the preconditioner context
406: - x - input vector
408: Output Parameter:
409: . y - output vector
411: Level: developer
413: Notes:
414: Currently, this routine is implemented only for PCICC and PCJACOBI preconditioners.
416: .keywords: PC, apply, symmetric, right
418: .seealso: PCApply(), PCApplySymmetricLeft()
419: @*/
420: int PCApplySymmetricRight(PC pc,Vec x,Vec y)
421: {
428: if (!pc->ops->applysymmetricright) SETERRQ(1,"PC does not have left symmetric apply");
430: if (pc->setupcalled < 2) {
431: PCSetUp(pc);
432: }
434: PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);
435: (*pc->ops->applysymmetricright)(pc,x,y);
436: PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);
437: return(0);
438: }
440: /*@
441: PCApplyTranspose - Applies the transpose of preconditioner to a vector.
443: Collective on PC and Vec
445: Input Parameters:
446: + pc - the preconditioner context
447: - x - input vector
449: Output Parameter:
450: . y - output vector
452: Level: developer
454: .keywords: PC, apply, transpose
456: .seealso: PCApplyTranspose(), PCApplyBAorAB(), PCApplyBAorABTranspose()
457: @*/
458: int PCApplyTranspose(PC pc,Vec x,Vec y)
459: {
466: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
467: if (!pc->ops->applytranspose) SETERRQ(PETSC_ERR_SUP,"");
469: if (pc->setupcalled < 2) {
470: PCSetUp(pc);
471: }
473: PetscLogEventBegin(PC_Apply,pc,x,y,0);
474: (*pc->ops->applytranspose)(pc,x,y);
475: PetscLogEventEnd(PC_Apply,pc,x,y,0);
476: return(0);
477: }
479: /*@
480: PCApplyBAorAB - Applies the preconditioner and operator to a vector.
482: Collective on PC and Vec
484: Input Parameters:
485: + pc - the preconditioner context
486: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
487: . x - input vector
488: - work - work vector
490: Output Parameter:
491: . y - output vector
493: Level: developer
495: .keywords: PC, apply, operator
497: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorABTranspose()
498: @*/
499: int PCApplyBAorAB(PC pc,PCSide side,Vec x,Vec y,Vec work)
500: {
501: int ierr;
508: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
509: if (side != PC_LEFT && side != PC_SYMMETRIC && side != PC_RIGHT) {
510: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right, left, or symmetric");
511: }
512: if (pc->diagonalscale && side == PC_SYMMETRIC) {
513: SETERRQ(1,"Cannot include diagonal scaling with symmetric preconditioner application");
514: }
516: if (pc->setupcalled < 2) {
517: PCSetUp(pc);
518: }
520: if (pc->diagonalscale) {
521: if (pc->ops->applyBA) {
522: Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */
523: VecDuplicate(x,&work2);
524: PCDiagonalScaleRight(pc,x,work2);
525: (*pc->ops->applyBA)(pc,side,work2,y,work);
526: /* Remove null space from preconditioned vector y */
527: if (pc->nullsp) {
528: MatNullSpaceRemove(pc->nullsp,y,PETSC_NULL);
529: }
530: PCDiagonalScaleLeft(pc,y,y);
531: VecDestroy(work2);
532: } else if (side == PC_RIGHT) {
533: PCDiagonalScaleRight(pc,x,y);
534: PCApply(pc,y,work);
535: MatMult(pc->mat,work,y);
536: PCDiagonalScaleLeft(pc,y,y);
537: } else if (side == PC_LEFT) {
538: PCDiagonalScaleRight(pc,x,y);
539: MatMult(pc->mat,y,work);
540: PCApply(pc,work,y);
541: PCDiagonalScaleLeft(pc,y,y);
542: } else if (side == PC_SYMMETRIC) {
543: SETERRQ(1,"Cannot provide diagonal scaling with symmetric application of preconditioner");
544: }
545: } else {
546: if (pc->ops->applyBA) {
547: (*pc->ops->applyBA)(pc,side,x,y,work);
548: /* Remove null space from preconditioned vector y */
549: if (pc->nullsp) {
550: MatNullSpaceRemove(pc->nullsp,y,PETSC_NULL);
551: }
552: } else if (side == PC_RIGHT) {
553: PCApply(pc,x,work);
554: MatMult(pc->mat,work,y);
555: } else if (side == PC_LEFT) {
556: MatMult(pc->mat,x,work);
557: PCApply(pc,work,y);
558: } else if (side == PC_SYMMETRIC) {
559: /* There's an extra copy here; maybe should provide 2 work vectors instead? */
560: PCApplySymmetricRight(pc,x,work);
561: MatMult(pc->mat,work,y);
562: VecCopy(y,work);
563: PCApplySymmetricLeft(pc,work,y);
564: }
565: }
566: return(0);
567: }
569: /*@
570: PCApplyBAorABTranspose - Applies the transpose of the preconditioner
571: and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning,
572: not tr(B*A) = tr(A)*tr(B).
574: Collective on PC and Vec
576: Input Parameters:
577: + pc - the preconditioner context
578: . side - indicates the preconditioner side, one of PC_LEFT, PC_RIGHT, or PC_SYMMETRIC
579: . x - input vector
580: - work - work vector
582: Output Parameter:
583: . y - output vector
585: Level: developer
587: .keywords: PC, apply, operator, transpose
589: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
590: @*/
591: int PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
592: {
600: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
601: if (pc->ops->applyBAtranspose) {
602: (*pc->ops->applyBAtranspose)(pc,side,x,y,work);
603: return(0);
604: }
605: if (side != PC_LEFT && side != PC_RIGHT) {
606: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
607: }
609: if (pc->setupcalled < 2) {
610: PCSetUp(pc);
611: }
613: if (side == PC_RIGHT) {
614: PCApplyTranspose(pc,x,work);
615: MatMultTranspose(pc->mat,work,y);
616: } else if (side == PC_LEFT) {
617: MatMultTranspose(pc->mat,x,work);
618: PCApplyTranspose(pc,work,y);
619: }
620: /* add support for PC_SYMMETRIC */
621: return(0); /* actually will never get here */
622: }
624: /* -------------------------------------------------------------------------------*/
626: /*@
627: PCApplyRichardsonExists - Determines whether a particular preconditioner has a
628: built-in fast application of Richardson's method.
630: Not Collective
632: Input Parameter:
633: . pc - the preconditioner
635: Output Parameter:
636: . exists - PETSC_TRUE or PETSC_FALSE
638: Level: developer
640: .keywords: PC, apply, Richardson, exists
642: .seealso: PCApplyRichardson()
643: @*/
644: int PCApplyRichardsonExists(PC pc,PetscTruth *exists)
645: {
649: if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
650: else *exists = PETSC_FALSE;
651: return(0);
652: }
654: /*@
655: PCApplyRichardson - Applies several steps of Richardson iteration with
656: the particular preconditioner. This routine is usually used by the
657: Krylov solvers and not the application code directly.
659: Collective on PC
661: Input Parameters:
662: + pc - the preconditioner context
663: . x - the initial guess
664: . w - one work vector
665: - its - the number of iterations to apply.
667: Output Parameter:
668: . y - the solution
670: Notes:
671: Most preconditioners do not support this function. Use the command
672: PCApplyRichardsonExists() to determine if one does.
674: Level: developer
676: .keywords: PC, apply, Richardson
678: .seealso: PCApplyRichardsonExists()
679: @*/
680: int PCApplyRichardson(PC pc,Vec x,Vec y,Vec w,int its)
681: {
689: if (!pc->ops->applyrichardson) SETERRQ(PETSC_ERR_SUP,"");
691: if (pc->setupcalled < 2) {
692: PCSetUp(pc);
693: }
695: (*pc->ops->applyrichardson)(pc,x,y,w,its);
696: return(0);
697: }
699: /*
700: a setupcall of 0 indicates never setup,
701: 1 needs to be resetup,
702: 2 does not need any changes.
703: */
704: /*@
705: PCSetUp - Prepares for the use of a preconditioner.
707: Collective on PC
709: Input Parameter:
710: . pc - the preconditioner context
712: Level: developer
714: .keywords: PC, setup
716: .seealso: PCCreate(), PCApply(), PCDestroy()
717: @*/
718: int PCSetUp(PC pc)
719: {
725: if (pc->setupcalled > 1) {
726: PetscLogInfo(pc,"PCSetUp:Setting PC with identical preconditionern");
727: } else if (pc->setupcalled == 0) {
728: PetscLogInfo(pc,"PCSetUp:Setting up new PCn");
729: } else if (pc->flag == SAME_NONZERO_PATTERN) {
730: PetscLogInfo(pc,"PCSetUp:Setting up PC with same nonzero patternn");
731: } else {
732: PetscLogInfo(pc,"PCSetUp:Setting up PC with different nonzero patternn");
733: }
734: if (pc->setupcalled > 1) return(0);
735: PetscLogEventBegin(PC_SetUp,pc,0,0,0);
736: if (!pc->vec) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Vector must be set first");}
737: if (!pc->mat) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
738: if (!pc->type_name) {
739: int size;
741: MPI_Comm_size(pc->comm,&size);
742: if (size == 1) {
743: PCSetType(pc,PCILU);
744: } else {
745: PCSetType(pc,PCBJACOBI);
746: }
747: }
748: if (pc->ops->setup) {
749: (*pc->ops->setup)(pc);
750: }
751: pc->setupcalled = 2;
752: if (pc->nullsp) {
753: PetscTruth test;
754: PetscOptionsHasName(pc->prefix,"-pc_test_null_space",&test);
755: if (test) {
756: MatNullSpaceTest(pc->nullsp,pc->mat);
757: }
758: }
759: PetscLogEventEnd(PC_SetUp,pc,0,0,0);
760: return(0);
761: }
763: /*@
764: PCSetUpOnBlocks - Sets up the preconditioner for each block in
765: the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
766: methods.
768: Collective on PC
770: Input Parameters:
771: . pc - the preconditioner context
773: Level: developer
775: .keywords: PC, setup, blocks
777: .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
778: @*/
779: int PCSetUpOnBlocks(PC pc)
780: {
785: if (!pc->ops->setuponblocks) return(0);
786: PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);
787: (*pc->ops->setuponblocks)(pc);
788: PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);
789: return(0);
790: }
792: /*@
793: PCSetModifySubMatrices - Sets a user-defined routine for modifying the
794: submatrices that arise within certain subdomain-based preconditioners.
795: The basic submatrices are extracted from the preconditioner matrix as
796: usual; the user can then alter these (for example, to set different boundary
797: conditions for each submatrix) before they are used for the local solves.
799: Collective on PC
801: Input Parameters:
802: + pc - the preconditioner context
803: . func - routine for modifying the submatrices
804: - ctx - optional user-defined context (may be null)
806: Calling sequence of func:
807: $ func (PC pc,int nsub,IS *row,IS *col,Mat *submat,void *ctx);
809: . row - an array of index sets that contain the global row numbers
810: that comprise each local submatrix
811: . col - an array of index sets that contain the global column numbers
812: that comprise each local submatrix
813: . submat - array of local submatrices
814: - ctx - optional user-defined context for private data for the
815: user-defined func routine (may be null)
817: Notes:
818: PCSetModifySubMatrices() MUST be called before SLESSetUp() and
819: SLESSolve().
821: A routine set by PCSetModifySubMatrices() is currently called within
822: the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
823: preconditioners. All other preconditioners ignore this routine.
825: Level: advanced
827: .keywords: PC, set, modify, submatrices
829: .seealso: PCModifySubMatrices()
830: @*/
831: int PCSetModifySubMatrices(PC pc,int(*func)(PC,int,IS*,IS*,Mat*,void*),void *ctx)
832: {
835: pc->modifysubmatrices = func;
836: pc->modifysubmatricesP = ctx;
837: return(0);
838: }
840: /*@
841: PCModifySubMatrices - Calls an optional user-defined routine within
842: certain preconditioners if one has been set with PCSetModifySubMarices().
844: Collective on PC
846: Input Parameters:
847: + pc - the preconditioner context
848: . nsub - the number of local submatrices
849: . row - an array of index sets that contain the global row numbers
850: that comprise each local submatrix
851: . col - an array of index sets that contain the global column numbers
852: that comprise each local submatrix
853: . submat - array of local submatrices
854: - ctx - optional user-defined context for private data for the
855: user-defined routine (may be null)
857: Output Parameter:
858: . submat - array of local submatrices (the entries of which may
859: have been modified)
861: Notes:
862: The user should NOT generally call this routine, as it will
863: automatically be called within certain preconditioners (currently
864: block Jacobi, additive Schwarz) if set.
866: The basic submatrices are extracted from the preconditioner matrix
867: as usual; the user can then alter these (for example, to set different
868: boundary conditions for each submatrix) before they are used for the
869: local solves.
871: Level: developer
873: .keywords: PC, modify, submatrices
875: .seealso: PCSetModifySubMatrices()
876: @*/
877: int PCModifySubMatrices(PC pc,int nsub,IS *row,IS *col,Mat *submat,void *ctx)
878: {
882: if (!pc->modifysubmatrices) return(0);
883: PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);
884: (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);
885: PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);
886: return(0);
887: }
889: /*@
890: PCSetOperators - Sets the matrix associated with the linear system and
891: a (possibly) different one associated with the preconditioner.
893: Collective on PC and Mat
895: Input Parameters:
896: + pc - the preconditioner context
897: . Amat - the matrix associated with the linear system
898: . Pmat - the matrix to be used in constructing the preconditioner, usually
899: the same as Amat.
900: - flag - flag indicating information about the preconditioner matrix structure
901: during successive linear solves. This flag is ignored the first time a
902: linear system is solved, and thus is irrelevant when solving just one linear
903: system.
905: Notes:
906: The flag can be used to eliminate unnecessary work in the preconditioner
907: during the repeated solution of linear systems of the same size. The
908: available options are
909: + SAME_PRECONDITIONER -
910: Pmat is identical during successive linear solves.
911: This option is intended for folks who are using
912: different Amat and Pmat matrices and wish to reuse the
913: same preconditioner matrix. For example, this option
914: saves work by not recomputing incomplete factorization
915: for ILU/ICC preconditioners.
916: . SAME_NONZERO_PATTERN -
917: Pmat has the same nonzero structure during
918: successive linear solves.
919: - DIFFERENT_NONZERO_PATTERN -
920: Pmat does not have the same nonzero structure.
922: Caution:
923: If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
924: and does not check the structure of the matrix. If you erroneously
925: claim that the structure is the same when it actually is not, the new
926: preconditioner will not function correctly. Thus, use this optimization
927: feature carefully!
929: If in doubt about whether your preconditioner matrix has changed
930: structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
932: More Notes about Repeated Solution of Linear Systems:
933: PETSc does NOT reset the matrix entries of either Amat or Pmat
934: to zero after a linear solve; the user is completely responsible for
935: matrix assembly. See the routine MatZeroEntries() if desiring to
936: zero all elements of a matrix.
938: Level: developer
940: .keywords: PC, set, operators, matrix, linear system
942: .seealso: PCGetOperators(), MatZeroEntries()
943: @*/
944: int PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
945: {
946: int ierr;
947: PetscTruth isbjacobi,isrowbs;
954: /*
955: BlockSolve95 cannot use default BJacobi preconditioning
956: */
957: PetscTypeCompare((PetscObject)Amat,MATMPIROWBS,&isrowbs);
958: if (isrowbs) {
959: PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);
960: if (isbjacobi) {
961: PCSetType(pc,PCILU);
962: PetscLogInfo(pc,"PCSetOperators:Switching default PC to PCILU since BS95 doesn't support PCBJACOBIn");
963: }
964: }
966: pc->mat = Amat;
967: pc->pmat = Pmat;
968: if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
969: pc->setupcalled = 1;
970: }
971: pc->flag = flag;
972: return(0);
973: }
975: /*@C
976: PCGetOperators - Gets the matrix associated with the linear system and
977: possibly a different one associated with the preconditioner.
979: Not collective, though parallel Mats are returned if the PC is parallel
981: Input Parameter:
982: . pc - the preconditioner context
984: Output Parameters:
985: + mat - the matrix associated with the linear system
986: . pmat - matrix associated with the preconditioner, usually the same
987: as mat.
988: - flag - flag indicating information about the preconditioner
989: matrix structure. See PCSetOperators() for details.
991: Level: developer
993: .keywords: PC, get, operators, matrix, linear system
995: .seealso: PCSetOperators()
996: @*/
997: int PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
998: {
1001: if (mat) *mat = pc->mat;
1002: if (pmat) *pmat = pc->pmat;
1003: if (flag) *flag = pc->flag;
1004: return(0);
1005: }
1007: /*@
1008: PCSetVector - Sets a vector associated with the preconditioner.
1010: Collective on PC and Vec
1012: Input Parameters:
1013: + pc - the preconditioner context
1014: - vec - the vector
1016: Notes:
1017: The vector must be set so that the preconditioner knows what type
1018: of vector to allocate if necessary.
1020: Level: developer
1022: .keywords: PC, set, vector
1024: .seealso: PCGetVector()
1026: @*/
1027: int PCSetVector(PC pc,Vec vec)
1028: {
1033: pc->vec = vec;
1034: return(0);
1035: }
1037: /*@
1038: PCGetVector - Gets a vector associated with the preconditioner; if the
1039: vector was not get set it will return a 0 pointer.
1041: Not collective, but vector is shared by all processors that share the PC
1043: Input Parameter:
1044: . pc - the preconditioner context
1046: Output Parameter:
1047: . vec - the vector
1049: Level: developer
1051: .keywords: PC, get, vector
1053: .seealso: PCSetVector()
1055: @*/
1056: int PCGetVector(PC pc,Vec *vec)
1057: {
1060: *vec = pc->vec;
1061: return(0);
1062: }
1064: /*@C
1065: PCGetFactoredMatrix - Gets the factored matrix from the
1066: preconditioner context. This routine is valid only for the LU,
1067: incomplete LU, Cholesky, and incomplete Cholesky methods.
1069: Not Collective on PC though Mat is parallel if PC is parallel
1071: Input Parameters:
1072: . pc - the preconditioner context
1074: Output parameters:
1075: . mat - the factored matrix
1077: Level: advanced
1079: .keywords: PC, get, factored, matrix
1080: @*/
1081: int PCGetFactoredMatrix(PC pc,Mat *mat)
1082: {
1087: if (pc->ops->getfactoredmatrix) {
1088: (*pc->ops->getfactoredmatrix)(pc,mat);
1089: }
1090: return(0);
1091: }
1093: /*@C
1094: PCSetOptionsPrefix - Sets the prefix used for searching for all
1095: PC options in the database.
1097: Collective on PC
1099: Input Parameters:
1100: + pc - the preconditioner context
1101: - prefix - the prefix string to prepend to all PC option requests
1103: Notes:
1104: A hyphen (-) must NOT be given at the beginning of the prefix name.
1105: The first character of all runtime options is AUTOMATICALLY the
1106: hyphen.
1108: Level: advanced
1110: .keywords: PC, set, options, prefix, database
1112: .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1113: @*/
1114: int PCSetOptionsPrefix(PC pc,char *prefix)
1115: {
1120: PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);
1121: return(0);
1122: }
1124: /*@C
1125: PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1126: PC options in the database.
1128: Collective on PC
1130: Input Parameters:
1131: + pc - the preconditioner context
1132: - prefix - the prefix string to prepend to all PC option requests
1134: Notes:
1135: A hyphen (-) must NOT be given at the beginning of the prefix name.
1136: The first character of all runtime options is AUTOMATICALLY the
1137: hyphen.
1139: Level: advanced
1141: .keywords: PC, append, options, prefix, database
1143: .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1144: @*/
1145: int PCAppendOptionsPrefix(PC pc,char *prefix)
1146: {
1151: PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);
1152: return(0);
1153: }
1155: /*@C
1156: PCGetOptionsPrefix - Gets the prefix used for searching for all
1157: PC options in the database.
1159: Not Collective
1161: Input Parameters:
1162: . pc - the preconditioner context
1164: Output Parameters:
1165: . prefix - pointer to the prefix string used, is returned
1167: Notes: On the fortran side, the user should pass in a string 'prifix' of
1168: sufficient length to hold the prefix.
1170: Level: advanced
1172: .keywords: PC, get, options, prefix, database
1174: .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1175: @*/
1176: int PCGetOptionsPrefix(PC pc,char **prefix)
1177: {
1182: PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);
1183: return(0);
1184: }
1186: /*@
1187: PCPreSolve - Optional pre-solve phase, intended for any
1188: preconditioner-specific actions that must be performed before
1189: the iterative solve itself.
1191: Collective on PC
1193: Input Parameters:
1194: + pc - the preconditioner context
1195: - ksp - the Krylov subspace context
1197: Level: developer
1199: Sample of Usage:
1200: .vb
1201: PCPreSolve(pc,ksp);
1202: KSPSolve(ksp,its);
1203: PCPostSolve(pc,ksp);
1204: .ve
1206: Notes:
1207: The pre-solve phase is distinct from the PCSetUp() phase.
1209: SLESSolve() calls this directly, so is rarely called by the user.
1211: .keywords: PC, pre-solve
1213: .seealso: PCPostSolve()
1214: @*/
1215: int PCPreSolve(PC pc,KSP ksp)
1216: {
1218: Vec x,rhs;
1219: Mat A,B;
1224: KSPGetSolution(ksp,&x);
1225: KSPGetRhs(ksp,&rhs);
1226: /*
1227: Scale the system and have the matrices use the scaled form
1228: only if the two matrices are actually the same (and hence
1229: have the same scaling
1230: */
1231: PCGetOperators(pc,&A,&B,PETSC_NULL);
1232: if (A == B) {
1233: MatScaleSystem(pc->mat,x,rhs);
1234: MatUseScaledForm(pc->mat,PETSC_TRUE);
1235: }
1237: if (pc->ops->presolve) {
1238: (*pc->ops->presolve)(pc,ksp,x,rhs);
1239: }
1240: return(0);
1241: }
1243: /*@
1244: PCPostSolve - Optional post-solve phase, intended for any
1245: preconditioner-specific actions that must be performed after
1246: the iterative solve itself.
1248: Collective on PC
1250: Input Parameters:
1251: + pc - the preconditioner context
1252: - ksp - the Krylov subspace context
1254: Sample of Usage:
1255: .vb
1256: PCPreSolve(pc,ksp);
1257: KSPSolve(ksp,its);
1258: PCPostSolve(pc,ksp);
1259: .ve
1261: Note:
1262: SLESSolve() calls this routine directly, so it is rarely called by the user.
1264: Level: developer
1266: .keywords: PC, post-solve
1268: .seealso: PCPreSolve(), SLESSolve()
1269: @*/
1270: int PCPostSolve(PC pc,KSP ksp)
1271: {
1273: Vec x,rhs;
1274: Mat A,B;
1278: KSPGetSolution(ksp,&x);
1279: KSPGetRhs(ksp,&rhs);
1280: if (pc->ops->postsolve) {
1281: (*pc->ops->postsolve)(pc,ksp,x,rhs);
1282: }
1284: /*
1285: Scale the system and have the matrices use the scaled form
1286: only if the two matrices are actually the same (and hence
1287: have the same scaling
1288: */
1289: PCGetOperators(pc,&A,&B,PETSC_NULL);
1290: if (A == B) {
1291: MatUnScaleSystem(pc->mat,x,rhs);
1292: MatUseScaledForm(pc->mat,PETSC_FALSE);
1293: }
1294: return(0);
1295: }
1297: /*@C
1298: PCView - Prints the PC data structure.
1300: Collective on PC
1302: Input Parameters:
1303: + PC - the PC context
1304: - viewer - optional visualization context
1306: Note:
1307: The available visualization contexts include
1308: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
1309: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1310: output where only the first processor opens
1311: the file. All other processors send their
1312: data to the first processor to print.
1314: The user can open an alternative visualization contexts with
1315: PetscViewerASCIIOpen() (output to a specified file).
1317: Level: developer
1319: .keywords: PC, view
1321: .seealso: KSPView(), PetscViewerASCIIOpen()
1322: @*/
1323: int PCView(PC pc,PetscViewer viewer)
1324: {
1325: PCType cstr;
1326: int ierr;
1327: PetscTruth mat_exists,isascii,isstring;
1328: PetscViewerFormat format;
1332: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(pc->comm);
1336: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
1337: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
1338: if (isascii) {
1339: PetscViewerGetFormat(viewer,&format);
1340: PetscViewerASCIIPrintf(viewer,"PC Object:n");
1341: PCGetType(pc,&cstr);
1342: if (cstr) {
1343: PetscViewerASCIIPrintf(viewer," type: %sn",cstr);
1344: } else {
1345: PetscViewerASCIIPrintf(viewer," type: not yet setn");
1346: }
1347: if (pc->ops->view) {
1348: PetscViewerASCIIPushTab(viewer);
1349: (*pc->ops->view)(pc,viewer);
1350: PetscViewerASCIIPopTab(viewer);
1351: }
1352: PetscObjectExists((PetscObject)pc->mat,&mat_exists);
1353: if (mat_exists) {
1354: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1355: if (pc->pmat == pc->mat) {
1356: PetscViewerASCIIPrintf(viewer," linear system matrix = precond matrix:n");
1357: PetscViewerASCIIPushTab(viewer);
1358: MatView(pc->mat,viewer);
1359: PetscViewerASCIIPopTab(viewer);
1360: } else {
1361: PetscObjectExists((PetscObject)pc->pmat,&mat_exists);
1362: if (mat_exists) {
1363: PetscViewerASCIIPrintf(viewer," linear system matrix followed by preconditioner matrix:n");
1364: } else {
1365: PetscViewerASCIIPrintf(viewer," linear system matrix:n");
1366: }
1367: PetscViewerASCIIPushTab(viewer);
1368: MatView(pc->mat,viewer);
1369: if (mat_exists) {MatView(pc->pmat,viewer);}
1370: PetscViewerASCIIPopTab(viewer);
1371: }
1372: PetscViewerPopFormat(viewer);
1373: }
1374: } else if (isstring) {
1375: PCGetType(pc,&cstr);
1376: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
1377: if (pc->ops->view) {(*pc->ops->view)(pc,viewer);}
1378: } else {
1379: SETERRQ1(1,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1380: }
1381: return(0);
1382: }
1384: /*MC
1385: PCRegisterDynamic - Adds a method to the preconditioner package.
1387: Synopsis:
1388: int PCRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(PC))
1390: Not collective
1392: Input Parameters:
1393: + name_solver - name of a new user-defined solver
1394: . path - path (either absolute or relative) the library containing this solver
1395: . name_create - name of routine to create method context
1396: - routine_create - routine to create method context
1398: Notes:
1399: PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
1401: If dynamic libraries are used, then the fourth input argument (routine_create)
1402: is ignored.
1404: Sample usage:
1405: .vb
1406: PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
1407: "MySolverCreate",MySolverCreate);
1408: .ve
1410: Then, your solver can be chosen with the procedural interface via
1411: $ PCSetType(pc,"my_solver")
1412: or at runtime via the option
1413: $ -pc_type my_solver
1415: Level: advanced
1417: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT}, or ${any environmental variable}
1418: occuring in pathname will be replaced with appropriate values.
1420: .keywords: PC, register
1422: .seealso: PCRegisterAll(), PCRegisterDestroy(), PCRegister()
1423: M*/
1425: int PCRegister(char *sname,char *path,char *name,int (*function)(PC))
1426: {
1427: int ierr;
1428: char fullname[256];
1432: PetscFListConcat(path,name,fullname);
1433: PetscFListAdd(&PCList,sname,fullname,(void (*)())function);
1434: return(0);
1435: }
1437: /*@
1438: PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1440: Collective on PC
1442: Input Parameter:
1443: . pc - the preconditioner object
1445: Output Parameter:
1446: . mat - the explict preconditioned operator
1448: Notes:
1449: This computation is done by applying the operators to columns of the
1450: identity matrix.
1452: Currently, this routine uses a dense matrix format when 1 processor
1453: is used and a sparse format otherwise. This routine is costly in general,
1454: and is recommended for use only with relatively small systems.
1456: Level: advanced
1457:
1458: .keywords: PC, compute, explicit, operator
1460: @*/
1461: int PCComputeExplicitOperator(PC pc,Mat *mat)
1462: {
1463: Vec in,out;
1464: int ierr,i,M,m,size,*rows,start,end;
1465: MPI_Comm comm;
1466: Scalar *array,zero = 0.0,one = 1.0;
1472: comm = pc->comm;
1473: MPI_Comm_size(comm,&size);
1475: PCGetVector(pc,&in);
1476: VecDuplicate(in,&out);
1477: VecGetOwnershipRange(in,&start,&end);
1478: VecGetSize(in,&M);
1479: VecGetLocalSize(in,&m);
1480: PetscMalloc((m+1)*sizeof(int),&rows);
1481: for (i=0; i<m; i++) {rows[i] = start + i;}
1483: if (size == 1) {
1484: MatCreateSeqDense(comm,M,M,PETSC_NULL,mat);
1485: } else {
1486: MatCreateMPIAIJ(comm,m,m,M,M,0,0,0,0,mat);
1487: }
1489: for (i=0; i<M; i++) {
1491: VecSet(&zero,in);
1492: VecSetValues(in,1,&i,&one,INSERT_VALUES);
1493: VecAssemblyBegin(in);
1494: VecAssemblyEnd(in);
1496: PCApply(pc,in,out);
1497:
1498: VecGetArray(out,&array);
1499: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
1500: VecRestoreArray(out,&array);
1502: }
1503: PetscFree(rows);
1504: VecDestroy(out);
1505: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
1506: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
1507: return(0);
1508: }