Actual source code: snes.c

petsc-3.12.5 2020-03-29
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: PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES,Vec,Mat,Mat,void*);

313: /*@C
314:    SNESView - Prints the SNES data structure.

316:    Collective on SNES

318:    Input Parameters:
319: +  SNES - the SNES context
320: -  viewer - visualization context

322:    Options Database Key:
323: .  -snes_view - Calls SNESView() at end of SNESSolve()

325:    Notes:
326:    The available visualization contexts include
327: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
328: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
329:          output where only the first processor opens
330:          the file.  All other processors send their
331:          data to the first processor to print.

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

336:    Level: beginner

338: .seealso: PetscViewerASCIIOpen()
339: @*/
340: PetscErrorCode  SNESView(SNES snes,PetscViewer viewer)
341: {
342:   SNESKSPEW      *kctx;
344:   KSP            ksp;
345:   SNESLineSearch linesearch;
346:   PetscBool      iascii,isstring,isbinary,isdraw;
347:   DMSNES         dmsnes;
348: #if defined(PETSC_HAVE_SAWS)
349:   PetscBool      issaws;
350: #endif

354:   if (!viewer) {
355:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
356:   }

360:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
361:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
362:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
363:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
364: #if defined(PETSC_HAVE_SAWS)
365:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
366: #endif
367:   if (iascii) {
368:     SNESNormSchedule normschedule;
369:     DM               dm;
370:     PetscErrorCode   (*cJ)(SNES,Vec,Mat,Mat,void*);
371:     void             *ctx;
372:     const char       *pre = "";

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

445:     PetscObjectGetComm((PetscObject)snes,&comm);
446:     MPI_Comm_rank(comm,&rank);
447:     if (!rank) {
448:       PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
449:       PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));
450:       PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR,PETSC_FALSE);
451:     }
452:     if (snes->ops->view) {
453:       (*snes->ops->view)(snes,viewer);
454:     }
455:   } else if (isdraw) {
456:     PetscDraw draw;
457:     char      str[36];
458:     PetscReal x,y,bottom,h;

460:     PetscViewerDrawGetDraw(viewer,0,&draw);
461:     PetscDrawGetCurrentPoint(draw,&x,&y);
462:     PetscStrncpy(str,"SNES: ",sizeof(str));
463:     PetscStrlcat(str,((PetscObject)snes)->type_name,sizeof(str));
464:     PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);
465:     bottom = y - h;
466:     PetscDrawPushCurrentPoint(draw,x,bottom);
467:     if (snes->ops->view) {
468:       (*snes->ops->view)(snes,viewer);
469:     }
470: #if defined(PETSC_HAVE_SAWS)
471:   } else if (issaws) {
472:     PetscMPIInt rank;
473:     const char *name;

475:     PetscObjectGetName((PetscObject)snes,&name);
476:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
477:     if (!((PetscObject)snes)->amsmem && !rank) {
478:       char       dir[1024];

480:       PetscObjectViewSAWs((PetscObject)snes,viewer);
481:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
482:       PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT));
483:       if (!snes->conv_hist) {
484:         SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);
485:       }
486:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);
487:       PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE));
488:     }
489: #endif
490:   }
491:   if (snes->linesearch) {
492:     SNESGetLineSearch(snes, &linesearch);
493:     PetscViewerASCIIPushTab(viewer);
494:     SNESLineSearchView(linesearch, viewer);
495:     PetscViewerASCIIPopTab(viewer);
496:   }
497:   if (snes->npc && snes->usesnpc) {
498:     PetscViewerASCIIPushTab(viewer);
499:     SNESView(snes->npc, viewer);
500:     PetscViewerASCIIPopTab(viewer);
501:   }
502:   PetscViewerASCIIPushTab(viewer);
503:   DMGetDMSNES(snes->dm,&dmsnes);
504:   DMSNESView(dmsnes, viewer);
505:   PetscViewerASCIIPopTab(viewer);
506:   if (snes->usesksp) {
507:     SNESGetKSP(snes,&ksp);
508:     PetscViewerASCIIPushTab(viewer);
509:     KSPView(ksp,viewer);
510:     PetscViewerASCIIPopTab(viewer);
511:   }
512:   if (isdraw) {
513:     PetscDraw draw;
514:     PetscViewerDrawGetDraw(viewer,0,&draw);
515:     PetscDrawPopCurrentPoint(draw);
516:   }
517:   return(0);
518: }

520: /*
521:   We retain a list of functions that also take SNES command
522:   line options. These are called at the end SNESSetFromOptions()
523: */
524: #define MAXSETFROMOPTIONS 5
525: static PetscInt numberofsetfromoptions;
526: static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);

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

531:   Not Collective

533:   Input Parameter:
534: . snescheck - function that checks for options

536:   Level: developer

538: .seealso: SNESSetFromOptions()
539: @*/
540: PetscErrorCode  SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
541: {
543:   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
544:   othersetfromoptions[numberofsetfromoptions++] = snescheck;
545:   return(0);
546: }

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

550: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
551: {
552:   Mat            J;
554:   MatNullSpace   nullsp;


559:   if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
560:     Mat A = snes->jacobian, B = snes->jacobian_pre;
561:     MatCreateVecs(A ? A : B, NULL,&snes->vec_func);
562:   }

564:   if (version == 1) {
565:     MatCreateSNESMF(snes,&J);
566:     MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
567:     MatSetFromOptions(J);
568:   } else if (version == 2) {
569:     if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
570: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
571:     SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);
572: #else
573:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)");
574: #endif
575:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2");

577:   /* attach any user provided null space that was on Amat to the newly created matrix free matrix */
578:   if (snes->jacobian) {
579:     MatGetNullSpace(snes->jacobian,&nullsp);
580:     if (nullsp) {
581:       MatSetNullSpace(J,nullsp);
582:     }
583:   }

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

588:     /* This version replaces the user provided Jacobian matrix with a
589:        matrix-free version but still employs the user-provided preconditioner matrix. */
590:     SNESSetJacobian(snes,J,0,0,0);
591:   } else {
592:     /* This version replaces both the user-provided Jacobian and the user-
593:      provided preconditioner Jacobian with the default matrix free version. */
594:     if ((snes->npcside== PC_LEFT) && snes->npc) {
595:       if (!snes->jacobian){SNESSetJacobian(snes,J,0,0,0);}
596:     } else {
597:       KSP       ksp;
598:       PC        pc;
599:       PetscBool match;

601:       SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,0);
602:       /* Force no preconditioner */
603:       SNESGetKSP(snes,&ksp);
604:       KSPGetPC(ksp,&pc);
605:       PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);
606:       if (!match) {
607:         PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");
608:         PCSetType(pc,PCNONE);
609:       }
610:     }
611:   }
612:   MatDestroy(&J);
613:   return(0);
614: }

616: static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx)
617: {
618:   SNES           snes = (SNES)ctx;
620:   Vec            Xfine,Xfine_named = NULL,Xcoarse;

623:   if (PetscLogPrintInfo) {
624:     PetscInt finelevel,coarselevel,fineclevel,coarseclevel;
625:     DMGetRefineLevel(dmfine,&finelevel);
626:     DMGetCoarsenLevel(dmfine,&fineclevel);
627:     DMGetRefineLevel(dmcoarse,&coarselevel);
628:     DMGetCoarsenLevel(dmcoarse,&coarseclevel);
629:     PetscInfo4(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);
630:   }
631:   if (dmfine == snes->dm) Xfine = snes->vec_sol;
632:   else {
633:     DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
634:     Xfine = Xfine_named;
635:   }
636:   DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
637:   if (Inject) {
638:     MatRestrict(Inject,Xfine,Xcoarse);
639:   } else {
640:     MatRestrict(Restrict,Xfine,Xcoarse);
641:     VecPointwiseMult(Xcoarse,Xcoarse,Rscale);
642:   }
643:   DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
644:   if (Xfine_named) {DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);}
645:   return(0);
646: }

648: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
649: {

653:   DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
654:   return(0);
655: }

657: /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can
658:  * safely call SNESGetDM() in their residual evaluation routine. */
659: static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx)
660: {
661:   SNES           snes = (SNES)ctx;
663:   Vec            X,Xnamed = NULL;
664:   DM             dmsave;
665:   void           *ctxsave;
666:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL;

669:   dmsave = snes->dm;
670:   KSPGetDM(ksp,&snes->dm);
671:   if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */
672:   else {                                     /* We are on a coarser level, this vec was initialized using a DM restrict hook */
673:     DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
674:     X    = Xnamed;
675:     SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);
676:     /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */
677:     if (jac == SNESComputeJacobianDefaultColor) {
678:       SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,0);
679:     }
680:   }
681:   /* Make sure KSP DM has the Jacobian computation routine */
682:   {
683:     DMSNES sdm;

685:     DMGetDMSNES(snes->dm, &sdm);
686:     if (!sdm->ops->computejacobian) {
687:       DMCopyDMSNES(dmsave, snes->dm);
688:     }
689:   }
690:   /* Compute the operators */
691:   SNESComputeJacobian(snes,X,A,B);
692:   /* Put the previous context back */
693:   if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
694:     SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
695:   }

697:   if (Xnamed) {DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);}
698:   snes->dm = dmsave;
699:   return(0);
700: }

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

705:    Collective

707:    Input Arguments:
708: .  snes - snes to configure

710:    Level: developer

712: .seealso: SNESSetUp()
713: @*/
714: PetscErrorCode SNESSetUpMatrices(SNES snes)
715: {
717:   DM             dm;
718:   DMSNES         sdm;

721:   SNESGetDM(snes,&dm);
722:   DMGetDMSNES(dm,&sdm);
723:   if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"DMSNES not properly configured");
724:   else if (!snes->jacobian && snes->mf) {
725:     Mat  J;
726:     void *functx;
727:     MatCreateSNESMF(snes,&J);
728:     MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
729:     MatSetFromOptions(J);
730:     SNESGetFunction(snes,NULL,NULL,&functx);
731:     SNESSetJacobian(snes,J,J,0,0);
732:     MatDestroy(&J);
733:   } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
734:     Mat J,B;
735:     MatCreateSNESMF(snes,&J);
736:     MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
737:     MatSetFromOptions(J);
738:     DMCreateMatrix(snes->dm,&B);
739:     /* sdm->computejacobian was already set to reach here */
740:     SNESSetJacobian(snes,J,B,NULL,NULL);
741:     MatDestroy(&J);
742:     MatDestroy(&B);
743:   } else if (!snes->jacobian_pre) {
744:     PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *);
745:     PetscDS          prob;
746:     Mat              J, B;
747:     MatNullSpace     nullspace = NULL;
748:     PetscBool        hasPrec   = PETSC_FALSE;
749:     PetscInt         Nf;

751:     J    = snes->jacobian;
752:     DMGetDS(dm, &prob);
753:     if (prob) {PetscDSHasJacobianPreconditioner(prob, &hasPrec);}
754:     if (J)            {PetscObjectReference((PetscObject) J);}
755:     else if (hasPrec) {DMCreateMatrix(snes->dm, &J);}
756:     DMCreateMatrix(snes->dm, &B);
757:     PetscDSGetNumFields(prob, &Nf);
758:     DMGetNullSpaceConstructor(snes->dm, Nf, &nspconstr);
759:     if (nspconstr) (*nspconstr)(snes->dm, -1, &nullspace);
760:     MatSetNullSpace(B, nullspace);
761:     MatNullSpaceDestroy(&nullspace);
762:     SNESSetJacobian(snes, J ? J : B, B, NULL, NULL);
763:     MatDestroy(&J);
764:     MatDestroy(&B);
765:   }
766:   {
767:     KSP ksp;
768:     SNESGetKSP(snes,&ksp);
769:     KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);
770:     DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
771:   }
772:   return(0);
773: }

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

778:    Collective on SNES

780:    Input Parameters:
781: +  snes - SNES object you wish to monitor
782: .  name - the monitor type one is seeking
783: .  help - message indicating what monitoring is done
784: .  manual - manual page for the monitor
785: .  monitor - the monitor function
786: -  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

788:    Level: developer

790: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
791:           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
792:           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
793:           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
794:           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
795:           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
796:           PetscOptionsFList(), PetscOptionsEList()
797: @*/
798: PetscErrorCode  SNESMonitorSetFromOptions(SNES snes,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNES,PetscViewerAndFormat*))
799: {
800:   PetscErrorCode    ierr;
801:   PetscViewer       viewer;
802:   PetscViewerFormat format;
803:   PetscBool         flg;

806:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,name,&viewer,&format,&flg);
807:   if (flg) {
808:     PetscViewerAndFormat *vf;
809:     PetscViewerAndFormatCreate(viewer,format,&vf);
810:     PetscObjectDereference((PetscObject)viewer);
811:     if (monitorsetup) {
812:       (*monitorsetup)(snes,vf);
813:     }
814:     SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
815:   }
816:   return(0);
817: }

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

822:    Collective on SNES

824:    Input Parameter:
825: .  snes - the SNES context

827:    Options Database Keys:
828: +  -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
829: .  -snes_stol - convergence tolerance in terms of the norm
830:                 of the change in the solution between steps
831: .  -snes_atol <abstol> - absolute tolerance of residual norm
832: .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
833: .  -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
834: .  -snes_force_iteration <force> - force SNESSolve() to take at least one iteration
835: .  -snes_max_it <max_it> - maximum number of iterations
836: .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
837: .  -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
838: .  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
839: .  -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
840: .  -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
841: .  -snes_trtol <trtol> - trust region tolerance
842: .  -snes_no_convergence_test - skip convergence test in nonlinear
843:                                solver; hence iterations will continue until max_it
844:                                or some other criterion is reached. Saves expense
845:                                of convergence test
846: .  -snes_monitor [ascii][:filename][:viewer format] - prints residual norm at each iteration. if no filename given prints to stdout
847: .  -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration
848: .  -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration
849: .  -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration
850: .  -snes_monitor_lg_residualnorm - plots residual norm at each iteration
851: .  -snes_monitor_lg_range - plots residual norm at each iteration
852: .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
853: .  -snes_fd_color - use finite differences with coloring to compute Jacobian
854: .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
855: .  -snes_converged_reason - print the reason for convergence/divergence after each solve
856: -  -npc_snes_type <type> - the SNES type to use as a nonlinear preconditioner

858:     Options Database for Eisenstat-Walker method:
859: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
860: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
861: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
862: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
863: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
864: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
865: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
866: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

868:    Notes:
869:    To see all options, run your program with the -help option or consult
870:    Users-Manual: ch_snes

872:    Level: beginner

874: .seealso: SNESSetOptionsPrefix(), SNESResetFromOptions()
875: @*/
876: PetscErrorCode  SNESSetFromOptions(SNES snes)
877: {
878:   PetscBool      flg,pcset,persist,set;
879:   PetscInt       i,indx,lag,grids;
880:   const char     *deft        = SNESNEWTONLS;
881:   const char     *convtests[] = {"default","skip"};
882:   SNESKSPEW      *kctx        = NULL;
883:   char           type[256], monfilename[PETSC_MAX_PATH_LEN];
885:   PCSide         pcside;
886:   const char     *optionsprefix;

890:   SNESRegisterAll();
891:   PetscObjectOptionsBegin((PetscObject)snes);
892:   if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
893:   PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
894:   if (flg) {
895:     SNESSetType(snes,type);
896:   } else if (!((PetscObject)snes)->type_name) {
897:     SNESSetType(snes,deft);
898:   }
899:   PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);
900:   PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);

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

912:   PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
913:   if (flg) {
914:     SNESSetLagPreconditioner(snes,lag);
915:   }
916:   PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
917:   if (flg) {
918:     SNESSetLagPreconditionerPersists(snes,persist);
919:   }
920:   PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
921:   if (flg) {
922:     SNESSetLagJacobian(snes,lag);
923:   }
924:   PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
925:   if (flg) {
926:     SNESSetLagJacobianPersists(snes,persist);
927:   }

