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;
10: PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
11: PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, 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: MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ICC,&flg1);
28: MatGetFactorAvailable(pc->pmat,"petsc",MAT_FACTOR_ILU,&flg2);
29: MatIsSymmetricKnown(pc->pmat,&set,&flg3);
30: if (flg1 && (!flg2 || (set && flg3))) {
31: *type = PCICC;
32: } else if (flg2) {
33: *type = PCILU;
34: } else if (f) { /* likely is a parallel matrix run on one processor */
35: *type = PCBJACOBI;
36: } else {
37: *type = PCNONE;
38: }
39: } else {
40: if (f) {
41: *type = PCBJACOBI;
42: } else {
43: *type = PCNONE;
44: }
45: }
46: } else {
47: if (size == 1) {
48: *type = PCILU;
49: } else {
50: *type = PCBJACOBI;
51: }
52: }
53: return(0);
54: }
58: /*@
59: 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 for sparse matrices is PCILU or PCICC with 0 fill on one process and block Jacobi with PCILU or ICC
282: in parallel. For dense matrices it is always PCNONE.
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");
356: PetscLogEventBegin(PC_Apply,pc,x,y,0);
357: (*pc->ops->apply)(pc,x,y);
358: PetscLogEventEnd(PC_Apply,pc,x,y,0);
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");
398: PetscLogEventBegin(PC_ApplySymmetricLeft,pc,x,y,0);
399: (*pc->ops->applysymmetricleft)(pc,x,y);
400: PetscLogEventEnd(PC_ApplySymmetricLeft,pc,x,y,0);
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");
440: PetscLogEventBegin(PC_ApplySymmetricRight,pc,x,y,0);
441: (*pc->ops->applysymmetricright)(pc,x,y);
442: PetscLogEventEnd(PC_ApplySymmetricRight,pc,x,y,0);
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");
479: PetscLogEventBegin(PC_Apply,pc,x,y,0);
480: (*pc->ops->applytranspose)(pc,x,y);
481: PetscLogEventEnd(PC_Apply,pc,x,y,0);
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. y = B*A*x or y = A*B*x.
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
617: Notes: this routine is used internally so that the same Krylov code can be used to solve A x = b and A' x = b, with a preconditioner
618: defined by B'. This is why this has the funny form that it computes tr(B) * tr(A)
619:
620: Level: developer
622: .keywords: PC, apply, operator, transpose
624: .seealso: PCApply(), PCApplyTranspose(), PCApplyBAorAB()
625: @*/
626: PetscErrorCode PCApplyBAorABTranspose(PC pc,PCSide side,Vec x,Vec y,Vec work)
627: {
635: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
636: if (pc->ops->applyBAtranspose) {
637: (*pc->ops->applyBAtranspose)(pc,side,x,y,work);
638: return(0);
639: }
640: if (side != PC_LEFT && side != PC_RIGHT) {
641: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Side must be right or left");
642: }
644: if (pc->setupcalled < 2) {
645: PCSetUp(pc);
646: }
648: if (side == PC_RIGHT) {
649: PCApplyTranspose(pc,x,work);
650: MatMultTranspose(pc->mat,work,y);
651: } else if (side == PC_LEFT) {
652: MatMultTranspose(pc->mat,x,work);
653: PCApplyTranspose(pc,work,y);
654: }
655: /* add support for PC_SYMMETRIC */
656: return(0); /* actually will never get here */
657: }
659: /* -------------------------------------------------------------------------------*/
663: /*@
664: PCApplyRichardsonExists - Determines whether a particular preconditioner has a
665: built-in fast application of Richardson's method.
667: Not Collective
669: Input Parameter:
670: . pc - the preconditioner
672: Output Parameter:
673: . exists - PETSC_TRUE or PETSC_FALSE
675: Level: developer
677: .keywords: PC, apply, Richardson, exists
679: .seealso: PCApplyRichardson()
680: @*/
681: PetscErrorCode PCApplyRichardsonExists(PC pc,PetscTruth *exists)
682: {
686: if (pc->ops->applyrichardson) *exists = PETSC_TRUE;
687: else *exists = PETSC_FALSE;
688: return(0);
689: }
693: /*@
694: PCApplyRichardson - Applies several steps of Richardson iteration with
695: the particular preconditioner. This routine is usually used by the
696: Krylov solvers and not the application code directly.
698: Collective on PC
700: Input Parameters:
701: + pc - the preconditioner context
702: . x - the initial guess
703: . w - one work vector
704: . rtol - relative decrease in residual norm convergence criteria
705: . abstol - absolute residual norm convergence criteria
706: . dtol - divergence residual norm increase criteria
707: - its - the number of iterations to apply.
709: Output Parameter:
710: + outits - number of iterations actually used (for SOR this always equals its)
711: . reason - the reason the apply terminated
712: - y - the solution
714: Notes:
715: Most preconditioners do not support this function. Use the command
716: PCApplyRichardsonExists() to determine if one does.
718: Except for the multigrid PC this routine ignores the convergence tolerances
719: and always runs for the number of iterations
720:
721: Level: developer
723: .keywords: PC, apply, Richardson
725: .seealso: PCApplyRichardsonExists()
726: @*/
727: PetscErrorCode PCApplyRichardson(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscInt *outits,PCRichardsonConvergedReason *reason)
728: {
736: if (x == y) SETERRQ(PETSC_ERR_ARG_IDN,"x and y must be different vectors");
737: if (pc->setupcalled < 2) {
738: PCSetUp(pc);
739: }
740: if (!pc->ops->applyrichardson) SETERRQ(PETSC_ERR_SUP,"PC does not have apply richardson");
741: (*pc->ops->applyrichardson)(pc,x,y,w,rtol,abstol,dtol,its,outits,reason);
742: return(0);
743: }
745: /*
746: a setupcall of 0 indicates never setup,
747: 1 needs to be resetup,
748: 2 does not need any changes.
749: */
752: /*@
753: PCSetUp - Prepares for the use of a preconditioner.
755: Collective on PC
757: Input Parameter:
758: . pc - the preconditioner context
760: Level: developer
762: .keywords: PC, setup
764: .seealso: PCCreate(), PCApply(), PCDestroy()
765: @*/
766: PetscErrorCode PCSetUp(PC pc)
767: {
769: const char *def;
774: if (!pc->mat) {SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be set first");}
776: if (pc->setupcalled > 1) {
777: PetscInfo(pc,"Setting PC with identical preconditioner\n");
778: return(0);
779: } else if (!pc->setupcalled) {
780: PetscInfo(pc,"Setting up new PC\n");
781: } else if (pc->flag == SAME_NONZERO_PATTERN) {
782: PetscInfo(pc,"Setting up PC with same nonzero pattern\n");
783: } else {
784: PetscInfo(pc,"Setting up PC with different nonzero pattern\n");
785: }
787: if (!((PetscObject)pc)->type_name) {
788: PCGetDefaultType_Private(pc,&def);
789: PCSetType(pc,def);
790: }
792: PetscLogEventBegin(PC_SetUp,pc,0,0,0);
793: if (pc->ops->setup) {
794: (*pc->ops->setup)(pc);
795: }
796: pc->setupcalled = 2;
797: PetscLogEventEnd(PC_SetUp,pc,0,0,0);
798: return(0);
799: }
803: /*@
804: PCSetUpOnBlocks - Sets up the preconditioner for each block in
805: the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
806: methods.
808: Collective on PC
810: Input Parameters:
811: . pc - the preconditioner context
813: Level: developer
815: .keywords: PC, setup, blocks
817: .seealso: PCCreate(), PCApply(), PCDestroy(), PCSetUp()
818: @*/
819: PetscErrorCode PCSetUpOnBlocks(PC pc)
820: {
825: if (!pc->ops->setuponblocks) return(0);
826: PetscLogEventBegin(PC_SetUpOnBlocks,pc,0,0,0);
827: (*pc->ops->setuponblocks)(pc);
828: PetscLogEventEnd(PC_SetUpOnBlocks,pc,0,0,0);
829: return(0);
830: }
834: /*@C
835: PCSetModifySubMatrices - Sets a user-defined routine for modifying the
836: submatrices that arise within certain subdomain-based preconditioners.
837: The basic submatrices are extracted from the preconditioner matrix as
838: usual; the user can then alter these (for example, to set different boundary
839: conditions for each submatrix) before they are used for the local solves.
841: Collective on PC
843: Input Parameters:
844: + pc - the preconditioner context
845: . func - routine for modifying the submatrices
846: - ctx - optional user-defined context (may be null)
848: Calling sequence of func:
849: $ func (PC pc,PetscInt nsub,IS *row,IS *col,Mat *submat,void *ctx);
851: . row - an array of index sets that contain the global row numbers
852: that comprise each local submatrix
853: . col - an array of index sets that contain the global column numbers
854: that comprise each local submatrix
855: . submat - array of local submatrices
856: - ctx - optional user-defined context for private data for the
857: user-defined func routine (may be null)
859: Notes:
860: PCSetModifySubMatrices() MUST be called before KSPSetUp() and
861: KSPSolve().
863: A routine set by PCSetModifySubMatrices() is currently called within
864: the block Jacobi (PCBJACOBI) and additive Schwarz (PCASM)
865: preconditioners. All other preconditioners ignore this routine.
867: Level: advanced
869: .keywords: PC, set, modify, submatrices
871: .seealso: PCModifySubMatrices()
872: @*/
873: PetscErrorCode PCSetModifySubMatrices(PC pc,PetscErrorCode (*func)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void *ctx)
874: {
877: pc->modifysubmatrices = func;
878: pc->modifysubmatricesP = ctx;
879: return(0);
880: }
884: /*@C
885: PCModifySubMatrices - Calls an optional user-defined routine within
886: certain preconditioners if one has been set with PCSetModifySubMarices().
888: Collective on PC
890: Input Parameters:
891: + pc - the preconditioner context
892: . nsub - the number of local submatrices
893: . row - an array of index sets that contain the global row numbers
894: that comprise each local submatrix
895: . col - an array of index sets that contain the global column numbers
896: that comprise each local submatrix
897: . submat - array of local submatrices
898: - ctx - optional user-defined context for private data for the
899: user-defined routine (may be null)
901: Output Parameter:
902: . submat - array of local submatrices (the entries of which may
903: have been modified)
905: Notes:
906: The user should NOT generally call this routine, as it will
907: automatically be called within certain preconditioners (currently
908: block Jacobi, additive Schwarz) if set.
910: The basic submatrices are extracted from the preconditioner matrix
911: as usual; the user can then alter these (for example, to set different
912: boundary conditions for each submatrix) before they are used for the
913: local solves.
915: Level: developer
917: .keywords: PC, modify, submatrices
919: .seealso: PCSetModifySubMatrices()
920: @*/
921: PetscErrorCode PCModifySubMatrices(PC pc,PetscInt nsub,const IS row[],const IS col[],Mat submat[],void *ctx)
922: {
927: if (!pc->modifysubmatrices) return(0);
928: PetscLogEventBegin(PC_ModifySubMatrices,pc,0,0,0);
929: (*pc->modifysubmatrices)(pc,nsub,row,col,submat,ctx);
930: PetscLogEventEnd(PC_ModifySubMatrices,pc,0,0,0);
931: return(0);
932: }
936: /*@
937: PCSetOperators - Sets the matrix associated with the linear system and
938: a (possibly) different one associated with the preconditioner.
940: Collective on PC and Mat
942: Input Parameters:
943: + pc - the preconditioner context
944: . Amat - the matrix associated with the linear system
945: . Pmat - the matrix to be used in constructing the preconditioner, usually
946: the same as Amat.
947: - flag - flag indicating information about the preconditioner matrix structure
948: during successive linear solves. This flag is ignored the first time a
949: linear system is solved, and thus is irrelevant when solving just one linear
950: system.
952: Notes:
953: The flag can be used to eliminate unnecessary work in the preconditioner
954: during the repeated solution of linear systems of the same size. The
955: available options are
956: + SAME_PRECONDITIONER -
957: Pmat is identical during successive linear solves.
958: This option is intended for folks who are using
959: different Amat and Pmat matrices and wish to reuse the
960: same preconditioner matrix. For example, this option
961: saves work by not recomputing incomplete factorization
962: for ILU/ICC preconditioners.
963: . SAME_NONZERO_PATTERN -
964: Pmat has the same nonzero structure during
965: successive linear solves.
966: - DIFFERENT_NONZERO_PATTERN -
967: Pmat does not have the same nonzero structure.
969: Passing a PETSC_NULL for Amat or Pmat removes the matrix that is currently used.
971: If you wish to replace either Amat or Pmat but leave the other one untouched then
972: first call KSPGetOperators() to get the one you wish to keep, call PetscObjectReference()
973: on it and then pass it back in in your call to KSPSetOperators().
975: Caution:
976: If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion
977: and does not check the structure of the matrix. If you erroneously
978: claim that the structure is the same when it actually is not, the new
979: preconditioner will not function correctly. Thus, use this optimization
980: feature carefully!
982: If in doubt about whether your preconditioner matrix has changed
983: structure or not, use the flag DIFFERENT_NONZERO_PATTERN.
985: More Notes about Repeated Solution of Linear Systems:
986: PETSc does NOT reset the matrix entries of either Amat or Pmat
987: to zero after a linear solve; the user is completely responsible for
988: matrix assembly. See the routine MatZeroEntries() if desiring to
989: zero all elements of a matrix.
991: Level: intermediate
993: .keywords: PC, set, operators, matrix, linear system
995: .seealso: PCGetOperators(), MatZeroEntries()
996: @*/
997: PetscErrorCode PCSetOperators(PC pc,Mat Amat,Mat Pmat,MatStructure flag)
998: {
1000: PetscTruth isbjacobi,isrowbs;
1009: /*
1010: BlockSolve95 cannot use default BJacobi preconditioning
1011: */
1012: if (Amat) {
1013: PetscTypeCompare((PetscObject)Amat,MATMPIROWBS,&isrowbs);
1014: if (isrowbs) {
1015: PetscTypeCompare((PetscObject)pc,PCBJACOBI,&isbjacobi);
1016: if (isbjacobi) {
1017: PCSetType(pc,PCILU);
1018: PetscInfo(pc,"Switching default PC to PCILU since BS95 doesn't support PCBJACOBI\n");
1019: }
1020: }
1021: }
1023: /* reference first in case the matrices are the same */
1024: if (Amat) {PetscObjectReference((PetscObject)Amat);}
1025: if (pc->mat) {MatDestroy(pc->mat);}
1026: if (Pmat) {PetscObjectReference((PetscObject)Pmat);}
1027: if (pc->pmat) {MatDestroy(pc->pmat);}
1028: pc->mat = Amat;
1029: pc->pmat = Pmat;
1031: if (pc->setupcalled == 2 && flag != SAME_PRECONDITIONER) {
1032: pc->setupcalled = 1;
1033: }
1034: pc->flag = flag;
1035: return(0);
1036: }
1040: /*@C
1041: PCGetOperators - Gets the matrix associated with the linear system and
1042: possibly a different one associated with the preconditioner.
1044: Not collective, though parallel Mats are returned if the PC is parallel
1046: Input Parameter:
1047: . pc - the preconditioner context
1049: Output Parameters:
1050: + mat - the matrix associated with the linear system
1051: . pmat - matrix associated with the preconditioner, usually the same
1052: as mat.
1053: - flag - flag indicating information about the preconditioner
1054: matrix structure. See PCSetOperators() for details.
1056: Level: intermediate
1058: Alternative usage: If the operators have NOT been set with KSP/PCSetOperators() then the operators
1059: are created in PC and returned to the user. In this case, if both operators
1060: mat and pmat are requested, two DIFFERENT operators will be returned. If
1061: only one is requested both operators in the PC will be the same (i.e. as
1062: if one had called KSP/PCSetOperators() with the same argument for both Mats).
1063: The user must set the sizes of the returned matrices and their type etc just
1064: as if the user created them with MatCreate(). For example,
1066: $ KSP/PCGetOperators(ksp/pc,&mat,PETSC_NULL,PETSC_NULL); is equivalent to
1067: $ set size, type, etc of mat
1069: $ MatCreate(comm,&mat);
1070: $ KSP/PCSetOperators(ksp/pc,mat,mat,SAME_NONZERO_PATTERN);
1071: $ PetscObjectDereference((PetscObject)mat);
1072: $ set size, type, etc of mat
1074: and
1076: $ KSP/PCGetOperators(ksp/pc,&mat,&pmat,PETSC_NULL); is equivalent to
1077: $ set size, type, etc of mat and pmat
1079: $ MatCreate(comm,&mat);
1080: $ MatCreate(comm,&pmat);
1081: $ KSP/PCSetOperators(ksp/pc,mat,pmat,SAME_NONZERO_PATTERN);
1082: $ PetscObjectDereference((PetscObject)mat);
1083: $ PetscObjectDereference((PetscObject)pmat);
1084: $ set size, type, etc of mat and pmat
1086: The rational for this support is so that when creating a TS, SNES, or KSP the hierarchy
1087: of underlying objects (i.e. SNES, KSP, PC, Mat) and their livespans can be completely
1088: managed by the top most level object (i.e. the TS, SNES, or KSP). Another way to look
1089: at this is when you create a SNES you do not NEED to create a KSP and attach it to
1090: the SNES object (the SNES object manages it for you). Similarly when you create a KSP
1091: you do not need to attach a PC to it (the KSP object manages the PC object for you).
1092: Thus, why should YOU have to create the Mat and attach it to the SNES/KSP/PC, when
1093: it can be created for you?
1094:
1096: .keywords: PC, get, operators, matrix, linear system
1098: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperatorsSet()
1099: @*/
1100: PetscErrorCode PCGetOperators(PC pc,Mat *mat,Mat *pmat,MatStructure *flag)
1101: {
1106: if (mat) {
1107: if (!pc->mat) {
1108: MatCreate(((PetscObject)pc)->comm,&pc->mat);
1109: if (!pc->pmat && !pmat) { /* user did NOT request pmat, so make same as mat */
1110: pc->pmat = pc->mat;
1111: PetscObjectReference((PetscObject)pc->pmat);
1112: }
1113: }
1114: *mat = pc->mat;
1115: }
1116: if (pmat) {
1117: if (!pc->pmat) {
1118: MatCreate(((PetscObject)pc)->comm,&pc->mat);
1119: if (!pc->mat && !mat) { /* user did NOT request mat, so make same as pmat */
1120: pc->mat = pc->pmat;
1121: PetscObjectReference((PetscObject)pc->mat);
1122: }
1123: }
1124: *pmat = pc->pmat;
1125: }
1126: if (flag) *flag = pc->flag;
1127: return(0);
1128: }
1132: /*@C
1133: PCGetOperatorsSet - Determines if the matrix associated with the linear system and
1134: possibly a different one associated with the preconditioner have been set in the PC.
1136: Not collective, though the results on all processes should be the same
1138: Input Parameter:
1139: . pc - the preconditioner context
1141: Output Parameters:
1142: + mat - the matrix associated with the linear system was set
1143: - pmat - matrix associated with the preconditioner was set, usually the same
1145: Level: intermediate
1147: .keywords: PC, get, operators, matrix, linear system
1149: .seealso: PCSetOperators(), KSPGetOperators(), KSPSetOperators(), PCGetOperators()
1150: @*/
1151: PetscErrorCode PCGetOperatorsSet(PC pc,PetscTruth *mat,PetscTruth *pmat)
1152: {
1155: if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE;
1156: if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE;
1157: return(0);
1158: }
1162: /*@
1163: PCFactorGetMatrix - Gets the factored matrix from the
1164: preconditioner context. This routine is valid only for the LU,
1165: incomplete LU, Cholesky, and incomplete Cholesky methods.
1167: Not Collective on PC though Mat is parallel if PC is parallel
1169: Input Parameters:
1170: . pc - the preconditioner context
1172: Output parameters:
1173: . mat - the factored matrix
1175: Level: advanced
1177: Notes: Does not increase the reference count for the matrix so DO NOT destroy it
1179: .keywords: PC, get, factored, matrix
1180: @*/
1181: PetscErrorCode PCFactorGetMatrix(PC pc,Mat *mat)
1182: {
1188: if (pc->ops->getfactoredmatrix) {
1189: (*pc->ops->getfactoredmatrix)(pc,mat);
1190: }
1191: return(0);
1192: }
1196: /*@C
1197: PCSetOptionsPrefix - Sets the prefix used for searching for all
1198: PC options in the database.
1200: Collective on PC
1202: Input Parameters:
1203: + pc - the preconditioner context
1204: - prefix - the prefix string to prepend to all PC option requests
1206: Notes:
1207: A hyphen (-) must NOT be given at the beginning of the prefix name.
1208: The first character of all runtime options is AUTOMATICALLY the
1209: hyphen.
1211: Level: advanced
1213: .keywords: PC, set, options, prefix, database
1215: .seealso: PCAppendOptionsPrefix(), PCGetOptionsPrefix()
1216: @*/
1217: PetscErrorCode PCSetOptionsPrefix(PC pc,const char prefix[])
1218: {
1223: PetscObjectSetOptionsPrefix((PetscObject)pc,prefix);
1224: return(0);
1225: }
1229: /*@C
1230: PCAppendOptionsPrefix - Appends to the prefix used for searching for all
1231: PC options in the database.
1233: Collective on PC
1235: Input Parameters:
1236: + pc - the preconditioner context
1237: - prefix - the prefix string to prepend to all PC option requests
1239: Notes:
1240: A hyphen (-) must NOT be given at the beginning of the prefix name.
1241: The first character of all runtime options is AUTOMATICALLY the
1242: hyphen.
1244: Level: advanced
1246: .keywords: PC, append, options, prefix, database
1248: .seealso: PCSetOptionsPrefix(), PCGetOptionsPrefix()
1249: @*/
1250: PetscErrorCode PCAppendOptionsPrefix(PC pc,const char prefix[])
1251: {
1256: PetscObjectAppendOptionsPrefix((PetscObject)pc,prefix);
1257: return(0);
1258: }
1262: /*@C
1263: PCGetOptionsPrefix - Gets the prefix used for searching for all
1264: PC options in the database.
1266: Not Collective
1268: Input Parameters:
1269: . pc - the preconditioner context
1271: Output Parameters:
1272: . prefix - pointer to the prefix string used, is returned
1274: Notes: On the fortran side, the user should pass in a string 'prifix' of
1275: sufficient length to hold the prefix.
1277: Level: advanced
1279: .keywords: PC, get, options, prefix, database
1281: .seealso: PCSetOptionsPrefix(), PCAppendOptionsPrefix()
1282: @*/
1283: PetscErrorCode PCGetOptionsPrefix(PC pc,const char *prefix[])
1284: {
1290: PetscObjectGetOptionsPrefix((PetscObject)pc,prefix);
1291: return(0);
1292: }
1296: /*@
1297: PCPreSolve - Optional pre-solve phase, intended for any
1298: preconditioner-specific actions that must be performed before
1299: the iterative solve itself.
1301: Collective on PC
1303: Input Parameters:
1304: + pc - the preconditioner context
1305: - ksp - the Krylov subspace context
1307: Level: developer
1309: Sample of Usage:
1310: .vb
1311: PCPreSolve(pc,ksp);
1312: KSPSolve(ksp,b,x);
1313: PCPostSolve(pc,ksp);
1314: .ve
1316: Notes:
1317: The pre-solve phase is distinct from the PCSetUp() phase.
1319: KSPSolve() calls this directly, so is rarely called by the user.
1321: .keywords: PC, pre-solve
1323: .seealso: PCPostSolve()
1324: @*/
1325: PetscErrorCode PCPreSolve(PC pc,KSP ksp)
1326: {
1328: Vec x,rhs;
1329: Mat A,B;
1334: KSPGetSolution(ksp,&x);
1335: KSPGetRhs(ksp,&rhs);
1336: /*
1337: Scale the system and have the matrices use the scaled form
1338: only if the two matrices are actually the same (and hence
1339: have the same scaling
1340: */
1341: PCGetOperators(pc,&A,&B,PETSC_NULL);
1342: if (A == B) {
1343: MatScaleSystem(pc->mat,rhs,x);
1344: MatUseScaledForm(pc->mat,PETSC_TRUE);
1345: }
1347: if (pc->ops->presolve) {
1348: (*pc->ops->presolve)(pc,ksp,rhs,x);
1349: }
1350: return(0);
1351: }
1355: /*@
1356: PCPostSolve - Optional post-solve phase, intended for any
1357: preconditioner-specific actions that must be performed after
1358: the iterative solve itself.
1360: Collective on PC
1362: Input Parameters:
1363: + pc - the preconditioner context
1364: - ksp - the Krylov subspace context
1366: Sample of Usage:
1367: .vb
1368: PCPreSolve(pc,ksp);
1369: KSPSolve(ksp,b,x);
1370: PCPostSolve(pc,ksp);
1371: .ve
1373: Note:
1374: KSPSolve() calls this routine directly, so it is rarely called by the user.
1376: Level: developer
1378: .keywords: PC, post-solve
1380: .seealso: PCPreSolve(), KSPSolve()
1381: @*/
1382: PetscErrorCode PCPostSolve(PC pc,KSP ksp)
1383: {
1385: Vec x,rhs;
1386: Mat A,B;
1391: KSPGetSolution(ksp,&x);
1392: KSPGetRhs(ksp,&rhs);
1393: if (pc->ops->postsolve) {
1394: (*pc->ops->postsolve)(pc,ksp,rhs,x);
1395: }
1396: /*
1397: Scale the system and have the matrices use the scaled form
1398: only if the two matrices are actually the same (and hence
1399: have the same scaling
1400: */
1401: PCGetOperators(pc,&A,&B,PETSC_NULL);
1402: if (A == B) {
1403: MatUnScaleSystem(pc->mat,rhs,x);
1404: MatUseScaledForm(pc->mat,PETSC_FALSE);
1405: }
1406: return(0);
1407: }
1411: /*@C
1412: PCView - Prints the PC data structure.
1414: Collective on PC
1416: Input Parameters:
1417: + PC - the PC context
1418: - viewer - optional visualization context
1420: Note:
1421: The available visualization contexts include
1422: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
1423: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1424: output where only the first processor opens
1425: the file. All other processors send their
1426: data to the first processor to print.
1428: The user can open an alternative visualization contexts with
1429: PetscViewerASCIIOpen() (output to a specified file).
1431: Level: developer
1433: .keywords: PC, view
1435: .seealso: KSPView(), PetscViewerASCIIOpen()
1436: @*/
1437: PetscErrorCode PCView(PC pc,PetscViewer viewer)
1438: {
1439: const PCType cstr;
1440: PetscErrorCode ierr;
1441: PetscTruth mat_exists,iascii,isstring;
1442: PetscViewerFormat format;
1446: if (!viewer) {
1447: PetscViewerASCIIGetStdout(((PetscObject)pc)->comm,&viewer);
1448: }
1452: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
1453: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
1454: if (iascii) {
1455: PetscViewerGetFormat(viewer,&format);
1456: if (((PetscObject)pc)->prefix) {
1457: PetscViewerASCIIPrintf(viewer,"PC Object:(%s)\n",((PetscObject)pc)->prefix);
1458: } else {
1459: PetscViewerASCIIPrintf(viewer,"PC Object:\n");
1460: }
1461: PCGetType(pc,&cstr);
1462: if (cstr) {
1463: PetscViewerASCIIPrintf(viewer," type: %s\n",cstr);
1464: } else {
1465: PetscViewerASCIIPrintf(viewer," type: not yet set\n");
1466: }
1467: if (pc->ops->view) {
1468: PetscViewerASCIIPushTab(viewer);
1469: (*pc->ops->view)(pc,viewer);
1470: PetscViewerASCIIPopTab(viewer);
1471: }
1472: PetscObjectExists((PetscObject)pc->mat,&mat_exists);
1473: if (mat_exists) {
1474: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1475: if (pc->pmat == pc->mat) {
1476: PetscViewerASCIIPrintf(viewer," linear system matrix = precond matrix:\n");
1477: PetscViewerASCIIPushTab(viewer);
1478: MatView(pc->mat,viewer);
1479: PetscViewerASCIIPopTab(viewer);
1480: } else {
1481: PetscObjectExists((PetscObject)pc->pmat,&mat_exists);
1482: if (mat_exists) {
1483: PetscViewerASCIIPrintf(viewer," linear system matrix followed by preconditioner matrix:\n");
1484: } else {
1485: PetscViewerASCIIPrintf(viewer," linear system matrix:\n");
1486: }
1487: PetscViewerASCIIPushTab(viewer);
1488: MatView(pc->mat,viewer);
1489: if (mat_exists) {MatView(pc->pmat,viewer);}
1490: PetscViewerASCIIPopTab(viewer);
1491: }
1492: PetscViewerPopFormat(viewer);
1493: }
1494: } else if (isstring) {
1495: PCGetType(pc,&cstr);
1496: PetscViewerStringSPrintf(viewer," %-7.7s",cstr);
1497: if (pc->ops->view) {(*pc->ops->view)(pc,viewer);}
1498: } else {
1499: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by PC",((PetscObject)viewer)->type_name);
1500: }
1501: return(0);
1502: }
1507: /*@
1508: PCSetInitialGuessNonzero - Tells the iterative solver that the
1509: initial guess is nonzero; otherwise PC assumes the initial guess
1510: is to be zero (and thus zeros it out before solving).
1512: Collective on PC
1514: Input Parameters:
1515: + pc - iterative context obtained from PCCreate()
1516: - flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero
1518: Level: Developer
1520: Notes:
1521: This is a weird function. Since PC's are linear operators on the right hand side they
1522: CANNOT use an initial guess. This function is for the "pass-through" preconditioners
1523: PCKSP, PCREDUNDANT and PCOPENMP and causes the inner KSP object to use the nonzero
1524: initial guess. Not currently working for PCREDUNDANT, that has to be rewritten to use KSP.
1527: .keywords: PC, set, initial guess, nonzero
1529: .seealso: PCGetInitialGuessNonzero(), PCSetInitialGuessKnoll(), PCGetInitialGuessKnoll()
1530: @*/
1531: PetscErrorCode PCSetInitialGuessNonzero(PC pc,PetscTruth flg)
1532: {
1534: pc->nonzero_guess = flg;
1535: return(0);
1536: }
1540: /*@C
1541: PCRegister - See PCRegisterDynamic()
1543: Level: advanced
1544: @*/
1545: PetscErrorCode PCRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PC))
1546: {
1548: char fullname[PETSC_MAX_PATH_LEN];
1552: PetscFListConcat(path,name,fullname);
1553: PetscFListAdd(&PCList,sname,fullname,(void (*)(void))function);
1554: return(0);
1555: }
1559: /*@
1560: PCComputeExplicitOperator - Computes the explicit preconditioned operator.
1562: Collective on PC
1564: Input Parameter:
1565: . pc - the preconditioner object
1567: Output Parameter:
1568: . mat - the explict preconditioned operator
1570: Notes:
1571: This computation is done by applying the operators to columns of the
1572: identity matrix.
1574: Currently, this routine uses a dense matrix format when 1 processor
1575: is used and a sparse format otherwise. This routine is costly in general,
1576: and is recommended for use only with relatively small systems.
1578: Level: advanced
1579:
1580: .keywords: PC, compute, explicit, operator
1582: .seealso: KSPComputeExplicitOperator()
1584: @*/
1585: PetscErrorCode PCComputeExplicitOperator(PC pc,Mat *mat)
1586: {
1587: Vec in,out;
1589: PetscInt i,M,m,*rows,start,end;
1590: PetscMPIInt size;
1591: MPI_Comm comm;
1592: PetscScalar *array,one = 1.0;
1593:
1598: comm = ((PetscObject)pc)->comm;
1599: MPI_Comm_size(comm,&size);
1601: if (!pc->pmat) SETERRQ(PETSC_ERR_ORDER,"You must call KSPSetOperators() or PCSetOperators() before this call");
1602: MatGetVecs(pc->pmat,&in,0);
1603: VecDuplicate(in,&out);
1604: VecGetOwnershipRange(in,&start,&end);
1605: VecGetSize(in,&M);
1606: VecGetLocalSize(in,&m);
1607: PetscMalloc((m+1)*sizeof(PetscInt),&rows);
1608: for (i=0; i<m; i++) {rows[i] = start + i;}
1610: MatCreate(comm,mat);
1611: MatSetSizes(*mat,m,m,M,M);
1612: if (size == 1) {
1613: MatSetType(*mat,MATSEQDENSE);
1614: MatSeqDenseSetPreallocation(*mat,PETSC_NULL);
1615: } else {
1616: MatSetType(*mat,MATMPIAIJ);
1617: MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);
1618: }
1620: for (i=0; i<M; i++) {
1622: VecSet(in,0.0);
1623: VecSetValues(in,1,&i,&one,INSERT_VALUES);
1624: VecAssemblyBegin(in);
1625: VecAssemblyEnd(in);
1627: /* should fix, allowing user to choose side */
1628: PCApply(pc,in,out);
1629:
1630: VecGetArray(out,&array);
1631: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
1632: VecRestoreArray(out,&array);
1634: }
1635: PetscFree(rows);
1636: VecDestroy(out);
1637: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
1638: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
1639: return(0);
1640: }