Actual source code: mg.c


  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,PetscBool matapp,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:     if (matapp) {
 26:       KSPMatSolve(mglevels->smoothd,mglevels->B,mglevels->X);  /* pre-smooth */
 27:       KSPCheckSolve(mglevels->smoothd,pc,NULL);
 28:     } else {
 29:       KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);  /* pre-smooth */
 30:       KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);
 31:     }
 32:   } else {
 33:     if (matapp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported");
 34:     KSPSolveTranspose(mglevels->smoothu,mglevels->b,mglevels->x); /* transpose of post-smooth */
 35:     KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);
 36:   }
 37:   if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
 38:   if (mglevels->level) {  /* not the coarsest grid */
 39:     if (mglevels->eventresidual) {PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);}
 40:     if (matapp && !mglevels->R) {
 41:       MatDuplicate(mglevels->B,MAT_DO_NOT_COPY_VALUES,&mglevels->R);
 42:     }
 43:     if (!transpose) {
 44:       if (matapp) { (*mglevels->matresidual)(mglevels->A,mglevels->B,mglevels->X,mglevels->R); }
 45:       else { (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r); }
 46:     } else {
 47:       if (matapp) { (*mglevels->matresidualtranspose)(mglevels->A,mglevels->B,mglevels->X,mglevels->R); }
 48:       else { (*mglevels->residualtranspose)(mglevels->A,mglevels->b,mglevels->x,mglevels->r); }
 49:     }
 50:     if (mglevels->eventresidual) {PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);}

 52:     /* if on finest level and have convergence criteria set */
 53:     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
 54:       PetscReal rnorm;
 55:       VecNorm(mglevels->r,NORM_2,&rnorm);
 56:       if (rnorm <= mg->ttol) {
 57:         if (rnorm < mg->abstol) {
 58:           *reason = PCRICHARDSON_CONVERGED_ATOL;
 59:           PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);
 60:         } else {
 61:           *reason = PCRICHARDSON_CONVERGED_RTOL;
 62:           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);
 63:         }
 64:         return(0);
 65:       }
 66:     }

 68:     mgc = *(mglevelsin - 1);
 69:     if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
 70:     if (!transpose) {
 71:       if (matapp) { MatMatRestrict(mglevels->restrct,mglevels->R,&mgc->B); }
 72:       else { MatRestrict(mglevels->restrct,mglevels->r,mgc->b); }
 73:     } else {
 74:       if (matapp) { MatMatRestrict(mglevels->interpolate,mglevels->R,&mgc->B); }
 75:       else { MatRestrict(mglevels->interpolate,mglevels->r,mgc->b); }
 76:     }
 77:     if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
 78:     if (matapp) {
 79:       if (!mgc->X) {
 80:         MatDuplicate(mgc->B,MAT_DO_NOT_COPY_VALUES,&mgc->X);
 81:       } else {
 82:         MatZeroEntries(mgc->X);
 83:       }
 84:     } else {
 85:       VecZeroEntries(mgc->x);
 86:     }
 87:     while (cycles--) {
 88:       PCMGMCycle_Private(pc,mglevelsin-1,transpose,matapp,reason);
 89:     }
 90:     if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
 91:     if (!transpose) {
 92:       if (matapp) { MatMatInterpolateAdd(mglevels->interpolate,mgc->X,mglevels->X,&mglevels->X); }
 93:       else { MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x); }
 94:     } else {
 95:       MatInterpolateAdd(mglevels->restrct,mgc->x,mglevels->x,mglevels->x);
 96:     }
 97:     if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
 98:     if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
 99:     if (!transpose) {
100:       if (matapp) {
101:         KSPMatSolve(mglevels->smoothu,mglevels->B,mglevels->X);    /* post smooth */
102:         KSPCheckSolve(mglevels->smoothu,pc,NULL);
103:       } else {
104:         KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);    /* post smooth */
105:         KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);
106:       }
107:     } else {
108:       if (matapp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported");
109:       KSPSolveTranspose(mglevels->smoothd,mglevels->b,mglevels->x);    /* post smooth */
110:       KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);
111:     }
112:     if (mglevels->cr) {
113:       if (matapp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported");
114:       /* TODO Turn on copy and turn off noisy if we have an exact solution
115:       VecCopy(mglevels->x, mglevels->crx);
116:       VecCopy(mglevels->b, mglevels->crb); */
117:       KSPSetNoisy_Private(mglevels->crx);
118:       KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx);    /* compatible relaxation */
119:       KSPCheckSolve(mglevels->cr,pc,mglevels->crx);
120:     }
121:     if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
122:   }
123:   return(0);
124: }

126: 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)
127: {
128:   PC_MG          *mg        = (PC_MG*)pc->data;
129:   PC_MG_Levels   **mglevels = mg->levels;
131:   PC             tpc;
132:   PetscBool      changeu,changed;
133:   PetscInt       levels = mglevels[0]->levels,i;

136:   /* When the DM is supplying the matrix then it will not exist until here */
137:   for (i=0; i<levels; i++) {
138:     if (!mglevels[i]->A) {
139:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
140:       PetscObjectReference((PetscObject)mglevels[i]->A);
141:     }
142:   }

144:   KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
145:   PCPreSolveChangeRHS(tpc,&changed);
146:   KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
147:   PCPreSolveChangeRHS(tpc,&changeu);
148:   if (!changed && !changeu) {
149:     VecDestroy(&mglevels[levels-1]->b);
150:     mglevels[levels-1]->b = b;
151:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
152:     if (!mglevels[levels-1]->b) {
153:       Vec *vec;

155:       KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
156:       mglevels[levels-1]->b = *vec;
157:       PetscFree(vec);
158:     }
159:     VecCopy(b,mglevels[levels-1]->b);
160:   }
161:   mglevels[levels-1]->x = x;

163:   mg->rtol   = rtol;
164:   mg->abstol = abstol;
165:   mg->dtol   = dtol;
166:   if (rtol) {
167:     /* compute initial residual norm for relative convergence test */
168:     PetscReal rnorm;
169:     if (zeroguess) {
170:       VecNorm(b,NORM_2,&rnorm);
171:     } else {
172:       (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);
173:       VecNorm(w,NORM_2,&rnorm);
174:     }
175:     mg->ttol = PetscMax(rtol*rnorm,abstol);
176:   } else if (abstol) mg->ttol = abstol;
177:   else mg->ttol = 0.0;

179:   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
180:      stop prematurely due to small residual */
181:   for (i=1; i<levels; i++) {
182:     KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
183:     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
184:       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
185:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
186:       KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
187:     }
188:   }

190:   *reason = (PCRichardsonConvergedReason)0;
191:   for (i=0; i<its; i++) {
192:     PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_FALSE,PETSC_FALSE,reason);
193:     if (*reason) break;
194:   }
195:   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
196:   *outits = i;
197:   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
198:   return(0);
199: }

