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:   PetscOptionsGetInt(0,"-dmmg_nlevels",&nlevels,PETSC_IGNORE);
 54:   PetscOptionsHasName(0,"-dmmg_galerkin",&galerkin);

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

 69:   PetscOptionsGetString(PETSC_NULL,"-dmmg_mat_type",mtype,256,&ftype);
 70:   if (ftype) {
 71:     DMMGSetMatType(*dmmg,mtype);
 72:   }
 73:   return(0);
 74: }

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

 81:     Collective on MPI_Comm 

 83:     Input Parameters:
 84: +    dmmg - the DMMG object created with DMMGCreate()
 85: -    mtype - the matrix type, defaults to MATAIJ

 87:     Level: intermediate

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

 91: @*/
 92: PetscErrorCode  DMMGSetMatType(DMMG *dmmg,MatType mtype)
 93: {
 94:   PetscInt i;
 95: 
 97:   for (i=0; i<dmmg[0]->nlevels; i++) {
 98:     dmmg[i]->mtype  = mtype;
 99:   }
100:   return(0);
101: }

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

108:     Collective on MPI_Comm 

110:     Input Parameters:
111: +    dmmg - the DMMG object created with DMMGCreate()
112: -    prefix - the prefix string

114:     Level: intermediate

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

118: @*/
119: PetscErrorCode  DMMGSetPrefix(DMMG *dmmg,const char* prefix)
120: {
121:   PetscInt       i;
123: 
125:   for (i=0; i<dmmg[0]->nlevels; i++) {
126:     PetscStrallocpy(prefix,&dmmg[i]->prefix);
127:   }
128:   return(0);
129: }

133: /*@C
134:     DMMGSetUseGalerkinCoarse - Courses the DMMG to use R*A_f*R^T to form
135:        the coarser matrices from finest 

137:     Collective on DMMG

139:     Input Parameter:
140: .    dmmg - the context

142:     Options Database Keys:
143: .    -dmmg_galerkin 

145:     Level: advanced

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

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

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

157: @*/
158: PetscErrorCode  DMMGSetUseGalerkinCoarse(DMMG* dmmg)
159: {
160:   PetscInt  i,nlevels = dmmg[0]->nlevels;

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

165:   for (i=0; i<nlevels-1; i++) {
166:     dmmg[i]->galerkin = PETSC_TRUE;
167:   }
168:   return(0);
169: }

173: /*@C
174:     DMMGDestroy - Destroys a DA based multigrid solver object. 

176:     Collective on DMMG

178:     Input Parameter:
179: .    - the context

181:     Level: advanced

183: .seealso DMMGCreate()

185: @*/
186: PetscErrorCode  DMMGDestroy(DMMG *dmmg)
187: {
189:   PetscInt       i,nlevels = dmmg[0]->nlevels;

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

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

222: /*@C
223:     DMMGSetDM - Sets the coarse grid information for the grids

225:     Collective on DMMG

227:     Input Parameter:
228: +   dmmg - the context
229: -   dm - the DA or DMComposite object

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

235:     Level: advanced

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

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

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

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

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

275:       DMCoarsenHierarchy(dm, nlevels-1, &hierarchy);
276:       for(i = 0; i < nlevels-1; ++i) {
277:         dmmg[nlevels-2-i]->dm = hierarchy[i];
278:       }
279:     } else {
280: /*       for(i = nlevels-2; i >= 0; --i) { */
281: /*         DMCoarsen(dmmg[i+1]->dm, dmmg[i]->comm, &dmmg[i]->dm); */
282: /*       } */
283:       SETERRQ(PETSC_ERR_SUP, "Sequential coarsening not yet implemented");
284:     }
285:   }
286:   DMMGSetUp(dmmg);
287:   return(0);
288: }

292: /*@C
293:     DMMGSetUp - Prepares the DMMG to solve a system

295:     Collective on DMMG

297:     Input Parameter:
298: .   dmmg - the context

300:     Level: advanced

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

304: @*/
305: PetscErrorCode  DMMGSetUp(DMMG *dmmg)
306: {
308:   PetscInt       i,nlevels = dmmg[0]->nlevels;


312:   /* Create work vectors and matrix for each level */
313:   for (i=0; i<nlevels; i++) {
314:     DMCreateGlobalVector(dmmg[i]->dm,&dmmg[i]->x);
315:     VecDuplicate(dmmg[i]->x,&dmmg[i]->b);
316:     VecDuplicate(dmmg[i]->x,&dmmg[i]->r);
317:   }

319:   /* Create interpolation/restriction between levels */
320:   for (i=1; i<nlevels; i++) {
321:     DMGetInterpolation(dmmg[i-1]->dm,dmmg[i]->dm,&dmmg[i]->R,PETSC_NULL);
322:   }

324:   return(0);
325: }

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

