Actual source code: mgfunc.c

petsc-3.12.5 2020-03-29
Report Typos and Errors

  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: /*@
 32:    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.

 34:    Not Collective

 36:    Input Parameter:
 37: .  pc - the multigrid context

 39:    Output Parameter:
 40: .  ksp - the coarse grid solver context

 42:    Level: advanced

 44: .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetSmoother()
 45: @*/
 46: PetscErrorCode  PCMGGetCoarseSolve(PC pc,KSP *ksp)
 47: {
 48:   PC_MG        *mg        = (PC_MG*)pc->data;
 49:   PC_MG_Levels **mglevels = mg->levels;

 53:   *ksp =  mglevels[0]->smoothd;
 54:   return(0);
 55: }

 57: /*@C
 58:    PCMGSetResidual - Sets the function to be used to calculate the residual
 59:    on the lth level.

 61:    Logically Collective on PC

 63:    Input Parameters:
 64: +  pc       - the multigrid context
 65: .  l        - the level (0 is coarsest) to supply
 66: .  residual - function used to form residual, if none is provided the previously provide one is used, if no
 67:               previous one were provided then a default is used
 68: -  mat      - matrix associated with residual

 70:    Level: advanced

 72: .seealso: PCMGResidualDefault()
 73: @*/
 74: PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
 75: {
 76:   PC_MG          *mg        = (PC_MG*)pc->data;
 77:   PC_MG_Levels   **mglevels = mg->levels;

 82:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
 83:   if (residual) mglevels[l]->residual = residual;
 84:   if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault;
 85:   if (mat) {PetscObjectReference((PetscObject)mat);}
 86:   MatDestroy(&mglevels[l]->A);
 87:   mglevels[l]->A = mat;
 88:   return(0);
 89: }

 91: /*@
 92:    PCMGSetInterpolation - Sets the function to be used to calculate the
 93:    interpolation from l-1 to the lth level

 95:    Logically Collective on PC

 97:    Input Parameters:
 98: +  pc  - the multigrid context
 99: .  mat - the interpolation operator
100: -  l   - the level (0 is coarsest) to supply [do not supply 0]

102:    Level: advanced

104:    Notes:
105:           Usually this is the same matrix used also to set the restriction
106:     for the same level.

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

111: .seealso: PCMGSetRestriction()
112: @*/
113: PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
114: {
115:   PC_MG          *mg        = (PC_MG*)pc->data;
116:   PC_MG_Levels   **mglevels = mg->levels;

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

126:   mglevels[l]->interpolate = mat;
127:   return(0);
128: }

130: /*@
131:    PCMGSetOperators - Sets operator and preconditioning matrix for lth level

133:    Logically Collective on PC

135:    Input Parameters:
136: +  pc  - the multigrid context
137: .  Amat - the operator
138: .  pmat - the preconditioning operator
139: -  l   - the level (0 is the coarsest) to supply

141:    Level: advanced

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

145: .seealso: PCMGSetRestriction(), PCMGSetInterpolation()
146: @*/
147: PetscErrorCode  PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat)
148: {
149:   PC_MG          *mg        = (PC_MG*)pc->data;
150:   PC_MG_Levels   **mglevels = mg->levels;

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

162: /*@
163:    PCMGGetInterpolation - Gets the function to be used to calculate the
164:    interpolation from l-1 to the lth level

166:    Logically Collective on PC

168:    Input Parameters:
169: +  pc - the multigrid context
170: -  l - the level (0 is coarsest) to supply [Do not supply 0]

172:    Output Parameter:
173: .  mat - the interpolation matrix, can be NULL

175:    Level: advanced

177: .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale()
178: @*/
179: PetscErrorCode  PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat)
180: {
181:   PC_MG          *mg        = (PC_MG*)pc->data;
182:   PC_MG_Levels   **mglevels = mg->levels;

188:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
189:   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);
190:   if (!mglevels[l]->interpolate) {
191:     if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()");
192:     PCMGSetInterpolation(pc,l,mglevels[l]->restrct);
193:   }
194:   if (mat) *mat = mglevels[l]->interpolate;
195:   return(0);
196: }

