Actual source code: mgfunc.c


  2: #include <petsc/private/pcmgimpl.h>

  4: /* ---------------------------------------------------------------------------*/
  5: /*@C
  6:    PCMGResidualDefault - Default routine to calculate the residual.

  8:    Collective on Mat

 10:    Input Parameters:
 11: +  mat - the matrix
 12: .  b   - the right-hand-side
 13: -  x   - the approximate solution

 15:    Output Parameter:
 16: .  r - location to store the residual

 18:    Level: developer

 20: .seealso: PCMGSetResidual()
 21: @*/
 22: PetscErrorCode  PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r)
 23: {

 27:   MatResidual(mat,b,x,r);
 28:   return(0);
 29: }

 31: /*@C
 32:    PCMGResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system

 34:    Collective on Mat

 36:    Input Parameters:
 37: +  mat - the matrix
 38: .  b   - the right-hand-side
 39: -  x   - the approximate solution

 41:    Output Parameter:
 42: .  r - location to store the residual

 44:    Level: developer

 46: .seealso: PCMGSetResidualTranspose()
 47: @*/
 48: PetscErrorCode PCMGResidualTransposeDefault(Mat mat,Vec b,Vec x,Vec r)
 49: {

 53:   MatMultTranspose(mat,x,r);
 54:   VecAYPX(r,-1.0,b);
 55:   return(0);
 56: }

 58: /*@C
 59:    PCMGMatResidualDefault - Default routine to calculate the residual.

 61:    Collective on Mat

 63:    Input Parameters:
 64: +  mat - the matrix
 65: .  b   - the right-hand-side
 66: -  x   - the approximate solution

 68:    Output Parameter:
 69: .  r - location to store the residual

 71:    Level: developer

 73: .seealso: PCMGSetMatResidual()
 74: @*/
 75: PetscErrorCode  PCMGMatResidualDefault(Mat mat,Mat b,Mat x,Mat r)
 76: {

 80:   MatMatMult(mat,x,MAT_REUSE_MATRIX,PETSC_DEFAULT,&r);
 81:   MatAYPX(r,-1.0,b,UNKNOWN_NONZERO_PATTERN);
 82:   return(0);
 83: }

 85: /*@C
 86:    PCMGMatResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system

 88:    Collective on Mat

 90:    Input Parameters:
 91: +  mat - the matrix
 92: .  b   - the right-hand-side
 93: -  x   - the approximate solution

 95:    Output Parameter:
 96: .  r - location to store the residual

 98:    Level: developer

100: .seealso: PCMGSetMatResidualTranspose()
101: @*/
102: PetscErrorCode PCMGMatResidualTransposeDefault(Mat mat,Mat b,Mat x,Mat r)
103: {

107:   MatTransposeMatMult(mat,x,MAT_REUSE_MATRIX,PETSC_DEFAULT,&r);
108:   MatAYPX(r,-1.0,b,UNKNOWN_NONZERO_PATTERN);
109:   return(0);
110: }
111: /*@
112:    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.

114:    Not Collective

116:    Input Parameter:
117: .  pc - the multigrid context

119:    Output Parameter:
120: .  ksp - the coarse grid solver context

122:    Level: advanced

124: .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetSmoother()
125: @*/
126: PetscErrorCode  PCMGGetCoarseSolve(PC pc,KSP *ksp)
127: {
128:   PC_MG        *mg        = (PC_MG*)pc->data;
129:   PC_MG_Levels **mglevels = mg->levels;

133:   *ksp =  mglevels[0]->smoothd;
134:   return(0);
135: }

137: /*@C
138:    PCMGSetResidual - Sets the function to be used to calculate the residual
139:    on the lth level.

141:    Logically Collective on PC

143:    Input Parameters:
144: +  pc       - the multigrid context
145: .  l        - the level (0 is coarsest) to supply
146: .  residual - function used to form residual, if none is provided the previously provide one is used, if no
147:               previous one were provided then a default is used
148: -  mat      - matrix associated with residual

150:    Level: advanced

152: .seealso: PCMGResidualDefault()
153: @*/
154: PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
155: {
156:   PC_MG          *mg        = (PC_MG*)pc->data;
157:   PC_MG_Levels   **mglevels = mg->levels;

162:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
163:   if (residual) mglevels[l]->residual = residual;
164:   if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault;
165:   mglevels[l]->matresidual = PCMGMatResidualDefault;
166:   if (mat) {PetscObjectReference((PetscObject)mat);}
167:   MatDestroy(&mglevels[l]->A);
168:   mglevels[l]->A = mat;
169:   return(0);
170: }

