Actual source code: mg.c

petsc-3.14.6 2021-03-30
Report Typos and Errors

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

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

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

 22:   if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
 23:   KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);  /* pre-smooth */
 24:   KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);
 25:   if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
 26:   if (mglevels->level) {  /* not the coarsest grid */
 27:     if (mglevels->eventresidual) {PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);}
 28:     (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);
 29:     if (mglevels->eventresidual) {PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);}

 31:     /* if on finest level and have convergence criteria set */
 32:     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
 33:       PetscReal rnorm;
 34:       VecNorm(mglevels->r,NORM_2,&rnorm);
 35:       if (rnorm <= mg->ttol) {
 36:         if (rnorm < mg->abstol) {
 37:           *reason = PCRICHARDSON_CONVERGED_ATOL;
 38:           PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);
 39:         } else {
 40:           *reason = PCRICHARDSON_CONVERGED_RTOL;
 41:           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);
 42:         }
 43:         return(0);
 44:       }
 45:     }

 47:     mgc = *(mglevelsin - 1);
 48:     if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
 49:     MatRestrict(mglevels->restrct,mglevels->r,mgc->b);
 50:     if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
 51:     VecSet(mgc->x,0.0);
 52:     while (cycles--) {
 53:       PCMGMCycle_Private(pc,mglevelsin-1,reason);
 54:     }
 55:     if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
 56:     MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);
 57:     if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
 58:     if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
 59:     KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);    /* post smooth */
 60:     KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);
 61:     if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
 62:   }
 63:   return(0);
 64: }

 66: 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)
 67: {
 68:   PC_MG          *mg        = (PC_MG*)pc->data;
 69:   PC_MG_Levels   **mglevels = mg->levels;
 71:   PC             tpc;
 72:   PetscBool      changeu,changed;
 73:   PetscInt       levels = mglevels[0]->levels,i;

 76:   /* When the DM is supplying the matrix then it will not exist until here */
 77:   for (i=0; i<levels; i++) {
 78:     if (!mglevels[i]->A) {
 79:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
 80:       PetscObjectReference((PetscObject)mglevels[i]->A);
 81:     }
 82:   }

 84:   KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
 85:   PCPreSolveChangeRHS(tpc,&changed);
 86:   KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
 87:   PCPreSolveChangeRHS(tpc,&changeu);
 88:   if (!changed && !changeu) {
 89:     VecDestroy(&mglevels[levels-1]->b);
 90:     mglevels[levels-1]->b = b;
 91:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
 92:     if (!mglevels[levels-1]->b) {
 93:       Vec *vec;

 95:       KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
 96:       mglevels[levels-1]->b = *vec;
 97:       PetscFree(vec);
 98:     }
 99:     VecCopy(b,mglevels[levels-1]->b);
100:   }
101:   mglevels[levels-1]->x = x;

103:   mg->rtol   = rtol;
104:   mg->abstol = abstol;
105:   mg->dtol   = dtol;
106:   if (rtol) {
107:     /* compute initial residual norm for relative convergence test */
108:     PetscReal rnorm;
109:     if (zeroguess) {
110:       VecNorm(b,NORM_2,&rnorm);
111:     } else {
112:       (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);
113:       VecNorm(w,NORM_2,&rnorm);
114:     }
115:     mg->ttol = PetscMax(rtol*rnorm,abstol);
116:   } else if (abstol) mg->ttol = abstol;
117:   else mg->ttol = 0.0;

119:   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
120:      stop prematurely due to small residual */
121:   for (i=1; i<levels; i++) {
122:     KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
123:     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
124:       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
125:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
126:       KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
127:     }
128:   }

130:   *reason = (PCRichardsonConvergedReason)0;
131:   for (i=0; i<its; i++) {
132:     PCMGMCycle_Private(pc,mglevels+levels-1,reason);
133:     if (*reason) break;
134:   }
135:   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
136:   *outits = i;
137:   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
138:   return(0);
139: }

141: PetscErrorCode PCReset_MG(PC pc)
142: {
143:   PC_MG          *mg        = (PC_MG*)pc->data;
144:   PC_MG_Levels   **mglevels = mg->levels;
146:   PetscInt       i,c,n;

149:   if (mglevels) {
150:     n = mglevels[0]->levels;
151:     for (i=0; i<n-1; i++) {
152:       VecDestroy(&mglevels[i+1]->r);
153:       VecDestroy(&mglevels[i]->b);
154:       VecDestroy(&mglevels[i]->x);
155:       MatDestroy(&mglevels[i+1]->restrct);
156:       MatDestroy(&mglevels[i+1]->interpolate);
157:       MatDestroy(&mglevels[i+1]->inject);
158:       VecDestroy(&mglevels[i+1]->rscale);
159:     }
160:     /* this is not null only if the smoother on the finest level
161:        changes the rhs during PreSolve */
162:     VecDestroy(&mglevels[n-1]->b);

164:     for (i=0; i<n; i++) {
165:       if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {VecDestroy(&mglevels[i]->coarseSpace[c]);}
166:       PetscFree(mglevels[i]->coarseSpace);
167:       mglevels[i]->coarseSpace = NULL;
168:       MatDestroy(&mglevels[i]->A);
169:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
170:         KSPReset(mglevels[i]->smoothd);
171:       }
172:       KSPReset(mglevels[i]->smoothu);
173:     }
174:     mg->Nc = 0;
175:   }
176:   return(0);
177: }

