Actual source code: snes.c

petsc-master 2020-01-21
Report Typos and Errors
  1:  #include <petsc/private/snesimpl.h>
  2:  #include <petscdmshell.h>
  3:  #include <petscdraw.h>
  4:  #include <petscds.h>
  5:  #include <petscdmadaptor.h>
  6:  #include <petscconvest.h>

  8: PetscBool         SNESRegisterAllCalled = PETSC_FALSE;
  9: PetscFunctionList SNESList              = NULL;

 11: /* Logging support */
 12: PetscClassId  SNES_CLASSID, DMSNES_CLASSID;
 13: PetscLogEvent SNES_Solve, SNES_Setup, SNES_FunctionEval, SNES_JacobianEval, SNES_NGSEval, SNES_NGSFuncEval, SNES_NPCSolve, SNES_ObjectiveEval;

 15: /*@
 16:    SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged.

 18:    Logically Collective on SNES

 20:    Input Parameters:
 21: +  snes - iterative context obtained from SNESCreate()
 22: -  flg - PETSC_TRUE indicates you want the error generated

 24:    Options database keys:
 25: .  -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)

 27:    Level: intermediate

 29:    Notes:
 30:     Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve()
 31:     to determine if it has converged.

 33: .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
 34: @*/
 35: PetscErrorCode  SNESSetErrorIfNotConverged(SNES snes,PetscBool flg)
 36: {
 40:   snes->errorifnotconverged = flg;
 41:   return(0);
 42: }

 44: /*@
 45:    SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge?

 47:    Not Collective

 49:    Input Parameter:
 50: .  snes - iterative context obtained from SNESCreate()

 52:    Output Parameter:
 53: .  flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE

 55:    Level: intermediate

 57: .seealso:  SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
 58: @*/
 59: PetscErrorCode  SNESGetErrorIfNotConverged(SNES snes,PetscBool  *flag)
 60: {
 64:   *flag = snes->errorifnotconverged;
 65:   return(0);
 66: }

 68: /*@
 69:     SNESSetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?

 71:    Logically Collective on SNES

 73:     Input Parameters:
 74: +   snes - the shell SNES
 75: -   flg - is the residual computed?

 77:    Level: advanced

 79: .seealso: SNESGetAlwaysComputesFinalResidual()
 80: @*/
 81: PetscErrorCode  SNESSetAlwaysComputesFinalResidual(SNES snes, PetscBool flg)
 82: {
 85:   snes->alwayscomputesfinalresidual = flg;
 86:   return(0);
 87: }

 89: /*@
 90:     SNESGetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?

 92:    Logically Collective on SNES

 94:     Input Parameter:
 95: .   snes - the shell SNES

 97:     Output Parameter:
 98: .   flg - is the residual computed?

100:    Level: advanced

102: .seealso: SNESSetAlwaysComputesFinalResidual()
103: @*/
104: PetscErrorCode  SNESGetAlwaysComputesFinalResidual(SNES snes, PetscBool *flg)
105: {
108:   *flg = snes->alwayscomputesfinalresidual;
109:   return(0);
110: }

112: /*@
113:    SNESSetFunctionDomainError - tells SNES that the input vector to your SNESFunction is not
114:      in the functions domain. For example, negative pressure.

116:    Logically Collective on SNES

118:    Input Parameters:
119: .  snes - the SNES context

121:    Level: advanced

123: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction
124: @*/
125: PetscErrorCode  SNESSetFunctionDomainError(SNES snes)
126: {
129:   if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates input vector is not in the function domain");
130:   snes->domainerror = PETSC_TRUE;
131:   return(0);
132: }

134: /*@
135:    SNESSetJacobianDomainError - tells SNES that computeJacobian does not make sense any more. For example there is a negative element transformation.

137:    Logically Collective on SNES

139:    Input Parameters:
140: .  snes - the SNES context

142:    Level: advanced

144: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError()
145: @*/
146: PetscErrorCode SNESSetJacobianDomainError(SNES snes)
147: {
150:   if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates computeJacobian does not make sense");
151:   snes->jacobiandomainerror = PETSC_TRUE;
152:   return(0);
153: }

155: /*@
156:    SNESSetCheckJacobianDomainError - if or not to check jacobian domain error after each Jacobian evaluation. By default, we check Jacobian domain error
157:    in the debug mode, and do not check it in the optimized mode.

159:    Logically Collective on SNES

161:    Input Parameters:
162: +  snes - the SNES context
163: -  flg  - indicates if or not to check jacobian domain error after each Jacobian evaluation

165:    Level: advanced

167: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESGetCheckJacobianDomainError()
168: @*/
169: PetscErrorCode SNESSetCheckJacobianDomainError(SNES snes, PetscBool flg)
170: {
173:   snes->checkjacdomainerror = flg;
174:   return(0);
175: }

177: /*@
178:    SNESGetCheckJacobianDomainError - Get an indicator whether or not we are checking Jacobian domain errors after each Jacobian evaluation.

180:    Logically Collective on SNES

182:    Input Parameters:
183: .  snes - the SNES context

185:    Output Parameters:
186: .  flg  - PETSC_FALSE indicates that we don't check jacobian domain errors after each Jacobian evaluation

188:    Level: advanced

190: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESSetCheckJacobianDomainError()
191: @*/
192: PetscErrorCode SNESGetCheckJacobianDomainError(SNES snes, PetscBool *flg)
193: {
197:   *flg = snes->checkjacdomainerror;
198:   return(0);
199: }

201: /*@
202:    SNESGetFunctionDomainError - Gets the status of the domain error after a call to SNESComputeFunction;

204:    Logically Collective on SNES

206:    Input Parameters:
207: .  snes - the SNES context

209:    Output Parameters:
210: .  domainerror - Set to PETSC_TRUE if there's a domain error; PETSC_FALSE otherwise.

212:    Level: advanced

214: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction()
215: @*/
216: PetscErrorCode  SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror)
217: {
221:   *domainerror = snes->domainerror;
222:   return(0);
223: }

225: /*@
226:    SNESGetJacobianDomainError - Gets the status of the Jacobian domain error after a call to SNESComputeJacobian;

228:    Logically Collective on SNES

230:    Input Parameters:
231: .  snes - the SNES context

233:    Output Parameters:
234: .  domainerror - Set to PETSC_TRUE if there's a jacobian domain error; PETSC_FALSE otherwise.

236:    Level: advanced

238: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction(),SNESGetFunctionDomainError()
239: @*/
240: PetscErrorCode SNESGetJacobianDomainError(SNES snes, PetscBool *domainerror)
241: {
245:   *domainerror = snes->jacobiandomainerror;
246:   return(0);
247: }

249: /*@C
250:   SNESLoad - Loads a SNES that has been stored in binary  with SNESView().

252:   Collective on PetscViewer

254:   Input Parameters:
255: + newdm - the newly loaded SNES, this needs to have been created with SNESCreate() or
256:            some related function before a call to SNESLoad().
257: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()

259:    Level: intermediate

261:   Notes:
262:    The type is determined by the data in the file, any type set into the SNES before this call is ignored.

264:   Notes for advanced users:
265:   Most users should not need to know the details of the binary storage
266:   format, since SNESLoad() and TSView() completely hide these details.
267:   But for anyone who's interested, the standard binary matrix storage
268:   format is
269: .vb
270:      has not yet been determined
271: .ve

273: .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad()
274: @*/
275: PetscErrorCode  SNESLoad(SNES snes, PetscViewer viewer)
276: {
278:   PetscBool      isbinary;
279:   PetscInt       classid;
280:   char           type[256];
281:   KSP            ksp;
282:   DM             dm;
283:   DMSNES         dmsnes;

288:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
289:   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");

291:   PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
292:   if (classid != SNES_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONG,"Not SNES next in file");
293:   PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
294:   SNESSetType(snes, type);
295:   if (snes->ops->load) {
296:     (*snes->ops->load)(snes,viewer);
297:   }
298:   SNESGetDM(snes,&dm);
299:   DMGetDMSNES(dm,&dmsnes);
300:   DMSNESLoad(dmsnes,viewer);
301:   SNESGetKSP(snes,&ksp);
302:   KSPLoad(ksp,viewer);
303:   return(0);
304: }

306:  #include <petscdraw.h>
307: #if defined(PETSC_HAVE_SAWS)
308:  #include <petscviewersaws.h>
309: #endif

311: /*@C
312:    SNESViewFromOptions - View from Options

314:    Collective on SNES

316:    Input Parameters:
317: +  A - the application ordering context
318: .  obj - Optional object
319: -  name - command line option

321:    Level: intermediate
322: .seealso:  SNES, SNESView, PetscObjectViewFromOptions(), SNESCreate()
323: @*/
324: PetscErrorCode  SNESViewFromOptions(SNES A,PetscObject obj,const char name[])
325: {

330:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
331:   return(0);
332: }

334: PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES,Vec,Mat,Mat,void*);

336: /*@C
337:    SNESView - Prints the SNES data structure.

339:    Collective on SNES

341:    Input Parameters:
342: +  SNES - the SNES context
343: -  viewer - visualization context

345:    Options Database Key:
346: .  -snes_view - Calls SNESView() at end of SNESSolve()

348:    Notes:
349:    The available visualization contexts include
350: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
351: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
352:          output where only the first processor opens
353:          the file.  All other processors send their
354:          data to the first processor to print.

356:    The user can open an alternative visualization context with
357:    PetscViewerASCIIOpen() - output to a specified file.

359:    Level: beginner

361: .seealso: PetscViewerASCIIOpen()
362: @*/
363: PetscErrorCode  SNESView(SNES snes,PetscViewer viewer)
364: {
365:   SNESKSPEW      *kctx;
367:   KSP            ksp;
368:   SNESLineSearch linesearch;
369:   PetscBool      iascii,isstring,isbinary,isdraw;
370:   DMSNES         dmsnes;
371: #if defined(PETSC_HAVE_SAWS)
372:   PetscBool      issaws;
373: #endif

377:   if (!viewer) {
378:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
379:   }

383:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
384:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
385:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
386:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
387: #if defined(PETSC_HAVE_SAWS)
388:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
389: #endif
390:   if (iascii) {
391:     SNESNormSchedule normschedule;
392:     DM               dm;
393:     PetscErrorCode   (*cJ)(SNES,Vec,Mat,Mat,void*);
394:     void             *ctx;
395:     const char       *pre = "";

397:     PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);
398:     if (!snes->setupcalled) {
399:       PetscViewerASCIIPrintf(viewer,"  SNES has not been set up so information may be incomplete\n");
400:     }
401:     if (snes->ops->view) {
402:       PetscViewerASCIIPushTab(viewer);
403:       (*snes->ops->view)(snes,viewer);
404:       PetscViewerASCIIPopTab(viewer);
405:     }
406:     PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);
407:     PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);
408:     if (snes->usesksp) {
409:       PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);
410:     }
411:     PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);
412:     SNESGetNormSchedule(snes, &normschedule);
413:     if (normschedule > 0) {PetscViewerASCIIPrintf(viewer,"  norm schedule %s\n",SNESNormSchedules[normschedule]);}
414:     if (snes->gridsequence) {
415:       PetscViewerASCIIPrintf(viewer,"  total number of grid sequence refinements=%D\n",snes->gridsequence);
416:     }
417:     if (snes->ksp_ewconv) {
418:       kctx = (SNESKSPEW*)snes->kspconvctx;
419:       if (kctx) {
420:         PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);
421:         PetscViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);
422:         PetscViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);
423:       }
424:     }
425:     if (snes->lagpreconditioner == -1) {
426:       PetscViewerASCIIPrintf(viewer,"  Preconditioned is never rebuilt\n");
427:     } else if (snes->lagpreconditioner > 1) {
428:       PetscViewerASCIIPrintf(viewer,"  Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);
429:     }
430:     if (snes->lagjacobian == -1) {
431:       PetscViewerASCIIPrintf(viewer,"  Jacobian is never rebuilt\n");
432:     } else if (snes->lagjacobian > 1) {
433:       PetscViewerASCIIPrintf(viewer,"  Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);
434:     }
435:     SNESGetDM(snes,&dm);
436:     DMSNESGetJacobian(dm,&cJ,&ctx);
437:     if (snes->mf_operator) {
438:       PetscViewerASCIIPrintf(viewer,"  Jacobian is applied matrix-free with differencing\n");
439:       pre  = "Preconditioning ";
440:     }
441:     if (cJ == SNESComputeJacobianDefault) {
442:       PetscViewerASCIIPrintf(viewer,"  %sJacobian is built using finite differences one column at a time\n",pre);
443:     } else if (cJ == SNESComputeJacobianDefaultColor) {
444:       PetscViewerASCIIPrintf(viewer,"  %sJacobian is built using finite differences with coloring\n",pre);
445:     /* it slightly breaks data encapsulation for access the DMDA information directly */
446:     } else if (cJ == SNESComputeJacobian_DMDA) {
447:       MatFDColoring fdcoloring;
448:       PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);
449:       if (fdcoloring) {
450:         PetscViewerASCIIPrintf(viewer,"  %sJacobian is built using colored finite differences on a DMDA\n",pre);
451:       } else {
452:         PetscViewerASCIIPrintf(viewer,"  %sJacobian is built using a DMDA local Jacobian\n",pre);
453:       }
454:     } else if (snes->mf) {
455:       PetscViewerASCIIPrintf(viewer,"  Jacobian is applied matrix-free with differencing, no explict Jacobian\n");
456:     }
457:   } else if (isstring) {
458:     const char *type;
459:     SNESGetType(snes,&type);
460:     PetscViewerStringSPrintf(viewer," SNESType: %-7.7s",type);
461:     if (snes->ops->view) {(*snes->ops->view)(snes,viewer);}
462:   } else if (isbinary) {
463:     PetscInt    classid = SNES_FILE_CLASSID;
464:     MPI_Comm    comm;
465:     PetscMPIInt rank;
466:     char        type[256];

468:     PetscObjectGetComm((PetscObject)snes,&comm);
469:     MPI_Comm_rank(comm,&rank);
470:     if (!rank) {
471:       PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
472:       PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));
473:       PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR,PETSC_FALSE);
474:     }
475:     if (snes->ops->view) {
476:       (*snes->ops->view)(snes,viewer);
477:     }
478:   } else if (isdraw) {
479:     PetscDraw draw;
480:     char      str[36];
481:     PetscReal x,y,bottom,h;

483:     PetscViewerDrawGetDraw(viewer,0,&draw);
484:     PetscDrawGetCurrentPoint(draw,&x,&y);
485:     PetscStrncpy(str,"SNES: ",sizeof(str));
486:     PetscStrlcat(str,((PetscObject)snes)->type_name,sizeof(str));
487:     PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);
488:     bottom = y - h;
489:     PetscDrawPushCurrentPoint(draw,x,bottom);
490:     if (snes->ops->view) {
491:       (*snes->ops->view)(snes,viewer);
492:     }
493: #if defined(PETSC_HAVE_SAWS)
494:   } else if (issaws) {
495:     PetscMPIInt rank;
496:     const char *name;

498:     PetscObjectGetName((PetscObject)snes,&name);
499:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
500:     if (!((PetscObject)snes)->amsmem && !rank) {
501:       char       dir[1024];

503:       PetscObjectViewSAWs((PetscObject)snes,viewer);
504:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
505:       PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT));
506:       if (!snes->conv_hist) {
507:         SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);
508:       }
509:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);
510:       PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE));
511:     }
512: #endif
513:   }
514:   if (snes->linesearch) {
515:     SNESGetLineSearch(snes, &linesearch);
516:     PetscViewerASCIIPushTab(viewer);
517:     SNESLineSearchView(linesearch, viewer);
518:     PetscViewerASCIIPopTab(viewer);
519:   }
520:   if (snes->npc && snes->usesnpc) {
521:     PetscViewerASCIIPushTab(viewer);
522:     SNESView(snes->npc, viewer);
523:     PetscViewerASCIIPopTab(viewer);
524:   }
525:   PetscViewerASCIIPushTab(viewer);
526:   DMGetDMSNES(snes->dm,&dmsnes);
527:   DMSNESView(dmsnes, viewer);
528:   PetscViewerASCIIPopTab(viewer);
529:   if (snes->usesksp) {
530:     SNESGetKSP(snes,&ksp);
531:     PetscViewerASCIIPushTab(viewer);
532:     KSPView(ksp,viewer);
533:     PetscViewerASCIIPopTab(viewer);
534:   }
535:   if (isdraw) {
536:     PetscDraw draw;
537:     PetscViewerDrawGetDraw(viewer,0,&draw);
538:     PetscDrawPopCurrentPoint(draw);
539:   }
540:   return(0);
541: }

543: /*
544:   We retain a list of functions that also take SNES command
545:   line options. These are called at the end SNESSetFromOptions()
546: */
547: #define MAXSETFROMOPTIONS 5
548: static PetscInt numberofsetfromoptions;
549: static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);

551: /*@C
552:   SNESAddOptionsChecker - Adds an additional function to check for SNES options.

554:   Not Collective

556:   Input Parameter:
557: . snescheck - function that checks for options

559:   Level: developer

561: .seealso: SNESSetFromOptions()
562: @*/
563: PetscErrorCode  SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
564: {
566:   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
567:   othersetfromoptions[numberofsetfromoptions++] = snescheck;
568:   return(0);
569: }

571: PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);

573: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
574: {
575:   Mat            J;
577:   MatNullSpace   nullsp;


582:   if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
583:     Mat A = snes->jacobian, B = snes->jacobian_pre;
584:     MatCreateVecs(A ? A : B, NULL,&snes->vec_func);
585:   }

587:   if (version == 1) {
588:     MatCreateSNESMF(snes,&J);
589:     MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
590:     MatSetFromOptions(J);
591:   } else if (version == 2) {
592:     if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
593: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
594:     SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);
595: #else
596:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)");
597: #endif
598:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2");

600:   /* attach any user provided null space that was on Amat to the newly created matrix free matrix */
601:   if (snes->jacobian) {
602:     MatGetNullSpace(snes->jacobian,&nullsp);
603:     if (nullsp) {
604:       MatSetNullSpace(J,nullsp);
605:     }
606:   }

608:   PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);
609:   if (hasOperator) {

611:     /* This version replaces the user provided Jacobian matrix with a
612:        matrix-free version but still employs the user-provided preconditioner matrix. */
613:     SNESSetJacobian(snes,J,0,0,0);
614:   } else {
615:     /* This version replaces both the user-provided Jacobian and the user-
616:      provided preconditioner Jacobian with the default matrix free version. */
617:     if ((snes->npcside== PC_LEFT) && snes->npc) {
618:       if (!snes->jacobian){SNESSetJacobian(snes,J,0,0,0);}
619:     } else {
620:       KSP       ksp;
621:       PC        pc;
622:       PetscBool match;

624:       SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,0);
625:       /* Force no preconditioner */
626:       SNESGetKSP(snes,&ksp);
627:       KSPGetPC(ksp,&pc);
628:       PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);
629:       if (!match) {
630:         PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");
631:         PCSetType(pc,PCNONE);
632:       }
633:     }
634:   }
635:   MatDestroy(&J);
636:   return(0);
637: }

639: static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx)
640: {
641:   SNES           snes = (SNES)ctx;
643:   Vec            Xfine,Xfine_named = NULL,Xcoarse;

646:   if (PetscLogPrintInfo) {
647:     PetscInt finelevel,coarselevel,fineclevel,coarseclevel;
648:     DMGetRefineLevel(dmfine,&finelevel);
649:     DMGetCoarsenLevel(dmfine,&fineclevel);
650:     DMGetRefineLevel(dmcoarse,&coarselevel);
651:     DMGetCoarsenLevel(dmcoarse,&coarseclevel);
652:     PetscInfo4(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);
653:   }
654:   if (dmfine == snes->dm) Xfine = snes->vec_sol;
655:   else {
656:     DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
657:     Xfine = Xfine_named;
658:   }
659:   DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
660:   if (Inject) {
661:     MatRestrict(Inject,Xfine,Xcoarse);
662:   } else {
663:     MatRestrict(Restrict,Xfine,Xcoarse);
664:     VecPointwiseMult(Xcoarse,Xcoarse,Rscale);
665:   }
666:   DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
667:   if (Xfine_named) {DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);}
668:   return(0);
669: }

671: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
672: {

676:   DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
677:   return(0);
678: }

