Actual source code: damg.c

  1: #define PETSCSNES_DLL
  2: 
 3:  #include petscda.h
 4:  #include petscksp.h
 5:  #include petscmg.h
 6:  #include petscdmmg.h
 7:  #include private/pcimpl.h

  9: /*
 10:    Code for almost fully managing multigrid/multi-level linear solvers for DA grids
 11: */

 15: /*@C
 16:     DMMGCreate - Creates a DA based multigrid solver object. This allows one to 
 17:       easily implement MG methods on regular grids.

 19:     Collective on MPI_Comm

 21:     Input Parameter:
 22: +   comm - the processors that will share the grids and solution process
 23: .   nlevels - number of multigrid levels 
 24: -   user - an optional user context

 26:     Output Parameters:
 27: .    - the context

 29:     Options Database:
 30: +     -dmmg_nlevels <levels> - number of levels to use
 31: .     -dmmg_galerkin - use Galerkin approach to compute coarser matrices
 32: -     -dmmg_mat_type <type> - matrix type that DMMG should create, defaults to MATAIJ

 34:     Notes:
 35:       To provide a different user context for each level call DMMGSetUser() after calling
 36:       this routine

 38:     Level: advanced

 40: .seealso DMMGDestroy(), DMMGSetUser(), DMMGGetUser(), DMMGSetMatType(), DMMGSetUseGalerkin(), DMMGSetNullSpace(), DMMGSetInitialGuess(),
 41:          DMMGSetISColoringType()

 43: @*/
 44: PetscErrorCode  DMMGCreate(MPI_Comm comm,PetscInt nlevels,void *user,DMMG **dmmg)
 45: {
 47:   PetscInt       i;
 48:   DMMG           *p;
 49:   PetscTruth     galerkin,ftype;
 50:   char           mtype[256];

 53:   if (nlevels < 0) {
 54:     nlevels = -nlevels;
 55:   } else {
 56:     PetscOptionsGetInt(0,"-dmmg_nlevels",&nlevels,PETSC_IGNORE);
 57:   }
 58:   PetscOptionsHasName(0,"-dmmg_galerkin",&galerkin);

 60:   PetscMalloc(nlevels*sizeof(DMMG),&p);
 61:   for (i=0; i<nlevels; i++) {
 62:     PetscNew(struct _n_DMMG,&p[i]);
 63:     p[i]->nlevels  = nlevels - i;
 64:     p[i]->comm     = comm;
 65:     p[i]->user     = user;
 66:     p[i]->galerkin = galerkin;
 67:     p[i]->mtype    = MATAIJ;
 68:     p[i]->isctype  = IS_COLORING_GHOSTED;   /* default to faster version, requires DMMGSetSNESLocal() */
 69:   }
 70:   p[nlevels-1]->galerkin = PETSC_FALSE;
 71:   *dmmg = p;

 73:   PetscOptionsGetString(PETSC_NULL,"-dmmg_mat_type",mtype,256,&ftype);
 74:   if (ftype) {
 75:     DMMGSetMatType(*dmmg,mtype);
 76:   }
 77:   return(0);
 78: }

 82: /*@C
 83:     DMMGSetMatType - Sets the type of matrices that DMMG will create for its solvers.

 85:     Collective on MPI_Comm 

 87:     Input Parameters:
 88: +    dmmg - the DMMG object created with DMMGCreate()
 89: -    mtype - the matrix type, defaults to MATAIJ

 91:     Level: intermediate

 93: .seealso DMMGDestroy(), DMMGSetUser(), DMMGGetUser(), DMMGCreate(), DMMGSetNullSpace()

 95: @*/
 96: PetscErrorCode  DMMGSetMatType(DMMG *dmmg,MatType mtype)
 97: {
 98:   PetscInt i;
 99: 
101:   for (i=0; i<dmmg[0]->nlevels; i++) {
102:     dmmg[i]->mtype  = mtype;
103:   }
104:   return(0);
105: }