179: PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
180: {
182:   PC_MG          *mg        = (PC_MG*)pc->data;
183:   MPI_Comm       comm;
184:   PC_MG_Levels   **mglevels = mg->levels;
185:   PCMGType       mgtype     = mg->am;
186:   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
187:   PetscInt       i;
188:   PetscMPIInt    size;
189:   const char     *prefix;
190:   PC             ipc;
191:   PetscInt       n;

196:   if (mg->nlevels == levels) return(0);
197:   PetscObjectGetComm((PetscObject)pc,&comm);
198:   if (mglevels) {
199:     mgctype = mglevels[0]->cycles;
200:     /* changing the number of levels so free up the previous stuff */
201:     PCReset_MG(pc);
202:     n    = mglevels[0]->levels;
203:     for (i=0; i<n; i++) {
204:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
205:         KSPDestroy(&mglevels[i]->smoothd);
206:       }
207:       KSPDestroy(&mglevels[i]->smoothu);
208:       PetscFree(mglevels[i]);
209:     }
210:     PetscFree(mg->levels);
211:   }

213:   mg->nlevels = levels;

215:   PetscMalloc1(levels,&mglevels);
216:   PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));

218:   PCGetOptionsPrefix(pc,&prefix);

220:   mg->stageApply = 0;
221:   for (i=0; i<levels; i++) {
222:     PetscNewLog(pc,&mglevels[i]);

224:     mglevels[i]->level               = i;
225:     mglevels[i]->levels              = levels;
226:     mglevels[i]->cycles              = mgctype;
227:     mg->default_smoothu              = 2;
228:     mg->default_smoothd              = 2;
229:     mglevels[i]->eventsmoothsetup    = 0;
230:     mglevels[i]->eventsmoothsolve    = 0;
231:     mglevels[i]->eventresidual       = 0;
232:     mglevels[i]->eventinterprestrict = 0;

234:     if (comms) comm = comms[i];
235:     if (comm != MPI_COMM_NULL) {
236:       KSPCreate(comm,&mglevels[i]->smoothd);
237:       KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);
238:       PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);
239:       KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);
240:       PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);
241:       if (i || levels == 1) {
242:         char tprefix[128];

244:         KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);
245:         KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);
246:         KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);
247:         KSPGetPC(mglevels[i]->smoothd,&ipc);
248:         PCSetType(ipc,PCSOR);
249:         KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);

251:         PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);
252:         KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);
253:       } else {
254:         KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");

256:         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
257:         KSPSetType(mglevels[0]->smoothd,KSPPREONLY);
258:         KSPGetPC(mglevels[0]->smoothd,&ipc);
259:         MPI_Comm_size(comm,&size);
260:         if (size > 1) {
261:           PCSetType(ipc,PCREDUNDANT);
262:         } else {
263:           PCSetType(ipc,PCLU);
264:         }
265:         PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);
266:       }
267:       PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);
268:     }
269:     mglevels[i]->smoothu = mglevels[i]->smoothd;
270:     mg->rtol             = 0.0;
271:     mg->abstol           = 0.0;
272:     mg->dtol             = 0.0;
273:     mg->ttol             = 0.0;
274:     mg->cyclesperpcapply = 1;
275:   }
276:   mg->levels = mglevels;
277:   PCMGSetType(pc,mgtype);
278:   return(0);
279: }

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

285:    Logically Collective on PC

287:    Input Parameters:
288: +  pc - the preconditioner context
289: .  levels - the number of levels
290: -  comms - optional communicators for each level; this is to allow solving the coarser problems
291:            on smaller sets of processes. For processes that are not included in the computation
292:            you must pass MPI_COMM_NULL.

294:    Level: intermediate

296:    Notes:
297:      If the number of levels is one then the multigrid uses the -mg_levels prefix
298:      for setting the level options rather than the -mg_coarse prefix.

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

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

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


318: .seealso: PCMGSetType(), PCMGGetLevels()
319: @*/
320: PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
321: {

327:   PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));
328:   return(0);
329: }


332: PetscErrorCode PCDestroy_MG(PC pc)
333: {
335:   PC_MG          *mg        = (PC_MG*)pc->data;
336:   PC_MG_Levels   **mglevels = mg->levels;
337:   PetscInt       i,n;

340:   PCReset_MG(pc);
341:   if (mglevels) {
342:     n = mglevels[0]->levels;
343:     for (i=0; i<n; i++) {
344:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
345:         KSPDestroy(&mglevels[i]->smoothd);
346:       }
347:       KSPDestroy(&mglevels[i]->smoothu);
348:       PetscFree(mglevels[i]);
349:     }
350:     PetscFree(mg->levels);
351:   }
352:   PetscFree(pc->data);
353:   PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);
354:   PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);
355:   return(0);
356: }