680: /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can
681:  * safely call SNESGetDM() in their residual evaluation routine. */
682: static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx)
683: {
684:   SNES           snes = (SNES)ctx;
686:   Vec            X,Xnamed = NULL;
687:   DM             dmsave;
688:   void           *ctxsave;
689:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL;

692:   dmsave = snes->dm;
693:   KSPGetDM(ksp,&snes->dm);
694:   if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */
695:   else {                                     /* We are on a coarser level, this vec was initialized using a DM restrict hook */
696:     DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
697:     X    = Xnamed;
698:     SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);
699:     /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */
700:     if (jac == SNESComputeJacobianDefaultColor) {
701:       SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,0);
702:     }
703:   }
704:   /* Make sure KSP DM has the Jacobian computation routine */
705:   {
706:     DMSNES sdm;

708:     DMGetDMSNES(snes->dm, &sdm);
709:     if (!sdm->ops->computejacobian) {
710:       DMCopyDMSNES(dmsave, snes->dm);
711:     }
712:   }
713:   /* Compute the operators */
714:   SNESComputeJacobian(snes,X,A,B);
715:   /* Put the previous context back */
716:   if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
717:     SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
718:   }

720:   if (Xnamed) {DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);}
721:   snes->dm = dmsave;
722:   return(0);
723: }

725: /*@
726:    SNESSetUpMatrices - ensures that matrices are available for SNES, to be called by SNESSetUp_XXX()

728:    Collective

730:    Input Arguments:
731: .  snes - snes to configure

733:    Level: developer

735: .seealso: SNESSetUp()
736: @*/
737: PetscErrorCode SNESSetUpMatrices(SNES snes)
738: {
740:   DM             dm;
741:   DMSNES         sdm;

744:   SNESGetDM(snes,&dm);
745:   DMGetDMSNES(dm,&sdm);
746:   if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"DMSNES not properly configured");
747:   else if (!snes->jacobian && snes->mf) {
748:     Mat  J;
749:     void *functx;
750:     MatCreateSNESMF(snes,&J);
751:     MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
752:     MatSetFromOptions(J);
753:     SNESGetFunction(snes,NULL,NULL,&functx);
754:     SNESSetJacobian(snes,J,J,0,0);
755:     MatDestroy(&J);
756:   } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
757:     Mat J,B;
758:     MatCreateSNESMF(snes,&J);
759:     MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
760:     MatSetFromOptions(J);
761:     DMCreateMatrix(snes->dm,&B);
762:     /* sdm->computejacobian was already set to reach here */
763:     SNESSetJacobian(snes,J,B,NULL,NULL);
764:     MatDestroy(&J);
765:     MatDestroy(&B);
766:   } else if (!snes->jacobian_pre) {
767:     PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *);
768:     PetscDS          prob;
769:     Mat              J, B;
770:     MatNullSpace     nullspace = NULL;
771:     PetscBool        hasPrec   = PETSC_FALSE;
772:     PetscInt         Nf;

774:     J    = snes->jacobian;
775:     DMGetDS(dm, &prob);
776:     if (prob) {PetscDSHasJacobianPreconditioner(prob, &hasPrec);}
777:     if (J)            {PetscObjectReference((PetscObject) J);}
778:     else if (hasPrec) {DMCreateMatrix(snes->dm, &J);}
779:     DMCreateMatrix(snes->dm, &B);
780:     PetscDSGetNumFields(prob, &Nf);
781:     DMGetNullSpaceConstructor(snes->dm, Nf, &nspconstr);
782:     if (nspconstr) (*nspconstr)(snes->dm, -1, &nullspace);
783:     MatSetNullSpace(B, nullspace);
784:     MatNullSpaceDestroy(&nullspace);
785:     SNESSetJacobian(snes, J ? J : B, B, NULL, NULL);
786:     MatDestroy(&J);
787:     MatDestroy(&B);
788:   }
789:   {
790:     KSP ksp;
791:     SNESGetKSP(snes,&ksp);
792:     KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);
793:     DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
794:   }
795:   return(0);
796: }

798: /*@C
799:    SNESMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

801:    Collective on SNES

803:    Input Parameters:
804: +  snes - SNES object you wish to monitor
805: .  name - the monitor type one is seeking
806: .  help - message indicating what monitoring is done
807: .  manual - manual page for the monitor
808: .  monitor - the monitor function
809: -  monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the SNES or PetscViewer objects

811:    Level: developer

813: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
814:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
815:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
816:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
817:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
818:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
819:           PetscOptionsFList(), PetscOptionsEList()
820: @*/
821: PetscErrorCode  SNESMonitorSetFromOptions(SNES snes,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNES,PetscViewerAndFormat*))
822: {
823:   PetscErrorCode    ierr;
824:   PetscViewer       viewer;
825:   PetscViewerFormat format;
826:   PetscBool         flg;

829:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,name,&viewer,&format,&flg);
830:   if (flg) {
831:     PetscViewerAndFormat *vf;
832:     PetscViewerAndFormatCreate(viewer,format,&vf);
833:     PetscObjectDereference((PetscObject)viewer);
834:     if (monitorsetup) {
835:       (*monitorsetup)(snes,vf);
836:     }
837:     SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
838:   }
839:   return(0);
840: }

842: /*@
843:    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.

845:    Collective on SNES

847:    Input Parameter:
848: .  snes - the SNES context

850:    Options Database Keys:
851: +  -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
852: .  -snes_stol - convergence tolerance in terms of the norm
853:                 of the change in the solution between steps
854: .  -snes_atol <abstol> - absolute tolerance of residual norm
855: .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
856: .  -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
857: .  -snes_force_iteration <force> - force SNESSolve() to take at least one iteration
858: .  -snes_max_it <max_it> - maximum number of iterations
859: .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
860: .  -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
861: .  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
862: .  -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
863: .  -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
864: .  -snes_trtol <trtol> - trust region tolerance
865: .  -snes_no_convergence_test - skip convergence test in nonlinear
866:                                solver; hence iterations will continue until max_it
867:                                or some other criterion is reached. Saves expense
868:                                of convergence test
869: .  -snes_monitor [ascii][:filename][:viewer format] - prints residual norm at each iteration. if no filename given prints to stdout
870: .  -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration
871: .  -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration
872: .  -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration
873: .  -snes_monitor_lg_residualnorm - plots residual norm at each iteration
874: .  -snes_monitor_lg_range - plots residual norm at each iteration
875: .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
876: .  -snes_fd_color - use finite differences with coloring to compute Jacobian
877: .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
878: .  -snes_converged_reason - print the reason for convergence/divergence after each solve
879: -  -npc_snes_type <type> - the SNES type to use as a nonlinear preconditioner

881:     Options Database for Eisenstat-Walker method:
882: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
883: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
884: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
885: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
886: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
887: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
888: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
889: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

891:    Notes:
892:    To see all options, run your program with the -help option or consult the users manual

894:    Notes:
895:       SNES supports three approaches for computing (approximate) Jacobians: user provided via SNESSetJacobian(), matrix free, and computing explictly with
896:       finite differences and coloring using MatFDColoring. It is also possible to use automatic differentiation and the MatFDColoring object.

898:    Level: beginner

900: .seealso: SNESSetOptionsPrefix(), SNESResetFromOptions(), SNES, SNESCreate()
901: @*/
902: PetscErrorCode  SNESSetFromOptions(SNES snes)
903: {
904:   PetscBool      flg,pcset,persist,set;
905:   PetscInt       i,indx,lag,grids;
906:   const char     *deft        = SNESNEWTONLS;
907:   const char     *convtests[] = {"default","skip"};
908:   SNESKSPEW      *kctx        = NULL;
909:   char           type[256], monfilename[PETSC_MAX_PATH_LEN];
911:   PCSide         pcside;
912:   const char     *optionsprefix;

916:   SNESRegisterAll();
917:   PetscObjectOptionsBegin((PetscObject)snes);
918:   if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
919:   PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
920:   if (flg) {
921:     SNESSetType(snes,type);
922:   } else if (!((PetscObject)snes)->type_name) {
923:     SNESSetType(snes,deft);
924:   }
925:   PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);
926:   PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);

928:   PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,NULL);
929:   PetscOptionsReal("-snes_divergence_tolerance","Stop if residual norm increases by this factor","SNESSetDivergenceTolerance",snes->divtol,&snes->divtol,NULL);
930:   PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,NULL);
931:   PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,NULL);
932:   PetscOptionsInt("-snes_max_fail","Maximum nonlinear step failures","SNESSetMaxNonlinearStepFailures",snes->maxFailures,&snes->maxFailures,NULL);
933:   PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,NULL);
934:   PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,NULL);
935:   PetscOptionsBool("-snes_force_iteration","Force SNESSolve() to take at least one iteration","SNESSetForceIteration",snes->forceiteration,&snes->forceiteration,NULL);
936:   PetscOptionsBool("-snes_check_jacobian_domain_error","Check Jacobian domain error after Jacobian evaluation","SNESCheckJacobianDomainError",snes->checkjacdomainerror,&snes->checkjacdomainerror,NULL);

938:   PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
939:   if (flg) {
940:     SNESSetLagPreconditioner(snes,lag);
941:   }
942:   PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple SNES solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
943:   if (flg) {
944:     SNESSetLagPreconditionerPersists(snes,persist);
945:   }
946:   PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
947:   if (flg) {
948:     SNESSetLagJacobian(snes,lag);
949:   }
950:   PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple SNES solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
951:   if (flg) {
952:     SNESSetLagJacobianPersists(snes,persist);
953:   }

955:   PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
956:   if (flg) {
957:     SNESSetGridSequence(snes,grids);
958:   }

960:   PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);
961:   if (flg) {
962:     switch (indx) {
963:     case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
964:     case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL);    break;
965:     }
966:   }

968:   PetscOptionsEList("-snes_norm_schedule","SNES Norm schedule","SNESSetNormSchedule",SNESNormSchedules,5,"function",&indx,&flg);
969:   if (flg) { SNESSetNormSchedule(snes,(SNESNormSchedule)indx); }

971:   PetscOptionsEList("-snes_function_type","SNES Norm schedule","SNESSetFunctionType",SNESFunctionTypes,2,"unpreconditioned",&indx,&flg);
972:   if (flg) { SNESSetFunctionType(snes,(SNESFunctionType)indx); }

974:   kctx = (SNESKSPEW*)snes->kspconvctx;

976:   PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,NULL);

978:   PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);
979:   PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);
980:   PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);
981:   PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);
982:   PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);
983:   PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);
984:   PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);

986:   flg  = PETSC_FALSE;
987:   PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);
988:   if (set && flg) {SNESMonitorCancel(snes);}

990:   SNESMonitorSetFromOptions(snes,"-snes_monitor","Monitor norm of function","SNESMonitorDefault",SNESMonitorDefault,NULL);
991:   SNESMonitorSetFromOptions(snes,"-snes_monitor_short","Monitor norm of function with fewer digits","SNESMonitorDefaultShort",SNESMonitorDefaultShort,NULL);
992:   SNESMonitorSetFromOptions(snes,"-snes_monitor_range","Monitor range of elements of function","SNESMonitorRange",SNESMonitorRange,NULL);

994:   SNESMonitorSetFromOptions(snes,"-snes_monitor_ratio","Monitor ratios of the norm of function for consecutive steps","SNESMonitorRatio",SNESMonitorRatio,SNESMonitorRatioSetUp);
995:   SNESMonitorSetFromOptions(snes,"-snes_monitor_field","Monitor norm of function (split into fields)","SNESMonitorDefaultField",SNESMonitorDefaultField,NULL);
996:   SNESMonitorSetFromOptions(snes,"-snes_monitor_solution","View solution at each iteration","SNESMonitorSolution",SNESMonitorSolution,NULL);
997:   SNESMonitorSetFromOptions(snes,"-snes_monitor_solution_update","View correction at each iteration","SNESMonitorSolutionUpdate",SNESMonitorSolutionUpdate,NULL);
998:   SNESMonitorSetFromOptions(snes,"-snes_monitor_residual","View residual at each iteration","SNESMonitorResidual",SNESMonitorResidual,NULL);
999:   SNESMonitorSetFromOptions(snes,"-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",SNESMonitorJacUpdateSpectrum,NULL);
1000:   SNESMonitorSetFromOptions(snes,"-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet",SNESMonitorFields,NULL);

1002:   PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
1003:   if (flg) {PetscPythonMonitorSet((PetscObject)snes,monfilename);}

1005:   flg  = PETSC_FALSE;
1006:   PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);
1007:   if (flg) {
1008:     PetscDrawLG ctx;

1010:     SNESMonitorLGCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
1011:     SNESMonitorSet(snes,SNESMonitorLGResidualNorm,ctx,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
1012:   }
1013:   flg  = PETSC_FALSE;
1014:   PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
1015:   if (flg) {
1016:     PetscViewer ctx;

1018:     PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
1019:     SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
1020:   }

1022:   flg  = PETSC_FALSE;
1023:   PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
1024:   if (flg) {
1025:     void    *functx;
1026:     DM      dm;
1027:     DMSNES  sdm;
1028:     SNESGetDM(snes,&dm);
1029:     DMGetDMSNES(dm,&sdm);
1030:     sdm->jacobianctx = NULL;
1031:     SNESGetFunction(snes,NULL,NULL,&functx);
1032:     SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
1033:     PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
1034:   }

1036:   flg  = PETSC_FALSE;
1037:   PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
1038:   if (flg) {
1039:     SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
1040:   }

1042:   flg  = PETSC_FALSE;
1043:   PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);
1044:   if (flg) {
1045:     DM             dm;
1046:     DMSNES         sdm;
1047:     SNESGetDM(snes,&dm);
1048:     DMGetDMSNES(dm,&sdm);
1049:     sdm->jacobianctx = NULL;
1050:     SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,0);
1051:     PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
1052:   }

1054:   flg  = PETSC_FALSE;
1055:   PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf_operator,&flg);
1056:   if (flg && snes->mf_operator) {
1057:     snes->mf_operator = PETSC_TRUE;
1058:     snes->mf          = PETSC_TRUE;
1059:   }
1060:   flg  = PETSC_FALSE;
1061:   PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf,&flg);
1062:   if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE;
1063:   PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,0);

1065:   flg  = PETSC_FALSE;
1066:   SNESGetNPCSide(snes,&pcside);
1067:   PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
1068:   if (flg) {SNESSetNPCSide(snes,pcside);}

1070: #if defined(PETSC_HAVE_SAWS)
1071:   /*
1072:     Publish convergence information using SAWs
1073:   */
1074:   flg  = PETSC_FALSE;
1075:   PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);
1076:   if (flg) {
1077:     void *ctx;
1078:     SNESMonitorSAWsCreate(snes,&ctx);
1079:     SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);
1080:   }
1081: #endif
1082: #if defined(PETSC_HAVE_SAWS)
1083:   {
1084:   PetscBool set;
1085:   flg  = PETSC_FALSE;
1086:   PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);
1087:   if (set) {
1088:     PetscObjectSAWsSetBlock((PetscObject)snes,flg);
1089:   }
1090:   }
1091: #endif

1093:   for (i = 0; i < numberofsetfromoptions; i++) {
1094:     (*othersetfromoptions[i])(snes);
1095:   }

1097:   if (snes->ops->setfromoptions) {
1098:     (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
1099:   }

1101:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
1102:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)snes);
1103:   PetscOptionsEnd();

1105:   if (snes->linesearch) {
1106:     SNESGetLineSearch(snes, &snes->linesearch);
1107:     SNESLineSearchSetFromOptions(snes->linesearch);
1108:   }

1110:   if (snes->usesksp) {
1111:     if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
1112:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
1113:     KSPSetFromOptions(snes->ksp);
1114:   }

1116:   /* if user has set the SNES NPC type via options database, create it. */
1117:   SNESGetOptionsPrefix(snes, &optionsprefix);
1118:   PetscOptionsHasName(((PetscObject)snes)->options,optionsprefix, "-npc_snes_type", &pcset);
1119:   if (pcset && (!snes->npc)) {
1120:     SNESGetNPC(snes, &snes->npc);
1121:   }
1122:   if (snes->npc) {
1123:     SNESSetFromOptions(snes->npc);
1124:   }
1125:   snes->setfromoptionscalled++;
1126:   return(0);
1127: }

1129: /*@
1130:    SNESResetFromOptions - Sets various SNES and KSP parameters from user options ONLY if the SNES was previously set from options

1132:    Collective on SNES

1134:    Input Parameter:
1135: .  snes - the SNES context

1137:    Level: beginner

1139: .seealso: SNESSetFromOptions(), SNESSetOptionsPrefix()
1140: @*/
1141: PetscErrorCode SNESResetFromOptions(SNES snes)
1142: {

1146:   if (snes->setfromoptionscalled) {SNESSetFromOptions(snes);}
1147:   return(0);
1148: }

1150: /*@C
1151:    SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
1152:    the nonlinear solvers.

1154:    Logically Collective on SNES

1156:    Input Parameters:
1157: +  snes - the SNES context
1158: .  compute - function to compute the context
1159: -  destroy - function to destroy the context

1161:    Level: intermediate

1163:    Notes:
1164:    This function is currently not available from Fortran.

1166: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
1167: @*/
1168: PetscErrorCode  SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
1169: {
1172:   snes->ops->usercompute = compute;
1173:   snes->ops->userdestroy = destroy;
1174:   return(0);
1175: }

1177: /*@
1178:    SNESSetApplicationContext - Sets the optional user-defined context for
1179:    the nonlinear solvers.

1181:    Logically Collective on SNES

1183:    Input Parameters:
1184: +  snes - the SNES context
1185: -  usrP - optional user context

1187:    Level: intermediate

1189:    Fortran Notes:
1190:     To use this from Fortran you must write a Fortran interface definition for this
1191:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

1193: .seealso: SNESGetApplicationContext()
1194: @*/
1195: PetscErrorCode  SNESSetApplicationContext(SNES snes,void *usrP)
1196: {
1198:   KSP            ksp;

1202:   SNESGetKSP(snes,&ksp);
1203:   KSPSetApplicationContext(ksp,usrP);
1204:   snes->user = usrP;
1205:   return(0);
1206: }

1208: /*@
1209:    SNESGetApplicationContext - Gets the user-defined context for the
1210:    nonlinear solvers.

1212:    Not Collective

1214:    Input Parameter:
1215: .  snes - SNES context

1217:    Output Parameter:
1218: .  usrP - user context

1220:    Fortran Notes:
1221:     To use this from Fortran you must write a Fortran interface definition for this
1222:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

1224:    Level: intermediate

1226: .seealso: SNESSetApplicationContext()
1227: @*/
1228: PetscErrorCode  SNESGetApplicationContext(SNES snes,void *usrP)
1229: {
1232:   *(void**)usrP = snes->user;
1233:   return(0);
1234: }

1236: /*@
1237:    SNESSetUseMatrixFree - indicates that SNES should use matrix free finite difference matrix vector products internally to apply the Jacobian.

1239:    Collective on SNES

1241:    Input Parameters:
1242: +  snes - SNES context
1243: .  mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1244: -  mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used

1246:    Options Database:
1247: + -snes_mf - use matrix free for both the mat and pmat operator
1248: . -snes_mf_operator - use matrix free only for the mat operator
1249: . -snes_fd_color - compute the Jacobian via coloring and finite differences.
1250: - -snes_fd - compute the Jacobian via finite differences (slow)

1252:    Level: intermediate

1254:    Notes:
1255:       SNES supports three approaches for computing (approximate) Jacobians: user provided via SNESSetJacobian(), matrix free, and computing explictly with
1256:       finite differences and coloring using MatFDColoring. It is also possible to use automatic differentiation and the MatFDColoring object.

1258: .seealso:   SNESGetUseMatrixFree(), MatCreateSNESMF(), SNESComputeJacobianDefaultColor()
1259: @*/
1260: PetscErrorCode  SNESSetUseMatrixFree(SNES snes,PetscBool mf_operator,PetscBool mf)
1261: {
1266:   if (mf && !mf_operator) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"If using mf must also use mf_operator");
1267:   snes->mf          = mf;
1268:   snes->mf_operator = mf_operator;
1269:   return(0);
1270: }

1272: /*@
1273:    SNESGetUseMatrixFree - indicates if the SNES uses matrix free finite difference matrix vector products to apply the Jacobian.

1275:    Collective on SNES

1277:    Input Parameter:
1278: .  snes - SNES context

1280:    Output Parameters:
1281: +  mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1282: -  mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used

1284:    Options Database:
1285: + -snes_mf - use matrix free for both the mat and pmat operator
1286: - -snes_mf_operator - use matrix free only for the mat operator

1288:    Level: intermediate

1290: .seealso:   SNESSetUseMatrixFree(), MatCreateSNESMF()
1291: @*/
1292: PetscErrorCode  SNESGetUseMatrixFree(SNES snes,PetscBool *mf_operator,PetscBool *mf)
1293: {
1296:   if (mf)          *mf          = snes->mf;
1297:   if (mf_operator) *mf_operator = snes->mf_operator;
1298:   return(0);
1299: }