109: /*@C
110:     DMMGSetPrefix - Sets the prefix used for the solvers inside a DMMG

112:     Collective on MPI_Comm 

114:     Input Parameters:
115: +    dmmg - the DMMG object created with DMMGCreate()
116: -    prefix - the prefix string

118:     Level: intermediate

120: .seealso DMMGDestroy(), DMMGSetUser(), DMMGGetUser(), DMMGCreate(), DMMGSetNullSpace()

122: @*/
123: PetscErrorCode  DMMGSetPrefix(DMMG *dmmg,const char* prefix)
124: {
125:   PetscInt       i;
127: 
129:   for (i=0; i<dmmg[0]->nlevels; i++) {
130:     PetscStrallocpy(prefix,&dmmg[i]->prefix);
131:   }
132:   return(0);
133: }

137: /*@C
138:     DMMGSetUseGalerkinCoarse - Courses the DMMG to use R*A_f*R^T to form
139:        the coarser matrices from finest 

141:     Collective on DMMG

143:     Input Parameter:
144: .    dmmg - the context

146:     Options Database Keys:
147: .    -dmmg_galerkin 

149:     Level: advanced

151:     Notes: After you have called this you can manually set dmmg[0]->galerkin = PETSC_FALSE
152:        to have the coarsest grid not compute via Galerkin but still have the intermediate
153:        grids computed via Galerkin.

155:        The default behavior of this should be idential to using -pc_mg_galerkin; this offers
156:        more potential flexibility since you can select exactly which levels are done via
157:        Galerkin and which are done via user provided function.

159: .seealso DMMGCreate(), PCMGSetUseGalerkin(), DMMGSetMatType(), DMMGSetNullSpace()

161: @*/
162: PetscErrorCode  DMMGSetUseGalerkinCoarse(DMMG* dmmg)
163: {
164:   PetscInt  i,nlevels = dmmg[0]->nlevels;

167:   if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");

169:   for (i=0; i<nlevels-1; i++) {
170:     dmmg[i]->galerkin = PETSC_TRUE;
171:   }
172:   return(0);
173: }

177: /*@C
178:     DMMGDestroy - Destroys a DA based multigrid solver object. 

180:     Collective on DMMG

182:     Input Parameter:
183: .    - the context

185:     Level: advanced

187: .seealso DMMGCreate()

189: @*/
190: PetscErrorCode  DMMGDestroy(DMMG *dmmg)
191: {
193:   PetscInt       i,nlevels = dmmg[0]->nlevels;

196:   if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");

198:   for (i=1; i<nlevels; i++) {
199:     if (dmmg[i]->R) {MatDestroy(dmmg[i]->R);}
200:   }
201:   for (i=0; i<nlevels; i++) {
202:     PetscStrfree(dmmg[i]->prefix);
203:     if (dmmg[i]->dm)      {DMDestroy(dmmg[i]->dm);}
204:     if (dmmg[i]->x)       {VecDestroy(dmmg[i]->x);}
205:     if (dmmg[i]->b)       {VecDestroy(dmmg[i]->b);}
206:     if (dmmg[i]->r)       {VecDestroy(dmmg[i]->r);}
207:     if (dmmg[i]->work1)   {VecDestroy(dmmg[i]->work1);}
208:     if (dmmg[i]->w)       {VecDestroy(dmmg[i]->w);}
209:     if (dmmg[i]->work2)   {VecDestroy(dmmg[i]->work2);}
210:     if (dmmg[i]->lwork1)  {VecDestroy(dmmg[i]->lwork1);}
211:     if (dmmg[i]->B)         {MatDestroy(dmmg[i]->B);}
212:     if (dmmg[i]->J)         {MatDestroy(dmmg[i]->J);}
213:     if (dmmg[i]->Rscale)    {VecDestroy(dmmg[i]->Rscale);}
214:     if (dmmg[i]->fdcoloring){MatFDColoringDestroy(dmmg[i]->fdcoloring);}
215:     if (dmmg[i]->ksp && !dmmg[i]->snes) {KSPDestroy(dmmg[i]->ksp);}
216:     if (dmmg[i]->snes)      {PetscObjectDestroy((PetscObject)dmmg[i]->snes);}
217:     if (dmmg[i]->inject)    {VecScatterDestroy(dmmg[i]->inject);}
218:     PetscFree(dmmg[i]);
219:   }
220:   PetscFree(dmmg);
221:   return(0);
222: }