198: /*@
199:    PCMGSetRestriction - Sets the function to be used to restrict dual vectors
200:    from level l to l-1.

202:    Logically Collective on PC

204:    Input Parameters:
205: +  pc - the multigrid context
206: .  l - the level (0 is coarsest) to supply [Do not supply 0]
207: -  mat - the restriction matrix

209:    Level: advanced

211:    Notes:
212:           Usually this is the same matrix used also to set the interpolation
213:     for the same level.

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

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

221: .seealso: PCMGSetInterpolation()
222: @*/
223: PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
224: {
226:   PC_MG          *mg        = (PC_MG*)pc->data;
227:   PC_MG_Levels   **mglevels = mg->levels;

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

237:   mglevels[l]->restrct = mat;
238:   return(0);
239: }

241: /*@
242:    PCMGGetRestriction - Gets the function to be used to restrict dual vectors
243:    from level l to l-1.

245:    Logically Collective on PC

247:    Input Parameters:
248: +  pc - the multigrid context
249: -  l - the level (0 is coarsest) to supply [Do not supply 0]

251:    Output Parameter:
252: .  mat - the restriction matrix

254:    Level: advanced

256: .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection()
257: @*/
258: PetscErrorCode  PCMGGetRestriction(PC pc,PetscInt l,Mat *mat)
259: {
260:   PC_MG          *mg        = (PC_MG*)pc->data;
261:   PC_MG_Levels   **mglevels = mg->levels;

267:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
268:   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);
269:   if (!mglevels[l]->restrct) {
270:     if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()");
271:     PCMGSetRestriction(pc,l,mglevels[l]->interpolate);
272:   }
273:   if (mat) *mat = mglevels[l]->restrct;
274:   return(0);
275: }

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

280:    Logically Collective on PC

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

287:    Level: advanced

289:    Notes:
290:        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.

292: .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection()
293: @*/
294: PetscErrorCode  PCMGSetRScale(PC pc,PetscInt l,Vec rscale)
295: {
297:   PC_MG          *mg        = (PC_MG*)pc->data;
298:   PC_MG_Levels   **mglevels = mg->levels;

302:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
303:   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);
304:   PetscObjectReference((PetscObject)rscale);
305:   VecDestroy(&mglevels[l]->rscale);

307:   mglevels[l]->rscale = rscale;
308:   return(0);
309: }

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

314:    Collective on PC

316:    Input Parameters:
317: +  pc - the multigrid context
318: .  rscale - the scaling
319: -  l - the level (0 is coarsest) to supply [Do not supply 0]

321:    Level: advanced

323:    Notes:
324:        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.

326: .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection()
327: @*/
328: PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale)
329: {
331:   PC_MG          *mg        = (PC_MG*)pc->data;
332:   PC_MG_Levels   **mglevels = mg->levels;

336:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
337:   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);
338:   if (!mglevels[l]->rscale) {
339:     Mat      R;
340:     Vec      X,Y,coarse,fine;
341:     PetscInt M,N;
342:     PCMGGetRestriction(pc,l,&R);
343:     MatCreateVecs(R,&X,&Y);
344:     MatGetSize(R,&M,&N);
345:     if (M < N) {
346:       fine = X;
347:       coarse = Y;
348:     } else if (N < M) {
349:       fine = Y; coarse = X;
350:     } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser");
351:     VecSet(fine,1.);
352:     MatRestrict(R,fine,coarse);
353:     VecDestroy(&fine);
354:     VecReciprocal(coarse);
355:     mglevels[l]->rscale = coarse;
356:   }
357:   *rscale = mglevels[l]->rscale;
358:   return(0);
359: }