1301: /*@
1302:    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1303:    at this time.

1305:    Not Collective

1307:    Input Parameter:
1308: .  snes - SNES context

1310:    Output Parameter:
1311: .  iter - iteration number

1313:    Notes:
1314:    For example, during the computation of iteration 2 this would return 1.

1316:    This is useful for using lagged Jacobians (where one does not recompute the
1317:    Jacobian at each SNES iteration). For example, the code
1318: .vb
1319:       SNESGetIterationNumber(snes,&it);
1320:       if (!(it % 2)) {
1321:         [compute Jacobian here]
1322:       }
1323: .ve
1324:    can be used in your ComputeJacobian() function to cause the Jacobian to be
1325:    recomputed every second SNES iteration.

1327:    After the SNES solve is complete this will return the number of nonlinear iterations used.

1329:    Level: intermediate

1331: .seealso:   SNESGetLinearSolveIterations()
1332: @*/
1333: PetscErrorCode  SNESGetIterationNumber(SNES snes,PetscInt *iter)
1334: {
1338:   *iter = snes->iter;
1339:   return(0);
1340: }

1342: /*@
1343:    SNESSetIterationNumber - Sets the current iteration number.

1345:    Not Collective

1347:    Input Parameter:
1348: +  snes - SNES context
1349: -  iter - iteration number

1351:    Level: developer

1353: .seealso:   SNESGetLinearSolveIterations()
1354: @*/
1355: PetscErrorCode  SNESSetIterationNumber(SNES snes,PetscInt iter)
1356: {

1361:   PetscObjectSAWsTakeAccess((PetscObject)snes);
1362:   snes->iter = iter;
1363:   PetscObjectSAWsGrantAccess((PetscObject)snes);
1364:   return(0);
1365: }

1367: /*@
1368:    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1369:    attempted by the nonlinear solver.

1371:    Not Collective

1373:    Input Parameter:
1374: .  snes - SNES context

1376:    Output Parameter:
1377: .  nfails - number of unsuccessful steps attempted

1379:    Notes:
1380:    This counter is reset to zero for each successive call to SNESSolve().

1382:    Level: intermediate

1384: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1385:           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1386: @*/
1387: PetscErrorCode  SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1388: {
1392:   *nfails = snes->numFailures;
1393:   return(0);
1394: }

1396: /*@
1397:    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1398:    attempted by the nonlinear solver before it gives up.

1400:    Not Collective

1402:    Input Parameters:
1403: +  snes     - SNES context
1404: -  maxFails - maximum of unsuccessful steps

1406:    Level: intermediate

1408: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1409:           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1410: @*/
1411: PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1412: {
1415:   snes->maxFailures = maxFails;
1416:   return(0);
1417: }

1419: /*@
1420:    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1421:    attempted by the nonlinear solver before it gives up.

1423:    Not Collective

1425:    Input Parameter:
1426: .  snes     - SNES context

1428:    Output Parameter:
1429: .  maxFails - maximum of unsuccessful steps

1431:    Level: intermediate

1433: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1434:           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()

1436: @*/
1437: PetscErrorCode  SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1438: {
1442:   *maxFails = snes->maxFailures;
1443:   return(0);
1444: }

1446: /*@
1447:    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1448:      done by SNES.

1450:    Not Collective

1452:    Input Parameter:
1453: .  snes     - SNES context

1455:    Output Parameter:
1456: .  nfuncs - number of evaluations

1458:    Level: intermediate

1460:    Notes:
1461:     Reset every time SNESSolve is called unless SNESSetCountersReset() is used.

1463: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1464: @*/
1465: PetscErrorCode  SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1466: {
1470:   *nfuncs = snes->nfuncs;
1471:   return(0);
1472: }

1474: /*@
1475:    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1476:    linear solvers.

1478:    Not Collective

1480:    Input Parameter:
1481: .  snes - SNES context

1483:    Output Parameter:
1484: .  nfails - number of failed solves

1486:    Level: intermediate

1488:    Options Database Keys:
1489: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated

1491:    Notes:
1492:    This counter is reset to zero for each successive call to SNESSolve().

1494: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1495: @*/
1496: PetscErrorCode  SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1497: {
1501:   *nfails = snes->numLinearSolveFailures;
1502:   return(0);
1503: }

1505: /*@
1506:    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1507:    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE

1509:    Logically Collective on SNES

1511:    Input Parameters:
1512: +  snes     - SNES context
1513: -  maxFails - maximum allowed linear solve failures

1515:    Level: intermediate

1517:    Options Database Keys:
1518: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated

1520:    Notes:
1521:     By default this is 0; that is SNES returns on the first failed linear solve

1523: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1524: @*/
1525: PetscErrorCode  SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1526: {
1530:   snes->maxLinearSolveFailures = maxFails;
1531:   return(0);
1532: }

1534: /*@
1535:    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1536:      are allowed before SNES terminates

1538:    Not Collective

1540:    Input Parameter:
1541: .  snes     - SNES context

1543:    Output Parameter:
1544: .  maxFails - maximum of unsuccessful solves allowed

1546:    Level: intermediate

1548:    Notes:
1549:     By default this is 1; that is SNES returns on the first failed linear solve

1551: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1552: @*/
1553: PetscErrorCode  SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1554: {
1558:   *maxFails = snes->maxLinearSolveFailures;
1559:   return(0);
1560: }

1562: /*@
1563:    SNESGetLinearSolveIterations - Gets the total number of linear iterations
1564:    used by the nonlinear solver.

1566:    Not Collective

1568:    Input Parameter:
1569: .  snes - SNES context

1571:    Output Parameter:
1572: .  lits - number of linear iterations

1574:    Notes:
1575:    This counter is reset to zero for each successive call to SNESSolve() unless SNESSetCountersReset() is used.

1577:    If the linear solver fails inside the SNESSolve() the iterations for that call to the linear solver are not included. If you wish to count them
1578:    then call KSPGetIterationNumber() after the failed solve.

1580:    Level: intermediate

1582: .seealso:  SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1583: @*/
1584: PetscErrorCode  SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1585: {
1589:   *lits = snes->linear_its;
1590:   return(0);
1591: }

1593: /*@
1594:    SNESSetCountersReset - Sets whether or not the counters for linear iterations and function evaluations
1595:    are reset every time SNESSolve() is called.

1597:    Logically Collective on SNES

1599:    Input Parameter:
1600: +  snes - SNES context
1601: -  reset - whether to reset the counters or not

1603:    Notes:
1604:    This defaults to PETSC_TRUE

1606:    Level: developer

1608: .seealso:  SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC()
1609: @*/
1610: PetscErrorCode  SNESSetCountersReset(SNES snes,PetscBool reset)
1611: {
1615:   snes->counters_reset = reset;
1616:   return(0);
1617: }


1620: /*@
1621:    SNESSetKSP - Sets a KSP context for the SNES object to use

1623:    Not Collective, but the SNES and KSP objects must live on the same MPI_Comm

1625:    Input Parameters:
1626: +  snes - the SNES context
1627: -  ksp - the KSP context

1629:    Notes:
1630:    The SNES object already has its KSP object, you can obtain with SNESGetKSP()
1631:    so this routine is rarely needed.

1633:    The KSP object that is already in the SNES object has its reference count
1634:    decreased by one.

1636:    Level: developer

1638: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1639: @*/
1640: PetscErrorCode  SNESSetKSP(SNES snes,KSP ksp)
1641: {

1648:   PetscObjectReference((PetscObject)ksp);
1649:   if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1650:   snes->ksp = ksp;
1651:   return(0);
1652: }

1654: /* -----------------------------------------------------------*/
1655: /*@
1656:    SNESCreate - Creates a nonlinear solver context.

1658:    Collective

1660:    Input Parameters:
1661: .  comm - MPI communicator

1663:    Output Parameter:
1664: .  outsnes - the new SNES context

1666:    Options Database Keys:
1667: +   -snes_mf - Activates default matrix-free Jacobian-vector products,
1668:                and no preconditioning matrix
1669: .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
1670:                products, and a user-provided preconditioning matrix
1671:                as set by SNESSetJacobian()
1672: -   -snes_fd - Uses (slow!) finite differences to compute Jacobian

1674:    Level: beginner

1676:    Developer Notes:
1677:     SNES always creates a KSP object even though many SNES methods do not use it. This is
1678:                     unfortunate and should be fixed at some point. The flag snes->usesksp indicates if the
1679:                     particular method does use KSP and regulates if the information about the KSP is printed
1680:                     in SNESView(). TSSetFromOptions() does call SNESSetFromOptions() which can lead to users being confused
1681:                     by help messages about meaningless SNES options.

1683:                     SNES always creates the snes->kspconvctx even though it is used by only one type. This should
1684:                     be fixed.

1686: .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()

1688: @*/
1689: PetscErrorCode  SNESCreate(MPI_Comm comm,SNES *outsnes)
1690: {
1692:   SNES           snes;
1693:   SNESKSPEW      *kctx;

1697:   *outsnes = NULL;
1698:   SNESInitializePackage();

1700:   PetscHeaderCreate(snes,SNES_CLASSID,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);

1702:   snes->ops->converged    = SNESConvergedDefault;
1703:   snes->usesksp           = PETSC_TRUE;
1704:   snes->tolerancesset     = PETSC_FALSE;
1705:   snes->max_its           = 50;
1706:   snes->max_funcs         = 10000;
1707:   snes->norm              = 0.0;
1708:   snes->xnorm             = 0.0;
1709:   snes->ynorm             = 0.0;
1710:   snes->normschedule      = SNES_NORM_ALWAYS;
1711:   snes->functype          = SNES_FUNCTION_DEFAULT;
1712: #if defined(PETSC_USE_REAL_SINGLE)
1713:   snes->rtol              = 1.e-5;
1714: #else
1715:   snes->rtol              = 1.e-8;
1716: #endif
1717:   snes->ttol              = 0.0;
1718: #if defined(PETSC_USE_REAL_SINGLE)
1719:   snes->abstol            = 1.e-25;
1720: #else
1721:   snes->abstol            = 1.e-50;
1722: #endif
1723: #if defined(PETSC_USE_REAL_SINGLE)
1724:   snes->stol              = 1.e-5;
1725: #else
1726:   snes->stol              = 1.e-8;
1727: #endif
1728: #if defined(PETSC_USE_REAL_SINGLE)
1729:   snes->deltatol          = 1.e-6;
1730: #else
1731:   snes->deltatol          = 1.e-12;
1732: #endif
1733:   snes->divtol            = 1.e4;
1734:   snes->rnorm0            = 0;
1735:   snes->nfuncs            = 0;
1736:   snes->numFailures       = 0;
1737:   snes->maxFailures       = 1;
1738:   snes->linear_its        = 0;
1739:   snes->lagjacobian       = 1;
1740:   snes->jac_iter          = 0;
1741:   snes->lagjac_persist    = PETSC_FALSE;
1742:   snes->lagpreconditioner = 1;
1743:   snes->pre_iter          = 0;
1744:   snes->lagpre_persist    = PETSC_FALSE;
1745:   snes->numbermonitors    = 0;
1746:   snes->data              = 0;
1747:   snes->setupcalled       = PETSC_FALSE;
1748:   snes->ksp_ewconv        = PETSC_FALSE;
1749:   snes->nwork             = 0;
1750:   snes->work              = 0;
1751:   snes->nvwork            = 0;
1752:   snes->vwork             = 0;
1753:   snes->conv_hist_len     = 0;
1754:   snes->conv_hist_max     = 0;
1755:   snes->conv_hist         = NULL;
1756:   snes->conv_hist_its     = NULL;
1757:   snes->conv_hist_reset   = PETSC_TRUE;
1758:   snes->counters_reset    = PETSC_TRUE;
1759:   snes->vec_func_init_set = PETSC_FALSE;
1760:   snes->reason            = SNES_CONVERGED_ITERATING;
1761:   snes->npcside           = PC_RIGHT;
1762:   snes->setfromoptionscalled = 0;

1764:   snes->mf          = PETSC_FALSE;
1765:   snes->mf_operator = PETSC_FALSE;
1766:   snes->mf_version  = 1;

1768:   snes->numLinearSolveFailures = 0;
1769:   snes->maxLinearSolveFailures = 1;

1771:   snes->vizerotolerance = 1.e-8;
1772: #if defined(PETSC_USE_DEBUG)
1773:   snes->checkjacdomainerror = PETSC_TRUE;
1774: #else
1775:   snes->checkjacdomainerror = PETSC_FALSE;
1776: #endif

1778:   /* Set this to true if the implementation of SNESSolve_XXX does compute the residual at the final solution. */
1779:   snes->alwayscomputesfinalresidual = PETSC_FALSE;

1781:   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
1782:   PetscNewLog(snes,&kctx);

1784:   snes->kspconvctx  = (void*)kctx;
1785:   kctx->version     = 2;
1786:   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1787:                              this was too large for some test cases */
1788:   kctx->rtol_last   = 0.0;
1789:   kctx->rtol_max    = .9;
1790:   kctx->gamma       = 1.0;
1791:   kctx->alpha       = .5*(1.0 + PetscSqrtReal(5.0));
1792:   kctx->alpha2      = kctx->alpha;
1793:   kctx->threshold   = .1;
1794:   kctx->lresid_last = 0.0;
1795:   kctx->norm_last   = 0.0;

1797:   *outsnes = snes;
1798:   return(0);
1799: }

1801: /*MC
1802:     SNESFunction - Functional form used to convey the nonlinear function to be solved by SNES

1804:      Synopsis:
1805:      #include "petscsnes.h"
1806:      PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);

1808:      Input Parameters:
1809: +     snes - the SNES context
1810: .     x    - state at which to evaluate residual
1811: -     ctx     - optional user-defined function context, passed in with SNESSetFunction()

1813:      Output Parameter:
1814: .     f  - vector to put residual (function value)

1816:    Level: intermediate

1818: .seealso:   SNESSetFunction(), SNESGetFunction()
1819: M*/

1821: /*@C
1822:    SNESSetFunction - Sets the function evaluation routine and function
1823:    vector for use by the SNES routines in solving systems of nonlinear
1824:    equations.

1826:    Logically Collective on SNES

1828:    Input Parameters:
1829: +  snes - the SNES context
1830: .  r - vector to store function value
1831: .  f - function evaluation routine; see SNESFunction for calling sequence details
1832: -  ctx - [optional] user-defined context for private data for the
1833:          function evaluation routine (may be NULL)

1835:    Notes:
1836:    The Newton-like methods typically solve linear systems of the form
1837: $      f'(x) x = -f(x),
1838:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

1840:    Level: beginner

1842: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1843: @*/
1844: PetscErrorCode  SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1845: {
1847:   DM             dm;

1851:   if (r) {
1854:     PetscObjectReference((PetscObject)r);
1855:     VecDestroy(&snes->vec_func);

1857:     snes->vec_func = r;
1858:   }
1859:   SNESGetDM(snes,&dm);
1860:   DMSNESSetFunction(dm,f,ctx);
1861:   return(0);
1862: }


1865: /*@C
1866:    SNESSetInitialFunction - Sets the function vector to be used as the
1867:    function norm at the initialization of the method.  In some
1868:    instances, the user has precomputed the function before calling
1869:    SNESSolve.  This function allows one to avoid a redundant call
1870:    to SNESComputeFunction in that case.

1872:    Logically Collective on SNES

1874:    Input Parameters:
1875: +  snes - the SNES context
1876: -  f - vector to store function value

1878:    Notes:
1879:    This should not be modified during the solution procedure.

1881:    This is used extensively in the SNESFAS hierarchy and in nonlinear preconditioning.

1883:    Level: developer

1885: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1886: @*/
1887: PetscErrorCode  SNESSetInitialFunction(SNES snes, Vec f)
1888: {
1890:   Vec            vec_func;

1896:   if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1897:     snes->vec_func_init_set = PETSC_FALSE;
1898:     return(0);
1899:   }
1900:   SNESGetFunction(snes,&vec_func,NULL,NULL);
1901:   VecCopy(f, vec_func);

1903:   snes->vec_func_init_set = PETSC_TRUE;
1904:   return(0);
1905: }

1907: /*@
1908:    SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1909:    of the SNES method.

1911:    Logically Collective on SNES

1913:    Input Parameters:
1914: +  snes - the SNES context
1915: -  normschedule - the frequency of norm computation

1917:    Options Database Key:
1918: .  -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly>

1920:    Notes:
1921:    Only certain SNES methods support certain SNESNormSchedules.  Most require evaluation
1922:    of the nonlinear function and the taking of its norm at every iteration to
1923:    even ensure convergence at all.  However, methods such as custom Gauss-Seidel methods
1924:    (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
1925:    may either be monitored for convergence or not.  As these are often used as nonlinear
1926:    preconditioners, monitoring the norm of their error is not a useful enterprise within
1927:    their solution.

1929:    Level: developer

1931: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1932: @*/
1933: PetscErrorCode  SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1934: {
1937:   snes->normschedule = normschedule;
1938:   return(0);
1939: }


1942: /*@
1943:    SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1944:    of the SNES method.

1946:    Logically Collective on SNES

1948:    Input Parameters:
1949: +  snes - the SNES context
1950: -  normschedule - the type of the norm used

1952:    Level: advanced

1954: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1955: @*/
1956: PetscErrorCode  SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1957: {
1960:   *normschedule = snes->normschedule;
1961:   return(0);
1962: }


1965: /*@
1966:   SNESSetFunctionNorm - Sets the last computed residual norm.

1968:   Logically Collective on SNES

1970:   Input Parameters:
1971: + snes - the SNES context

1973: - normschedule - the frequency of norm computation

1975:   Level: developer

1977: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1978: @*/
1979: PetscErrorCode SNESSetFunctionNorm(SNES snes, PetscReal norm)
1980: {
1983:   snes->norm = norm;
1984:   return(0);
1985: }

1987: /*@
1988:   SNESGetFunctionNorm - Gets the last computed norm of the residual

1990:   Not Collective

1992:   Input Parameter:
1993: . snes - the SNES context

1995:   Output Parameter:
1996: . norm - the last computed residual norm

1998:   Level: developer

2000: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2001: @*/
2002: PetscErrorCode SNESGetFunctionNorm(SNES snes, PetscReal *norm)
2003: {
2007:   *norm = snes->norm;
2008:   return(0);
2009: }

2011: /*@
2012:   SNESGetUpdateNorm - Gets the last computed norm of the Newton update

2014:   Not Collective

2016:   Input Parameter:
2017: . snes - the SNES context

2019:   Output Parameter:
2020: . ynorm - the last computed update norm

2022:   Level: developer

2024: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm()
2025: @*/
2026: PetscErrorCode SNESGetUpdateNorm(SNES snes, PetscReal *ynorm)
2027: {
2031:   *ynorm = snes->ynorm;
2032:   return(0);
2033: }

2035: /*@
2036:   SNESGetSolutionNorm - Gets the last computed norm of the solution

2038:   Not Collective

2040:   Input Parameter:
2041: . snes - the SNES context

2043:   Output Parameter:
2044: . xnorm - the last computed solution norm

2046:   Level: developer

2048: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm(), SNESGetUpdateNorm()
2049: @*/
2050: PetscErrorCode SNESGetSolutionNorm(SNES snes, PetscReal *xnorm)
2051: {
2055:   *xnorm = snes->xnorm;
2056:   return(0);
2057: }

2059: /*@C
2060:    SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
2061:    of the SNES method.

2063:    Logically Collective on SNES

2065:    Input Parameters:
2066: +  snes - the SNES context
2067: -  normschedule - the frequency of norm computation

2069:    Notes:
2070:    Only certain SNES methods support certain SNESNormSchedules.  Most require evaluation
2071:    of the nonlinear function and the taking of its norm at every iteration to
2072:    even ensure convergence at all.  However, methods such as custom Gauss-Seidel methods
2073:    (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
2074:    may either be monitored for convergence or not.  As these are often used as nonlinear
2075:    preconditioners, monitoring the norm of their error is not a useful enterprise within
2076:    their solution.

2078:    Level: developer

2080: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2081: @*/
2082: PetscErrorCode  SNESSetFunctionType(SNES snes, SNESFunctionType type)
2083: {
2086:   snes->functype = type;
2087:   return(0);
2088: }


2091: /*@C
2092:    SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
2093:    of the SNES method.

2095:    Logically Collective on SNES

2097:    Input Parameters:
2098: +  snes - the SNES context
2099: -  normschedule - the type of the norm used

2101:    Level: advanced

2103: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2104: @*/
2105: PetscErrorCode  SNESGetFunctionType(SNES snes, SNESFunctionType *type)
2106: {
2109:   *type = snes->functype;
2110:   return(0);
2111: }

