Actual source code: mg.c

petsc-main 2021-03-05
Report Typos and Errors

  2: /*
  3:     Defines the multigrid preconditioner interface.
  4: */
  5: #include <petsc/private/pcmgimpl.h>
  6: #include <petsc/private/kspimpl.h>
  7: #include <petscdm.h>
  8: PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*);

 10: /*
 11:    Contains the list of registered coarse space construction routines
 12: */
 13: PetscFunctionList PCMGCoarseList = NULL;

 15: PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PetscBool transpose,PCRichardsonConvergedReason *reason)
 16: {
 17:   PC_MG          *mg = (PC_MG*)pc->data;
 18:   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
 20:   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;

 23:   if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
 24:   if (!transpose) {
 25:     KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);  /* pre-smooth */
 26:     KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);
 27:   } else {
 28:     KSPSolveTranspose(mglevels->smoothu,mglevels->b,mglevels->x); /* transpose of post-smooth */
 29:     KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);
 30:   }
 31:   if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
 32:   if (mglevels->level) {  /* not the coarsest grid */
 33:     if (mglevels->eventresidual) {PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);}
 34:     if (!transpose) {
 35:       (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);
 36:     } else {
 37:       (*mglevels->residualtranspose)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);
 38:     }
 39:     if (mglevels->eventresidual) {PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);}

 41:     /* if on finest level and have convergence criteria set */
 42:     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
 43:       PetscReal rnorm;
 44:       VecNorm(mglevels->r,NORM_2,&rnorm);
 45:       if (rnorm <= mg->ttol) {
 46:         if (rnorm < mg->abstol) {
 47:           *reason = PCRICHARDSON_CONVERGED_ATOL;
 48:           PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);
 49:         } else {
 50:           *reason = PCRICHARDSON_CONVERGED_RTOL;
 51:           PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);
 52:         }
 53:         return(0);
 54:       }
 55:     }

 57:     mgc = *(mglevelsin - 1);
 58:     if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
 59:     if (!transpose) {
 60:       MatRestrict(mglevels->restrct,mglevels->r,mgc->b);
 61:     } else {
 62:       MatRestrict(mglevels->interpolate,mglevels->r,mgc->b);
 63:     }
 64:     if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
 65:     VecSet(mgc->x,0.0);
 66:     while (cycles--) {
 67:       PCMGMCycle_Private(pc,mglevelsin-1,transpose,reason);
 68:     }
 69:     if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
 70:     if (!transpose) {
 71:       MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);
 72:     } else {
 73:       MatInterpolateAdd(mglevels->restrct,mgc->x,mglevels->x,mglevels->x);
 74:     }
 75:     if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
 76:     if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
 77:     if (!transpose) {
 78:       KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);    /* post smooth */
 79:       KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);
 80:     } else {
 81:       KSPSolveTranspose(mglevels->smoothd,mglevels->b,mglevels->x);    /* post smooth */
 82:       KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);
 83:     }
 84:     if (mglevels->cr) {
 85:       /* TODO Turn on copy and turn off noisy if we have an exact solution
 86:       VecCopy(mglevels->x, mglevels->crx);
 87:       VecCopy(mglevels->b, mglevels->crb); */
 88:       KSPSetNoisy_Private(mglevels->crx);
 89:       KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx);    /* compatible relaxation */
 90:       KSPCheckSolve(mglevels->cr,pc,mglevels->crx);
 91:     }
 92:     if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
 93:   }
 94:   return(0);
 95: }

 97: static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
 98: {
 99:   PC_MG          *mg        = (PC_MG*)pc->data;
100:   PC_MG_Levels   **mglevels = mg->levels;
102:   PC             tpc;
103:   PetscBool      changeu,changed;
104:   PetscInt       levels = mglevels[0]->levels,i;

107:   /* When the DM is supplying the matrix then it will not exist until here */
108:   for (i=0; i<levels; i++) {
109:     if (!mglevels[i]->A) {
110:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
111:       PetscObjectReference((PetscObject)mglevels[i]->A);
112:     }
113:   }

115:   KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
116:   PCPreSolveChangeRHS(tpc,&changed);
117:   KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
118:   PCPreSolveChangeRHS(tpc,&changeu);
119:   if (!changed && !changeu) {
120:     VecDestroy(&mglevels[levels-1]->b);
121:     mglevels[levels-1]->b = b;
122:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
123:     if (!mglevels[levels-1]->b) {
124:       Vec *vec;

126:       KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
127:       mglevels[levels-1]->b = *vec;
128:       PetscFree(vec);
129:     }
130:     VecCopy(b,mglevels[levels-1]->b);
131:   }
132:   mglevels[levels-1]->x = x;

134:   mg->rtol   = rtol;
135:   mg->abstol = abstol;
136:   mg->dtol   = dtol;
137:   if (rtol) {
138:     /* compute initial residual norm for relative convergence test */
139:     PetscReal rnorm;
140:     if (zeroguess) {
141:       VecNorm(b,NORM_2,&rnorm);
142:     } else {
143:       (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);
144:       VecNorm(w,NORM_2,&rnorm);
145:     }
146:     mg->ttol = PetscMax(rtol*rnorm,abstol);
147:   } else if (abstol) mg->ttol = abstol;
148:   else mg->ttol = 0.0;

150:   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
151:      stop prematurely due to small residual */
152:   for (i=1; i<levels; i++) {
153:     KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
154:     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
155:       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
156:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
157:       KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
158:     }
159:   }

161:   *reason = (PCRichardsonConvergedReason)0;
162:   for (i=0; i<its; i++) {
163:     PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_FALSE,reason);
164:     if (*reason) break;
165:   }
166:   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
167:   *outits = i;
168:   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
169:   return(0);
170: }

172: PetscErrorCode PCReset_MG(PC pc)
173: {
174:   PC_MG          *mg        = (PC_MG*)pc->data;
175:   PC_MG_Levels   **mglevels = mg->levels;
177:   PetscInt       i,c,n;

180:   if (mglevels) {
181:     n = mglevels[0]->levels;
182:     for (i=0; i<n-1; i++) {
183:       VecDestroy(&mglevels[i+1]->r);
184:       VecDestroy(&mglevels[i]->b);
185:       VecDestroy(&mglevels[i]->x);
186:       VecDestroy(&mglevels[i]->crx);
187:       VecDestroy(&mglevels[i]->crb);
188:       MatDestroy(&mglevels[i+1]->restrct);
189:       MatDestroy(&mglevels[i+1]->interpolate);
190:       MatDestroy(&mglevels[i+1]->inject);
191:       VecDestroy(&mglevels[i+1]->rscale);
192:     }
193:     VecDestroy(&mglevels[n-1]->crx);
194:     VecDestroy(&mglevels[n-1]->crb);
195:     /* this is not null only if the smoother on the finest level
196:        changes the rhs during PreSolve */
197:     VecDestroy(&mglevels[n-1]->b);

199:     for (i=0; i<n; i++) {
200:       if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {VecDestroy(&mglevels[i]->coarseSpace[c]);}
201:       PetscFree(mglevels[i]->coarseSpace);
202:       mglevels[i]->coarseSpace = NULL;
203:       MatDestroy(&mglevels[i]->A);
204:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
205:         KSPReset(mglevels[i]->smoothd);
206:       }
207:       KSPReset(mglevels[i]->smoothu);
208:       if (mglevels[i]->cr) {KSPReset(mglevels[i]->cr);}
209:     }
210:     mg->Nc = 0;
211:   }
212:   return(0);
213: }

