Actual source code: schur.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: schur.c,v 1.15 2000/10/04 02:16:21 knepley Exp $";
3: #endif
4: /*
5: Defines the Schur preconditioner due to Ahmed Sameh for AIJ matrices
6: */
7: #include src/sles/pc/pcimpl.h
8: #include schur.h
10: extern int GridResetConstrainedMultiply_Private(Grid, GMat);
12: /*
13: Math:
15: There are several possible options for this type of preconditioner. I call it
16: a Schur complement preconditioner because we split the problem into blocks
18: / A B
19: B^T 0 /
21: and generally must provide a mechanism for solving (preconditioning) the
22: Schur complement B^T A^{-1} B at each itertion of the overall preconditioner.
23: We could include Uzawa type preconditioners, but since they are so simple, I
24: have already abstracted them using GMatCreateUzawa() and related functions.
25: Here I am concerned with preconditioners for the entire saddle point problem
26: above, the simplest being a block diagonal (Jacobi)
28: / hat A 0
29: 0 -hat S /
31: or the slightly more sophisticated block triangular
33: / hat A B
34: 0 -hat S /
36: culminating in the full block LU factorization
38: / hat A 0 / I hat A^{-1} B
39: B^T I / 0 -hat S /
41: These all depend on solving the Schur complement S, so we need an effective
42: preconditioner for it. Elman, Wathen and Silvester propose the "BFBT"
43: preconditioner which essentially projects the problem onto the pressure space
44: and then applies the convection diffusion operator. So we have
46: (B^T A^{-1} B)^{-1} = (B^T B)^{-1} B^T A B (B^T B)^{-1}
47: = A_p (B^T B)^{-1}
49: where in the last step they assume a commutation relation of the form
51: A B = B A_p
53: so that the momentum equations on the pressure mesh are equivalent in some
54: sense to those on the velocity mesh. In addition, they recommend several cycles
55: of multigrid for hat A^{-1}.
56: Sameh instead proposes to precondition S with the pressure lapalcian L_p and
57: perform fast solves with hat A using the balanced method. I do not understand
58: why Wathen, et.al. do not do the full block LU instead of the triangular thing,
59: but we can test this. We can afford to do an iterative method for S since
60: applications of hat A^{-1} are cheap, but this may be compared with the sparse
61: factorizations of Saad.
62: Lastly, it must be noted that you can solve a preconditioned system for hat A
63: or hat S, or merely apply the preconditioners. The cost of solving must be
64: weighed against the extra outer iterations incured by the weaker approximation.
66: Programming:
68: We allow several modes for this preconditioner since parts of the matrix may
69: not be available, e.g. if the matrix is given as implicitly P^T A P. The flag
70: explicitConstraints signals that a constrained matrix has been explicitly formed.
71: Otherwise, we use GridGetConstraintMatrix() to retrieve the projector P onto the
72: constraint space.
73: Also, we may have new variables introduced by constraints. These will be present
74: already in explicit matrices, but are not yet included in the implicit formulation.
75: Thus, we have functions which add this contribution to a given vector, with functions
76: like GVecEvaluateJacobianConstrained(). This function operates in the projected space
77: so that we get an action like
79: P^T A P + G
81: where G incorporates the contribution of the new fields.
82: */
84: /*--------------------------------------------- Public Functions ----------------------------------------------------*/
85: /*@C PCSchurSetGradientOperator
86: This function sets the operator and its transpose, which define the
87: constraint manifold in the saddle-point problem. For instance, in
88: the Stokes problem, this operator is the gradient, and its transpose
89: is the divergence. It must be called prior to PCSetUp().
91: Collective on PC
93: Input Parameters:
94: + pc - The preconditioner context
95: . gradOp - The gradient operator
96: - divOp - The divergence operator
98: Level: intermediate
100: .keywords schur, set, gradient
101: .seealso PCSchurGetIterationNumber()
102: @*/
103: int PCSchurSetGradientOperator(PC pc, int gradOp, int divOp)
104: {
105: PC_Schur *s;
109: s = (PC_Schur *) pc->data;
110: s->gradOp = gradOp;
111: s->divOp = divOp;
112: return(0);
113: }
114: /*@C
115: PCSchurGetIterationNumber - Get the iteration numbers since the solve was started.
117: Not collective
119: Input Parameter:
120: . pc - The preconditioner context
122: Output
123: + its - The total number of iterations used to solve A
124: - schurIts - The number of iterations used to solve the Schur complement
126: Level: intermediate
128: .keywords schur, iteration
129: .seealso PCSchurSetGradientOperator()
130: @*/
131: int PCSchurGetIterationNumber(PC pc, int *its, int *schurIts)
132: {
133: PC_Schur *s;
139: s = (PC_Schur *) pc->data;
140: *its = s->iter;
141: *schurIts = s->schurIter;
142: return(0);
143: }
145: /*@C
146: PCSchurInnerMonitor - Print the residual norm of the solver for the inner system A.
148: Collective on KSP
150: Input Parameters:
151: . ksp - The KSP context
152: . n - The iteration number
153: . rnorm - The 2-norm of residual
154: . ctx - The Schur context
156: Level: intermediate
158: .keywords: KSP, schur, monitor, norm
159: .seealso: KSPSetMonitor(), PCSchurMonitor()
160: @*/
161: int PCSchurInnerMonitor(KSP ksp, int n, PetscReal rnorm, void *ctx)
162: {
163: MPI_Comm comm;
164: int ierr;
168: PetscObjectGetComm((PetscObject) ksp, &comm);
169: PetscPrintf(comm, " inner iter = %d, residual norm %g n", n, rnorm);
170: return(0);
171: }
173: /*@C
174: PCSchurMonitor - Print the residual norm of the solver for the Schur complement system S.
176: Collective on KSP
178: Input Parameters:
179: . ksp - The KSP context
180: . n - The iteration number
181: . rnorm - The 2-norm of residual
182: . ctx - The Schur context
184: Level: intermediate
186: .keywords: KSP, schur, monitor, norm
187: .seealso: KSPSetMonitor(), PCSchurInnerMonitor()
188: @*/
189: int PCSchurMonitor(KSP ksp, int n, PetscReal rnorm, void *ctx)
190: {
191: MPI_Comm comm;
192: int ierr;
196: PetscObjectGetComm((PetscObject) ksp, &comm);
197: PetscPrintf(comm, " schur iter = %d, residual norm %g n", n, rnorm);
198: return(0);
199: }
202: /*@C
203: PCSchurSolveMonitor - Print the number of iterates for each inner solve in the Schur complement system.
205: Collective on KSP
207: Input Parameters:
208: . ksp - The KSP context
209: . n - The iteration number
210: . rnorm - The 2-norm of residual
211: . ctx - The Schur context
213: Level: intermediate
215: .keywords: KSP, schur, monitor, norm
216: .seealso: KSPSetMonitor(), PCSchurMonitor(), PCSchurInnerMonitor()
217: @*/
218: int PCSchurSolveMonitor(KSP ksp, int n, PetscReal rnorm, void *ctx)
219: {
220: PC_Schur *s = (PC_Schur *) ctx;
221: KSP innerKSP;
222: int its;
223: int ierr;
227: SLESGetKSP(s->sles, &innerKSP);
228: KSPGetIterationNumber(innerKSP, &its);
229: if (n > 0)
230: PetscLogInfo(ksp, "PCSchur: %d iterations in A solve %dn", its, n);
231: return(0);
232: }
234: /*--------------------------------------------- Private Functions ---------------------------------------------------*/
235: static int PCDestroyStructures_Schur(PC pc)
236: {
237: PC_Schur *s = (PC_Schur *) pc->data;
238: int ierr;
241: /* Destroy operator structures */
242: PetscFree(s->momOps);
243: PetscFree(s->momOpIsALE);
244: PetscFree(s->momOpAlphas);
246: /* Destroy variable orderings */
247: VarOrderingDestroy(s->sOrder);
248: VarOrderingDestroy(s->tOrder);
249: LocalVarOrderingDestroy(s->sLocOrder);
250: LocalVarOrderingDestroy(s->tLocOrder);
252: /* Destroy matrices and reorderings */
253: GMatDestroyUzawa(s->S);
254: GMatDestroy(s->A);
255: if (s->sparseA != PETSC_NULL) {
256: GMatDestroy(s->sparseA);
257: }
258: GMatDestroy(s->B);
259: if (s->lap != PETSC_NULL) {
260: GMatDestroy(s->lap);
261: }
262: if (s->rowPerm != PETSC_NULL) {
263: ISDestroy(s->rowPerm);
264: ISDestroy(s->colPerm);
265: }
266: /* Destroy field structures */
267: VecDestroy(s->u);
268: VecDestroy(s->uNew);
269: VecDestroy(s->uNew2);
270: VecDestroy(s->p);
271: VecDestroy(s->pNew);
272: VecScatterDestroy(s->uScatter);
273: VecScatterDestroy(s->pScatter);
274: if ((s->isConstrained == PETSC_TRUE) && (s->explicitConstraints == PETSC_FALSE)) {
275: VecDestroy(s->projX);
276: VecDestroy(s->projY);
277: VecDestroy(s->projZ);
278: }
280: return(0);
281: }
283: int PCValidQ_Schur(PC pc)
284: {
285: PC_Schur *s = (PC_Schur *) pc->data;
288: if (pc->setupcalled < 2)
289: PetscFunctionReturn(1);
292: /* Check dimensions */
294: /* Check fields */
300: if ((s->isConstrained == PETSC_TRUE) && (s->explicitConstraints == PETSC_FALSE)) {
304: }
308: /* Check inner solvers */
312: /* Check matrices */
316: if (s->sparseA != PETSC_NULL) {
318: }
320: return(0);
321: }
323: static int PCDebug_Schur(PC pc)
324: {
325: PC_Schur *s = (PC_Schur *) pc->data;
326: Vec x, y;
327: Mat P;
328: Vec xx;
329: Vec yy;
330: PetscScalar *xArray;
331: PetscScalar zero = 0.0, one = 1.0, minusOne = -1.0;
332: PetscReal uNorm, pNorm;
333: int row, size;
334: int ierr;
337: VecDuplicate(pc->vec, &x);
338: VecDuplicate(pc->vec, &y);
339: if (s->explicitConstraints == PETSC_FALSE) {
340: GridGetConstraintMatrix(s->grid, &P);
341: xx = s->projX;
342: yy = s->projY;
343: } else {
344: xx = x;
345: yy = y;
346: }
348: /* Check the matrix-vector product */
349: VecGetSize(x, &size);
350: VecGetArray(x, &xArray);
351: for(row = 0; row < size; row++) {
352: /* Start with e_row */
353: VecSet(&zero, x);
354: xArray[row] = 1.0;
355: /* hat A x */
356: MatMult(pc->mat, x, y);
357: if (s->explicitConstraints == PETSC_FALSE) {
358: /* Project vectors */
359: MatMult(P, x, xx);
360: MatMult(P, y, yy);
361: }
362: /* New hat A x */
363: VecScatterBegin(xx, s->u, INSERT_VALUES, SCATTER_FORWARD, s->uScatter);
364: VecScatterEnd(xx, s->u, INSERT_VALUES, SCATTER_FORWARD, s->uScatter);
365: VecScatterBegin(xx, s->p, INSERT_VALUES, SCATTER_FORWARD, s->pScatter);
366: VecScatterEnd(xx, s->p, INSERT_VALUES, SCATTER_FORWARD, s->pScatter);
367: if (s->colPerm != PETSC_NULL) {
368: VecPermute(s->u, s->colPerm, PETSC_FALSE);
369: }
370: /* uNew2 = hat A u + B p, pNew = B u */
371: MatMult(s->A, s->u, s->uNew);
372: MatMult(s->B, s->p, s->uNew2);
373: VecAXPY(&one, s->uNew2, s->uNew);
374: MatMultTranspose(s->B, s->u, s->pNew);
375: if (s->rowPerm != PETSC_NULL) {
376: VecPermute(s->uNew, s->rowPerm, PETSC_TRUE);
377: }
378: if (s->explicitConstraints == PETSC_FALSE) {
379: /* Add in contribution from new variables */
380: /* GVecEvaluateJacobianConstrained(xx, xx); */
381: VecScatterBegin(xx, s->uNew, ADD_VALUES, SCATTER_FORWARD, s->uScatter);
382: VecScatterEnd(xx, s->uNew, ADD_VALUES, SCATTER_FORWARD, s->uScatter);
383: VecScatterBegin(xx, s->pNew, ADD_VALUES, SCATTER_FORWARD, s->pScatter);
384: VecScatterEnd(xx, s->pNew, ADD_VALUES, SCATTER_FORWARD, s->pScatter);
385: }
386: /* Compare field vectors */
387: VecScatterBegin(yy, s->u, INSERT_VALUES, SCATTER_FORWARD, s->uScatter);
388: VecScatterEnd(yy, s->u, INSERT_VALUES, SCATTER_FORWARD, s->uScatter);
389: VecScatterBegin(yy, s->p, INSERT_VALUES, SCATTER_FORWARD, s->pScatter);
390: VecScatterEnd(yy, s->p, INSERT_VALUES, SCATTER_FORWARD, s->pScatter);
391: VecAXPY(&minusOne, s->u, s->uNew);
392: VecAXPY(&minusOne, s->p, s->pNew);
393: VecNorm(s->uNew, NORM_2, &uNorm);
394: VecNorm(s->pNew, NORM_2, &pNorm);
395: if ((uNorm > 1.0e-8) || (pNorm > 1.0e-8)) {
396: SETERRQ(PETSC_ERR_PLIB, "Invalid field matrix");
397: }
398: }
400: VecRestoreArray(x, &xArray);
401: VecDestroy(x);
402: VecDestroy(y);
403: return(0);
404: }
406: int PCSetFromOptions_Inner(PC pc, SLES sles, SLES schurSles)
407: {
408: PC_Schur *s = (PC_Schur *) pc->data;
409: KSP ksp;
410: char *prefix;
411: PetscTruth opt;
412: int ierr;
415: SLESGetOptionsPrefix(sles, &prefix);
416: SLESGetKSP(sles, &ksp);
417: PetscOptionsHasName(prefix, "-ksp_monitor", &opt);
418: if (opt == PETSC_TRUE) {
419: KSPClearMonitor(ksp);
420: KSPSetMonitor(ksp, PCSchurInnerMonitor, s, PETSC_NULL);
421: }
422: SLESGetOptionsPrefix(schurSles, &prefix);
423: SLESGetKSP(schurSles, &ksp);
424: PetscOptionsHasName(prefix, "-ksp_monitor", &opt);
425: if (opt == PETSC_TRUE) {
426: KSPClearMonitor(ksp);
427: KSPSetMonitor(ksp, PCSchurMonitor, s, PETSC_NULL);
428: }
429: PetscOptionsHasName(prefix, "-ksp_solve_monitor", &opt);
430: if (opt == PETSC_TRUE) {
431: KSPClearMonitor(ksp);
432: KSPSetMonitor(ksp, PCSchurSolveMonitor, s, PETSC_NULL);
433: }
434: return(0);
435: }
437: static int PCSetFromOptions_Schur(PC pc)
438: {
439: PC_Schur *s = (PC_Schur *) pc->data;
440: PetscTruth opt;
441: int ierr;
444: /* Enable explicit constraints */
445: PetscOptionsHasName(pc->prefix, "-pc_schur_explicit", &opt);
446: if (opt == PETSC_TRUE)
447: s->explicitConstraints = PETSC_TRUE;
449: /* Enable laplacian preconditioning of the Schur complement */
450: PetscOptionsHasName(pc->prefix, "-pc_schur_lap", &opt);
451: if (opt == PETSC_TRUE)
452: s->useLaplacian = PETSC_TRUE;
454: /* Enable Mathematica support */
455: PetscOptionsHasName(pc->prefix, "-pc_schur_math", &opt);
456: if (opt == PETSC_TRUE) {
457: s->useMath = PETSC_TRUE;
458: PetscViewerMathematicaOpen(pc->comm, PETSC_DECIDE, PETSC_NULL, PETSC_NULL, &s->mathViewer);
459: }
461: /* Setup inner solvers */
462: PCSetFromOptions_Inner(pc, s->sles, s->schurSles);
463: return(0);
464: }
466: static int PCSchurGetMomentumOperators(PC pc)
467: {
468: PC_Schur *s = (PC_Schur *) pc->data;
469: int numOps; /* Grid op query: The number of grid matrix operators used for A */
470: int sField, tField; /* Grid op query: The shape and test function fields */
471: PetscScalar alpha; /* Grid op query: The scalar multiplier */
472: PetscTruth isALE; /* Grid op query: The flag for ALE operators */
473: int op;
474: int ierr;
477: GridGetNumMatOperators(s->grid, &s->numMomOps);
478: s->numMomOps -= 2;
479: PetscMalloc(s->numMomOps * sizeof(int), &s->momOps);
480: PetscMalloc(s->numMomOps * sizeof(PetscTruth), &s->momOpIsALE);
481: PetscMalloc(s->numMomOps * sizeof(double), &s->momOpAlphas);
482: /* Find gradient fields */
483: ierr = GridGetMatOperatorStart(s->grid, &op, &sField, &tField, &alpha, &isALE);
484: numOps = 0;
485: while(op >= 0) {
486: if (op == s->gradOp) {
487: s->sField = sField;
488: s->tField = tField;
489: s->gradOpAlpha = alpha;
490: } else if (op != s->divOp) {
491: s->momOps[numOps] = op;
492: s->momOpIsALE[numOps] = isALE;
493: s->momOpAlphas[numOps] = alpha;
494: numOps++;
495: }
496: GridGetMatOperatorNext(s->grid, &op, &sField, &tField, &alpha, &isALE);
497: }
498: if ((s->sField < 0) || (s->tField < 0)) {
499: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Gradient operator not found in grid.");
500: }
501: if (numOps != s->numMomOps) {
502: SETERRQ(PETSC_ERR_PLIB, "Inconsistent operator setup in grid");
503: }
505: return(0);
506: }
508: static int PCSchurReorderMatrices(PC pc)
509: {
510: PC_Schur *s = (PC_Schur *) pc->data;
511: Mat newA, newB; /* The reordering matrices */
512: char orderType[256]; /* The type of reordering */
513: IS tempIS; /* The identity ordering for B */
514: int bw; /* The bandwidth limitation on the reordered matrix */
515: PetscReal frac; /* The fractional bandwidth limitation on the reordered matrix */
516: double tol; /* The drop tolerance */
517: int m, n;
518: PetscTruth opt, issparse;
519: int ierr;
522: /* Reduce the bandwidth of A */
523: PetscOptionsHasName(pc->prefix, "-pc_schur_reorder", &opt);
524: if (opt == PETSC_TRUE) {
525: MatGetOrdering(s->A, MATORDERING_RCM, &s->rowPerm, &s->colPerm);
526: PetscOptionsGetString(pc->prefix, "-pc_schur_reorder", orderType, 255, &opt);
527: if (opt == PETSC_TRUE) {
528: PetscStrcasecmp(orderType, "sparse", &issparse);
529: if (issparse) {
530: bw = PETSC_DECIDE;
531: PetscOptionsGetInt(pc->prefix, "-pc_schur_sparsify_bw", &bw, &opt);
532: frac = 0.0;
533: PetscOptionsGetReal(pc->prefix, "-pc_schur_sparsify_frac", &frac, &opt);
534: tol = 0.0;
535: PetscOptionsGetReal(pc->prefix, "-pc_schur_sparsify_tol", &frac, &opt);
536: MatPermuteSparsify(s->A, bw, frac, tol, s->rowPerm, s->colPerm, &s->sparseA);
537: }
538: MatPermute(s->A, s->rowPerm, s->colPerm, &newA);
539: }
540: MatDestroy(s->A);
541: MatGetSize(s->B, &m, &n);
542: ISCreateStride(pc->comm, n, 0, 1, &tempIS);
543: ISSetPermutation(tempIS);
544: MatPermute(s->B, s->rowPerm, tempIS, &newB);
545: ISDestroy(tempIS);
546: MatDestroy(s->B);
547: s->A = newA;
548: s->B = newB;
549: PetscObjectCompose((PetscObject) s->A, "Grid", (PetscObject) s->grid);
550: PetscObjectCompose((PetscObject) s->B, "Grid", (PetscObject) s->grid);
551: }
553: /* View the new matrices */
554: PetscOptionsHasName(PETSC_NULL, "-pc_schur_view", &opt);
555: if (opt == PETSC_TRUE) {
556: PetscViewer viewer = PETSC_VIEWER_DRAW_(pc->comm);
557: PetscDraw draw;
559: PetscViewerDrawGetDraw(viewer, 0, &draw);
560: PetscDrawSetPause(draw, -1);
561: MatView(s->A, viewer);
562: if (s->sparseA != PETSC_NULL) {
563: MatView(s->sparseA, viewer);
564: }
565: MatView(s->B, viewer);
566: }
568: return(0);
569: }
571: static int PCSchurCreateMatrices(PC pc)
572: {
573: PC_Schur *s = (PC_Schur *) pc->data;
574: Grid grid = s->grid;
575: VarOrdering tempOrder;
576: int op;
577: int ierr;
580: GridIsConstrained(grid, &s->isConstrained);
581: /* Create test function variable orderings */
582: VarOrderingCreateGeneral(grid, 1, &s->tField, &s->tOrder);
583: LocalVarOrderingCreate(grid, 1, &s->tField, &s->tLocOrder);
585: /* Create matrices */
586: if ((s->isConstrained == PETSC_FALSE) || (s->explicitConstraints == PETSC_FALSE)) {
587: /* Create shape function variable orderings */
588: VarOrderingCreateGeneral(grid, 1, &s->sField, &s->sOrder);
589: LocalVarOrderingCreate(grid, 1, &s->sField, &s->sLocOrder);
590: GMatCreateRectangular(grid, s->tOrder, s->tOrder, &s->A);
591: GMatCreateRectangular(grid, s->sOrder, s->tOrder, &s->B);
592: if (s->useLaplacian == PETSC_TRUE) {
593: GMatCreateRectangular(grid, s->sOrder, s->sOrder, &s->lap);
594: }
595: } else {
596: /* Constrain test function variable ordering */
597: VarOrderingConstrain(grid, s->tOrder, &tempOrder);
598: VarOrderingDestroy(s->tOrder);
599: s->tOrder = tempOrder;
600: /* Create shape function variable orderings */
601: VarOrderingCreateGeneral(grid, 1, &s->sField, &s->sOrder);
602: LocalVarOrderingCreate(grid, 1, &s->sField, &s->sLocOrder);
603: GMatCreateRectangular(grid, s->tOrder, s->tOrder, &s->A);
604: GMatCreateRectangular(grid, s->sOrder, s->tOrder, &s->B);
605: if (s->useLaplacian == PETSC_TRUE) {
606: GMatCreateRectangular(grid, s->sOrder, s->sOrder, &s->lap);
607: }
608: }
609: MatSetOption(s->A, MAT_NEW_NONZERO_ALLOCATION_ERR);
610: MatSetOption(s->B, MAT_NEW_NONZERO_ALLOCATION_ERR);
611: if (s->useLaplacian == PETSC_TRUE) {
612: MatSetOption(s->lap, MAT_NEW_NONZERO_ALLOCATION_ERR);
613: }
615: /* Construct the operators */
616: if ((s->isConstrained == PETSC_FALSE) || (s->explicitConstraints == PETSC_FALSE)) {
617: for(op = 0; op < s->numMomOps; op++) {
618: if (s->momOpIsALE[op] == PETSC_FALSE) {
619: GMatEvaluateOperatorGalerkin(s->A, PETSC_NULL, 1, &s->tField, &s->tField, s->momOps[op], s->momOpAlphas[op],
620: MAT_FINAL_ASSEMBLY, s);
621:
622: } else {
623: GMatEvaluateALEOperatorGalerkin(s->A, 1, &s->tField, s->tOrder, s->tLocOrder, &s->tField, s->tOrder,
624: s->tLocOrder, s->momOps[op], s->momOpAlphas[op],
625: MAT_FINAL_ASSEMBLY, s);
626:
627: }
628: }
629: GMatEvaluateOperatorGalerkin(s->B, PETSC_NULL, 1, &s->sField, &s->tField, s->gradOp, s->gradOpAlpha, MAT_FINAL_ASSEMBLY, s);
630:
631: if (s->useLaplacian == PETSC_TRUE) {
632: GMatEvaluateOperatorGalerkin(s->lap, PETSC_NULL, 1, &s->sField, &s->sField, LAPLACIAN, 1.0, MAT_FINAL_ASSEMBLY, s);
633:
634: }
635: } else {
636: for(op = 0; op < s->numMomOps; op++) {
637: if (s->momOpIsALE[op] == PETSC_FALSE) {
638: GMatEvaluateOperatorGalerkin(s->A, PETSC_NULL, 1, &s->tField, &s->tField, s->momOps[op], s->momOpAlphas[op],
639: MAT_FINAL_ASSEMBLY, s);
640:
641: } else {
642: GMatEvaluateALEConstrainedOperatorGalerkin(s->A, 1, &s->tField, s->tOrder, s->tLocOrder, &s->tField, s->tOrder,
643: s->tLocOrder, s->momOps[op], s->momOpAlphas[op], MAT_FINAL_ASSEMBLY, s);
644:
645: }
646: }
647: GMatEvaluateNewFields(s->A, 1, &s->tField, s->tOrder, s->tLocOrder, &s->tField, s->tOrder, s->tLocOrder, 1.0,
648: MAT_FINAL_ASSEMBLY);
649:
650: GMatEvaluateOperatorGalerkin(s->B, PETSC_NULL, 1, &s->sField, &s->tField, s->gradOp, s->gradOpAlpha, MAT_FINAL_ASSEMBLY, s);
651:
652: if (s->useLaplacian == PETSC_TRUE) {
653: GMatEvaluateOperatorGalerkin(s->lap, PETSC_NULL, 1, &s->sField, &s->sField, LAPLACIAN, 1.0, MAT_FINAL_ASSEMBLY, s);
654:
655: }
656: /* Reset matrix multiplication */
657: GridResetConstrainedMultiply_Private(grid, s->A);
658: GridResetConstrainedMultiply_Private(grid, s->B);
659: if (s->useLaplacian == PETSC_TRUE) {
660: GridResetConstrainedMultiply_Private(grid, s->lap);
661: }
662: }
664: /* Reduce the bandwidth of A */
665: PCSchurReorderMatrices(pc);
667: return(0);
668: }
670: static int PCSchurCreateInnerSLES(PC pc)
671: {
672: PC_Schur *s = (PC_Schur *) pc->data;
673: int ierr;
676: /* Create the momentum solver */
677: SLESCreate(pc->comm, &s->sles);
678: SLESAppendOptionsPrefix(s->sles, "pc_schur_inner_");
680: /* Create the Schur complement solver */
681: SLESCreate(pc->comm, &s->schurSles);
682: SLESAppendOptionsPrefix(s->schurSles, "pc_schur_");
683: return(0);
684: }
686: static int PCSchurSetupInnerSLES(PC pc)
687: {
688: PC_Schur *s = (PC_Schur *) pc->data;
689: KSP innerKsp;
690: PC innerPc;
691: int ierr;
694: /* Create the momentum solver */
695: SLESCreate(pc->comm, &s->sles);
696: if (s->sparseA != PETSC_NULL) {
697: SLESSetOperators(s->sles, s->A, s->sparseA, DIFFERENT_NONZERO_PATTERN);
698: } else {
699: SLESSetOperators(s->sles, s->A, s->A, SAME_NONZERO_PATTERN);
700: }
701: SLESAppendOptionsPrefix(s->sles, "pc_schur_inner_");
702: SLESGetKSP(s->sles, &innerKsp);
703: SLESGetPC(s->sles, &innerPc);
704: KSPSetType(innerKsp, KSPPREONLY);
705: PCSetType(innerPc, PCLU);
706: SLESSetFromOptions(s->sles);
708: /* Create the Schur complement solver */
709: SLESCreate(pc->comm, &s->schurSles);
710: GMatCreateUzawa(s->A, s->B, s->sles, &s->S);
711: SLESSetOperators(s->schurSles, s->S, s->S, SAME_NONZERO_PATTERN);
712: SLESAppendOptionsPrefix(s->schurSles, "pc_schur_");
713: SLESGetKSP(s->schurSles, &innerKsp);
714: SLESGetPC(s->schurSles, &innerPc);
715: KSPSetType(innerKsp, KSPGMRES);
716: if (s->useLaplacian == PETSC_TRUE) {
717: PCSetOperators(innerPc, s->S, s->lap, DIFFERENT_NONZERO_PATTERN);
718: PCSetType(innerPc, PCLU);
719: } else {
720: PCSetType(innerPc, PCNONE);
721: }
722: SLESSetFromOptions(s->schurSles);
724: /* Setup monitors */
725: PCSetFromOptions_Inner(pc, s->sles, s->schurSles);
726: return(0);
727: }
729: static int PCCreateStructures_Schur(PC pc)
730: {
731: PC_Schur *s = (PC_Schur *) pc->data;
732: Grid grid = s->grid;
733: int ierr;
736: /* Create matrices */
737: PCSchurCreateMatrices(pc);
739: /* Create the work vectors and field scatters */
740: GVecCreateRectangular(grid, s->tOrder, &s->u);
741: GVecDuplicate(s->u, &s->uNew);
742: GVecDuplicate(s->u, &s->uNew2);
743: GVecCreateRectangular(grid, s->sOrder, &s->p);
744: GVecDuplicate(s->p, &s->pNew);
745: if ((s->isConstrained == PETSC_TRUE) && (s->explicitConstraints == PETSC_FALSE)) {
746: GVecCreate(grid, &s->projX);
747: GVecDuplicate(s->projX, &s->projY);
748: GVecCreateConstrained(grid, &s->projZ);
749: GVecScatterCreate(s->projX, s->u, &s->uScatter);
750: GVecScatterCreate(s->projX, s->p, &s->pScatter);
751: } else {
752: GVecScatterCreate(pc->vec, s->u, &s->uScatter);
753: GVecScatterCreate(pc->vec, s->p, &s->pScatter);
754: }
756: /* Setup Inner Solvers */
757: PCSchurSetupInnerSLES(pc);
759: return(0);
760: }
762: static int PCSetUp_Schur(PC pc)
763: {
764: PC_Schur *s = (PC_Schur *) pc->data;
765: PetscTruth opt;
766: int ierr;
769: /* This lets the factorization become stale, that is uses the same preconditioner for several Newton steps */
770: if (pc->setupcalled > 0) {
771: return(0);
772: }
773: pc->setupcalled = 2;
775: /* Get the grid */
776: GVecGetGrid(pc->vec, &s->grid);
778: /* Divide up the problem */
779: PCSchurGetMomentumOperators(pc);
781: /* Create matrices, vector storage, scatters, and inner solvers */
782: PCCreateStructures_Schur(pc);
784: #ifdef PETSC_USE_BOPT_g
785: #endif
786: PetscOptionsHasName(PETSC_NULL, "-pc_schur_debug", &opt);
787: if (opt == PETSC_TRUE) {
788: PCDebug_Schur(pc);
789: }
791: return(0);
792: }
794: static int PCPreSolve_Schur(PC pc, KSP ksp, Vec x, Vec rhs)
795: {
796: PC_Schur *s = (PC_Schur *) pc->data;
799: s->iter = 0;
800: s->schurIter = 0;
801: return(0);
802: }
804: static int PCPostSolve_Schur(PC pc, KSP ksp, Vec x, Vec rhs)
805: {
807: return(0);
808: }
810: static int PCApply_Schur(PC pc, Vec x, Vec y)
811: {
812: PC_Schur *s = (PC_Schur *) pc->data;
813: Vec xx = x;
814: Vec yy = y;
815: Vec z = s->projZ;
816: Vec tempSol;
817: Mat P, invP;
818: PetscScalar minusOne = -1.0;
819: int its, schurIts;
820: PetscTruth opt;
821: int ierr;
824: PetscOptionsHasName(PETSC_NULL, "-pc_schur_test", &opt);
825: if (opt == PETSC_TRUE) {
826: /* Set the input vector to Ax */
827: VecDuplicate(x, &tempSol);
828: VecCopy(x, tempSol);
829: MatMult(pc->mat, x, y);
830: VecCopy(y, x);
831: }
833: /* Project: Apply (P^+)^T = P (P^T P)^{-1} */
834: if ((s->isConstrained == PETSC_TRUE) && (s->explicitConstraints == PETSC_FALSE)) {
835: xx = s->projX;
836: yy = s->projY;
837: GridGetConstraintMatrix(s->grid, &P);
838: GridGetConstraintInverse(s->grid, &invP);
839: MatMult(invP, x, z);
840: MatMult(P, z, xx);
841: }
843: /* / A^{-1} 0 0 / u / A^{-1} u / w_u
844: L^{-1} x = | -B^T A^{-1} I 0 | | p | = | p - B^T A^{-1} u | = | w_p |
845: 0 0 G^{-1} / U / G^{-1} U / w_U /
846: */
847: VecScatterBegin(xx, s->u, INSERT_VALUES, SCATTER_FORWARD, s->uScatter);
848: VecScatterEnd(xx, s->u, INSERT_VALUES, SCATTER_FORWARD, s->uScatter);
849: /* Permute u */
850: if (s->colPerm != PETSC_NULL) {
851: VecPermute(s->u, s->colPerm, PETSC_FALSE);
852: }
853: /* Velocity solve: Now w_u = s->uNew */
854: SLESSolve(s->sles, s->u, s->uNew, &its);
855: s->iter += its;
856: PetscLogInfo(pc, "PCSchur: %d iterations in A solve %dn", its, 0);
857: /* Pressure solve: Now w_p = s->p */
858: MatMultTranspose(s->B, s->uNew, s->p);
859: VecScale(&minusOne, s->p);
860: VecScatterBegin(xx, s->p, ADD_VALUES, SCATTER_FORWARD, s->pScatter);
861: VecScatterEnd(xx, s->p, ADD_VALUES, SCATTER_FORWARD, s->pScatter);
863: /* / I A^{-1} B (B^T A^{-1} B)^{-1} 0 / w_u / w_u + A^{-1} B (B^T A^{-1} B)^{-1} w_p
864: U^{-1} w = | 0 -(B^T A^{-1} B)^{-1} 0 | | w_p | = | -(B^T A^{-1} B)^{-1} w_p |
865: 0 0 I / w_U / w_U /
866: */
867: /* Schur solve: Now s->pNew = -S^{-1} w_p */
868: SLESSolve(s->schurSles, s->p, s->pNew, &schurIts);
869: VecScale(&minusOne, s->pNew);
870: s->schurIter += schurIts;
871: PetscLogInfo(pc, "PCSchur: %d iterations in B^T A^{-1} B solven", schurIts);
872: /* Velocity solve: Now s->uNew2 = -A^{-1} B S^{-1} w_p */
873: MatMult(s->B, s->pNew, s->u);
874: SLESSolve(s->sles, s->u, s->uNew2, &its);
875: s->iter += its;
876: PetscLogInfo(pc, "PCSchur: %d iterations in A solve %dn", its, schurIts+1);
877: /* Now s->uNew = w_u - A^{-1} B S^{-1} w_p */
878: VecAXPY(&minusOne, s->uNew2, s->uNew);
879: /* Permute uNew */
880: if (s->rowPerm != PETSC_NULL) {
881: VecPermute(s->uNew, s->rowPerm, PETSC_TRUE);
882: }
883: /* Put u and p in y */
884: VecScatterBegin(s->uNew, yy, INSERT_VALUES, SCATTER_REVERSE, s->uScatter);
885: VecScatterEnd(s->uNew, yy, INSERT_VALUES, SCATTER_REVERSE, s->uScatter);
886: VecScatterBegin(s->pNew, yy, INSERT_VALUES, SCATTER_REVERSE, s->pScatter);
887: VecScatterEnd(s->pNew, yy, INSERT_VALUES, SCATTER_REVERSE, s->pScatter);
889: /* Project back: Apply P^+ = (P^T P)^{-1} P^T */
890: if ((s->isConstrained == PETSC_TRUE) && (s->explicitConstraints == PETSC_FALSE)) {
891: /* Project solution back onto constrained unknowns */
892: MatMultTranspose(P, yy, z);
893: /* Particle solve: Add G^{-1} U = w_U to z */
894: GVecSolveJacobianConstrained(z, x);
895: /* Projection Solve: Now y = (P^T P)^{-1} z */
896: MatMult(invP, z, y);
897: }
899: PetscOptionsHasName(PETSC_NULL, "-pc_schur_test", &opt);
900: if (opt == PETSC_TRUE) {
901: PetscReal solNorm;
903: /* Check that we get A^{-1} A x = x */
904: VecAXPY(&minusOne, y, tempSol);
905: VecNorm(tempSol, NORM_2, &solNorm);
906: VecDestroy(tempSol);
907: PetscPrintf(pc->comm, "Error in test solution: %gn", solNorm);
908: }
909: return(0);
910: }
912: static int PCApplyTrans_Schur(PC pc, Vec x, Vec y)
913: {
915: return(0);
916: }
918: static int PCApplySymmetricLeft_Schur(PC pc, Vec x, Vec y)
919: {
921: return(0);
922: }
924: static int PCApplySymmetricRight_Schur(PC pc, Vec x, Vec y)
925: {
927: return(0);
928: }
930: static int PCDestroy_Schur(PC pc)
931: {
932: PC_Schur *s = (PC_Schur *) pc->data;
933: int ierr;
936: /* Destroy inner solvers */
937: SLESDestroy(s->sles);
938: SLESDestroy(s->schurSles);
939: if (pc->setupcalled) {
940: PCDestroyStructures_Schur(pc);
941: }
942: if (s->mathViewer != PETSC_NULL) {
943: PetscViewerDestroy(s->mathViewer);
944: }
945: PetscFree(s);
946: return(0);
947: }
949: EXTERN_C_BEGIN
950: int PCCreate_Schur(PC pc)
951: {
952: PC_Schur *s;
953: int ierr;
956: PetscNew(PC_Schur, &s);
957: PetscLogObjectMemory(pc, sizeof(PC_Schur));
958: pc->data = (void *) s;
960: pc->ops->setup = PCSetUp_Schur;
961: pc->ops->apply = PCApply_Schur;
962: pc->ops->applyrichardson = PETSC_NULL;
963: pc->ops->applyBA = PETSC_NULL;
964: pc->ops->applytranspose = PCApplyTrans_Schur;
965: pc->ops->applyBAtranspose = PETSC_NULL;
966: pc->ops->setfromoptions = PCSetFromOptions_Schur;
967: pc->ops->presolve = PCPreSolve_Schur;
968: pc->ops->postsolve = PCPostSolve_Schur;
969: pc->ops->getfactoredmatrix = PETSC_NULL;
970: pc->ops->applysymmetricleft = PCApplySymmetricLeft_Schur;
971: pc->ops->applysymmetricright = PCApplySymmetricRight_Schur;
972: pc->ops->setuponblocks = PETSC_NULL;
973: pc->ops->destroy = PCDestroy_Schur;
974: pc->ops->view = PETSC_NULL;
975: pc->modifysubmatrices = PETSC_NULL;
977: s->grid = PETSC_NULL;
978: s->isConstrained = PETSC_FALSE;
979: s->explicitConstraints = PETSC_FALSE;
980: s->useMath = PETSC_FALSE;
981: s->mathViewer = PETSC_NULL;
983: s->u = PETSC_NULL;
984: s->uNew = PETSC_NULL;
985: s->uNew2 = PETSC_NULL;
986: s->p = PETSC_NULL;
987: s->pNew = PETSC_NULL;
988: s->uScatter = PETSC_NULL;
989: s->pScatter = PETSC_NULL;
990: s->projX = PETSC_NULL;
991: s->projY = PETSC_NULL;
992: s->projZ = PETSC_NULL;
994: s->sles = PETSC_NULL;
995: s->schurSles = PETSC_NULL;
996: s->iter = 0;
997: s->schurIter = 0;
999: s->numMomOps = -1;
1000: s->gradOp = GRADIENT;
1001: s->divOp = DIVERGENCE;
1002: s->gradOpAlpha = 0.0;
1003: s->momOps = PETSC_NULL;
1004: s->momOpIsALE = PETSC_NULL;
1005: s->momOpAlphas = PETSC_NULL;
1006: s->sField = -1;
1007: s->tField = -1;
1008: s->sOrder = PETSC_NULL;
1009: s->sLocOrder = PETSC_NULL;
1010: s->tOrder = PETSC_NULL;
1011: s->tLocOrder = PETSC_NULL;
1012: s->A = PETSC_NULL;
1013: s->sparseA = PETSC_NULL;
1014: s->B = PETSC_NULL;
1015: s->S = PETSC_NULL;
1016: s->lap = PETSC_NULL;
1017: s->rowPerm = PETSC_NULL;
1018: s->colPerm = PETSC_NULL;
1020: /* Create the inner solvers */
1021: PCSchurCreateInnerSLES(pc);
1023: return(0);
1024: }
1025: EXTERN_C_END