361: /*@
362:    PCMGSetInjection - Sets the function to be used to inject primal vectors
363:    from level l to l-1.

365:    Logically Collective on PC

367:    Input Parameters:
368: +  pc - the multigrid context
369: .  l - the level (0 is coarsest) to supply [Do not supply 0]
370: -  mat - the injection matrix

372:    Level: advanced

374: .seealso: PCMGSetRestriction()
375: @*/
376: PetscErrorCode  PCMGSetInjection(PC pc,PetscInt l,Mat mat)
377: {
379:   PC_MG          *mg        = (PC_MG*)pc->data;
380:   PC_MG_Levels   **mglevels = mg->levels;

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

390:   mglevels[l]->inject = mat;
391:   return(0);
392: }

394: /*@
395:    PCMGGetInjection - Gets the function to be used to inject primal vectors
396:    from level l to l-1.

398:    Logically Collective on PC

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

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

407:    Level: advanced

409: .seealso: PCMGSetInjection(), PCMGetGetRestriction()
410: @*/
411: PetscErrorCode  PCMGGetInjection(PC pc,PetscInt l,Mat *mat)
412: {
413:   PC_MG          *mg        = (PC_MG*)pc->data;
414:   PC_MG_Levels   **mglevels = mg->levels;

419:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
420:   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);
421:   if (mat) *mat = mglevels[l]->inject;
422:   return(0);
423: }

425: /*@
426:    PCMGGetSmoother - Gets the KSP context to be used as smoother for
427:    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
428:    PCMGGetSmootherDown() to use different functions for pre- and
429:    post-smoothing.

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

433:    Input Parameters:
434: +  pc - the multigrid context
435: -  l - the level (0 is coarsest) to supply

437:    Ouput Parameters:
438: .  ksp - the smoother

440:    Notes:
441:    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.
442:    You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc)
443:    and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing.

445:    Level: advanced

447: .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve()
448: @*/
449: PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
450: {
451:   PC_MG        *mg        = (PC_MG*)pc->data;
452:   PC_MG_Levels **mglevels = mg->levels;

456:   *ksp = mglevels[l]->smoothd;
457:   return(0);
458: }

460: /*@
461:    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
462:    coarse grid correction (post-smoother).

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

466:    Input Parameters:
467: +  pc - the multigrid context
468: -  l  - the level (0 is coarsest) to supply

470:    Ouput Parameters:
471: .  ksp - the smoother

473:    Level: advanced

475:    Notes:
476:     calling this will result in a different pre and post smoother so you may need to
477:          set options on the pre smoother also

479: .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
480: @*/
481: PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
482: {
483:   PC_MG          *mg        = (PC_MG*)pc->data;
484:   PC_MG_Levels   **mglevels = mg->levels;
486:   const char     *prefix;
487:   MPI_Comm       comm;

491:   /*
492:      This is called only if user wants a different pre-smoother from post.
493:      Thus we check if a different one has already been allocated,
494:      if not we allocate it.
495:   */
496:   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
497:   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
498:     KSPType     ksptype;
499:     PCType      pctype;
500:     PC          ipc;
501:     PetscReal   rtol,abstol,dtol;
502:     PetscInt    maxits;
503:     KSPNormType normtype;
504:     PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);
505:     KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);
506:     KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);
507:     KSPGetType(mglevels[l]->smoothd,&ksptype);
508:     KSPGetNormType(mglevels[l]->smoothd,&normtype);
509:     KSPGetPC(mglevels[l]->smoothd,&ipc);
510:     PCGetType(ipc,&pctype);

512:     KSPCreate(comm,&mglevels[l]->smoothu);
513:     KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);
514:     PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);
515:     KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);
516:     KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);
517:     KSPSetType(mglevels[l]->smoothu,ksptype);
518:     KSPSetNormType(mglevels[l]->smoothu,normtype);
519:     KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);
520:     KSPGetPC(mglevels[l]->smoothu,&ipc);
521:     PCSetType(ipc,pctype);
522:     PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);
523:     PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);
524:   }
525:   if (ksp) *ksp = mglevels[l]->smoothu;
526:   return(0);
527: }