215: /* Implementing CR

217: We only want to make corrections that ``do not change'' the coarse solution. What we mean by not changing is that if I prolong my coarse solution to the fine grid and then inject that fine solution back to the coarse grid, I get the same answer. Injection is what Brannick calls R. We want the complementary projector to Inj, which we will call S, after Brannick, so that Inj S = 0. Now the orthogonal projector onto the range of Inj^T is

219:   Inj^T (Inj Inj^T)^{-1} Inj

221: and if Inj is a VecScatter, as it is now in PETSc, we have

223:   Inj^T Inj

225: and

227:   S = I - Inj^T Inj

229: since

231:   Inj S = Inj - (Inj Inj^T) Inj = 0.

233: Brannick suggests

235:   A \to S^T A S  \qquad\mathrm{and}\qquad M \to S^T M S

237: but I do not think his :math:`S^T S = I` is correct. Our S is an orthogonal projector, so :math:`S^T S = S^2 = S`. We will use

239:   M^{-1} A \to S M^{-1} A S

241: In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left.

243:   Check: || Inj P - I ||_F < tol
244:   Check: In general, Inj Inj^T = I
245: */

247: typedef struct {
248:   PC       mg;  /* The PCMG object */
249:   PetscInt l;   /* The multigrid level for this solver */
250:   Mat      Inj; /* The injection matrix */
251:   Mat      S;   /* I - Inj^T Inj */
252: } CRContext;

254: static PetscErrorCode CRSetup_Private(PC pc)
255: {
256:   CRContext     *ctx;
257:   Mat            It;

261:   PCShellGetContext(pc, (void **) &ctx);
262:   PCMGGetInjection(ctx->mg, ctx->l, &It);
263:   if (!It) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG");
264:   MatCreateTranspose(It, &ctx->Inj);
265:   MatCreateNormal(ctx->Inj, &ctx->S);
266:   MatScale(ctx->S, -1.0);
267:   MatShift(ctx->S,  1.0);
268:   return(0);
269: }

271: static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y)
272: {
273:   CRContext     *ctx;

277:   PCShellGetContext(pc, (void **) &ctx);
278:   MatMult(ctx->S, x, y);
279:   return(0);
280: }

282: static PetscErrorCode CRDestroy_Private(PC pc)
283: {
284:   CRContext     *ctx;

288:   PCShellGetContext(pc, (void **) &ctx);
289:   MatDestroy(&ctx->Inj);
290:   MatDestroy(&ctx->S);
291:   PetscFree(ctx);
292:   PCShellSetContext(pc, NULL);
293:   return(0);
294: }

296: static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr)
297: {
298:   CRContext     *ctx;

302:   PCCreate(PetscObjectComm((PetscObject) pc), cr);
303:   PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)");
304:   PetscCalloc1(1, &ctx);
305:   ctx->mg = pc;
306:   ctx->l  = l;
307:   PCSetType(*cr, PCSHELL);
308:   PCShellSetContext(*cr, ctx);
309:   PCShellSetApply(*cr, CRApply_Private);
310:   PCShellSetSetUp(*cr, CRSetup_Private);
311:   PCShellSetDestroy(*cr, CRDestroy_Private);
312:   return(0);
313: }

315: PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
316: {
318:   PC_MG          *mg        = (PC_MG*)pc->data;
319:   MPI_Comm       comm;
320:   PC_MG_Levels   **mglevels = mg->levels;
321:   PCMGType       mgtype     = mg->am;
322:   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
323:   PetscInt       i;
324:   PetscMPIInt    size;
325:   const char     *prefix;
326:   PC             ipc;
327:   PetscInt       n;

332:   if (mg->nlevels == levels) return(0);
333:   PetscObjectGetComm((PetscObject)pc,&comm);
334:   if (mglevels) {
335:     mgctype = mglevels[0]->cycles;
336:     /* changing the number of levels so free up the previous stuff */
337:     PCReset_MG(pc);
338:     n    = mglevels[0]->levels;
339:     for (i=0; i<n; i++) {
340:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
341:         KSPDestroy(&mglevels[i]->smoothd);
342:       }
343:       KSPDestroy(&mglevels[i]->smoothu);
344:       KSPDestroy(&mglevels[i]->cr);
345:       PetscFree(mglevels[i]);
346:     }
347:     PetscFree(mg->levels);
348:   }

350:   mg->nlevels = levels;

352:   PetscMalloc1(levels,&mglevels);
353:   PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));

355:   PCGetOptionsPrefix(pc,&prefix);

357:   mg->stageApply = 0;
358:   for (i=0; i<levels; i++) {
359:     PetscNewLog(pc,&mglevels[i]);

361:     mglevels[i]->level               = i;
362:     mglevels[i]->levels              = levels;
363:     mglevels[i]->cycles              = mgctype;
364:     mg->default_smoothu              = 2;
365:     mg->default_smoothd              = 2;
366:     mglevels[i]->eventsmoothsetup    = 0;
367:     mglevels[i]->eventsmoothsolve    = 0;
368:     mglevels[i]->eventresidual       = 0;
369:     mglevels[i]->eventinterprestrict = 0;

371:     if (comms) comm = comms[i];
372:     if (comm != MPI_COMM_NULL) {
373:       KSPCreate(comm,&mglevels[i]->smoothd);
374:       KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);
375:       PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);
376:       KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);
377:       PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);
378:       if (i || levels == 1) {
379:         char tprefix[128];

381:         KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);
382:         KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);
383:         KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);
384:         KSPGetPC(mglevels[i]->smoothd,&ipc);
385:         PCSetType(ipc,PCSOR);
386:         KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);

388:         PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);
389:         KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);
390:       } else {
391:         KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");

393:         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
394:         KSPSetType(mglevels[0]->smoothd,KSPPREONLY);
395:         KSPGetPC(mglevels[0]->smoothd,&ipc);
396:         MPI_Comm_size(comm,&size);
397:         if (size > 1) {
398:           PCSetType(ipc,PCREDUNDANT);
399:         } else {
400:           PCSetType(ipc,PCLU);
401:         }
402:         PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);
403:       }
404:       PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);
405:     }
406:     mglevels[i]->smoothu = mglevels[i]->smoothd;
407:     mg->rtol             = 0.0;
408:     mg->abstol           = 0.0;
409:     mg->dtol             = 0.0;
410:     mg->ttol             = 0.0;
411:     mg->cyclesperpcapply = 1;
412:   }
413:   mg->levels = mglevels;
414:   PCMGSetType(pc,mgtype);
415:   return(0);
416: }