360: extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
361: extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
362: extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);

364: /*
365:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
366:              or full cycle of multigrid.

368:   Note:
369:   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
370: */
371: static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
372: {
373:   PC_MG          *mg        = (PC_MG*)pc->data;
374:   PC_MG_Levels   **mglevels = mg->levels;
376:   PC             tpc;
377:   PetscInt       levels = mglevels[0]->levels,i;
378:   PetscBool      changeu,changed;

381:   if (mg->stageApply) {PetscLogStagePush(mg->stageApply);}
382:   /* When the DM is supplying the matrix then it will not exist until here */
383:   for (i=0; i<levels; i++) {
384:     if (!mglevels[i]->A) {
385:       KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
386:       PetscObjectReference((PetscObject)mglevels[i]->A);
387:     }
388:   }

390:   KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
391:   PCPreSolveChangeRHS(tpc,&changed);
392:   KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
393:   PCPreSolveChangeRHS(tpc,&changeu);
394:   if (!changeu && !changed) {
395:     VecDestroy(&mglevels[levels-1]->b);
396:     mglevels[levels-1]->b = b;
397:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
398:     if (!mglevels[levels-1]->b) {
399:       Vec *vec;

401:       KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
402:       mglevels[levels-1]->b = *vec;
403:       PetscFree(vec);
404:     }
405:     VecCopy(b,mglevels[levels-1]->b);
406:   }
407:   mglevels[levels-1]->x = x;

409:   if (mg->am == PC_MG_MULTIPLICATIVE) {
410:     VecSet(x,0.0);
411:     for (i=0; i<mg->cyclesperpcapply; i++) {
412:       PCMGMCycle_Private(pc,mglevels+levels-1,NULL);
413:     }
414:   } else if (mg->am == PC_MG_ADDITIVE) {
415:     PCMGACycle_Private(pc,mglevels);
416:   } else if (mg->am == PC_MG_KASKADE) {
417:     PCMGKCycle_Private(pc,mglevels);
418:   } else {
419:     PCMGFCycle_Private(pc,mglevels);
420:   }
421:   if (mg->stageApply) {PetscLogStagePop();}
422:   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
423:   return(0);
424: }


427: PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
428: {
429:   PetscErrorCode   ierr;
430:   PetscInt         levels,cycles;
431:   PetscBool        flg, flg2;
432:   PC_MG            *mg = (PC_MG*)pc->data;
433:   PC_MG_Levels     **mglevels;
434:   PCMGType         mgtype;
435:   PCMGCycleType    mgctype;
436:   PCMGGalerkinType gtype;

439:   levels = PetscMax(mg->nlevels,1);
440:   PetscOptionsHead(PetscOptionsObject,"Multigrid options");
441:   PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
442:   if (!flg && !mg->levels && pc->dm) {
443:     DMGetRefineLevel(pc->dm,&levels);
444:     levels++;
445:     mg->usedmfornumberoflevels = PETSC_TRUE;
446:   }
447:   PCMGSetLevels(pc,levels,NULL);
448:   mglevels = mg->levels;

450:   mgctype = (PCMGCycleType) mglevels[0]->cycles;
451:   PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);
452:   if (flg) {
453:     PCMGSetCycleType(pc,mgctype);
454:   }
455:   gtype = mg->galerkin;
456:   PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);
457:   if (flg) {
458:     PCMGSetGalerkin(pc,gtype);
459:   }
460:   flg2 = PETSC_FALSE;
461:   PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);
462:   if (flg) {PCMGSetAdaptInterpolation(pc, flg2);}
463:   PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg);
464:   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);
465:   PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);
466:   flg = PETSC_FALSE;
467:   PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);
468:   if (flg) {
469:     PCMGSetDistinctSmoothUp(pc);
470:   }
471:   mgtype = mg->am;
472:   PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
473:   if (flg) {
474:     PCMGSetType(pc,mgtype);
475:   }
476:   if (mg->am == PC_MG_MULTIPLICATIVE) {
477:     PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);
478:     if (flg) {
479:       PCMGMultiplicativeSetCycles(pc,cycles);
480:     }
481:   }
482:   flg  = PETSC_FALSE;
483:   PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);
484:   if (flg) {
485:     PetscInt i;
486:     char     eventname[128];

488:     levels = mglevels[0]->levels;
489:     for (i=0; i<levels; i++) {
490:       sprintf(eventname,"MGSetup Level %d",(int)i);
491:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);
492:       sprintf(eventname,"MGSmooth Level %d",(int)i);
493:       PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);
494:       if (i) {
495:         sprintf(eventname,"MGResid Level %d",(int)i);
496:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);
497:         sprintf(eventname,"MGInterp Level %d",(int)i);
498:         PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);
499:       }
500:     }