201: PetscErrorCode PCReset_MG(PC pc)
202: {
203:   PC_MG          *mg        = (PC_MG*)pc->data;
204:   PC_MG_Levels   **mglevels = mg->levels;
206:   PetscInt       i,c,n;

209:   if (mglevels) {
210:     n = mglevels[0]->levels;
211:     for (i=0; i<n-1; i++) {
212:       VecDestroy(&mglevels[i+1]->r);
213:       VecDestroy(&mglevels[i]->b);
214:       VecDestroy(&mglevels[i]->x);
215:       MatDestroy(&mglevels[i+1]->R);
216:       MatDestroy(&mglevels[i]->B);
217:       MatDestroy(&mglevels[i]->X);
218:       VecDestroy(&mglevels[i]->crx);
219:       VecDestroy(&mglevels[i]->crb);
220:       MatDestroy(&mglevels[i+1]->restrct);
221:       MatDestroy(&mglevels[i+1]->interpolate);
222:       MatDestroy(&mglevels[i+1]->inject);
223:       VecDestroy(&mglevels[i+1]->rscale);
224:     }
225:     VecDestroy(&mglevels[n-1]->crx);
226:     VecDestroy(&mglevels[n-1]->crb);
227:     /* this is not null only if the smoother on the finest level
228:        changes the rhs during PreSolve */
229:     VecDestroy(&mglevels[n-1]->b);
230:     MatDestroy(&mglevels[n-1]->B);

232:     for (i=0; i<n; i++) {
233:       if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {VecDestroy(&mglevels[i]->coarseSpace[c]);}
234:       PetscFree(mglevels[i]->coarseSpace);
235:       mglevels[i]->coarseSpace = NULL;
236:       MatDestroy(&mglevels[i]->A);
237:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
238:         KSPReset(mglevels[i]->smoothd);
239:       }
240:       KSPReset(mglevels[i]->smoothu);
241:       if (mglevels[i]->cr) {KSPReset(mglevels[i]->cr);}
242:     }
243:     mg->Nc = 0;
244:   }
245:   return(0);
246: }

248: /* Implementing CR

250: 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

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

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

256:   Inj^T Inj

258: and

260:   S = I - Inj^T Inj

262: since

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

266: Brannick suggests

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

270: 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

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

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

276:   Check: || Inj P - I ||_F < tol
277:   Check: In general, Inj Inj^T = I
278: */

280: typedef struct {
281:   PC       mg;  /* The PCMG object */
282:   PetscInt l;   /* The multigrid level for this solver */
283:   Mat      Inj; /* The injection matrix */
284:   Mat      S;   /* I - Inj^T Inj */
285: } CRContext;

287: static PetscErrorCode CRSetup_Private(PC pc)
288: {
289:   CRContext     *ctx;
290:   Mat            It;

294:   PCShellGetContext(pc, &ctx);
295:   PCMGGetInjection(ctx->mg, ctx->l, &It);
296:   if (!It) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG");
297:   MatCreateTranspose(It, &ctx->Inj);
298:   MatCreateNormal(ctx->Inj, &ctx->S);
299:   MatScale(ctx->S, -1.0);
300:   MatShift(ctx->S,  1.0);
301:   return(0);
302: }

304: static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y)
305: {
306:   CRContext     *ctx;

310:   PCShellGetContext(pc, &ctx);
311:   MatMult(ctx->S, x, y);
312:   return(0);
313: }

315: static PetscErrorCode CRDestroy_Private(PC pc)
316: {
317:   CRContext     *ctx;

321:   PCShellGetContext(pc, &ctx);
322:   MatDestroy(&ctx->Inj);
323:   MatDestroy(&ctx->S);
324:   PetscFree(ctx);
325:   PCShellSetContext(pc, NULL);
326:   return(0);
327: }

329: static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr)
330: {
331:   CRContext     *ctx;

335:   PCCreate(PetscObjectComm((PetscObject) pc), cr);
336:   PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)");
337:   PetscCalloc1(1, &ctx);
338:   ctx->mg = pc;
339:   ctx->l  = l;
340:   PCSetType(*cr, PCSHELL);
341:   PCShellSetContext(*cr, ctx);
342:   PCShellSetApply(*cr, CRApply_Private);
343:   PCShellSetSetUp(*cr, CRSetup_Private);
344:   PCShellSetDestroy(*cr, CRDestroy_Private);
345:   return(0);
346: }

348: PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
349: {
351:   PC_MG          *mg        = (PC_MG*)pc->data;
352:   MPI_Comm       comm;
353:   PC_MG_Levels   **mglevels = mg->levels;
354:   PCMGType       mgtype     = mg->am;
355:   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
356:   PetscInt       i;
357:   PetscMPIInt    size;
358:   const char     *prefix;
359:   PC             ipc;
360:   PetscInt       n;

365:   if (mg->nlevels == levels) return(0);
366:   PetscObjectGetComm((PetscObject)pc,&comm);
367:   if (mglevels) {
368:     mgctype = mglevels[0]->cycles;
369:     /* changing the number of levels so free up the previous stuff */
370:     PCReset_MG(pc);
371:     n    = mglevels[0]->levels;
372:     for (i=0; i<n; i++) {
373:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
374:         KSPDestroy(&mglevels[i]->smoothd);
375:       }
376:       KSPDestroy(&mglevels[i]->smoothu);
377:       KSPDestroy(&mglevels[i]->cr);
378:       PetscFree(mglevels[i]);
379:     }
380:     PetscFree(mg->levels);
381:   }

383:   mg->nlevels = levels;

385:   PetscMalloc1(levels,&mglevels);
386:   PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));

388:   PCGetOptionsPrefix(pc,&prefix);

390:   mg->stageApply = 0;
391:   for (i=0; i<levels; i++) {
392:     PetscNewLog(pc,&mglevels[i]);

394:     mglevels[i]->level               = i;
395:     mglevels[i]->levels              = levels;
396:     mglevels[i]->cycles              = mgctype;
397:     mg->default_smoothu              = 2;
398:     mg->default_smoothd              = 2;
399:     mglevels[i]->eventsmoothsetup    = 0;
400:     mglevels[i]->eventsmoothsolve    = 0;
401:     mglevels[i]->eventresidual       = 0;
402:     mglevels[i]->eventinterprestrict = 0;

404:     if (comms) comm = comms[i];
405:     if (comm != MPI_COMM_NULL) {
406:       KSPCreate(comm,&mglevels[i]->smoothd);
407:       KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);
408:       PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);
409:       KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);
410:       PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);
411:       if (i || levels == 1) {
412:         char tprefix[128];

414:         KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);
415:         KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);
416:         KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);
417:         KSPGetPC(mglevels[i]->smoothd,&ipc);
418:         PCSetType(ipc,PCSOR);
419:         KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);