418: /*@C
419:    PCMGSetLevels - Sets the number of levels to use with MG.
420:    Must be called before any other MG routine.

422:    Logically Collective on PC

424:    Input Parameters:
425: +  pc - the preconditioner context
426: .  levels - the number of levels
427: -  comms - optional communicators for each level; this is to allow solving the coarser problems
428:            on smaller sets of processes. For processes that are not included in the computation
429:            you must pass MPI_COMM_NULL.

431:    Level: intermediate

433:    Notes:
434:      If the number of levels is one then the multigrid uses the -mg_levels prefix
435:      for setting the level options rather than the -mg_coarse prefix.

437:      You can free the information in comms after this routine is called.

439:      The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level
440:      are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on
441:      the two levels, rank 0 in the original communicator will pass in an array of 2 communicators
442:      of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators
443:      the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate
444:      in the coarse grid solve.

446:      Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one
447:      must take special care in providing the restriction and interpolation operation. We recommend
448:      providing these as two step operations; first perform a standard restriction or interpolation on
449:      the full number of ranks for that level and then use an MPI call to copy the resulting vector
450:      array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, not in both
451:      cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and
452:      recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries.


455: .seealso: PCMGSetType(), PCMGGetLevels()
456: @*/
457: PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
458: {

464:   PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));
465:   return(0);
466: }


469: PetscErrorCode PCDestroy_MG(PC pc)
470: {
472:   PC_MG          *mg        = (PC_MG*)pc->data;
473:   PC_MG_Levels   **mglevels = mg->levels;
474:   PetscInt       i,n;

477:   PCReset_MG(pc);
478:   if (mglevels) {
479:     n = mglevels[0]->levels;
480:     for (i=0; i<n; i++) {
481:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
482:         KSPDestroy(&mglevels[i]->smoothd);
483:       }
484:       KSPDestroy(&mglevels[i]->smoothu);
485:       KSPDestroy(&mglevels[i]->cr);
486:       PetscFree(mglevels[i]);
487:     }
488:     PetscFree(mg->levels);
489:   }
490:   PetscFree(pc->data);
491:   PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);
492:   PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);
493:   return(0);
494: }


497: /*
498:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
499:              or full cycle of multigrid.

501:   Note:
502:   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
503: */
504: static PetscErrorCode PCApply_MG_Internal(PC pc,Vec b,Vec x,PetscBool transpose)
505: {
506:   PC_MG          *mg        = (PC_MG*)pc->data;
507:   PC_MG_Levels   **mglevels = mg->levels;
509:   PC             tpc;
510:   PetscInt       levels = mglevels[0]->levels,i;
511:   PetscBool      changeu,changed;

514:   if (mg->stageApply) {PetscLogStagePush(mg->stageApply);}
515:   /* When the DM is supplying the matrix then it will not exist until here */
516:   for (i=0; i<levels; i++) {
517:     if (!mglevels[i]->A) {
518:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
519:       PetscObjectReference((PetscObject)mglevels[i]->A);
520:     }
521:   }

523:   KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
524:   PCPreSolveChangeRHS(tpc,&changed);
525:   KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
526:   PCPreSolveChangeRHS(tpc,&changeu);
527:   if (!changeu && !changed) {
528:     VecDestroy(&mglevels[levels-1]->b);
529:     mglevels[levels-1]->b = b;
530:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
531:     if (!mglevels[levels-1]->b) {
532:       Vec *vec;

534:       KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
535:       mglevels[levels-1]->b = *vec;
536:       PetscFree(vec);
537:     }
538:     VecCopy(b,mglevels[levels-1]->b);
539:   }
540:   mglevels[levels-1]->x = x;

542:   if (mg->am == PC_MG_MULTIPLICATIVE) {
543:     VecZeroEntries(x);
544:     for (i=0; i<mg->cyclesperpcapply; i++) {
545:       PCMGMCycle_Private(pc,mglevels+levels-1,transpose,NULL);
546:     }
547:   } else if (mg->am == PC_MG_ADDITIVE) {
548:     PCMGACycle_Private(pc,mglevels,transpose);
549:   } else if (mg->am == PC_MG_KASKADE) {
550:     PCMGKCycle_Private(pc,mglevels,transpose);
551:   } else {
552:     PCMGFCycle_Private(pc,mglevels,transpose);
553:   }
554:   if (mg->stageApply) {PetscLogStagePop();}
555:   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
556:   return(0);
557: }

559: static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
560: {

564:     PCApply_MG_Internal(pc,b,x,PETSC_FALSE);
565:     return(0);
566: }

568: static PetscErrorCode PCApplyTranspose_MG(PC pc,Vec b,Vec x)
569: {

573:     PCApply_MG_Internal(pc,b,x,PETSC_TRUE);
574:     return(0);
575: }

577: PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
578: {
579:   PetscErrorCode   ierr;
580:   PetscInt         levels,cycles;
581:   PetscBool        flg, flg2;
582:   PC_MG            *mg = (PC_MG*)pc->data;
583:   PC_MG_Levels     **mglevels;
584:   PCMGType         mgtype;
585:   PCMGCycleType    mgctype;
586:   PCMGGalerkinType gtype;

589:   levels = PetscMax(mg->nlevels,1);
590:   PetscOptionsHead(PetscOptionsObject,"Multigrid options");
591:   PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
592:   if (!flg && !mg->levels && pc->dm) {
593:     DMGetRefineLevel(pc->dm,&levels);
594:     levels++;
595:     mg->usedmfornumberoflevels = PETSC_TRUE;
596:   }
597:   PCMGSetLevels(pc,levels,NULL);
598:   mglevels = mg->levels;

600:   mgctype = (PCMGCycleType) mglevels[0]->cycles;
601:   PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);
602:   if (flg) {
603:     PCMGSetCycleType(pc,mgctype);
604:   }
605:   gtype = mg->galerkin;
606:   PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);
607:   if (flg) {
608:     PCMGSetGalerkin(pc,gtype);
609:   }
610:   flg2 = PETSC_FALSE;
611:   PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);
612:   if (flg) {PCMGSetAdaptInterpolation(pc, flg2);}
613:   PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg);
614:   PetscOptionsEnum("-pc_mg_adapt_interp_coarse_space","Type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector","PCMGSetAdaptCoarseSpaceType",PCMGCoarseSpaceTypes,(PetscEnum)mg->coarseSpaceType,(PetscEnum*)&mg->coarseSpaceType,&flg);
615:   PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);
616:   flg2 = PETSC_FALSE;
617:   PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg);
618:   if (flg) {PCMGSetAdaptCR(pc, flg2);}
619:   flg = PETSC_FALSE;
620:   PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);
621:   if (flg) {
622:     PCMGSetDistinctSmoothUp(pc);
623:   }
624:   mgtype = mg->am;
625:   PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
626:   if (flg) {
627:     PCMGSetType(pc,mgtype);
628:   }
629:   if (mg->am == PC_MG_MULTIPLICATIVE) {
630:     PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);
631:     if (flg) {
632:       PCMGMultiplicativeSetCycles(pc,cycles);
633:     }
634:   }
635:   flg  = PETSC_FALSE;
636:   PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);
637:   if (flg) {
638:     PetscInt i;
639:     char     eventname[128];

641:     levels = mglevels[0]->levels;
642:     for (i=0; i<levels; i++) {
643:       sprintf(eventname,"MGSetup Level %d",(int)i);
644:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);
645:       sprintf(eventname,"MGSmooth Level %d",(int)i);
646:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);
647:       if (i) {
648:         sprintf(eventname,"MGResid Level %d",(int)i);
649:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);
650:         sprintf(eventname,"MGInterp Level %d",(int)i);
651:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);
652:       }
653:     }