929:   PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
930:   if (flg) {
931:     SNESSetGridSequence(snes,grids);
932:   }

934:   PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);
935:   if (flg) {
936:     switch (indx) {
937:     case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
938:     case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL);    break;
939:     }
940:   }

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

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

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

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

952:   PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);
953:   PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);
954:   PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);
955:   PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);
956:   PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);
957:   PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);
958:   PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);

960:   flg  = PETSC_FALSE;
961:   PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);
962:   if (set && flg) {SNESMonitorCancel(snes);}

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

968:   SNESMonitorSetFromOptions(snes,"-snes_monitor_ratio","Monitor ratios of the norm of function for consecutive steps","SNESMonitorRatio",SNESMonitorRatio,SNESMonitorRatioSetUp);
969:   SNESMonitorSetFromOptions(snes,"-snes_monitor_field","Monitor norm of function (split into fields)","SNESMonitorDefaultField",SNESMonitorDefaultField,NULL);
970:   SNESMonitorSetFromOptions(snes,"-snes_monitor_solution","View solution at each iteration","SNESMonitorSolution",SNESMonitorSolution,NULL);
971:   SNESMonitorSetFromOptions(snes,"-snes_monitor_solution_update","View correction at each iteration","SNESMonitorSolutionUpdate",SNESMonitorSolutionUpdate,NULL);
972:   SNESMonitorSetFromOptions(snes,"-snes_monitor_residual","View residual at each iteration","SNESMonitorResidual",SNESMonitorResidual,NULL);
973:   SNESMonitorSetFromOptions(snes,"-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",SNESMonitorJacUpdateSpectrum,NULL);
974:   SNESMonitorSetFromOptions(snes,"-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet",SNESMonitorFields,NULL);

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

979:   flg  = PETSC_FALSE;
980:   PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);
981:   if (flg) {
982:     PetscDrawLG ctx;

984:     SNESMonitorLGCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
985:     SNESMonitorSet(snes,SNESMonitorLGResidualNorm,ctx,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
986:   }
987:   flg  = PETSC_FALSE;
988:   PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
989:   if (flg) {
990:     PetscViewer ctx;

992:     PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
993:     SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
994:   }

996:   flg  = PETSC_FALSE;
997:   PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
998:   if (flg) {
999:     void    *functx;
1000:     DM      dm;
1001:     DMSNES  sdm;
1002:     SNESGetDM(snes,&dm);
1003:     DMGetDMSNES(dm,&sdm);
1004:     sdm->jacobianctx = NULL;
1005:     SNESGetFunction(snes,NULL,NULL,&functx);
1006:     SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
1007:     PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
1008:   }

1010:   flg  = PETSC_FALSE;
1011:   PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
1012:   if (flg) {
1013:     SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
1014:   }

1016:   flg  = PETSC_FALSE;
1017:   PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);
1018:   if (flg) {
1019:     DM             dm;
1020:     DMSNES         sdm;
1021:     SNESGetDM(snes,&dm);
1022:     DMGetDMSNES(dm,&sdm);
1023:     sdm->jacobianctx = NULL;
1024:     SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,0);
1025:     PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
1026:   }

1028:   flg  = PETSC_FALSE;
1029:   PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf_operator,&flg);
1030:   if (flg && snes->mf_operator) {
1031:     snes->mf_operator = PETSC_TRUE;
1032:     snes->mf          = PETSC_TRUE;
1033:   }
1034:   flg  = PETSC_FALSE;
1035:   PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf,&flg);
1036:   if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE;
1037:   PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,0);

1039:   flg  = PETSC_FALSE;
1040:   SNESGetNPCSide(snes,&pcside);
1041:   PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
1042:   if (flg) {SNESSetNPCSide(snes,pcside);}

1044: #if defined(PETSC_HAVE_SAWS)
1045:   /*
1046:     Publish convergence information using SAWs
1047:   */
1048:   flg  = PETSC_FALSE;
1049:   PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);
1050:   if (flg) {
1051:     void *ctx;
1052:     SNESMonitorSAWsCreate(snes,&ctx);
1053:     SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);
1054:   }
1055: #endif
1056: #if defined(PETSC_HAVE_SAWS)
1057:   {
1058:   PetscBool set;
1059:   flg  = PETSC_FALSE;
1060:   PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);
1061:   if (set) {
1062:     PetscObjectSAWsSetBlock((PetscObject)snes,flg);
1063:   }
1064:   }
1065: #endif

1067:   for (i = 0; i < numberofsetfromoptions; i++) {
1068:     (*othersetfromoptions[i])(snes);
1069:   }

1071:   if (snes->ops->setfromoptions) {
1072:     (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
1073:   }

1075:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
1076:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)snes);
1077:   PetscOptionsEnd();

1079:   if (snes->linesearch) {
1080:     SNESGetLineSearch(snes, &snes->linesearch);
1081:     SNESLineSearchSetFromOptions(snes->linesearch);
1082:   }

1084:   if (snes->usesksp) {
1085:     if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
1086:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
1087:     KSPSetFromOptions(snes->ksp);
1088:   }

1090:   /* if user has set the SNES NPC type via options database, create it. */
1091:   SNESGetOptionsPrefix(snes, &optionsprefix);
1092:   PetscOptionsHasName(((PetscObject)snes)->options,optionsprefix, "-npc_snes_type", &pcset);
1093:   if (pcset && (!snes->npc)) {
1094:     SNESGetNPC(snes, &snes->npc);
1095:   }
1096:   if (snes->npc) {
1097:     SNESSetFromOptions(snes->npc);
1098:   }
1099:   snes->setfromoptionscalled++;
1100:   return(0);
1101: }

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

1106:    Collective on SNES

1108:    Input Parameter:
1109: .  snes - the SNES context

1111:    Level: beginner

1113: .seealso: SNESSetFromOptions(), SNESSetOptionsPrefix()
1114: @*/
1115: PetscErrorCode SNESResetFromOptions(SNES snes)
1116: {

1120:   if (snes->setfromoptionscalled) {SNESSetFromOptions(snes);}
1121:   return(0);
1122: }

1124: /*@C
1125:    SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
1126:    the nonlinear solvers.

1128:    Logically Collective on SNES

1130:    Input Parameters:
1131: +  snes - the SNES context
1132: .  compute - function to compute the context
1133: -  destroy - function to destroy the context

1135:    Level: intermediate

1137:    Notes:
1138:    This function is currently not available from Fortran.

1140: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
1141: @*/
1142: PetscErrorCode  SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
1143: {
1146:   snes->ops->usercompute = compute;
1147:   snes->ops->userdestroy = destroy;
1148:   return(0);
1149: }

1151: /*@
1152:    SNESSetApplicationContext - Sets the optional user-defined context for
1153:    the nonlinear solvers.

1155:    Logically Collective on SNES

1157:    Input Parameters:
1158: +  snes - the SNES context
1159: -  usrP - optional user context

1161:    Level: intermediate

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

1167: .seealso: SNESGetApplicationContext()
1168: @*/
1169: PetscErrorCode  SNESSetApplicationContext(SNES snes,void *usrP)
1170: {
1172:   KSP            ksp;

1176:   SNESGetKSP(snes,&ksp);
1177:   KSPSetApplicationContext(ksp,usrP);
1178:   snes->user = usrP;
1179:   return(0);
1180: }

1182: /*@
1183:    SNESGetApplicationContext - Gets the user-defined context for the
1184:    nonlinear solvers.

1186:    Not Collective

1188:    Input Parameter:
1189: .  snes - SNES context

1191:    Output Parameter:
1192: .  usrP - user context

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

1198:    Level: intermediate

1200: .seealso: SNESSetApplicationContext()
1201: @*/
1202: PetscErrorCode  SNESGetApplicationContext(SNES snes,void *usrP)
1203: {
1206:   *(void**)usrP = snes->user;
1207:   return(0);
1208: }

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

1214:    Collective on SNES

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

1221:    Options Database:
1222: + -snes_mf - use matrix free for both the mat and pmat operator
1223: - -snes_mf_operator - use matrix free only for the mat operator

1225:    Level: intermediate

1227: .seealso:   SNESGetUseMatrixFree(), MatCreateSNESMF()
1228: @*/
1229: PetscErrorCode  SNESSetUseMatrixFree(SNES snes,PetscBool mf_operator,PetscBool mf)
1230: {
1235:   if (mf && !mf_operator) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"If using mf must also use mf_operator");
1236:   snes->mf          = mf;
1237:   snes->mf_operator = mf_operator;
1238:   return(0);
1239: }

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

1245:    Collective on SNES

1247:    Input Parameter:
1248: .  snes - SNES context

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

1254:    Options Database:
1255: + -snes_mf - use matrix free for both the mat and pmat operator
1256: - -snes_mf_operator - use matrix free only for the mat operator

1258:    Level: intermediate

1260: .seealso:   SNESSetUseMatrixFree(), MatCreateSNESMF()
1261: @*/
1262: PetscErrorCode  SNESGetUseMatrixFree(SNES snes,PetscBool *mf_operator,PetscBool *mf)
1263: {
1266:   if (mf)          *mf          = snes->mf;
1267:   if (mf_operator) *mf_operator = snes->mf_operator;
1268:   return(0);
1269: }

1271: /*@
1272:    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1273:    at this time.

1275:    Not Collective

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

1280:    Output Parameter:
1281: .  iter - iteration number

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

1286:    This is useful for using lagged Jacobians (where one does not recompute the
1287:    Jacobian at each SNES iteration). For example, the code
1288: .vb
1289:       SNESGetIterationNumber(snes,&it);
1290:       if (!(it % 2)) {
1291:         [compute Jacobian here]
1292:       }
1293: .ve
1294:    can be used in your ComputeJacobian() function to cause the Jacobian to be
1295:    recomputed every second SNES iteration.

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

1299:    Level: intermediate

1301: .seealso:   SNESGetLinearSolveIterations()
1302: @*/
1303: PetscErrorCode  SNESGetIterationNumber(SNES snes,PetscInt *iter)
1304: {
1308:   *iter = snes->iter;
1309:   return(0);
1310: }

1312: /*@
1313:    SNESSetIterationNumber - Sets the current iteration number.

1315:    Not Collective

1317:    Input Parameter:
1318: +  snes - SNES context
1319: -  iter - iteration number

1321:    Level: developer

1323: .seealso:   SNESGetLinearSolveIterations()
1324: @*/
1325: PetscErrorCode  SNESSetIterationNumber(SNES snes,PetscInt iter)
1326: {

1331:   PetscObjectSAWsTakeAccess((PetscObject)snes);
1332:   snes->iter = iter;
1333:   PetscObjectSAWsGrantAccess((PetscObject)snes);
1334:   return(0);
1335: }

1337: /*@
1338:    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1339:    attempted by the nonlinear solver.

1341:    Not Collective

1343:    Input Parameter:
1344: .  snes - SNES context

1346:    Output Parameter:
1347: .  nfails - number of unsuccessful steps attempted

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

1352:    Level: intermediate

1354: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1355:           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1356: @*/
1357: PetscErrorCode  SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1358: {
1362:   *nfails = snes->numFailures;
1363:   return(0);
1364: }

1366: /*@
1367:    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1368:    attempted by the nonlinear solver before it gives up.

1370:    Not Collective

1372:    Input Parameters:
1373: +  snes     - SNES context
1374: -  maxFails - maximum of unsuccessful steps

1376:    Level: intermediate

1378: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1379:           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1380: @*/
1381: PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1382: {
1385:   snes->maxFailures = maxFails;
1386:   return(0);
1387: }

1389: /*@
1390:    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1391:    attempted by the nonlinear solver before it gives up.

1393:    Not Collective

1395:    Input Parameter:
1396: .  snes     - SNES context

1398:    Output Parameter:
1399: .  maxFails - maximum of unsuccessful steps

1401:    Level: intermediate

1403: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1404:           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()

1406: @*/
1407: PetscErrorCode  SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1408: {
1412:   *maxFails = snes->maxFailures;
1413:   return(0);
1414: }

1416: /*@
1417:    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1418:      done by SNES.

1420:    Not Collective

1422:    Input Parameter:
1423: .  snes     - SNES context

1425:    Output Parameter:
1426: .  nfuncs - number of evaluations

1428:    Level: intermediate

1430:    Notes:
1431:     Reset every time SNESSolve is called unless SNESSetCountersReset() is used.

1433: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1434: @*/
1435: PetscErrorCode  SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1436: {
1440:   *nfuncs = snes->nfuncs;
1441:   return(0);
1442: }

1444: /*@
1445:    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1446:    linear solvers.

1448:    Not Collective

1450:    Input Parameter:
1451: .  snes - SNES context

1453:    Output Parameter:
1454: .  nfails - number of failed solves

1456:    Level: intermediate

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

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

1464: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1465: @*/
1466: PetscErrorCode  SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1467: {
1471:   *nfails = snes->numLinearSolveFailures;
1472:   return(0);
1473: }

1475: /*@
1476:    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1477:    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE

1479:    Logically Collective on SNES

1481:    Input Parameters:
1482: +  snes     - SNES context
1483: -  maxFails - maximum allowed linear solve failures

1485:    Level: intermediate

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

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

1493: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1494: @*/
1495: PetscErrorCode  SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1496: {
1500:   snes->maxLinearSolveFailures = maxFails;
1501:   return(0);
1502: }

1504: /*@
1505:    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1506:      are allowed before SNES terminates

1508:    Not Collective

1510:    Input Parameter:
1511: .  snes     - SNES context

1513:    Output Parameter:
1514: .  maxFails - maximum of unsuccessful solves allowed

1516:    Level: intermediate

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

1521: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1522: @*/
1523: PetscErrorCode  SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1524: {
1528:   *maxFails = snes->maxLinearSolveFailures;
1529:   return(0);
1530: }

1532: /*@
1533:    SNESGetLinearSolveIterations - Gets the total number of linear iterations
1534:    used by the nonlinear solver.

1536:    Not Collective

1538:    Input Parameter:
1539: .  snes - SNES context

1541:    Output Parameter:
1542: .  lits - number of linear iterations

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

1547:    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
1548:    then call KSPGetIterationNumber() after the failed solve.

1550:    Level: intermediate

1552: .seealso:  SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1553: @*/
1554: PetscErrorCode  SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1555: {
1559:   *lits = snes->linear_its;
1560:   return(0);
1561: }

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

1567:    Logically Collective on SNES

1569:    Input Parameter:
1570: +  snes - SNES context
1571: -  reset - whether to reset the counters or not

1573:    Notes:
1574:    This defaults to PETSC_TRUE

1576:    Level: developer

1578: .seealso:  SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC()
1579: @*/
1580: PetscErrorCode  SNESSetCountersReset(SNES snes,PetscBool reset)
1581: {
1585:   snes->counters_reset = reset;
1586:   return(0);
1587: }


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

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

1595:    Input Parameters:
1596: +  snes - the SNES context
1597: -  ksp - the KSP context

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

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

1606:    Level: developer

1608: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1609: @*/
1610: PetscErrorCode  SNESSetKSP(SNES snes,KSP ksp)
1611: {

1618:   PetscObjectReference((PetscObject)ksp);
1619:   if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1620:   snes->ksp = ksp;
1621:   return(0);
1622: }