421:         PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);
422:         KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);
423:       } else {
424:         KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");

426:         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
427:         KSPSetType(mglevels[0]->smoothd,KSPPREONLY);
428:         KSPGetPC(mglevels[0]->smoothd,&ipc);
429:         MPI_Comm_size(comm,&size);
430:         if (size > 1) {
431:           PCSetType(ipc,PCREDUNDANT);
432:         } else {
433:           PCSetType(ipc,PCLU);
434:         }
435:         PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);
436:       }
437:       PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);
438:     }
439:     mglevels[i]->smoothu = mglevels[i]->smoothd;
440:     mg->rtol             = 0.0;
441:     mg->abstol           = 0.0;
442:     mg->dtol             = 0.0;
443:     mg->ttol             = 0.0;
444:     mg->cyclesperpcapply = 1;
445:   }
446:   mg->levels = mglevels;
447:   PCMGSetType(pc,mgtype);
448:   return(0);
449: }

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

455:    Logically Collective on PC

457:    Input Parameters:
458: +  pc - the preconditioner context
459: .  levels - the number of levels
460: -  comms - optional communicators for each level; this is to allow solving the coarser problems
461:            on smaller sets of processes. For processes that are not included in the computation
462:            you must pass MPI_COMM_NULL. Use comms = NULL to specify that all processes
463:            should participate in each level of problem.

465:    Level: intermediate

467:    Notes:
468:      If the number of levels is one then the multigrid uses the -mg_levels prefix
469:      for setting the level options rather than the -mg_coarse prefix.

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

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

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

488:    Fortran Notes:
489:      Use comms = PETSC_NULL_MPI_COMM as the equivalent of NULL in the C interface. Note PETSC_NULL_MPI_COMM
490:      is not MPI_COMM_NULL. It is more like PETSC_NULL_INTEGER, PETSC_NULL_REAL etc.

492: .seealso: PCMGSetType(), PCMGGetLevels()
493: @*/
494: PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
495: {

501:   PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));
502:   return(0);
503: }

505: PetscErrorCode PCDestroy_MG(PC pc)
506: {
508:   PC_MG          *mg        = (PC_MG*)pc->data;
509:   PC_MG_Levels   **mglevels = mg->levels;
510:   PetscInt       i,n;

513:   PCReset_MG(pc);
514:   if (mglevels) {
515:     n = mglevels[0]->levels;
516:     for (i=0; i<n; i++) {
517:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
518:         KSPDestroy(&mglevels[i]->smoothd);
519:       }
520:       KSPDestroy(&mglevels[i]->smoothu);
521:       KSPDestroy(&mglevels[i]->cr);
522:       PetscFree(mglevels[i]);
523:     }
524:     PetscFree(mg->levels);
525:   }
526:   PetscFree(pc->data);
527:   PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);
528:   PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);
529:   return(0);
530: }

532: /*
533:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
534:              or full cycle of multigrid.

536:   Note:
537:   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
538: */
539: static PetscErrorCode PCApply_MG_Internal(PC pc,Vec b,Vec x,Mat B,Mat X,PetscBool transpose)
540: {
541:   PC_MG          *mg        = (PC_MG*)pc->data;
542:   PC_MG_Levels   **mglevels = mg->levels;
544:   PC             tpc;
545:   PetscInt       levels = mglevels[0]->levels,i;
546:   PetscBool      changeu,changed,matapp;

549:   matapp = (PetscBool)(B && X);
550:   if (mg->stageApply) {PetscLogStagePush(mg->stageApply);}
551:   /* When the DM is supplying the matrix then it will not exist until here */
552:   for (i=0; i<levels; i++) {
553:     if (!mglevels[i]->A) {
554:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
555:       PetscObjectReference((PetscObject)mglevels[i]->A);
556:     }
557:   }

559:   KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
560:   PCPreSolveChangeRHS(tpc,&changed);
561:   KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
562:   PCPreSolveChangeRHS(tpc,&changeu);
563:   if (!changeu && !changed) {
564:     if (matapp) {
565:       MatDestroy(&mglevels[levels-1]->B);
566:       mglevels[levels-1]->B = B;
567:     } else {
568:       VecDestroy(&mglevels[levels-1]->b);
569:       mglevels[levels-1]->b = b;
570:     }
571:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
572:     if (matapp) {
573:       if (mglevels[levels-1]->B) {
574:         PetscInt  N1,N2;
575:         PetscBool flg;

577:         MatGetSize(mglevels[levels-1]->B,NULL,&N1);
578:         MatGetSize(B,NULL,&N2);
579:         PetscObjectTypeCompare((PetscObject)mglevels[levels-1]->B,((PetscObject)B)->type_name,&flg);
580:         if (N1 != N2 || !flg) {
581:           MatDestroy(&mglevels[levels-1]->B);
582:         }
583:       }
584:       if (!mglevels[levels-1]->B) {
585:         MatDuplicate(B,MAT_COPY_VALUES,&mglevels[levels-1]->B);
586:       } else {
587:         MatCopy(B,mglevels[levels-1]->B,SAME_NONZERO_PATTERN);
588:       }
589:     } else {
590:       if (!mglevels[levels-1]->b) {
591:         Vec *vec;

593:         KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
594:         mglevels[levels-1]->b = *vec;
595:         PetscFree(vec);
596:       }
597:       VecCopy(b,mglevels[levels-1]->b);
598:     }
599:   }
600:   if (matapp) { mglevels[levels-1]->X = X; }
601:   else { mglevels[levels-1]->x = x; }

603:   /* If coarser Xs are present, it means we have already block applied the PC at least once
604:      Reset operators if sizes/type do no match */
605:   if (matapp && levels > 1 && mglevels[levels-2]->X) {
606:     PetscInt  Xc,Bc;
607:     PetscBool flg;

609:     MatGetSize(mglevels[levels-2]->X,NULL,&Xc);
610:     MatGetSize(mglevels[levels-1]->B,NULL,&Bc);
611:     PetscObjectTypeCompare((PetscObject)mglevels[levels-2]->X,((PetscObject)mglevels[levels-1]->X)->type_name,&flg);
612:     if (Xc != Bc || !flg) {
613:       MatDestroy(&mglevels[levels-1]->R);
614:       for (i=0;i<levels-1;i++) {
615:         MatDestroy(&mglevels[i]->R);
616:         MatDestroy(&mglevels[i]->B);
617:         MatDestroy(&mglevels[i]->X);
618:       }
619:     }
620:   }

622:   if (mg->am == PC_MG_MULTIPLICATIVE) {
623:     if (matapp) { MatZeroEntries(X); }
624:     else { VecZeroEntries(x); }
625:     for (i=0; i<mg->cyclesperpcapply; i++) {
626:       PCMGMCycle_Private(pc,mglevels+levels-1,transpose,matapp,NULL);
627:     }
628:   } else if (mg->am == PC_MG_ADDITIVE) {
629:     PCMGACycle_Private(pc,mglevels,transpose,matapp);
630:   } else if (mg->am == PC_MG_KASKADE) {
631:     PCMGKCycle_Private(pc,mglevels,transpose,matapp);
632:   } else {
633:     PCMGFCycle_Private(pc,mglevels,transpose,matapp);
634:   }
635:   if (mg->stageApply) {PetscLogStagePop();}
636:   if (!changeu && !changed) {
637:     if (matapp) { mglevels[levels-1]->B = NULL; }
638:     else { mglevels[levels-1]->b = NULL; }
639:   }
640:   return(0);
641: }