655: #if defined(PETSC_USE_LOG)
656:     {
657:       const char    *sname = "MG Apply";
658:       PetscStageLog stageLog;
659:       PetscInt      st;

661:       PetscLogGetStageLog(&stageLog);
662:       for (st = 0; st < stageLog->numStages; ++st) {
663:         PetscBool same;

665:         PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
666:         if (same) mg->stageApply = st;
667:       }
668:       if (!mg->stageApply) {
669:         PetscLogStageRegister(sname, &mg->stageApply);
670:       }
671:     }
672: #endif
673:   }
674:   PetscOptionsTail();
675:   /* Check option consistency */
676:   PCMGGetGalerkin(pc, &gtype);
677:   PCMGGetAdaptInterpolation(pc, &flg);
678:   if (flg && (gtype >= PC_MG_GALERKIN_NONE)) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
679:   return(0);
680: }

682: const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL};
683: const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL};
684: const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL};
685: const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL};

687: #include <petscdraw.h>
688: PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
689: {
690:   PC_MG          *mg        = (PC_MG*)pc->data;
691:   PC_MG_Levels   **mglevels = mg->levels;
693:   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
694:   PetscBool      iascii,isbinary,isdraw;

697:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
698:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
699:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
700:   if (iascii) {
701:     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
702:     PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);
703:     if (mg->am == PC_MG_MULTIPLICATIVE) {
704:       PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);
705:     }
706:     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
707:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");
708:     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
709:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");
710:     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
711:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");
712:     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
713:       PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");
714:     } else {
715:       PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");
716:     }
717:     if (mg->view){
718:       (*mg->view)(pc,viewer);
719:     }
720:     for (i=0; i<levels; i++) {
721:       if (!i) {
722:         PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);
723:       } else {
724:         PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);
725:       }
726:       PetscViewerASCIIPushTab(viewer);
727:       KSPView(mglevels[i]->smoothd,viewer);
728:       PetscViewerASCIIPopTab(viewer);
729:       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
730:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
731:       } else if (i) {
732:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);
733:         PetscViewerASCIIPushTab(viewer);
734:         KSPView(mglevels[i]->smoothu,viewer);
735:         PetscViewerASCIIPopTab(viewer);
736:       }
737:       if (i && mglevels[i]->cr) {
738:         PetscViewerASCIIPrintf(viewer,"CR solver on level %D -------------------------------\n",i);
739:         PetscViewerASCIIPushTab(viewer);
740:         KSPView(mglevels[i]->cr,viewer);
741:         PetscViewerASCIIPopTab(viewer);
742:       }
743:     }
744:   } else if (isbinary) {
745:     for (i=levels-1; i>=0; i--) {
746:       KSPView(mglevels[i]->smoothd,viewer);
747:       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
748:         KSPView(mglevels[i]->smoothu,viewer);
749:       }
750:     }
751:   } else if (isdraw) {
752:     PetscDraw draw;
753:     PetscReal x,w,y,bottom,th;
754:     PetscViewerDrawGetDraw(viewer,0,&draw);
755:     PetscDrawGetCurrentPoint(draw,&x,&y);
756:     PetscDrawStringGetSize(draw,NULL,&th);
757:     bottom = y - th;
758:     for (i=levels-1; i>=0; i--) {
759:       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
760:         PetscDrawPushCurrentPoint(draw,x,bottom);
761:         KSPView(mglevels[i]->smoothd,viewer);
762:         PetscDrawPopCurrentPoint(draw);
763:       } else {
764:         w    = 0.5*PetscMin(1.0-x,x);
765:         PetscDrawPushCurrentPoint(draw,x+w,bottom);
766:         KSPView(mglevels[i]->smoothd,viewer);
767:         PetscDrawPopCurrentPoint(draw);
768:         PetscDrawPushCurrentPoint(draw,x-w,bottom);
769:         KSPView(mglevels[i]->smoothu,viewer);
770:         PetscDrawPopCurrentPoint(draw);
771:       }
772:       PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);
773:       bottom -= th;
774:     }
775:   }
776:   return(0);
777: }

779: #include <petsc/private/kspimpl.h>

781: /*
782:     Calls setup for the KSP on each level
783: */
784: PetscErrorCode PCSetUp_MG(PC pc)
785: {
786:   PC_MG          *mg        = (PC_MG*)pc->data;
787:   PC_MG_Levels   **mglevels = mg->levels;
789:   PetscInt       i,n;
790:   PC             cpc;
791:   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
792:   Mat            dA,dB;
793:   Vec            tvec;
794:   DM             *dms;
795:   PetscViewer    viewer = NULL;
796:   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE;

799:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
800:   n = mglevels[0]->levels;
801:   /* FIX: Move this to PCSetFromOptions_MG? */
802:   if (mg->usedmfornumberoflevels) {
803:     PetscInt levels;
804:     DMGetRefineLevel(pc->dm,&levels);
805:     levels++;
806:     if (levels > n) { /* the problem is now being solved on a finer grid */
807:       PCMGSetLevels(pc,levels,NULL);
808:       n        = levels;
809:       PCSetFromOptions(pc); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
810:       mglevels = mg->levels;
811:     }
812:   }
813:   KSPGetPC(mglevels[0]->smoothd,&cpc);


816:   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
817:   /* so use those from global PC */
818:   /* Is this what we always want? What if user wants to keep old one? */
819:   KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);
820:   if (opsset) {
821:     Mat mmat;
822:     KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);
823:     if (mmat == pc->pmat) opsset = PETSC_FALSE;
824:   }