502: #if defined(PETSC_USE_LOG)
503:     {
504:       const char    *sname = "MG Apply";
505:       PetscStageLog stageLog;
506:       PetscInt      st;

508:       PetscLogGetStageLog(&stageLog);
509:       for (st = 0; st < stageLog->numStages; ++st) {
510:         PetscBool same;

512:         PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
513:         if (same) mg->stageApply = st;
514:       }
515:       if (!mg->stageApply) {
516:         PetscLogStageRegister(sname, &mg->stageApply);
517:       }
518:     }
519: #endif
520:   }
521:   PetscOptionsTail();
522:   /* Check option consistency */
523:   PCMGGetGalerkin(pc, &gtype);
524:   PCMGGetAdaptInterpolation(pc, &flg);
525:   if (flg && (gtype >= PC_MG_GALERKIN_NONE)) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
526:   return(0);
527: }

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

534: #include <petscdraw.h>
535: PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
536: {
537:   PC_MG          *mg        = (PC_MG*)pc->data;
538:   PC_MG_Levels   **mglevels = mg->levels;
540:   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
541:   PetscBool      iascii,isbinary,isdraw;

544:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
545:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
546:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
547:   if (iascii) {
548:     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
549:     PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);
550:     if (mg->am == PC_MG_MULTIPLICATIVE) {
551:       PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);
552:     }
553:     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
554:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");
555:     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
556:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");
557:     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
558:       PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");
559:     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
560:       PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");
561:     } else {
562:       PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");
563:     }
564:     if (mg->view){
565:       (*mg->view)(pc,viewer);
566:     }
567:     for (i=0; i<levels; i++) {
568:       if (!i) {
569:         PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);
570:       } else {
571:         PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);
572:       }
573:       PetscViewerASCIIPushTab(viewer);
574:       KSPView(mglevels[i]->smoothd,viewer);
575:       PetscViewerASCIIPopTab(viewer);
576:       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
577:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
578:       } else if (i) {
579:         PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);
580:         PetscViewerASCIIPushTab(viewer);
581:         KSPView(mglevels[i]->smoothu,viewer);
582:         PetscViewerASCIIPopTab(viewer);
583:       }
584:     }
585:   } else if (isbinary) {
586:     for (i=levels-1; i>=0; i--) {
587:       KSPView(mglevels[i]->smoothd,viewer);
588:       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
589:         KSPView(mglevels[i]->smoothu,viewer);
590:       }
591:     }
592:   } else if (isdraw) {
593:     PetscDraw draw;
594:     PetscReal x,w,y,bottom,th;
595:     PetscViewerDrawGetDraw(viewer,0,&draw);
596:     PetscDrawGetCurrentPoint(draw,&x,&y);
597:     PetscDrawStringGetSize(draw,NULL,&th);
598:     bottom = y - th;
599:     for (i=levels-1; i>=0; i--) {
600:       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
601:         PetscDrawPushCurrentPoint(draw,x,bottom);
602:         KSPView(mglevels[i]->smoothd,viewer);
603:         PetscDrawPopCurrentPoint(draw);
604:       } else {
605:         w    = 0.5*PetscMin(1.0-x,x);
606:         PetscDrawPushCurrentPoint(draw,x+w,bottom);
607:         KSPView(mglevels[i]->smoothd,viewer);
608:         PetscDrawPopCurrentPoint(draw);
609:         PetscDrawPushCurrentPoint(draw,x-w,bottom);
610:         KSPView(mglevels[i]->smoothu,viewer);
611:         PetscDrawPopCurrentPoint(draw);
612:       }
613:       PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);
614:       bottom -= th;
615:     }
616:   }
617:   return(0);
618: }

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

622: /*
623:     Calls setup for the KSP on each level
624: */
625: PetscErrorCode PCSetUp_MG(PC pc)
626: {
627:   PC_MG          *mg        = (PC_MG*)pc->data;
628:   PC_MG_Levels   **mglevels = mg->levels;
630:   PetscInt       i,n;
631:   PC             cpc;
632:   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
633:   Mat            dA,dB;
634:   Vec            tvec;
635:   DM             *dms;
636:   PetscViewer    viewer = NULL;
637:   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;

640:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
641:   n = mglevels[0]->levels;
642:   /* FIX: Move this to PCSetFromOptions_MG? */
643:   if (mg->usedmfornumberoflevels) {
644:     PetscInt levels;
645:     DMGetRefineLevel(pc->dm,&levels);
646:     levels++;
647:     if (levels > n) { /* the problem is now being solved on a finer grid */
648:       PCMGSetLevels(pc,levels,NULL);
649:       n        = levels;
650:       PCSetFromOptions(pc); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
651:       mglevels = mg->levels;
652:     }
653:   }
654:   KSPGetPC(mglevels[0]->smoothd,&cpc);


657:   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
658:   /* so use those from global PC */
659:   /* Is this what we always want? What if user wants to keep old one? */
660:   KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);
661:   if (opsset) {
662:     Mat mmat;
663:     KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);
664:     if (mmat == pc->pmat) opsset = PETSC_FALSE;
665:   }

667:   if (!opsset) {
668:     PCGetUseAmat(pc,&use_amat);
669:     if (use_amat) {
670:       PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
671:       KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);
672:     } else {
673:       PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
674:       KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);
675:     }
676:   }