2113: /*MC
2114:     SNESNGSFunction - function used to convey a Gauss-Seidel sweep on the nonlinear function

2116:      Synopsis:
2117:      #include <petscsnes.h>
2118: $    SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);

2120: +  X   - solution vector
2121: .  B   - RHS vector
2122: -  ctx - optional user-defined Gauss-Seidel context

2124:    Level: intermediate

2126: .seealso:   SNESSetNGS(), SNESGetNGS()
2127: M*/

2129: /*@C
2130:    SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
2131:    use with composed nonlinear solvers.

2133:    Input Parameters:
2134: +  snes   - the SNES context
2135: .  f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
2136: -  ctx    - [optional] user-defined context for private data for the
2137:             smoother evaluation routine (may be NULL)

2139:    Notes:
2140:    The NGS routines are used by the composed nonlinear solver to generate
2141:     a problem appropriate update to the solution, particularly FAS.

2143:    Level: intermediate

2145: .seealso: SNESGetFunction(), SNESComputeNGS()
2146: @*/
2147: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
2148: {
2150:   DM             dm;

2154:   SNESGetDM(snes,&dm);
2155:   DMSNESSetNGS(dm,f,ctx);
2156:   return(0);
2157: }

2159: PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
2160: {
2162:   DM             dm;
2163:   DMSNES         sdm;

2166:   SNESGetDM(snes,&dm);
2167:   DMGetDMSNES(dm,&sdm);
2168:   if (!sdm->ops->computepfunction) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");
2169:   if (!sdm->ops->computepjacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard Jacobian.");
2170:   /*  A(x)*x - b(x) */
2171:   PetscStackPush("SNES Picard user function");
2172:   (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
2173:   PetscStackPop;
2174:   PetscStackPush("SNES Picard user Jacobian");
2175:   (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
2176:   PetscStackPop;
2177:   VecScale(f,-1.0);
2178:   MatMultAdd(snes->jacobian,x,f,f);
2179:   return(0);
2180: }

2182: PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
2183: {
2185:   /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
2186:   return(0);
2187: }

2189: /*@C
2190:    SNESSetPicard - Use SNES to solve the semilinear-system A(x) x = b(x) via a Picard type iteration (Picard linearization)

2192:    Logically Collective on SNES

2194:    Input Parameters:
2195: +  snes - the SNES context
2196: .  r - vector to store function value
2197: .  b - function evaluation routine
2198: .  Amat - matrix with which A(x) x - b(x) is to be computed
2199: .  Pmat - matrix from which preconditioner is computed (usually the same as Amat)
2200: .  J  - function to compute matrix value, see SNESJacobianFunction for details on its calling sequence
2201: -  ctx - [optional] user-defined context for private data for the
2202:          function evaluation routine (may be NULL)

2204:    Notes:
2205:     We do not recomemend using this routine. It is far better to provide the nonlinear function F() and some approximation to the Jacobian and use
2206:     an approximate Newton solver. This interface is provided to allow porting/testing a previous Picard based code in PETSc before converting it to approximate Newton.

2208:     One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both

2210: $     Solves the equation A(x) x = b(x) via the defect correction algorithm A(x^{n}) (x^{n+1} - x^{n}) = b(x^{n}) - A(x^{n})x^{n}
2211: $     Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration.

2213:      Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner.

2215:    We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then
2216:    the direct Picard iteration A(x^n) x^{n+1} = b(x^n)

2218:    There is some controversity over the definition of a Picard iteration for nonlinear systems but almost everyone agrees that it involves a linear solve and some
2219:    believe it is the iteration  A(x^{n}) x^{n+1} = b(x^{n}) hence we use the name Picard. If anyone has an authoritative  reference that defines the Picard iteration
2220:    different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).

2222:    Level: intermediate

2224: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
2225: @*/
2226: PetscErrorCode  SNESSetPicard(SNES snes,Vec r,PetscErrorCode (*b)(SNES,Vec,Vec,void*),Mat Amat, Mat Pmat, PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2227: {
2229:   DM             dm;

2233:   SNESGetDM(snes, &dm);
2234:   DMSNESSetPicard(dm,b,J,ctx);
2235:   SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
2236:   SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
2237:   return(0);
2238: }

2240: /*@C
2241:    SNESGetPicard - Returns the context for the Picard iteration

2243:    Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.

2245:    Input Parameter:
2246: .  snes - the SNES context

2248:    Output Parameter:
2249: +  r - the function (or NULL)
2250: .  f - the function (or NULL); see SNESFunction for calling sequence details
2251: .  Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
2252: .  Pmat  - the matrix from which the preconditioner will be constructed (or NULL)
2253: .  J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
2254: -  ctx - the function context (or NULL)

2256:    Level: advanced

2258: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
2259: @*/
2260: PetscErrorCode  SNESGetPicard(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),Mat *Amat, Mat *Pmat, PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2261: {
2263:   DM             dm;

2267:   SNESGetFunction(snes,r,NULL,NULL);
2268:   SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
2269:   SNESGetDM(snes,&dm);
2270:   DMSNESGetPicard(dm,f,J,ctx);
2271:   return(0);
2272: }

2274: /*@C
2275:    SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem

2277:    Logically Collective on SNES

2279:    Input Parameters:
2280: +  snes - the SNES context
2281: .  func - function evaluation routine
2282: -  ctx - [optional] user-defined context for private data for the
2283:          function evaluation routine (may be NULL)

2285:    Calling sequence of func:
2286: $    func (SNES snes,Vec x,void *ctx);

2288: .  f - function vector
2289: -  ctx - optional user-defined function context

2291:    Level: intermediate

2293: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
2294: @*/
2295: PetscErrorCode  SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
2296: {
2299:   if (func) snes->ops->computeinitialguess = func;
2300:   if (ctx)  snes->initialguessP            = ctx;
2301:   return(0);
2302: }

2304: /* --------------------------------------------------------------- */
2305: /*@C
2306:    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
2307:    it assumes a zero right hand side.

2309:    Logically Collective on SNES

2311:    Input Parameter:
2312: .  snes - the SNES context

2314:    Output Parameter:
2315: .  rhs - the right hand side vector or NULL if the right hand side vector is null

2317:    Level: intermediate

2319: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
2320: @*/
2321: PetscErrorCode  SNESGetRhs(SNES snes,Vec *rhs)
2322: {
2326:   *rhs = snes->vec_rhs;
2327:   return(0);
2328: }

2330: /*@
2331:    SNESComputeFunction - Calls the function that has been set with SNESSetFunction().

2333:    Collective on SNES

2335:    Input Parameters:
2336: +  snes - the SNES context
2337: -  x - input vector

2339:    Output Parameter:
2340: .  y - function vector, as set by SNESSetFunction()

2342:    Notes:
2343:    SNESComputeFunction() is typically used within nonlinear solvers
2344:    implementations, so most users would not generally call this routine
2345:    themselves.

2347:    Level: developer

2349: .seealso: SNESSetFunction(), SNESGetFunction()
2350: @*/
2351: PetscErrorCode  SNESComputeFunction(SNES snes,Vec x,Vec y)
2352: {
2354:   DM             dm;
2355:   DMSNES         sdm;

2363:   VecValidValues(x,2,PETSC_TRUE);

2365:   SNESGetDM(snes,&dm);
2366:   DMGetDMSNES(dm,&sdm);
2367:   if (sdm->ops->computefunction) {
2368:     if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2369:       PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2370:     }
2371:     VecLockReadPush(x);
2372:     PetscStackPush("SNES user function");
2373:     /* ensure domainerror is false prior to computefunction evaluation (may not have been reset) */
2374:     snes->domainerror = PETSC_FALSE;
2375:     (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2376:     PetscStackPop;
2377:     VecLockReadPop(x);
2378:     if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2379:       PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2380:     }
2381:   } else if (snes->vec_rhs) {
2382:     MatMult(snes->jacobian, x, y);
2383:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2384:   if (snes->vec_rhs) {
2385:     VecAXPY(y,-1.0,snes->vec_rhs);
2386:   }
2387:   snes->nfuncs++;
2388:   /*
2389:      domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2390:      propagate the value to all processes
2391:   */
2392:   if (snes->domainerror) {
2393:     VecSetInf(y);
2394:   }
2395:   return(0);
2396: }

2398: /*@
2399:    SNESComputeNGS - Calls the Gauss-Seidel function that has been set with  SNESSetNGS().

2401:    Collective on SNES

2403:    Input Parameters:
2404: +  snes - the SNES context
2405: .  x - input vector
2406: -  b - rhs vector

2408:    Output Parameter:
2409: .  x - new solution vector

2411:    Notes:
2412:    SNESComputeNGS() is typically used within composed nonlinear solver
2413:    implementations, so most users would not generally call this routine
2414:    themselves.

2416:    Level: developer

2418: .seealso: SNESSetNGS(), SNESComputeFunction()
2419: @*/
2420: PetscErrorCode  SNESComputeNGS(SNES snes,Vec b,Vec x)
2421: {
2423:   DM             dm;
2424:   DMSNES         sdm;

2432:   if (b) {VecValidValues(b,2,PETSC_TRUE);}
2433:   PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2434:   SNESGetDM(snes,&dm);
2435:   DMGetDMSNES(dm,&sdm);
2436:   if (sdm->ops->computegs) {
2437:     if (b) {VecLockReadPush(b);}
2438:     PetscStackPush("SNES user NGS");
2439:     (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2440:     PetscStackPop;
2441:     if (b) {VecLockReadPop(b);}
2442:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2443:   PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2444:   return(0);
2445: }

2447: PetscErrorCode SNESTestJacobian(SNES snes)
2448: {
2449:   Mat               A,B,C,D,jacobian;
2450:   Vec               x = snes->vec_sol,f = snes->vec_func;
2451:   PetscErrorCode    ierr;
2452:   PetscReal         nrm,gnorm;
2453:   PetscReal         threshold = 1.e-5;
2454:   MatType           mattype;
2455:   PetscInt          m,n,M,N;
2456:   void              *functx;
2457:   PetscBool         complete_print = PETSC_FALSE,threshold_print = PETSC_FALSE,test = PETSC_FALSE,flg;
2458:   PetscViewer       viewer,mviewer;
2459:   MPI_Comm          comm;
2460:   PetscInt          tabs;
2461:   static PetscBool  directionsprinted = PETSC_FALSE;
2462:   PetscViewerFormat format;

2465:   PetscObjectOptionsBegin((PetscObject)snes);
2466:   PetscOptionsName("-snes_test_jacobian","Compare hand-coded and finite difference Jacobians","None",&test);
2467:   PetscOptionsReal("-snes_test_jacobian", "Threshold for element difference between hand-coded and finite difference being meaningful", "None", threshold, &threshold,NULL);
2468:   PetscOptionsViewer("-snes_test_jacobian_view","View difference between hand-coded and finite difference Jacobians element entries","None",&mviewer,&format,&complete_print);
2469:   if (!complete_print) {
2470:     PetscOptionsViewer("-snes_test_jacobian_display","Display difference between hand-coded and finite difference Jacobians","None",&mviewer,&format,&complete_print);
2471:   }
2472:   /* for compatibility with PETSc 3.9 and older. */
2473:   PetscOptionsReal("-snes_test_jacobian_display_threshold", "Display difference between hand-coded and finite difference Jacobians which exceed input threshold", "None", threshold, &threshold, &threshold_print);
2474:   PetscOptionsEnd();
2475:   if (!test) return(0);

2477:   PetscObjectGetComm((PetscObject)snes,&comm);
2478:   PetscViewerASCIIGetStdout(comm,&viewer);
2479:   PetscViewerASCIIGetTab(viewer, &tabs);
2480:   PetscViewerASCIISetTab(viewer, ((PetscObject)snes)->tablevel);
2481:   PetscViewerASCIIPrintf(viewer,"  ---------- Testing Jacobian -------------\n");
2482:   if (!complete_print && !directionsprinted) {
2483:     PetscViewerASCIIPrintf(viewer,"  Run with -snes_test_jacobian_view and optionally -snes_test_jacobian <threshold> to show difference\n");
2484:     PetscViewerASCIIPrintf(viewer,"    of hand-coded and finite difference Jacobian entries greater than <threshold>.\n");
2485:   }
2486:   if (!directionsprinted) {
2487:     PetscViewerASCIIPrintf(viewer,"  Testing hand-coded Jacobian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");
2488:     PetscViewerASCIIPrintf(viewer,"    O(1.e-8), the hand-coded Jacobian is probably correct.\n");
2489:     directionsprinted = PETSC_TRUE;
2490:   }
2491:   if (complete_print) {
2492:     PetscViewerPushFormat(mviewer,format);
2493:   }

2495:   PetscObjectTypeCompare((PetscObject)snes->jacobian,MATMFFD,&flg);
2496:   if (!flg) jacobian = snes->jacobian;
2497:   else jacobian = snes->jacobian_pre;

2499:   if (!x) {
2500:     MatCreateVecs(jacobian, &x, NULL);
2501:   } else {
2502:     PetscObjectReference((PetscObject) x);
2503:   }
2504:   if (!f) {
2505:     VecDuplicate(x, &f);
2506:   } else {
2507:     PetscObjectReference((PetscObject) f);
2508:   }
2509:   /* evaluate the function at this point because SNESComputeJacobianDefault() assumes that the function has been evaluated and put into snes->vec_func */
2510:   SNESComputeFunction(snes,x,f);
2511:   VecDestroy(&f);

2513:   while (jacobian) {
2514:     PetscObjectBaseTypeCompareAny((PetscObject)jacobian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPISBAIJ,"");
2515:     if (flg) {
2516:       A    = jacobian;
2517:       PetscObjectReference((PetscObject)A);
2518:     } else {
2519:       MatComputeOperator(jacobian,MATAIJ,&A);
2520:     }

2522:     MatGetType(A,&mattype);
2523:     MatGetSize(A,&M,&N);
2524:     MatGetLocalSize(A,&m,&n);

2526:     MatCreate(PetscObjectComm((PetscObject)A),&B);
2527:     MatSetType(B,mattype);
2528:     MatSetSizes(B,m,n,M,N);
2529:     MatSetBlockSizesFromMats(B,A,A);
2530:     MatSetUp(B);
2531:     MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);

2533:     SNESGetFunction(snes,NULL,NULL,&functx);
2534:     SNESComputeJacobianDefault(snes,x,B,B,functx);

2536:     MatDuplicate(B,MAT_COPY_VALUES,&D);
2537:     MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2538:     MatNorm(D,NORM_FROBENIUS,&nrm);
2539:     MatNorm(A,NORM_FROBENIUS,&gnorm);
2540:     MatDestroy(&D);
2541:     if (!gnorm) gnorm = 1; /* just in case */
2542:     PetscViewerASCIIPrintf(viewer,"  ||J - Jfd||_F/||J||_F = %g, ||J - Jfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);

2544:     if (complete_print) {
2545:       PetscViewerASCIIPrintf(viewer,"  Hand-coded Jacobian ----------\n");
2546:       MatView(jacobian,mviewer);
2547:       PetscViewerASCIIPrintf(viewer,"  Finite difference Jacobian ----------\n");
2548:       MatView(B,mviewer);
2549:     }

2551:     if (threshold_print || complete_print) {
2552:       PetscInt          Istart, Iend, *ccols, bncols, cncols, j, row;
2553:       PetscScalar       *cvals;
2554:       const PetscInt    *bcols;
2555:       const PetscScalar *bvals;

2557:       MatCreate(PetscObjectComm((PetscObject)A),&C);
2558:       MatSetType(C,mattype);
2559:       MatSetSizes(C,m,n,M,N);
2560:       MatSetBlockSizesFromMats(C,A,A);
2561:       MatSetUp(C);
2562:       MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);

2564:       MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2565:       MatGetOwnershipRange(B,&Istart,&Iend);

2567:       for (row = Istart; row < Iend; row++) {
2568:         MatGetRow(B,row,&bncols,&bcols,&bvals);
2569:         PetscMalloc2(bncols,&ccols,bncols,&cvals);
2570:         for (j = 0, cncols = 0; j < bncols; j++) {
2571:           if (PetscAbsScalar(bvals[j]) > threshold) {
2572:             ccols[cncols] = bcols[j];
2573:             cvals[cncols] = bvals[j];
2574:             cncols += 1;
2575:           }
2576:         }
2577:         if (cncols) {
2578:           MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);
2579:         }
2580:         MatRestoreRow(B,row,&bncols,&bcols,&bvals);
2581:         PetscFree2(ccols,cvals);
2582:       }
2583:       MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2584:       MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2585:       PetscViewerASCIIPrintf(viewer,"  Hand-coded minus finite-difference Jacobian with tolerance %g ----------\n",(double)threshold);
2586:       MatView(C,complete_print ? mviewer : viewer);
2587:       MatDestroy(&C);
2588:     }
2589:     MatDestroy(&A);
2590:     MatDestroy(&B);

2592:     if (jacobian != snes->jacobian_pre) {
2593:       jacobian = snes->jacobian_pre;
2594:       PetscViewerASCIIPrintf(viewer,"  ---------- Testing Jacobian for preconditioner -------------\n");
2595:     }
2596:     else jacobian = NULL;
2597:   }
2598:   VecDestroy(&x);
2599:   if (complete_print) {
2600:     PetscViewerPopFormat(mviewer);
2601:   }
2602:   if (mviewer) { PetscViewerDestroy(&mviewer); }
2603:   PetscViewerASCIISetTab(viewer,tabs);
2604:   return(0);
2605: }

2607: /*@
2608:    SNESComputeJacobian - Computes the Jacobian matrix that has been set with SNESSetJacobian().

2610:    Collective on SNES

2612:    Input Parameters:
2613: +  snes - the SNES context
2614: -  x - input vector

2616:    Output Parameters:
2617: +  A - Jacobian matrix
2618: -  B - optional preconditioning matrix

2620:   Options Database Keys:
2621: +    -snes_lag_preconditioner <lag>
2622: .    -snes_lag_jacobian <lag>
2623: .    -snes_test_jacobian - compare the user provided Jacobian with one compute via finite differences to check for errors
2624: .    -snes_test_jacobian_display - display the user provided Jacobian, the finite difference Jacobian and the difference between them to help users detect the location of errors in the user provided Jacobian
2625: .    -snes_test_jacobian_display_threshold <numerical value>  - display entries in the difference between the user provided Jacobian and finite difference Jacobian that are greater than a certain value to help users detect errors
2626: .    -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2627: .    -snes_compare_explicit_draw  - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2628: .    -snes_compare_explicit_contour  - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2629: .    -snes_compare_operator  - Make the comparison options above use the operator instead of the preconditioning matrix
2630: .    -snes_compare_coloring - Compute the finite difference Jacobian using coloring and display norms of difference
2631: .    -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
2632: .    -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2633: .    -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2634: .    -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2635: .    -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
2636: -    -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences


2639:    Notes:
2640:    Most users should not need to explicitly call this routine, as it
2641:    is used internally within the nonlinear solvers.

2643:    Developer Notes:
2644:     This has duplicative ways of checking the accuracy of the user provided Jacobian (see the options above). This is for historical reasons, the routine SNESTestJacobian() use to used
2645:       for with the SNESType of test that has been removed.

2647:    Level: developer

2649: .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2650: @*/
2651: PetscErrorCode  SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2652: {
2654:   PetscBool      flag;
2655:   DM             dm;
2656:   DMSNES         sdm;
2657:   KSP            ksp;

2663:   VecValidValues(X,2,PETSC_TRUE);
2664:   SNESGetDM(snes,&dm);
2665:   DMGetDMSNES(dm,&sdm);

2667:   if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_USER,"Must call SNESSetJacobian(), DMSNESSetJacobian(), DMDASNESSetJacobianLocal(), etc");

2669:   /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */

2671:   if (snes->lagjacobian == -2) {
2672:     snes->lagjacobian = -1;

2674:     PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2675:   } else if (snes->lagjacobian == -1) {
2676:     PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2677:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2678:     if (flag) {
2679:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2680:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2681:     }
2682:     return(0);
2683:   } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2684:     PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2685:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2686:     if (flag) {
2687:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2688:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2689:     }
2690:     return(0);
2691:   }
2692:   if (snes->npc && snes->npcside== PC_LEFT) {
2693:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2694:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2695:     return(0);
2696:   }

2698:   PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2699:   VecLockReadPush(X);
2700:   PetscStackPush("SNES user Jacobian function");
2701:   (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2702:   PetscStackPop;
2703:   VecLockReadPop(X);
2704:   PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);

2706:   /* attach latest linearization point to the preconditioning matrix */
2707:   PetscObjectCompose((PetscObject)B,"__SNES_latest_X",(PetscObject)X);

2709:   /* the next line ensures that snes->ksp exists */
2710:   SNESGetKSP(snes,&ksp);
2711:   if (snes->lagpreconditioner == -2) {
2712:     PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2713:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2714:     snes->lagpreconditioner = -1;
2715:   } else if (snes->lagpreconditioner == -1) {
2716:     PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2717:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2718:   } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2719:     PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2720:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2721:   } else {
2722:     PetscInfo(snes,"Rebuilding preconditioner\n");
2723:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2724:   }

2726:   SNESTestJacobian(snes);
2727:   /* make sure user returned a correct Jacobian and preconditioner */
2730:   {
2731:     PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
2732:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit",NULL,NULL,&flag);
2733:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",NULL,NULL,&flag_draw);
2734:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",NULL,NULL,&flag_contour);
2735:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_operator",NULL,NULL,&flag_operator);
2736:     if (flag || flag_draw || flag_contour) {
2737:       Mat          Bexp_mine = NULL,Bexp,FDexp;
2738:       PetscViewer  vdraw,vstdout;
2739:       PetscBool    flg;
2740:       if (flag_operator) {
2741:         MatComputeOperator(A,MATAIJ,&Bexp_mine);
2742:         Bexp = Bexp_mine;
2743:       } else {
2744:         /* See if the preconditioning matrix can be viewed and added directly */
2745:         PetscObjectBaseTypeCompareAny((PetscObject)B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2746:         if (flg) Bexp = B;
2747:         else {
2748:           /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
2749:           MatComputeOperator(B,MATAIJ,&Bexp_mine);
2750:           Bexp = Bexp_mine;
2751:         }
2752:       }
2753:       MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);
2754:       SNESComputeJacobianDefault(snes,X,FDexp,FDexp,NULL);
2755:       PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2756:       if (flag_draw || flag_contour) {
2757:         PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2758:         if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2759:       } else vdraw = NULL;
2760:       PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator ? "Jacobian" : "preconditioning Jacobian");
2761:       if (flag) {MatView(Bexp,vstdout);}
2762:       if (vdraw) {MatView(Bexp,vdraw);}
2763:       PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");
2764:       if (flag) {MatView(FDexp,vstdout);}
2765:       if (vdraw) {MatView(FDexp,vdraw);}
2766:       MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);
2767:       PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");
2768:       if (flag) {MatView(FDexp,vstdout);}
2769:       if (vdraw) {              /* Always use contour for the difference */
2770:         PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2771:         MatView(FDexp,vdraw);
2772:         PetscViewerPopFormat(vdraw);
2773:       }
2774:       if (flag_contour) {PetscViewerPopFormat(vdraw);}
2775:       PetscViewerDestroy(&vdraw);
2776:       MatDestroy(&Bexp_mine);
2777:       MatDestroy(&FDexp);
2778:     }
2779:   }
2780:   {
2781:     PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
2782:     PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
2783:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring",NULL,NULL,&flag);
2784:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_display",NULL,NULL,&flag_display);
2785:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",NULL,NULL,&flag_draw);
2786:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",NULL,NULL,&flag_contour);
2787:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",NULL,NULL,&flag_threshold);
2788:     if (flag_threshold) {
2789:       PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,NULL);
2790:       PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,NULL);
2791:     }
2792:     if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
2793:       Mat            Bfd;
2794:       PetscViewer    vdraw,vstdout;
2795:       MatColoring    coloring;
2796:       ISColoring     iscoloring;
2797:       MatFDColoring  matfdcoloring;
2798:       PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2799:       void           *funcctx;
2800:       PetscReal      norm1,norm2,normmax;