529: /*@
530:    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
531:    coarse grid correction (pre-smoother).

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

535:    Input Parameters:
536: +  pc - the multigrid context
537: -  l  - the level (0 is coarsest) to supply

539:    Ouput Parameters:
540: .  ksp - the smoother

542:    Level: advanced

544:    Notes:
545:     calling this will result in a different pre and post smoother so you may need to
546:          set options on the post smoother also

548: .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
549: @*/
550: PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
551: {
553:   PC_MG          *mg        = (PC_MG*)pc->data;
554:   PC_MG_Levels   **mglevels = mg->levels;

558:   /* make sure smoother up and down are different */
559:   if (l) {
560:     PCMGGetSmootherUp(pc,l,NULL);
561:   }
562:   *ksp = mglevels[l]->smoothd;
563:   return(0);
564: }

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

569:    Logically Collective on PC

571:    Input Parameters:
572: +  pc - the multigrid context
573: .  l  - the level (0 is coarsest)
574: -  c  - either PC_MG_CYCLE_V or PC_MG_CYCLE_W

576:    Level: advanced

578: .seealso: PCMGSetCycleType()
579: @*/
580: PetscErrorCode  PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c)
581: {
582:   PC_MG        *mg        = (PC_MG*)pc->data;
583:   PC_MG_Levels **mglevels = mg->levels;

587:   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
590:   mglevels[l]->cycles = c;
591:   return(0);
592: }

594: /*@
595:    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
596:    on a particular level.

598:    Logically Collective on PC

600:    Input Parameters:
601: +  pc - the multigrid context
602: .  l  - the level (0 is coarsest) this is to be used for
603: -  c  - the space

605:    Level: advanced

607:    Notes:
608:     If this is not provided PETSc will automatically generate one.

610:           You do not need to keep a reference to this vector if you do
611:           not need it PCDestroy() will properly free it.

613: .seealso: PCMGSetX(), PCMGSetR()
614: @*/
615: PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
616: {
618:   PC_MG          *mg        = (PC_MG*)pc->data;
619:   PC_MG_Levels   **mglevels = mg->levels;

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

628:   mglevels[l]->b = c;
629:   return(0);
630: }

632: /*@
633:    PCMGSetX - Sets the vector space to be used to store the solution on a
634:    particular level.

636:    Logically Collective on PC

638:    Input Parameters:
639: +  pc - the multigrid context
640: .  l - the level (0 is coarsest) this is to be used for (do not supply the finest level)
641: -  c - the space

643:    Level: advanced

645:    Notes:
646:     If this is not provided PETSc will automatically generate one.

648:           You do not need to keep a reference to this vector if you do
649:           not need it PCDestroy() will properly free it.

651: .seealso: PCMGSetRhs(), PCMGSetR()
652: @*/
653: PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
654: {
656:   PC_MG          *mg        = (PC_MG*)pc->data;
657:   PC_MG_Levels   **mglevels = mg->levels;

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

666:   mglevels[l]->x = c;
667:   return(0);
668: }

670: /*@
671:    PCMGSetR - Sets the vector space to be used to store the residual on a
672:    particular level.

674:    Logically Collective on PC

676:    Input Parameters:
677: +  pc - the multigrid context
678: .  l - the level (0 is coarsest) this is to be used for
679: -  c - the space

681:    Level: advanced

683:    Notes:
684:     If this is not provided PETSc will automatically generate one.

686:           You do not need to keep a reference to this vector if you do
687:           not need it PCDestroy() will properly free it.

689: @*/
690: PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
691: {
693:   PC_MG          *mg        = (PC_MG*)pc->data;
694:   PC_MG_Levels   **mglevels = mg->levels;

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

703:   mglevels[l]->r = c;
704:   return(0);
705: }