643: static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
644: {

648:   PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_FALSE);
649:   return(0);
650: }

652: static PetscErrorCode PCApplyTranspose_MG(PC pc,Vec b,Vec x)
653: {

657:   PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_TRUE);
658:   return(0);
659: }

661: static PetscErrorCode PCMatApply_MG(PC pc,Mat b,Mat x)
662: {

666:   PCApply_MG_Internal(pc,NULL,NULL,b,x,PETSC_FALSE);
667:   return(0);
668: }

670: PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
671: {
672:   PetscErrorCode   ierr;
673:   PetscInt         levels,cycles;
674:   PetscBool        flg, flg2;
675:   PC_MG            *mg = (PC_MG*)pc->data;
676:   PC_MG_Levels     **mglevels;
677:   PCMGType         mgtype;
678:   PCMGCycleType    mgctype;
679:   PCMGGalerkinType gtype;

682:   levels = PetscMax(mg->nlevels,1);
683:   PetscOptionsHead(PetscOptionsObject,"Multigrid options");
684:   PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
685:   if (!flg && !mg->levels && pc->dm) {
686:     DMGetRefineLevel(pc->dm,&levels);
687:     levels++;
688:     mg->usedmfornumberoflevels = PETSC_TRUE;
689:   }
690:   PCMGSetLevels(pc,levels,NULL);
691:   mglevels = mg->levels;

693:   mgctype = (PCMGCycleType) mglevels[0]->cycles;
694:   PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);
695:   if (flg) {
696:     PCMGSetCycleType(pc,mgctype);
697:   }
698:   gtype = mg->galerkin;
699:   PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);
700:   if (flg) {
701:     PCMGSetGalerkin(pc,gtype);
702:   }
703:   flg2 = PETSC_FALSE;
704:   PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);
705:   if (flg) {PCMGSetAdaptInterpolation(pc, flg2);}
706:   PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg);
707:   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);
708:   PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);
709:   flg2 = PETSC_FALSE;
710:   PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg);
711:   if (flg) {PCMGSetAdaptCR(pc, flg2);}
712:   flg = PETSC_FALSE;
713:   PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);
714:   if (flg) {
715:     PCMGSetDistinctSmoothUp(pc);
716:   }
717:   mgtype = mg->am;
718:   PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
719:   if (flg) {
720:     PCMGSetType(pc,mgtype);
721:   }
722:   if (mg->am == PC_MG_MULTIPLICATIVE) {
723:     PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);
724:     if (flg) {
725:       PCMGMultiplicativeSetCycles(pc,cycles);
726:     }
727:   }
728:   flg  = PETSC_FALSE;
729:   PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);
730:   if (flg) {
731:     PetscInt i;
732:     char     eventname[128];

734:     levels = mglevels[0]->levels;
735:     for (i=0; i<levels; i++) {
736:       sprintf(eventname,"MGSetup Level %d",(int)i);
737:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);
738:       sprintf(eventname,"MGSmooth Level %d",(int)i);
739:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);
740:       if (i) {
741:         sprintf(eventname,"MGResid Level %d",(int)i);
742:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);
743:         sprintf(eventname,"MGInterp Level %d",(int)i);
744:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);
745:       }
746:     }

748: #if defined(PETSC_USE_LOG)
749:     {
750:       const char    *sname = "MG Apply";
751:       PetscStageLog stageLog;
752:       PetscInt      st;

754:       PetscLogGetStageLog(&stageLog);
755:       for (st = 0; st < stageLog->numStages; ++st) {
756:         PetscBool same;

758:         PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
759:         if (same) mg->stageApply = st;
760:       }
761:       if (!mg->stageApply) {
762:         PetscLogStageRegister(sname, &mg->stageApply);
763:       }
764:     }
765: #endif
766:   }
767:   PetscOptionsTail();
768:   /* Check option consistency */
769:   PCMGGetGalerkin(pc, &gtype);
770:   PCMGGetAdaptInterpolation(pc, &flg);
771:   if (flg && (gtype >= PC_MG_GALERKIN_NONE)) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
772:   return(0);
773: }

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

780: #include <petscdraw.h>
781: PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
782: {
783:   PC_MG          *mg        = (PC_MG*)pc->data;
784:   PC_MG_Levels   **mglevels = mg->levels;
786:   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
787:   PetscBool      iascii,isbinary,isdraw;

790:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
791:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
792:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
793:   if (iascii) {
794:     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
795:     PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);
796:     if (mg->am == PC_MG_MULTIPLICATIVE) {
797:       PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);
798:     }
799:     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
800:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");
801:     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
802:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");
803:     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
804:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");
805:     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
806:       PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");
807:     } else {
808:       PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");
809:     }
810:     if (mg->view) {
811:       (*mg->view)(pc,viewer);
812:     }
813:     for (i=0; i<levels; i++) {
814:       if (!i) {
815:         PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);
816:       } else {
817:         PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);
818:       }
819:       PetscViewerASCIIPushTab(viewer);
820:       KSPView(mglevels[i]->smoothd,viewer);
821:       PetscViewerASCIIPopTab(viewer);
822:       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
823:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
824:       } else if (i) {
825:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);
826:         PetscViewerASCIIPushTab(viewer);
827:         KSPView(mglevels[i]->smoothu,viewer);
828:         PetscViewerASCIIPopTab(viewer);
829:       }
830:       if (i && mglevels[i]->cr) {
831:         PetscViewerASCIIPrintf(viewer,"CR solver on level %D -------------------------------\n",i);
832:         PetscViewerASCIIPushTab(viewer);
833:         KSPView(mglevels[i]->cr,viewer);
834:         PetscViewerASCIIPopTab(viewer);
835:       }
836:     }
837:   } else if (isbinary) {
838:     for (i=levels-1; i>=0; i--) {
839:       KSPView(mglevels[i]->smoothd,viewer);
840:       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
841:         KSPView(mglevels[i]->smoothu,viewer);
842:       }
843:     }
844:   } else if (isdraw) {
845:     PetscDraw draw;
846:     PetscReal x,w,y,bottom,th;
847:     PetscViewerDrawGetDraw(viewer,0,&draw);
848:     PetscDrawGetCurrentPoint(draw,&x,&y);
849:     PetscDrawStringGetSize(draw,NULL,&th);
850:     bottom = y - th;
851:     for (i=levels-1; i>=0; i--) {
852:       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
853:         PetscDrawPushCurrentPoint(draw,x,bottom);
854:         KSPView(mglevels[i]->smoothd,viewer);
855:         PetscDrawPopCurrentPoint(draw);
856:       } else {
857:         w    = 0.5*PetscMin(1.0-x,x);
858:         PetscDrawPushCurrentPoint(draw,x+w,bottom);
859:         KSPView(mglevels[i]->smoothd,viewer);
860:         PetscDrawPopCurrentPoint(draw);
861:         PetscDrawPushCurrentPoint(draw,x-w,bottom);
862:         KSPView(mglevels[i]->smoothu,viewer);
863:         PetscDrawPopCurrentPoint(draw);
864:       }
865:       PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);
866:       bottom -= th;
867:     }
868:   }
869:   return(0);
870: }

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