826:   /* Create CR solvers */
827:   PCMGGetAdaptCR(pc, &doCR);
828:   if (doCR) {
829:     const char *prefix;

831:     PCGetOptionsPrefix(pc, &prefix);
832:     for (i = 1; i < n; ++i) {
833:       PC   ipc, cr;
834:       char crprefix[128];

836:       KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr);
837:       KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE);
838:       PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i);
839:       KSPSetOptionsPrefix(mglevels[i]->cr, prefix);
840:       PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level);
841:       KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV);
842:       KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL);
843:       KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED);
844:       KSPGetPC(mglevels[i]->cr, &ipc);

846:       PCSetType(ipc, PCCOMPOSITE);
847:       PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE);
848:       PCCompositeAddPCType(ipc, PCSOR);
849:       CreateCR_Private(pc, i, &cr);
850:       PCCompositeAddPC(ipc, cr);
851:       PCDestroy(&cr);

853:       KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd);
854:       KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE);
855:       PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i);
856:       KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix);
857:       PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr);
858:     }
859:   }

861:   if (!opsset) {
862:     PCGetUseAmat(pc,&use_amat);
863:     if (use_amat) {
864:       PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
865:       KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);
866:     } else {
867:       PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
868:       KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);
869:     }
870:   }

872:   for (i=n-1; i>0; i--) {
873:     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
874:       missinginterpolate = PETSC_TRUE;
875:       continue;
876:     }
877:   }

879:   KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);
880:   if (dA == dB) dAeqdB = PETSC_TRUE;
881:   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
882:     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
883:   }


886:   /*
887:    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
888:    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
889:   */
890:   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
891:         /* construct the interpolation from the DMs */
892:     Mat p;
893:     Vec rscale;
894:     PetscMalloc1(n,&dms);
895:     dms[n-1] = pc->dm;
896:     /* Separately create them so we do not get DMKSP interference between levels */
897:     for (i=n-2; i>-1; i--) {DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);}
898:         /*
899:            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
900:            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
901:            But it is safe to use -dm_mat_type.

903:            The mat type should not be hardcoded like this, we need to find a better way.
904:     DMSetMatType(dms[0],MATAIJ);
905:     */
906:     for (i=n-2; i>-1; i--) {
907:       DMKSP     kdm;
908:       PetscBool dmhasrestrict, dmhasinject;
909:       KSPSetDM(mglevels[i]->smoothd,dms[i]);
910:       if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);}
911:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
912:         KSPSetDM(mglevels[i]->smoothu,dms[i]);
913:         if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);}
914:       }
915:       if (mglevels[i]->cr) {
916:         KSPSetDM(mglevels[i]->cr,dms[i]);
917:         if (!needRestricts) {KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE);}
918:       }
919:       DMGetDMKSPWrite(dms[i],&kdm);
920:       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
921:        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
922:       kdm->ops->computerhs = NULL;
923:       kdm->rhsctx          = NULL;
924:       if (!mglevels[i+1]->interpolate) {
925:         DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);
926:         PCMGSetInterpolation(pc,i+1,p);
927:         if (rscale) {PCMGSetRScale(pc,i+1,rscale);}
928:         VecDestroy(&rscale);
929:         MatDestroy(&p);
930:       }
931:       DMHasCreateRestriction(dms[i],&dmhasrestrict);
932:       if (dmhasrestrict && !mglevels[i+1]->restrct){
933:         DMCreateRestriction(dms[i],dms[i+1],&p);
934:         PCMGSetRestriction(pc,i+1,p);
935:         MatDestroy(&p);
936:       }
937:       DMHasCreateInjection(dms[i],&dmhasinject);
938:       if (dmhasinject && !mglevels[i+1]->inject){
939:         DMCreateInjection(dms[i],dms[i+1],&p);
940:         PCMGSetInjection(pc,i+1,p);
941:         MatDestroy(&p);
942:       }
943:     }

945:     for (i=n-2; i>-1; i--) {DMDestroy(&dms[i]);}
946:     PetscFree(dms);
947:   }

949:   if (pc->dm && !pc->setupcalled) {
950:     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
951:     KSPSetDM(mglevels[n-1]->smoothd,pc->dm);
952:     KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);
953:     if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
954:       KSPSetDM(mglevels[n-1]->smoothu,pc->dm);
955:       KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);
956:     }
957:     if (mglevels[n-1]->cr) {
958:       KSPSetDM(mglevels[n-1]->cr,pc->dm);
959:       KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE);
960:     }
961:   }

963:   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
964:     Mat       A,B;
965:     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
966:     MatReuse  reuse = MAT_INITIAL_MATRIX;

968:     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
969:     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
970:     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
971:     for (i=n-2; i>-1; i--) {
972:       if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0");
973:       if (!mglevels[i+1]->interpolate) {
974:         PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);
975:       }
976:       if (!mglevels[i+1]->restrct) {
977:         PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);
978:       }
979:       if (reuse == MAT_REUSE_MATRIX) {
980:         KSPGetOperators(mglevels[i]->smoothd,&A,&B);
981:       }
982:       if (doA) {
983:         MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);
984:       }
985:       if (doB) {
986:         MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);
987:       }
988:       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
989:       if (!doA && dAeqdB) {
990:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)B);}
991:         A = B;
992:       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
993:         KSPGetOperators(mglevels[i]->smoothd,&A,NULL);
994:         PetscObjectReference((PetscObject)A);
995:       }
996:       if (!doB && dAeqdB) {
997:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)A);}
998:         B = A;
999:       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
1000:         KSPGetOperators(mglevels[i]->smoothd,NULL,&B);
1001:         PetscObjectReference((PetscObject)B);
1002:       }
1003:       if (reuse == MAT_INITIAL_MATRIX) {
1004:         KSPSetOperators(mglevels[i]->smoothd,A,B);
1005:         PetscObjectDereference((PetscObject)A);
1006:         PetscObjectDereference((PetscObject)B);
1007:       }
1008:       dA = A;
1009:       dB = B;
1010:     }
1011:   }


1014:   /* Adapt interpolation matrices */
1015:   if (mg->adaptInterpolation) {
1016:     mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */
1017:     for (i = 0; i < n; ++i) {
1018:       PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);
1019:       if (i) {PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mg->Nc, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace);}
1020:     }
1021:     for (i = n-2; i > -1; --i) {
1022:       PCMGRecomputeLevelOperators_Internal(pc, i);
1023:     }
1024:   }