332:     Collective on DMMG

334:     Input Parameter:
335: .   dmmg - the context

337:     Level: advanced

339:     Options Database:
340: +   -dmmg_grid_sequence - use grid sequencing to get the initial solution for each level from the previous
341: -   -dmmg_monitor_solution - display the solution at each iteration

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

349: @*/
350: PetscErrorCode  DMMGSolve(DMMG *dmmg)
351: {
353:   PetscInt       i,nlevels = dmmg[0]->nlevels;
354:   PetscTruth     gridseq,vecmonitor,flg;

357:   PetscOptionsHasName(0,"-dmmg_grid_sequence",&gridseq);
358:   PetscOptionsHasName(0,"-dmmg_monitor_solution",&vecmonitor);
359:   if (gridseq) {
360:     if (dmmg[0]->initialguess) {
361:       (*dmmg[0]->initialguess)(dmmg[0],dmmg[0]->x);
362:       if (dmmg[0]->ksp && !dmmg[0]->snes) {
363:         KSPSetInitialGuessNonzero(dmmg[0]->ksp,PETSC_TRUE);
364:       }
365:     }
366:     for (i=0; i<nlevels-1; i++) {
367:       (*dmmg[i]->solve)(dmmg,i);
368:       if (vecmonitor) {
369:         VecView(dmmg[i]->x,PETSC_VIEWER_DRAW_(dmmg[i]->comm));
370:       }
371:       MatInterpolate(dmmg[i+1]->R,dmmg[i]->x,dmmg[i+1]->x);
372:       if (dmmg[i+1]->ksp && !dmmg[i+1]->snes) {
373:         KSPSetInitialGuessNonzero(dmmg[i+1]->ksp,PETSC_TRUE);
374:       }
375:     }
376:   } else {
377:     if (dmmg[nlevels-1]->initialguess) {
378:       (*dmmg[nlevels-1]->initialguess)(dmmg[nlevels-1],dmmg[nlevels-1]->x);
379:     }
380:   }

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

384:   (*DMMGGetFine(dmmg)->solve)(dmmg,nlevels-1);
385:   if (vecmonitor) {
386:      VecView(dmmg[nlevels-1]->x,PETSC_VIEWER_DRAW_(dmmg[nlevels-1]->comm));
387:   }

389:   PetscOptionsHasName(PETSC_NULL,"-dmmg_view",&flg);
390:   if (flg && !PetscPreLoadingOn) {
391:     PetscViewer viewer;
392:     PetscViewerASCIIGetStdout(dmmg[0]->comm,&viewer);
393:     DMMGView(dmmg,viewer);
394:   }
395:   PetscOptionsHasName(PETSC_NULL,"-dmmg_view_binary",&flg);
396:   if (flg && !PetscPreLoadingOn) {
397:     DMMGView(dmmg,PETSC_VIEWER_BINARY_(dmmg[0]->comm));
398:   }
399:   return(0);
400: }

404: PetscErrorCode  DMMGSolveKSP(DMMG *dmmg,PetscInt level)
405: {

409:   if (dmmg[level]->rhs) {
410:     CHKMEMQ;
411:     (*dmmg[level]->rhs)(dmmg[level],dmmg[level]->b);
412:     CHKMEMQ;
413:   }
414:   KSPSolve(dmmg[level]->ksp,dmmg[level]->b,dmmg[level]->x);
415:   return(0);
416: }

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

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