1624: /* -----------------------------------------------------------*/
1625: /*@
1626:    SNESCreate - Creates a nonlinear solver context.

1628:    Collective

1630:    Input Parameters:
1631: .  comm - MPI communicator

1633:    Output Parameter:
1634: .  outsnes - the new SNES context

1636:    Options Database Keys:
1637: +   -snes_mf - Activates default matrix-free Jacobian-vector products,
1638:                and no preconditioning matrix
1639: .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
1640:                products, and a user-provided preconditioning matrix
1641:                as set by SNESSetJacobian()
1642: -   -snes_fd - Uses (slow!) finite differences to compute Jacobian

1644:    Level: beginner

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

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

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

1658: @*/
1659: PetscErrorCode  SNESCreate(MPI_Comm comm,SNES *outsnes)
1660: {
1662:   SNES           snes;
1663:   SNESKSPEW      *kctx;

1667:   *outsnes = NULL;
1668:   SNESInitializePackage();

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

1672:   snes->ops->converged    = SNESConvergedDefault;
1673:   snes->usesksp           = PETSC_TRUE;
1674:   snes->tolerancesset     = PETSC_FALSE;
1675:   snes->max_its           = 50;
1676:   snes->max_funcs         = 10000;
1677:   snes->norm              = 0.0;
1678:   snes->xnorm             = 0.0;
1679:   snes->ynorm             = 0.0;
1680:   snes->normschedule      = SNES_NORM_ALWAYS;
1681:   snes->functype          = SNES_FUNCTION_DEFAULT;
1682: #if defined(PETSC_USE_REAL_SINGLE)
1683:   snes->rtol              = 1.e-5;
1684: #else
1685:   snes->rtol              = 1.e-8;
1686: #endif
1687:   snes->ttol              = 0.0;
1688: #if defined(PETSC_USE_REAL_SINGLE)
1689:   snes->abstol            = 1.e-25;
1690: #else
1691:   snes->abstol            = 1.e-50;
1692: #endif
1693: #if defined(PETSC_USE_REAL_SINGLE)
1694:   snes->stol              = 1.e-5;
1695: #else
1696:   snes->stol              = 1.e-8;
1697: #endif
1698: #if defined(PETSC_USE_REAL_SINGLE)
1699:   snes->deltatol          = 1.e-6;
1700: #else
1701:   snes->deltatol          = 1.e-12;
1702: #endif
1703:   snes->divtol            = 1.e4;
1704:   snes->rnorm0            = 0;
1705:   snes->nfuncs            = 0;
1706:   snes->numFailures       = 0;
1707:   snes->maxFailures       = 1;
1708:   snes->linear_its        = 0;
1709:   snes->lagjacobian       = 1;
1710:   snes->jac_iter          = 0;
1711:   snes->lagjac_persist    = PETSC_FALSE;
1712:   snes->lagpreconditioner = 1;
1713:   snes->pre_iter          = 0;
1714:   snes->lagpre_persist    = PETSC_FALSE;
1715:   snes->numbermonitors    = 0;
1716:   snes->data              = 0;
1717:   snes->setupcalled       = PETSC_FALSE;
1718:   snes->ksp_ewconv        = PETSC_FALSE;
1719:   snes->nwork             = 0;
1720:   snes->work              = 0;
1721:   snes->nvwork            = 0;
1722:   snes->vwork             = 0;
1723:   snes->conv_hist_len     = 0;
1724:   snes->conv_hist_max     = 0;
1725:   snes->conv_hist         = NULL;
1726:   snes->conv_hist_its     = NULL;
1727:   snes->conv_hist_reset   = PETSC_TRUE;
1728:   snes->counters_reset    = PETSC_TRUE;
1729:   snes->vec_func_init_set = PETSC_FALSE;
1730:   snes->reason            = SNES_CONVERGED_ITERATING;
1731:   snes->npcside           = PC_RIGHT;
1732:   snes->setfromoptionscalled = 0;

1734:   snes->mf          = PETSC_FALSE;
1735:   snes->mf_operator = PETSC_FALSE;
1736:   snes->mf_version  = 1;

1738:   snes->numLinearSolveFailures = 0;
1739:   snes->maxLinearSolveFailures = 1;

1741:   snes->vizerotolerance = 1.e-8;
1742: #if defined(PETSC_USE_DEBUG)
1743:   snes->checkjacdomainerror = PETSC_TRUE;
1744: #else
1745:   snes->checkjacdomainerror = PETSC_FALSE;
1746: #endif

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

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

1754:   snes->kspconvctx  = (void*)kctx;
1755:   kctx->version     = 2;
1756:   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1757:                              this was too large for some test cases */
1758:   kctx->rtol_last   = 0.0;
1759:   kctx->rtol_max    = .9;
1760:   kctx->gamma       = 1.0;
1761:   kctx->alpha       = .5*(1.0 + PetscSqrtReal(5.0));
1762:   kctx->alpha2      = kctx->alpha;
1763:   kctx->threshold   = .1;
1764:   kctx->lresid_last = 0.0;
1765:   kctx->norm_last   = 0.0;

1767:   *outsnes = snes;
1768:   return(0);
1769: }

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

1774:      Synopsis:
1775:      #include "petscsnes.h"
1776:      PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);

1778:      Collective on snes

1780:      Input Parameters:
1781: +     snes - the SNES context
1782: .     x    - state at which to evaluate residual
1783: -     ctx     - optional user-defined function context, passed in with SNESSetFunction()

1785:      Output Parameter:
1786: .     f  - vector to put residual (function value)

1788:    Level: intermediate

1790: .seealso:   SNESSetFunction(), SNESGetFunction()
1791: M*/

1793: /*@C
1794:    SNESSetFunction - Sets the function evaluation routine and function
1795:    vector for use by the SNES routines in solving systems of nonlinear
1796:    equations.

1798:    Logically Collective on SNES

1800:    Input Parameters:
1801: +  snes - the SNES context
1802: .  r - vector to store function value
1803: .  f - function evaluation routine; see SNESFunction for calling sequence details
1804: -  ctx - [optional] user-defined context for private data for the
1805:          function evaluation routine (may be NULL)

1807:    Notes:
1808:    The Newton-like methods typically solve linear systems of the form
1809: $      f'(x) x = -f(x),
1810:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

1812:    Level: beginner

1814: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1815: @*/
1816: PetscErrorCode  SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1817: {
1819:   DM             dm;

1823:   if (r) {
1826:     PetscObjectReference((PetscObject)r);
1827:     VecDestroy(&snes->vec_func);

1829:     snes->vec_func = r;
1830:   }
1831:   SNESGetDM(snes,&dm);
1832:   DMSNESSetFunction(dm,f,ctx);
1833:   return(0);
1834: }


1837: /*@C
1838:    SNESSetInitialFunction - Sets the function vector to be used as the
1839:    function norm at the initialization of the method.  In some
1840:    instances, the user has precomputed the function before calling
1841:    SNESSolve.  This function allows one to avoid a redundant call
1842:    to SNESComputeFunction in that case.

1844:    Logically Collective on SNES

1846:    Input Parameters:
1847: +  snes - the SNES context
1848: -  f - vector to store function value

1850:    Notes:
1851:    This should not be modified during the solution procedure.

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

1855:    Level: developer

1857: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1858: @*/
1859: PetscErrorCode  SNESSetInitialFunction(SNES snes, Vec f)
1860: {
1862:   Vec            vec_func;

1868:   if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1869:     snes->vec_func_init_set = PETSC_FALSE;
1870:     return(0);
1871:   }
1872:   SNESGetFunction(snes,&vec_func,NULL,NULL);
1873:   VecCopy(f, vec_func);

1875:   snes->vec_func_init_set = PETSC_TRUE;
1876:   return(0);
1877: }

1879: /*@
1880:    SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1881:    of the SNES method.

1883:    Logically Collective on SNES

1885:    Input Parameters:
1886: +  snes - the SNES context
1887: -  normschedule - the frequency of norm computation

1889:    Options Database Key:
1890: .  -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly>

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

1901:    Level: developer

1903: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1904: @*/
1905: PetscErrorCode  SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1906: {
1909:   snes->normschedule = normschedule;
1910:   return(0);
1911: }


1914: /*@
1915:    SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1916:    of the SNES method.

1918:    Logically Collective on SNES

1920:    Input Parameters:
1921: +  snes - the SNES context
1922: -  normschedule - the type of the norm used

1924:    Level: advanced

1926: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1927: @*/
1928: PetscErrorCode  SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1929: {
1932:   *normschedule = snes->normschedule;
1933:   return(0);
1934: }


1937: /*@
1938:   SNESSetFunctionNorm - Sets the last computed residual norm.

1940:   Logically Collective on SNES

1942:   Input Parameters:
1943: + snes - the SNES context

1945: - normschedule - the frequency of norm computation

1947:   Level: developer

1949: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1950: @*/
1951: PetscErrorCode SNESSetFunctionNorm(SNES snes, PetscReal norm)
1952: {
1955:   snes->norm = norm;
1956:   return(0);
1957: }

1959: /*@
1960:   SNESGetFunctionNorm - Gets the last computed norm of the residual

1962:   Not Collective

1964:   Input Parameter:
1965: . snes - the SNES context

1967:   Output Parameter:
1968: . norm - the last computed residual norm

1970:   Level: developer

1972: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1973: @*/
1974: PetscErrorCode SNESGetFunctionNorm(SNES snes, PetscReal *norm)
1975: {
1979:   *norm = snes->norm;
1980:   return(0);
1981: }

1983: /*@
1984:   SNESGetUpdateNorm - Gets the last computed norm of the Newton update

1986:   Not Collective

1988:   Input Parameter:
1989: . snes - the SNES context

1991:   Output Parameter:
1992: . ynorm - the last computed update norm

1994:   Level: developer

1996: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm()
1997: @*/
1998: PetscErrorCode SNESGetUpdateNorm(SNES snes, PetscReal *ynorm)
1999: {
2003:   *ynorm = snes->ynorm;
2004:   return(0);
2005: }

2007: /*@
2008:   SNESGetSolutionNorm - Gets the last computed norm of the solution

2010:   Not Collective

2012:   Input Parameter:
2013: . snes - the SNES context

2015:   Output Parameter:
2016: . xnorm - the last computed solution norm

2018:   Level: developer

2020: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm(), SNESGetUpdateNorm()
2021: @*/
2022: PetscErrorCode SNESGetSolutionNorm(SNES snes, PetscReal *xnorm)
2023: {
2027:   *xnorm = snes->xnorm;
2028:   return(0);
2029: }

2031: /*@C
2032:    SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
2033:    of the SNES method.

2035:    Logically Collective on SNES

2037:    Input Parameters:
2038: +  snes - the SNES context
2039: -  normschedule - the frequency of norm computation

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

2050:    Level: developer

2052: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2053: @*/
2054: PetscErrorCode  SNESSetFunctionType(SNES snes, SNESFunctionType type)
2055: {
2058:   snes->functype = type;
2059:   return(0);
2060: }


2063: /*@C
2064:    SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
2065:    of the SNES method.

2067:    Logically Collective on SNES

2069:    Input Parameters:
2070: +  snes - the SNES context
2071: -  normschedule - the type of the norm used

2073:    Level: advanced

2075: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2076: @*/
2077: PetscErrorCode  SNESGetFunctionType(SNES snes, SNESFunctionType *type)
2078: {
2081:   *type = snes->functype;
2082:   return(0);
2083: }

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

2088:      Synopsis:
2089:      #include <petscsnes.h>
2090: $    SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);

2092:      Collective on snes

2094:      Input Parameters:
2095: +  X   - solution vector
2096: .  B   - RHS vector
2097: -  ctx - optional user-defined Gauss-Seidel context

2099:      Output Parameter:
2100: .  X   - solution vector

2102:    Level: intermediate

2104: .seealso:   SNESSetNGS(), SNESGetNGS()
2105: M*/

2107: /*@C
2108:    SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
2109:    use with composed nonlinear solvers.

2111:    Input Parameters:
2112: +  snes   - the SNES context
2113: .  f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
2114: -  ctx    - [optional] user-defined context for private data for the
2115:             smoother evaluation routine (may be NULL)

2117:    Notes:
2118:    The NGS routines are used by the composed nonlinear solver to generate
2119:     a problem appropriate update to the solution, particularly FAS.

2121:    Level: intermediate

2123: .seealso: SNESGetFunction(), SNESComputeNGS()
2124: @*/
2125: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
2126: {
2128:   DM             dm;

2132:   SNESGetDM(snes,&dm);
2133:   DMSNESSetNGS(dm,f,ctx);
2134:   return(0);
2135: }

2137: PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
2138: {
2140:   DM             dm;
2141:   DMSNES         sdm;

2144:   SNESGetDM(snes,&dm);
2145:   DMGetDMSNES(dm,&sdm);
2146:   if (!sdm->ops->computepfunction) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");
2147:   if (!sdm->ops->computepjacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard Jacobian.");
2148:   /*  A(x)*x - b(x) */
2149:   PetscStackPush("SNES Picard user function");
2150:   (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
2151:   PetscStackPop;
2152:   PetscStackPush("SNES Picard user Jacobian");
2153:   (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
2154:   PetscStackPop;
2155:   VecScale(f,-1.0);
2156:   MatMultAdd(snes->jacobian,x,f,f);
2157:   return(0);
2158: }

2160: PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
2161: {
2163:   /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
2164:   return(0);
2165: }

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

2170:    Logically Collective on SNES

2172:    Input Parameters:
2173: +  snes - the SNES context
2174: .  r - vector to store function value
2175: .  b - function evaluation routine
2176: .  Amat - matrix with which A(x) x - b(x) is to be computed
2177: .  Pmat - matrix from which preconditioner is computed (usually the same as Amat)
2178: .  J  - function to compute matrix value, see SNESJacobianFunction for details on its calling sequence
2179: -  ctx - [optional] user-defined context for private data for the
2180:          function evaluation routine (may be NULL)

2182:    Notes:
2183:     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
2184:     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.

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

2188: $     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}
2189: $     Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration.

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

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

2196:    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
2197:    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
2198:    different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).

2200:    Level: intermediate

2202: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
2203: @*/
2204: 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)
2205: {
2207:   DM             dm;

2211:   SNESGetDM(snes, &dm);
2212:   DMSNESSetPicard(dm,b,J,ctx);
2213:   SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
2214:   SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
2215:   return(0);
2216: }

2218: /*@C
2219:    SNESGetPicard - Returns the context for the Picard iteration

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

2223:    Input Parameter:
2224: .  snes - the SNES context

2226:    Output Parameter:
2227: +  r - the function (or NULL)
2228: .  f - the function (or NULL); see SNESFunction for calling sequence details
2229: .  Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
2230: .  Pmat  - the matrix from which the preconditioner will be constructed (or NULL)
2231: .  J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
2232: -  ctx - the function context (or NULL)

2234:    Level: advanced

2236: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
2237: @*/
2238: 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)
2239: {
2241:   DM             dm;

2245:   SNESGetFunction(snes,r,NULL,NULL);
2246:   SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
2247:   SNESGetDM(snes,&dm);
2248:   DMSNESGetPicard(dm,f,J,ctx);
2249:   return(0);
2250: }

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

2255:    Logically Collective on SNES

2257:    Input Parameters:
2258: +  snes - the SNES context
2259: .  func - function evaluation routine
2260: -  ctx - [optional] user-defined context for private data for the
2261:          function evaluation routine (may be NULL)

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

2266: .  f - function vector
2267: -  ctx - optional user-defined function context

2269:    Level: intermediate

2271: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
2272: @*/
2273: PetscErrorCode  SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
2274: {
2277:   if (func) snes->ops->computeinitialguess = func;
2278:   if (ctx)  snes->initialguessP            = ctx;
2279:   return(0);
2280: }

2282: /* --------------------------------------------------------------- */
2283: /*@C
2284:    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
2285:    it assumes a zero right hand side.

2287:    Logically Collective on SNES

2289:    Input Parameter:
2290: .  snes - the SNES context

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

2295:    Level: intermediate

2297: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
2298: @*/
2299: PetscErrorCode  SNESGetRhs(SNES snes,Vec *rhs)
2300: {
2304:   *rhs = snes->vec_rhs;
2305:   return(0);
2306: }

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

2311:    Collective on SNES

2313:    Input Parameters:
2314: +  snes - the SNES context
2315: -  x - input vector

2317:    Output Parameter:
2318: .  y - function vector, as set by SNESSetFunction()

2320:    Notes:
2321:    SNESComputeFunction() is typically used within nonlinear solvers
2322:    implementations, so most users would not generally call this routine
2323:    themselves.