874: /*
875:     Calls setup for the KSP on each level
876: */
877: PetscErrorCode PCSetUp_MG(PC pc)
878: {
879:   PC_MG          *mg        = (PC_MG*)pc->data;
880:   PC_MG_Levels   **mglevels = mg->levels;
882:   PetscInt       i,n;
883:   PC             cpc;
884:   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
885:   Mat            dA,dB;
886:   Vec            tvec;
887:   DM             *dms;
888:   PetscViewer    viewer = NULL;
889:   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE;

892:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
893:   n = mglevels[0]->levels;
894:   /* FIX: Move this to PCSetFromOptions_MG? */
895:   if (mg->usedmfornumberoflevels) {
896:     PetscInt levels;
897:     DMGetRefineLevel(pc->dm,&levels);
898:     levels++;
899:     if (levels > n) { /* the problem is now being solved on a finer grid */
900:       PCMGSetLevels(pc,levels,NULL);
901:       n        = levels;
902:       PCSetFromOptions(pc); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
903:       mglevels = mg->levels;
904:     }
905:   }
906:   KSPGetPC(mglevels[0]->smoothd,&cpc);

908:   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
909:   /* so use those from global PC */
910:   /* Is this what we always want? What if user wants to keep old one? */
911:   KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);
912:   if (opsset) {
913:     Mat mmat;
914:     KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);
915:     if (mmat == pc->pmat) opsset = PETSC_FALSE;
916:   }

918:   /* Create CR solvers */
919:   PCMGGetAdaptCR(pc, &doCR);
920:   if (doCR) {
921:     const char *prefix;

923:     PCGetOptionsPrefix(pc, &prefix);
924:     for (i = 1; i < n; ++i) {
925:       PC   ipc, cr;
926:       char crprefix[128];

928:       KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr);
929:       KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE);
930:       PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i);
931:       KSPSetOptionsPrefix(mglevels[i]->cr, prefix);
932:       PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level);
933:       KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV);
934:       KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL);
935:       KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED);
936:       KSPGetPC(mglevels[i]->cr, &ipc);

938:       PCSetType(ipc, PCCOMPOSITE);
939:       PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE);
940:       PCCompositeAddPCType(ipc, PCSOR);
941:       CreateCR_Private(pc, i, &cr);
942:       PCCompositeAddPC(ipc, cr);
943:       PCDestroy(&cr);

945:       KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd);
946:       KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE);
947:       PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i);
948:       KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix);
949:       PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr);
950:     }
951:   }

953:   if (!opsset) {
954:     PCGetUseAmat(pc,&use_amat);
955:     if (use_amat) {
956:       PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
957:       KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);
958:     } else {
959:       PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
960:       KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);
961:     }
962:   }

964:   for (i=n-1; i>0; i--) {
965:     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
966:       missinginterpolate = PETSC_TRUE;
967:       continue;
968:     }
969:   }

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

977:   /*
978:    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
979:    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
980:   */
981:   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
982:         /* construct the interpolation from the DMs */
983:     Mat p;
984:     Vec rscale;
985:     PetscMalloc1(n,&dms);
986:     dms[n-1] = pc->dm;
987:     /* Separately create them so we do not get DMKSP interference between levels */
988:     for (i=n-2; i>-1; i--) {DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);}
989:         /*
990:            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
991:            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
992:            But it is safe to use -dm_mat_type.

994:            The mat type should not be hardcoded like this, we need to find a better way.
995:     DMSetMatType(dms[0],MATAIJ);
996:     */
997:     for (i=n-2; i>-1; i--) {
998:       DMKSP     kdm;
999:       PetscBool dmhasrestrict, dmhasinject;
1000:       KSPSetDM(mglevels[i]->smoothd,dms[i]);
1001:       if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);}
1002:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1003:         KSPSetDM(mglevels[i]->smoothu,dms[i]);
1004:         if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);}
1005:       }
1006:       if (mglevels[i]->cr) {
1007:         KSPSetDM(mglevels[i]->cr,dms[i]);
1008:         if (!needRestricts) {KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE);}
1009:       }
1010:       DMGetDMKSPWrite(dms[i],&kdm);
1011:       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
1012:        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
1013:       kdm->ops->computerhs = NULL;
1014:       kdm->rhsctx          = NULL;
1015:       if (!mglevels[i+1]->interpolate) {
1016:         DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);
1017:         PCMGSetInterpolation(pc,i+1,p);
1018:         if (rscale) {PCMGSetRScale(pc,i+1,rscale);}
1019:         VecDestroy(&rscale);
1020:         MatDestroy(&p);
1021:       }
1022:       DMHasCreateRestriction(dms[i],&dmhasrestrict);
1023:       if (dmhasrestrict && !mglevels[i+1]->restrct) {
1024:         DMCreateRestriction(dms[i],dms[i+1],&p);
1025:         PCMGSetRestriction(pc,i+1,p);
1026:         MatDestroy(&p);
1027:       }
1028:       DMHasCreateInjection(dms[i],&dmhasinject);
1029:       if (dmhasinject && !mglevels[i+1]->inject) {
1030:         DMCreateInjection(dms[i],dms[i+1],&p);
1031:         PCMGSetInjection(pc,i+1,p);
1032:         MatDestroy(&p);
1033:       }
1034:     }

1036:     for (i=n-2; i>-1; i--) {DMDestroy(&dms[i]);}
1037:     PetscFree(dms);
1038:   }

1040:   if (pc->dm && !pc->setupcalled) {
1041:     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
1042:     KSPSetDM(mglevels[n-1]->smoothd,pc->dm);
1043:     KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);
1044:     if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
1045:       KSPSetDM(mglevels[n-1]->smoothu,pc->dm);
1046:       KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);
1047:     }
1048:     if (mglevels[n-1]->cr) {
1049:       KSPSetDM(mglevels[n-1]->cr,pc->dm);
1050:       KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE);
1051:     }
1052:   }