226: /*@C
227:     DMMGSetDM - Sets the coarse grid information for the grids

229:     Collective on DMMG

231:     Input Parameter:
232: +   dmmg - the context
233: -   dm - the DA or DMComposite object

235:     Options Database Keys:
236: +   -dmmg_refine: Use the input problem as the coarse level and refine. Otherwise, use it as the fine level and coarsen.
237: -   -dmmg_hierarchy: Construct all grids at once

239:     Level: advanced

241: .seealso DMMGCreate(), DMMGDestroy(), DMMGSetUseGalerkin(), DMMGSetMatType()

243: @*/
244: PetscErrorCode  DMMGSetDM(DMMG *dmmg, DM dm)
245: {
246:   PetscInt       nlevels     = dmmg[0]->nlevels;
247:   PetscTruth     doRefine    = PETSC_TRUE;
248:   PetscTruth     doHierarchy = PETSC_FALSE;
249:   PetscInt       i;

253:   if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");

255:   /* Create DM data structure for all the levels */
256:   PetscOptionsGetTruth(PETSC_NULL, "-dmmg_refine", &doRefine, PETSC_IGNORE);
257:   PetscOptionsHasName(PETSC_NULL, "-dmmg_hierarchy", &doHierarchy);
258:   PetscObjectReference((PetscObject) dm);
259:   if (doRefine) {
260:     dmmg[0]->dm = dm;
261:     if (doHierarchy) {
262: /*       DM *hierarchy; */

264: /*       DMRefineHierarchy(dm, nlevels-1, &hierarchy); */
265: /*       for(i = 1; i < nlevels; ++i) { */
266: /*         dmmg[i]->dm = hierarchy[i-1]; */
267: /*       } */
268:       SETERRQ(PETSC_ERR_SUP, "Refinement hierarchy not yet implemented");
269:     } else {
270:       for(i = 1; i < nlevels; ++i) {
271:         DMRefine(dmmg[i-1]->dm, dmmg[i]->comm, &dmmg[i]->dm);
272:       }
273:     }
274:   } else {
275:     dmmg[nlevels-1]->dm = dm;
276:     if (doHierarchy) {
277:       DM *hierarchy;

279:       DMCoarsenHierarchy(dm, nlevels-1, &hierarchy);
280:       for(i = 0; i < nlevels-1; ++i) {
281:         dmmg[nlevels-2-i]->dm = hierarchy[i];
282:       }
283:     } else {
284: /*       for(i = nlevels-2; i >= 0; --i) { */
285: /*         DMCoarsen(dmmg[i+1]->dm, dmmg[i]->comm, &dmmg[i]->dm); */
286: /*       } */
287:       SETERRQ(PETSC_ERR_SUP, "Sequential coarsening not yet implemented");
288:     }
289:   }
290:   /* Cleanup old structures (should use some private Destroy() instead) */
291:   for(i = 0; i < nlevels; ++i) {
292:     if (dmmg[i]->B) {MatDestroy(dmmg[i]->B); dmmg[i]->B = PETSC_NULL;}
293:     if (dmmg[i]->J) {MatDestroy(dmmg[i]->J); dmmg[i]->J = PETSC_NULL;}
294:   }
295:   DMMGSetUp(dmmg);
296:   return(0);
297: }

301: /*@C
302:     DMMGSetUp - Prepares the DMMG to solve a system

304:     Collective on DMMG

306:     Input Parameter:
307: .   dmmg - the context

309:     Level: advanced

311: .seealso DMMGCreate(), DMMGDestroy(), DMMG, DMMGSetSNES(), DMMGSetKSP(), DMMGSolve(), DMMGSetMatType()

313: @*/
314: PetscErrorCode  DMMGSetUp(DMMG *dmmg)
315: {
317:   PetscInt       i,nlevels = dmmg[0]->nlevels;


321:   /* Create work vectors and matrix for each level */
322:   for (i=0; i<nlevels; i++) {
323:     DMCreateGlobalVector(dmmg[i]->dm,&dmmg[i]->x);
324:     VecDuplicate(dmmg[i]->x,&dmmg[i]->b);
325:     VecDuplicate(dmmg[i]->x,&dmmg[i]->r);
326:   }

328:   /* Create interpolation/restriction between levels */
329:   for (i=1; i<nlevels; i++) {
330:     DMGetInterpolation(dmmg[i-1]->dm,dmmg[i]->dm,&dmmg[i]->R,PETSC_NULL);
331:   }

333:   return(0);
334: }

