Actual source code: lsc.c

  2: #include <private/pcimpl.h>   /*I "petscpc.h" I*/

  4: typedef struct {
  5:   PetscBool  allocated;
  6:   PetscBool  scalediag;
  7:   KSP        kspL;
  8:   Vec        scale;
  9:   Vec        x0,y0,x1;
 10:   Mat        L;            /* keep a copy to reuse when obtained with L = A10*A01 */
 11: } PC_LSC;

 15: static PetscErrorCode PCLSCAllocate_Private(PC pc)
 16: {
 17:   PC_LSC         *lsc = (PC_LSC*)pc->data;
 18:   Mat             A;
 19:   PetscErrorCode  ierr;

 22:   if (lsc->allocated) return(0);
 23:   KSPCreate(((PetscObject)pc)->comm,&lsc->kspL);
 24:   PetscObjectIncrementTabLevel((PetscObject)lsc->kspL,(PetscObject)pc,1);
 25:   KSPSetType(lsc->kspL,KSPPREONLY);
 26:   KSPSetOptionsPrefix(lsc->kspL,((PetscObject)pc)->prefix);
 27:   KSPAppendOptionsPrefix(lsc->kspL,"lsc_");
 28:   KSPSetFromOptions(lsc->kspL);
 29:   MatSchurComplementGetSubmatrices(pc->mat,&A,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 30:   MatGetVecs(A,&lsc->x0,&lsc->y0);
 31:   MatGetVecs(pc->pmat,&lsc->x1,PETSC_NULL);
 32:   if (lsc->scalediag) {
 33:     VecDuplicate(lsc->x0,&lsc->scale);
 34:   }
 35:   lsc->allocated = PETSC_TRUE;
 36:   return(0);
 37: }

 41: static PetscErrorCode PCSetUp_LSC(PC pc)
 42: {
 43:   PC_LSC         *lsc = (PC_LSC*)pc->data;
 44:   Mat             L,Lp,B,C;
 45:   PetscErrorCode  ierr;

 48:   PCLSCAllocate_Private(pc);
 49:   PetscObjectQuery((PetscObject)pc->pmat,"LSC_L",(PetscObject*)&L);
 50:   PetscObjectQuery((PetscObject)pc->pmat,"LSC_Lp",(PetscObject*)&Lp);
 51:   if (!L) {
 52:     MatSchurComplementGetSubmatrices(pc->mat,PETSC_NULL,PETSC_NULL,&B,&C,PETSC_NULL);
 53:     if (!lsc->L) {
 54:       MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L);
 55:     } else {
 56:       MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L);
 57:     }
 58:     Lp = L = lsc->L;
 59:   }
 60:   if (lsc->scale) {
 61:     Mat Ap;
 62:     MatSchurComplementGetSubmatrices(pc->mat,PETSC_NULL,&Ap,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 63:     MatGetDiagonal(Ap,lsc->scale); /* Should be the mass matrix, but we don't have plumbing for that yet */
 64:     VecReciprocal(lsc->scale);
 65:   }
 66:   KSPSetOperators(lsc->kspL,L,Lp,SAME_NONZERO_PATTERN);
 67:   return(0);
 68: }

 72: static PetscErrorCode PCApply_LSC(PC pc,Vec x,Vec y)
 73: {
 74:   PC_LSC        *lsc = (PC_LSC*)pc->data;
 75:   Mat            A,B,C;

 79:   MatSchurComplementGetSubmatrices(pc->mat,&A,PETSC_NULL,&B,&C,PETSC_NULL);
 80:   KSPSolve(lsc->kspL,x,lsc->x1);
 81:   MatMult(B,lsc->x1,lsc->x0);
 82:   if (lsc->scale) {
 83:     VecPointwiseMult(lsc->x0,lsc->x0,lsc->scale);
 84:   }
 85:   MatMult(A,lsc->x0,lsc->y0);
 86:   if (lsc->scale) {
 87:     VecPointwiseMult(lsc->y0,lsc->y0,lsc->scale);
 88:   }
 89:   MatMult(C,lsc->y0,lsc->x1);
 90:   KSPSolve(lsc->kspL,lsc->x1,y);
 91:   return(0);
 92: }

 96: static PetscErrorCode PCReset_LSC(PC pc)
 97: {
 98:   PC_LSC         *lsc = (PC_LSC*)pc->data;

102:   VecDestroy(&lsc->x0);
103:   VecDestroy(&lsc->y0);
104:   VecDestroy(&lsc->x1);
105:   VecDestroy(&lsc->scale);
106:   KSPDestroy(&lsc->kspL);
107:   MatDestroy(&lsc->L);
108:   return(0);
109: }