2325:    Level: developer

2327: .seealso: SNESSetFunction(), SNESGetFunction()
2328: @*/
2329: PetscErrorCode  SNESComputeFunction(SNES snes,Vec x,Vec y)
2330: {
2332:   DM             dm;
2333:   DMSNES         sdm;

2341:   VecValidValues(x,2,PETSC_TRUE);

2343:   SNESGetDM(snes,&dm);
2344:   DMGetDMSNES(dm,&sdm);
2345:   if (sdm->ops->computefunction) {
2346:     if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2347:       PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2348:     }
2349:     VecLockReadPush(x);
2350:     PetscStackPush("SNES user function");
2351:     /* ensure domainerror is false prior to computefunction evaluation (may not have been reset) */
2352:     snes->domainerror = PETSC_FALSE;
2353:     (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2354:     PetscStackPop;
2355:     VecLockReadPop(x);
2356:     if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2357:       PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2358:     }
2359:   } else if (snes->vec_rhs) {
2360:     MatMult(snes->jacobian, x, y);
2361:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2362:   if (snes->vec_rhs) {
2363:     VecAXPY(y,-1.0,snes->vec_rhs);
2364:   }
2365:   snes->nfuncs++;
2366:   /*
2367:      domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2368:      propagate the value to all processes
2369:   */
2370:   if (snes->domainerror) {
2371:     VecSetInf(y);
2372:   }
2373:   return(0);
2374: }

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

2379:    Collective on SNES

2381:    Input Parameters:
2382: +  snes - the SNES context
2383: .  x - input vector
2384: -  b - rhs vector

2386:    Output Parameter:
2387: .  x - new solution vector

2389:    Notes:
2390:    SNESComputeNGS() is typically used within composed nonlinear solver
2391:    implementations, so most users would not generally call this routine
2392:    themselves.

2394:    Level: developer

2396: .seealso: SNESSetNGS(), SNESComputeFunction()
2397: @*/
2398: PetscErrorCode  SNESComputeNGS(SNES snes,Vec b,Vec x)
2399: {
2401:   DM             dm;
2402:   DMSNES         sdm;

2410:   if (b) {VecValidValues(b,2,PETSC_TRUE);}
2411:   PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2412:   SNESGetDM(snes,&dm);
2413:   DMGetDMSNES(dm,&sdm);
2414:   if (sdm->ops->computegs) {
2415:     if (b) {VecLockReadPush(b);}
2416:     PetscStackPush("SNES user NGS");
2417:     (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2418:     PetscStackPop;
2419:     if (b) {VecLockReadPop(b);}
2420:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2421:   PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2422:   return(0);
2423: }

2425: PetscErrorCode SNESTestJacobian(SNES snes)
2426: {
2427:   Mat               A,B,C,D,jacobian;
2428:   Vec               x = snes->vec_sol,f = snes->vec_func;
2429:   PetscErrorCode    ierr;
2430:   PetscReal         nrm,gnorm;
2431:   PetscReal         threshold = 1.e-5;
2432:   MatType           mattype;
2433:   PetscInt          m,n,M,N;
2434:   void              *functx;
2435:   PetscBool         complete_print = PETSC_FALSE,threshold_print = PETSC_FALSE,test = PETSC_FALSE,flg;
2436:   PetscViewer       viewer,mviewer;
2437:   MPI_Comm          comm;
2438:   PetscInt          tabs;
2439:   static PetscBool  directionsprinted = PETSC_FALSE;
2440:   PetscViewerFormat format;

2443:   PetscObjectOptionsBegin((PetscObject)snes);
2444:   PetscOptionsName("-snes_test_jacobian","Compare hand-coded and finite difference Jacobians","None",&test);
2445:   PetscOptionsReal("-snes_test_jacobian", "Threshold for element difference between hand-coded and finite difference being meaningful", "None", threshold, &threshold,NULL);
2446:   PetscOptionsViewer("-snes_test_jacobian_view","View difference between hand-coded and finite difference Jacobians element entries","None",&mviewer,&format,&complete_print);
2447:   if (!complete_print) {
2448:     PetscOptionsViewer("-snes_test_jacobian_display","Display difference between hand-coded and finite difference Jacobians","None",&mviewer,&format,&complete_print);
2449:   }
2450:   /* for compatibility with PETSc 3.9 and older. */
2451:   PetscOptionsReal("-snes_test_jacobian_display_threshold", "Display difference between hand-coded and finite difference Jacobians which exceed input threshold", "None", threshold, &threshold, &threshold_print);
2452:   PetscOptionsEnd();
2453:   if (!test) return(0);

2455:   PetscObjectGetComm((PetscObject)snes,&comm);
2456:   PetscViewerASCIIGetStdout(comm,&viewer);
2457:   PetscViewerASCIIGetTab(viewer, &tabs);
2458:   PetscViewerASCIISetTab(viewer, ((PetscObject)snes)->tablevel);
2459:   PetscViewerASCIIPrintf(viewer,"  ---------- Testing Jacobian -------------\n");
2460:   if (!complete_print && !directionsprinted) {
2461:     PetscViewerASCIIPrintf(viewer,"  Run with -snes_test_jacobian_view and optionally -snes_test_jacobian <threshold> to show difference\n");
2462:     PetscViewerASCIIPrintf(viewer,"    of hand-coded and finite difference Jacobian entries greater than <threshold>.\n");
2463:   }
2464:   if (!directionsprinted) {
2465:     PetscViewerASCIIPrintf(viewer,"  Testing hand-coded Jacobian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");
2466:     PetscViewerASCIIPrintf(viewer,"    O(1.e-8), the hand-coded Jacobian is probably correct.\n");
2467:     directionsprinted = PETSC_TRUE;
2468:   }
2469:   if (complete_print) {
2470:     PetscViewerPushFormat(mviewer,format);
2471:   }

2473:   PetscObjectTypeCompare((PetscObject)snes->jacobian,MATMFFD,&flg);
2474:   if (!flg) jacobian = snes->jacobian;
2475:   else jacobian = snes->jacobian_pre;

2477:   if (!x) {
2478:     MatCreateVecs(jacobian, &x, NULL);
2479:   } else {
2480:     PetscObjectReference((PetscObject) x);
2481:   }
2482:   if (!f) {
2483:     VecDuplicate(x, &f);
2484:   } else {
2485:     PetscObjectReference((PetscObject) f);
2486:   }
2487:   /* evaluate the function at this point because SNESComputeJacobianDefault() assumes that the function has been evaluated and put into snes->vec_func */
2488:   SNESComputeFunction(snes,x,f);
2489:   VecDestroy(&f);

2491:   while (jacobian) {
2492:     PetscObjectBaseTypeCompareAny((PetscObject)jacobian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPISBAIJ,"");
2493:     if (flg) {
2494:       A    = jacobian;
2495:       PetscObjectReference((PetscObject)A);
2496:     } else {
2497:       MatComputeOperator(jacobian,MATAIJ,&A);
2498:     }

2500:     MatGetType(A,&mattype);
2501:     MatGetSize(A,&M,&N);
2502:     MatGetLocalSize(A,&m,&n);

2504:     MatCreate(PetscObjectComm((PetscObject)A),&B);
2505:     MatSetType(B,mattype);
2506:     MatSetSizes(B,m,n,M,N);
2507:     MatSetBlockSizesFromMats(B,A,A);
2508:     MatSetUp(B);
2509:     MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);

2511:     SNESGetFunction(snes,NULL,NULL,&functx);
2512:     SNESComputeJacobianDefault(snes,x,B,B,functx);

2514:     MatDuplicate(B,MAT_COPY_VALUES,&D);
2515:     MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2516:     MatNorm(D,NORM_FROBENIUS,&nrm);
2517:     MatNorm(A,NORM_FROBENIUS,&gnorm);
2518:     MatDestroy(&D);
2519:     if (!gnorm) gnorm = 1; /* just in case */
2520:     PetscViewerASCIIPrintf(viewer,"  ||J - Jfd||_F/||J||_F = %g, ||J - Jfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);

2522:     if (complete_print) {
2523:       PetscViewerASCIIPrintf(viewer,"  Hand-coded Jacobian ----------\n");
2524:       MatView(A,mviewer);
2525:       PetscViewerASCIIPrintf(viewer,"  Finite difference Jacobian ----------\n");
2526:       MatView(B,mviewer);
2527:     }

2529:     if (threshold_print || complete_print) {
2530:       PetscInt          Istart, Iend, *ccols, bncols, cncols, j, row;
2531:       PetscScalar       *cvals;
2532:       const PetscInt    *bcols;
2533:       const PetscScalar *bvals;

2535:       MatCreate(PetscObjectComm((PetscObject)A),&C);
2536:       MatSetType(C,mattype);
2537:       MatSetSizes(C,m,n,M,N);
2538:       MatSetBlockSizesFromMats(C,A,A);
2539:       MatSetUp(C);
2540:       MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);

2542:       MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2543:       MatGetOwnershipRange(B,&Istart,&Iend);

2545:       for (row = Istart; row < Iend; row++) {
2546:         MatGetRow(B,row,&bncols,&bcols,&bvals);
2547:         PetscMalloc2(bncols,&ccols,bncols,&cvals);
2548:         for (j = 0, cncols = 0; j < bncols; j++) {
2549:           if (PetscAbsScalar(bvals[j]) > threshold) {
2550:             ccols[cncols] = bcols[j];
2551:             cvals[cncols] = bvals[j];
2552:             cncols += 1;
2553:           }
2554:         }
2555:         if (cncols) {
2556:           MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);
2557:         }
2558:         MatRestoreRow(B,row,&bncols,&bcols,&bvals);
2559:         PetscFree2(ccols,cvals);
2560:       }
2561:       MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2562:       MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2563:       PetscViewerASCIIPrintf(viewer,"  Hand-coded minus finite-difference Jacobian with tolerance %g ----------\n",(double)threshold);
2564:       MatView(C,complete_print ? mviewer : viewer);
2565:       MatDestroy(&C);
2566:     }
2567:     MatDestroy(&A);
2568:     MatDestroy(&B);

2570:     if (jacobian != snes->jacobian_pre) {
2571:       jacobian = snes->jacobian_pre;
2572:       PetscViewerASCIIPrintf(viewer,"  ---------- Testing Jacobian for preconditioner -------------\n");
2573:     }
2574:     else jacobian = NULL;
2575:   }
2576:   VecDestroy(&x);
2577:   if (complete_print) {
2578:     PetscViewerPopFormat(mviewer);
2579:   }
2580:   if (mviewer) { PetscViewerDestroy(&mviewer); }
2581:   PetscViewerASCIISetTab(viewer,tabs);
2582:   return(0);
2583: }

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

2588:    Collective on SNES

2590:    Input Parameters:
2591: +  snes - the SNES context
2592: -  x - input vector

2594:    Output Parameters:
2595: +  A - Jacobian matrix
2596: -  B - optional preconditioning matrix

2598:   Options Database Keys:
2599: +    -snes_lag_preconditioner <lag>
2600: .    -snes_lag_jacobian <lag>
2601: .    -snes_test_jacobian - compare the user provided Jacobian with one compute via finite differences to check for errors
2602: .    -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
2603: .    -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
2604: .    -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2605: .    -snes_compare_explicit_draw  - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2606: .    -snes_compare_explicit_contour  - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2607: .    -snes_compare_operator  - Make the comparison options above use the operator instead of the preconditioning matrix
2608: .    -snes_compare_coloring - Compute the finite difference Jacobian using coloring and display norms of difference
2609: .    -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
2610: .    -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2611: .    -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2612: .    -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2613: .    -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
2614: -    -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences


2617:    Notes:
2618:    Most users should not need to explicitly call this routine, as it
2619:    is used internally within the nonlinear solvers.

2621:    Developer Notes:
2622:     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
2623:       for with the SNESType of test that has been removed.

2625:    Level: developer

2627: .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2628: @*/
2629: PetscErrorCode  SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2630: {
2632:   PetscBool      flag;
2633:   DM             dm;
2634:   DMSNES         sdm;
2635:   KSP            ksp;

2641:   VecValidValues(X,2,PETSC_TRUE);
2642:   SNESGetDM(snes,&dm);
2643:   DMGetDMSNES(dm,&sdm);

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

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

2649:   if (snes->lagjacobian == -2) {
2650:     snes->lagjacobian = -1;

2652:     PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2653:   } else if (snes->lagjacobian == -1) {
2654:     PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2655:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2656:     if (flag) {
2657:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2658:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2659:     }
2660:     return(0);
2661:   } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2662:     PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2663:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2664:     if (flag) {
2665:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2666:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2667:     }
2668:     return(0);
2669:   }
2670:   if (snes->npc && snes->npcside== PC_LEFT) {
2671:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2672:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2673:     return(0);
2674:   }

2676:   PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2677:   VecLockReadPush(X);
2678:   PetscStackPush("SNES user Jacobian function");
2679:   (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2680:   PetscStackPop;
2681:   VecLockReadPop(X);
2682:   PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);

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

2687:   /* the next line ensures that snes->ksp exists */
2688:   SNESGetKSP(snes,&ksp);
2689:   if (snes->lagpreconditioner == -2) {
2690:     PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2691:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2692:     snes->lagpreconditioner = -1;
2693:   } else if (snes->lagpreconditioner == -1) {
2694:     PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2695:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2696:   } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2697:     PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2698:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2699:   } else {
2700:     PetscInfo(snes,"Rebuilding preconditioner\n");
2701:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2702:   }

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

2780:       MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2781:       MatColoringCreate(Bfd,&coloring);
2782:       MatColoringSetType(coloring,MATCOLORINGSL);
2783:       MatColoringSetFromOptions(coloring);
2784:       MatColoringApply(coloring,&iscoloring);
2785:       MatColoringDestroy(&coloring);
2786:       MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2787:       MatFDColoringSetFromOptions(matfdcoloring);
2788:       MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2789:       ISColoringDestroy(&iscoloring);

2791:       /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2792:       SNESGetFunction(snes,NULL,&func,&funcctx);
2793:       MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2794:       PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2795:       PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2796:       MatFDColoringSetFromOptions(matfdcoloring);
2797:       MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2798:       MatFDColoringDestroy(&matfdcoloring);

2800:       PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2801:       if (flag_draw || flag_contour) {
2802:         PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2803:         if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2804:       } else vdraw = NULL;
2805:       PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2806:       if (flag_display) {MatView(B,vstdout);}
2807:       if (vdraw) {MatView(B,vdraw);}
2808:       PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2809:       if (flag_display) {MatView(Bfd,vstdout);}
2810:       if (vdraw) {MatView(Bfd,vdraw);}
2811:       MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2812:       MatNorm(Bfd,NORM_1,&norm1);
2813:       MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2814:       MatNorm(Bfd,NORM_MAX,&normmax);
2815:       PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2816:       if (flag_display) {MatView(Bfd,vstdout);}
2817:       if (vdraw) {              /* Always use contour for the difference */
2818:         PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2819:         MatView(Bfd,vdraw);
2820:         PetscViewerPopFormat(vdraw);
2821:       }
2822:       if (flag_contour) {PetscViewerPopFormat(vdraw);}