678:   for (i=n-1; i>0; i--) {
679:     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
680:       missinginterpolate = PETSC_TRUE;
681:       continue;
682:     }
683:   }

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


692:   /*
693:    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
694:    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
695:   */
696:   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
697:         /* construct the interpolation from the DMs */
698:     Mat p;
699:     Vec rscale;
700:     PetscMalloc1(n,&dms);
701:     dms[n-1] = pc->dm;
702:     /* Separately create them so we do not get DMKSP interference between levels */
703:     for (i=n-2; i>-1; i--) {DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);}
704:         /*
705:            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
706:            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
707:            But it is safe to use -dm_mat_type.

709:            The mat type should not be hardcoded like this, we need to find a better way.
710:     DMSetMatType(dms[0],MATAIJ);
711:     */
712:     for (i=n-2; i>-1; i--) {
713:       DMKSP     kdm;
714:       PetscBool dmhasrestrict, dmhasinject;
715:       KSPSetDM(mglevels[i]->smoothd,dms[i]);
716:       if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);}
717:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
718:         KSPSetDM(mglevels[i]->smoothu,dms[i]);
719:         if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);}
720:       }
721:       DMGetDMKSPWrite(dms[i],&kdm);
722:       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
723:        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
724:       kdm->ops->computerhs = NULL;
725:       kdm->rhsctx          = NULL;
726:       if (!mglevels[i+1]->interpolate) {
727:         DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);
728:         PCMGSetInterpolation(pc,i+1,p);
729:         if (rscale) {PCMGSetRScale(pc,i+1,rscale);}
730:         VecDestroy(&rscale);
731:         MatDestroy(&p);
732:       }
733:       DMHasCreateRestriction(dms[i],&dmhasrestrict);
734:       if (dmhasrestrict && !mglevels[i+1]->restrct){
735:         DMCreateRestriction(dms[i],dms[i+1],&p);
736:         PCMGSetRestriction(pc,i+1,p);
737:         MatDestroy(&p);
738:       }
739:       DMHasCreateInjection(dms[i],&dmhasinject);
740:       if (dmhasinject && !mglevels[i+1]->inject){
741:         DMCreateInjection(dms[i],dms[i+1],&p);
742:         PCMGSetInjection(pc,i+1,p);
743:         MatDestroy(&p);
744:       }
745:     }

747:     for (i=n-2; i>-1; i--) {DMDestroy(&dms[i]);}
748:     PetscFree(dms);
749:   }

751:   if (pc->dm && !pc->setupcalled) {
752:     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
753:     KSPSetDM(mglevels[n-1]->smoothd,pc->dm);
754:     KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);
755:     if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
756:       KSPSetDM(mglevels[n-1]->smoothu,pc->dm);
757:       KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);
758:     }
759:   }

761:   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
762:     Mat       A,B;
763:     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
764:     MatReuse  reuse = MAT_INITIAL_MATRIX;

766:     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
767:     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
768:     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
769:     for (i=n-2; i>-1; i--) {
770:       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");
771:       if (!mglevels[i+1]->interpolate) {
772:         PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);
773:       }
774:       if (!mglevels[i+1]->restrct) {
775:         PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);
776:       }
777:       if (reuse == MAT_REUSE_MATRIX) {
778:         KSPGetOperators(mglevels[i]->smoothd,&A,&B);
779:       }
780:       if (doA) {
781:         MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);
782:       }
783:       if (doB) {
784:         MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);
785:       }
786:       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
787:       if (!doA && dAeqdB) {
788:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)B);}
789:         A = B;
790:       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
791:         KSPGetOperators(mglevels[i]->smoothd,&A,NULL);
792:         PetscObjectReference((PetscObject)A);
793:       }
794:       if (!doB && dAeqdB) {
795:         if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)A);}
796:         B = A;
797:       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
798:         KSPGetOperators(mglevels[i]->smoothd,NULL,&B);
799:         PetscObjectReference((PetscObject)B);
800:       }
801:       if (reuse == MAT_INITIAL_MATRIX) {
802:         KSPSetOperators(mglevels[i]->smoothd,A,B);
803:         PetscObjectDereference((PetscObject)A);
804:         PetscObjectDereference((PetscObject)B);
805:       }
806:       dA = A;
807:       dB = B;
808:     }
809:   }


812:   /* Adapt interpolation matrices */
813:   if (mg->adaptInterpolation) {
814:     mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */
815:     for (i = 0; i < n; ++i) {
816:       PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);
817:       if (i) {PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mg->Nc, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace);}
818:     }
819:     for (i = n-2; i > -1; --i) {
820:       PCMGRecomputeLevelOperators_Internal(pc, i);
821:     }
822:   }