2802:       MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2803:       MatColoringCreate(Bfd,&coloring);
2804:       MatColoringSetType(coloring,MATCOLORINGSL);
2805:       MatColoringSetFromOptions(coloring);
2806:       MatColoringApply(coloring,&iscoloring);
2807:       MatColoringDestroy(&coloring);
2808:       MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2809:       MatFDColoringSetFromOptions(matfdcoloring);
2810:       MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2811:       ISColoringDestroy(&iscoloring);

2813:       /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2814:       SNESGetFunction(snes,NULL,&func,&funcctx);
2815:       MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2816:       PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2817:       PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2818:       MatFDColoringSetFromOptions(matfdcoloring);
2819:       MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2820:       MatFDColoringDestroy(&matfdcoloring);

2822:       PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2823:       if (flag_draw || flag_contour) {
2824:         PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2825:         if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2826:       } else vdraw = NULL;
2827:       PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2828:       if (flag_display) {MatView(B,vstdout);}
2829:       if (vdraw) {MatView(B,vdraw);}
2830:       PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2831:       if (flag_display) {MatView(Bfd,vstdout);}
2832:       if (vdraw) {MatView(Bfd,vdraw);}
2833:       MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2834:       MatNorm(Bfd,NORM_1,&norm1);
2835:       MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2836:       MatNorm(Bfd,NORM_MAX,&normmax);
2837:       PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2838:       if (flag_display) {MatView(Bfd,vstdout);}
2839:       if (vdraw) {              /* Always use contour for the difference */
2840:         PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2841:         MatView(Bfd,vdraw);
2842:         PetscViewerPopFormat(vdraw);
2843:       }
2844:       if (flag_contour) {PetscViewerPopFormat(vdraw);}

2846:       if (flag_threshold) {
2847:         PetscInt bs,rstart,rend,i;
2848:         MatGetBlockSize(B,&bs);
2849:         MatGetOwnershipRange(B,&rstart,&rend);
2850:         for (i=rstart; i<rend; i++) {
2851:           const PetscScalar *ba,*ca;
2852:           const PetscInt    *bj,*cj;
2853:           PetscInt          bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
2854:           PetscReal         maxentry = 0,maxdiff = 0,maxrdiff = 0;
2855:           MatGetRow(B,i,&bn,&bj,&ba);
2856:           MatGetRow(Bfd,i,&cn,&cj,&ca);
2857:           if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
2858:           for (j=0; j<bn; j++) {
2859:             PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2860:             if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
2861:               maxentrycol = bj[j];
2862:               maxentry    = PetscRealPart(ba[j]);
2863:             }
2864:             if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
2865:               maxdiffcol = bj[j];
2866:               maxdiff    = PetscRealPart(ca[j]);
2867:             }
2868:             if (rdiff > maxrdiff) {
2869:               maxrdiffcol = bj[j];
2870:               maxrdiff    = rdiff;
2871:             }
2872:           }
2873:           if (maxrdiff > 1) {
2874:             PetscViewerASCIIPrintf(vstdout,"row %D (maxentry=%g at %D, maxdiff=%g at %D, maxrdiff=%g at %D):",i,(double)maxentry,maxentrycol,(double)maxdiff,maxdiffcol,(double)maxrdiff,maxrdiffcol);
2875:             for (j=0; j<bn; j++) {
2876:               PetscReal rdiff;
2877:               rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2878:               if (rdiff > 1) {
2879:                 PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
2880:               }
2881:             }
2882:             PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
2883:           }
2884:           MatRestoreRow(B,i,&bn,&bj,&ba);
2885:           MatRestoreRow(Bfd,i,&cn,&cj,&ca);
2886:         }
2887:       }
2888:       PetscViewerDestroy(&vdraw);
2889:       MatDestroy(&Bfd);
2890:     }
2891:   }
2892:   return(0);
2893: }

2895: /*MC
2896:     SNESJacobianFunction - Function used to convey the nonlinear Jacobian of the function to be solved by SNES

2898:      Synopsis:
2899:      #include "petscsnes.h"
2900:      PetscErrorCode SNESJacobianFunction(SNES snes,Vec x,Mat Amat,Mat Pmat,void *ctx);

2902: +  x - input vector
2903: .  Amat - the matrix that defines the (approximate) Jacobian
2904: .  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2905: -  ctx - [optional] user-defined Jacobian context

2907:    Level: intermediate

2909: .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2910: M*/

2912: /*@C
2913:    SNESSetJacobian - Sets the function to compute Jacobian as well as the
2914:    location to store the matrix.

2916:    Logically Collective on SNES

2918:    Input Parameters:
2919: +  snes - the SNES context
2920: .  Amat - the matrix that defines the (approximate) Jacobian
2921: .  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2922: .  J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
2923: -  ctx - [optional] user-defined context for private data for the
2924:          Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)

2926:    Notes:
2927:    If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2928:    each matrix.

2930:    If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
2931:    space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.

2933:    If using SNESComputeJacobianDefaultColor() to assemble a Jacobian, the ctx argument
2934:    must be a MatFDColoring.

2936:    Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian.  One common
2937:    example is to use the "Picard linearization" which only differentiates through the highest order parts of each term.

2939:    Level: beginner

2941: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J,
2942:           SNESSetPicard(), SNESJacobianFunction
2943: @*/
2944: PetscErrorCode  SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2945: {
2947:   DM             dm;

2955:   SNESGetDM(snes,&dm);
2956:   DMSNESSetJacobian(dm,J,ctx);
2957:   if (Amat) {
2958:     PetscObjectReference((PetscObject)Amat);
2959:     MatDestroy(&snes->jacobian);

2961:     snes->jacobian = Amat;
2962:   }
2963:   if (Pmat) {
2964:     PetscObjectReference((PetscObject)Pmat);
2965:     MatDestroy(&snes->jacobian_pre);

2967:     snes->jacobian_pre = Pmat;
2968:   }
2969:   return(0);
2970: }

2972: /*@C
2973:    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2974:    provided context for evaluating the Jacobian.

2976:    Not Collective, but Mat object will be parallel if SNES object is

2978:    Input Parameter:
2979: .  snes - the nonlinear solver context

2981:    Output Parameters:
2982: +  Amat - location to stash (approximate) Jacobian matrix (or NULL)
2983: .  Pmat - location to stash matrix used to compute the preconditioner (or NULL)
2984: .  J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
2985: -  ctx - location to stash Jacobian ctx (or NULL)

2987:    Level: advanced

2989: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
2990: @*/
2991: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2992: {
2994:   DM             dm;
2995:   DMSNES         sdm;

2999:   if (Amat) *Amat = snes->jacobian;
3000:   if (Pmat) *Pmat = snes->jacobian_pre;
3001:   SNESGetDM(snes,&dm);
3002:   DMGetDMSNES(dm,&sdm);
3003:   if (J) *J = sdm->ops->computejacobian;
3004:   if (ctx) *ctx = sdm->jacobianctx;
3005:   return(0);
3006: }

3008: /*@
3009:    SNESSetUp - Sets up the internal data structures for the later use
3010:    of a nonlinear solver.

3012:    Collective on SNES

3014:    Input Parameters:
3015: .  snes - the SNES context

3017:    Notes:
3018:    For basic use of the SNES solvers the user need not explicitly call
3019:    SNESSetUp(), since these actions will automatically occur during
3020:    the call to SNESSolve().  However, if one wishes to control this
3021:    phase separately, SNESSetUp() should be called after SNESCreate()
3022:    and optional routines of the form SNESSetXXX(), but before SNESSolve().

3024:    Level: advanced

3026: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
3027: @*/
3028: PetscErrorCode  SNESSetUp(SNES snes)
3029: {
3031:   DM             dm;
3032:   DMSNES         sdm;
3033:   SNESLineSearch linesearch, pclinesearch;
3034:   void           *lsprectx,*lspostctx;
3035:   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
3036:   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
3037:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
3038:   Vec            f,fpc;
3039:   void           *funcctx;
3040:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
3041:   void           *jacctx,*appctx;
3042:   Mat            j,jpre;

3046:   if (snes->setupcalled) return(0);
3047:   PetscLogEventBegin(SNES_Setup,snes,0,0,0);

3049:   if (!((PetscObject)snes)->type_name) {
3050:     SNESSetType(snes,SNESNEWTONLS);
3051:   }

3053:   SNESGetFunction(snes,&snes->vec_func,NULL,NULL);

3055:   SNESGetDM(snes,&dm);
3056:   DMGetDMSNES(dm,&sdm);
3057:   if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
3058:   if (!sdm->ops->computejacobian) {
3059:     DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
3060:   }
3061:   if (!snes->vec_func) {
3062:     DMCreateGlobalVector(dm,&snes->vec_func);
3063:   }

3065:   if (!snes->ksp) {
3066:     SNESGetKSP(snes, &snes->ksp);
3067:   }

3069:   if (snes->linesearch) {
3070:     SNESGetLineSearch(snes, &snes->linesearch);
3071:     SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);
3072:   }

3074:   if (snes->npc && (snes->npcside== PC_LEFT)) {
3075:     snes->mf          = PETSC_TRUE;
3076:     snes->mf_operator = PETSC_FALSE;
3077:   }

3079:   if (snes->npc) {
3080:     /* copy the DM over */
3081:     SNESGetDM(snes,&dm);
3082:     SNESSetDM(snes->npc,dm);

3084:     SNESGetFunction(snes,&f,&func,&funcctx);
3085:     VecDuplicate(f,&fpc);
3086:     SNESSetFunction(snes->npc,fpc,func,funcctx);
3087:     SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
3088:     SNESSetJacobian(snes->npc,j,jpre,jac,jacctx);
3089:     SNESGetApplicationContext(snes,&appctx);
3090:     SNESSetApplicationContext(snes->npc,appctx);
3091:     VecDestroy(&fpc);

3093:     /* copy the function pointers over */
3094:     PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)snes->npc);

3096:     /* default to 1 iteration */
3097:     SNESSetTolerances(snes->npc,0.0,0.0,0.0,1,snes->npc->max_funcs);
3098:     if (snes->npcside==PC_RIGHT) {
3099:       SNESSetNormSchedule(snes->npc,SNES_NORM_FINAL_ONLY);
3100:     } else {
3101:       SNESSetNormSchedule(snes->npc,SNES_NORM_NONE);
3102:     }
3103:     SNESSetFromOptions(snes->npc);

3105:     /* copy the line search context over */
3106:     if (snes->linesearch && snes->npc->linesearch) {
3107:       SNESGetLineSearch(snes,&linesearch);
3108:       SNESGetLineSearch(snes->npc,&pclinesearch);
3109:       SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
3110:       SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
3111:       SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
3112:       SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
3113:       PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
3114:     }
3115:   }
3116:   if (snes->mf) {
3117:     SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
3118:   }
3119:   if (snes->ops->usercompute && !snes->user) {
3120:     (*snes->ops->usercompute)(snes,(void**)&snes->user);
3121:   }

3123:   snes->jac_iter = 0;
3124:   snes->pre_iter = 0;

3126:   if (snes->ops->setup) {
3127:     (*snes->ops->setup)(snes);
3128:   }

3130:   if (snes->npc && (snes->npcside== PC_LEFT)) {
3131:     if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3132:       if (snes->linesearch){
3133:         SNESGetLineSearch(snes,&linesearch);
3134:         SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
3135:       }
3136:     }
3137:   }
3138:   PetscLogEventEnd(SNES_Setup,snes,0,0,0);
3139:   snes->setupcalled = PETSC_TRUE;
3140:   return(0);
3141: }

3143: /*@
3144:    SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats

3146:    Collective on SNES

3148:    Input Parameter:
3149: .  snes - iterative context obtained from SNESCreate()

3151:    Level: intermediate

3153:    Notes:
3154:     Also calls the application context destroy routine set with SNESSetComputeApplicationContext()

3156: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
3157: @*/
3158: PetscErrorCode  SNESReset(SNES snes)
3159: {

3164:   if (snes->ops->userdestroy && snes->user) {
3165:     (*snes->ops->userdestroy)((void**)&snes->user);
3166:     snes->user = NULL;
3167:   }
3168:   if (snes->npc) {
3169:     SNESReset(snes->npc);
3170:   }

3172:   if (snes->ops->reset) {
3173:     (*snes->ops->reset)(snes);
3174:   }
3175:   if (snes->ksp) {
3176:     KSPReset(snes->ksp);
3177:   }

3179:   if (snes->linesearch) {
3180:     SNESLineSearchReset(snes->linesearch);
3181:   }

3183:   VecDestroy(&snes->vec_rhs);
3184:   VecDestroy(&snes->vec_sol);
3185:   VecDestroy(&snes->vec_sol_update);
3186:   VecDestroy(&snes->vec_func);
3187:   MatDestroy(&snes->jacobian);
3188:   MatDestroy(&snes->jacobian_pre);
3189:   VecDestroyVecs(snes->nwork,&snes->work);
3190:   VecDestroyVecs(snes->nvwork,&snes->vwork);

3192:   snes->alwayscomputesfinalresidual = PETSC_FALSE;

3194:   snes->nwork       = snes->nvwork = 0;
3195:   snes->setupcalled = PETSC_FALSE;
3196:   return(0);
3197: }

3199: /*@
3200:    SNESDestroy - Destroys the nonlinear solver context that was created
3201:    with SNESCreate().

3203:    Collective on SNES

3205:    Input Parameter:
3206: .  snes - the SNES context

3208:    Level: beginner

3210: .seealso: SNESCreate(), SNESSolve()
3211: @*/
3212: PetscErrorCode  SNESDestroy(SNES *snes)
3213: {

3217:   if (!*snes) return(0);
3219:   if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; return(0);}

3221:   SNESReset((*snes));
3222:   SNESDestroy(&(*snes)->npc);

3224:   /* if memory was published with SAWs then destroy it */
3225:   PetscObjectSAWsViewOff((PetscObject)*snes);
3226:   if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}

3228:   if ((*snes)->dm) {DMCoarsenHookRemove((*snes)->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,*snes);}
3229:   DMDestroy(&(*snes)->dm);
3230:   KSPDestroy(&(*snes)->ksp);
3231:   SNESLineSearchDestroy(&(*snes)->linesearch);

3233:   PetscFree((*snes)->kspconvctx);
3234:   if ((*snes)->ops->convergeddestroy) {
3235:     (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
3236:   }
3237:   if ((*snes)->conv_hist_alloc) {
3238:     PetscFree2((*snes)->conv_hist,(*snes)->conv_hist_its);
3239:   }
3240:   SNESMonitorCancel((*snes));
3241:   PetscHeaderDestroy(snes);
3242:   return(0);
3243: }

3245: /* ----------- Routines to set solver parameters ---------- */

3247: /*@
3248:    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.

3250:    Logically Collective on SNES

3252:    Input Parameters:
3253: +  snes - the SNES context
3254: -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3255:          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that

3257:    Options Database Keys:
3258: .    -snes_lag_preconditioner <lag>

3260:    Notes:
3261:    The default is 1
3262:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3263:    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use

3265:    Level: intermediate

3267: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()

3269: @*/
3270: PetscErrorCode  SNESSetLagPreconditioner(SNES snes,PetscInt lag)
3271: {
3274:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3275:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3277:   snes->lagpreconditioner = lag;
3278:   return(0);
3279: }

3281: /*@
3282:    SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does

3284:    Logically Collective on SNES

3286:    Input Parameters:
3287: +  snes - the SNES context
3288: -  steps - the number of refinements to do, defaults to 0

3290:    Options Database Keys:
3291: .    -snes_grid_sequence <steps>

3293:    Level: intermediate

3295:    Notes:
3296:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

3298: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetGridSequence()

3300: @*/
3301: PetscErrorCode  SNESSetGridSequence(SNES snes,PetscInt steps)
3302: {
3306:   snes->gridsequence = steps;
3307:   return(0);
3308: }

3310: /*@
3311:    SNESGetGridSequence - gets the number of steps of grid sequencing that SNES does

3313:    Logically Collective on SNES

3315:    Input Parameter:
3316: .  snes - the SNES context

3318:    Output Parameter:
3319: .  steps - the number of refinements to do, defaults to 0

3321:    Options Database Keys:
3322: .    -snes_grid_sequence <steps>

3324:    Level: intermediate

3326:    Notes:
3327:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

3329: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESSetGridSequence()

3331: @*/
3332: PetscErrorCode  SNESGetGridSequence(SNES snes,PetscInt *steps)
3333: {
3336:   *steps = snes->gridsequence;
3337:   return(0);
3338: }

3340: /*@
3341:    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt

3343:    Not Collective

3345:    Input Parameter:
3346: .  snes - the SNES context

3348:    Output Parameter:
3349: .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3350:          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that

3352:    Options Database Keys:
3353: .    -snes_lag_preconditioner <lag>

3355:    Notes:
3356:    The default is 1
3357:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

3359:    Level: intermediate

3361: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()

3363: @*/
3364: PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
3365: {
3368:   *lag = snes->lagpreconditioner;
3369:   return(0);
3370: }