172: /*@C
173:    PCMGSetResidualTranspose - Sets the function to be used to calculate the residual of the transposed linear system
174:    on the lth level.

176:    Logically Collective on PC

178:    Input Parameters:
179: +  pc        - the multigrid context
180: .  l         - the level (0 is coarsest) to supply
181: .  residualt - function used to form transpose of residual, if none is provided the previously provide one is used, if no
182:                previous one were provided then a default is used
183: -  mat       - matrix associated with residual

185:    Level: advanced

187: .seealso: PCMGResidualTransposeDefault()
188: @*/
189: PetscErrorCode  PCMGSetResidualTranspose(PC pc,PetscInt l,PetscErrorCode (*residualt)(Mat,Vec,Vec,Vec),Mat mat)
190: {
191:   PC_MG          *mg        = (PC_MG*)pc->data;
192:   PC_MG_Levels   **mglevels = mg->levels;

197:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
198:   if (residualt) mglevels[l]->residualtranspose = residualt;
199:   if (!mglevels[l]->residualtranspose) mglevels[l]->residualtranspose = PCMGResidualTransposeDefault;
200:   mglevels[l]->matresidualtranspose = PCMGMatResidualTransposeDefault;
201:   if (mat) {PetscObjectReference((PetscObject)mat);}
202:   MatDestroy(&mglevels[l]->A);
203:   mglevels[l]->A = mat;
204:   return(0);
205: }

207: /*@
208:    PCMGSetInterpolation - Sets the function to be used to calculate the
209:    interpolation from l-1 to the lth level

211:    Logically Collective on PC

213:    Input Parameters:
214: +  pc  - the multigrid context
215: .  mat - the interpolation operator
216: -  l   - the level (0 is coarsest) to supply [do not supply 0]

218:    Level: advanced

220:    Notes:
221:           Usually this is the same matrix used also to set the restriction
222:     for the same level.

224:           One can pass in the interpolation matrix or its transpose; PETSc figures
225:     out from the matrix size which one it is.

227: .seealso: PCMGSetRestriction()
228: @*/
229: PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
230: {
231:   PC_MG          *mg        = (PC_MG*)pc->data;
232:   PC_MG_Levels   **mglevels = mg->levels;

237:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
238:   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
239:   PetscObjectReference((PetscObject)mat);
240:   MatDestroy(&mglevels[l]->interpolate);

242:   mglevels[l]->interpolate = mat;
243:   return(0);
244: }

246: /*@
247:    PCMGSetOperators - Sets operator and preconditioning matrix for lth level

249:    Logically Collective on PC

251:    Input Parameters:
252: +  pc  - the multigrid context
253: .  Amat - the operator
254: .  pmat - the preconditioning operator
255: -  l   - the level (0 is the coarsest) to supply

257:    Level: advanced

259: .keywords:  multigrid, set, interpolate, level

261: .seealso: PCMGSetRestriction(), PCMGSetInterpolation()
262: @*/
263: PetscErrorCode  PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat)
264: {
265:   PC_MG          *mg        = (PC_MG*)pc->data;
266:   PC_MG_Levels   **mglevels = mg->levels;

273:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
274:   KSPSetOperators(mglevels[l]->smoothd,Amat,Pmat);
275:   return(0);
276: }

278: /*@
279:    PCMGGetInterpolation - Gets the function to be used to calculate the
280:    interpolation from l-1 to the lth level

282:    Logically Collective on PC

284:    Input Parameters:
285: +  pc - the multigrid context
286: -  l - the level (0 is coarsest) to supply [Do not supply 0]

288:    Output Parameter:
289: .  mat - the interpolation matrix, can be NULL

291:    Level: advanced

293: .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale()
294: @*/
295: PetscErrorCode  PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat)
296: {
297:   PC_MG          *mg        = (PC_MG*)pc->data;
298:   PC_MG_Levels   **mglevels = mg->levels;

304:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
305:   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
306:   if (!mglevels[l]->interpolate) {
307:     if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()");
308:     PCMGSetInterpolation(pc,l,mglevels[l]->restrct);
309:   }
310:   if (mat) *mat = mglevels[l]->interpolate;
311:   return(0);
312: }