2824:       if (flag_threshold) {
2825:         PetscInt bs,rstart,rend,i;
2826:         MatGetBlockSize(B,&bs);
2827:         MatGetOwnershipRange(B,&rstart,&rend);
2828:         for (i=rstart; i<rend; i++) {
2829:           const PetscScalar *ba,*ca;
2830:           const PetscInt    *bj,*cj;
2831:           PetscInt          bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
2832:           PetscReal         maxentry = 0,maxdiff = 0,maxrdiff = 0;
2833:           MatGetRow(B,i,&bn,&bj,&ba);
2834:           MatGetRow(Bfd,i,&cn,&cj,&ca);
2835:           if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
2836:           for (j=0; j<bn; j++) {
2837:             PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2838:             if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
2839:               maxentrycol = bj[j];
2840:               maxentry    = PetscRealPart(ba[j]);
2841:             }
2842:             if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
2843:               maxdiffcol = bj[j];
2844:               maxdiff    = PetscRealPart(ca[j]);
2845:             }
2846:             if (rdiff > maxrdiff) {
2847:               maxrdiffcol = bj[j];
2848:               maxrdiff    = rdiff;
2849:             }
2850:           }
2851:           if (maxrdiff > 1) {
2852:             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);
2853:             for (j=0; j<bn; j++) {
2854:               PetscReal rdiff;
2855:               rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2856:               if (rdiff > 1) {
2857:                 PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
2858:               }
2859:             }
2860:             PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
2861:           }
2862:           MatRestoreRow(B,i,&bn,&bj,&ba);
2863:           MatRestoreRow(Bfd,i,&cn,&cj,&ca);
2864:         }
2865:       }
2866:       PetscViewerDestroy(&vdraw);
2867:       MatDestroy(&Bfd);
2868:     }
2869:   }
2870:   return(0);
2871: }

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

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

2880:      Collective on snes

2882:     Input Parameters:
2883: +  x - input vector, the Jacobian is to be computed at this value
2884: -  ctx - [optional] user-defined Jacobian context

2886:     Output Parameters:
2887: +  Amat - the matrix that defines the (approximate) Jacobian
2888: -  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.

2890:    Level: intermediate

2892: .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2893: M*/

2895: /*@C
2896:    SNESSetJacobian - Sets the function to compute Jacobian as well as the
2897:    location to store the matrix.

2899:    Logically Collective on SNES

2901:    Input Parameters:
2902: +  snes - the SNES context
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: .  J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
2906: -  ctx - [optional] user-defined context for private data for the
2907:          Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)

2909:    Notes:
2910:    If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2911:    each matrix.

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

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

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

2922:    Level: beginner

2924: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J,
2925:           SNESSetPicard(), SNESJacobianFunction
2926: @*/
2927: PetscErrorCode  SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2928: {
2930:   DM             dm;

2938:   SNESGetDM(snes,&dm);
2939:   DMSNESSetJacobian(dm,J,ctx);
2940:   if (Amat) {
2941:     PetscObjectReference((PetscObject)Amat);
2942:     MatDestroy(&snes->jacobian);

2944:     snes->jacobian = Amat;
2945:   }
2946:   if (Pmat) {
2947:     PetscObjectReference((PetscObject)Pmat);
2948:     MatDestroy(&snes->jacobian_pre);

2950:     snes->jacobian_pre = Pmat;
2951:   }
2952:   return(0);
2953: }

2955: /*@C
2956:    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2957:    provided context for evaluating the Jacobian.

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

2961:    Input Parameter:
2962: .  snes - the nonlinear solver context

2964:    Output Parameters:
2965: +  Amat - location to stash (approximate) Jacobian matrix (or NULL)
2966: .  Pmat - location to stash matrix used to compute the preconditioner (or NULL)
2967: .  J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
2968: -  ctx - location to stash Jacobian ctx (or NULL)

2970:    Level: advanced

2972: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
2973: @*/
2974: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2975: {
2977:   DM             dm;
2978:   DMSNES         sdm;

2982:   if (Amat) *Amat = snes->jacobian;
2983:   if (Pmat) *Pmat = snes->jacobian_pre;
2984:   SNESGetDM(snes,&dm);
2985:   DMGetDMSNES(dm,&sdm);
2986:   if (J) *J = sdm->ops->computejacobian;
2987:   if (ctx) *ctx = sdm->jacobianctx;
2988:   return(0);
2989: }

2991: /*@
2992:    SNESSetUp - Sets up the internal data structures for the later use
2993:    of a nonlinear solver.

2995:    Collective on SNES

2997:    Input Parameters:
2998: .  snes - the SNES context

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

3007:    Level: advanced

3009: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
3010: @*/
3011: PetscErrorCode  SNESSetUp(SNES snes)
3012: {
3014:   DM             dm;
3015:   DMSNES         sdm;
3016:   SNESLineSearch linesearch, pclinesearch;
3017:   void           *lsprectx,*lspostctx;
3018:   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
3019:   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
3020:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
3021:   Vec            f,fpc;
3022:   void           *funcctx;
3023:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
3024:   void           *jacctx,*appctx;
3025:   Mat            j,jpre;

3029:   if (snes->setupcalled) return(0);
3030:   PetscLogEventBegin(SNES_Setup,snes,0,0,0);

3032:   if (!((PetscObject)snes)->type_name) {
3033:     SNESSetType(snes,SNESNEWTONLS);
3034:   }

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

3038:   SNESGetDM(snes,&dm);
3039:   DMGetDMSNES(dm,&sdm);
3040:   if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
3041:   if (!sdm->ops->computejacobian) {
3042:     DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
3043:   }
3044:   if (!snes->vec_func) {
3045:     DMCreateGlobalVector(dm,&snes->vec_func);
3046:   }

3048:   if (!snes->ksp) {
3049:     SNESGetKSP(snes, &snes->ksp);
3050:   }

3052:   if (snes->linesearch) {
3053:     SNESGetLineSearch(snes, &snes->linesearch);
3054:     SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);
3055:   }

3057:   if (snes->npc && (snes->npcside== PC_LEFT)) {
3058:     snes->mf          = PETSC_TRUE;
3059:     snes->mf_operator = PETSC_FALSE;
3060:   }

3062:   if (snes->npc) {
3063:     /* copy the DM over */
3064:     SNESGetDM(snes,&dm);
3065:     SNESSetDM(snes->npc,dm);

3067:     SNESGetFunction(snes,&f,&func,&funcctx);
3068:     VecDuplicate(f,&fpc);
3069:     SNESSetFunction(snes->npc,fpc,func,funcctx);
3070:     SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
3071:     SNESSetJacobian(snes->npc,j,jpre,jac,jacctx);
3072:     SNESGetApplicationContext(snes,&appctx);
3073:     SNESSetApplicationContext(snes->npc,appctx);
3074:     VecDestroy(&fpc);

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

3079:     /* default to 1 iteration */
3080:     SNESSetTolerances(snes->npc,0.0,0.0,0.0,1,snes->npc->max_funcs);
3081:     if (snes->npcside==PC_RIGHT) {
3082:       SNESSetNormSchedule(snes->npc,SNES_NORM_FINAL_ONLY);
3083:     } else {
3084:       SNESSetNormSchedule(snes->npc,SNES_NORM_NONE);
3085:     }
3086:     SNESSetFromOptions(snes->npc);

3088:     /* copy the line search context over */
3089:     if (snes->linesearch && snes->npc->linesearch) {
3090:       SNESGetLineSearch(snes,&linesearch);
3091:       SNESGetLineSearch(snes->npc,&pclinesearch);
3092:       SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
3093:       SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
3094:       SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
3095:       SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
3096:       PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
3097:     }
3098:   }
3099:   if (snes->mf) {
3100:     SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
3101:   }
3102:   if (snes->ops->usercompute && !snes->user) {
3103:     (*snes->ops->usercompute)(snes,(void**)&snes->user);
3104:   }

3106:   snes->jac_iter = 0;
3107:   snes->pre_iter = 0;

3109:   if (snes->ops->setup) {
3110:     (*snes->ops->setup)(snes);
3111:   }

3113:   if (snes->npc && (snes->npcside== PC_LEFT)) {
3114:     if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3115:       if (snes->linesearch){
3116:         SNESGetLineSearch(snes,&linesearch);
3117:         SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
3118:       }
3119:     }
3120:   }
3121:   PetscLogEventEnd(SNES_Setup,snes,0,0,0);
3122:   snes->setupcalled = PETSC_TRUE;
3123:   return(0);
3124: }

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

3129:    Collective on SNES

3131:    Input Parameter:
3132: .  snes - iterative context obtained from SNESCreate()

3134:    Level: intermediate

3136:    Notes:
3137:     Also calls the application context destroy routine set with SNESSetComputeApplicationContext()

3139: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
3140: @*/
3141: PetscErrorCode  SNESReset(SNES snes)
3142: {

3147:   if (snes->ops->userdestroy && snes->user) {
3148:     (*snes->ops->userdestroy)((void**)&snes->user);
3149:     snes->user = NULL;
3150:   }
3151:   if (snes->npc) {
3152:     SNESReset(snes->npc);
3153:   }

3155:   if (snes->ops->reset) {
3156:     (*snes->ops->reset)(snes);
3157:   }
3158:   if (snes->ksp) {
3159:     KSPReset(snes->ksp);
3160:   }

3162:   if (snes->linesearch) {
3163:     SNESLineSearchReset(snes->linesearch);
3164:   }

3166:   VecDestroy(&snes->vec_rhs);
3167:   VecDestroy(&snes->vec_sol);
3168:   VecDestroy(&snes->vec_sol_update);
3169:   VecDestroy(&snes->vec_func);
3170:   MatDestroy(&snes->jacobian);
3171:   MatDestroy(&snes->jacobian_pre);
3172:   VecDestroyVecs(snes->nwork,&snes->work);
3173:   VecDestroyVecs(snes->nvwork,&snes->vwork);

3175:   snes->alwayscomputesfinalresidual = PETSC_FALSE;

3177:   snes->nwork       = snes->nvwork = 0;
3178:   snes->setupcalled = PETSC_FALSE;
3179:   return(0);
3180: }

3182: /*@
3183:    SNESDestroy - Destroys the nonlinear solver context that was created
3184:    with SNESCreate().

3186:    Collective on SNES

3188:    Input Parameter:
3189: .  snes - the SNES context

3191:    Level: beginner

3193: .seealso: SNESCreate(), SNESSolve()
3194: @*/
3195: PetscErrorCode  SNESDestroy(SNES *snes)
3196: {

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

3204:   SNESReset((*snes));
3205:   SNESDestroy(&(*snes)->npc);

3207:   /* if memory was published with SAWs then destroy it */
3208:   PetscObjectSAWsViewOff((PetscObject)*snes);
3209:   if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}

3211:   if ((*snes)->dm) {DMCoarsenHookRemove((*snes)->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,*snes);}
3212:   DMDestroy(&(*snes)->dm);
3213:   KSPDestroy(&(*snes)->ksp);
3214:   SNESLineSearchDestroy(&(*snes)->linesearch);

3216:   PetscFree((*snes)->kspconvctx);
3217:   if ((*snes)->ops->convergeddestroy) {
3218:     (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
3219:   }
3220:   if ((*snes)->conv_hist_alloc) {
3221:     PetscFree2((*snes)->conv_hist,(*snes)->conv_hist_its);
3222:   }
3223:   SNESMonitorCancel((*snes));
3224:   PetscHeaderDestroy(snes);
3225:   return(0);
3226: }

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

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

3233:    Logically Collective on SNES

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

3240:    Options Database Keys:
3241: .    -snes_lag_preconditioner <lag>

3243:    Notes:
3244:    The default is 1
3245:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3246:    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use

3248:    Level: intermediate

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

3252: @*/
3253: PetscErrorCode  SNESSetLagPreconditioner(SNES snes,PetscInt lag)
3254: {
3257:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3258:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3260:   snes->lagpreconditioner = lag;
3261:   return(0);
3262: }

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

3267:    Logically Collective on SNES

3269:    Input Parameters:
3270: +  snes - the SNES context
3271: -  steps - the number of refinements to do, defaults to 0

3273:    Options Database Keys:
3274: .    -snes_grid_sequence <steps>

3276:    Level: intermediate

3278:    Notes:
3279:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

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

3283: @*/
3284: PetscErrorCode  SNESSetGridSequence(SNES snes,PetscInt steps)
3285: {
3289:   snes->gridsequence = steps;
3290:   return(0);
3291: }

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

3296:    Logically Collective on SNES

3298:    Input Parameter:
3299: .  snes - the SNES context

3301:    Output Parameter:
3302: .  steps - the number of refinements to do, defaults to 0

3304:    Options Database Keys:
3305: .    -snes_grid_sequence <steps>

3307:    Level: intermediate

3309:    Notes:
3310:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

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

3314: @*/
3315: PetscErrorCode  SNESGetGridSequence(SNES snes,PetscInt *steps)
3316: {
3319:   *steps = snes->gridsequence;
3320:   return(0);
3321: }

3323: /*@
3324:    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt

3326:    Not Collective

3328:    Input Parameter:
3329: .  snes - the SNES context

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

3335:    Options Database Keys:
3336: .    -snes_lag_preconditioner <lag>

3338:    Notes:
3339:    The default is 1
3340:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

3342:    Level: intermediate

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

3346: @*/
3347: PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
3348: {
3351:   *lag = snes->lagpreconditioner;
3352:   return(0);
3353: }

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

3359:    Logically Collective on SNES

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

3366:    Options Database Keys:
3367: .    -snes_lag_jacobian <lag>

3369:    Notes:
3370:    The default is 1
3371:    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3372:    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
3373:    at the next Newton step but never again (unless it is reset to another value)

3375:    Level: intermediate

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

3379: @*/
3380: PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
3381: {
3384:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3385:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3387:   snes->lagjacobian = lag;
3388:   return(0);
3389: }

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

3394:    Not Collective

3396:    Input Parameter:
3397: .  snes - the SNES context

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

3403:    Options Database Keys:
3404: .    -snes_lag_jacobian <lag>

3406:    Notes:
3407:    The default is 1
3408:    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

3410:    Level: intermediate

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

3414: @*/
3415: PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
3416: {
3419:   *lag = snes->lagjacobian;
3420:   return(0);
3421: }

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

3426:    Logically collective on SNES

3428:    Input Parameter:
3429: +  snes - the SNES context
3430: -   flg - jacobian lagging persists if true

3432:    Options Database Keys:
3433: .    -snes_lag_jacobian_persists <flg>

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

3440:    Level: developer

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

3444: @*/
3445: PetscErrorCode  SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
3446: {
3450:   snes->lagjac_persist = flg;
3451:   return(0);
3452: }

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

3457:    Logically Collective on SNES

3459:    Input Parameter:
3460: +  snes - the SNES context
3461: -   flg - preconditioner lagging persists if true

3463:    Options Database Keys:
3464: .    -snes_lag_jacobian_persists <flg>

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

3471:    Level: developer

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

3475: @*/
3476: PetscErrorCode  SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3477: {
3481:   snes->lagpre_persist = flg;
3482:   return(0);
3483: }

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

3488:    Logically Collective on SNES

3490:    Input Parameters:
3491: +  snes - the SNES context
3492: -  force - PETSC_TRUE require at least one iteration

3494:    Options Database Keys:
3495: .    -snes_force_iteration <force> - Sets forcing an iteration

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

3500:    Level: intermediate

3502: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3503: @*/
3504: PetscErrorCode  SNESSetForceIteration(SNES snes,PetscBool force)
3505: {
3508:   snes->forceiteration = force;
3509:   return(0);
3510: }

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

3515:    Logically Collective on SNES

3517:    Input Parameters:
3518: .  snes - the SNES context

3520:    Output Parameter:
3521: .  force - PETSC_TRUE requires at least one iteration.

3523:    Level: intermediate

3525: .seealso: SNESSetForceIteration(), SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3526: @*/
3527: PetscErrorCode  SNESGetForceIteration(SNES snes,PetscBool *force)
3528: {
3531:   *force = snes->forceiteration;
3532:   return(0);
3533: }