425: */
428: PetscErrorCode  DMMGSetUpLevel(DMMG *dmmg,KSP ksp,PetscInt nlevels)
429: {
430:   PetscErrorCode          ierr;
431:   PetscInt                i;
432:   PC                      pc;
433:   PetscTruth              ismg,monitor,monitor_short,ismf,isshell,ismffd;
434:   KSP                     lksp; /* solver internal to the multigrid preconditioner */
435:   MPI_Comm                *comms,comm;
436:   PetscViewerASCIIMonitor ascii;

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

441:   PetscOptionsHasName(PETSC_NULL,"-dmmg_ksp_monitor",&monitor);
442:   PetscOptionsHasName(PETSC_NULL,"-dmmg_ksp_monitor_short",&monitor_short);
443:   if (monitor || monitor_short) {
444:     PetscObjectGetComm((PetscObject)ksp,&comm);
445:     PetscViewerASCIIMonitorCreate(comm,"stdout",1+dmmg[0]->nlevels-nlevels,&ascii);
446:     if (monitor) {
447:       KSPMonitorSet(ksp,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
448:     } else {
449:       KSPMonitorSet(ksp,KSPMonitorDefaultShort,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);
450:     }
451:   }

453:   /* use fgmres on outer iteration by default */
454:   KSPSetType(ksp,KSPFGMRES);
455:   KSPGetPC(ksp,&pc);
456:   PCSetType(pc,PCMG);
457:   PetscMalloc(nlevels*sizeof(MPI_Comm),&comms);
458:   for (i=0; i<nlevels; i++) {
459:     comms[i] = dmmg[i]->comm;
460:   }
461:   PCMGSetLevels(pc,nlevels,comms);
462:   PetscFree(comms);
463:    PCMGSetType(pc,PC_MG_FULL);

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

500:     /* Set interpolation/restriction between levels */
501:     for (i=1; i<nlevels; i++) {
502:       PCMGSetInterpolation(pc,i,dmmg[i]->R);
503:       PCMGSetRestriction(pc,i,dmmg[i]->R);
504:     }
505:   }
506:   return(0);
507: }

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

514:     Collective on DMMG

516:     Input Parameter:
517: +   dmmg - the context
518: .   func - function to compute linear system matrix on each grid level
519: -   rhs - function to compute right hand side on each level (need only work on the finest grid
520:           if you do not use grid sequencing)

522:     Level: advanced

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

530: @*/
531: PetscErrorCode  DMMGSetKSP(DMMG *dmmg,PetscErrorCode (*rhs)(DMMG,Vec),PetscErrorCode (*func)(DMMG,Mat,Mat))
532: {
534:   PetscInt       i,nlevels = dmmg[0]->nlevels,level;
535:   PetscTruth     galerkin,ismg;
536:   PC             pc;
537:   KSP            lksp;

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

543:   if (galerkin) {
544:     DMGetMatrix(dmmg[nlevels-1]->dm,dmmg[nlevels-1]->mtype,&dmmg[nlevels-1]->B);
545:     if (!dmmg[nlevels-1]->J) {
546:       dmmg[nlevels-1]->J = dmmg[nlevels-1]->B;
547:     }
548:     if (func) {
549:       (*func)(dmmg[nlevels-1],dmmg[nlevels-1]->J,dmmg[nlevels-1]->B);
550:     }
551:     for (i=nlevels-2; i>-1; i--) {
552:       if (dmmg[i]->galerkin) {
553:         MatPtAP(dmmg[i+1]->B,dmmg[i+1]->R,MAT_INITIAL_MATRIX,1.0,&dmmg[i]->B);
554:         if (!dmmg[i]->J) {
555:           dmmg[i]->J = dmmg[i]->B;
556:         }
557:       }
558:     }
559:   }

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

565:       if (!dmmg[i]->B && !dmmg[i]->galerkin) {
566:         DMGetMatrix(dmmg[i]->dm,dmmg[nlevels-1]->mtype,&dmmg[i]->B);
567:       }
568:       if (!dmmg[i]->J) {
569:         dmmg[i]->J = dmmg[i]->B;
570:       }

572:       KSPCreate(dmmg[i]->comm,&dmmg[i]->ksp);
573:       KSPSetOptionsPrefix(dmmg[i]->ksp,dmmg[i]->prefix);
574:       DMMGSetUpLevel(dmmg,dmmg[i]->ksp,i+1);
575:       KSPSetFromOptions(dmmg[i]->ksp);
576:       dmmg[i]->solve = DMMGSolveKSP;
577:       dmmg[i]->rhs   = rhs;
578:     }
579:   }

581:   /* evalute matrix on each level */
582:   for (i=0; i<nlevels; i++) {
583:     if (!dmmg[i]->galerkin) {
584:       if (func) {
585:         (*func)(dmmg[i],dmmg[i]->J,dmmg[i]->B);
586:       }
587:     }
588:   }

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

594:   for (level=0; level<nlevels; level++) {
595:     KSPSetOperators(dmmg[level]->ksp,dmmg[level]->J,dmmg[level]->B,SAME_NONZERO_PATTERN);
596:     KSPGetPC(dmmg[level]->ksp,&pc);
597:     PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
598:     if (ismg) {
599:       for (i=0; i<=level; i++) {
600:         PCMGGetSmoother(pc,i,&lksp);
601:         KSPSetOperators(lksp,dmmg[i]->J,dmmg[i]->B,SAME_NONZERO_PATTERN);
602:       }
603:     }
604:   }

606:   return(0);
607: }

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

614:     Collective on DMMG and PetscViewer

616:     Input Parameter:
617: +   dmmg - the context
618: -   viewer - the viewer

