Actual source code: damgsnes.c

  1: 
 2:  #include petscda.h
 3:  #include petscmg.h


  6: /*
  7:       period of -1 indicates update only on zeroth iteration of SNES
  8: */
  9: #define ShouldUpdate(l,it) (((dmmg[l-1]->updatejacobianperiod == -1) && (it == 0)) || \
 10:                             ((dmmg[l-1]->updatejacobianperiod >   0) && !(it % dmmg[l-1]->updatejacobianperiod)))
 11: /*
 12:    Evaluates the Jacobian on all of the grids. It is used by DMMG to provide the 
 13:    ComputeJacobian() function that SNESSetJacobian() requires.
 14: */
 17: PetscErrorCode DMMGComputeJacobian_Multigrid(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
 18: {
 19:   DMMG           *dmmg = (DMMG*)ptr;
 21:   PetscInt       i,nlevels = dmmg[0]->nlevels,it;
 22:   KSP            ksp,lksp;
 23:   PC             pc;
 24:   PetscTruth     ismg;
 25:   Vec            W;
 26:   MatStructure   flg;

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

 32:   /* compute Jacobian on finest grid */
 33:   if (dmmg[nlevels-1]->updatejacobian && ShouldUpdate(nlevels,it)) {
 34:     (*DMMGGetFine(dmmg)->computejacobian)(snes,X,J,B,flag,DMMGGetFine(dmmg));
 35:   } else {
 36:     PetscLogInfo(0,"DMMGComputeJacobian_Multigrid:Skipping Jacobian, SNES iteration %D frequence %D level %D\n",it,dmmg[nlevels-1]->updatejacobianperiod,nlevels-1);
 37:     *flag = SAME_PRECONDITIONER;
 38:   }
 39:   MatSNESMFSetBase(DMMGGetFine(dmmg)->J,X);

 41:   /* create coarser grid Jacobians for preconditioner if multigrid is the preconditioner */
 42:   SNESGetKSP(snes,&ksp);
 43:   KSPGetPC(ksp,&pc);
 44:   PetscTypeCompare((PetscObject)pc,PCMG,&ismg);
 45:   if (ismg) {

 47:     MGGetSmoother(pc,nlevels-1,&lksp);
 48:     KSPSetOperators(lksp,DMMGGetFine(dmmg)->J,DMMGGetFine(dmmg)->B,*flag);

 50:     if (dmmg[0]->galerkin) {
 51:       for (i=nlevels-2; i>-1; i--) {
 52:         PetscTruth JeqB = (PetscTruth)( dmmg[i]->B == dmmg[i]->J);
 53:         MatDestroy(dmmg[i]->B);
 54:         MatPtAP(dmmg[i+1]->B,dmmg[i+1]->R,MAT_INITIAL_MATRIX,1.0,&dmmg[i]->B);
 55:         if (JeqB) dmmg[i]->J = dmmg[i]->B;
 56:         MGGetSmoother(pc,i,&lksp);
 57:         KSPSetOperators(lksp,dmmg[i]->J,dmmg[i]->B,*flag);
 58:       }
 59:     } else {
 60:       for (i=nlevels-1; i>0; i--) {
 61:         if (!dmmg[i-1]->w) {
 62:           VecDuplicate(dmmg[i-1]->x,&dmmg[i-1]->w);
 63:         }
 64:         W    = dmmg[i-1]->w;
 65:         /* restrict X to coarser grid */
 66:         MatRestrict(dmmg[i]->R,X,W);
 67:         X    = W;
 68:         /* scale to "natural" scaling for that grid */
 69:         VecPointwiseMult(dmmg[i]->Rscale,X,X);
 70:         /* tell the base vector for matrix free multiplies */
 71:         MatSNESMFSetBase(dmmg[i-1]->J,X);
 72:         /* compute Jacobian on coarse grid */
 73:         if (dmmg[i-1]->updatejacobian && ShouldUpdate(i,it)) {
 74:           (*dmmg[i-1]->computejacobian)(snes,X,&dmmg[i-1]->J,&dmmg[i-1]->B,&flg,dmmg[i-1]);
 75:           flg = SAME_NONZERO_PATTERN;
 76:         } else {
 77:           PetscLogInfo(0,"DMMGComputeJacobian_Multigrid:Skipping Jacobian, SNES iteration %D frequence %D level %D\n",it,dmmg[i-1]->updatejacobianperiod,i-1);
 78:           flg = SAME_PRECONDITIONER;
 79:         }
 80:         MGGetSmoother(pc,i-1,&lksp);
 81:         KSPSetOperators(lksp,dmmg[i-1]->J,dmmg[i-1]->B,flg);
 82:       }
 83:     }
 84:   }
 85:   return(0);
 86: }

 88: /* ---------------------------------------------------------------------------*/

 92: /* 
 93:    DMMGFormFunction - This is a universal global FormFunction used by the DMMG code
 94:    when the user provides a local function.

 96:    Input Parameters:
 97: +  snes - the SNES context
 98: .  X - input vector
 99: -  ptr - optional user-defined context, as set by SNESSetFunction()

101:    Output Parameter:
102: .  F - function vector

104:  */
105: PetscErrorCode DMMGFormFunction(SNES snes,Vec X,Vec F,void *ptr)
106: {
107:   DMMG           dmmg = (DMMG)ptr;
109:   Vec            localX;
110:   DA             da = (DA)dmmg->dm;

113:   DAGetLocalVector(da,&localX);
114:   /*
115:      Scatter ghost points to local vector, using the 2-step process
116:         DAGlobalToLocalBegin(), DAGlobalToLocalEnd().
117:   */
118:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
119:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
120:   DAFormFunction1(da,localX,F,dmmg->user);
121:   DARestoreLocalVector(da,&localX);
122:   return(0);
123: }

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

132:    Collective on SNES

134:    Input Parameters:
135: +  snes - the SNES context
136: .  X - input vector
137: .  F - function vector
138: -  ptr - pointer to a structure that must have a DA as its first entry. For example this 
139:          could be a DMMG

141:    Level: intermediate

143: .seealso: DASetLocalFunction(), DASetLocalJacobian(), DASetLocalAdicFunction(), DASetLocalAdicMFFunction(),
144:           SNESSetFunction(), SNESSetJacobian()

146: @*/
147: PetscErrorCode SNESDAFormFunction(SNES snes,Vec X,Vec F,void *ptr)
148: {
150:   Vec            localX;
151:   DA             da = *(DA*)ptr;

154:   DAGetLocalVector(da,&localX);
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:   DAFormFunction1(da,localX,F,ptr);
162:   DARestoreLocalVector(da,&localX);
163:   return(0);
164: }

166: /* ---------------------------------------------------------------------------------------------------------------------------*/

170: PetscErrorCode DMMGComputeJacobianWithFD(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
171: {
173:   DMMG           dmmg = (DMMG)ctx;
174: 
176:   SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,dmmg->fdcoloring);
177:   return(0);
178: }

182: PetscErrorCode DMMGComputeJacobianWithMF(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
183: {
185: 
187:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
188:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
189:   return(0);
190: }

194: /*
195:     DMMGComputeJacobianWithAdic - Evaluates the Jacobian via Adic when the user has provided
196:     a local function evaluation routine.
197: */
198: PetscErrorCode DMMGComputeJacobianWithAdic(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
199: {
200:   DMMG           dmmg = (DMMG) ptr;
202:   Vec            localX;
203:   DA             da = (DA) dmmg->dm;

206:   DAGetLocalVector(da,&localX);
207:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
208:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
209:   DAComputeJacobian1WithAdic(da,localX,*B,dmmg->user);
210:   DARestoreLocalVector(da,&localX);
211:   /* Assemble true Jacobian; if it is different */
212:   if (*J != *B) {
213:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
214:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
215:   }
216:   MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);
217:   *flag = SAME_NONZERO_PATTERN;
218:   return(0);
219: }

223: /*
224:     DMMGComputeJacobian - Evaluates the Jacobian when the user has provided
225:     a local function evaluation routine.
226: */
227: PetscErrorCode DMMGComputeJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
228: {
229:   DMMG           dmmg = (DMMG) ptr;
231:   Vec            localX;
232:   DA             da = (DA) dmmg->dm;

235:   DAGetLocalVector(da,&localX);
236:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
237:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
238:   DAComputeJacobian1(da,localX,*B,dmmg->user);
239:   DARestoreLocalVector(da,&localX);
240:   /* Assemble true Jacobian; if it is different */
241:   if (*J != *B) {
242:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
243:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
244:   }
245:   MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);
246:   *flag = SAME_NONZERO_PATTERN;
247:   return(0);
248: }

252: /*@
253:     SNESDAComputeJacobianWithAdic - This is a universal Jacobian evaluation routine
254:     that may be used with SNESSetJacobian() as long as the user context has a DA as
255:     its first record and DASetLocalAdicFunction() has been called.  

257:    Collective on SNES

259:    Input Parameters:
260: +  snes - the SNES context
261: .  X - input vector
262: .  J - Jacobian
263: .  B - Jacobian used in preconditioner (usally same as J)
264: .  flag - indicates if the matrix changed its structure
265: -  ptr - optional user-defined context, as set by SNESSetFunction()

267:    Level: intermediate

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

271: @*/
272: PetscErrorCode SNESDAComputeJacobianWithAdic(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
273: {
274:   DA             da = *(DA*) ptr;
276:   Vec            localX;

279:   DAGetLocalVector(da,&localX);
280:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
281:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
282:   DAComputeJacobian1WithAdic(da,localX,*B,ptr);
283:   DARestoreLocalVector(da,&localX);
284:   /* Assemble true Jacobian; if it is different */
285:   if (*J != *B) {
286:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
287:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
288:   }
289:   MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);
290:   *flag = SAME_NONZERO_PATTERN;
291:   return(0);
292: }

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

301:    Collective on SNES

303:    Input Parameters:
304: +  snes - the SNES context
305: .  X - input vector
306: .  J - Jacobian
307: .  B - Jacobian used in preconditioner (usally same as J)
308: .  flag - indicates if the matrix changed its structure
309: -  ptr - optional user-defined context, as set by SNESSetFunction()

311:    Level: intermediate

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

315: */
316: PetscErrorCode SNESDAComputeJacobianWithAdifor(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
317: {
318:   DA             da = *(DA*) ptr;
320:   Vec            localX;

323:   DAGetLocalVector(da,&localX);
324:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
325:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
326:   DAComputeJacobian1WithAdifor(da,localX,*B,ptr);
327:   DARestoreLocalVector(da,&localX);
328:   /* Assemble true Jacobian; if it is different */
329:   if (*J != *B) {
330:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
331:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
332:   }
333:   MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);
334:   *flag = SAME_NONZERO_PATTERN;
335:   return(0);
336: }

340: /*
341:    SNESDAComputeJacobian - This is a universal Jacobian evaluation routine for a
342:    locally provided Jacobian.

344:    Collective on SNES

346:    Input Parameters:
347: +  snes - the SNES context
348: .  X - input vector
349: .  J - Jacobian
350: .  B - Jacobian used in preconditioner (usally same as J)
351: .  flag - indicates if the matrix changed its structure
352: -  ptr - optional user-defined context, as set by SNESSetFunction()

354:    Level: intermediate

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

358: */
359: PetscErrorCode SNESDAComputeJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
360: {
361:   DA             da = *(DA*) ptr;
363:   Vec            localX;

366:   DAGetLocalVector(da,&localX);
367:   DAGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
368:   DAGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
369:   DAComputeJacobian1(da,localX,*B,ptr);
370:   DARestoreLocalVector(da,&localX);
371:   /* Assemble true Jacobian; if it is different */
372:   if (*J != *B) {
373:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
374:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
375:   }
376:   MatSetOption(*B,MAT_NEW_NONZERO_LOCATION_ERR);
377:   *flag = SAME_NONZERO_PATTERN;
378:   return(0);
379: }

383: PetscErrorCode DMMGSolveSNES(DMMG *dmmg,PetscInt level)
384: {
386:   PetscInt       nlevels = dmmg[0]->nlevels;

389:   dmmg[0]->nlevels = level+1;
390:   SNESSolve(dmmg[level]->snes,dmmg[level]->x);
391:   dmmg[0]->nlevels = nlevels;
392:   return(0);
393: }

396: EXTERN PetscErrorCode NLFCreate_DAAD(NLF*);
397: EXTERN PetscErrorCode NLFRelax_DAAD(NLF,MatSORType,PetscInt,Vec);
398: EXTERN PetscErrorCode NLFDAADSetDA_DAAD(NLF,DA);
399: EXTERN PetscErrorCode NLFDAADSetCtx_DAAD(NLF,void*);
400: EXTERN PetscErrorCode NLFDAADSetResidual_DAAD(NLF,Vec);
401: EXTERN PetscErrorCode NLFDAADSetNewtonIterations_DAAD(NLF,PetscInt);

404: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
405:  #include src/ksp/pc/impls/mg/mgimpl.h
406: /*
407:           This is pre-beta FAS code. It's design should not be taken seriously!
408: */
411: PetscErrorCode DMMGSolveFAS(DMMG *dmmg,PetscInt level)
412: {
414:   PetscInt       i,j,k;
415:   PetscReal      norm;
416:   PetscScalar    zero = 0.0,mone = -1.0,one = 1.0;
417:   MG             *mg;
418:   PC             pc;

421:   VecSet(&zero,dmmg[level]->r);
422:   for (j=1; j<=level; j++) {
423:     if (!dmmg[j]->inject) {
424:       DMGetInjection(dmmg[j-1]->dm,dmmg[j]->dm,&dmmg[j]->inject);
425:     }
426:   }

428:   KSPGetPC(dmmg[level]->ksp,&pc);
429:   mg   = ((MG*)pc->data);

431:   for (i=0; i<100; i++) {

433:     for (j=level; j>0; j--) {

435:       /* Relax residual_fine - F(x_fine) = 0 */
436:       for (k=0; k<dmmg[j]->presmooth; k++) {
437:         NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
438:       }

440:       /* R*(residual_fine - F(x_fine)) */
441:       DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
442:       VecAYPX(&mone,dmmg[j]->r,dmmg[j]->w);

444:       if (j == level || dmmg[j]->monitorall) {
445:         /* norm( residual_fine - f(x_fine) ) */
446:         VecNorm(dmmg[j]->w,NORM_2,&norm);
447:         if (j == level) {
448:           if (norm < dmmg[level]->abstol) goto theend;
449:           if (i == 0) {
450:             dmmg[level]->rrtol = norm*dmmg[level]->rtol;
451:           } else {
452:             if (norm < dmmg[level]->rrtol) goto theend;
453:           }
454:         }
455:       }

457:       if (dmmg[j]->monitorall) {
458:         for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm,"  ");}
459:         PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
460:       }
461:       MatRestrict(mg[j]->restrct,dmmg[j]->w,dmmg[j-1]->r);
462: 
463:       /* F(R*x_fine) */
464:       VecScatterBegin(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);
465:       VecScatterEnd(dmmg[j]->x,dmmg[j-1]->x,INSERT_VALUES,SCATTER_FORWARD,dmmg[j]->inject);
466:       DMMGFormFunction(0,dmmg[j-1]->x,dmmg[j-1]->w,dmmg[j-1]);

468:       /* residual_coarse = F(R*x_fine) + R*(residual_fine - F(x_fine)) */
469:       VecAYPX(&one,dmmg[j-1]->w,dmmg[j-1]->r);

471:       /* save R*x_fine into b (needed when interpolating compute x back up */
472:       VecCopy(dmmg[j-1]->x,dmmg[j-1]->b);
473:     }

475:     for (j=0; j<dmmg[0]->presmooth; j++) {
476:       NLFRelax_DAAD(dmmg[0]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[0]->x);
477:     }
478:     if (dmmg[0]->monitorall){
479:       DMMGFormFunction(0,dmmg[0]->x,dmmg[0]->w,dmmg[0]);
480:       VecAXPY(&mone,dmmg[0]->r,dmmg[0]->w);
481:       VecNorm(dmmg[0]->w,NORM_2,&norm);
482:       for (k=0; k<level+1; k++) {PetscPrintf(dmmg[0]->comm,"  ");}
483:       PetscPrintf(dmmg[0]->comm,"FAS coarse grid function norm %g\n",norm);
484:     }

486:     for (j=1; j<=level; j++) {
487:       /* x_fine = x_fine + R'*(x_coarse - R*x_fine) */
488:       VecAXPY(&mone,dmmg[j-1]->b,dmmg[j-1]->x);
489:       MatInterpolateAdd(mg[j]->restrct,dmmg[j-1]->x,dmmg[j]->x,dmmg[j]->x);

491:       if (dmmg[j]->monitorall) {
492:         /* norm( F(x_fine) - residual_fine ) */
493:         DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
494:         VecAXPY(&mone,dmmg[j]->r,dmmg[j]->w);
495:         VecNorm(dmmg[j]->w,NORM_2,&norm);
496:         for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm,"  ");}
497:         PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
498:       }