1026:   if (needRestricts && pc->dm) {
1027:     for (i=n-2; i>=0; i--) {
1028:       DM  dmfine,dmcoarse;
1029:       Mat Restrict,Inject;
1030:       Vec rscale;
1031:       KSPGetDM(mglevels[i+1]->smoothd,&dmfine);
1032:       KSPGetDM(mglevels[i]->smoothd,&dmcoarse);
1033:       PCMGGetRestriction(pc,i+1,&Restrict);
1034:       PCMGGetRScale(pc,i+1,&rscale);
1035:       PCMGGetInjection(pc,i+1,&Inject);
1036:       DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);
1037:     }
1038:   }

1040:   if (!pc->setupcalled) {
1041:     for (i=0; i<n; i++) {
1042:       KSPSetFromOptions(mglevels[i]->smoothd);
1043:     }
1044:     for (i=1; i<n; i++) {
1045:       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
1046:         KSPSetFromOptions(mglevels[i]->smoothu);
1047:       }
1048:       if (mglevels[i]->cr) {
1049:         KSPSetFromOptions(mglevels[i]->cr);
1050:       }
1051:     }
1052:     /* insure that if either interpolation or restriction is set the other other one is set */
1053:     for (i=1; i<n; i++) {
1054:       PCMGGetInterpolation(pc,i,NULL);
1055:       PCMGGetRestriction(pc,i,NULL);
1056:     }
1057:     for (i=0; i<n-1; i++) {
1058:       if (!mglevels[i]->b) {
1059:         Vec *vec;
1060:         KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);
1061:         PCMGSetRhs(pc,i,*vec);
1062:         VecDestroy(vec);
1063:         PetscFree(vec);
1064:       }
1065:       if (!mglevels[i]->r && i) {
1066:         VecDuplicate(mglevels[i]->b,&tvec);
1067:         PCMGSetR(pc,i,tvec);
1068:         VecDestroy(&tvec);
1069:       }
1070:       if (!mglevels[i]->x) {
1071:         VecDuplicate(mglevels[i]->b,&tvec);
1072:         PCMGSetX(pc,i,tvec);
1073:         VecDestroy(&tvec);
1074:       }
1075:       if (doCR) {
1076:         VecDuplicate(mglevels[i]->b,&mglevels[i]->crx);
1077:         VecDuplicate(mglevels[i]->b,&mglevels[i]->crb);
1078:         VecZeroEntries(mglevels[i]->crb);
1079:       }
1080:     }
1081:     if (n != 1 && !mglevels[n-1]->r) {
1082:       /* PCMGSetR() on the finest level if user did not supply it */
1083:       Vec *vec;
1084:       KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);
1085:       PCMGSetR(pc,n-1,*vec);
1086:       VecDestroy(vec);
1087:       PetscFree(vec);
1088:     }
1089:     if (doCR) {
1090:       VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx);
1091:       VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb);
1092:       VecZeroEntries(mglevels[n-1]->crb);
1093:     }
1094:   }

1096:   if (pc->dm) {
1097:     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
1098:     for (i=0; i<n-1; i++) {
1099:       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
1100:     }
1101:   }

1103:   for (i=1; i<n; i++) {
1104:     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
1105:       /* if doing only down then initial guess is zero */
1106:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
1107:     }
1108:     if (mglevels[i]->cr) {KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);}
1109:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1110:     KSPSetUp(mglevels[i]->smoothd);
1111:     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1112:       pc->failedreason = PC_SUBPC_ERROR;
1113:     }
1114:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1115:     if (!mglevels[i]->residual) {
1116:       Mat mat;
1117:       KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);
1118:       PCMGSetResidual(pc,i,PCMGResidualDefault,mat);
1119:     }
1120:     if (!mglevels[i]->residualtranspose) {
1121:       Mat mat;
1122:       KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);
1123:       PCMGSetResidualTranspose(pc,i,PCMGResidualTransposeDefault,mat);
1124:     }
1125:   }
1126:   for (i=1; i<n; i++) {
1127:     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
1128:       Mat downmat,downpmat;

1130:       /* check if operators have been set for up, if not use down operators to set them */
1131:       KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);
1132:       if (!opsset) {
1133:         KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
1134:         KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);
1135:       }

1137:       KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);
1138:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1139:       KSPSetUp(mglevels[i]->smoothu);
1140:       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
1141:         pc->failedreason = PC_SUBPC_ERROR;
1142:       }
1143:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1144:     }
1145:     if (mglevels[i]->cr) {
1146:       Mat downmat,downpmat;

1148:       /* check if operators have been set for up, if not use down operators to set them */
1149:       KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL);
1150:       if (!opsset) {
1151:         KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
1152:         KSPSetOperators(mglevels[i]->cr,downmat,downpmat);
1153:       }

1155:       KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);
1156:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1157:       KSPSetUp(mglevels[i]->cr);
1158:       if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) {
1159:         pc->failedreason = PC_SUBPC_ERROR;
1160:       }
1161:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1162:     }
1163:   }

1165:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
1166:   KSPSetUp(mglevels[0]->smoothd);
1167:   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1168:     pc->failedreason = PC_SUBPC_ERROR;
1169:   }
1170:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}

1172:   /*
1173:      Dump the interpolation/restriction matrices plus the
1174:    Jacobian/stiffness on each level. This allows MATLAB users to
1175:    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.

1177:    Only support one or the other at the same time.
1178:   */
1179: #if defined(PETSC_USE_SOCKET_VIEWER)
1180:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);
1181:   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
1182:   dump = PETSC_FALSE;
1183: #endif
1184:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);
1185:   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));

1187:   if (viewer) {
1188:     for (i=1; i<n; i++) {
1189:       MatView(mglevels[i]->restrct,viewer);
1190:     }
1191:     for (i=0; i<n; i++) {
1192:       KSPGetPC(mglevels[i]->smoothd,&pc);
1193:       MatView(pc->mat,viewer);
1194:     }
1195:   }
1196:   return(0);
1197: }

1199: /* -------------------------------------------------------------------------------------*/

1201: PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
1202: {
1203:   PC_MG *mg = (PC_MG *) pc->data;

1206:   *levels = mg->nlevels;
1207:   return(0);
1208: }

1210: /*@
1211:    PCMGGetLevels - Gets the number of levels to use with MG.

1213:    Not Collective

1215:    Input Parameter:
1216: .  pc - the preconditioner context

1218:    Output parameter:
1219: .  levels - the number of levels

1221:    Level: advanced

1223: .seealso: PCMGSetLevels()
1224: @*/
1225: PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
1226: {

1232:   *levels = 0;
1233:   PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));
1234:   return(0);
1235: }

1237: /*@
1238:    PCMGSetType - Determines the form of multigrid to use:
1239:    multiplicative, additive, full, or the Kaskade algorithm.

1241:    Logically Collective on PC

1243:    Input Parameters:
1244: +  pc - the preconditioner context
1245: -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
1246:    PC_MG_FULL, PC_MG_KASKADE

1248:    Options Database Key:
1249: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
1250:    additive, full, kaskade

1252:    Level: advanced

1254: .seealso: PCMGSetLevels()
1255: @*/
1256: PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
1257: {
1258:   PC_MG *mg = (PC_MG*)pc->data;

1263:   mg->am = form;
1264:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1265:   else pc->ops->applyrichardson = NULL;
1266:   return(0);
1267: }