620:     Level: advanced

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

624: @*/
625: PetscErrorCode  DMMGView(DMMG *dmmg,PetscViewer viewer)
626: {
628:   PetscInt       i,nlevels = dmmg[0]->nlevels;
629:   PetscMPIInt    flag;
630:   MPI_Comm       comm;
631:   PetscTruth     iascii,isbinary;

636:   PetscObjectGetComm((PetscObject)viewer,&comm);
637:   MPI_Comm_compare(comm,dmmg[0]->comm,&flag);
638:   if (flag != MPI_CONGRUENT && flag != MPI_IDENT) {
639:     SETERRQ(PETSC_ERR_ARG_NOTSAMECOMM,"Different communicators in the DMMG and the PetscViewer");
640:   }

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

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

682:     Collective on DMMG

684:     Input Parameter:
685: +   dmmg - the context
686: .   has_cnst - is the constant vector in the null space
687: .   n - number of null vectors (excluding the possible constant vector)
688: -   func - a function that fills an array of vectors with the null vectors (must be orthonormal), may be PETSC_NULL

690:     Level: advanced

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

694: @*/
695: PetscErrorCode  DMMGSetNullSpace(DMMG *dmmg,PetscTruth has_cnst,PetscInt n,PetscErrorCode (*func)(DMMG,Vec[]))
696: {
698:   PetscInt       i,j,nlevels = dmmg[0]->nlevels;
699:   Vec            *nulls = 0;
700:   MatNullSpace   nullsp;
701:   KSP            iksp;
702:   PC             pc,ipc;
703:   PetscTruth     ismg,isred;

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

711:   for (i=0; i<nlevels; i++) {
712:     if (n) {
713:       VecDuplicateVecs(dmmg[i]->b,n,&nulls);
714:       (*func)(dmmg[i],nulls);
715:     }
716:     MatNullSpaceCreate(dmmg[i]->comm,has_cnst,n,nulls,&nullsp);
717:     KSPSetNullSpace(dmmg[i]->ksp,nullsp);
718:     for (j=i; j<nlevels; j++) {
719:       KSPGetPC(dmmg[j]->ksp,&pc);
720:       PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
721:       if (ismg) {
722:         PCMGGetSmoother(pc,i,&iksp);
723:         KSPSetNullSpace(iksp, nullsp);
724:       }
725:     }
726:     MatNullSpaceDestroy(nullsp);
727:     if (n) {
728:       PetscFree(nulls);
729:     }
730:   }
731:   /* make all the coarse grid solvers have LU shift since they are singular */
732:   for (i=0; i<nlevels; i++) {
733:     KSPGetPC(dmmg[i]->ksp,&pc);
734:     PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
735:     if (ismg) {
736:       PCMGGetSmoother(pc,0,&iksp);
737:       KSPGetPC(iksp,&ipc);
738:       PetscTypeCompare((PetscObject)ipc,PCREDUNDANT,&isred);
739:       if (isred) {
740:         PCRedundantGetPC(ipc,&ipc);
741:       }
742:       PCFactorSetShiftPd(ipc,PETSC_TRUE);
743:     }
744:   }
745:   return(0);
746: }

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

756:     Collective on DMMG

758:     Input Parameter:
759: +   dmmg - the context
760: -   vec - dummy argument

762:     Level: intermediate

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

766: @*/
767: PetscErrorCode  DMMGInitialGuessCurrent(DMMG dmmg,Vec vec)
768: {
770:   return(0);
771: }

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

778:     Collective on DMMG

780:     Input Parameter:
781: +   dmmg - the context
782: -   guess - the function

784:     Notes: For nonlinear problems, if this is not set, then the current value in the 
785:              solution vector (obtained with DMMGGetX()) is used. Thus is if you doing 'time
786:              stepping' it will use your current solution as the guess for the next timestep.
787:            If grid sequencing is used (via -dmmg_grid_sequence) then the "guess" function
788:              is used only on the coarsest grid.
789:            For linear problems, if this is not set, then 0 is used as an initial guess.
790:              If you would like the linear solver to also (like the nonlinear solver) use
791:              the current solution vector as the initial guess then use DMMGInitialGuessCurrent()
792:              as the function you pass in

794:     Level: intermediate


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

799: @*/
800: PetscErrorCode  DMMGSetInitialGuess(DMMG *dmmg,PetscErrorCode (*guess)(DMMG,Vec))
801: {
802:   PetscInt i,nlevels = dmmg[0]->nlevels;

805:   for (i=0; i<nlevels; i++) {
806:     dmmg[i]->initialguess = guess;
807:   }
808:   return(0);
809: }