824:   if (needRestricts && pc->dm) {
825:     for (i=n-2; i>=0; i--) {
826:       DM  dmfine,dmcoarse;
827:       Mat Restrict,Inject;
828:       Vec rscale;
829:       KSPGetDM(mglevels[i+1]->smoothd,&dmfine);
830:       KSPGetDM(mglevels[i]->smoothd,&dmcoarse);
831:       PCMGGetRestriction(pc,i+1,&Restrict);
832:       PCMGGetRScale(pc,i+1,&rscale);
833:       PCMGGetInjection(pc,i+1,&Inject);
834:       DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);
835:     }
836:   }

838:   if (!pc->setupcalled) {
839:     for (i=0; i<n; i++) {
840:       KSPSetFromOptions(mglevels[i]->smoothd);
841:     }
842:     for (i=1; i<n; i++) {
843:       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
844:         KSPSetFromOptions(mglevels[i]->smoothu);
845:       }
846:     }
847:     /* insure that if either interpolation or restriction is set the other other one is set */
848:     for (i=1; i<n; i++) {
849:       PCMGGetInterpolation(pc,i,NULL);
850:       PCMGGetRestriction(pc,i,NULL);
851:     }
852:     for (i=0; i<n-1; i++) {
853:       if (!mglevels[i]->b) {
854:         Vec *vec;
855:         KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);
856:         PCMGSetRhs(pc,i,*vec);
857:         VecDestroy(vec);
858:         PetscFree(vec);
859:       }
860:       if (!mglevels[i]->r && i) {
861:         VecDuplicate(mglevels[i]->b,&tvec);
862:         PCMGSetR(pc,i,tvec);
863:         VecDestroy(&tvec);
864:       }
865:       if (!mglevels[i]->x) {
866:         VecDuplicate(mglevels[i]->b,&tvec);
867:         PCMGSetX(pc,i,tvec);
868:         VecDestroy(&tvec);
869:       }
870:     }
871:     if (n != 1 && !mglevels[n-1]->r) {
872:       /* PCMGSetR() on the finest level if user did not supply it */
873:       Vec *vec;
874:       KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);
875:       PCMGSetR(pc,n-1,*vec);
876:       VecDestroy(vec);
877:       PetscFree(vec);
878:     }
879:   }

881:   if (pc->dm) {
882:     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
883:     for (i=0; i<n-1; i++) {
884:       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
885:     }
886:   }

888:   for (i=1; i<n; i++) {
889:     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
890:       /* if doing only down then initial guess is zero */
891:       KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
892:     }
893:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
894:     KSPSetUp(mglevels[i]->smoothd);
895:     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
896:       pc->failedreason = PC_SUBPC_ERROR;
897:     }
898:     if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
899:     if (!mglevels[i]->residual) {
900:       Mat mat;
901:       KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);
902:       PCMGSetResidual(pc,i,PCMGResidualDefault,mat);
903:     }
904:   }
905:   for (i=1; i<n; i++) {
906:     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
907:       Mat downmat,downpmat;

909:       /* check if operators have been set for up, if not use down operators to set them */
910:       KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);
911:       if (!opsset) {
912:         KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
913:         KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);
914:       }

916:       KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);
917:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
918:       KSPSetUp(mglevels[i]->smoothu);
919:       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
920:         pc->failedreason = PC_SUBPC_ERROR;
921:       }
922:       if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
923:     }
924:   }

926:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
927:   KSPSetUp(mglevels[0]->smoothd);
928:   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
929:     pc->failedreason = PC_SUBPC_ERROR;
930:   }
931:   if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}

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

938:    Only support one or the other at the same time.
939:   */
940: #if defined(PETSC_USE_SOCKET_VIEWER)
941:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);
942:   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
943:   dump = PETSC_FALSE;
944: #endif
945:   PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);
946:   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));

948:   if (viewer) {
949:     for (i=1; i<n; i++) {
950:       MatView(mglevels[i]->restrct,viewer);
951:     }
952:     for (i=0; i<n; i++) {
953:       KSPGetPC(mglevels[i]->smoothd,&pc);
954:       MatView(pc->mat,viewer);
955:     }
956:   }
957:   return(0);
958: }

960: /* -------------------------------------------------------------------------------------*/

962: PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
963: {
964:   PC_MG *mg = (PC_MG *) pc->data;

967:   *levels = mg->nlevels;
968:   return(0);
969: }

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

974:    Not Collective

976:    Input Parameter:
977: .  pc - the preconditioner context

979:    Output parameter:
980: .  levels - the number of levels

982:    Level: advanced

984: .seealso: PCMGSetLevels()
985: @*/
986: PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
987: {

993:   *levels = 0;
994:   PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));
995:   return(0);
996: }

998: /*@
999:    PCMGSetType - Determines the form of multigrid to use:
1000:    multiplicative, additive, full, or the Kaskade algorithm.

1002:    Logically Collective on PC

1004:    Input Parameters:
1005: +  pc - the preconditioner context
1006: -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
1007:    PC_MG_FULL, PC_MG_KASKADE

1009:    Options Database Key:
1010: .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
1011:    additive, full, kaskade

1013:    Level: advanced

1015: .seealso: PCMGSetLevels()
1016: @*/
1017: PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
1018: {
1019:   PC_MG *mg = (PC_MG*)pc->data;

1024:   mg->am = form;
1025:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1026:   else pc->ops->applyrichardson = NULL;
1027:   return(0);
1028: }