1054:   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
1055:     Mat       A,B;
1056:     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
1057:     MatReuse  reuse = MAT_INITIAL_MATRIX;

1059:     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
1060:     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
1061:     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
1062:     for (i=n-2; i>-1; i--) {
1063:       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");
1064:       if (!mglevels[i+1]->interpolate) {
1065:         PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);
1066:       }
1067:       if (!mglevels[i+1]->restrct) {
1068:         PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);
1069:       }
1070:       if (reuse == MAT_REUSE_MATRIX) {
1071:         KSPGetOperators(mglevels[i]->smoothd,&A,&B);
1072:       }
1073:       if (doA) {
1074:         MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);
1075:       }
1076:       if (doB) {
1077:         MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);
1078:       }
1079:       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
1080:       if (!doA && dAeqdB) {
1081:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)B);}
1082:         A = B;
1083:       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
1084:         KSPGetOperators(mglevels[i]->smoothd,&A,NULL);
1085:         PetscObjectReference((PetscObject)A);
1086:       }
1087:       if (!doB && dAeqdB) {
1088:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)A);}
1089:         B = A;
1090:       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
1091:         KSPGetOperators(mglevels[i]->smoothd,NULL,&B);
1092:         PetscObjectReference((PetscObject)B);
1093:       }
1094:       if (reuse == MAT_INITIAL_MATRIX) {
1095:         KSPSetOperators(mglevels[i]->smoothd,A,B);
1096:         PetscObjectDereference((PetscObject)A);
1097:         PetscObjectDereference((PetscObject)B);
1098:       }
1099:       dA = A;
1100:       dB = B;
1101:     }
1102:   }

1104:   /* Adapt interpolation matrices */
1105:   if (mg->adaptInterpolation) {
1106:     mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */
1107:     for (i = 0; i < n; ++i) {
1108:       PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);
1109:       if (i) {PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mg->Nc, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace);}
1110:     }
1111:     for (i = n-2; i > -1; --i) {
1112:       PCMGRecomputeLevelOperators_Internal(pc, i);
1113:     }
1114:   }

1116:   if (needRestricts && pc->dm) {
1117:     for (i=n-2; i>=0; i--) {
1118:       DM  dmfine,dmcoarse;
1119:       Mat Restrict,Inject;
1120:       Vec rscale;
1121:       KSPGetDM(mglevels[i+1]->smoothd,&dmfine);
1122:       KSPGetDM(mglevels[i]->smoothd,&dmcoarse);
1123:       PCMGGetRestriction(pc,i+1,&Restrict);
1124:       PCMGGetRScale(pc,i+1,&rscale);
1125:       PCMGGetInjection(pc,i+1,&Inject);
1126:       DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);
1127:     }
1128:   }

1130:   if (!pc->setupcalled) {
1131:     for (i=0; i<n; i++) {
1132:       KSPSetFromOptions(mglevels[i]->smoothd);
1133:     }
1134:     for (i=1; i<n; i++) {
1135:       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
1136:         KSPSetFromOptions(mglevels[i]->smoothu);
1137:       }
1138:       if (mglevels[i]->cr) {
1139:         KSPSetFromOptions(mglevels[i]->cr);
1140:       }
1141:     }
1142:     /* insure that if either interpolation or restriction is set the other other one is set */
1143:     for (i=1; i<n; i++) {
1144:       PCMGGetInterpolation(pc,i,NULL);
1145:       PCMGGetRestriction(pc,i,NULL);
1146:     }
1147:     for (i=0; i<n-1; i++) {
1148:       if (!mglevels[i]->b) {
1149:         Vec *vec;
1150:         KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);
1151:         PCMGSetRhs(pc,i,*vec);
1152:         VecDestroy(vec);
1153:         PetscFree(vec);
1154:       }
1155:       if (!mglevels[i]->r && i) {
1156:         VecDuplicate(mglevels[i]->b,&tvec);
1157:         PCMGSetR(pc,i,tvec);
1158:         VecDestroy(&tvec);
1159:       }
1160:       if (!mglevels[i]->x) {
1161:         VecDuplicate(mglevels[i]->b,&tvec);
1162:         PCMGSetX(pc,i,tvec);
1163:         VecDestroy(&tvec);
1164:       }
1165:       if (doCR) {
1166:         VecDuplicate(mglevels[i]->b,&mglevels[i]->crx);
1167:         VecDuplicate(mglevels[i]->b,&mglevels[i]->crb);
1168:         VecZeroEntries(mglevels[i]->crb);
1169:       }
1170:     }
1171:     if (n != 1 && !mglevels[n-1]->r) {
1172:       /* PCMGSetR() on the finest level if user did not supply it */
1173:       Vec *vec;
1174:       KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);
1175:       PCMGSetR(pc,n-1,*vec);
1176:       VecDestroy(vec);
1177:       PetscFree(vec);
1178:     }
1179:     if (doCR) {
1180:       VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx);
1181:       VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb);
1182:       VecZeroEntries(mglevels[n-1]->crb);
1183:     }
1184:   }

1186:   if (pc->dm) {
1187:     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
1188:     for (i=0; i<n-1; i++) {
1189:       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
1190:     }
1191:   }

1193:   for (i=1; i<n; i++) {
1194:     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) {
1195:       /* if doing only down then initial guess is zero */
1196:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
1197:     }
1198:     if (mglevels[i]->cr) {KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);}
1199:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1200:     KSPSetUp(mglevels[i]->smoothd);
1201:     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1202:       pc->failedreason = PC_SUBPC_ERROR;
1203:     }
1204:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1205:     if (!mglevels[i]->residual) {
1206:       Mat mat;
1207:       KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);
1208:       PCMGSetResidual(pc,i,PCMGResidualDefault,mat);
1209:     }
1210:     if (!mglevels[i]->residualtranspose) {
1211:       Mat mat;
1212:       KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);
1213:       PCMGSetResidualTranspose(pc,i,PCMGResidualTransposeDefault,mat);
1214:     }
1215:   }
1216:   for (i=1; i<n; i++) {
1217:     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
1218:       Mat downmat,downpmat;

1220:       /* check if operators have been set for up, if not use down operators to set them */
1221:       KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);
1222:       if (!opsset) {
1223:         KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
1224:         KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);
1225:       }

1227:       KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);
1228:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1229:       KSPSetUp(mglevels[i]->smoothu);
1230:       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
1231:         pc->failedreason = PC_SUBPC_ERROR;
1232:       }
1233:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1234:     }
1235:     if (mglevels[i]->cr) {
1236:       Mat downmat,downpmat;

1238:       /* check if operators have been set for up, if not use down operators to set them */
1239:       KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL);
1240:       if (!opsset) {
1241:         KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
1242:         KSPSetOperators(mglevels[i]->cr,downmat,downpmat);
1243:       }