500:       /* Relax residual_fine - F(x_fine)  = 0 */
501:       for (k=0; k<dmmg[j]->postsmooth; k++) {
502:         NLFRelax_DAAD(dmmg[j]->nlf,SOR_SYMMETRIC_SWEEP,1,dmmg[j]->x);
503:       }

505:       if (dmmg[j]->monitorall) {
506:         /* norm( F(x_fine) - residual_fine ) */
507:         DMMGFormFunction(0,dmmg[j]->x,dmmg[j]->w,dmmg[j]);
508:         VecAXPY(&mone,dmmg[j]->r,dmmg[j]->w);
509:         VecNorm(dmmg[j]->w,NORM_2,&norm);
510:         for (k=0; k<level-j+1; k++) {PetscPrintf(dmmg[j]->comm,"  ");}
511:         PetscPrintf(dmmg[j]->comm,"FAS function norm %g\n",norm);
512:       }
513:     }

515:     if (dmmg[level]->monitor){
516:       DMMGFormFunction(0,dmmg[level]->x,dmmg[level]->w,dmmg[level]);
517:       VecNorm(dmmg[level]->w,NORM_2,&norm);
518:       PetscPrintf(dmmg[level]->comm,"%D FAS function norm %g\n",i,norm);
519:     }
520:   }
521:   theend:
522:   return(0);
523: }
524: #endif