113: static PetscErrorCode PCDestroy_LSC(PC pc)
114: {

118:   PCReset_LSC(pc);
119:   PetscFree(pc->data);
120:   return(0);
121: }

125: static PetscErrorCode PCSetFromOptions_LSC(PC pc)
126: {
127:   PC_LSC         *lsc = (PC_LSC*)pc->data;
128:   PetscErrorCode  ierr;

131:   PetscOptionsHead("LSC options");
132:   {
133:     PetscOptionsBool("-pc_lsc_scale_diag","Use diagonal of velocity block (A) for scaling","None",lsc->scalediag,&lsc->scalediag,PETSC_NULL);
134:   }
135:   PetscOptionsTail();
136:   return(0);
137: }

141: static PetscErrorCode PCView_LSC(PC pc,PetscViewer viewer)
142: {
143:   PC_LSC           *jac = (PC_LSC*)pc->data;
144:   PetscErrorCode   ierr;
145:   PetscBool        iascii;

148:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
149:   if (iascii) {
150:     PetscViewerASCIIPushTab(viewer);
151:     KSPView(jac->kspL,viewer);
152:     PetscViewerASCIIPopTab(viewer);
153:   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for LSC",((PetscObject)viewer)->type_name);
154:   return(0);
155: }

157: /*MC
158:      PCLSC - Preconditioning for Schur complements, based on Least Squares Commutators

160:    Options Database Key:
161: .    -pc_lsc_scale_diag - Use the diagonal of A for scaling

163:    Level: intermediate

165:    Notes:
166:    This preconditioner will normally be used with PCFieldSplit to precondition the Schur complement, but
167:    it can be used for any Schur complement system.  Consider the Schur complement

169: .vb
170:    S = A11 - A10 inv(A00) A01
171: .ve

173:    PCLSC currently doesn't do anything with A11, so let's assume it is 0.  The idea is that a good approximation to
174:    inv(S) is given by

176: .vb
177:    inv(A10 A01) A10 A00 A01 inv(A10 A01)
178: .ve

180:    The product A10 A01 can be computed for you, but you can provide it (this is
181:    usually more efficient anyway).  In the case of incompressible flow, A10 A10 is a Laplacian, call it L.  The current
182:    interface is to hang L and a preconditioning matrix Lp on the preconditioning matrix.

184:    If you had called KSPSetOperators(ksp,S,Sp,flg), S should have type MATSCHURCOMPLEMENT and Sp can be any type you
185:    like (PCLSC doesn't use it directly) but should have matrices composed with it, under the names "LSC_L" and "LSC_Lp".
186:    For example, you might have setup code like this

188: .vb
189:    PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L);
190:    PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp);
191: .ve

193:    And then your Jacobian assembly would look like

195: .vb
196:    PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L);
197:    PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp);
198:    if (L) { assembly L }
199:    if (Lp) { assemble Lp }
200: .ve

202:    With this, you should be able to choose LSC preconditioning, using e.g. ML's algebraic multigrid to solve with L

204: .vb
205:    -fieldsplit_1_pc_type lsc -fieldsplit_1_lsc_pc_type ml
206: .ve

208:    Since we do not use the values in Sp, you can still put an assembled matrix there to use normal preconditioners.

210:    Concepts: physics based preconditioners, block preconditioners

212: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners, PCFIELDSPLIT,
213:            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition(),
214:            MatCreateSchurComplement()
215: M*/

220: PetscErrorCode  PCCreate_LSC(PC pc)
221: {
222:   PC_LSC         *lsc;

226:   PetscNewLog(pc,PC_LSC,&lsc);
227:   pc->data  = (void*)lsc;

229:   pc->ops->apply               = PCApply_LSC;
230:   pc->ops->applytranspose      = 0;
231:   pc->ops->setup               = PCSetUp_LSC;
232:   pc->ops->reset               = PCReset_LSC;
233:   pc->ops->destroy             = PCDestroy_LSC;
234:   pc->ops->setfromoptions      = PCSetFromOptions_LSC;
235:   pc->ops->view                = PCView_LSC;
236:   pc->ops->applyrichardson     = 0;
237:   return(0);
238: }