1030: /*@
1031:    PCMGGetType - Determines the form of multigrid to use:
1032:    multiplicative, additive, full, or the Kaskade algorithm.

1034:    Logically Collective on PC

1036:    Input Parameter:
1037: .  pc - the preconditioner context

1039:    Output Parameter:
1040: .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE


1043:    Level: advanced

1045: .seealso: PCMGSetLevels()
1046: @*/
1047: PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1048: {
1049:   PC_MG *mg = (PC_MG*)pc->data;

1053:   *type = mg->am;
1054:   return(0);
1055: }

1057: /*@
1058:    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
1059:    complicated cycling.

1061:    Logically Collective on PC

1063:    Input Parameters:
1064: +  pc - the multigrid context
1065: -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W

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

1070:    Level: advanced

1072: .seealso: PCMGSetCycleTypeOnLevel()
1073: @*/
1074: PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
1075: {
1076:   PC_MG        *mg        = (PC_MG*)pc->data;
1077:   PC_MG_Levels **mglevels = mg->levels;
1078:   PetscInt     i,levels;

1083:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1084:   levels = mglevels[0]->levels;
1085:   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
1086:   return(0);
1087: }

1089: /*@
1090:    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1091:          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used

1093:    Logically Collective on PC

1095:    Input Parameters:
1096: +  pc - the multigrid context
1097: -  n - number of cycles (default is 1)

1099:    Options Database Key:
1100: .  -pc_mg_multiplicative_cycles n

1102:    Level: advanced

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

1107: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1108: @*/
1109: PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1110: {
1111:   PC_MG        *mg        = (PC_MG*)pc->data;

1116:   mg->cyclesperpcapply = n;
1117:   return(0);
1118: }

1120: PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1121: {
1122:   PC_MG *mg = (PC_MG*)pc->data;

1125:   mg->galerkin = use;
1126:   return(0);
1127: }

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

1133:    Logically Collective on PC

1135:    Input Parameters:
1136: +  pc - the multigrid context
1137: -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE

1139:    Options Database Key:
1140: .  -pc_mg_galerkin <both,pmat,mat,none>

1142:    Level: intermediate

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

1148: .seealso: PCMGGetGalerkin(), PCMGGalerkinType

1150: @*/
1151: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1152: {

1157:   PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));
1158:   return(0);
1159: }

1161: /*@
1162:    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1163:       A_i-1 = r_i * A_i * p_i

1165:    Not Collective

1167:    Input Parameter:
1168: .  pc - the multigrid context

1170:    Output Parameter:
1171: .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL

1173:    Level: intermediate

1175: .seealso: PCMGSetGalerkin(), PCMGGalerkinType

1177: @*/
1178: PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1179: {
1180:   PC_MG *mg = (PC_MG*)pc->data;

1184:   *galerkin = mg->galerkin;
1185:   return(0);
1186: }

1188: PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1189: {
1190:   PC_MG *mg = (PC_MG *) pc->data;

1193:   mg->adaptInterpolation = adapt;
1194:   return(0);
1195: }

1197: PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1198: {
1199:   PC_MG *mg = (PC_MG *) pc->data;

1202:   *adapt = mg->adaptInterpolation;
1203:   return(0);
1204: }

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

1209:   Logically Collective on PC

1211:   Input Parameters:
1212: + pc    - the multigrid context
1213: - adapt - flag for adaptation of the interpolator

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

1220:   Level: intermediate

1222: .keywords: MG, set, Galerkin
1223: .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1224: @*/
1225: PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1226: {

1231:   PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));
1232:   return(0);
1233: }

1235: /*@
1236:   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.

1238:   Not collective

1240:   Input Parameter:
1241: . pc    - the multigrid context

1243:   Output Parameter:
1244: . adapt - flag for adaptation of the interpolator

1246:   Level: intermediate

1248: .keywords: MG, set, Galerkin
1249: .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1250: @*/
1251: PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1252: {

1258:   PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));
1259:   return(0);
1260: }

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

1267:    Logically Collective on PC

1269:    Input Parameters:
1270: +  mg - the multigrid context
1271: -  n - the number of smoothing steps

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

1276:    Level: advanced

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

1282: .seealso: PCMGSetDistinctSmoothUp()
1283: @*/
1284: PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1285: {
1286:   PC_MG          *mg        = (PC_MG*)pc->data;
1287:   PC_MG_Levels   **mglevels = mg->levels;
1289:   PetscInt       i,levels;

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

1297:   for (i=1; i<levels; i++) {
1298:     KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1299:     KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1300:     mg->default_smoothu = n;
1301:     mg->default_smoothd = n;
1302:   }
1303:   return(0);
1304: }

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

1310:    Logically Collective on PC

1312:    Input Parameters:
1313: .  pc - the preconditioner context

1315:    Options Database Key:
1316: .  -pc_mg_distinct_smoothup

1318:    Level: advanced

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