314: /*@
315:    PCMGSetRestriction - Sets the function to be used to restrict dual vectors
316:    from level l to l-1.

318:    Logically Collective on PC

320:    Input Parameters:
321: +  pc - the multigrid context
322: .  l - the level (0 is coarsest) to supply [Do not supply 0]
323: -  mat - the restriction matrix

325:    Level: advanced

327:    Notes:
328:           Usually this is the same matrix used also to set the interpolation
329:     for the same level.

331:           One can pass in the interpolation matrix or its transpose; PETSc figures
332:     out from the matrix size which one it is.

334:          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
335:     is used.

337: .seealso: PCMGSetInterpolation()
338: @*/
339: PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
340: {
342:   PC_MG          *mg        = (PC_MG*)pc->data;
343:   PC_MG_Levels   **mglevels = mg->levels;

348:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
349:   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
350:   PetscObjectReference((PetscObject)mat);
351:   MatDestroy(&mglevels[l]->restrct);

353:   mglevels[l]->restrct = mat;
354:   return(0);
355: }

357: /*@
358:    PCMGGetRestriction - Gets the function to be used to restrict dual vectors
359:    from level l to l-1.

361:    Logically Collective on PC

363:    Input Parameters:
364: +  pc - the multigrid context
365: -  l - the level (0 is coarsest) to supply [Do not supply 0]

367:    Output Parameter:
368: .  mat - the restriction matrix

370:    Level: advanced

372: .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection()
373: @*/
374: PetscErrorCode  PCMGGetRestriction(PC pc,PetscInt l,Mat *mat)
375: {
376:   PC_MG          *mg        = (PC_MG*)pc->data;
377:   PC_MG_Levels   **mglevels = mg->levels;

383:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
384:   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
385:   if (!mglevels[l]->restrct) {
386:     if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()");
387:     PCMGSetRestriction(pc,l,mglevels[l]->interpolate);
388:   }
389:   if (mat) *mat = mglevels[l]->restrct;
390:   return(0);
391: }

393: /*@
394:    PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.

396:    Logically Collective on PC

398:    Input Parameters:
399: +  pc - the multigrid context
400: -  l - the level (0 is coarsest) to supply [Do not supply 0]
401: .  rscale - the scaling

403:    Level: advanced

405:    Notes:
406:        When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R.  It is preferable to use PCMGSetInjection() to control moving primal vectors.

408: .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection()
409: @*/
410: PetscErrorCode  PCMGSetRScale(PC pc,PetscInt l,Vec rscale)
411: {
413:   PC_MG          *mg        = (PC_MG*)pc->data;
414:   PC_MG_Levels   **mglevels = mg->levels;

418:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
419:   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
420:   PetscObjectReference((PetscObject)rscale);
421:   VecDestroy(&mglevels[l]->rscale);

423:   mglevels[l]->rscale = rscale;
424:   return(0);
425: }

427: /*@
428:    PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1.

430:    Collective on PC

432:    Input Parameters:
433: +  pc - the multigrid context
434: .  rscale - the scaling
435: -  l - the level (0 is coarsest) to supply [Do not supply 0]

437:    Level: advanced

439:    Notes:
440:        When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R.  It is preferable to use PCMGGetInjection() to control moving primal vectors.

442: .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection()
443: @*/
444: PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale)
445: {
447:   PC_MG          *mg        = (PC_MG*)pc->data;
448:   PC_MG_Levels   **mglevels = mg->levels;

452:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
453:   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
454:   if (!mglevels[l]->rscale) {
455:     Mat      R;
456:     Vec      X,Y,coarse,fine;
457:     PetscInt M,N;
458:     PCMGGetRestriction(pc,l,&R);
459:     MatCreateVecs(R,&X,&Y);
460:     MatGetSize(R,&M,&N);
461:     if (M < N) {
462:       fine = X;
463:       coarse = Y;
464:     } else if (N < M) {
465:       fine = Y; coarse = X;
466:     } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser");
467:     VecSet(fine,1.);
468:     MatRestrict(R,fine,coarse);
469:     VecDestroy(&fine);
470:     VecReciprocal(coarse);
471:     mglevels[l]->rscale = coarse;
472:   }
473:   *rscale = mglevels[l]->rscale;
474:   return(0);
475: }