1269: /*@
1270:    PCMGGetType - Determines the form of multigrid to use:
1271:    multiplicative, additive, full, or the Kaskade algorithm.

1273:    Logically Collective on PC

1275:    Input Parameter:
1276: .  pc - the preconditioner context

1278:    Output Parameter:
1279: .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE


1282:    Level: advanced

1284: .seealso: PCMGSetLevels()
1285: @*/
1286: PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1287: {
1288:   PC_MG *mg = (PC_MG*)pc->data;

1292:   *type = mg->am;
1293:   return(0);
1294: }

1296: /*@
1297:    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
1298:    complicated cycling.

1300:    Logically Collective on PC

1302:    Input Parameters:
1303: +  pc - the multigrid context
1304: -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W

1306:    Options Database Key:
1307: .  -pc_mg_cycle_type <v,w> - provide the cycle desired

1309:    Level: advanced

1311: .seealso: PCMGSetCycleTypeOnLevel()
1312: @*/
1313: PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
1314: {
1315:   PC_MG        *mg        = (PC_MG*)pc->data;
1316:   PC_MG_Levels **mglevels = mg->levels;
1317:   PetscInt     i,levels;

1322:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1323:   levels = mglevels[0]->levels;
1324:   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
1325:   return(0);
1326: }

1328: /*@
1329:    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1330:          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used

1332:    Logically Collective on PC

1334:    Input Parameters:
1335: +  pc - the multigrid context
1336: -  n - number of cycles (default is 1)

1338:    Options Database Key:
1339: .  -pc_mg_multiplicative_cycles n

1341:    Level: advanced

1343:    Notes:
1344:     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()

1346: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1347: @*/
1348: PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1349: {
1350:   PC_MG        *mg        = (PC_MG*)pc->data;

1355:   mg->cyclesperpcapply = n;
1356:   return(0);
1357: }

1359: PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1360: {
1361:   PC_MG *mg = (PC_MG*)pc->data;

1364:   mg->galerkin = use;
1365:   return(0);
1366: }

1368: /*@
1369:    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
1370:       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i

1372:    Logically Collective on PC

1374:    Input Parameters:
1375: +  pc - the multigrid context
1376: -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE

1378:    Options Database Key:
1379: .  -pc_mg_galerkin <both,pmat,mat,none>

1381:    Level: intermediate

1383:    Notes:
1384:     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
1385:      use the PCMG construction of the coarser grids.

1387: .seealso: PCMGGetGalerkin(), PCMGGalerkinType

1389: @*/
1390: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1391: {

1396:   PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));
1397:   return(0);
1398: }

1400: /*@
1401:    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1402:       A_i-1 = r_i * A_i * p_i

1404:    Not Collective

1406:    Input Parameter:
1407: .  pc - the multigrid context

1409:    Output Parameter:
1410: .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL

1412:    Level: intermediate

1414: .seealso: PCMGSetGalerkin(), PCMGGalerkinType

1416: @*/
1417: PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1418: {
1419:   PC_MG *mg = (PC_MG*)pc->data;

1423:   *galerkin = mg->galerkin;
1424:   return(0);
1425: }

1427: PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1428: {
1429:   PC_MG *mg = (PC_MG *) pc->data;

1432:   mg->adaptInterpolation = adapt;
1433:   return(0);
1434: }

1436: PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1437: {
1438:   PC_MG *mg = (PC_MG *) pc->data;

1441:   *adapt = mg->adaptInterpolation;
1442:   return(0);
1443: }

1445: PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr)
1446: {
1447:   PC_MG *mg = (PC_MG *) pc->data;

1450:   mg->compatibleRelaxation = cr;
1451:   return(0);
1452: }

1454: PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr)
1455: {
1456:   PC_MG *mg = (PC_MG *) pc->data;

1459:   *cr = mg->compatibleRelaxation;
1460:   return(0);
1461: }

1463: /*@
1464:   PCMGSetAdaptInterpolation - Adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.

1466:   Logically Collective on PC

1468:   Input Parameters:
1469: + pc    - the multigrid context
1470: - adapt - flag for adaptation of the interpolator

1472:   Options Database Keys:
1473: + -pc_mg_adapt_interp                     - Turn on adaptation
1474: . -pc_mg_adapt_interp_n <int>             - The number of modes to use, should be divisible by dimension
1475: - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector

1477:   Level: intermediate

1479: .keywords: MG, set, Galerkin
1480: .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1481: @*/
1482: PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1483: {

1488:   PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));
1489:   return(0);
1490: }

1492: /*@
1493:   PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.

1495:   Not collective

1497:   Input Parameter:
1498: . pc    - the multigrid context

1500:   Output Parameter:
1501: . adapt - flag for adaptation of the interpolator

1503:   Level: intermediate

1505: .keywords: MG, set, Galerkin
1506: .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1507: @*/
1508: PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1509: {

1515:   PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));
1516:   return(0);
1517: }

1519: /*@
1520:   PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation.

1522:   Logically Collective on PC

1524:   Input Parameters:
1525: + pc - the multigrid context
1526: - cr - flag for compatible relaxation

1528:   Options Database Keys:
1529: . -pc_mg_adapt_cr - Turn on compatible relaxation

1531:   Level: intermediate

1533: .keywords: MG, set, Galerkin
1534: .seealso: PCMGGetAdaptCR(), PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1535: @*/
1536: PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr)
1537: {

1542:   PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr));
1543:   return(0);
1544: }

1546: /*@
1547:   PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation.

1549:   Not collective

1551:   Input Parameter:
1552: . pc    - the multigrid context

1554:   Output Parameter:
1555: . cr - flag for compatible relaxaion

1557:   Level: intermediate

1559: .keywords: MG, set, Galerkin
1560: .seealso: PCMGSetAdaptCR(), PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1561: @*/
1562: PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr)
1563: {

1569:   PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr));
1570:   return(0);
1571: }

1573: /*@
1574:    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1575:    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1576:    pre- and post-smoothing steps.

1578:    Logically Collective on PC

1580:    Input Parameters:
1581: +  mg - the multigrid context
1582: -  n - the number of smoothing steps

1584:    Options Database Key:
1585: .  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps

1587:    Level: advanced

1589:    Notes:
1590:     this does not set a value on the coarsest grid, since we assume that
1591:     there is no separate smooth up on the coarsest grid.

1593: .seealso: PCMGSetDistinctSmoothUp()
1594: @*/
1595: PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1596: {
1597:   PC_MG          *mg        = (PC_MG*)pc->data;
1598:   PC_MG_Levels   **mglevels = mg->levels;
1600:   PetscInt       i,levels;

1605:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1606:   levels = mglevels[0]->levels;

1608:   for (i=1; i<levels; i++) {
1609:     KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1610:     KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1611:     mg->default_smoothu = n;
1612:     mg->default_smoothd = n;
1613:   }
1614:   return(0);
1615: }