338: /*@C
339:     DMMGSolve - Actually solves the (non)linear system defined with the DMMG

341:     Collective on DMMG

343:     Input Parameter:
344: .   dmmg - the context

346:     Level: advanced

348:     Options Database:
349: +   -dmmg_grid_sequence - use grid sequencing to get the initial solution for each level from the previous
350: -   -dmmg_monitor_solution - display the solution at each iteration

352:      Notes: For linear (KSP) problems may be called more than once, uses the same 
353:     matrices but recomputes the right hand side for each new solve. Call DMMGSetKSP()
354:     to generate new matrices.
355:  
356: .seealso DMMGCreate(), DMMGDestroy(), DMMG, DMMGSetSNES(), DMMGSetKSP(), DMMGSetUp(), DMMGSetMatType()

358: @*/
359: PetscErrorCode  DMMGSolve(DMMG *dmmg)
360: {
362:   PetscInt       i,nlevels = dmmg[0]->nlevels;
363:   PetscTruth     gridseq,vecmonitor,flg;

366:   PetscOptionsHasName(0,"-dmmg_grid_sequence",&gridseq);
367:   PetscOptionsHasName(0,"-dmmg_monitor_solution",&vecmonitor);
368:   if (gridseq) {
369:     if (dmmg[0]->initialguess) {
370:       (*dmmg[0]->initialguess)(dmmg[0],dmmg[0]->x);
371:       if (dmmg[0]->ksp && !dmmg[0]->snes) {
372:         KSPSetInitialGuessNonzero(dmmg[0]->ksp,PETSC_TRUE);
373:       }
374:     }
375:     for (i=0; i<nlevels-1; i++) {
376:       (*dmmg[i]->solve)(dmmg,i);
377:       if (vecmonitor) {
378:         VecView(dmmg[i]->x,PETSC_VIEWER_DRAW_(dmmg[i]->comm));
379:       }
380:       MatInterpolate(dmmg[i+1]->R,dmmg[i]->x,dmmg[i+1]->x);
381:       if (dmmg[i+1]->ksp && !dmmg[i+1]->snes) {
382:         KSPSetInitialGuessNonzero(dmmg[i+1]->ksp,PETSC_TRUE);
383:       }
384:     }
385:   } else {
386:     if (dmmg[nlevels-1]->initialguess) {
387:       (*dmmg[nlevels-1]->initialguess)(dmmg[nlevels-1],dmmg[nlevels-1]->x);
388:     }
389:   }

391:   /*VecView(dmmg[nlevels-1]->x,PETSC_VIEWER_DRAW_WORLD);*/

393:   (*DMMGGetFine(dmmg)->solve)(dmmg,nlevels-1);
394:   if (vecmonitor) {
395:      VecView(dmmg[nlevels-1]->x,PETSC_VIEWER_DRAW_(dmmg[nlevels-1]->comm));
396:   }

398:   PetscOptionsHasName(PETSC_NULL,"-dmmg_view",&flg);
399:   if (flg && !PetscPreLoadingOn) {
400:     PetscViewer viewer;
401:     PetscViewerASCIIGetStdout(dmmg[0]->comm,&viewer);
402:     DMMGView(dmmg,viewer);
403:   }
404:   PetscOptionsHasName(PETSC_NULL,"-dmmg_view_binary",&flg);
405:   if (flg && !PetscPreLoadingOn) {
406:     DMMGView(dmmg,PETSC_VIEWER_BINARY_(dmmg[0]->comm));
407:   }
408:   return(0);
409: }

413: PetscErrorCode  DMMGSolveKSP(DMMG *dmmg,PetscInt level)
414: {

418:   if (dmmg[level]->rhs) {
419:     CHKMEMQ;
420:     (*dmmg[level]->rhs)(dmmg[level],dmmg[level]->b);
421:     CHKMEMQ;
422:   }
423:   KSPSolve(dmmg[level]->ksp,dmmg[level]->b,dmmg[level]->x);
424:   return(0);
425: }