3372: /*@
3373:    SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
3374:      often the preconditioner is rebuilt.

3376:    Logically Collective on SNES

3378:    Input Parameters:
3379: +  snes - the SNES context
3380: -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3381:          the Jacobian is built etc. -2 means rebuild at next chance but then never again

3383:    Options Database Keys:
3384: .    -snes_lag_jacobian <lag>

3386:    Notes:
3387:    The default is 1
3388:    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3389:    If  -1 is used before the very first nonlinear solve the CODE WILL FAIL! because no Jacobian is used, use -2 to indicate you want it recomputed
3390:    at the next Newton step but never again (unless it is reset to another value)

3392:    Level: intermediate

3394: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()

3396: @*/
3397: PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
3398: {
3401:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3402:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3404:   snes->lagjacobian = lag;
3405:   return(0);
3406: }

3408: /*@
3409:    SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt

3411:    Not Collective

3413:    Input Parameter:
3414: .  snes - the SNES context

3416:    Output Parameter:
3417: .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3418:          the Jacobian is built etc.

3420:    Options Database Keys:
3421: .    -snes_lag_jacobian <lag>

3423:    Notes:
3424:    The default is 1
3425:    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

3427:    Level: intermediate

3429: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()

3431: @*/
3432: PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
3433: {
3436:   *lag = snes->lagjacobian;
3437:   return(0);
3438: }

3440: /*@
3441:    SNESSetLagJacobianPersists - Set whether or not the Jacobian lagging persists through multiple solves

3443:    Logically collective on SNES

3445:    Input Parameter:
3446: +  snes - the SNES context
3447: -   flg - jacobian lagging persists if true

3449:    Options Database Keys:
3450: .    -snes_lag_jacobian_persists <flg>

3452:    Notes:
3453:     This is useful both for nonlinear preconditioning, where it's appropriate to have the Jacobian be stale by
3454:    several solves, and for implicit time-stepping, where Jacobian lagging in the inner nonlinear solve over several
3455:    timesteps may present huge efficiency gains.

3457:    Level: developer

3459: .seealso: SNESSetLagPreconditionerPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()

3461: @*/
3462: PetscErrorCode  SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
3463: {
3467:   snes->lagjac_persist = flg;
3468:   return(0);
3469: }

3471: /*@
3472:    SNESSetLagPreconditionerPersists - Set whether or not the preconditioner lagging persists through multiple solves

3474:    Logically Collective on SNES

3476:    Input Parameter:
3477: +  snes - the SNES context
3478: -   flg - preconditioner lagging persists if true

3480:    Options Database Keys:
3481: .    -snes_lag_jacobian_persists <flg>

3483:    Notes:
3484:     This is useful both for nonlinear preconditioning, where it's appropriate to have the preconditioner be stale
3485:    by several solves, and for implicit time-stepping, where preconditioner lagging in the inner nonlinear solve over
3486:    several timesteps may present huge efficiency gains.

3488:    Level: developer

3490: .seealso: SNESSetLagJacobianPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()

3492: @*/
3493: PetscErrorCode  SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3494: {
3498:   snes->lagpre_persist = flg;
3499:   return(0);
3500: }

3502: /*@
3503:    SNESSetForceIteration - force SNESSolve() to take at least one iteration regardless of the initial residual norm

3505:    Logically Collective on SNES

3507:    Input Parameters:
3508: +  snes - the SNES context
3509: -  force - PETSC_TRUE require at least one iteration

3511:    Options Database Keys:
3512: .    -snes_force_iteration <force> - Sets forcing an iteration

3514:    Notes:
3515:    This is used sometimes with TS to prevent TS from detecting a false steady state solution

3517:    Level: intermediate

3519: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3520: @*/
3521: PetscErrorCode  SNESSetForceIteration(SNES snes,PetscBool force)
3522: {
3525:   snes->forceiteration = force;
3526:   return(0);
3527: }

3529: /*@
3530:    SNESGetForceIteration - Whether or not to force SNESSolve() take at least one iteration regardless of the initial residual norm

3532:    Logically Collective on SNES

3534:    Input Parameters:
3535: .  snes - the SNES context

3537:    Output Parameter:
3538: .  force - PETSC_TRUE requires at least one iteration.

3540:    Level: intermediate

3542: .seealso: SNESSetForceIteration(), SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3543: @*/
3544: PetscErrorCode  SNESGetForceIteration(SNES snes,PetscBool *force)
3545: {
3548:   *force = snes->forceiteration;
3549:   return(0);
3550: }

3552: /*@
3553:    SNESSetTolerances - Sets various parameters used in convergence tests.

3555:    Logically Collective on SNES

3557:    Input Parameters:
3558: +  snes - the SNES context
3559: .  abstol - absolute convergence tolerance
3560: .  rtol - relative convergence tolerance
3561: .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
3562: .  maxit - maximum number of iterations
3563: -  maxf - maximum number of function evaluations (-1 indicates no limit)

3565:    Options Database Keys:
3566: +    -snes_atol <abstol> - Sets abstol
3567: .    -snes_rtol <rtol> - Sets rtol
3568: .    -snes_stol <stol> - Sets stol
3569: .    -snes_max_it <maxit> - Sets maxit
3570: -    -snes_max_funcs <maxf> - Sets maxf

3572:    Notes:
3573:    The default maximum number of iterations is 50.
3574:    The default maximum number of function evaluations is 1000.

3576:    Level: intermediate

3578: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance(), SNESSetForceIteration()
3579: @*/
3580: PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3581: {

3590:   if (abstol != PETSC_DEFAULT) {
3591:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3592:     snes->abstol = abstol;
3593:   }
3594:   if (rtol != PETSC_DEFAULT) {
3595:     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %g must be non-negative and less than 1.0",(double)rtol);
3596:     snes->rtol = rtol;
3597:   }
3598:   if (stol != PETSC_DEFAULT) {
3599:     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3600:     snes->stol = stol;
3601:   }
3602:   if (maxit != PETSC_DEFAULT) {
3603:     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3604:     snes->max_its = maxit;
3605:   }
3606:   if (maxf != PETSC_DEFAULT) {
3607:     if (maxf < -1) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be -1 or nonnegative",maxf);
3608:     snes->max_funcs = maxf;
3609:   }
3610:   snes->tolerancesset = PETSC_TRUE;
3611:   return(0);
3612: }

3614: /*@
3615:    SNESSetDivergenceTolerance - Sets the divergence tolerance used for the SNES divergence test.

3617:    Logically Collective on SNES

3619:    Input Parameters:
3620: +  snes - the SNES context
3621: -  divtol - the divergence tolerance. Use -1 to deactivate the test.

3623:    Options Database Keys:
3624: .    -snes_divergence_tolerance <divtol> - Sets divtol

3626:    Notes:
3627:    The default divergence tolerance is 1e4.

3629:    Level: intermediate

3631: .seealso: SNESSetTolerances(), SNESGetDivergenceTolerance
3632: @*/
3633: PetscErrorCode  SNESSetDivergenceTolerance(SNES snes,PetscReal divtol)
3634: {

3639:   if (divtol != PETSC_DEFAULT) {
3640:     snes->divtol = divtol;
3641:   }
3642:   else {
3643:     snes->divtol = 1.0e4;
3644:   }
3645:   return(0);
3646: }

3648: /*@
3649:    SNESGetTolerances - Gets various parameters used in convergence tests.

3651:    Not Collective

3653:    Input Parameters:
3654: +  snes - the SNES context
3655: .  atol - absolute convergence tolerance
3656: .  rtol - relative convergence tolerance
3657: .  stol -  convergence tolerance in terms of the norm
3658:            of the change in the solution between steps
3659: .  maxit - maximum number of iterations
3660: -  maxf - maximum number of function evaluations

3662:    Notes:
3663:    The user can specify NULL for any parameter that is not needed.

3665:    Level: intermediate

3667: .seealso: SNESSetTolerances()
3668: @*/
3669: PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3670: {
3673:   if (atol)  *atol  = snes->abstol;
3674:   if (rtol)  *rtol  = snes->rtol;
3675:   if (stol)  *stol  = snes->stol;
3676:   if (maxit) *maxit = snes->max_its;
3677:   if (maxf)  *maxf  = snes->max_funcs;
3678:   return(0);
3679: }

3681: /*@
3682:    SNESGetDivergenceTolerance - Gets divergence tolerance used in divergence test.

3684:    Not Collective

3686:    Input Parameters:
3687: +  snes - the SNES context
3688: -  divtol - divergence tolerance

3690:    Level: intermediate

3692: .seealso: SNESSetDivergenceTolerance()
3693: @*/
3694: PetscErrorCode  SNESGetDivergenceTolerance(SNES snes,PetscReal *divtol)
3695: {
3698:   if (divtol) *divtol = snes->divtol;
3699:   return(0);
3700: }

3702: /*@
3703:    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.

3705:    Logically Collective on SNES

3707:    Input Parameters:
3708: +  snes - the SNES context
3709: -  tol - tolerance

3711:    Options Database Key:
3712: .  -snes_trtol <tol> - Sets tol

3714:    Level: intermediate

3716: .seealso: SNESSetTolerances()
3717: @*/
3718: PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3719: {
3723:   snes->deltatol = tol;
3724:   return(0);
3725: }

3727: /*
3728:    Duplicate the lg monitors for SNES from KSP; for some reason with
3729:    dynamic libraries things don't work under Sun4 if we just use
3730:    macros instead of functions
3731: */
3732: PetscErrorCode  SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,void *ctx)
3733: {

3738:   KSPMonitorLGResidualNorm((KSP)snes,it,norm,ctx);
3739:   return(0);
3740: }

3742: PetscErrorCode  SNESMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)
3743: {

3747:   KSPMonitorLGResidualNormCreate(comm,host,label,x,y,m,n,lgctx);
3748:   return(0);
3749: }

3751: PETSC_INTERN PetscErrorCode  SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);

3753: PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3754: {
3755:   PetscDrawLG      lg;
3756:   PetscErrorCode   ierr;
3757:   PetscReal        x,y,per;
3758:   PetscViewer      v = (PetscViewer)monctx;
3759:   static PetscReal prev; /* should be in the context */
3760:   PetscDraw        draw;

3764:   PetscViewerDrawGetDrawLG(v,0,&lg);
3765:   if (!n) {PetscDrawLGReset(lg);}
3766:   PetscDrawLGGetDraw(lg,&draw);
3767:   PetscDrawSetTitle(draw,"Residual norm");
3768:   x    = (PetscReal)n;
3769:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3770:   else y = -15.0;
3771:   PetscDrawLGAddPoint(lg,&x,&y);
3772:   if (n < 20 || !(n % 5) || snes->reason) {
3773:     PetscDrawLGDraw(lg);
3774:     PetscDrawLGSave(lg);
3775:   }

3777:   PetscViewerDrawGetDrawLG(v,1,&lg);
3778:   if (!n) {PetscDrawLGReset(lg);}
3779:   PetscDrawLGGetDraw(lg,&draw);
3780:   PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3781:    SNESMonitorRange_Private(snes,n,&per);
3782:   x    = (PetscReal)n;
3783:   y    = 100.0*per;
3784:   PetscDrawLGAddPoint(lg,&x,&y);
3785:   if (n < 20 || !(n % 5) || snes->reason) {
3786:     PetscDrawLGDraw(lg);
3787:     PetscDrawLGSave(lg);
3788:   }

3790:   PetscViewerDrawGetDrawLG(v,2,&lg);
3791:   if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3792:   PetscDrawLGGetDraw(lg,&draw);
3793:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3794:   x    = (PetscReal)n;
3795:   y    = (prev - rnorm)/prev;
3796:   PetscDrawLGAddPoint(lg,&x,&y);
3797:   if (n < 20 || !(n % 5) || snes->reason) {
3798:     PetscDrawLGDraw(lg);
3799:     PetscDrawLGSave(lg);
3800:   }

3802:   PetscViewerDrawGetDrawLG(v,3,&lg);
3803:   if (!n) {PetscDrawLGReset(lg);}
3804:   PetscDrawLGGetDraw(lg,&draw);
3805:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3806:   x    = (PetscReal)n;
3807:   y    = (prev - rnorm)/(prev*per);
3808:   if (n > 2) { /*skip initial crazy value */
3809:     PetscDrawLGAddPoint(lg,&x,&y);
3810:   }
3811:   if (n < 20 || !(n % 5) || snes->reason) {
3812:     PetscDrawLGDraw(lg);
3813:     PetscDrawLGSave(lg);
3814:   }
3815:   prev = rnorm;
3816:   return(0);
3817: }

3819: /*@
3820:    SNESMonitor - runs the user provided monitor routines, if they exist

3822:    Collective on SNES

3824:    Input Parameters:
3825: +  snes - nonlinear solver context obtained from SNESCreate()
3826: .  iter - iteration number
3827: -  rnorm - relative norm of the residual

3829:    Notes:
3830:    This routine is called by the SNES implementations.
3831:    It does not typically need to be called by the user.

3833:    Level: developer

3835: .seealso: SNESMonitorSet()
3836: @*/
3837: PetscErrorCode  SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3838: {
3840:   PetscInt       i,n = snes->numbermonitors;

3843:   VecLockReadPush(snes->vec_sol);
3844:   for (i=0; i<n; i++) {
3845:     (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3846:   }
3847:   VecLockReadPop(snes->vec_sol);
3848:   return(0);
3849: }

3851: /* ------------ Routines to set performance monitoring options ----------- */

3853: /*MC
3854:     SNESMonitorFunction - functional form passed to SNESMonitorSet() to monitor convergence of nonlinear solver

3856:      Synopsis:
3857:      #include <petscsnes.h>
3858: $    PetscErrorCode SNESMonitorFunction(SNES snes,PetscInt its, PetscReal norm,void *mctx)

3860: +    snes - the SNES context
3861: .    its - iteration number
3862: .    norm - 2-norm function value (may be estimated)
3863: -    mctx - [optional] monitoring context

3865:    Level: advanced

3867: .seealso:   SNESMonitorSet(), SNESMonitorGet()
3868: M*/

3870: /*@C
3871:    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3872:    iteration of the nonlinear solver to display the iteration's
3873:    progress.

3875:    Logically Collective on SNES

3877:    Input Parameters:
3878: +  snes - the SNES context
3879: .  f - the monitor function, see SNESMonitorFunction for the calling sequence
3880: .  mctx - [optional] user-defined context for private data for the
3881:           monitor routine (use NULL if no context is desired)
3882: -  monitordestroy - [optional] routine that frees monitor context
3883:           (may be NULL)

3885:    Options Database Keys:
3886: +    -snes_monitor        - sets SNESMonitorDefault()
3887: .    -snes_monitor_lg_residualnorm    - sets line graph monitor,
3888:                             uses SNESMonitorLGCreate()
3889: -    -snes_monitor_cancel - cancels all monitors that have
3890:                             been hardwired into a code by
3891:                             calls to SNESMonitorSet(), but
3892:                             does not cancel those set via
3893:                             the options database.

3895:    Notes:
3896:    Several different monitoring routines may be set by calling
3897:    SNESMonitorSet() multiple times; all will be called in the
3898:    order in which they were set.

3900:    Fortran Notes:
3901:     Only a single monitor function can be set for each SNES object

3903:    Level: intermediate

3905: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3906: @*/
3907: PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3908: {
3909:   PetscInt       i;
3911:   PetscBool      identical;

3915:   for (i=0; i<snes->numbermonitors;i++) {
3916:     PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))snes->monitor[i],snes->monitorcontext[i],snes->monitordestroy[i],&identical);
3917:     if (identical) return(0);
3918:   }
3919:   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3920:   snes->monitor[snes->numbermonitors]          = f;
3921:   snes->monitordestroy[snes->numbermonitors]   = monitordestroy;
3922:   snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3923:   return(0);
3924: }