3535: /*@
3536:    SNESSetTolerances - Sets various parameters used in convergence tests.

3538:    Logically Collective on SNES

3540:    Input Parameters:
3541: +  snes - the SNES context
3542: .  abstol - absolute convergence tolerance
3543: .  rtol - relative convergence tolerance
3544: .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
3545: .  maxit - maximum number of iterations
3546: -  maxf - maximum number of function evaluations (-1 indicates no limit)

3548:    Options Database Keys:
3549: +    -snes_atol <abstol> - Sets abstol
3550: .    -snes_rtol <rtol> - Sets rtol
3551: .    -snes_stol <stol> - Sets stol
3552: .    -snes_max_it <maxit> - Sets maxit
3553: -    -snes_max_funcs <maxf> - Sets maxf

3555:    Notes:
3556:    The default maximum number of iterations is 50.
3557:    The default maximum number of function evaluations is 1000.

3559:    Level: intermediate

3561: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance(), SNESSetForceIteration()
3562: @*/
3563: PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3564: {

3573:   if (abstol != PETSC_DEFAULT) {
3574:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3575:     snes->abstol = abstol;
3576:   }
3577:   if (rtol != PETSC_DEFAULT) {
3578:     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);
3579:     snes->rtol = rtol;
3580:   }
3581:   if (stol != PETSC_DEFAULT) {
3582:     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3583:     snes->stol = stol;
3584:   }
3585:   if (maxit != PETSC_DEFAULT) {
3586:     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3587:     snes->max_its = maxit;
3588:   }
3589:   if (maxf != PETSC_DEFAULT) {
3590:     if (maxf < -1) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be -1 or nonnegative",maxf);
3591:     snes->max_funcs = maxf;
3592:   }
3593:   snes->tolerancesset = PETSC_TRUE;
3594:   return(0);
3595: }

3597: /*@
3598:    SNESSetDivergenceTolerance - Sets the divergence tolerance used for the SNES divergence test.

3600:    Logically Collective on SNES

3602:    Input Parameters:
3603: +  snes - the SNES context
3604: -  divtol - the divergence tolerance. Use -1 to deactivate the test.

3606:    Options Database Keys:
3607: .    -snes_divergence_tolerance <divtol> - Sets divtol

3609:    Notes:
3610:    The default divergence tolerance is 1e4.

3612:    Level: intermediate

3614: .seealso: SNESSetTolerances(), SNESGetDivergenceTolerance
3615: @*/
3616: PetscErrorCode  SNESSetDivergenceTolerance(SNES snes,PetscReal divtol)
3617: {

3622:   if (divtol != PETSC_DEFAULT) {
3623:     snes->divtol = divtol;
3624:   }
3625:   else {
3626:     snes->divtol = 1.0e4;
3627:   }
3628:   return(0);
3629: }

3631: /*@
3632:    SNESGetTolerances - Gets various parameters used in convergence tests.

3634:    Not Collective

3636:    Input Parameters:
3637: +  snes - the SNES context
3638: .  atol - absolute convergence tolerance
3639: .  rtol - relative convergence tolerance
3640: .  stol -  convergence tolerance in terms of the norm
3641:            of the change in the solution between steps
3642: .  maxit - maximum number of iterations
3643: -  maxf - maximum number of function evaluations

3645:    Notes:
3646:    The user can specify NULL for any parameter that is not needed.

3648:    Level: intermediate

3650: .seealso: SNESSetTolerances()
3651: @*/
3652: PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3653: {
3656:   if (atol)  *atol  = snes->abstol;
3657:   if (rtol)  *rtol  = snes->rtol;
3658:   if (stol)  *stol  = snes->stol;
3659:   if (maxit) *maxit = snes->max_its;
3660:   if (maxf)  *maxf  = snes->max_funcs;
3661:   return(0);
3662: }

3664: /*@
3665:    SNESGetDivergenceTolerance - Gets divergence tolerance used in divergence test.

3667:    Not Collective

3669:    Input Parameters:
3670: +  snes - the SNES context
3671: -  divtol - divergence tolerance

3673:    Level: intermediate

3675: .seealso: SNESSetDivergenceTolerance()
3676: @*/
3677: PetscErrorCode  SNESGetDivergenceTolerance(SNES snes,PetscReal *divtol)
3678: {
3681:   if (divtol) *divtol = snes->divtol;
3682:   return(0);
3683: }

3685: /*@
3686:    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.

3688:    Logically Collective on SNES

3690:    Input Parameters:
3691: +  snes - the SNES context
3692: -  tol - tolerance

3694:    Options Database Key:
3695: .  -snes_trtol <tol> - Sets tol

3697:    Level: intermediate

3699: .seealso: SNESSetTolerances()
3700: @*/
3701: PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3702: {
3706:   snes->deltatol = tol;
3707:   return(0);
3708: }

3710: /*
3711:    Duplicate the lg monitors for SNES from KSP; for some reason with
3712:    dynamic libraries things don't work under Sun4 if we just use
3713:    macros instead of functions
3714: */
3715: PetscErrorCode  SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,void *ctx)
3716: {

3721:   KSPMonitorLGResidualNorm((KSP)snes,it,norm,ctx);
3722:   return(0);
3723: }

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

3730:   KSPMonitorLGResidualNormCreate(comm,host,label,x,y,m,n,lgctx);
3731:   return(0);
3732: }

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

3736: PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3737: {
3738:   PetscDrawLG      lg;
3739:   PetscErrorCode   ierr;
3740:   PetscReal        x,y,per;
3741:   PetscViewer      v = (PetscViewer)monctx;
3742:   static PetscReal prev; /* should be in the context */
3743:   PetscDraw        draw;

3747:   PetscViewerDrawGetDrawLG(v,0,&lg);
3748:   if (!n) {PetscDrawLGReset(lg);}
3749:   PetscDrawLGGetDraw(lg,&draw);
3750:   PetscDrawSetTitle(draw,"Residual norm");
3751:   x    = (PetscReal)n;
3752:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3753:   else y = -15.0;
3754:   PetscDrawLGAddPoint(lg,&x,&y);
3755:   if (n < 20 || !(n % 5) || snes->reason) {
3756:     PetscDrawLGDraw(lg);
3757:     PetscDrawLGSave(lg);
3758:   }

3760:   PetscViewerDrawGetDrawLG(v,1,&lg);
3761:   if (!n) {PetscDrawLGReset(lg);}
3762:   PetscDrawLGGetDraw(lg,&draw);
3763:   PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3764:    SNESMonitorRange_Private(snes,n,&per);
3765:   x    = (PetscReal)n;
3766:   y    = 100.0*per;
3767:   PetscDrawLGAddPoint(lg,&x,&y);
3768:   if (n < 20 || !(n % 5) || snes->reason) {
3769:     PetscDrawLGDraw(lg);
3770:     PetscDrawLGSave(lg);
3771:   }

3773:   PetscViewerDrawGetDrawLG(v,2,&lg);
3774:   if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3775:   PetscDrawLGGetDraw(lg,&draw);
3776:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3777:   x    = (PetscReal)n;
3778:   y    = (prev - rnorm)/prev;
3779:   PetscDrawLGAddPoint(lg,&x,&y);
3780:   if (n < 20 || !(n % 5) || snes->reason) {
3781:     PetscDrawLGDraw(lg);
3782:     PetscDrawLGSave(lg);
3783:   }

3785:   PetscViewerDrawGetDrawLG(v,3,&lg);
3786:   if (!n) {PetscDrawLGReset(lg);}
3787:   PetscDrawLGGetDraw(lg,&draw);
3788:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3789:   x    = (PetscReal)n;
3790:   y    = (prev - rnorm)/(prev*per);
3791:   if (n > 2) { /*skip initial crazy value */
3792:     PetscDrawLGAddPoint(lg,&x,&y);
3793:   }
3794:   if (n < 20 || !(n % 5) || snes->reason) {
3795:     PetscDrawLGDraw(lg);
3796:     PetscDrawLGSave(lg);
3797:   }
3798:   prev = rnorm;
3799:   return(0);
3800: }

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

3805:    Collective on SNES

3807:    Input Parameters:
3808: +  snes - nonlinear solver context obtained from SNESCreate()
3809: .  iter - iteration number
3810: -  rnorm - relative norm of the residual

3812:    Notes:
3813:    This routine is called by the SNES implementations.
3814:    It does not typically need to be called by the user.

3816:    Level: developer

3818: .seealso: SNESMonitorSet()
3819: @*/
3820: PetscErrorCode  SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3821: {
3823:   PetscInt       i,n = snes->numbermonitors;

3826:   VecLockReadPush(snes->vec_sol);
3827:   for (i=0; i<n; i++) {
3828:     (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3829:   }
3830:   VecLockReadPop(snes->vec_sol);
3831:   return(0);
3832: }

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

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

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

3843:      Collective on snes

3845:     Input Parameters:
3846: +    snes - the SNES context
3847: .    its - iteration number
3848: .    norm - 2-norm function value (may be estimated)
3849: -    mctx - [optional] monitoring context

3851:    Level: advanced

3853: .seealso:   SNESMonitorSet(), SNESMonitorGet()
3854: M*/

3856: /*@C
3857:    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3858:    iteration of the nonlinear solver to display the iteration's
3859:    progress.

3861:    Logically Collective on SNES

3863:    Input Parameters:
3864: +  snes - the SNES context
3865: .  f - the monitor function, see SNESMonitorFunction for the calling sequence
3866: .  mctx - [optional] user-defined context for private data for the
3867:           monitor routine (use NULL if no context is desired)
3868: -  monitordestroy - [optional] routine that frees monitor context
3869:           (may be NULL)

3871:    Options Database Keys:
3872: +    -snes_monitor        - sets SNESMonitorDefault()
3873: .    -snes_monitor_lg_residualnorm    - sets line graph monitor,
3874:                             uses SNESMonitorLGCreate()
3875: -    -snes_monitor_cancel - cancels all monitors that have
3876:                             been hardwired into a code by
3877:                             calls to SNESMonitorSet(), but
3878:                             does not cancel those set via
3879:                             the options database.

3881:    Notes:
3882:    Several different monitoring routines may be set by calling
3883:    SNESMonitorSet() multiple times; all will be called in the
3884:    order in which they were set.

3886:    Fortran Notes:
3887:     Only a single monitor function can be set for each SNES object

3889:    Level: intermediate

3891: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3892: @*/
3893: PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3894: {
3895:   PetscInt       i;
3897:   PetscBool      identical;

3901:   for (i=0; i<snes->numbermonitors;i++) {
3902:     PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))snes->monitor[i],snes->monitorcontext[i],snes->monitordestroy[i],&identical);
3903:     if (identical) return(0);
3904:   }
3905:   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3906:   snes->monitor[snes->numbermonitors]          = f;
3907:   snes->monitordestroy[snes->numbermonitors]   = monitordestroy;
3908:   snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3909:   return(0);
3910: }

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

3915:    Logically Collective on SNES

3917:    Input Parameters:
3918: .  snes - the SNES context

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

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

3928:    Level: intermediate

3930: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3931: @*/
3932: PetscErrorCode  SNESMonitorCancel(SNES snes)
3933: {
3935:   PetscInt       i;

3939:   for (i=0; i<snes->numbermonitors; i++) {
3940:     if (snes->monitordestroy[i]) {
3941:       (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3942:     }
3943:   }
3944:   snes->numbermonitors = 0;
3945:   return(0);
3946: }

3948: /*MC
3949:     SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver

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

3955:      Collective on snes

3957:     Input Parameters:
3958: +    snes - the SNES context
3959: .    it - current iteration (0 is the first and is before any Newton step)
3960: .    xnorm - 2-norm of current iterate
3961: .    gnorm - 2-norm of current step
3962: .    f - 2-norm of function
3963: -    cctx - [optional] convergence context

3965:     Output Parameter:
3966: .    reason - reason for convergence/divergence, only needs to be set when convergence or divergence is detected

3968:    Level: intermediate

3970: .seealso:   SNESSetConvergenceTest(), SNESGetConvergenceTest()
3971: M*/

3973: /*@C
3974:    SNESSetConvergenceTest - Sets the function that is to be used
3975:    to test for convergence of the nonlinear iterative solution.

3977:    Logically Collective on SNES

3979:    Input Parameters:
3980: +  snes - the SNES context
3981: .  SNESConvergenceTestFunction - routine to test for convergence
3982: .  cctx - [optional] context for private data for the convergence routine  (may be NULL)
3983: -  destroy - [optional] destructor for the context (may be NULL; PETSC_NULL_FUNCTION in Fortran)

3985:    Level: advanced

3987: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3988: @*/
3989: PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3990: {

3995:   if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
3996:   if (snes->ops->convergeddestroy) {
3997:     (*snes->ops->convergeddestroy)(snes->cnvP);
3998:   }
3999:   snes->ops->converged        = SNESConvergenceTestFunction;
4000:   snes->ops->convergeddestroy = destroy;
4001:   snes->cnvP                  = cctx;
4002:   return(0);
4003: }

4005: /*@
4006:    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.

4008:    Not Collective

4010:    Input Parameter:
4011: .  snes - the SNES context

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

4017:    Options Database:
4018: .   -snes_converged_reason - prints the reason to standard out

4020:    Level: intermediate

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

4025: .seealso: SNESSetConvergenceTest(), SNESSetConvergedReason(), SNESConvergedReason
4026: @*/
4027: PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
4028: {
4032:   *reason = snes->reason;
4033:   return(0);
4034: }

4036: /*@
4037:    SNESSetConvergedReason - Sets the reason the SNES iteration was stopped.

4039:    Not Collective

4041:    Input Parameters:
4042: +  snes - the SNES context
4043: -  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
4044:             manual pages for the individual convergence tests for complete lists

4046:    Level: intermediate

4048: .seealso: SNESGetConvergedReason(), SNESSetConvergenceTest(), SNESConvergedReason
4049: @*/
4050: PetscErrorCode SNESSetConvergedReason(SNES snes,SNESConvergedReason reason)
4051: {
4054:   snes->reason = reason;
4055:   return(0);
4056: }

4058: /*@
4059:    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.

4061:    Logically Collective on SNES

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

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

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

4079:    Level: intermediate

4081: .seealso: SNESGetConvergenceHistory()

4083: @*/
4084: PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
4085: {

4092:   if (!a) {
4093:     if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
4094:     PetscCalloc2(na,&a,na,&its);
4095:     snes->conv_hist_alloc = PETSC_TRUE;
4096:   }
4097:   snes->conv_hist       = a;
4098:   snes->conv_hist_its   = its;
4099:   snes->conv_hist_max   = na;
4100:   snes->conv_hist_len   = 0;
4101:   snes->conv_hist_reset = reset;
4102:   return(0);
4103: }

4105: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4106: #include <engine.h>   /* MATLAB include file */
4107: #include <mex.h>      /* MATLAB include file */

4109: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
4110: {
4111:   mxArray   *mat;
4112:   PetscInt  i;
4113:   PetscReal *ar;

4116:   mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
4117:   ar  = (PetscReal*) mxGetData(mat);
4118:   for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
4119:   PetscFunctionReturn(mat);
4120: }
4121: #endif

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

4126:    Not Collective

4128:    Input Parameter:
4129: .  snes - iterative context obtained from SNESCreate()

4131:    Output Parameters:
4132: +  a   - array to hold history
4133: .  its - integer array holds the number of linear iterations (or
4134:          negative if not converged) for each solve.
4135: -  na  - size of a and its

4137:    Notes:
4138:     The calling sequence for this routine in Fortran is
4139: $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)

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

4145:    Level: intermediate

4147: .seealso: SNESSetConvergencHistory()

4149: @*/
4150: PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
4151: {
4154:   if (a)   *a   = snes->conv_hist;
4155:   if (its) *its = snes->conv_hist_its;
4156:   if (na)  *na  = snes->conv_hist_len;
4157:   return(0);
4158: }

4160: /*@C
4161:   SNESSetUpdate - Sets the general-purpose update function called
4162:   at the beginning of every iteration of the nonlinear solve. Specifically
4163:   it is called just before the Jacobian is "evaluated".

4165:   Logically Collective on SNES

4167:   Input Parameters:
4168: + snes - The nonlinear solver context
4169: - func - The function

4171:   Calling sequence of func:
4172: $ func (SNES snes, PetscInt step);

4174: . step - The current step of the iteration

4176:   Level: advanced

4178:   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()
4179:         This is not used by most users.

4181: .seealso SNESSetJacobian(), SNESSolve()
4182: @*/
4183: PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
4184: {
4187:   snes->ops->update = func;
4188:   return(0);
4189: }