427: /*
428:     For each level (of grid sequencing) this sets the interpolation/restriction and 
429:     work vectors needed by the multigrid preconditioner within the KSP 
430:     (for nonlinear problems the KSP inside the SNES) of that level.

432:     Also sets the KSP monitoring on all the levels if requested by user.

434: */
437: PetscErrorCode  DMMGSetUpLevel(DMMG *dmmg,KSP ksp,PetscInt nlevels)
438: {
439:   PetscErrorCode          ierr;
440:   PetscInt                i;
441:   PC                      pc;
442:   PetscTruth              ismg,monitor,monitor_short,ismf,isshell,ismffd;
443:   KSP                     lksp; /* solver internal to the multigrid preconditioner */
444:   MPI_Comm                *comms,comm;
445:   PetscViewerASCIIMonitor ascii;

448:   if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");

450:   PetscOptionsHasName(PETSC_NULL,"-dmmg_ksp_monitor",&monitor);
451:   PetscOptionsHasName(PETSC_NULL,"-dmmg_ksp_monitor_short",&monitor_short);
452:   if (monitor || monitor_short) {
453:     PetscObjectGetComm((PetscObject)ksp,&comm);
454:     PetscViewerASCIIMonitorCreate(comm,"stdout",1+dmmg[0]->nlevels-nlevels,&ascii);
455:     if (monitor) {
456:       KSPMonitorSet(ksp,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
457:     } else {
458:       KSPMonitorSet(ksp,KSPMonitorDefaultShort,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
459:     }
460:   }

462:   /* use fgmres on outer iteration by default */
463:   KSPSetType(ksp,KSPFGMRES);
464:   KSPGetPC(ksp,&pc);
465:   PCSetType(pc,PCMG);
466:   PetscMalloc(nlevels*sizeof(MPI_Comm),&comms);
467:   for (i=0; i<nlevels; i++) {
468:     comms[i] = dmmg[i]->comm;
469:   }
470:   PCMGSetLevels(pc,nlevels,comms);
471:   PetscFree(comms);
472:    PCMGSetType(pc,PC_MG_FULL);

474:   PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
475:   if (ismg) {
476:     /* set solvers for each level */
477:     for (i=0; i<nlevels; i++) {
478:       if (i < nlevels-1) { /* don't set for finest level, they are set in PCApply_MG()*/
479:         PCMGSetX(pc,i,dmmg[i]->x);
480:         PCMGSetRhs(pc,i,dmmg[i]->b);
481:       }
482:       if (i > 0) {
483:         PCMGSetR(pc,i,dmmg[i]->r);
484:       }
485:       if (monitor || monitor_short) {
486:         PCMGGetSmoother(pc,i,&lksp);
487:         PetscObjectGetComm((PetscObject)lksp,&comm);
488:         PetscViewerASCIIMonitorCreate(comm,"stdout",1+dmmg[0]->nlevels-i,&ascii);
489:         if (monitor) {
490:           KSPMonitorSet(lksp,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
491:         } else {
492:           KSPMonitorSet(lksp,KSPMonitorDefaultShort,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
493:         }
494:       }
495:       /* If using a matrix free multiply and did not provide an explicit matrix to build
496:          the preconditioner then must use no preconditioner 
497:       */
498:       PetscTypeCompare((PetscObject)dmmg[i]->B,MATSHELL,&isshell);
499:       PetscTypeCompare((PetscObject)dmmg[i]->B,MATDAAD,&ismf);
500:       PetscTypeCompare((PetscObject)dmmg[i]->B,MATMFFD,&ismffd);
501:       if (isshell || ismf || ismffd) {
502:         PC  lpc;
503:         PCMGGetSmoother(pc,i,&lksp);
504:         KSPGetPC(lksp,&lpc);
505:         PCSetType(lpc,PCNONE);
506:       }
507:     }

509:     /* Set interpolation/restriction between levels */
510:     for (i=1; i<nlevels; i++) {
511:       PCMGSetInterpolation(pc,i,dmmg[i]->R);
512:       PCMGSetRestriction(pc,i,dmmg[i]->R);
513:     }
514:   }
515:   return(0);
516: }

520: /*@C
521:     DMMGSetKSP - Sets the linear solver object that will use the grid hierarchy

523:     Collective on DMMG

525:     Input Parameter:
526: +   dmmg - the context
527: .   func - function to compute linear system matrix on each grid level
528: -   rhs - function to compute right hand side on each level (need only work on the finest grid
529:           if you do not use grid sequencing)

531:     Level: advanced

533:     Notes: For linear problems my be called more than once, reevaluates the matrices if it is called more
534:        than once. Call DMMGSolve() directly several times to solve with the same matrix but different 
535:        right hand sides.
536:    
537: .seealso DMMGCreate(), DMMGDestroy, DMMGSetDM(), DMMGSolve(), DMMGSetMatType()

539: @*/
540: PetscErrorCode  DMMGSetKSP(DMMG *dmmg,PetscErrorCode (*rhs)(DMMG,Vec),PetscErrorCode (*func)(DMMG,Mat,Mat))
541: {
543:   PetscInt       i,nlevels = dmmg[0]->nlevels,level;
544:   PetscTruth     galerkin,ismg;
545:   PC             pc;
546:   KSP            lksp;

549:   if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
550:   galerkin = dmmg[nlevels - 2 > 0 ? nlevels - 2 : 0]->galerkin;

552:   if (galerkin) {
553:     DMGetMatrix(dmmg[nlevels-1]->dm,dmmg[nlevels-1]->mtype,&dmmg[nlevels-1]->B);
554:     if (!dmmg[nlevels-1]->J) {
555:       dmmg[nlevels-1]->J = dmmg[nlevels-1]->B;
556:       PetscObjectReference((PetscObject) dmmg[nlevels-1]->J);
557:     }
558:     if (func) {
559:       (*func)(dmmg[nlevels-1],dmmg[nlevels-1]->J,dmmg[nlevels-1]->B);
560:     }
561:     for (i=nlevels-2; i>-1; i--) {
562:       if (dmmg[i]->galerkin) {
563:         MatPtAP(dmmg[i+1]->B,dmmg[i+1]->R,MAT_INITIAL_MATRIX,1.0,&dmmg[i]->B);
564:         if (!dmmg[i]->J) {
565:           dmmg[i]->J = dmmg[i]->B;
566:           PetscObjectReference((PetscObject) dmmg[i]->J);
567:         }
568:       }
569:     }
570:   }

572:   if (!dmmg[0]->ksp) {
573:     /* create solvers for each level if they don't already exist*/
574:     for (i=0; i<nlevels; i++) {

576:       if (!dmmg[i]->B && !dmmg[i]->galerkin) {
577:         DMGetMatrix(dmmg[i]->dm,dmmg[nlevels-1]->mtype,&dmmg[i]->B);
578:       }
579:       if (!dmmg[i]->J) {
580:         dmmg[i]->J = dmmg[i]->B;
581:         PetscObjectReference((PetscObject) dmmg[i]->J);
582:       }

584:       KSPCreate(dmmg[i]->comm,&dmmg[i]->ksp);
585:       KSPSetOptionsPrefix(dmmg[i]->ksp,dmmg[i]->prefix);
586:       DMMGSetUpLevel(dmmg,dmmg[i]->ksp,i+1);
587:       KSPSetFromOptions(dmmg[i]->ksp);
588:       dmmg[i]->solve = DMMGSolveKSP;
589:       dmmg[i]->rhs   = rhs;
590:     }
591:   }

593:   /* evalute matrix on each level */
594:   for (i=0; i<nlevels; i++) {
595:     if (!dmmg[i]->galerkin) {
596:       if (func) {
597:         (*func)(dmmg[i],dmmg[i]->J,dmmg[i]->B);
598:       }
599:     }
600:   }

602:   for (i=0; i<nlevels-1; i++) {
603:     KSPSetOptionsPrefix(dmmg[i]->ksp,"dmmg_");
604:   }

606:   for (level=0; level<nlevels; level++) {
607:     KSPSetOperators(dmmg[level]->ksp,dmmg[level]->J,dmmg[level]->B,SAME_NONZERO_PATTERN);
608:     KSPGetPC(dmmg[level]->ksp,&pc);
609:     PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
610:     if (ismg) {
611:       for (i=0; i<=level; i++) {
612:         PCMGGetSmoother(pc,i,&lksp);
613:         KSPSetOperators(lksp,dmmg[i]->J,dmmg[i]->B,SAME_NONZERO_PATTERN);
614:       }
615:     }
616:   }

618:   return(0);
619: }

623: /*@C
624:     DMMGView - prints information on a DA based multi-level preconditioner

626:     Collective on DMMG and PetscViewer

628:     Input Parameter:
629: +   dmmg - the context
630: -   viewer - the viewer

632:     Level: advanced

634: .seealso DMMGCreate(), DMMGDestroy(), DMMGSetMatType(), DMMGSetUseGalerkin()

636: @*/
637: PetscErrorCode  DMMGView(DMMG *dmmg,PetscViewer viewer)
638: {
640:   PetscInt       i,nlevels = dmmg[0]->nlevels;
641:   PetscMPIInt    flag;
642:   MPI_Comm       comm;
643:   PetscTruth     iascii,isbinary;

648:   PetscObjectGetComm((PetscObject)viewer,&comm);
649:   MPI_Comm_compare(comm,dmmg[0]->comm,&flag);
650:   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) {
651:     SETERRQ(PETSC_ERR_ARG_NOTSAMECOMM,"Different communicators in the DMMG and the PetscViewer");
652:   }

654:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
655:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
656:   if (isbinary) {
657:     for (i=0; i<nlevels; i++) {
658:       MatView(dmmg[i]->J,viewer);
659:     }
660:     for (i=1; i<nlevels; i++) {
661:       MatView(dmmg[i]->R,viewer);
662:     }
663:   } else {
664:     if (iascii) {
665:       PetscViewerASCIIPrintf(viewer,"DMMG Object with %D levels\n",nlevels);
666:     }
667:     for (i=0; i<nlevels; i++) {
668:       PetscViewerASCIIPushTab(viewer);
669:       DMView(dmmg[i]->dm,viewer);
670:       PetscViewerASCIIPopTab(viewer);
671:     }
672:     if (iascii) {
673:       PetscViewerASCIIPrintf(viewer,"%s Object on finest level\n",dmmg[nlevels-1]->ksp ? "KSP" : "SNES");
674:       if (dmmg[nlevels-2 > 0 ? nlevels-2 : 0]->galerkin) {
675:         PetscViewerASCIIPrintf(viewer,"Using Galerkin R^T*A*R process to compute coarser matrices\n");
676:       }
677:       PetscViewerASCIIPrintf(viewer,"Using matrix type %s\n",dmmg[nlevels-1]->mtype);
678:     }
679:     if (dmmg[nlevels-1]->ksp) {
680:       KSPView(dmmg[nlevels-1]->ksp,viewer);
681:     } else {
682:       /* use of PetscObjectView() means we do not have to link with libpetscsnes if SNES is not being used */
683:       PetscObjectView((PetscObject)dmmg[nlevels-1]->snes,viewer);
684:     }
685:   }
686:   return(0);
687: }

691: /*@C
692:     DMMGSetNullSpace - Indicates the null space in the linear operator (this is needed by the linear solver)

694:     Collective on DMMG

696:     Input Parameter:
697: +   dmmg - the context
698: .   has_cnst - is the constant vector in the null space
699: .   n - number of null vectors (excluding the possible constant vector)
700: -   func - a function that fills an array of vectors with the null vectors (must be orthonormal), may be PETSC_NULL

702:     Level: advanced

704: .seealso DMMGCreate(), DMMGDestroy, DMMGSetDM(), DMMGSolve(), MatNullSpaceCreate(), KSPSetNullSpace(), DMMGSetUseGalerkin(), DMMGSetMatType()

706: @*/
707: PetscErrorCode  DMMGSetNullSpace(DMMG *dmmg,PetscTruth has_cnst,PetscInt n,PetscErrorCode (*func)(DMMG,Vec[]))
708: {
710:   PetscInt       i,j,nlevels = dmmg[0]->nlevels;
711:   Vec            *nulls = 0;
712:   MatNullSpace   nullsp;
713:   KSP            iksp;
714:   PC             pc,ipc;
715:   PetscTruth     ismg,isred;

718:   if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
719:   if (!dmmg[0]->ksp) SETERRQ(PETSC_ERR_ORDER,"Must call AFTER DMMGSetKSP() or DMMGSetSNES()");
720:   if ((n && !func) || (!n && func)) SETERRQ(PETSC_ERR_ARG_INCOMP,"Both n and func() must be set together");
721:   if (n < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative number of vectors in null space n = %D",n)

723:   for (i=0; i<nlevels; i++) {
724:     if (n) {
725:       VecDuplicateVecs(dmmg[i]->b,n,&nulls);
726:       (*func)(dmmg[i],nulls);
727:     }
728:     MatNullSpaceCreate(dmmg[i]->comm,has_cnst,n,nulls,&nullsp);
729:     KSPSetNullSpace(dmmg[i]->ksp,nullsp);
730:     for (j=i; j<nlevels; j++) {
731:       KSPGetPC(dmmg[j]->ksp,&pc);
732:       PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
733:       if (ismg) {
734:         PCMGGetSmoother(pc,i,&iksp);
735:         KSPSetNullSpace(iksp, nullsp);
736:       }
737:     }
738:     MatNullSpaceDestroy(nullsp);
739:     if (n) {
740:       PetscFree(nulls);
741:     }
742:   }
743:   /* make all the coarse grid solvers have LU shift since they are singular */
744:   for (i=0; i<nlevels; i++) {
745:     KSPGetPC(dmmg[i]->ksp,&pc);
746:     PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
747:     if (ismg) {
748:       PCMGGetSmoother(pc,0,&iksp);
749:       KSPGetPC(iksp,&ipc);
750:       PetscTypeCompare((PetscObject)ipc,PCREDUNDANT,&isred);
751:       if (isred) {
752:         PCRedundantGetPC(ipc,&ipc);
753:       }
754:       PCFactorSetShiftPd(ipc,PETSC_TRUE);
755:     }
756:   }
757:   return(0);
758: }

762: /*@C
763:     DMMGInitialGuessCurrent - Use with DMMGSetInitialGuess() to use the current value in the 
764:        solution vector (obtainable with DMMGGetx()) as the initial guess. Otherwise for linear
765:        problems zero is used for the initial guess (unless grid sequencing is used). For nonlinear 
766:        problems this is not needed; it always uses the previous solution as the initial guess.

768:     Collective on DMMG

770:     Input Parameter:
771: +   dmmg - the context
772: -   vec - dummy argument

774:     Level: intermediate

776: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess()

778: @*/
779: PetscErrorCode  DMMGInitialGuessCurrent(DMMG dmmg,Vec vec)
780: {
782:   return(0);
783: }

787: /*@C
788:     DMMGSetInitialGuess - Sets the function that computes an initial guess.

790:     Collective on DMMG

792:     Input Parameter:
793: +   dmmg - the context
794: -   guess - the function

796:     Notes: For nonlinear problems, if this is not set, then the current value in the 
797:              solution vector (obtained with DMMGGetX()) is used. Thus is if you doing 'time
798:              stepping' it will use your current solution as the guess for the next timestep.
799:            If grid sequencing is used (via -dmmg_grid_sequence) then the "guess" function
800:              is used only on the coarsest grid.
801:            For linear problems, if this is not set, then 0 is used as an initial guess.
802:              If you would like the linear solver to also (like the nonlinear solver) use
803:              the current solution vector as the initial guess then use DMMGInitialGuessCurrent()
804:              as the function you pass in

806:     Level: intermediate


809: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGInitialGuessCurrent(), DMMGSetGalekin(), DMMGSetMatType(), DMMGSetNullSpace()

811: @*/
812: PetscErrorCode  DMMGSetInitialGuess(DMMG *dmmg,PetscErrorCode (*guess)(DMMG,Vec))
813: {
814:   PetscInt i,nlevels = dmmg[0]->nlevels;

817:   for (i=0; i<nlevels; i++) {
818:     dmmg[i]->initialguess = guess;
819:   }
820:   return(0);
821: }