477: /*@
478:    PCMGSetInjection - Sets the function to be used to inject primal vectors
479:    from level l to l-1.

481:    Logically Collective on PC

483:    Input Parameters:
484: +  pc - the multigrid context
485: .  l - the level (0 is coarsest) to supply [Do not supply 0]
486: -  mat - the injection matrix

488:    Level: advanced

490: .seealso: PCMGSetRestriction()
491: @*/
492: PetscErrorCode  PCMGSetInjection(PC pc,PetscInt l,Mat mat)
493: {
495:   PC_MG          *mg        = (PC_MG*)pc->data;
496:   PC_MG_Levels   **mglevels = mg->levels;

501:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
502:   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
503:   PetscObjectReference((PetscObject)mat);
504:   MatDestroy(&mglevels[l]->inject);

506:   mglevels[l]->inject = mat;
507:   return(0);
508: }

510: /*@
511:    PCMGGetInjection - Gets the function to be used to inject primal vectors
512:    from level l to l-1.

514:    Logically Collective on PC

516:    Input Parameters:
517: +  pc - the multigrid context
518: -  l - the level (0 is coarsest) to supply [Do not supply 0]

520:    Output Parameter:
521: .  mat - the restriction matrix (may be NULL if no injection is available).

523:    Level: advanced

525: .seealso: PCMGSetInjection(), PCMGetGetRestriction()
526: @*/
527: PetscErrorCode  PCMGGetInjection(PC pc,PetscInt l,Mat *mat)
528: {
529:   PC_MG          *mg        = (PC_MG*)pc->data;
530:   PC_MG_Levels   **mglevels = mg->levels;

535:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
536:   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
537:   if (mat) *mat = mglevels[l]->inject;
538:   return(0);
539: }

541: /*@
542:    PCMGGetSmoother - Gets the KSP context to be used as smoother for
543:    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
544:    PCMGGetSmootherDown() to use different functions for pre- and
545:    post-smoothing.

547:    Not Collective, KSP returned is parallel if PC is

549:    Input Parameters:
550: +  pc - the multigrid context
551: -  l - the level (0 is coarsest) to supply

553:    Output Parameter:
554: .  ksp - the smoother

556:    Notes:
557:    Once you have called this routine, you can call KSPSetOperators(ksp,...) on the resulting ksp to provide the operators for the smoother for this level.
558:    You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc)
559:    and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing.

561:    Level: advanced

563: .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve()
564: @*/
565: PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
566: {
567:   PC_MG        *mg        = (PC_MG*)pc->data;
568:   PC_MG_Levels **mglevels = mg->levels;

572:   *ksp = mglevels[l]->smoothd;
573:   return(0);
574: }

576: /*@
577:    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
578:    coarse grid correction (post-smoother).

580:    Not Collective, KSP returned is parallel if PC is

582:    Input Parameters:
583: +  pc - the multigrid context
584: -  l  - the level (0 is coarsest) to supply

586:    Output Parameter:
587: .  ksp - the smoother

589:    Level: advanced

591:    Notes:
592:     calling this will result in a different pre and post smoother so you may need to
593:          set options on the pre smoother also

595: .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
596: @*/
597: PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
598: {
599:   PC_MG          *mg        = (PC_MG*)pc->data;
600:   PC_MG_Levels   **mglevels = mg->levels;
602:   const char     *prefix;
603:   MPI_Comm       comm;

607:   /*
608:      This is called only if user wants a different pre-smoother from post.
609:      Thus we check if a different one has already been allocated,
610:      if not we allocate it.
611:   */
612:   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
613:   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
614:     KSPType     ksptype;
615:     PCType      pctype;
616:     PC          ipc;
617:     PetscReal   rtol,abstol,dtol;
618:     PetscInt    maxits;
619:     KSPNormType normtype;
620:     PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);
621:     KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);
622:     KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);
623:     KSPGetType(mglevels[l]->smoothd,&ksptype);
624:     KSPGetNormType(mglevels[l]->smoothd,&normtype);
625:     KSPGetPC(mglevels[l]->smoothd,&ipc);
626:     PCGetType(ipc,&pctype);

