Actual source code: damgsnes.c

  1: #define PETSCSNES_DLL
  2: 
 3:  #include petscda.h
 4:  #include ../src/dm/da/daimpl.h
  5: /* It appears that preprocessor directives are not respected by bfort */
  6: #ifdef PETSC_HAVE_SIEVE
 7:  #include petscmesh.h
  8: #endif
 9:  #include petscmg.h
 10:  #include petscdmmg.h

 12: #if defined(PETSC_HAVE_ADIC)
 19: #endif
 20: #if defined(PETSC_HAVE_SIEVE)
 22: #endif

 25: EXTERN PetscErrorCode  NLFCreate_DAAD(NLF*);
 26: EXTERN PetscErrorCode  NLFDAADSetDA_DAAD(NLF,DA);
 27: EXTERN PetscErrorCode  NLFDAADSetCtx_DAAD(NLF,void*);
 28: EXTERN PetscErrorCode  NLFDAADSetResidual_DAAD(NLF,Vec);
 29: EXTERN PetscErrorCode  NLFDAADSetNewtonIterations_DAAD(NLF,PetscInt);

 32: /*
 33:       period of -1 indicates update only on zeroth iteration of SNES
 34: */
 35: #define ShouldUpdate(l,it) (((dmmg[l-1]->updatejacobianperiod == -1) && (it == 0)) || \
 36:                             ((dmmg[l-1]->updatejacobianperiod >   0) && !(it % dmmg[l-1]->updatejacobianperiod)))
 37: /*
 38:    Evaluates the Jacobian on all of the grids. It is used by DMMG to provide the 
 39:    ComputeJacobian() function that SNESSetJacobian() requires.
 40: */
 43: PetscErrorCode DMMGComputeJacobian_Multigrid(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
 44: {
 45:   DMMG           *dmmg = (DMMG*)ptr;
 47:   PetscInt       i,nlevels = dmmg[0]->nlevels,it;
 48:   KSP            ksp,lksp;
 49:   PC             pc;
 50:   PetscTruth     ismg,galerkin = PETSC_FALSE;
 51:   Vec            W;
 52:   MatStructure   flg;

 55:   if (!dmmg) SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as user context which should contain DMMG");
 56:   SNESGetIterationNumber(snes,&it);

 58:   /* compute Jacobian on finest grid */
 59:   if (dmmg[nlevels-1]->updatejacobian && ShouldUpdate(nlevels,it)) {
 60:     (*DMMGGetFine(dmmg)->computejacobian)(snes,X,J,B,flag,DMMGGetFine(dmmg));
 61:   } else {
 62:     PetscInfo3(0,"Skipping Jacobian, SNES iteration %D frequence %D level %D\n",it,dmmg[nlevels-1]->updatejacobianperiod,nlevels-1);
 63:     *flag = SAME_PRECONDITIONER;
 64:   }
 65:   MatMFFDSetBase(DMMGGetFine(dmmg)->J,X,PETSC_NULL);

 67:   /* create coarser grid Jacobians for preconditioner if multigrid is the preconditioner */
 68:   SNESGetKSP(snes,&ksp);
 69:   KSPGetPC(ksp,&pc);
 70:   PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
 71:   if (ismg) {
 72:     PCMGGetGalerkin(pc,&galerkin);
 73:   }

 75:   if (!galerkin) {
 76:     for (i=nlevels-1; i>0; i--) {
 77:       if (!dmmg[i-1]->w) {
 78:         VecDuplicate(dmmg[i-1]->x,&dmmg[i-1]->w);
 79:       }
 80:       W    = dmmg[i-1]->w;
 81:       /* restrict X to coarser grid */
 82:       MatRestrict(dmmg[i]->R,X,W);
 83:       X    = W;
 84:       /* scale to "natural" scaling for that grid */
 85:       VecPointwiseMult(X,X,dmmg[i]->Rscale);
 86:       /* tell the base vector for matrix free multiplies */
 87:       MatMFFDSetBase(dmmg[i-1]->J,X,PETSC_NULL);
 88:       /* compute Jacobian on coarse grid */
 89:       if (dmmg[i-1]->updatejacobian && ShouldUpdate(i,it)) {
 90:         (*dmmg[i-1]->computejacobian)(snes,X,&dmmg[i-1]->J,&dmmg[i-1]->B,&flg,dmmg[i-1]);
 91:         flg = SAME_NONZERO_PATTERN;
 92:       } else {
 93:         PetscInfo3(0,"Skipping Jacobian, SNES iteration %D frequence %D level %D\n",it,dmmg[i-1]->updatejacobianperiod,i-1);
 94:         flg = SAME_PRECONDITIONER;
 95:       }
 96:       if (ismg) {
 97:         PCMGGetSmoother(pc,i-1,&lksp);
 98:         KSPSetOperators(lksp,dmmg[i-1]->J,dmmg[i-1]->B,flg);
 99:       }
100:     }
101:   }
102:   return(0);
103: }

105: /* ---------------------------------------------------------------------------*/


110: /* 
111:    DMMGFormFunction - This is a universal global FormFunction used by the DMMG code
112:    when the user provides a local function.

114:    Input Parameters:
115: +  snes - the SNES context
116: .  X - input vector
117: -  ptr - optional user-defined context, as set by SNESSetFunction()

119:    Output Parameter:
120: .  F - function vector

122:  */
123: PetscErrorCode DMMGFormFunction(SNES snes,Vec X,Vec F,void *ptr)
124: {
125:   DMMG           dmmg = (DMMG)ptr;
127:   Vec            localX;
128:   DA             da = (DA)dmmg->dm;

131:   DAGetLocalVector(da,&localX);
132:   /*
133:      Scatter ghost points to local vector, using the 2-step process
134:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
135:   */
136:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
137:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
138:   DAFormFunction1(da,localX,F,dmmg->user);
139:   DARestoreLocalVector(da,&localX);
140:   return(0);
141: }

145: PetscErrorCode DMMGFormFunctionGhost(SNES snes,Vec X,Vec F,void *ptr)
146: {
147:   DMMG           dmmg = (DMMG)ptr;
149:   Vec            localX, localF;
150:   DA             da = (DA)dmmg->dm;

153:   DAGetLocalVector(da,&localX);
154:   DAGetLocalVector(da,&localF);
155:   /*
156:      Scatter ghost points to local vector, using the 2-step process
157:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
158:   */
159:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
160:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
161:   VecSet(F, 0.0);
162:   VecSet(localF, 0.0);
163:   DAFormFunction1(da,localX,localF,dmmg->user);
164:   DALocalToGlobalBegin(da,localF,F);
165:   DALocalToGlobalEnd(da,localF,F);
166:   DARestoreLocalVector(da,&localX);
167:   DARestoreLocalVector(da,&localF);
168:   return(0);
169: }

171: #ifdef PETSC_HAVE_SIEVE
174: /* 
175:    DMMGFormFunctionMesh - This is a universal global FormFunction used by the DMMG code
176:    when the user provides a local function.

178:    Input Parameters:
179: +  snes - the SNES context
180: .  X - input vector
181: -  ptr - This is the DMMG object

183:    Output Parameter:
184: .  F - function vector

186:  */
187: PetscErrorCode DMMGFormFunctionMesh(SNES snes, Vec X, Vec F, void *ptr)
188: {
189:   DMMG           dmmg = (DMMG) ptr;
190:   Mesh           mesh = (Mesh) dmmg->dm;
191:   SectionReal    sectionF, section;

195:   MeshGetSectionReal(mesh, "default", &section);
196:   SectionRealDuplicate(section, &sectionF);
197:   SectionRealToVec(section, mesh, SCATTER_REVERSE, X);
198:   MeshFormFunction(mesh, section, sectionF, dmmg->user);
199:   SectionRealToVec(sectionF, mesh, SCATTER_FORWARD, F);
200:   SectionRealDestroy(sectionF);
201:   SectionRealDestroy(section);
202:   return(0);
203: }

207: /* 
208:    DMMGComputeJacobianMesh - This is a universal global FormJacobian used by the DMMG code
209:    when the user provides a local function.

211:    Input Parameters:
212: +  snes - the SNES context
213: .  X - input vector
214: -  ptr - This is the DMMG object

216:    Output Parameter:
217: .  F - function vector

219:  */
220: PetscErrorCode DMMGComputeJacobianMesh(SNES snes, Vec X, Mat *J, Mat *B, MatStructure *flag, void *ptr)
221: {
222:   DMMG           dmmg = (DMMG) ptr;
223:   Mesh           mesh = (Mesh) dmmg->dm;
224:   SectionReal    sectionX;

228:   MeshGetSectionReal(mesh, "default", &sectionX);
229:   SectionRealToVec(sectionX, mesh, SCATTER_REVERSE, X);
230:   MeshFormJacobian(mesh, sectionX, *B, dmmg->user);
231:   /* Assemble true Jacobian; if it is different */
232:   if (*J != *B) {
233:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
234:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
235:   }
236:   MatSetOption(*B, MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
237:   *flag = SAME_NONZERO_PATTERN;
238:   SectionRealDestroy(sectionX);
239:   return(0);
240: }
241: #endif

245: /* 
246:    DMMGFormFunctionFD - This is a universal global FormFunction used by the DMMG code
247:    when the user provides a local function used to compute the Jacobian via FD.

249:    Input Parameters:
250: +  snes - the SNES context
251: .  X - input vector
252: -  ptr - optional user-defined context, as set by SNESSetFunction()

254:    Output Parameter:
255: .  F - function vector

257:  */
258: PetscErrorCode DMMGFormFunctionFD(SNES snes,Vec X,Vec F,void *ptr)
259: {
260:   DMMG           dmmg = (DMMG)ptr;
262:   Vec            localX;
263:   DA             da = (DA)dmmg->dm;
264:   PetscInt       N,n;
265: 
267:   /* determine whether X=localX */
268:   DAGetLocalVector(da,&localX);
269:   VecGetSize(X,&N);
270:   VecGetSize(localX,&n);

272:   if (n != N){ /* X != localX */
273:     /* Scatter ghost points to local vector, using the 2-step process
274:        DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
275:     */
276:     DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
277:     DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
278:   } else {
279:     DARestoreLocalVector(da,&localX);
280:     localX = X;
281:   }
282:   DAFormFunction(da,dmmg->lfj,localX,F,dmmg->user);
283:   if (n != N){
284:     DARestoreLocalVector(da,&localX);
285:   }
286:   return(0);
287: }

291: /*@C 
292:    SNESDAFormFunction - This is a universal function evaluation routine that
293:    may be used with SNESSetFunction() as long as the user context has a DA
294:    as its first record and the user has called DASetLocalFunction().

296:    Collective on SNES

298:    Input Parameters:
299: +  snes - the SNES context
300: .  X - input vector
301: .  F - function vector
302: -  ptr - pointer to a structure that must have a DA as its first entry. For example this 
303:          could be a DMMG, this ptr must have been passed into SNESDAFormFunction() as the context

305:    Level: intermediate

307: .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
308:           SNESSetFunction(), SNESSetJacobian()

310: @*/
311: PetscErrorCode  SNESDAFormFunction(SNES snes,Vec X,Vec F,void *ptr)
312: {
314:   Vec            localX;
315:   DA             da = *(DA*)ptr;
316:   PetscInt       N,n;
317: 
319:   if (!da) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Looks like you called SNESSetFromFuntion(snes,SNESDAFormFunction,) without the DA context");

321:   /* determine whether X=localX */
322:   DAGetLocalVector(da,&localX);
323:   VecGetSize(X,&N);
324:   VecGetSize(localX,&n);
325: 
326: 
327:   if (n != N){ /* X != localX */
328:     /* Scatter ghost points to local vector, using the 2-step process
329:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
330:     */
331:     DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
332:     DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
333:   } else {
334:     DARestoreLocalVector(da,&localX);
335:     localX = X;
336:   }
337:   DAFormFunction1(da,localX,F,ptr);
338:   if (n != N){
339:     if (PetscExceptionValue(ierr)) {
340:       PetscErrorCode pDARestoreLocalVector(da,&localX);CHKERRQ(pierr);
341:     }
342: 
343:     DARestoreLocalVector(da,&localX);
344:   }
345:   return(0);
346: }

348: /* ------------------------------------------------------------------------------*/
349:  #include private/matimpl.h
352: PetscErrorCode DMMGComputeJacobianWithFD(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
353: {
355:   DMMG           dmmg = (DMMG)ctx;
356:   MatFDColoring  color = (MatFDColoring)dmmg->fdcoloring;
357: 
359:   if (color->ctype == IS_COLORING_GHOSTED){
360:     DA            da=(DA)dmmg->dm;
361:     Vec           x1_loc;
362:     DAGetLocalVector(da,&x1_loc);
363:     DAGlobalToLocalBegin(da,x1,INSERT_VALUES,x1_loc);
364:     DAGlobalToLocalEnd(da,x1,INSERT_VALUES,x1_loc);
365:     SNESDefaultComputeJacobianColor(snes,x1_loc,J,B,flag,dmmg->fdcoloring);
366:     DARestoreLocalVector(da,&x1_loc);
367:   } else {
368:     SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,dmmg->fdcoloring);
369:   }
370:   return(0);
371: }

375: PetscErrorCode DMMGComputeJacobianWithMF(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
376: {
378: 
380:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
381:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
382:   return(0);
383: }

387: /*
388:     DMMGComputeJacobian - Evaluates the Jacobian when the user has provided
389:     a local function evaluation routine.
390: */
391: PetscErrorCode DMMGComputeJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
392: {
393:   DMMG           dmmg = (DMMG) ptr;
395:   Vec            localX;
396:   DA             da = (DA) dmmg->dm;

399:   DAGetLocalVector(da,&localX);
400:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
401:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
402:   DAComputeJacobian1(da,localX,*B,dmmg->user);
403:   DARestoreLocalVector(da,&localX);
404:   /* Assemble true Jacobian; if it is different */
405:   if (*J != *B) {
406:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
407:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
408:   }
409:   MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
410:   *flag = SAME_NONZERO_PATTERN;
411:   return(0);
412: }

416: /*
417:     SNESDAComputeJacobianWithAdifor - This is a universal Jacobian evaluation routine
418:     that may be used with SNESSetJacobian() from Fortran as long as the user context has 
419:     a DA as its first record and DASetLocalAdiforFunction() has been called.  

421:    Collective on SNES

423:    Input Parameters:
424: +  snes - the SNES context
425: .  X - input vector
426: .  J - Jacobian
427: .  B - Jacobian used in preconditioner (usally same as J)
428: .  flag - indicates if the matrix changed its structure
429: -  ptr - optional user-defined context, as set by SNESSetFunction()

431:    Level: intermediate

433: .seealso: DASetLocalFunction(), DASetLocalAdicFunction(), SNESSetFunction(), SNESSetJacobian()

435: */
436: PetscErrorCode  SNESDAComputeJacobianWithAdifor(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
437: {
438:   DA             da = *(DA*) ptr;
440:   Vec            localX;

443:   DAGetLocalVector(da,&localX);
444:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
445:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
446:   DAComputeJacobian1WithAdifor(da,localX,*B,ptr);
447:   DARestoreLocalVector(da,&localX);
448:   /* Assemble true Jacobian; if it is different */
449:   if (*J != *B) {
450:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
451:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
452:   }
453:   MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
454:   *flag = SAME_NONZERO_PATTERN;
455:   return(0);
456: }

460: /*
461:    SNESDAComputeJacobian - This is a universal Jacobian evaluation routine for a
462:    locally provided Jacobian.

464:    Collective on SNES

466:    Input Parameters:
467: +  snes - the SNES context
468: .  X - input vector
469: .  J - Jacobian
470: .  B - Jacobian used in preconditioner (usally same as J)
471: .  flag - indicates if the matrix changed its structure
472: -  ptr - optional user-defined context, as set by SNESSetFunction()

474:    Level: intermediate

476: .seealso: DASetLocalFunction(), DASetLocalJacobian(), SNESSetFunction(), SNESSetJacobian()

478: */
479: PetscErrorCode  SNESDAComputeJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
480: {
481:   DA             da = *(DA*) ptr;
483:   Vec            localX;

486:   DAGetLocalVector(da,&localX);
487:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
488:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
489:   DAComputeJacobian1(da,localX,*B,ptr);
490:   DARestoreLocalVector(da,&localX);
491:   /* Assemble true Jacobian; if it is different */
492:   if (*J != *B) {
493:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
494:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
495:   }
496:   MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
497:   *flag = SAME_NONZERO_PATTERN;
498:   return(0);
499: }

503: PetscErrorCode DMMGSolveSNES(DMMG *dmmg,PetscInt level)
504: {
506:   PetscInt       nlevels = dmmg[0]->nlevels;

509:   dmmg[0]->nlevels = level+1;
510:   SNESSolve(dmmg[level]->snes,PETSC_NULL,dmmg[level]->x);
511:   dmmg[0]->nlevels = nlevels;
512:   return(0);
513: }

515: /* ===========================================================================================================*/

519: /*@C
520:     DMMGSetSNES - Sets the nonlinear function that defines the nonlinear set of equations
521:     to be solved using the grid hierarchy.

523:     Collective on DMMG

525:     Input Parameter:
526: +   dmmg - the context
527: .   function - the function that defines the nonlinear system
528: -   jacobian - optional function to compute Jacobian

530:     Options Database Keys:
531: +    -snes_monitor
532: .    -dmmg_coloring_from_mat - use graph coloring on the actual matrix nonzero structure instead of getting the coloring from the DM
533: .    -dmmg_jacobian_fd
534: .    -dmmg_jacobian_ad
535: .    -dmmg_jacobian_mf_fd_operator
536: .    -dmmg_jacobian_mf_fd
537: .    -dmmg_jacobian_mf_ad_operator
538: .    -dmmg_jacobian_mf_ad
539: .    -dmmg_iscoloring_type
540: -
541:                                  The period at which the Jacobian is recomputed can be set differently for different levels
542:                                  of the Jacobian (for example lag all Jacobians except on the finest level).
543:                                  There is no user interface currently for setting a different period on the different levels, one must set the
544:                                  fields dmmg[i]->updatejacobian and dmmg[i]->updatejacobianperiod directly in the DMMG data structure.
545:                                  

547:     Level: advanced

549: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNESLocal(), DMMGSetFromOptions()

551: @*/
552: PetscErrorCode  DMMGSetSNES(DMMG *dmmg,PetscErrorCode (*function)(SNES,Vec,Vec,void*),PetscErrorCode (*jacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
553: {
554:   PetscErrorCode          ierr;
555:   PetscInt                i,nlevels = dmmg[0]->nlevels;
556:   PetscTruth              mffdoperator,mffd,fdjacobian;
557:   PetscTruth              useFAS = PETSC_FALSE, fasBlock, fasGMRES;
558:   PetscTruth              monitor, monitorAll;
559:   PetscInt                fasPresmooth = 1, fasPostsmooth = 1, fasCoarsesmooth = 1, fasMaxIter = 2;
560:   PetscReal               fasRtol = 1.0e-8, fasAbstol = 1.0e-50;
561: #if defined(PETSC_HAVE_ADIC)
562:   PetscTruth              mfadoperator,mfad,adjacobian;
563: #endif
564:   PetscCookie             cookie;

567:   if (!dmmg)     SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
568:   if (!jacobian) jacobian = DMMGComputeJacobianWithFD;
569:   PetscObjectGetCookie((PetscObject) dmmg[0]->dm, &cookie);

571:   PetscOptionsBegin(dmmg[0]->comm,dmmg[0]->prefix,"DMMG Options","SNES");
572:     PetscOptionsName("-dmmg_monitor","Monitor DMMG iterations","DMMG",&monitor);
573:     PetscOptionsName("-dmmg_monitor_all","Monitor DMMG iterations","DMMG",&monitorAll);
574:     PetscOptionsTruth("-dmmg_fas","Use the Full Approximation Scheme","DMMGSetSNES",useFAS,&useFAS,PETSC_NULL);
575:     PetscOptionsName("-dmmg_fas_block","Use point-block smoothing","DMMG",&fasBlock);
576:     PetscOptionsName("-dmmg_fas_ngmres","Use Nonlinear GMRES","DMMG",&fasGMRES);
577:     PetscOptionsInt("-dmmg_fas_presmooth","Number of downward smoother iterates","DMMG",fasPresmooth,&fasPresmooth,PETSC_NULL);
578:     PetscOptionsInt("-dmmg_fas_postsmooth","Number of upward smoother iterates","DMMG",fasPostsmooth,&fasPostsmooth,PETSC_NULL);
579:     PetscOptionsInt("-dmmg_fas_coarsesmooth","Number of coarse smoother iterates","DMMG",fasCoarsesmooth,&fasCoarsesmooth,PETSC_NULL);
580:     PetscOptionsReal("-dmmg_fas_rtol","Relative tolerance for FAS","DMMG",fasRtol,&fasRtol,PETSC_NULL);
581:     PetscOptionsReal("-dmmg_fas_atol","Absolute tolerance for FAS","DMMG",fasAbstol,&fasAbstol,PETSC_NULL);
582:     PetscOptionsInt("-dmmg_fas_max_its","Maximum number of iterates per smoother","DMMG",fasMaxIter,&fasMaxIter,PETSC_NULL);

584:     PetscOptionsTruth("-dmmg_coloring_from_mat","Compute the coloring directly from the matrix nonzero structure","DMMGSetSNES",dmmg[0]->getcoloringfrommat,&dmmg[0]->getcoloringfrommat,PETSC_NULL);

586:     PetscOptionsName("-dmmg_jacobian_fd","Compute sparse Jacobian explicitly with finite differencing","DMMGSetSNES",&fdjacobian);
587:     if (fdjacobian) jacobian = DMMGComputeJacobianWithFD;
588: #if defined(PETSC_HAVE_ADIC)
589:     PetscOptionsName("-dmmg_jacobian_ad","Compute sparse Jacobian explicitly with ADIC (automatic differentiation)","DMMGSetSNES",&adjacobian);
590:     if (adjacobian) jacobian = DMMGComputeJacobianWithAdic;
591: #endif

593:     PetscOptionsTruthGroupBegin("-dmmg_jacobian_mf_fd_operator","Apply Jacobian via matrix free finite differencing","DMMGSetSNES",&mffdoperator);
594:     PetscOptionsTruthGroupEnd("-dmmg_jacobian_mf_fd","Apply Jacobian via matrix free finite differencing even in computing preconditioner","DMMGSetSNES",&mffd);
595:     if (mffd) mffdoperator = PETSC_TRUE;
596: #if defined(PETSC_HAVE_ADIC)
597:     PetscOptionsTruthGroupBegin("-dmmg_jacobian_mf_ad_operator","Apply Jacobian via matrix free ADIC (automatic differentiation)","DMMGSetSNES",&mfadoperator);
598:     PetscOptionsTruthGroupEnd("-dmmg_jacobian_mf_ad","Apply Jacobian via matrix free ADIC (automatic differentiation) even in computing preconditioner","DMMGSetSNES",&mfad);
599:     if (mfad) mfadoperator = PETSC_TRUE;
600: #endif
601:     PetscOptionsEnum("-dmmg_iscoloring_type","Type of ISColoring","None",ISColoringTypes,(PetscEnum)dmmg[0]->isctype,(PetscEnum*)&dmmg[0]->isctype,PETSC_NULL);
602: 
603:   PetscOptionsEnd();

605:   /* create solvers for each level */
606:   for (i=0; i<nlevels; i++) {
607:     SNESCreate(dmmg[i]->comm,&dmmg[i]->snes);
608:     PetscObjectIncrementTabLevel((PetscObject)dmmg[i]->snes,PETSC_NULL,nlevels - i - 1);
609:     SNESSetFunction(dmmg[i]->snes,dmmg[i]->b,function,dmmg[i]);
610:     SNESSetOptionsPrefix(dmmg[i]->snes,dmmg[i]->prefix);
611:     SNESGetKSP(dmmg[i]->snes,&dmmg[i]->ksp);

613:     if (mffdoperator) {
614:       MatCreateSNESMF(dmmg[i]->snes,&dmmg[i]->J);
615:       VecDuplicate(dmmg[i]->x,&dmmg[i]->work1);
616:       VecDuplicate(dmmg[i]->x,&dmmg[i]->work2);
617:       MatMFFDSetFunction(dmmg[i]->J,(PetscErrorCode (*)(void*, Vec,Vec))SNESComputeFunction,dmmg[i]->snes);
618:       if (mffd) {
619:         dmmg[i]->B = dmmg[i]->J;
620:         jacobian   = DMMGComputeJacobianWithMF;
621:         PetscObjectReference((PetscObject) dmmg[i]->B);
622:       }
623: #if defined(PETSC_HAVE_ADIC)
624:     } else if (mfadoperator) {
625:       MatRegisterDAAD();
626:       MatCreateDAAD((DA)dmmg[i]->dm,&dmmg[i]->J);
627:       MatDAADSetCtx(dmmg[i]->J,dmmg[i]->user);
628:       if (mfad) {
629:         dmmg[i]->B = dmmg[i]->J;
630:         jacobian   = DMMGComputeJacobianWithMF;
631:         PetscObjectReference((PetscObject) dmmg[i]->B);
632:       }
633: #endif
634:     }
635: 
636:     if (!useFAS) {
637:       if (!dmmg[i]->B) {
638:         DMGetMatrix(dmmg[i]->dm,dmmg[i]->mtype,&dmmg[i]->B);
639:       }
640:       if (!dmmg[i]->J) {
641:         dmmg[i]->J = dmmg[i]->B;
642:         PetscObjectReference((PetscObject) dmmg[i]->B);
643:       }
644:     }

646:     DMMGSetUpLevel(dmmg,dmmg[i]->ksp,i+1);
647: 
648:     /*
649:        if the number of levels is > 1 then we want the coarse solve in the grid sequencing to use LU
650:        when possible 
651:     */
652:     if (nlevels > 1 && i == 0) {
653:       PC         pc;
654:       KSP        cksp;
655:       PetscTruth flg1,flg2,flg3;

657:       KSPGetPC(dmmg[i]->ksp,&pc);
658:       PCMGGetCoarseSolve(pc,&cksp);
659:       KSPGetPC(cksp,&pc);
660:       PetscTypeCompare((PetscObject)pc,PCILU,&flg1);
661:       PetscTypeCompare((PetscObject)pc,PCSOR,&flg2);
662:       PetscTypeCompare((PetscObject)pc,PETSC_NULL,&flg3);
663:       if (flg1 || flg2 || flg3) {
664:         PCSetType(pc,PCLU);
665:       }
666:     }

668:     dmmg[i]->computejacobian = jacobian;
669:     dmmg[i]->computefunction = function;
670:     if (useFAS) {
671:       if (cookie == DM_COOKIE) {
672: #if defined(PETSC_HAVE_ADIC)
673:         if (fasBlock) {
674:           dmmg[i]->solve     = DMMGSolveFASb;
675:         } else if(fasGMRES) {
676:           dmmg[i]->solve     = DMMGSolveFAS_NGMRES;
677:         } else {
678:           dmmg[i]->solve     = DMMGSolveFAS4;
679:         }
680: #else
681:         SETERRQ(PETSC_ERR_SUP, "Must use ADIC for structured FAS.");
682: #endif
683:       } else {
684: #if defined(PETSC_HAVE_SIEVE)
685:         dmmg[i]->solve       = DMMGSolveFAS_Mesh;
686: #endif
687:       }
688:     } else {
689:       dmmg[i]->solve         = DMMGSolveSNES;
690:     }
691:   }

693:   if (jacobian == DMMGComputeJacobianWithFD) {
694:     ISColoring iscoloring;

696:     for (i=0; i<nlevels; i++) {
697:       if (dmmg[0]->getcoloringfrommat) {
698:         MatGetColoring(dmmg[i]->B,(MatColoringType)MATCOLORING_SL,&iscoloring);
699:       } else {
700:         DMGetColoring(dmmg[i]->dm,dmmg[0]->isctype,&iscoloring);
701:       }
702:       MatFDColoringCreate(dmmg[i]->B,iscoloring,&dmmg[i]->fdcoloring);
703:       ISColoringDestroy(iscoloring);
704:       if (function == DMMGFormFunction) function = DMMGFormFunctionFD;
705:       MatFDColoringSetFunction(dmmg[i]->fdcoloring,(PetscErrorCode(*)(void))function,dmmg[i]);
706:       MatFDColoringSetFromOptions(dmmg[i]->fdcoloring);
707:     }
708: #if defined(PETSC_HAVE_ADIC)
709:   } else if (jacobian == DMMGComputeJacobianWithAdic) {
710:     for (i=0; i<nlevels; i++) {
711:       ISColoring iscoloring;
712:       DMGetColoring(dmmg[i]->dm,IS_COLORING_GHOSTED,&iscoloring);
713:       MatSetColoring(dmmg[i]->B,iscoloring);
714:       ISColoringDestroy(iscoloring);
715:     }
716: #endif
717:   }
718:   if (!useFAS) {
719:     for (i=0; i<nlevels; i++) {
720:       SNESSetJacobian(dmmg[i]->snes,dmmg[i]->J,dmmg[i]->B,DMMGComputeJacobian_Multigrid,dmmg);
721:     }
722:   }

724:   /* Create interpolation scaling */
725:   for (i=1; i<nlevels; i++) {
726:     DMGetInterpolationScale(dmmg[i-1]->dm,dmmg[i]->dm,dmmg[i]->R,&dmmg[i]->Rscale);
727:   }

729:   if (useFAS) {
730:     PetscTruth flg;

732:     for(i = 0; i < nlevels; i++) {
733:       VecDuplicate(dmmg[i]->b,&dmmg[i]->w);
734:       if (cookie == DM_COOKIE) {
735: #if defined(PETSC_HAVE_ADIC)
736:         NLFCreate_DAAD(&dmmg[i]->nlf);
737:         NLFDAADSetDA_DAAD(dmmg[i]->nlf,(DA)dmmg[i]->dm);
738:         NLFDAADSetCtx_DAAD(dmmg[i]->nlf,dmmg[i]->user);
739:         NLFDAADSetResidual_DAAD(dmmg[i]->nlf,dmmg[i]->r);
740:         NLFDAADSetNewtonIterations_DAAD(dmmg[i]->nlf,fasMaxIter);
741: #endif
742:       } else {
743:       }

745:       dmmg[i]->monitor      = monitor;
746:       dmmg[i]->monitorall   = monitorAll;
747:       dmmg[i]->presmooth    = fasPresmooth;
748:       dmmg[i]->postsmooth   = fasPostsmooth;
749:       dmmg[i]->coarsesmooth = fasCoarsesmooth;
750:       dmmg[i]->rtol         = fasRtol;
751:       dmmg[i]->abstol       = fasAbstol;
752:     }

754:     PetscOptionsHasName(0, "-dmmg_fas_view", &flg);
755:     if (flg) {
756:       for(i = 0; i < nlevels; i++) {
757:         if (i == 0) {
758:           PetscPrintf(dmmg[i]->comm,"FAS Solver Parameters\n");
759:           PetscPrintf(dmmg[i]->comm,"  rtol %G atol %G\n",dmmg[i]->rtol,dmmg[i]->abstol);
760:           PetscPrintf(dmmg[i]->comm,"             coarsesmooths %D\n",dmmg[i]->coarsesmooth);
761:           PetscPrintf(dmmg[i]->comm,"             Newton iterations %D\n",fasMaxIter);
762:         } else {
763:           PetscPrintf(dmmg[i]->comm,"  level %D   presmooths    %D\n",i,dmmg[i]->presmooth);
764:           PetscPrintf(dmmg[i]->comm,"             postsmooths   %D\n",dmmg[i]->postsmooth);
765:           PetscPrintf(dmmg[i]->comm,"             Newton iterations %D\n",fasMaxIter);
766:         }
767:         if (fasBlock) {
768:           PetscPrintf(dmmg[i]->comm,"  using point-block smoothing\n");
769:         } else if(fasGMRES) {
770:           PetscPrintf(dmmg[i]->comm,"  using non-linear gmres\n");
771:         }
772:       }
773:     }
774:   }
775:   return(0);
776: }

780: /*@C
781:     DMMGSetFromOptions - Sets various options associated with the DMMG object

783:     Collective on DMMG

785:     Input Parameter:
786: .   dmmg - the context


789:     Notes: this is currently only supported for use with DMMGSetSNES() NOT DMMGSetKSP()

791:            Most options are checked in DMMGSetSNES(); this does call SNESSetFromOptions()

793:     Level: advanced

795: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNESLocal(), DMMGSetSNES()

797: @*/
798: PetscErrorCode  DMMGSetFromOptions(DMMG *dmmg)
799: {
800:   PetscErrorCode          ierr;
801:   PetscInt                i,nlevels = dmmg[0]->nlevels;

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

806:   for (i=0; i<nlevels; i++) {
807:     SNESSetFromOptions(dmmg[i]->snes);
808:   }
809:   return(0);
810: }

814: /*@C
815:     DMMGSetSNESLocalFD - Sets the local user function that is used to approximately compute the Jacobian
816:         via finite differences.

818:     Collective on DMMG

820:     Input Parameter:
821: +   dmmg - the context
822: -   function - the function that defines the nonlinear system

824:     Level: intermediate

826: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetSNESLocal()

828: @*/
829: PetscErrorCode DMMGSetSNESLocalFD(DMMG *dmmg,DALocalFunction1 function)
830: {
831:   PetscInt       i,nlevels = dmmg[0]->nlevels;

834:   for (i=0; i<nlevels; i++) {
835:     dmmg[i]->lfj = (PetscErrorCode (*)(void))function;
836:   }
837:   return(0);
838: }


843: /*@C
844:   DMMGGetSNESLocal - Returns the local functions for residual and Jacobian evaluation.

846:   Collective on DMMG

848:   Input Parameter:
849: . dmmg - the context

851:   Output Parameters:
852: + function - the function that defines the nonlinear system
853: - jacobian - function defines the local part of the Jacobian

855:   Level: intermediate

857: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetSNESLocal()
858: @*/
859: PetscErrorCode DMMGGetSNESLocal(DMMG *dmmg,DALocalFunction1 *function, DALocalFunction1 *jacobian)
860: {
861:   PetscCookie    cookie;

865:   PetscObjectGetCookie((PetscObject) dmmg[0]->dm, &cookie);
866:   if (cookie == DM_COOKIE) {
867:     DAGetLocalFunction((DA) dmmg[0]->dm, function);
868:     DAGetLocalJacobian((DA) dmmg[0]->dm, jacobian);
869:   } else {
870: #ifdef PETSC_HAVE_SIEVE
871:     MeshGetLocalFunction((Mesh) dmmg[0]->dm, (PetscErrorCode (**)(Mesh,SectionReal,SectionReal,void*)) function);
872:     MeshGetLocalJacobian((Mesh) dmmg[0]->dm, (PetscErrorCode (**)(Mesh,SectionReal,Mat,void*)) jacobian);
873: #else
874:     SETERRQ(PETSC_ERR_SUP, "Unstructured grids only supported when Sieve is enabled.\nReconfigure with --with-sieve.");
875: #endif
876:   }
877:   return(0);
878: }

880: /*MC
881:     DMMGSetSNESLocal - Sets the local user function that defines the nonlinear set of equations
882:     that will use the grid hierarchy and (optionally) its derivative.

884:     Collective on DMMG

886:    Synopsis:
887:    PetscErrorCode DMMGSetSNESLocal(DMMG *dmmg,DALocalFunction1 function, DALocalFunction1 jacobian,
888:                         DALocalFunction1 ad_function, DALocalFunction1 admf_function);

890:     Input Parameter:
891: +   dmmg - the context
892: .   function - the function that defines the nonlinear system
893: .   jacobian - function defines the local part of the Jacobian
894: .   ad_function - the name of the function with an ad_ prefix. This is ignored currently
895: -   admf_function - the name of the function with an ad_ prefix. This is ignored currently

897:     Options Database Keys:
898: +    -dmmg_jacobian_fd
899: .    -dmmg_jacobian_ad
900: .    -dmmg_jacobian_mf_fd_operator
901: .    -dmmg_jacobian_mf_fd
902: .    -dmmg_jacobian_mf_ad_operator
903: -    -dmmg_jacobian_mf_ad


906:     Level: intermediate

908:     Notes: 
909:     If the Jacobian is not provided, this routine
910:     uses finite differencing to approximate the Jacobian.

912: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES()

914: M*/

918: PetscErrorCode DMMGSetSNESLocal_Private(DMMG *dmmg,DALocalFunction1 function,DALocalFunction1 jacobian,DALocalFunction1 ad_function,DALocalFunction1 admf_function)
919: {
921:   PetscInt       i,nlevels = dmmg[0]->nlevels;
922:   PetscCookie    cookie;
923:   PetscErrorCode (*computejacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*) = 0;


927:   if (jacobian)         computejacobian = DMMGComputeJacobian;
928: #if defined(PETSC_HAVE_ADIC)
929:   else if (ad_function) computejacobian = DMMGComputeJacobianWithAdic;
930: #endif
931:   CHKMEMQ;
932:   PetscObjectGetCookie((PetscObject) dmmg[0]->dm,&cookie);
933:   if (cookie == DM_COOKIE) {
934:     PetscTruth flag;
935:     /* it makes no sense to use an option to decide on ghost, it depends on whether the 
936:        formfunctionlocal computes ghost values in F or not. */
937:     PetscOptionsHasName(PETSC_NULL, "-dmmg_form_function_ghost", &flag);
938:     if (flag) {
939:       DMMGSetSNES(dmmg,DMMGFormFunctionGhost,computejacobian);
940:     } else {
941:       DMMGSetSNES(dmmg,DMMGFormFunction,computejacobian);
942:     }
943:     for (i=0; i<nlevels; i++) {
944:       DASetLocalFunction((DA)dmmg[i]->dm,function);
945:       dmmg[i]->lfj = (PetscErrorCode (*)(void))function;
946:       DASetLocalJacobian((DA)dmmg[i]->dm,jacobian);
947:       DASetLocalAdicFunction((DA)dmmg[i]->dm,ad_function);
948:       DASetLocalAdicMFFunction((DA)dmmg[i]->dm,admf_function);
949:     }
950:     CHKMEMQ;
951:   } else {
952: #ifdef PETSC_HAVE_SIEVE
953:     DMMGSetSNES(dmmg, DMMGFormFunctionMesh, DMMGComputeJacobianMesh);
954:     for (i=0; i<nlevels; i++) {
955:       MeshSetLocalFunction((Mesh) dmmg[i]->dm, (PetscErrorCode (*)(Mesh,SectionReal,SectionReal,void*)) function);
956:       dmmg[i]->lfj = (PetscErrorCode (*)(void)) function;
957:       MeshSetLocalJacobian((Mesh) dmmg[i]->dm, (PetscErrorCode (*)(Mesh,SectionReal,Mat,void*)) jacobian);
958:       // Setup a work section
959:       SectionReal defaultSec, constantSec;
960:       PetscTruth  hasConstant;

962:       MeshGetSectionReal((Mesh) dmmg[i]->dm, "default", &defaultSec);
963:       MeshHasSectionReal((Mesh) dmmg[i]->dm, "constant", &hasConstant);
964:       if (!hasConstant) {
965:         SectionRealDuplicate(defaultSec, &constantSec);
966:         PetscObjectSetName((PetscObject) constantSec, "constant");
967:         MeshSetSectionReal((Mesh) dmmg[i]->dm, constantSec);
968:         SectionRealDestroy(constantSec);
969:       }
970:     }
971:     CHKMEMQ;
972: #else
973:     SETERRQ(PETSC_ERR_SUP, "Unstructured grids only supported when Sieve is enabled.\nReconfigure with --with-sieve.");
974: #endif
975:   }
976:   CHKMEMQ;
977:   return(0);
978: }

982: PetscErrorCode DMMGFunctioni(void* ctx,PetscInt i,Vec u,PetscScalar* r)
983: {
984:   DMMG           dmmg = (DMMG)ctx;
985:   Vec            U = dmmg->lwork1;
987:   VecScatter     gtol;

990:   /* copy u into interior part of U */
991:   DAGetScatter((DA)dmmg->dm,0,&gtol,0);
992:   VecScatterBegin(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
993:   VecScatterEnd(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
994:   DAFormFunctioni1((DA)dmmg->dm,i,U,r,dmmg->user);
995:   return(0);
996: }

1000: PetscErrorCode DMMGFunctionib(PetscInt i,Vec u,PetscScalar* r,void* ctx)
1001: {
1002:   DMMG           dmmg = (DMMG)ctx;
1003:   Vec            U = dmmg->lwork1;
1005:   VecScatter     gtol;

1008:   /* copy u into interior part of U */
1009:   DAGetScatter((DA)dmmg->dm,0,&gtol,0);
1010:   VecScatterBegin(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
1011:   VecScatterEnd(gtol,u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL);
1012:   DAFormFunctionib1((DA)dmmg->dm,i,U,r,dmmg->user);
1013:   return(0);
1014: }

1018: PetscErrorCode DMMGFunctioniBase(void* ctx,Vec u)
1019: {
1020:   DMMG           dmmg = (DMMG)ctx;
1021:   Vec            U = dmmg->lwork1;

1025:   DAGlobalToLocalBegin((DA)dmmg->dm,u,INSERT_VALUES,U);
1026:   DAGlobalToLocalEnd((DA)dmmg->dm,u,INSERT_VALUES,U);
1027:   return(0);
1028: }

1032: PetscErrorCode DMMGSetSNESLocali_Private(DMMG *dmmg,PetscErrorCode (*functioni)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*),PetscErrorCode (*adi)(DALocalInfo*,MatStencil*,void*,void*,void*),PetscErrorCode (*adimf)(DALocalInfo*,MatStencil*,void*,void*,void*))
1033: {
1035:   PetscInt       i,nlevels = dmmg[0]->nlevels;

1038:   for (i=0; i<nlevels; i++) {
1039:     DASetLocalFunctioni((DA)dmmg[i]->dm,functioni);
1040:     DASetLocalAdicFunctioni((DA)dmmg[i]->dm,adi);
1041:     DASetLocalAdicMFFunctioni((DA)dmmg[i]->dm,adimf);
1042:     MatMFFDSetFunctioni(dmmg[i]->J,DMMGFunctioni);
1043:     MatMFFDSetFunctioniBase(dmmg[i]->J,DMMGFunctioniBase);
1044:     if (!dmmg[i]->lwork1) {
1045:       DACreateLocalVector((DA)dmmg[i]->dm,&dmmg[i]->lwork1);
1046:     }
1047:   }
1048:   return(0);
1049: }

1053: PetscErrorCode DMMGSetSNESLocalib_Private(DMMG *dmmg,PetscErrorCode (*functioni)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*),PetscErrorCode (*adi)(DALocalInfo*,MatStencil*,void*,void*,void*),PetscErrorCode (*adimf)(DALocalInfo*,MatStencil*,void*,void*,void*))
1054: {
1056:   PetscInt       i,nlevels = dmmg[0]->nlevels;

1059:   for (i=0; i<nlevels; i++) {
1060:     DASetLocalFunctionib((DA)dmmg[i]->dm,functioni);
1061:     DASetLocalAdicFunctionib((DA)dmmg[i]->dm,adi);
1062:     DASetLocalAdicMFFunctionib((DA)dmmg[i]->dm,adimf);
1063:     if (!dmmg[i]->lwork1) {
1064:       DACreateLocalVector((DA)dmmg[i]->dm,&dmmg[i]->lwork1);
1065:     }
1066:   }
1067:   return(0);
1068: }

1070: static PetscErrorCode (*localfunc)(void) = 0;

1074: /*
1075:     Uses the DM object to call the user provided function with the correct calling
1076:   sequence.
1077: */
1078: PetscErrorCode  DMMGInitialGuess_Local(DMMG dmmg,Vec x)
1079: {

1083:   (*dmmg->dm->ops->forminitialguess)(dmmg->dm,localfunc,x,0);
1084:   return(0);
1085: }

1089: /*@C
1090:     DMMGSetInitialGuessLocal - sets code to compute the initial guess for each level

1092:     Collective on DMMG

1094:     Input Parameter:
1095: +   dmmg - the context
1096: -   localguess - the function that computes the initial guess

1098:     Level: intermediate

1100: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess(), DMMGSetSNESLocal()

1102: @*/
1103: PetscErrorCode DMMGSetInitialGuessLocal(DMMG *dmmg,PetscErrorCode (*localguess)(void))
1104: {

1108:   localfunc = localguess;  /* stash into ugly static for now */

1110:   DMMGSetInitialGuess(dmmg,DMMGInitialGuess_Local);
1111:   return(0);
1112: }

1116: /*@C
1117:     DMMGSetISColoringType - type of coloring used to compute Jacobian via finite differencing

1119:     Collective on DMMG

1121:     Input Parameter:
1122: +   dmmg - the context
1123: -   isctype - IS_COLORING_GHOSTED (default) or IS_COLORING_GLOBAL

1125:     Options Database:
1126: .   -dmmg_iscoloring_type <ghosted or global>

1128:     Notes: ghosted coloring requires using DMMGSetSNESLocal()

1130:     Level: intermediate

1132: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNES(), DMMGSetInitialGuess(), DMMGSetSNESLocal()

1134: @*/
1135: PetscErrorCode DMMGSetISColoringType(DMMG *dmmg,ISColoringType isctype)
1136: {
1138:   dmmg[0]->isctype = isctype;
1139:   return(0);
1140: }


1145: /*@C
1146:     DMMGSetUp - Called after DMMGSetSNES() and (optionally) DMMGSetFromOptions()

1148:     Collective on DMMG

1150:     Input Parameter:
1151: .   dmmg - the context

1153:     Notes: Currently this must be called by the user (if they want to). It checks to see if fieldsplit preconditioner
1154:            is being used and manages it.

1156:     Level: advanced

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

1160: @*/
1161: PetscErrorCode  DMMGSetUp(DMMG *dmmg)
1162: {
1164:   PetscInt       i,nDM;
1165:   PetscTruth     fieldsplit,dmcomposite;
1166:   KSP            ksp;
1167:   PC             pc;
1168:   IS             *fields;

1171:   SNESGetKSP(DMMGGetSNES(dmmg),&ksp);
1172:   KSPGetPC(ksp,&pc);
1173:   PetscTypeCompare((PetscObject)pc,PCFIELDSPLIT,&fieldsplit);
1174:   PetscTypeCompare((PetscObject)DMMGGetDM(dmmg),"DMComposite",&dmcomposite);
1175:   if (fieldsplit && dmcomposite) {
1176:     PetscInfo(ksp,"Setting up physics based fieldsplit preconditioner\n");
1177:     DMCompositeGetNumberDM((DMComposite)DMMGGetDM(dmmg),&nDM);
1178:     DMCompositeGetGlobalISs((DMComposite)DMMGGetDM(dmmg),&fields);
1179:     for (i=0; i<nDM; i++) {
1180:       PCFieldSplitSetIS(pc,fields[i]);
1181:       ISDestroy(fields[i]);
1182:     }
1183:     PetscFree(fields);
1184:   }

1186:   return(0);
1187: }