4191: /*
4192:    SNESScaleStep_Private - Scales a step so that its length is less than the
4193:    positive parameter delta.

4195:     Input Parameters:
4196: +   snes - the SNES context
4197: .   y - approximate solution of linear system
4198: .   fnorm - 2-norm of current function
4199: -   delta - trust region size

4201:     Output Parameters:
4202: +   gpnorm - predicted function norm at the new point, assuming local
4203:     linearization.  The value is zero if the step lies within the trust
4204:     region, and exceeds zero otherwise.
4205: -   ynorm - 2-norm of the step

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

4211: */
4212: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
4213: {
4214:   PetscReal      nrm;
4215:   PetscScalar    cnorm;


4223:   VecNorm(y,NORM_2,&nrm);
4224:   if (nrm > *delta) {
4225:     nrm     = *delta/nrm;
4226:     *gpnorm = (1.0 - nrm)*(*fnorm);
4227:     cnorm   = nrm;
4228:     VecScale(y,cnorm);
4229:     *ynorm  = *delta;
4230:   } else {
4231:     *gpnorm = 0.0;
4232:     *ynorm  = nrm;
4233:   }
4234:   return(0);
4235: }

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

4240:    Collective on SNES

4242:    Parameter:
4243: +  snes - iterative context obtained from SNESCreate()
4244: -  viewer - the viewer to display the reason


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

4250:    Level: beginner

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

4254: @*/
4255: PetscErrorCode  SNESReasonView(SNES snes,PetscViewer viewer)
4256: {
4257:   PetscViewerFormat format;
4258:   PetscBool         isAscii;
4259:   PetscErrorCode    ierr;

4262:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
4263:   if (isAscii) {
4264:     PetscViewerGetFormat(viewer, &format);
4265:     PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
4266:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
4267:       DM                dm;
4268:       Vec               u;
4269:       PetscDS           prob;
4270:       PetscInt          Nf, f;
4271:       PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
4272:       void            **exactCtx;
4273:       PetscReal         error;

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

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

4307:   Collective on SNES

4309:   Input Parameters:
4310: . snes   - the SNES object

4312:   Level: intermediate

4314: @*/
4315: PetscErrorCode SNESReasonViewFromOptions(SNES snes)
4316: {
4317:   PetscErrorCode    ierr;
4318:   PetscViewer       viewer;
4319:   PetscBool         flg;
4320:   static PetscBool  incall = PETSC_FALSE;
4321:   PetscViewerFormat format;

4324:   if (incall) return(0);
4325:   incall = PETSC_TRUE;
4326:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
4327:   if (flg) {
4328:     PetscViewerPushFormat(viewer,format);
4329:     SNESReasonView(snes,viewer);
4330:     PetscViewerPopFormat(viewer);
4331:     PetscViewerDestroy(&viewer);
4332:   }
4333:   incall = PETSC_FALSE;
4334:   return(0);
4335: }

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

4341:    Collective on SNES

4343:    Input Parameters:
4344: +  snes - the SNES context
4345: .  b - the constant part of the equation F(x) = b, or NULL to use zero.
4346: -  x - the solution vector.

4348:    Notes:
4349:    The user should initialize the vector,x, with the initial guess
4350:    for the nonlinear solve prior to calling SNESSolve.  In particular,
4351:    to employ an initial guess of zero, the user should explicitly set
4352:    this vector to zero by calling VecSet().

4354:    Level: beginner

4356: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
4357: @*/
4358: PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
4359: {
4360:   PetscErrorCode    ierr;
4361:   PetscBool         flg;
4362:   PetscInt          grid;
4363:   Vec               xcreated = NULL;
4364:   DM                dm;


4373:   /* High level operations using the nonlinear solver */
4374:   {
4375:     PetscViewer       viewer;
4376:     PetscViewerFormat format;
4377:     PetscInt          num;
4378:     PetscBool         flg;
4379:     static PetscBool  incall = PETSC_FALSE;

4381:     if (!incall) {
4382:       /* Estimate the convergence rate of the discretization */
4383:       PetscOptionsGetViewer(PetscObjectComm((PetscObject) snes),((PetscObject)snes)->options, ((PetscObject) snes)->prefix, "-snes_convergence_estimate", &viewer, &format, &flg);
4384:       if (flg) {
4385:         PetscConvEst conv;
4386:         DM           dm;
4387:         PetscReal   *alpha; /* Convergence rate of the solution error for each field in the L_2 norm */
4388:         PetscInt     Nf;

4390:         incall = PETSC_TRUE;
4391:         SNESGetDM(snes, &dm);
4392:         DMGetNumFields(dm, &Nf);
4393:         PetscCalloc1(Nf, &alpha);
4394:         PetscConvEstCreate(PetscObjectComm((PetscObject) snes), &conv);
4395:         PetscConvEstSetSolver(conv, snes);
4396:         PetscConvEstSetFromOptions(conv);
4397:         PetscConvEstSetUp(conv);
4398:         PetscConvEstGetConvRate(conv, alpha);
4399:         PetscViewerPushFormat(viewer, format);
4400:         PetscConvEstRateView(conv, alpha, viewer);
4401:         PetscViewerPopFormat(viewer);
4402:         PetscViewerDestroy(&viewer);
4403:         PetscConvEstDestroy(&conv);
4404:         PetscFree(alpha);
4405:         incall = PETSC_FALSE;
4406:       }
4407:       /* Adaptively refine the initial grid */
4408:       num  = 1;
4409:       PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_initial", &num, &flg);
4410:       if (flg) {
4411:         DMAdaptor adaptor;

4413:         incall = PETSC_TRUE;
4414:         DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);
4415:         DMAdaptorSetSolver(adaptor, snes);
4416:         DMAdaptorSetSequenceLength(adaptor, num);
4417:         DMAdaptorSetFromOptions(adaptor);
4418:         DMAdaptorSetUp(adaptor);
4419:         DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_INITIAL, &dm, &x);
4420:         DMAdaptorDestroy(&adaptor);
4421:         incall = PETSC_FALSE;
4422:       }
4423:       /* Use grid sequencing to adapt */
4424:       num  = 0;
4425:       PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_sequence", &num, NULL);
4426:       if (num) {
4427:         DMAdaptor adaptor;

4429:         incall = PETSC_TRUE;
4430:         DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);
4431:         DMAdaptorSetSolver(adaptor, snes);
4432:         DMAdaptorSetSequenceLength(adaptor, num);
4433:         DMAdaptorSetFromOptions(adaptor);
4434:         DMAdaptorSetUp(adaptor);
4435:         DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_SEQUENTIAL, &dm, &x);
4436:         DMAdaptorDestroy(&adaptor);
4437:         incall = PETSC_FALSE;
4438:       }
4439:     }
4440:   }
4441:   if (!x) {
4442:     SNESGetDM(snes,&dm);
4443:     DMCreateGlobalVector(dm,&xcreated);
4444:     x    = xcreated;
4445:   }
4446:   SNESViewFromOptions(snes,NULL,"-snes_view_pre");

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

4451:     /* set solution vector */
4452:     if (!grid) {PetscObjectReference((PetscObject)x);}
4453:     VecDestroy(&snes->vec_sol);
4454:     snes->vec_sol = x;
4455:     SNESGetDM(snes,&dm);

4457:     /* set affine vector if provided */
4458:     if (b) { PetscObjectReference((PetscObject)b); }
4459:     VecDestroy(&snes->vec_rhs);
4460:     snes->vec_rhs = b;

4462:     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");
4463:     if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
4464:     if (snes->vec_rhs  == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
4465:     if (!snes->vec_sol_update /* && snes->vec_sol */) {
4466:       VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
4467:       PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
4468:     }
4469:     DMShellSetGlobalVector(dm,snes->vec_sol);
4470:     SNESSetUp(snes);

4472:     if (!grid) {
4473:       if (snes->ops->computeinitialguess) {
4474:         (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
4475:       }
4476:     }

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

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

4487:     if (snes->lagjac_persist) snes->jac_iter += snes->iter;
4488:     if (snes->lagpre_persist) snes->pre_iter += snes->iter;

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

4494:     if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
4495:     if (snes->reason < 0) break;
4496:     if (grid <  snes->gridsequence) {
4497:       DM  fine;
4498:       Vec xnew;
4499:       Mat interp;

4501:       DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
4502:       if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
4503:       DMCreateInterpolation(snes->dm,fine,&interp,NULL);
4504:       DMCreateGlobalVector(fine,&xnew);
4505:       MatInterpolate(interp,x,xnew);
4506:       DMInterpolate(snes->dm,interp,fine);
4507:       MatDestroy(&interp);
4508:       x    = xnew;

4510:       SNESReset(snes);
4511:       SNESSetDM(snes,fine);
4512:       SNESResetFromOptions(snes);
4513:       DMDestroy(&fine);
4514:       PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
4515:     }
4516:   }
4517:   SNESViewFromOptions(snes,NULL,"-snes_view");
4518:   VecViewFromOptions(snes->vec_sol,(PetscObject)snes,"-snes_view_solution");

4520:   VecDestroy(&xcreated);
4521:   PetscObjectSAWsBlock((PetscObject)snes);
4522:   return(0);
4523: }

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

4527: /*@C
4528:    SNESSetType - Sets the method for the nonlinear solver.

4530:    Collective on SNES

4532:    Input Parameters:
4533: +  snes - the SNES context
4534: -  type - a known method

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

4540:    Notes:
4541:    See "petsc/include/petscsnes.h" for available methods (for instance)
4542: +    SNESNEWTONLS - Newton's method with line search
4543:      (systems of nonlinear equations)
4544: -    SNESNEWTONTR - Newton's method with trust region
4545:      (systems of nonlinear equations)

4547:   Normally, it is best to use the SNESSetFromOptions() command and then
4548:   set the SNES solver type from the options database rather than by using
4549:   this routine.  Using the options database provides the user with
4550:   maximum flexibility in evaluating the many nonlinear solvers.
4551:   The SNESSetType() routine is provided for those situations where it
4552:   is necessary to set the nonlinear solver independently of the command
4553:   line or options database.  This might be the case, for example, when
4554:   the choice of solver changes during the execution of the program,
4555:   and the user's application is taking responsibility for choosing the
4556:   appropriate method.

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

4562:   Level: intermediate

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

4566: @*/
4567: PetscErrorCode  SNESSetType(SNES snes,SNESType type)
4568: {
4569:   PetscErrorCode ierr,(*r)(SNES);
4570:   PetscBool      match;


4576:   PetscObjectTypeCompare((PetscObject)snes,type,&match);
4577:   if (match) return(0);

4579:   PetscFunctionListFind(SNESList,type,&r);
4580:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
4581:   /* Destroy the previous private SNES context */
4582:   if (snes->ops->destroy) {
4583:     (*(snes)->ops->destroy)(snes);
4584:     snes->ops->destroy = NULL;
4585:   }
4586:   /* Reinitialize function pointers in SNESOps structure */
4587:   snes->ops->setup          = 0;
4588:   snes->ops->solve          = 0;
4589:   snes->ops->view           = 0;
4590:   snes->ops->setfromoptions = 0;
4591:   snes->ops->destroy        = 0;

4593:   /* It may happen the user has customized the line search before calling SNESSetType */
4594:   if (((PetscObject)snes)->type_name) {
4595:     SNESLineSearchDestroy(&snes->linesearch);
4596:   }

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

4601:   PetscObjectChangeTypeName((PetscObject)snes,type);
4602:   (*r)(snes);
4603:   return(0);
4604: }

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

4609:    Not Collective

4611:    Input Parameter:
4612: .  snes - nonlinear solver context

4614:    Output Parameter:
4615: .  type - SNES method (a character string)

4617:    Level: intermediate

4619: @*/
4620: PetscErrorCode  SNESGetType(SNES snes,SNESType *type)
4621: {
4625:   *type = ((PetscObject)snes)->type_name;
4626:   return(0);
4627: }

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

4632:   Logically Collective on SNES

4634:   Input Parameters:
4635: + snes - the SNES context obtained from SNESCreate()
4636: - u    - the solution vector

4638:   Level: beginner

4640: @*/
4641: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4642: {
4643:   DM             dm;

4649:   PetscObjectReference((PetscObject) u);
4650:   VecDestroy(&snes->vec_sol);

4652:   snes->vec_sol = u;

4654:   SNESGetDM(snes, &dm);
4655:   DMShellSetGlobalVector(dm, u);
4656:   return(0);
4657: }

4659: /*@
4660:    SNESGetSolution - Returns the vector where the approximate solution is
4661:    stored. This is the fine grid solution when using SNESSetGridSequence().

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

4665:    Input Parameter:
4666: .  snes - the SNES context

4668:    Output Parameter:
4669: .  x - the solution

4671:    Level: intermediate

4673: .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
4674: @*/
4675: PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
4676: {
4680:   *x = snes->vec_sol;
4681:   return(0);
4682: }

4684: /*@
4685:    SNESGetSolutionUpdate - Returns the vector where the solution update is
4686:    stored.

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

4690:    Input Parameter:
4691: .  snes - the SNES context

4693:    Output Parameter:
4694: .  x - the solution update

4696:    Level: advanced

4698: .seealso: SNESGetSolution(), SNESGetFunction()
4699: @*/
4700: PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
4701: {
4705:   *x = snes->vec_sol_update;
4706:   return(0);
4707: }

4709: /*@C
4710:    SNESGetFunction - Returns the vector where the function is stored.

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

4714:    Input Parameter:
4715: .  snes - the SNES context

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

4722:    Level: advanced

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

4726: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4727: @*/
4728: PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4729: {
4731:   DM             dm;

4735:   if (r) {
4736:     if (!snes->vec_func) {
4737:       if (snes->vec_rhs) {
4738:         VecDuplicate(snes->vec_rhs,&snes->vec_func);
4739:       } else if (snes->vec_sol) {
4740:         VecDuplicate(snes->vec_sol,&snes->vec_func);
4741:       } else if (snes->dm) {
4742:         DMCreateGlobalVector(snes->dm,&snes->vec_func);
4743:       }
4744:     }
4745:     *r = snes->vec_func;
4746:   }
4747:   SNESGetDM(snes,&dm);
4748:   DMSNESGetFunction(dm,f,ctx);
4749:   return(0);
4750: }

4752: /*@C
4753:    SNESGetNGS - Returns the NGS function and context.

4755:    Input Parameter:
4756: .  snes - the SNES context

4758:    Output Parameter:
4759: +  f - the function (or NULL) see SNESNGSFunction for details
4760: -  ctx    - the function context (or NULL)

4762:    Level: advanced

4764: .seealso: SNESSetNGS(), SNESGetFunction()
4765: @*/

4767: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4768: {
4770:   DM             dm;

4774:   SNESGetDM(snes,&dm);
4775:   DMSNESGetNGS(dm,f,ctx);
4776:   return(0);
4777: }

4779: /*@C
4780:    SNESSetOptionsPrefix - Sets the prefix used for searching for all
4781:    SNES options in the database.

4783:    Logically Collective on SNES

4785:    Input Parameter:
4786: +  snes - the SNES context
4787: -  prefix - the prefix to prepend to all option names

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

4793:    Level: advanced

4795: .seealso: SNESSetFromOptions()
4796: @*/
4797: PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
4798: {

4803:   PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4804:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4805:   if (snes->linesearch) {
4806:     SNESGetLineSearch(snes,&snes->linesearch);
4807:     PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4808:   }
4809:   KSPSetOptionsPrefix(snes->ksp,prefix);
4810:   return(0);
4811: }