628:     KSPCreate(comm,&mglevels[l]->smoothu);
629:     KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);
630:     PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);
631:     KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);
632:     KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);
633:     KSPSetType(mglevels[l]->smoothu,ksptype);
634:     KSPSetNormType(mglevels[l]->smoothu,normtype);
635:     KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);
636:     KSPGetPC(mglevels[l]->smoothu,&ipc);
637:     PCSetType(ipc,pctype);
638:     PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);
639:     PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);
640:   }
641:   if (ksp) *ksp = mglevels[l]->smoothu;
642:   return(0);
643: }

645: /*@
646:    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
647:    coarse grid correction (pre-smoother).

649:    Not Collective, KSP returned is parallel if PC is

651:    Input Parameters:
652: +  pc - the multigrid context
653: -  l  - the level (0 is coarsest) to supply

655:    Output Parameter:
656: .  ksp - the smoother

658:    Level: advanced

660:    Notes:
661:     calling this will result in a different pre and post smoother so you may need to
662:          set options on the post smoother also

664: .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
665: @*/
666: PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
667: {
669:   PC_MG          *mg        = (PC_MG*)pc->data;
670:   PC_MG_Levels   **mglevels = mg->levels;

674:   /* make sure smoother up and down are different */
675:   if (l) {
676:     PCMGGetSmootherUp(pc,l,NULL);
677:   }
678:   *ksp = mglevels[l]->smoothd;
679:   return(0);
680: }

682: /*@
683:    PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level.

685:    Logically Collective on PC

687:    Input Parameters:
688: +  pc - the multigrid context
689: .  l  - the level (0 is coarsest)
690: -  c  - either PC_MG_CYCLE_V or PC_MG_CYCLE_W

692:    Level: advanced

694: .seealso: PCMGSetCycleType()
695: @*/
696: PetscErrorCode  PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c)
697: {
698:   PC_MG        *mg        = (PC_MG*)pc->data;
699:   PC_MG_Levels **mglevels = mg->levels;

703:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
706:   mglevels[l]->cycles = c;
707:   return(0);
708: }

710: /*@
711:   PCMGSetRhs - Sets the vector to be used to store the right-hand side on a particular level.

713:    Logically Collective on PC

715:   Input Parameters:
716: + pc - the multigrid context
717: . l  - the level (0 is coarsest) this is to be used for
718: - c  - the Vec

720:   Level: advanced

722:   Notes:
723:   If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it.

725: .keywords: MG, multigrid, set, right-hand-side, rhs, level
726: .seealso: PCMGSetX(), PCMGSetR()
727: @*/
728: PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
729: {
731:   PC_MG          *mg        = (PC_MG*)pc->data;
732:   PC_MG_Levels   **mglevels = mg->levels;

736:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
737:   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
738:   PetscObjectReference((PetscObject)c);
739:   VecDestroy(&mglevels[l]->b);

741:   mglevels[l]->b = c;
742:   return(0);
743: }

745: /*@
746:   PCMGSetX - Sets the vector to be used to store the solution on a particular level.

748:   Logically Collective on PC

750:   Input Parameters:
751: + pc - the multigrid context
752: . l - the level (0 is coarsest) this is to be used for (do not supply the finest level)
753: - c - the Vec

755:   Level: advanced

757:   Notes:
758:   If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it.

760: .keywords: MG, multigrid, set, solution, level
761: .seealso: PCMGSetRhs(), PCMGSetR()
762: @*/
763: PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
764: {
766:   PC_MG          *mg        = (PC_MG*)pc->data;
767:   PC_MG_Levels   **mglevels = mg->levels;

771:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
772:   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
773:   PetscObjectReference((PetscObject)c);
774:   VecDestroy(&mglevels[l]->x);

776:   mglevels[l]->x = c;
777:   return(0);
778: }

780: /*@
781:   PCMGSetR - Sets the vector to be used to store the residual on a particular level.

783:   Logically Collective on PC

785:   Input Parameters:
786: + pc - the multigrid context
787: . l - the level (0 is coarsest) this is to be used for
788: - c - the Vec

790:   Level: advanced

792:   Notes:
793:   If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it.

795: .keywords: MG, multigrid, set, residual, level
796: .seealso: PCMGSetRhs(), PCMGSetX()
797: @*/
798: PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
799: {
801:   PC_MG          *mg        = (PC_MG*)pc->data;
802:   PC_MG_Levels   **mglevels = mg->levels;

806:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
807:   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
808:   PetscObjectReference((PetscObject)c);
809:   VecDestroy(&mglevels[l]->r);

811:   mglevels[l]->r = c;
812:   return(0);
813: }