3926: /*@
3927:    SNESMonitorCancel - Clears all the monitor functions for a SNES object.

3929:    Logically Collective on SNES

3931:    Input Parameters:
3932: .  snes - the SNES context

3934:    Options Database Key:
3935: .  -snes_monitor_cancel - cancels all monitors that have been hardwired
3936:     into a code by calls to SNESMonitorSet(), but does not cancel those
3937:     set via the options database

3939:    Notes:
3940:    There is no way to clear one specific monitor from a SNES object.

3942:    Level: intermediate

3944: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3945: @*/
3946: PetscErrorCode  SNESMonitorCancel(SNES snes)
3947: {
3949:   PetscInt       i;

3953:   for (i=0; i<snes->numbermonitors; i++) {
3954:     if (snes->monitordestroy[i]) {
3955:       (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3956:     }
3957:   }
3958:   snes->numbermonitors = 0;
3959:   return(0);
3960: }

3962: /*MC
3963:     SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver

3965:      Synopsis:
3966:      #include <petscsnes.h>
3967: $     PetscErrorCode SNESConvergenceTest(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)

3969: +    snes - the SNES context
3970: .    it - current iteration (0 is the first and is before any Newton step)
3971: .    cctx - [optional] convergence context
3972: .    reason - reason for convergence/divergence
3973: .    xnorm - 2-norm of current iterate
3974: .    gnorm - 2-norm of current step
3975: -    f - 2-norm of function

3977:    Level: intermediate

3979: .seealso:   SNESSetConvergenceTest(), SNESGetConvergenceTest()
3980: M*/

3982: /*@C
3983:    SNESSetConvergenceTest - Sets the function that is to be used
3984:    to test for convergence of the nonlinear iterative solution.

3986:    Logically Collective on SNES

3988:    Input Parameters:
3989: +  snes - the SNES context
3990: .  SNESConvergenceTestFunction - routine to test for convergence
3991: .  cctx - [optional] context for private data for the convergence routine  (may be NULL)
3992: -  destroy - [optional] destructor for the context (may be NULL; PETSC_NULL_FUNCTION in Fortran)

3994:    Level: advanced

3996: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3997: @*/
3998: PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3999: {

4004:   if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
4005:   if (snes->ops->convergeddestroy) {
4006:     (*snes->ops->convergeddestroy)(snes->cnvP);
4007:   }
4008:   snes->ops->converged        = SNESConvergenceTestFunction;
4009:   snes->ops->convergeddestroy = destroy;
4010:   snes->cnvP                  = cctx;
4011:   return(0);
4012: }

4014: /*@
4015:    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.

4017:    Not Collective

4019:    Input Parameter:
4020: .  snes - the SNES context

4022:    Output Parameter:
4023: .  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
4024:             manual pages for the individual convergence tests for complete lists

4026:    Options Database:
4027: .   -snes_converged_reason - prints the reason to standard out

4029:    Level: intermediate

4031:    Notes:
4032:     Should only be called after the call the SNESSolve() is complete, if it is called earlier it returns the value SNES__CONVERGED_ITERATING.

4034: .seealso: SNESSetConvergenceTest(), SNESSetConvergedReason(), SNESConvergedReason
4035: @*/
4036: PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
4037: {
4041:   *reason = snes->reason;
4042:   return(0);
4043: }

4045: /*@
4046:    SNESSetConvergedReason - Sets the reason the SNES iteration was stopped.

4048:    Not Collective

4050:    Input Parameters:
4051: +  snes - the SNES context
4052: -  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
4053:             manual pages for the individual convergence tests for complete lists

4055:    Level: intermediate

4057: .seealso: SNESGetConvergedReason(), SNESSetConvergenceTest(), SNESConvergedReason
4058: @*/
4059: PetscErrorCode SNESSetConvergedReason(SNES snes,SNESConvergedReason reason)
4060: {
4063:   snes->reason = reason;
4064:   return(0);
4065: }

4067: /*@
4068:    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.

4070:    Logically Collective on SNES

4072:    Input Parameters:
4073: +  snes - iterative context obtained from SNESCreate()
4074: .  a   - array to hold history, this array will contain the function norms computed at each step
4075: .  its - integer array holds the number of linear iterations for each solve.
4076: .  na  - size of a and its
4077: -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
4078:            else it continues storing new values for new nonlinear solves after the old ones

4080:    Notes:
4081:    If 'a' and 'its' are NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
4082:    default array of length 10000 is allocated.

4084:    This routine is useful, e.g., when running a code for purposes
4085:    of accurate performance monitoring, when no I/O should be done
4086:    during the section of code that is being timed.

4088:    Level: intermediate

4090: .seealso: SNESGetConvergenceHistory()

4092: @*/
4093: PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
4094: {

4101:   if (!a) {
4102:     if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
4103:     PetscCalloc2(na,&a,na,&its);
4104:     snes->conv_hist_alloc = PETSC_TRUE;
4105:   }
4106:   snes->conv_hist       = a;
4107:   snes->conv_hist_its   = its;
4108:   snes->conv_hist_max   = na;
4109:   snes->conv_hist_len   = 0;
4110:   snes->conv_hist_reset = reset;
4111:   return(0);
4112: }

4114: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4115: #include <engine.h>   /* MATLAB include file */
4116: #include <mex.h>      /* MATLAB include file */

4118: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
4119: {
4120:   mxArray   *mat;
4121:   PetscInt  i;
4122:   PetscReal *ar;

4125:   mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
4126:   ar  = (PetscReal*) mxGetData(mat);
4127:   for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
4128:   PetscFunctionReturn(mat);
4129: }
4130: #endif

4132: /*@C
4133:    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.

4135:    Not Collective

4137:    Input Parameter:
4138: .  snes - iterative context obtained from SNESCreate()

4140:    Output Parameters:
4141: +  a   - array to hold history
4142: .  its - integer array holds the number of linear iterations (or
4143:          negative if not converged) for each solve.
4144: -  na  - size of a and its

4146:    Notes:
4147:     The calling sequence for this routine in Fortran is
4148: $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)

4150:    This routine is useful, e.g., when running a code for purposes
4151:    of accurate performance monitoring, when no I/O should be done
4152:    during the section of code that is being timed.

4154:    Level: intermediate

4156: .seealso: SNESSetConvergencHistory()

4158: @*/
4159: PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
4160: {
4163:   if (a)   *a   = snes->conv_hist;
4164:   if (its) *its = snes->conv_hist_its;
4165:   if (na)  *na  = snes->conv_hist_len;
4166:   return(0);
4167: }

4169: /*@C
4170:   SNESSetUpdate - Sets the general-purpose update function called
4171:   at the beginning of every iteration of the nonlinear solve. Specifically
4172:   it is called just before the Jacobian is "evaluated".

4174:   Logically Collective on SNES

4176:   Input Parameters:
4177: + snes - The nonlinear solver context
4178: - func - The function

4180:   Calling sequence of func:
4181: $ func (SNES snes, PetscInt step);

4183: . step - The current step of the iteration

4185:   Level: advanced

4187:   Note: This is NOT what one uses to update the ghost points before a function evaluation, that should be done at the beginning of your FormFunction()
4188:         This is not used by most users.

4190: .seealso SNESSetJacobian(), SNESSolve()
4191: @*/
4192: PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
4193: {
4196:   snes->ops->update = func;
4197:   return(0);
4198: }

4200: /*
4201:    SNESScaleStep_Private - Scales a step so that its length is less than the
4202:    positive parameter delta.

4204:     Input Parameters:
4205: +   snes - the SNES context
4206: .   y - approximate solution of linear system
4207: .   fnorm - 2-norm of current function
4208: -   delta - trust region size

4210:     Output Parameters:
4211: +   gpnorm - predicted function norm at the new point, assuming local
4212:     linearization.  The value is zero if the step lies within the trust
4213:     region, and exceeds zero otherwise.
4214: -   ynorm - 2-norm of the step

4216:     Note:
4217:     For non-trust region methods such as SNESNEWTONLS, the parameter delta
4218:     is set to be the maximum allowable step size.

4220: */
4221: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
4222: {
4223:   PetscReal      nrm;
4224:   PetscScalar    cnorm;


4232:   VecNorm(y,NORM_2,&nrm);
4233:   if (nrm > *delta) {
4234:     nrm     = *delta/nrm;
4235:     *gpnorm = (1.0 - nrm)*(*fnorm);
4236:     cnorm   = nrm;
4237:     VecScale(y,cnorm);
4238:     *ynorm  = *delta;
4239:   } else {
4240:     *gpnorm = 0.0;
4241:     *ynorm  = nrm;
4242:   }
4243:   return(0);
4244: }

4246: /*@
4247:    SNESReasonView - Displays the reason a SNES solve converged or diverged to a viewer

4249:    Collective on SNES

4251:    Parameter:
4252: +  snes - iterative context obtained from SNESCreate()
4253: -  viewer - the viewer to display the reason


4256:    Options Database Keys:
4257: .  -snes_converged_reason - print reason for converged or diverged, also prints number of iterations

4259:    Level: beginner

4261: .seealso: SNESCreate(), SNESSetUp(), SNESDestroy(), SNESSetTolerances(), SNESConvergedDefault()

4263: @*/
4264: PetscErrorCode  SNESReasonView(SNES snes,PetscViewer viewer)
4265: {
4266:   PetscViewerFormat format;
4267:   PetscBool         isAscii;
4268:   PetscErrorCode    ierr;

4271:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
4272:   if (isAscii) {
4273:     PetscViewerGetFormat(viewer, &format);
4274:     PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
4275:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
4276:       DM                dm;
4277:       Vec               u;
4278:       PetscDS           prob;
4279:       PetscInt          Nf, f;
4280:       PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
4281:       void            **exactCtx;
4282:       PetscReal         error;

4284:       SNESGetDM(snes, &dm);
4285:       SNESGetSolution(snes, &u);
4286:       DMGetDS(dm, &prob);
4287:       PetscDSGetNumFields(prob, &Nf);
4288:       PetscMalloc2(Nf, &exactSol, Nf, &exactCtx);
4289:       for (f = 0; f < Nf; ++f) {PetscDSGetExactSolution(prob, f, &exactSol[f], &exactCtx[f]);}
4290:       DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error);
4291:       PetscFree2(exactSol, exactCtx);
4292:       if (error < 1.0e-11) {PetscViewerASCIIPrintf(viewer, "L_2 Error: < 1.0e-11\n");}
4293:       else                 {PetscViewerASCIIPrintf(viewer, "L_2 Error: %g\n", error);}
4294:     }
4295:     if (snes->reason > 0) {
4296:       if (((PetscObject) snes)->prefix) {
4297:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4298:       } else {
4299:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4300:       }
4301:     } else {
4302:       if (((PetscObject) snes)->prefix) {
4303:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4304:       } else {
4305:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4306:       }
4307:     }
4308:     PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
4309:   }
4310:   return(0);
4311: }

4313: /*@C
4314:   SNESReasonViewFromOptions - Processes command line options to determine if/how a SNESReason is to be viewed.

4316:   Collective on SNES

4318:   Input Parameters:
4319: . snes   - the SNES object

4321:   Level: intermediate

4323: @*/
4324: PetscErrorCode SNESReasonViewFromOptions(SNES snes)
4325: {
4326:   PetscErrorCode    ierr;
4327:   PetscViewer       viewer;
4328:   PetscBool         flg;
4329:   static PetscBool  incall = PETSC_FALSE;
4330:   PetscViewerFormat format;

4333:   if (incall) return(0);
4334:   incall = PETSC_TRUE;
4335:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
4336:   if (flg) {
4337:     PetscViewerPushFormat(viewer,format);
4338:     SNESReasonView(snes,viewer);
4339:     PetscViewerPopFormat(viewer);
4340:     PetscViewerDestroy(&viewer);
4341:   }
4342:   incall = PETSC_FALSE;
4343:   return(0);
4344: }

4346: /*@
4347:    SNESSolve - Solves a nonlinear system F(x) = b.
4348:    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().

4350:    Collective on SNES

4352:    Input Parameters:
4353: +  snes - the SNES context
4354: .  b - the constant part of the equation F(x) = b, or NULL to use zero.
4355: -  x - the solution vector.

4357:    Notes:
4358:    The user should initialize the vector,x, with the initial guess
4359:    for the nonlinear solve prior to calling SNESSolve.  In particular,
4360:    to employ an initial guess of zero, the user should explicitly set
4361:    this vector to zero by calling VecSet().

4363:    Level: beginner

4365: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
4366: @*/
4367: PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
4368: {
4369:   PetscErrorCode    ierr;
4370:   PetscBool         flg;
4371:   PetscInt          grid;
4372:   Vec               xcreated = NULL;
4373:   DM                dm;


4382:   /* High level operations using the nonlinear solver */
4383:   {
4384:     PetscViewer       viewer;
4385:     PetscViewerFormat format;
4386:     PetscInt          num;
4387:     PetscBool         flg;
4388:     static PetscBool  incall = PETSC_FALSE;

4390:     if (!incall) {
4391:       /* Estimate the convergence rate of the discretization */
4392:       PetscOptionsGetViewer(PetscObjectComm((PetscObject) snes),((PetscObject)snes)->options, ((PetscObject) snes)->prefix, "-snes_convergence_estimate", &viewer, &format, &flg);
4393:       if (flg) {
4394:         PetscConvEst conv;
4395:         DM           dm;
4396:         PetscReal   *alpha; /* Convergence rate of the solution error for each field in the L_2 norm */
4397:         PetscInt     Nf;

4399:         incall = PETSC_TRUE;
4400:         SNESGetDM(snes, &dm);
4401:         DMGetNumFields(dm, &Nf);
4402:         PetscCalloc1(Nf, &alpha);
4403:         PetscConvEstCreate(PetscObjectComm((PetscObject) snes), &conv);
4404:         PetscConvEstSetSolver(conv, (PetscObject) snes);
4405:         PetscConvEstSetFromOptions(conv);
4406:         PetscConvEstSetUp(conv);
4407:         PetscConvEstGetConvRate(conv, alpha);
4408:         PetscViewerPushFormat(viewer, format);
4409:         PetscConvEstRateView(conv, alpha, viewer);
4410:         PetscViewerPopFormat(viewer);
4411:         PetscViewerDestroy(&viewer);
4412:         PetscConvEstDestroy(&conv);
4413:         PetscFree(alpha);
4414:         incall = PETSC_FALSE;
4415:       }
4416:       /* Adaptively refine the initial grid */
4417:       num  = 1;
4418:       PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_initial", &num, &flg);
4419:       if (flg) {
4420:         DMAdaptor adaptor;

4422:         incall = PETSC_TRUE;
4423:         DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);
4424:         DMAdaptorSetSolver(adaptor, snes);
4425:         DMAdaptorSetSequenceLength(adaptor, num);
4426:         DMAdaptorSetFromOptions(adaptor);
4427:         DMAdaptorSetUp(adaptor);
4428:         DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_INITIAL, &dm, &x);
4429:         DMAdaptorDestroy(&adaptor);
4430:         incall = PETSC_FALSE;
4431:       }
4432:       /* Use grid sequencing to adapt */
4433:       num  = 0;
4434:       PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_sequence", &num, NULL);
4435:       if (num) {
4436:         DMAdaptor adaptor;

4438:         incall = PETSC_TRUE;
4439:         DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);
4440:         DMAdaptorSetSolver(adaptor, snes);
4441:         DMAdaptorSetSequenceLength(adaptor, num);
4442:         DMAdaptorSetFromOptions(adaptor);
4443:         DMAdaptorSetUp(adaptor);
4444:         DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_SEQUENTIAL, &dm, &x);
4445:         DMAdaptorDestroy(&adaptor);
4446:         incall = PETSC_FALSE;
4447:       }
4448:     }
4449:   }
4450:   if (!x) {
4451:     SNESGetDM(snes,&dm);
4452:     DMCreateGlobalVector(dm,&xcreated);
4453:     x    = xcreated;
4454:   }
4455:   SNESViewFromOptions(snes,NULL,"-snes_view_pre");

4457:   for (grid=0; grid<snes->gridsequence; grid++) {PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));}
4458:   for (grid=0; grid<snes->gridsequence+1; grid++) {

4460:     /* set solution vector */
4461:     if (!grid) {PetscObjectReference((PetscObject)x);}
4462:     VecDestroy(&snes->vec_sol);
4463:     snes->vec_sol = x;
4464:     SNESGetDM(snes,&dm);

4466:     /* set affine vector if provided */
4467:     if (b) { PetscObjectReference((PetscObject)b); }
4468:     VecDestroy(&snes->vec_rhs);
4469:     snes->vec_rhs = b;

4471:     if (snes->vec_rhs && (snes->vec_func == snes->vec_rhs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Right hand side vector cannot be function vector");
4472:     if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
4473:     if (snes->vec_rhs  == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
4474:     if (!snes->vec_sol_update /* && snes->vec_sol */) {
4475:       VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
4476:       PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
4477:     }
4478:     DMShellSetGlobalVector(dm,snes->vec_sol);
4479:     SNESSetUp(snes);

4481:     if (!grid) {
4482:       if (snes->ops->computeinitialguess) {
4483:         (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
4484:       }
4485:     }

4487:     if (snes->conv_hist_reset) snes->conv_hist_len = 0;
4488:     if (snes->counters_reset) {snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;}

4490:     PetscLogEventBegin(SNES_Solve,snes,0,0,0);
4491:     (*snes->ops->solve)(snes);
4492:     PetscLogEventEnd(SNES_Solve,snes,0,0,0);
4493:     if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
4494:     snes->domainerror = PETSC_FALSE; /* clear the flag if it has been set */

4496:     if (snes->lagjac_persist) snes->jac_iter += snes->iter;
4497:     if (snes->lagpre_persist) snes->pre_iter += snes->iter;

4499:     PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_local_min",NULL,NULL,&flg);
4500:     if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
4501:     SNESReasonViewFromOptions(snes);

4503:     if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
4504:     if (snes->reason < 0) break;
4505:     if (grid <  snes->gridsequence) {
4506:       DM  fine;
4507:       Vec xnew;
4508:       Mat interp;

4510:       DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
4511:       if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
4512:       DMCreateInterpolation(snes->dm,fine,&interp,NULL);
4513:       DMCreateGlobalVector(fine,&xnew);
4514:       MatInterpolate(interp,x,xnew);
4515:       DMInterpolate(snes->dm,interp,fine);
4516:       MatDestroy(&interp);
4517:       x    = xnew;

4519:       SNESReset(snes);
4520:       SNESSetDM(snes,fine);
4521:       SNESResetFromOptions(snes);
4522:       DMDestroy(&fine);
4523:       PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
4524:     }
4525:   }
4526:   SNESViewFromOptions(snes,NULL,"-snes_view");
4527:   VecViewFromOptions(snes->vec_sol,(PetscObject)snes,"-snes_view_solution");
4528:   DMMonitor(snes->dm);

4530:   VecDestroy(&xcreated);
4531:   PetscObjectSAWsBlock((PetscObject)snes);
4532:   return(0);
4533: }

4535: /* --------- Internal routines for SNES Package --------- */

4537: /*@C
4538:    SNESSetType - Sets the method for the nonlinear solver.

4540:    Collective on SNES

4542:    Input Parameters:
4543: +  snes - the SNES context
4544: -  type - a known method

4546:    Options Database Key:
4547: .  -snes_type <type> - Sets the method; use -help for a list
4548:    of available methods (for instance, newtonls or newtontr)

4550:    Notes:
4551:    See "petsc/include/petscsnes.h" for available methods (for instance)
4552: +    SNESNEWTONLS - Newton's method with line search
4553:      (systems of nonlinear equations)
4554: -    SNESNEWTONTR - Newton's method with trust region
4555:      (systems of nonlinear equations)

4557:   Normally, it is best to use the SNESSetFromOptions() command and then
4558:   set the SNES solver type from the options database rather than by using
4559:   this routine.  Using the options database provides the user with
4560:   maximum flexibility in evaluating the many nonlinear solvers.
4561:   The SNESSetType() routine is provided for those situations where it
4562:   is necessary to set the nonlinear solver independently of the command
4563:   line or options database.  This might be the case, for example, when
4564:   the choice of solver changes during the execution of the program,
4565:   and the user's application is taking responsibility for choosing the
4566:   appropriate method.

4568:     Developer Notes:
4569:     SNESRegister() adds a constructor for a new SNESType to SNESList, SNESSetType() locates
4570:     the constructor in that list and calls it to create the spexific object.

4572:   Level: intermediate

4574: .seealso: SNESType, SNESCreate(), SNESDestroy(), SNESGetType(), SNESSetFromOptions()

4576: @*/
4577: PetscErrorCode  SNESSetType(SNES snes,SNESType type)
4578: {
4579:   PetscErrorCode ierr,(*r)(SNES);
4580:   PetscBool      match;


4586:   PetscObjectTypeCompare((PetscObject)snes,type,&match);
4587:   if (match) return(0);

4589:   PetscFunctionListFind(SNESList,type,&r);
4590:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
4591:   /* Destroy the previous private SNES context */
4592:   if (snes->ops->destroy) {
4593:     (*(snes)->ops->destroy)(snes);
4594:     snes->ops->destroy = NULL;
4595:   }
4596:   /* Reinitialize function pointers in SNESOps structure */
4597:   snes->ops->setup          = 0;
4598:   snes->ops->solve          = 0;
4599:   snes->ops->view           = 0;
4600:   snes->ops->setfromoptions = 0;
4601:   snes->ops->destroy        = 0;

4603:   /* It may happen the user has customized the line search before calling SNESSetType */
4604:   if (((PetscObject)snes)->type_name) {
4605:     SNESLineSearchDestroy(&snes->linesearch);
4606:   }

4608:   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
4609:   snes->setupcalled = PETSC_FALSE;

4611:   PetscObjectChangeTypeName((PetscObject)snes,type);
4612:   (*r)(snes);
4613:   return(0);
4614: }

4616: /*@C
4617:    SNESGetType - Gets the SNES method type and name (as a string).

4619:    Not Collective

4621:    Input Parameter:
4622: .  snes - nonlinear solver context

4624:    Output Parameter:
4625: .  type - SNES method (a character string)

4627:    Level: intermediate

4629: @*/
4630: PetscErrorCode  SNESGetType(SNES snes,SNESType *type)
4631: {
4635:   *type = ((PetscObject)snes)->type_name;
4636:   return(0);
4637: }

4639: /*@
4640:   SNESSetSolution - Sets the solution vector for use by the SNES routines.

4642:   Logically Collective on SNES

4644:   Input Parameters:
4645: + snes - the SNES context obtained from SNESCreate()
4646: - u    - the solution vector

4648:   Level: beginner

4650: @*/
4651: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4652: {
4653:   DM             dm;

4659:   PetscObjectReference((PetscObject) u);
4660:   VecDestroy(&snes->vec_sol);

4662:   snes->vec_sol = u;

4664:   SNESGetDM(snes, &dm);
4665:   DMShellSetGlobalVector(dm, u);
4666:   return(0);
4667: }

4669: /*@
4670:    SNESGetSolution - Returns the vector where the approximate solution is
4671:    stored. This is the fine grid solution when using SNESSetGridSequence().

4673:    Not Collective, but Vec is parallel if SNES is parallel

4675:    Input Parameter:
4676: .  snes - the SNES context

4678:    Output Parameter:
4679: .  x - the solution

4681:    Level: intermediate

4683: .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
4684: @*/
4685: PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
4686: {
4690:   *x = snes->vec_sol;
4691:   return(0);
4692: }

4694: /*@
4695:    SNESGetSolutionUpdate - Returns the vector where the solution update is
4696:    stored.

4698:    Not Collective, but Vec is parallel if SNES is parallel

4700:    Input Parameter:
4701: .  snes - the SNES context

4703:    Output Parameter:
4704: .  x - the solution update

4706:    Level: advanced

4708: .seealso: SNESGetSolution(), SNESGetFunction()
4709: @*/
4710: PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
4711: {
4715:   *x = snes->vec_sol_update;
4716:   return(0);
4717: }

4719: /*@C
4720:    SNESGetFunction - Returns the vector where the function is stored.

4722:    Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.

4724:    Input Parameter:
4725: .  snes - the SNES context

4727:    Output Parameter:
4728: +  r - the vector that is used to store residuals (or NULL if you don't want it)
4729: .  f - the function (or NULL if you don't want it); see SNESFunction for calling sequence details
4730: -  ctx - the function context (or NULL if you don't want it)

4732:    Level: advanced

4734:     Notes: The vector r DOES NOT, in general contain the current value of the SNES nonlinear function

4736: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4737: @*/
4738: PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4739: {
4741:   DM             dm;

4745:   if (r) {
4746:     if (!snes->vec_func) {
4747:       if (snes->vec_rhs) {
4748:         VecDuplicate(snes->vec_rhs,&snes->vec_func);
4749:       } else if (snes->vec_sol) {
4750:         VecDuplicate(snes->vec_sol,&snes->vec_func);
4751:       } else if (snes->dm) {
4752:         DMCreateGlobalVector(snes->dm,&snes->vec_func);
4753:       }
4754:     }
4755:     *r = snes->vec_func;
4756:   }
4757:   SNESGetDM(snes,&dm);
4758:   DMSNESGetFunction(dm,f,ctx);
4759:   return(0);
4760: }

4762: /*@C
4763:    SNESGetNGS - Returns the NGS function and context.

4765:    Input Parameter:
4766: .  snes - the SNES context

4768:    Output Parameter:
4769: +  f - the function (or NULL) see SNESNGSFunction for details
4770: -  ctx    - the function context (or NULL)

4772:    Level: advanced

4774: .seealso: SNESSetNGS(), SNESGetFunction()
4775: @*/

4777: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4778: {
4780:   DM             dm;

4784:   SNESGetDM(snes,&dm);
4785:   DMSNESGetNGS(dm,f,ctx);
4786:   return(0);
4787: }

4789: /*@C
4790:    SNESSetOptionsPrefix - Sets the prefix used for searching for all
4791:    SNES options in the database.

4793:    Logically Collective on SNES

4795:    Input Parameter:
4796: +  snes - the SNES context
4797: -  prefix - the prefix to prepend to all option names

4799:    Notes:
4800:    A hyphen (-) must NOT be given at the beginning of the prefix name.
4801:    The first character of all runtime options is AUTOMATICALLY the hyphen.

4803:    Level: advanced

4805: .seealso: SNESSetFromOptions()
4806: @*/
4807: PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
4808: {

4813:   PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4814:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4815:   if (snes->linesearch) {
4816:     SNESGetLineSearch(snes,&snes->linesearch);
4817:     PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4818:   }
4819:   KSPSetOptionsPrefix(snes->ksp,prefix);
4820:   return(0);
4821: }

4823: /*@C
4824:    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4825:    SNES options in the database.

4827:    Logically Collective on SNES

4829:    Input Parameters:
4830: +  snes - the SNES context
4831: -  prefix - the prefix to prepend to all option names

4833:    Notes:
4834:    A hyphen (-) must NOT be given at the beginning of the prefix name.
4835:    The first character of all runtime options is AUTOMATICALLY the hyphen.

4837:    Level: advanced

4839: .seealso: SNESGetOptionsPrefix()
4840: @*/
4841: PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4842: {

4847:   PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4848:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4849:   if (snes->linesearch) {
4850:     SNESGetLineSearch(snes,&snes->linesearch);
4851:     PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4852:   }
4853:   KSPAppendOptionsPrefix(snes->ksp,prefix);
4854:   return(0);
4855: }

4857: /*@C
4858:    SNESGetOptionsPrefix - Sets the prefix used for searching for all
4859:    SNES options in the database.

4861:    Not Collective

4863:    Input Parameter:
4864: .  snes - the SNES context

4866:    Output Parameter:
4867: .  prefix - pointer to the prefix string used

4869:    Notes:
4870:     On the fortran side, the user should pass in a string 'prefix' of
4871:    sufficient length to hold the prefix.

4873:    Level: advanced

4875: .seealso: SNESAppendOptionsPrefix()
4876: @*/
4877: PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4878: {

4883:   PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4884:   return(0);
4885: }


4888: /*@C
4889:   SNESRegister - Adds a method to the nonlinear solver package.

4891:    Not collective

4893:    Input Parameters:
4894: +  name_solver - name of a new user-defined solver
4895: -  routine_create - routine to create method context

4897:    Notes:
4898:    SNESRegister() may be called multiple times to add several user-defined solvers.

4900:    Sample usage:
4901: .vb
4902:    SNESRegister("my_solver",MySolverCreate);
4903: .ve

4905:    Then, your solver can be chosen with the procedural interface via
4906: $     SNESSetType(snes,"my_solver")
4907:    or at runtime via the option
4908: $     -snes_type my_solver

4910:    Level: advanced

4912:     Note: If your function is not being put into a shared library then use SNESRegister() instead

4914: .seealso: SNESRegisterAll(), SNESRegisterDestroy()

4916:   Level: advanced
4917: @*/
4918: PetscErrorCode  SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4919: {

4923:   SNESInitializePackage();
4924:   PetscFunctionListAdd(&SNESList,sname,function);
4925:   return(0);
4926: }

4928: PetscErrorCode  SNESTestLocalMin(SNES snes)
4929: {
4931:   PetscInt       N,i,j;
4932:   Vec            u,uh,fh;
4933:   PetscScalar    value;
4934:   PetscReal      norm;

4937:   SNESGetSolution(snes,&u);
4938:   VecDuplicate(u,&uh);
4939:   VecDuplicate(u,&fh);

4941:   /* currently only works for sequential */
4942:   PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4943:   VecGetSize(u,&N);
4944:   for (i=0; i<N; i++) {
4945:     VecCopy(u,uh);
4946:     PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4947:     for (j=-10; j<11; j++) {
4948:       value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4949:       VecSetValue(uh,i,value,ADD_VALUES);
4950:       SNESComputeFunction(snes,uh,fh);
4951:       VecNorm(fh,NORM_2,&norm);
4952:       PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);
4953:       value = -value;
4954:       VecSetValue(uh,i,value,ADD_VALUES);
4955:     }
4956:   }
4957:   VecDestroy(&uh);
4958:   VecDestroy(&fh);
4959:   return(0);
4960: }

4962: /*@
4963:    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4964:    computing relative tolerance for linear solvers within an inexact
4965:    Newton method.

4967:    Logically Collective on SNES

4969:    Input Parameters:
4970: +  snes - SNES context
4971: -  flag - PETSC_TRUE or PETSC_FALSE

4973:     Options Database:
4974: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4975: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
4976: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4977: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4978: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
4979: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
4980: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4981: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

4983:    Notes:
4984:    Currently, the default is to use a constant relative tolerance for
4985:    the inner linear solvers.  Alternatively, one can use the
4986:    Eisenstat-Walker method, where the relative convergence tolerance
4987:    is reset at each Newton iteration according progress of the nonlinear
4988:    solver.

4990:    Level: advanced

4992:    Reference:
4993:    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4994:    inexact Newton method", SISC 17 (1), pp.16-32, 1996.

4996: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4997: @*/
4998: PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool flag)
4999: {
5003:   snes->ksp_ewconv = flag;
5004:   return(0);
5005: }

5007: /*@
5008:    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
5009:    for computing relative tolerance for linear solvers within an
5010:    inexact Newton method.

5012:    Not Collective

5014:    Input Parameter:
5015: .  snes - SNES context

5017:    Output Parameter:
5018: .  flag - PETSC_TRUE or PETSC_FALSE

5020:    Notes:
5021:    Currently, the default is to use a constant relative tolerance for
5022:    the inner linear solvers.  Alternatively, one can use the
5023:    Eisenstat-Walker method, where the relative convergence tolerance
5024:    is reset at each Newton iteration according progress of the nonlinear
5025:    solver.

5027:    Level: advanced

5029:    Reference:
5030:    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
5031:    inexact Newton method", SISC 17 (1), pp.16-32, 1996.

5033: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
5034: @*/
5035: PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
5036: {
5040:   *flag = snes->ksp_ewconv;
5041:   return(0);
5042: }

5044: /*@
5045:    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
5046:    convergence criteria for the linear solvers within an inexact
5047:    Newton method.

5049:    Logically Collective on SNES

5051:    Input Parameters:
5052: +    snes - SNES context
5053: .    version - version 1, 2 (default is 2) or 3
5054: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5055: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5056: .    gamma - multiplicative factor for version 2 rtol computation
5057:              (0 <= gamma2 <= 1)
5058: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
5059: .    alpha2 - power for safeguard
5060: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

5062:    Note:
5063:    Version 3 was contributed by Luis Chacon, June 2006.

5065:    Use PETSC_DEFAULT to retain the default for any of the parameters.

5067:    Level: advanced

5069:    Reference:
5070:    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
5071:    inexact Newton method", Utah State University Math. Stat. Dept. Res.
5072:    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.

5074: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
5075: @*/
5076: PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
5077: {
5078:   SNESKSPEW *kctx;

5082:   kctx = (SNESKSPEW*)snes->kspconvctx;
5083:   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");

5092:   if (version != PETSC_DEFAULT)   kctx->version   = version;
5093:   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
5094:   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
5095:   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
5096:   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
5097:   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
5098:   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;

5100:   if (kctx->version < 1 || kctx->version > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
5101:   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %g",(double)kctx->rtol_0);
5102:   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%g) < 1.0\n",(double)kctx->rtol_max);
5103:   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%g) <= 1.0\n",(double)kctx->gamma);
5104:   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%g) <= 2.0\n",(double)kctx->alpha);
5105:   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%g) < 1.0\n",(double)kctx->threshold);
5106:   return(0);
5107: }