4813: /*@C
4814:    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4815:    SNES options in the database.

4817:    Logically Collective on SNES

4819:    Input Parameters:
4820: +  snes - the SNES context
4821: -  prefix - the prefix to prepend to all option names

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

4827:    Level: advanced

4829: .seealso: SNESGetOptionsPrefix()
4830: @*/
4831: PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4832: {

4837:   PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4838:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4839:   if (snes->linesearch) {
4840:     SNESGetLineSearch(snes,&snes->linesearch);
4841:     PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4842:   }
4843:   KSPAppendOptionsPrefix(snes->ksp,prefix);
4844:   return(0);
4845: }

4847: /*@C
4848:    SNESGetOptionsPrefix - Sets the prefix used for searching for all
4849:    SNES options in the database.

4851:    Not Collective

4853:    Input Parameter:
4854: .  snes - the SNES context

4856:    Output Parameter:
4857: .  prefix - pointer to the prefix string used

4859:    Notes:
4860:     On the fortran side, the user should pass in a string 'prefix' of
4861:    sufficient length to hold the prefix.

4863:    Level: advanced

4865: .seealso: SNESAppendOptionsPrefix()
4866: @*/
4867: PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4868: {

4873:   PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4874:   return(0);
4875: }


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

4881:    Not collective

4883:    Input Parameters:
4884: +  name_solver - name of a new user-defined solver
4885: -  routine_create - routine to create method context

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

4890:    Sample usage:
4891: .vb
4892:    SNESRegister("my_solver",MySolverCreate);
4893: .ve

4895:    Then, your solver can be chosen with the procedural interface via
4896: $     SNESSetType(snes,"my_solver")
4897:    or at runtime via the option
4898: $     -snes_type my_solver

4900:    Level: advanced

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

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

4906:   Level: advanced
4907: @*/
4908: PetscErrorCode  SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4909: {

4913:   SNESInitializePackage();
4914:   PetscFunctionListAdd(&SNESList,sname,function);
4915:   return(0);
4916: }

4918: PetscErrorCode  SNESTestLocalMin(SNES snes)
4919: {
4921:   PetscInt       N,i,j;
4922:   Vec            u,uh,fh;
4923:   PetscScalar    value;
4924:   PetscReal      norm;

4927:   SNESGetSolution(snes,&u);
4928:   VecDuplicate(u,&uh);
4929:   VecDuplicate(u,&fh);

4931:   /* currently only works for sequential */
4932:   PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4933:   VecGetSize(u,&N);
4934:   for (i=0; i<N; i++) {
4935:     VecCopy(u,uh);
4936:     PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4937:     for (j=-10; j<11; j++) {
4938:       value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4939:       VecSetValue(uh,i,value,ADD_VALUES);
4940:       SNESComputeFunction(snes,uh,fh);
4941:       VecNorm(fh,NORM_2,&norm);
4942:       PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);
4943:       value = -value;
4944:       VecSetValue(uh,i,value,ADD_VALUES);
4945:     }
4946:   }
4947:   VecDestroy(&uh);
4948:   VecDestroy(&fh);
4949:   return(0);
4950: }

4952: /*@
4953:    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4954:    computing relative tolerance for linear solvers within an inexact
4955:    Newton method.

4957:    Logically Collective on SNES

4959:    Input Parameters:
4960: +  snes - SNES context
4961: -  flag - PETSC_TRUE or PETSC_FALSE

4963:     Options Database:
4964: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4965: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
4966: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4967: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4968: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
4969: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
4970: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4971: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

4973:    Notes:
4974:    Currently, the default is to use a constant relative tolerance for
4975:    the inner linear solvers.  Alternatively, one can use the
4976:    Eisenstat-Walker method, where the relative convergence tolerance
4977:    is reset at each Newton iteration according progress of the nonlinear
4978:    solver.

4980:    Level: advanced

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

4986: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4987: @*/
4988: PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool flag)
4989: {
4993:   snes->ksp_ewconv = flag;
4994:   return(0);
4995: }

4997: /*@
4998:    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
4999:    for computing relative tolerance for linear solvers within an
5000:    inexact Newton method.

5002:    Not Collective

5004:    Input Parameter:
5005: .  snes - SNES context

5007:    Output Parameter:
5008: .  flag - PETSC_TRUE or PETSC_FALSE

5010:    Notes:
5011:    Currently, the default is to use a constant relative tolerance for
5012:    the inner linear solvers.  Alternatively, one can use the
5013:    Eisenstat-Walker method, where the relative convergence tolerance
5014:    is reset at each Newton iteration according progress of the nonlinear
5015:    solver.

5017:    Level: advanced

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

5023: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
5024: @*/
5025: PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
5026: {
5030:   *flag = snes->ksp_ewconv;
5031:   return(0);
5032: }

5034: /*@
5035:    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
5036:    convergence criteria for the linear solvers within an inexact
5037:    Newton method.

5039:    Logically Collective on SNES

5041:    Input Parameters:
5042: +    snes - SNES context
5043: .    version - version 1, 2 (default is 2) or 3
5044: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5045: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5046: .    gamma - multiplicative factor for version 2 rtol computation
5047:              (0 <= gamma2 <= 1)
5048: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
5049: .    alpha2 - power for safeguard
5050: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

5052:    Note:
5053:    Version 3 was contributed by Luis Chacon, June 2006.

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

5057:    Level: advanced

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

5064: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
5065: @*/
5066: PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
5067: {
5068:   SNESKSPEW *kctx;

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

5082:   if (version != PETSC_DEFAULT)   kctx->version   = version;
5083:   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
5084:   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
5085:   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
5086:   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
5087:   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
5088:   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;

5090:   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);
5091:   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);
5092:   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);
5093:   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);
5094:   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);
5095:   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);
5096:   return(0);
5097: }

5099: /*@
5100:    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
5101:    convergence criteria for the linear solvers within an inexact
5102:    Newton method.

5104:    Not Collective

5106:    Input Parameters:
5107:      snes - SNES context

5109:    Output Parameters:
5110: +    version - version 1, 2 (default is 2) or 3
5111: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5112: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5113: .    gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
5114: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
5115: .    alpha2 - power for safeguard
5116: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

5118:    Level: advanced

5120: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
5121: @*/
5122: PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
5123: {
5124:   SNESKSPEW *kctx;

5128:   kctx = (SNESKSPEW*)snes->kspconvctx;
5129:   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
5130:   if (version)   *version   = kctx->version;
5131:   if (rtol_0)    *rtol_0    = kctx->rtol_0;
5132:   if (rtol_max)  *rtol_max  = kctx->rtol_max;
5133:   if (gamma)     *gamma     = kctx->gamma;
5134:   if (alpha)     *alpha     = kctx->alpha;
5135:   if (alpha2)    *alpha2    = kctx->alpha2;
5136:   if (threshold) *threshold = kctx->threshold;
5137:   return(0);
5138: }

5140:  PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5141: {
5143:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
5144:   PetscReal      rtol  = PETSC_DEFAULT,stol;

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

5181: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5182: {
5184:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
5185:   PCSide         pcside;
5186:   Vec            lres;

5189:   if (!snes->ksp_ewconv) return(0);
5190:   KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
5191:   kctx->norm_last = snes->norm;
5192:   if (kctx->version == 1) {
5193:     PC        pc;
5194:     PetscBool isNone;

5196:     KSPGetPC(ksp, &pc);
5197:     PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
5198:     KSPGetPCSide(ksp,&pcside);
5199:      if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
5200:       /* KSP residual is true linear residual */
5201:       KSPGetResidualNorm(ksp,&kctx->lresid_last);
5202:     } else {
5203:       /* KSP residual is preconditioned residual */
5204:       /* compute true linear residual norm */
5205:       VecDuplicate(b,&lres);
5206:       MatMult(snes->jacobian,x,lres);
5207:       VecAYPX(lres,-1.0,b);
5208:       VecNorm(lres,NORM_2,&kctx->lresid_last);
5209:       VecDestroy(&lres);
5210:     }
5211:   }
5212:   return(0);
5213: }

5215: /*@
5216:    SNESGetKSP - Returns the KSP context for a SNES solver.

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

5220:    Input Parameter:
5221: .  snes - the SNES context

5223:    Output Parameter:
5224: .  ksp - the KSP context

5226:    Notes:
5227:    The user can then directly manipulate the KSP context to set various
5228:    options, etc.  Likewise, the user can then extract and manipulate the
5229:    PC contexts as well.

5231:    Level: beginner

5233: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
5234: @*/
5235: PetscErrorCode  SNESGetKSP(SNES snes,KSP *ksp)
5236: {


5243:   if (!snes->ksp) {
5244:     PetscBool monitor = PETSC_FALSE;

5246:     KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
5247:     PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
5248:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);

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

5253:     PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
5254:     if (monitor) {
5255:       KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
5256:     }
5257:     monitor = PETSC_FALSE;
5258:     PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
5259:     if (monitor) {
5260:       PetscObject *objs;
5261:       KSPMonitorSNESLGResidualNormCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
5262:       objs[0] = (PetscObject) snes;
5263:       KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
5264:     }
5265:     PetscObjectSetOptions((PetscObject)snes->ksp,((PetscObject)snes)->options);
5266:   }
5267:   *ksp = snes->ksp;
5268:   return(0);
5269: }


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

5276:    Logically Collective on SNES

5278:    Input Parameters:
5279: +  snes - the nonlinear solver context
5280: -  dm - the dm, cannot be NULL

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

5287:    Level: intermediate

5289: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
5290: @*/
5291: PetscErrorCode  SNESSetDM(SNES snes,DM dm)
5292: {
5294:   KSP            ksp;
5295:   DMSNES         sdm;

5300:   PetscObjectReference((PetscObject)dm);
5301:   if (snes->dm) {               /* Move the DMSNES context over to the new DM unless the new DM already has one */
5302:     if (snes->dm->dmsnes && !dm->dmsnes) {
5303:       DMCopyDMSNES(snes->dm,dm);
5304:       DMGetDMSNES(snes->dm,&sdm);
5305:       if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
5306:     }
5307:     DMCoarsenHookRemove(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
5308:     DMDestroy(&snes->dm);
5309:   }
5310:   snes->dm     = dm;
5311:   snes->dmAuto = PETSC_FALSE;

5313:   SNESGetKSP(snes,&ksp);
5314:   KSPSetDM(ksp,dm);
5315:   KSPSetDMActive(ksp,PETSC_FALSE);
5316:   if (snes->npc) {
5317:     SNESSetDM(snes->npc, snes->dm);
5318:     SNESSetNPCSide(snes,snes->npcside);
5319:   }
5320:   return(0);
5321: }

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

5326:    Not Collective but DM obtained is parallel on SNES

5328:    Input Parameter:
5329: . snes - the preconditioner context

5331:    Output Parameter:
5332: .  dm - the dm

5334:    Level: intermediate

5336: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
5337: @*/
5338: PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
5339: {

5344:   if (!snes->dm) {
5345:     DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
5346:     snes->dmAuto = PETSC_TRUE;
5347:   }
5348:   *dm = snes->dm;
5349:   return(0);
5350: }

5352: /*@
5353:   SNESSetNPC - Sets the nonlinear preconditioner to be used.

5355:   Collective on SNES

5357:   Input Parameters:
5358: + snes - iterative context obtained from SNESCreate()
5359: - pc   - the preconditioner object

5361:   Notes:
5362:   Use SNESGetNPC() to retrieve the preconditioner context (for example,
5363:   to configure it using the API).

5365:   Level: developer

5367: .seealso: SNESGetNPC(), SNESHasNPC()
5368: @*/
5369: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
5370: {

5377:   PetscObjectReference((PetscObject) pc);
5378:   SNESDestroy(&snes->npc);
5379:   snes->npc = pc;
5380:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->npc);
5381:   return(0);
5382: }

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

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

5389:   Input Parameter:
5390: . snes - iterative context obtained from SNESCreate()

5392:   Output Parameter:
5393: . pc - preconditioner context

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

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

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

5404:   Level: developer

5406: .seealso: SNESSetNPC(), SNESHasNPC(), SNES, SNESCreate()
5407: @*/
5408: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
5409: {
5411:   const char     *optionsprefix;

5416:   if (!snes->npc) {
5417:     SNESCreate(PetscObjectComm((PetscObject)snes),&snes->npc);
5418:     PetscObjectIncrementTabLevel((PetscObject)snes->npc,(PetscObject)snes,1);
5419:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->npc);
5420:     SNESGetOptionsPrefix(snes,&optionsprefix);
5421:     SNESSetOptionsPrefix(snes->npc,optionsprefix);
5422:     SNESAppendOptionsPrefix(snes->npc,"npc_");
5423:     SNESSetCountersReset(snes->npc,PETSC_FALSE);
5424:   }
5425:   *pc = snes->npc;
5426:   return(0);
5427: }

5429: /*@
5430:   SNESHasNPC - Returns whether a nonlinear preconditioner exists

5432:   Not Collective

5434:   Input Parameter:
5435: . snes - iterative context obtained from SNESCreate()

5437:   Output Parameter:
5438: . has_npc - whether the SNES has an NPC or not

5440:   Level: developer

5442: .seealso: SNESSetNPC(), SNESGetNPC()
5443: @*/
5444: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
5445: {
5448:   *has_npc = (PetscBool) (snes->npc ? PETSC_TRUE : PETSC_FALSE);
5449:   return(0);
5450: }

5452: /*@
5453:     SNESSetNPCSide - Sets the preconditioning side.

5455:     Logically Collective on SNES

5457:     Input Parameter:
5458: .   snes - iterative context obtained from SNESCreate()

5460:     Output Parameter:
5461: .   side - the preconditioning side, where side is one of
5462: .vb
5463:       PC_LEFT - left preconditioning
5464:       PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5465: .ve

5467:     Options Database Keys:
5468: .   -snes_pc_side <right,left>

5470:     Notes:
5471:     SNESNRICHARDSON and SNESNCG only support left preconditioning.

5473:     Level: intermediate

5475: .seealso: SNESGetNPCSide(), KSPSetPCSide()
5476: @*/
5477: PetscErrorCode  SNESSetNPCSide(SNES snes,PCSide side)
5478: {
5482:   snes->npcside= side;
5483:   return(0);
5484: }

5486: /*@
5487:     SNESGetNPCSide - Gets the preconditioning side.

5489:     Not Collective

5491:     Input Parameter:
5492: .   snes - iterative context obtained from SNESCreate()

5494:     Output Parameter:
5495: .   side - the preconditioning side, where side is one of
5496: .vb
5497:       PC_LEFT - left preconditioning
5498:       PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5499: .ve

5501:     Level: intermediate

5503: .seealso: SNESSetNPCSide(), KSPGetPCSide()
5504: @*/
5505: PetscErrorCode  SNESGetNPCSide(SNES snes,PCSide *side)
5506: {
5510:   *side = snes->npcside;
5511:   return(0);
5512: }

5514: /*@
5515:   SNESSetLineSearch - Sets the linesearch on the SNES instance.

5517:   Collective on SNES

5519:   Input Parameters:
5520: + snes - iterative context obtained from SNESCreate()
5521: - linesearch   - the linesearch object

5523:   Notes:
5524:   Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
5525:   to configure it using the API).

5527:   Level: developer

5529: .seealso: SNESGetLineSearch()
5530: @*/
5531: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
5532: {

5539:   PetscObjectReference((PetscObject) linesearch);
5540:   SNESLineSearchDestroy(&snes->linesearch);

5542:   snes->linesearch = linesearch;

5544:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5545:   return(0);
5546: }

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

5552:   Not Collective

5554:   Input Parameter:
5555: . snes - iterative context obtained from SNESCreate()

5557:   Output Parameter:
5558: . linesearch - linesearch context

5560:   Level: beginner

5562: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5563: @*/
5564: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5565: {
5567:   const char     *optionsprefix;

5572:   if (!snes->linesearch) {
5573:     SNESGetOptionsPrefix(snes, &optionsprefix);
5574:     SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5575:     SNESLineSearchSetSNES(snes->linesearch, snes);
5576:     SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5577:     PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5578:     PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5579:   }
5580:   *linesearch = snes->linesearch;
5581:   return(0);
5582: }