526: /* ===========================================================================================================*/

530: /*@C
531:     DMMGSetSNES - Sets the nonlinear function that defines the nonlinear set of equations
532:     to be solved using the grid hierarchy.

534:     Collective on DMMG

536:     Input Parameter:
537: +   dmmg - the context
538: .   function - the function that defines the nonlinear system
539: -   jacobian - optional function to compute Jacobian

541:     Options Database Keys:
542: +    -dmmg_snes_monitor
543: .    -dmmg_jacobian_fd
544: .    -dmmg_jacobian_ad
545: .    -dmmg_jacobian_mf_fd_operator
546: .    -dmmg_jacobian_mf_fd
547: .    -dmmg_jacobian_mf_ad_operator
548: .    -dmmg_jacobian_mf_ad
549: -    -dmmg_jacobian_period <p> - Indicates how often in the SNES solve the Jacobian is recomputed (on all levels)
550:                                  as suggested by Florin Dobrian if p is -1 then Jacobian is computed only on first
551:                                  SNES iteration (i.e. -1 is equivalent to infinity) 

553:     Level: advanced

555: .seealso DMMGCreate(), DMMGDestroy, DMMGSetKSP(), DMMGSetSNESLocal()

557: @*/
558: PetscErrorCode DMMGSetSNES(DMMG *dmmg,PetscErrorCode (*function)(SNES,Vec,Vec,void*),PetscErrorCode (*jacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
559: {
561:   PetscMPIInt    size;
562:   PetscInt       i,nlevels = dmmg[0]->nlevels,period = 1;
563:   PetscTruth     snesmonitor,mffdoperator,mffd,fdjacobian;
564: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
565:   PetscTruth     mfadoperator,mfad,adjacobian;
566: #endif
567:   PetscViewer    ascii;
568:   MPI_Comm       comm;

571:   if (!dmmg)     SETERRQ(PETSC_ERR_ARG_NULL,"Passing null as DMMG");
572:   if (!jacobian) jacobian = DMMGComputeJacobianWithFD;

574:   PetscOptionsBegin(dmmg[0]->comm,PETSC_NULL,"DMMG Options","SNES");
575:     PetscOptionsName("-dmmg_snes_monitor","Monitor nonlinear convergence","SNESSetMonitor",&snesmonitor);


578:     PetscOptionsName("-dmmg_jacobian_fd","Compute sparse Jacobian explicitly with finite differencing","DMMGSetSNES",&fdjacobian);
579:     if (fdjacobian) jacobian = DMMGComputeJacobianWithFD;
580: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
581:     PetscOptionsName("-dmmg_jacobian_ad","Compute sparse Jacobian explicitly with ADIC (automatic differentiation)","DMMGSetSNES",&adjacobian);
582:     if (adjacobian) jacobian = DMMGComputeJacobianWithAdic;
583: #endif

585:     PetscOptionsLogicalGroupBegin("-dmmg_jacobian_mf_fd_operator","Apply Jacobian via matrix free finite differencing","DMMGSetSNES",&mffdoperator);
586:     PetscOptionsLogicalGroupEnd("-dmmg_jacobian_mf_fd","Apply Jacobian via matrix free finite differencing even in computing preconditioner","DMMGSetSNES",&mffd);
587:     if (mffd) mffdoperator = PETSC_TRUE;
588: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
589:     PetscOptionsLogicalGroupBegin("-dmmg_jacobian_mf_ad_operator","Apply Jacobian via matrix free ADIC (automatic differentiation)","DMMGSetSNES",&mfadoperator);
590:     PetscOptionsLogicalGroupEnd("-dmmg_jacobian_mf_ad","Apply Jacobian via matrix free ADIC (automatic differentiation) even in computing preconditioner","DMMGSetSNES",&mfad);
591:     if (mfad) mfadoperator = PETSC_TRUE;
592: #endif
593:   PetscOptionsEnd();

595:   /* create solvers for each level */
596:   for (i=0; i<nlevels; i++) {
597:     SNESCreate(dmmg[i]->comm,&dmmg[i]->snes);
598:     SNESGetKSP(dmmg[i]->snes,&dmmg[i]->ksp);
599:     if (snesmonitor) {
600:       PetscObjectGetComm((PetscObject)dmmg[i]->snes,&comm);
601:       PetscViewerASCIIOpen(comm,"stdout",&ascii);
602:       PetscViewerASCIISetTab(ascii,nlevels-i);
603:       SNESSetMonitor(dmmg[i]->snes,SNESDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);
604:     }

606:     if (mffdoperator) {
607:       MatCreateSNESMF(dmmg[i]->snes,dmmg[i]->x,&dmmg[i]->J);
608:       VecDuplicate(dmmg[i]->x,&dmmg[i]->work1);
609:       VecDuplicate(dmmg[i]->x,&dmmg[i]->work2);
610:       MatSNESMFSetFunction(dmmg[i]->J,dmmg[i]->work1,function,dmmg[i]);
611:       if (mffd) {
612:         dmmg[i]->B = dmmg[i]->J;
613:         jacobian   = DMMGComputeJacobianWithMF;
614:       }
615: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
616:     } else if (mfadoperator) {
617:       MatRegisterDAAD();
618:       MatCreateDAAD((DA)dmmg[i]->dm,&dmmg[i]->J);
619:       MatDAADSetCtx(dmmg[i]->J,dmmg[i]->user);
620:       if (mfad) {
621:         dmmg[i]->B = dmmg[i]->J;
622:         jacobian   = DMMGComputeJacobianWithMF;
623:       }
624: #endif
625:     }
626: 
627:     if (!dmmg[i]->B) {
628:       MPI_Comm_size(dmmg[i]->comm,&size);
629:       DMGetMatrix(dmmg[i]->dm,MATAIJ,&dmmg[i]->B);
630:     }
631:     if (!dmmg[i]->J) {
632:       dmmg[i]->J = dmmg[i]->B;
633:     }

635:     DMMGSetUpLevel(dmmg,dmmg[i]->ksp,i+1);
636: 
637:     /*
638:        if the number of levels is > 1 then we want the coarse solve in the grid sequencing to use LU
639:        when possible 
640:     */
641:     if (nlevels > 1 && i == 0) {
642:       PC         pc;
643:       KSP        cksp;
644:       PetscTruth flg1,flg2,flg3;

646:       KSPGetPC(dmmg[i]->ksp,&pc);
647:       MGGetCoarseSolve(pc,&cksp);
648:       KSPGetPC(cksp,&pc);
649:       PetscTypeCompare((PetscObject)pc,PCILU,&flg1);
650:       PetscTypeCompare((PetscObject)pc,PCSOR,&flg2);
651:       PetscTypeCompare((PetscObject)pc,PETSC_NULL,&flg3);
652:       if (flg1 || flg2 || flg3) {
653:         PCSetType(pc,PCLU);
654:       }
655:     }

657:     dmmg[i]->solve           = DMMGSolveSNES;
658:     dmmg[i]->computejacobian = jacobian;
659:     dmmg[i]->computefunction = function;
660:   }

662:   if (jacobian == DMMGComputeJacobianWithFD) {
663:     ISColoring iscoloring;
664:     for (i=0; i<nlevels; i++) {
665:       DMGetColoring(dmmg[i]->dm,IS_COLORING_LOCAL,&iscoloring);
666:       MatFDColoringCreate(dmmg[i]->B,iscoloring,&dmmg[i]->fdcoloring);
667:       ISColoringDestroy(iscoloring);
668:       MatFDColoringSetFunction(dmmg[i]->fdcoloring,(PetscErrorCode(*)(void))function,dmmg[i]);
669:       MatFDColoringSetFromOptions(dmmg[i]->fdcoloring);
670:     }
671: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
672:   } else if (jacobian == DMMGComputeJacobianWithAdic) {
673:     for (i=0; i<nlevels; i++) {
674:       ISColoring iscoloring;
675:       DMGetColoring(dmmg[i]->dm,IS_COLORING_GHOSTED,&iscoloring);
676:       MatSetColoring(dmmg[i]->B,iscoloring);
677:       ISColoringDestroy(iscoloring);
678:     }
679: #endif
680:   }

682:   for (i=0; i<nlevels; i++) {
683:     SNESSetJacobian(dmmg[i]->snes,dmmg[i]->J,dmmg[i]->B,DMMGComputeJacobian_Multigrid,dmmg);
684:     SNESSetFunction(dmmg[i]->snes,dmmg[i]->b,function,dmmg[i]);
685:     SNESSetFromOptions(dmmg[i]->snes);
686:   }

688:   /* Create interpolation scaling */
689:   for (i=1; i<nlevels; i++) {
690:     DMGetInterpolationScale(dmmg[i-1]->dm,dmmg[i]->dm,dmmg[i]->R,&dmmg[i]->Rscale);
691:   }

693:   PetscOptionsGetInt(PETSC_NULL,"-dmmg_jacobian_period",&period,PETSC_NULL);
694:   for (i=0; i<nlevels; i++) {
695:     dmmg[i]->updatejacobian       = PETSC_TRUE;
696:     dmmg[i]->updatejacobianperiod = period;
697:   }

699: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
700:   {
701:     PetscTruth flg;
702:     PetscOptionsHasName(PETSC_NULL,"-dmmg_fas",&flg);
703:     if (flg) {
704:       PetscInt newton_its;
705:       PetscOptionsHasName(0,"-fas_view",&flg);
706:       for (i=0; i<nlevels; i++) {
707:         NLFCreate_DAAD(&dmmg[i]->nlf);
708:         NLFDAADSetDA_DAAD(dmmg[i]->nlf,(DA)dmmg[i]->dm);
709:         NLFDAADSetCtx_DAAD(dmmg[i]->nlf,dmmg[i]->user);
710:         NLFDAADSetResidual_DAAD(dmmg[i]->nlf,dmmg[i]->r);
711:         VecDuplicate(dmmg[i]->b,&dmmg[i]->w);

713:         dmmg[i]->monitor    = PETSC_FALSE;
714:         PetscOptionsHasName(0,"-dmmg_fas_monitor",&dmmg[i]->monitor);
715:         dmmg[i]->monitorall = PETSC_FALSE;
716:         PetscOptionsHasName(0,"-dmmg_fas_monitor_all",&dmmg[i]->monitorall);
717:         dmmg[i]->presmooth  = 2;
718:         PetscOptionsGetInt(0,"-dmmg_fas_presmooth",&dmmg[i]->presmooth,0);
719:         dmmg[i]->postsmooth = 2;
720:         PetscOptionsGetInt(0,"-dmmg_fas_postsmooth",&dmmg[i]->postsmooth,0);
721:         dmmg[i]->coarsesmooth = 2;
722:         PetscOptionsGetInt(0,"-dmmg_fas_coarsesmooth",&dmmg[i]->coarsesmooth,0);

724:         dmmg[i]->rtol = 1.e-8;
725:         PetscOptionsGetReal(0,"-dmmg_fas_rtol",&dmmg[i]->rtol,0);
726:         dmmg[i]->abstol = 1.e-50;
727:         PetscOptionsGetReal(0,"-dmmg_fas_atol",&dmmg[i]->abstol,0);

729:         newton_its = 2;
730:         PetscOptionsGetInt(0,"-dmmg_fas_newton_its",&newton_its,0);
731:         NLFDAADSetNewtonIterations_DAAD(dmmg[i]->nlf,newton_its);

733:         if (flg) {
734:           if (i == 0) {
735:             PetscPrintf(dmmg[i]->comm,"FAS Solver Parameters\n");
736:             PetscPrintf(dmmg[i]->comm,"  rtol %g atol %g\n",dmmg[i]->rtol,dmmg[i]->abstol);
737:             PetscPrintf(dmmg[i]->comm,"             coarsesmooths %D\n",dmmg[i]->coarsesmooth);
738:             PetscPrintf(dmmg[i]->comm,"             Newton iterations %D\n",newton_its);
739:           } else {
740:             PetscPrintf(dmmg[i]->comm,"  level %D   presmooths    %D\n",i,dmmg[i]->presmooth);
741:             PetscPrintf(dmmg[i]->comm,"             postsmooths   %D\n",dmmg[i]->postsmooth);
742:             PetscPrintf(dmmg[i]->comm,"             Newton iterations %D\n",newton_its);
743:           }
744:         }
745:         dmmg[i]->solve = DMMGSolveFAS;
746:       }
747:     }
748:   }
749: #endif
750: 
751:   return(0);
752: }

754: /*M
755:     DMMGSetSNESLocal - Sets the local user function that defines the nonlinear set of equations
756:     that will use the grid hierarchy and (optionally) its derivative.

758:     Collective on DMMG

760:    Synopsis:
761:    PetscErrorCode DMMGSetSNESLocal(DMMG *dmmg,DALocalFunction1 function, DALocalFunction1 jacobian,
762:                         DALocalFunction1 ad_function, DALocalFunction1 admf_function);

764:     Input Parameter:
765: +   dmmg - the context
766: .   function - the function that defines the nonlinear system
767: .   jacobian - function defines the local part of the Jacobian
768: .   ad_function - the name of the function with an ad_ prefix. This is ignored if ADIC is
769:                   not installed
770: -   admf_function - the name of the function with an ad_ prefix. This is ignored if ADIC is
771:                   not installed

773:     Options Database Keys:
774: +    -dmmg_snes_monitor
775: .    -dmmg_jacobian_fd
776: .    -dmmg_jacobian_ad
777: .    -dmmg_jacobian_mf_fd_operator
778: .    -dmmg_jacobian_mf_fd
779: .    -dmmg_jacobian_mf_ad_operator
780: .    -dmmg_jacobian_mf_ad
781: -    -dmmg_jacobian_period <p> - Indicates how often in the SNES solve the Jacobian is recomputed (on all levels)
782:                                  as suggested by Florin Dobrian if p is -1 then Jacobian is computed only on first
783:                                  SNES iteration (i.e. -1 is equivalent to infinity) 


786:     Level: intermediate

788:     Notes: 
789:     If ADIC or ADIFOR have been installed, this routine can use ADIC or ADIFOR to compute
790:     the derivative; however, that function cannot call other functions except those in
791:     standard C math libraries.

793:     If ADIC/ADIFOR have not been installed and the Jacobian is not provided, this routine
794:     uses finite differencing to approximate the Jacobian.

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

798: M*/

802: PetscErrorCode DMMGSetSNESLocal_Private(DMMG *dmmg,DALocalFunction1 function,DALocalFunction1 jacobian,DALocalFunction1 ad_function,DALocalFunction1 admf_function)
803: {
805:   PetscInt       i,nlevels = dmmg[0]->nlevels;
806:   PetscErrorCode (*computejacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*) = 0;


810:   if (jacobian)         computejacobian = DMMGComputeJacobian;
811: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
812:   else if (ad_function) computejacobian = DMMGComputeJacobianWithAdic;
813: #endif

815:   DMMGSetSNES(dmmg,DMMGFormFunction,computejacobian);
816:   for (i=0; i<nlevels; i++) {
817:     DASetLocalFunction((DA)dmmg[i]->dm,function);
818:     DASetLocalJacobian((DA)dmmg[i]->dm,jacobian);
819:     DASetLocalAdicFunction((DA)dmmg[i]->dm,ad_function);
820:     DASetLocalAdicMFFunction((DA)dmmg[i]->dm,admf_function);
821:   }
822:   return(0);
823: }

827: static PetscErrorCode DMMGFunctioni(PetscInt i,Vec u,PetscScalar* r,void* ctx)
828: {
829:   DMMG           dmmg = (DMMG)ctx;
830:   Vec            U = dmmg->lwork1;
832:   VecScatter     gtol;

835:   /* copy u into interior part of U */
836:   DAGetScatter((DA)dmmg->dm,0,&gtol,0);
837:   VecScatterBegin(u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL,gtol);
838:   VecScatterEnd(u,U,INSERT_VALUES,SCATTER_FORWARD_LOCAL,gtol);
839:   DAFormFunctioni1((DA)dmmg->dm,i,U,r,dmmg->user);
840:   return(0);
841: }

845: static PetscErrorCode DMMGFunctioniBase(Vec u,void* ctx)
846: {
847:   DMMG           dmmg = (DMMG)ctx;
848:   Vec            U = dmmg->lwork1;

852:   DAGlobalToLocalBegin((DA)dmmg->dm,u,INSERT_VALUES,U);
853:   DAGlobalToLocalEnd((DA)dmmg->dm,u,INSERT_VALUES,U);
854:   return(0);
855: }

859: 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*))
860: {
862:   PetscInt       i,nlevels = dmmg[0]->nlevels;

865:   for (i=0; i<nlevels; i++) {
866:     DASetLocalFunctioni((DA)dmmg[i]->dm,functioni);
867:     DASetLocalAdicFunctioni((DA)dmmg[i]->dm,adi);
868:     DASetLocalAdicMFFunctioni((DA)dmmg[i]->dm,adimf);
869:     MatSNESMFSetFunctioni(dmmg[i]->J,DMMGFunctioni);
870:     MatSNESMFSetFunctioniBase(dmmg[i]->J,DMMGFunctioniBase);
871:     DACreateLocalVector((DA)dmmg[i]->dm,&dmmg[i]->lwork1);
872:   }
873:   return(0);
874: }


877: #if defined(PETSC_HAVE_ADIC) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
879: #include "adic/ad_utils.h"

884: PetscErrorCode PetscADView(PetscInt N,PetscInt nc,double *ptr,PetscViewer viewer)
885: {
886:   PetscInt       i,j,nlen  = PetscADGetDerivTypeSize();
887:   char           *cptr = (char*)ptr;
888:   double         *values;

892:   for (i=0; i<N; i++) {
893:     PetscPrintf(PETSC_COMM_SELF,"Element %D value %g derivatives: ",i,*(double*)cptr);
894:     values = PetscADGetGradArray(cptr);
895:     for (j=0; j<nc; j++) {
896:       PetscPrintf(PETSC_COMM_SELF,"%g ",*values++);
897:     }
898:     PetscPrintf(PETSC_COMM_SELF,"\n");
899:     cptr += nlen;
900:   }

902:   return(0);
903: }

905: #endif