5109: /*@
5110:    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
5111:    convergence criteria for the linear solvers within an inexact
5112:    Newton method.

5114:    Not Collective

5116:    Input Parameters:
5117:      snes - SNES context

5119:    Output Parameters:
5120: +    version - version 1, 2 (default is 2) or 3
5121: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5122: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5123: .    gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
5124: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
5125: .    alpha2 - power for safeguard
5126: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

5128:    Level: advanced

5130: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
5131: @*/
5132: PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
5133: {
5134:   SNESKSPEW *kctx;

5138:   kctx = (SNESKSPEW*)snes->kspconvctx;
5139:   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
5140:   if (version)   *version   = kctx->version;
5141:   if (rtol_0)    *rtol_0    = kctx->rtol_0;
5142:   if (rtol_max)  *rtol_max  = kctx->rtol_max;
5143:   if (gamma)     *gamma     = kctx->gamma;
5144:   if (alpha)     *alpha     = kctx->alpha;
5145:   if (alpha2)    *alpha2    = kctx->alpha2;
5146:   if (threshold) *threshold = kctx->threshold;
5147:   return(0);
5148: }

5150:  PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5151: {
5153:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
5154:   PetscReal      rtol  = PETSC_DEFAULT,stol;

5157:   if (!snes->ksp_ewconv) return(0);
5158:   if (!snes->iter) {
5159:     rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
5160:     VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
5161:   }
5162:   else {
5163:     if (kctx->version == 1) {
5164:       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
5165:       if (rtol < 0.0) rtol = -rtol;
5166:       stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
5167:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5168:     } else if (kctx->version == 2) {
5169:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5170:       stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
5171:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5172:     } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
5173:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5174:       /* safeguard: avoid sharp decrease of rtol */
5175:       stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
5176:       stol = PetscMax(rtol,stol);
5177:       rtol = PetscMin(kctx->rtol_0,stol);
5178:       /* safeguard: avoid oversolving */
5179:       stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
5180:       stol = PetscMax(rtol,stol);
5181:       rtol = PetscMin(kctx->rtol_0,stol);
5182:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
5183:   }
5184:   /* safeguard: avoid rtol greater than one */
5185:   rtol = PetscMin(rtol,kctx->rtol_max);
5186:   KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
5187:   PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
5188:   return(0);
5189: }

5191: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5192: {
5194:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
5195:   PCSide         pcside;
5196:   Vec            lres;

5199:   if (!snes->ksp_ewconv) return(0);
5200:   KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
5201:   kctx->norm_last = snes->norm;
5202:   if (kctx->version == 1) {
5203:     PC        pc;
5204:     PetscBool isNone;

5206:     KSPGetPC(ksp, &pc);
5207:     PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
5208:     KSPGetPCSide(ksp,&pcside);
5209:      if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
5210:       /* KSP residual is true linear residual */
5211:       KSPGetResidualNorm(ksp,&kctx->lresid_last);
5212:     } else {
5213:       /* KSP residual is preconditioned residual */
5214:       /* compute true linear residual norm */
5215:       VecDuplicate(b,&lres);
5216:       MatMult(snes->jacobian,x,lres);
5217:       VecAYPX(lres,-1.0,b);
5218:       VecNorm(lres,NORM_2,&kctx->lresid_last);
5219:       VecDestroy(&lres);
5220:     }
5221:   }
5222:   return(0);
5223: }

5225: /*@
5226:    SNESGetKSP - Returns the KSP context for a SNES solver.

5228:    Not Collective, but if SNES object is parallel, then KSP object is parallel

5230:    Input Parameter:
5231: .  snes - the SNES context

5233:    Output Parameter:
5234: .  ksp - the KSP context

5236:    Notes:
5237:    The user can then directly manipulate the KSP context to set various
5238:    options, etc.  Likewise, the user can then extract and manipulate the
5239:    PC contexts as well.

5241:    Level: beginner

5243: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
5244: @*/
5245: PetscErrorCode  SNESGetKSP(SNES snes,KSP *ksp)
5246: {


5253:   if (!snes->ksp) {
5254:     PetscBool monitor = PETSC_FALSE;

5256:     KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
5257:     PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
5258:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);

5260:     KSPSetPreSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPreSolve_SNESEW,snes);
5261:     KSPSetPostSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPostSolve_SNESEW,snes);

5263:     PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
5264:     if (monitor) {
5265:       KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
5266:     }
5267:     monitor = PETSC_FALSE;
5268:     PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
5269:     if (monitor) {
5270:       PetscObject *objs;
5271:       KSPMonitorSNESLGResidualNormCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
5272:       objs[0] = (PetscObject) snes;
5273:       KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
5274:     }
5275:     PetscObjectSetOptions((PetscObject)snes->ksp,((PetscObject)snes)->options);
5276:   }
5277:   *ksp = snes->ksp;
5278:   return(0);
5279: }


5282:  #include <petsc/private/dmimpl.h>
5283: /*@
5284:    SNESSetDM - Sets the DM that may be used by some nonlinear solvers or their underlying preconditioners

5286:    Logically Collective on SNES

5288:    Input Parameters:
5289: +  snes - the nonlinear solver context
5290: -  dm - the dm, cannot be NULL

5292:    Notes:
5293:    A DM can only be used for solving one problem at a time because information about the problem is stored on the DM,
5294:    even when not using interfaces like DMSNESSetFunction().  Use DMClone() to get a distinct DM when solving different
5295:    problems using the same function space.

5297:    Level: intermediate

5299: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
5300: @*/
5301: PetscErrorCode  SNESSetDM(SNES snes,DM dm)
5302: {
5304:   KSP            ksp;
5305:   DMSNES         sdm;

5310:   PetscObjectReference((PetscObject)dm);
5311:   if (snes->dm) {               /* Move the DMSNES context over to the new DM unless the new DM already has one */
5312:     if (snes->dm->dmsnes && !dm->dmsnes) {
5313:       DMCopyDMSNES(snes->dm,dm);
5314:       DMGetDMSNES(snes->dm,&sdm);
5315:       if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
5316:     }
5317:     DMCoarsenHookRemove(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
5318:     DMDestroy(&snes->dm);
5319:   }
5320:   snes->dm     = dm;
5321:   snes->dmAuto = PETSC_FALSE;

5323:   SNESGetKSP(snes,&ksp);
5324:   KSPSetDM(ksp,dm);
5325:   KSPSetDMActive(ksp,PETSC_FALSE);
5326:   if (snes->npc) {
5327:     SNESSetDM(snes->npc, snes->dm);
5328:     SNESSetNPCSide(snes,snes->npcside);
5329:   }
5330:   return(0);
5331: }

5333: /*@
5334:    SNESGetDM - Gets the DM that may be used by some preconditioners

5336:    Not Collective but DM obtained is parallel on SNES

5338:    Input Parameter:
5339: . snes - the preconditioner context

5341:    Output Parameter:
5342: .  dm - the dm

5344:    Level: intermediate

5346: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
5347: @*/
5348: PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
5349: {

5354:   if (!snes->dm) {
5355:     DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
5356:     snes->dmAuto = PETSC_TRUE;
5357:   }
5358:   *dm = snes->dm;
5359:   return(0);
5360: }

5362: /*@
5363:   SNESSetNPC - Sets the nonlinear preconditioner to be used.

5365:   Collective on SNES

5367:   Input Parameters:
5368: + snes - iterative context obtained from SNESCreate()
5369: - pc   - the preconditioner object

5371:   Notes:
5372:   Use SNESGetNPC() to retrieve the preconditioner context (for example,
5373:   to configure it using the API).

5375:   Level: developer

5377: .seealso: SNESGetNPC(), SNESHasNPC()
5378: @*/
5379: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
5380: {

5387:   PetscObjectReference((PetscObject) pc);
5388:   SNESDestroy(&snes->npc);
5389:   snes->npc = pc;
5390:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->npc);
5391:   return(0);
5392: }

5394: /*@
5395:   SNESGetNPC - Creates a nonlinear preconditioning solver (SNES) to be used to precondition the nonlinear solver.

5397:   Not Collective; but any changes to the obtained SNES object must be applied collectively

5399:   Input Parameter:
5400: . snes - iterative context obtained from SNESCreate()

5402:   Output Parameter:
5403: . pc - preconditioner context

5405:   Options Database:
5406: . -npc_snes_type <type> - set the type of the SNES to use as the nonlinear preconditioner

5408:   Notes:
5409:     If a SNES was previously set with SNESSetNPC() then that SNES is returned, otherwise a new SNES object is created.

5411:     The (preconditioner) SNES returned automatically inherits the same nonlinear function and Jacobian supplied to the original
5412:     SNES during SNESSetUp()

5414:   Level: developer

5416: .seealso: SNESSetNPC(), SNESHasNPC(), SNES, SNESCreate()
5417: @*/
5418: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
5419: {
5421:   const char     *optionsprefix;

5426:   if (!snes->npc) {
5427:     SNESCreate(PetscObjectComm((PetscObject)snes),&snes->npc);
5428:     PetscObjectIncrementTabLevel((PetscObject)snes->npc,(PetscObject)snes,1);
5429:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->npc);
5430:     SNESGetOptionsPrefix(snes,&optionsprefix);
5431:     SNESSetOptionsPrefix(snes->npc,optionsprefix);
5432:     SNESAppendOptionsPrefix(snes->npc,"npc_");
5433:     SNESSetCountersReset(snes->npc,PETSC_FALSE);
5434:   }
5435:   *pc = snes->npc;
5436:   return(0);
5437: }

5439: /*@
5440:   SNESHasNPC - Returns whether a nonlinear preconditioner exists

5442:   Not Collective

5444:   Input Parameter:
5445: . snes - iterative context obtained from SNESCreate()

5447:   Output Parameter:
5448: . has_npc - whether the SNES has an NPC or not

5450:   Level: developer

5452: .seealso: SNESSetNPC(), SNESGetNPC()
5453: @*/
5454: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
5455: {
5458:   *has_npc = (PetscBool) (snes->npc ? PETSC_TRUE : PETSC_FALSE);
5459:   return(0);
5460: }

5462: /*@
5463:     SNESSetNPCSide - Sets the preconditioning side.

5465:     Logically Collective on SNES

5467:     Input Parameter:
5468: .   snes - iterative context obtained from SNESCreate()

5470:     Output Parameter:
5471: .   side - the preconditioning side, where side is one of
5472: .vb
5473:       PC_LEFT - left preconditioning
5474:       PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5475: .ve

5477:     Options Database Keys:
5478: .   -snes_pc_side <right,left>

5480:     Notes:
5481:     SNESNRICHARDSON and SNESNCG only support left preconditioning.

5483:     Level: intermediate

5485: .seealso: SNESGetNPCSide(), KSPSetPCSide()
5486: @*/
5487: PetscErrorCode  SNESSetNPCSide(SNES snes,PCSide side)
5488: {
5492:   snes->npcside= side;
5493:   return(0);
5494: }

5496: /*@
5497:     SNESGetNPCSide - Gets the preconditioning side.

5499:     Not Collective

5501:     Input Parameter:
5502: .   snes - iterative context obtained from SNESCreate()

5504:     Output Parameter:
5505: .   side - the preconditioning side, where side is one of
5506: .vb
5507:       PC_LEFT - left preconditioning
5508:       PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5509: .ve

5511:     Level: intermediate

5513: .seealso: SNESSetNPCSide(), KSPGetPCSide()
5514: @*/
5515: PetscErrorCode  SNESGetNPCSide(SNES snes,PCSide *side)
5516: {
5520:   *side = snes->npcside;
5521:   return(0);
5522: }

5524: /*@
5525:   SNESSetLineSearch - Sets the linesearch on the SNES instance.

5527:   Collective on SNES

5529:   Input Parameters:
5530: + snes - iterative context obtained from SNESCreate()
5531: - linesearch   - the linesearch object

5533:   Notes:
5534:   Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
5535:   to configure it using the API).

5537:   Level: developer

5539: .seealso: SNESGetLineSearch()
5540: @*/
5541: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
5542: {

5549:   PetscObjectReference((PetscObject) linesearch);
5550:   SNESLineSearchDestroy(&snes->linesearch);

5552:   snes->linesearch = linesearch;

5554:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5555:   return(0);
5556: }

5558: /*@
5559:   SNESGetLineSearch - Returns a pointer to the line search context set with SNESSetLineSearch()
5560:   or creates a default line search instance associated with the SNES and returns it.

5562:   Not Collective

5564:   Input Parameter:
5565: . snes - iterative context obtained from SNESCreate()

5567:   Output Parameter:
5568: . linesearch - linesearch context

5570:   Level: beginner

5572: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5573: @*/
5574: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5575: {
5577:   const char     *optionsprefix;

5582:   if (!snes->linesearch) {
5583:     SNESGetOptionsPrefix(snes, &optionsprefix);
5584:     SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5585:     SNESLineSearchSetSNES(snes->linesearch, snes);
5586:     SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5587:     PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5588:     PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5589:   }
5590:   *linesearch = snes->linesearch;
5591:   return(0);
5592: }