1324: .seealso: PCMGSetNumberSmooth()
1325: @*/
1326: PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1327: {
1328:   PC_MG          *mg        = (PC_MG*)pc->data;
1329:   PC_MG_Levels   **mglevels = mg->levels;
1331:   PetscInt       i,levels;
1332:   KSP            subksp;

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

1339:   for (i=1; i<levels; i++) {
1340:     const char *prefix = NULL;
1341:     /* make sure smoother up and down are different */
1342:     PCMGGetSmootherUp(pc,i,&subksp);
1343:     KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);
1344:     KSPSetOptionsPrefix(subksp,prefix);
1345:     KSPAppendOptionsPrefix(subksp,"up_");
1346:   }
1347:   return(0);
1348: }

1350: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1351: PetscErrorCode  PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[])
1352: {
1353:   PC_MG          *mg        = (PC_MG*)pc->data;
1354:   PC_MG_Levels   **mglevels = mg->levels;
1355:   Mat            *mat;
1356:   PetscInt       l;

1360:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1361:   PetscMalloc1(mg->nlevels,&mat);
1362:   for (l=1; l< mg->nlevels; l++) {
1363:     mat[l-1] = mglevels[l]->interpolate;
1364:     PetscObjectReference((PetscObject)mat[l-1]);
1365:   }
1366:   *num_levels = mg->nlevels;
1367:   *interpolations = mat;
1368:   return(0);
1369: }

1371: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1372: PetscErrorCode  PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
1373: {
1374:   PC_MG          *mg        = (PC_MG*)pc->data;
1375:   PC_MG_Levels   **mglevels = mg->levels;
1376:   PetscInt       l;
1377:   Mat            *mat;

1381:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1382:   PetscMalloc1(mg->nlevels,&mat);
1383:   for (l=0; l<mg->nlevels-1; l++) {
1384:     KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));
1385:     PetscObjectReference((PetscObject)mat[l]);
1386:   }
1387:   *num_levels = mg->nlevels;
1388:   *coarseOperators = mat;
1389:   return(0);
1390: }

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

1395:   Not collective

1397:   Input Parameters:
1398: + name     - name of the constructor
1399: - function - constructor routine

1401:   Notes:
1402:   Calling sequence for the routine:
1403: $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1404: $   pc        - The PC object
1405: $   l         - The multigrid level, 0 is the coarse level
1406: $   dm        - The DM for this level
1407: $   smooth    - The level smoother
1408: $   Nc        - The size of the coarse space
1409: $   initGuess - Basis for an initial guess for the space
1410: $   coarseSp  - A basis for the computed coarse space

1412:   Level: advanced

1414: .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister()
1415: @*/
1416: PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1417: {

1421:   PCInitializePackage();
1422:   PetscFunctionListAdd(&PCMGCoarseList,name,function);
1423:   return(0);
1424: }

1426: /*@C
1427:   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.

1429:   Not collective

1431:   Input Parameter:
1432: . name     - name of the constructor

1434:   Output Parameter:
1435: . function - constructor routine

1437:   Notes:
1438:   Calling sequence for the routine:
1439: $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1440: $   pc        - The PC object
1441: $   l         - The multigrid level, 0 is the coarse level
1442: $   dm        - The DM for this level
1443: $   smooth    - The level smoother
1444: $   Nc        - The size of the coarse space
1445: $   initGuess - Basis for an initial guess for the space
1446: $   coarseSp  - A basis for the computed coarse space

1448:   Level: advanced

1450: .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister()
1451: @*/
1452: PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1453: {

1457:   PetscFunctionListFind(PCMGCoarseList,name,function);
1458:   return(0);
1459: }

1461: /* ----------------------------------------------------------------------------------------*/

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

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

1480:    Notes:
1481:     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

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

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

1490:    Level: intermediate

1492: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1493:            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1494:            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1495:            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1496:            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1497: M*/

1499: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1500: {
1501:   PC_MG          *mg;

1505:   PetscNewLog(pc,&mg);
1506:   pc->data     = (void*)mg;
1507:   mg->nlevels  = -1;
1508:   mg->am       = PC_MG_MULTIPLICATIVE;
1509:   mg->galerkin = PC_MG_GALERKIN_NONE;
1510:   mg->adaptInterpolation = PETSC_FALSE;
1511:   mg->Nc                 = -1;
1512:   mg->eigenvalue         = -1;

1514:   pc->useAmat = PETSC_TRUE;

1516:   pc->ops->apply          = PCApply_MG;
1517:   pc->ops->setup          = PCSetUp_MG;
1518:   pc->ops->reset          = PCReset_MG;
1519:   pc->ops->destroy        = PCDestroy_MG;
1520:   pc->ops->setfromoptions = PCSetFromOptions_MG;
1521:   pc->ops->view           = PCView_MG;

1523:   PetscObjectComposedDataRegister(&mg->eigenvalue);
1524:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1525:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);
1526:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);
1527:   PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);
1528:   PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);
1529:   PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);
1530:   PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);
1531:   return(0);
1532: }