1617: /*@
1618:    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels
1619:        and adds the suffix _up to the options name

1621:    Logically Collective on PC

1623:    Input Parameters:
1624: .  pc - the preconditioner context

1626:    Options Database Key:
1627: .  -pc_mg_distinct_smoothup

1629:    Level: advanced

1631:    Notes:
1632:     this does not set a value on the coarsest grid, since we assume that
1633:     there is no separate smooth up on the coarsest grid.

1635: .seealso: PCMGSetNumberSmooth()
1636: @*/
1637: PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1638: {
1639:   PC_MG          *mg        = (PC_MG*)pc->data;
1640:   PC_MG_Levels   **mglevels = mg->levels;
1642:   PetscInt       i,levels;
1643:   KSP            subksp;

1647:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1648:   levels = mglevels[0]->levels;

1650:   for (i=1; i<levels; i++) {
1651:     const char *prefix = NULL;
1652:     /* make sure smoother up and down are different */
1653:     PCMGGetSmootherUp(pc,i,&subksp);
1654:     KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);
1655:     KSPSetOptionsPrefix(subksp,prefix);
1656:     KSPAppendOptionsPrefix(subksp,"up_");
1657:   }
1658:   return(0);
1659: }

1661: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1662: PetscErrorCode  PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[])
1663: {
1664:   PC_MG          *mg        = (PC_MG*)pc->data;
1665:   PC_MG_Levels   **mglevels = mg->levels;
1666:   Mat            *mat;
1667:   PetscInt       l;

1671:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1672:   PetscMalloc1(mg->nlevels,&mat);
1673:   for (l=1; l< mg->nlevels; l++) {
1674:     mat[l-1] = mglevels[l]->interpolate;
1675:     PetscObjectReference((PetscObject)mat[l-1]);
1676:   }
1677:   *num_levels = mg->nlevels;
1678:   *interpolations = mat;
1679:   return(0);
1680: }

1682: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1683: PetscErrorCode  PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
1684: {
1685:   PC_MG          *mg        = (PC_MG*)pc->data;
1686:   PC_MG_Levels   **mglevels = mg->levels;
1687:   PetscInt       l;
1688:   Mat            *mat;

1692:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1693:   PetscMalloc1(mg->nlevels,&mat);
1694:   for (l=0; l<mg->nlevels-1; l++) {
1695:     KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));
1696:     PetscObjectReference((PetscObject)mat[l]);
1697:   }
1698:   *num_levels = mg->nlevels;
1699:   *coarseOperators = mat;
1700:   return(0);
1701: }

1703: /*@C
1704:   PCMGRegisterCoarseSpaceConstructor -  Adds a method to the PCMG package for coarse space construction.

1706:   Not collective

1708:   Input Parameters:
1709: + name     - name of the constructor
1710: - function - constructor routine

1712:   Notes:
1713:   Calling sequence for the routine:
1714: $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1715: $   pc        - The PC object
1716: $   l         - The multigrid level, 0 is the coarse level
1717: $   dm        - The DM for this level
1718: $   smooth    - The level smoother
1719: $   Nc        - The size of the coarse space
1720: $   initGuess - Basis for an initial guess for the space
1721: $   coarseSp  - A basis for the computed coarse space

1723:   Level: advanced

1725: .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister()
1726: @*/
1727: PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1728: {

1732:   PCInitializePackage();
1733:   PetscFunctionListAdd(&PCMGCoarseList,name,function);
1734:   return(0);
1735: }

1737: /*@C
1738:   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.

1740:   Not collective

1742:   Input Parameter:
1743: . name     - name of the constructor

1745:   Output Parameter:
1746: . function - constructor routine

1748:   Notes:
1749:   Calling sequence for the routine:
1750: $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1751: $   pc        - The PC object
1752: $   l         - The multigrid level, 0 is the coarse level
1753: $   dm        - The DM for this level
1754: $   smooth    - The level smoother
1755: $   Nc        - The size of the coarse space
1756: $   initGuess - Basis for an initial guess for the space
1757: $   coarseSp  - A basis for the computed coarse space

1759:   Level: advanced

1761: .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister()
1762: @*/
1763: PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1764: {

1768:   PetscFunctionListFind(PCMGCoarseList,name,function);
1769:   return(0);
1770: }

1772: /* ----------------------------------------------------------------------------------------*/

1774: /*MC
1775:    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1776:     information about the coarser grid matrices and restriction/interpolation operators.

1778:    Options Database Keys:
1779: +  -pc_mg_levels <nlevels> - number of levels including finest
1780: .  -pc_mg_cycle_type <v,w> - provide the cycle desired
1781: .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1782: .  -pc_mg_log - log information about time spent on each level of the solver
1783: .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
1784: .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1785: .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1786: .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1787:                         to the Socket viewer for reading from MATLAB.
1788: -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1789:                         to the binary output file called binaryoutput

1791:    Notes:
1792:     If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method

1794:        When run with a single level the smoother options are used on that level NOT the coarse grid solver options

1796:        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
1797:        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
1798:        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
1799:        residual is computed at the end of each cycle.

1801:    Level: intermediate

1803: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1804:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1805:            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1806:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1807:            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1808: M*/

1810: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1811: {
1812:   PC_MG          *mg;

1816:   PetscNewLog(pc,&mg);
1817:   pc->data     = (void*)mg;
1818:   mg->nlevels  = -1;
1819:   mg->am       = PC_MG_MULTIPLICATIVE;
1820:   mg->galerkin = PC_MG_GALERKIN_NONE;
1821:   mg->adaptInterpolation = PETSC_FALSE;
1822:   mg->Nc                 = -1;
1823:   mg->eigenvalue         = -1;

1825:   pc->useAmat = PETSC_TRUE;

1827:   pc->ops->apply          = PCApply_MG;
1828:   pc->ops->applytranspose = PCApplyTranspose_MG;
1829:   pc->ops->setup          = PCSetUp_MG;
1830:   pc->ops->reset          = PCReset_MG;
1831:   pc->ops->destroy        = PCDestroy_MG;
1832:   pc->ops->setfromoptions = PCSetFromOptions_MG;
1833:   pc->ops->view           = PCView_MG;

1835:   PetscObjectComposedDataRegister(&mg->eigenvalue);
1836:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1837:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);
1838:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);
1839:   PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);
1840:   PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);
1841:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);
1842:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);
1843:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG);
1844:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG);
1845:   return(0);
1846: }