1245:       KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);
1246:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1247:       KSPSetUp(mglevels[i]->cr);
1248:       if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) {
1249:         pc->failedreason = PC_SUBPC_ERROR;
1250:       }
1251:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1252:     }
1253:   }

1255:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
1256:   KSPSetUp(mglevels[0]->smoothd);
1257:   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1258:     pc->failedreason = PC_SUBPC_ERROR;
1259:   }
1260:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}

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

1267:    Only support one or the other at the same time.
1268:   */
1269: #if defined(PETSC_USE_SOCKET_VIEWER)
1270:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);
1271:   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
1272:   dump = PETSC_FALSE;
1273: #endif
1274:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);
1275:   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));

1277:   if (viewer) {
1278:     for (i=1; i<n; i++) {
1279:       MatView(mglevels[i]->restrct,viewer);
1280:     }
1281:     for (i=0; i<n; i++) {
1282:       KSPGetPC(mglevels[i]->smoothd,&pc);
1283:       MatView(pc->mat,viewer);
1284:     }
1285:   }
1286:   return(0);
1287: }

1289: /* -------------------------------------------------------------------------------------*/

1291: PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
1292: {
1293:   PC_MG *mg = (PC_MG *) pc->data;

1296:   *levels = mg->nlevels;
1297:   return(0);
1298: }

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

1303:    Not Collective

1305:    Input Parameter:
1306: .  pc - the preconditioner context

1308:    Output parameter:
1309: .  levels - the number of levels

1311:    Level: advanced

1313: .seealso: PCMGSetLevels()
1314: @*/
1315: PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
1316: {

1322:   *levels = 0;
1323:   PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));
1324:   return(0);
1325: }

1327: /*@
1328:    PCMGSetType - Determines the form of multigrid to use:
1329:    multiplicative, additive, full, or the Kaskade algorithm.

1331:    Logically Collective on PC

1333:    Input Parameters:
1334: +  pc - the preconditioner context
1335: -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
1336:    PC_MG_FULL, PC_MG_KASKADE

1338:    Options Database Key:
1339: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
1340:    additive, full, kaskade

1342:    Level: advanced

1344: .seealso: PCMGSetLevels()
1345: @*/
1346: PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
1347: {
1348:   PC_MG *mg = (PC_MG*)pc->data;

1353:   mg->am = form;
1354:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1355:   else pc->ops->applyrichardson = NULL;
1356:   return(0);
1357: }

1359: /*@
1360:    PCMGGetType - Determines the form of multigrid to use:
1361:    multiplicative, additive, full, or the Kaskade algorithm.

1363:    Logically Collective on PC

1365:    Input Parameter:
1366: .  pc - the preconditioner context

1368:    Output Parameter:
1369: .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE

1371:    Level: advanced

1373: .seealso: PCMGSetLevels()
1374: @*/
1375: PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1376: {
1377:   PC_MG *mg = (PC_MG*)pc->data;

1381:   *type = mg->am;
1382:   return(0);
1383: }

1385: /*@
1386:    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
1387:    complicated cycling.

1389:    Logically Collective on PC

1391:    Input Parameters:
1392: +  pc - the multigrid context
1393: -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W

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

1398:    Level: advanced

1400: .seealso: PCMGSetCycleTypeOnLevel()
1401: @*/
1402: PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
1403: {
1404:   PC_MG        *mg        = (PC_MG*)pc->data;
1405:   PC_MG_Levels **mglevels = mg->levels;
1406:   PetscInt     i,levels;

1411:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1412:   levels = mglevels[0]->levels;
1413:   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
1414:   return(0);
1415: }

1417: /*@
1418:    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1419:          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used

1421:    Logically Collective on PC

1423:    Input Parameters:
1424: +  pc - the multigrid context
1425: -  n - number of cycles (default is 1)

1427:    Options Database Key:
1428: .  -pc_mg_multiplicative_cycles n

1430:    Level: advanced

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

1435: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1436: @*/
1437: PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1438: {
1439:   PC_MG        *mg        = (PC_MG*)pc->data;

1444:   mg->cyclesperpcapply = n;
1445:   return(0);
1446: }

1448: PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1449: {
1450:   PC_MG *mg = (PC_MG*)pc->data;

1453:   mg->galerkin = use;
1454:   return(0);
1455: }

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

1461:    Logically Collective on PC

1463:    Input Parameters:
1464: +  pc - the multigrid context
1465: -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE

1467:    Options Database Key:
1468: .  -pc_mg_galerkin <both,pmat,mat,none>

1470:    Level: intermediate

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

1476: .seealso: PCMGGetGalerkin(), PCMGGalerkinType

1478: @*/
1479: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1480: {

1485:   PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));
1486:   return(0);
1487: }

1489: /*@
1490:    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1491:       A_i-1 = r_i * A_i * p_i

1493:    Not Collective

1495:    Input Parameter:
1496: .  pc - the multigrid context

1498:    Output Parameter:
1499: .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL

1501:    Level: intermediate

1503: .seealso: PCMGSetGalerkin(), PCMGGalerkinType

1505: @*/
1506: PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1507: {
1508:   PC_MG *mg = (PC_MG*)pc->data;

1512:   *galerkin = mg->galerkin;
1513:   return(0);
1514: }

1516: PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1517: {
1518:   PC_MG *mg = (PC_MG *) pc->data;

1521:   mg->adaptInterpolation = adapt;
1522:   return(0);
1523: }

1525: PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1526: {
1527:   PC_MG *mg = (PC_MG *) pc->data;

1530:   *adapt = mg->adaptInterpolation;
1531:   return(0);
1532: }

1534: PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr)
1535: {
1536:   PC_MG *mg = (PC_MG *) pc->data;

1539:   mg->compatibleRelaxation = cr;
1540:   return(0);
1541: }

1543: PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr)
1544: {
1545:   PC_MG *mg = (PC_MG *) pc->data;

1548:   *cr = mg->compatibleRelaxation;
1549:   return(0);
1550: }

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

1555:   Logically Collective on PC

1557:   Input Parameters:
1558: + pc    - the multigrid context
1559: - adapt - flag for adaptation of the interpolator

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

1566:   Level: intermediate

1568: .keywords: MG, set, Galerkin
1569: .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1570: @*/
1571: PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1572: {

1577:   PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));
1578:   return(0);
1579: }

1581: /*@
1582:   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.

1584:   Not collective

1586:   Input Parameter:
1587: . pc    - the multigrid context

1589:   Output Parameter:
1590: . adapt - flag for adaptation of the interpolator

1592:   Level: intermediate

1594: .keywords: MG, set, Galerkin
1595: .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1596: @*/
1597: PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1598: {

1604:   PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));
1605:   return(0);
1606: }

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

1611:   Logically Collective on PC

1613:   Input Parameters:
1614: + pc - the multigrid context
1615: - cr - flag for compatible relaxation

1617:   Options Database Keys:
1618: . -pc_mg_adapt_cr - Turn on compatible relaxation

1620:   Level: intermediate

1622: .keywords: MG, set, Galerkin
1623: .seealso: PCMGGetAdaptCR(), PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1624: @*/
1625: PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr)
1626: {

1631:   PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr));
1632:   return(0);
1633: }

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

1638:   Not collective

1640:   Input Parameter:
1641: . pc    - the multigrid context

1643:   Output Parameter:
1644: . cr - flag for compatible relaxaion

1646:   Level: intermediate

1648: .keywords: MG, set, Galerkin
1649: .seealso: PCMGSetAdaptCR(), PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1650: @*/
1651: PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr)
1652: {

1658:   PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr));
1659:   return(0);
1660: }

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

1667:    Logically Collective on PC

1669:    Input Parameters:
1670: +  mg - the multigrid context
1671: -  n - the number of smoothing steps

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

1676:    Level: advanced

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

1682: .seealso: PCMGSetDistinctSmoothUp()
1683: @*/
1684: PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1685: {
1686:   PC_MG          *mg        = (PC_MG*)pc->data;
1687:   PC_MG_Levels   **mglevels = mg->levels;
1689:   PetscInt       i,levels;

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

1697:   for (i=1; i<levels; i++) {
1698:     KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1699:     KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1700:     mg->default_smoothu = n;
1701:     mg->default_smoothd = n;
1702:   }
1703:   return(0);
1704: }

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

1710:    Logically Collective on PC

1712:    Input Parameters:
1713: .  pc - the preconditioner context

1715:    Options Database Key:
1716: .  -pc_mg_distinct_smoothup

1718:    Level: advanced

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

1724: .seealso: PCMGSetNumberSmooth()
1725: @*/
1726: PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1727: {
1728:   PC_MG          *mg        = (PC_MG*)pc->data;
1729:   PC_MG_Levels   **mglevels = mg->levels;
1731:   PetscInt       i,levels;
1732:   KSP            subksp;

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

1739:   for (i=1; i<levels; i++) {
1740:     const char *prefix = NULL;
1741:     /* make sure smoother up and down are different */
1742:     PCMGGetSmootherUp(pc,i,&subksp);
1743:     KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);
1744:     KSPSetOptionsPrefix(subksp,prefix);
1745:     KSPAppendOptionsPrefix(subksp,"up_");
1746:   }
1747:   return(0);
1748: }

1750: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1751: PetscErrorCode  PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[])
1752: {
1753:   PC_MG          *mg        = (PC_MG*)pc->data;
1754:   PC_MG_Levels   **mglevels = mg->levels;
1755:   Mat            *mat;
1756:   PetscInt       l;

1760:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1761:   PetscMalloc1(mg->nlevels,&mat);
1762:   for (l=1; l< mg->nlevels; l++) {
1763:     mat[l-1] = mglevels[l]->interpolate;
1764:     PetscObjectReference((PetscObject)mat[l-1]);
1765:   }
1766:   *num_levels = mg->nlevels;
1767:   *interpolations = mat;
1768:   return(0);
1769: }

1771: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1772: PetscErrorCode  PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
1773: {
1774:   PC_MG          *mg        = (PC_MG*)pc->data;
1775:   PC_MG_Levels   **mglevels = mg->levels;
1776:   PetscInt       l;
1777:   Mat            *mat;

1781:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1782:   PetscMalloc1(mg->nlevels,&mat);
1783:   for (l=0; l<mg->nlevels-1; l++) {
1784:     KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));
1785:     PetscObjectReference((PetscObject)mat[l]);
1786:   }
1787:   *num_levels = mg->nlevels;
1788:   *coarseOperators = mat;
1789:   return(0);
1790: }

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

1795:   Not collective

1797:   Input Parameters:
1798: + name     - name of the constructor
1799: - function - constructor routine

1801:   Notes:
1802:   Calling sequence for the routine:
1803: $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1804: $   pc        - The PC object
1805: $   l         - The multigrid level, 0 is the coarse level
1806: $   dm        - The DM for this level
1807: $   smooth    - The level smoother
1808: $   Nc        - The size of the coarse space
1809: $   initGuess - Basis for an initial guess for the space
1810: $   coarseSp  - A basis for the computed coarse space

1812:   Level: advanced

1814: .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister()
1815: @*/
1816: PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1817: {

1821:   PCInitializePackage();
1822:   PetscFunctionListAdd(&PCMGCoarseList,name,function);
1823:   return(0);
1824: }

1826: /*@C
1827:   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.

1829:   Not collective

1831:   Input Parameter:
1832: . name     - name of the constructor

1834:   Output Parameter:
1835: . function - constructor routine

1837:   Notes:
1838:   Calling sequence for the routine:
1839: $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1840: $   pc        - The PC object
1841: $   l         - The multigrid level, 0 is the coarse level
1842: $   dm        - The DM for this level
1843: $   smooth    - The level smoother
1844: $   Nc        - The size of the coarse space
1845: $   initGuess - Basis for an initial guess for the space
1846: $   coarseSp  - A basis for the computed coarse space

1848:   Level: advanced

1850: .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister()
1851: @*/
1852: PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1853: {

1857:   PetscFunctionListFind(PCMGCoarseList,name,function);
1858:   return(0);
1859: }

1861: /* ----------------------------------------------------------------------------------------*/

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

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

1880:    Notes:
1881:     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

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

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

1890:    Level: intermediate

1892: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1893:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1894:            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1895:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1896:            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1897: M*/

1899: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1900: {
1901:   PC_MG          *mg;

1905:   PetscNewLog(pc,&mg);
1906:   pc->data     = mg;
1907:   mg->nlevels  = -1;
1908:   mg->am       = PC_MG_MULTIPLICATIVE;
1909:   mg->galerkin = PC_MG_GALERKIN_NONE;
1910:   mg->adaptInterpolation = PETSC_FALSE;
1911:   mg->Nc                 = -1;
1912:   mg->eigenvalue         = -1;

1914:   pc->useAmat = PETSC_TRUE;

1916:   pc->ops->apply          = PCApply_MG;
1917:   pc->ops->applytranspose = PCApplyTranspose_MG;
1918:   pc->ops->matapply       = PCMatApply_MG;
1919:   pc->ops->setup          = PCSetUp_MG;
1920:   pc->ops->reset          = PCReset_MG;
1921:   pc->ops->destroy        = PCDestroy_MG;
1922:   pc->ops->setfromoptions = PCSetFromOptions_MG;
1923:   pc->ops->view           = PCView_MG;

1925:   PetscObjectComposedDataRegister(&mg->eigenvalue);
1926:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1927:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);
1928:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);
1929:   PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);
1930:   PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);
1931:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);
1932:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);
1933:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG);
1934:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG);
1935:   return(0);
1936: }