Actual source code: snes.c

petsc-master 2019-05-25
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: .keywords: SNES, set, initial guess, nonzero

 35: .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
 36: @*/
 37: PetscErrorCode  SNESSetErrorIfNotConverged(SNES snes,PetscBool flg)
 38: {
 42:   snes->errorifnotconverged = flg;
 43:   return(0);
 44: }

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

 49:    Not Collective

 51:    Input Parameter:
 52: .  snes - iterative context obtained from SNESCreate()

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

 57:    Level: intermediate

 59: .keywords: SNES, set, initial guess, nonzero

 61: .seealso:  SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
 62: @*/
 63: PetscErrorCode  SNESGetErrorIfNotConverged(SNES snes,PetscBool  *flag)
 64: {
 68:   *flag = snes->errorifnotconverged;
 69:   return(0);
 70: }

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

 75:    Logically Collective on SNES

 77:     Input Parameters:
 78: +   snes - the shell SNES
 79: -   flg - is the residual computed?

 81:    Level: advanced

 83: .seealso: SNESGetAlwaysComputesFinalResidual()
 84: @*/
 85: PetscErrorCode  SNESSetAlwaysComputesFinalResidual(SNES snes, PetscBool flg)
 86: {
 89:   snes->alwayscomputesfinalresidual = flg;
 90:   return(0);
 91: }

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

 96:    Logically Collective on SNES

 98:     Input Parameter:
 99: .   snes - the shell SNES

101:     Output Parameter:
102: .   flg - is the residual computed?

104:    Level: advanced

106: .seealso: SNESSetAlwaysComputesFinalResidual()
107: @*/
108: PetscErrorCode  SNESGetAlwaysComputesFinalResidual(SNES snes, PetscBool *flg)
109: {
112:   *flg = snes->alwayscomputesfinalresidual;
113:   return(0);
114: }

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

120:    Logically Collective on SNES

122:    Input Parameters:
123: .  snes - the SNES context

125:    Level: advanced

127: .keywords: SNES, view

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

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

143:    Logically Collective on SNES

145:    Input Parameters:
146: .  snes - the SNES context

148:    Level: advanced

150: .keywords: SNES, view

152: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError()
153: @*/
154: PetscErrorCode SNESSetJacobianDomainError(SNES snes)
155: {
158:   if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates computeJacobian does not make sense");
159:   snes->jacobiandomainerror = PETSC_TRUE;
160:   return(0);
161: }

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

167:    Logically Collective on SNES

169:    Input Parameters:
170: .  snes - the SNES context
171: .  flg  - indicates if or not to check jacobian domain error after each Jacobian evaluation

173:    Level: advanced

175: .keywords: SNES, view

177: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESGetCheckJacobianDomainError()
178: @*/
179: PetscErrorCode SNESSetCheckJacobianDomainError(SNES snes, PetscBool flg)
180: {
183:   snes->checkjacdomainerror = flg;
184:   return(0);
185: }

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

190:    Logically Collective on SNES

192:    Input Parameters:
193: .  snes - the SNES context

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

198:    Level: advanced

200: .keywords: SNES, view

202: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESSetCheckJacobianDomainError()
203: @*/
204: PetscErrorCode SNESGetCheckJacobianDomainError(SNES snes, PetscBool *flg)
205: {
209:   *flg = snes->checkjacdomainerror;
210:   return(0);
211: }

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

216:    Logically Collective on SNES

218:    Input Parameters:
219: .  snes - the SNES context

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

224:    Level: advanced

226: .keywords: SNES, view

228: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction()
229: @*/
230: PetscErrorCode  SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror)
231: {
235:   *domainerror = snes->domainerror;
236:   return(0);
237: }

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

242:    Logically Collective on SNES

244:    Input Parameters:
245: .  snes - the SNES context

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

250:    Level: advanced

252: .keywords: SNES, view

254: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction(),SNESGetFunctionDomainError()
255: @*/
256: PetscErrorCode SNESGetJacobianDomainError(SNES snes, PetscBool *domainerror)
257: {
261:   *domainerror = snes->jacobiandomainerror;
262:   return(0);
263: }

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

268:   Collective on PetscViewer

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

275:    Level: intermediate

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

280:   Notes for advanced users:
281:   Most users should not need to know the details of the binary storage
282:   format, since SNESLoad() and TSView() completely hide these details.
283:   But for anyone who's interested, the standard binary matrix storage
284:   format is
285: .vb
286:      has not yet been determined
287: .ve

289: .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad()
290: @*/
291: PetscErrorCode  SNESLoad(SNES snes, PetscViewer viewer)
292: {
294:   PetscBool      isbinary;
295:   PetscInt       classid;
296:   char           type[256];
297:   KSP            ksp;
298:   DM             dm;
299:   DMSNES         dmsnes;

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

307:   PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
308:   if (classid != SNES_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONG,"Not SNES next in file");
309:   PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
310:   SNESSetType(snes, type);
311:   if (snes->ops->load) {
312:     (*snes->ops->load)(snes,viewer);
313:   }
314:   SNESGetDM(snes,&dm);
315:   DMGetDMSNES(dm,&dmsnes);
316:   DMSNESLoad(dmsnes,viewer);
317:   SNESGetKSP(snes,&ksp);
318:   KSPLoad(ksp,viewer);
319:   return(0);
320: }

322:  #include <petscdraw.h>
323: #if defined(PETSC_HAVE_SAWS)
324:  #include <petscviewersaws.h>
325: #endif

327: /*@C
328:    SNESView - Prints the SNES data structure.

330:    Collective on SNES

332:    Input Parameters:
333: +  SNES - the SNES context
334: -  viewer - visualization context

336:    Options Database Key:
337: .  -snes_view - Calls SNESView() at end of SNESSolve()

339:    Notes:
340:    The available visualization contexts include
341: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
342: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
343:          output where only the first processor opens
344:          the file.  All other processors send their
345:          data to the first processor to print.

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

350:    Level: beginner

352: .keywords: SNES, view

354: .seealso: PetscViewerASCIIOpen()
355: @*/
356: PetscErrorCode  SNESView(SNES snes,PetscViewer viewer)
357: {
358:   SNESKSPEW      *kctx;
360:   KSP            ksp;
361:   SNESLineSearch linesearch;
362:   PetscBool      iascii,isstring,isbinary,isdraw;
363:   DMSNES         dmsnes;
364: #if defined(PETSC_HAVE_SAWS)
365:   PetscBool      issaws;
366: #endif

370:   if (!viewer) {
371:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
372:   }

376:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
377:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
378:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
379:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
380: #if defined(PETSC_HAVE_SAWS)
381:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
382: #endif
383:   if (iascii) {
384:     SNESNormSchedule normschedule;
385:     DM               dm;
386:     PetscErrorCode   (*cJ)(SNES,Vec,Mat,Mat,void*);
387:     void             *ctx;

389:     PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);
390:     if (!snes->setupcalled) {
391:       PetscViewerASCIIPrintf(viewer,"  SNES has not been set up so information may be incomplete\n");
392:     }
393:     if (snes->ops->view) {
394:       PetscViewerASCIIPushTab(viewer);
395:       (*snes->ops->view)(snes,viewer);
396:       PetscViewerASCIIPopTab(viewer);
397:     }
398:     PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);
399:     PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);
400:     if (snes->usesksp) {
401:       PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);
402:     }
403:     PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);
404:     SNESGetNormSchedule(snes, &normschedule);
405:     if (normschedule > 0) {PetscViewerASCIIPrintf(viewer,"  norm schedule %s\n",SNESNormSchedules[normschedule]);}
406:     if (snes->gridsequence) {
407:       PetscViewerASCIIPrintf(viewer,"  total number of grid sequence refinements=%D\n",snes->gridsequence);
408:     }
409:     if (snes->ksp_ewconv) {
410:       kctx = (SNESKSPEW*)snes->kspconvctx;
411:       if (kctx) {
412:         PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);
413:         PetscViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);
414:         PetscViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);
415:       }
416:     }
417:     if (snes->lagpreconditioner == -1) {
418:       PetscViewerASCIIPrintf(viewer,"  Preconditioned is never rebuilt\n");
419:     } else if (snes->lagpreconditioner > 1) {
420:       PetscViewerASCIIPrintf(viewer,"  Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);
421:     }
422:     if (snes->lagjacobian == -1) {
423:       PetscViewerASCIIPrintf(viewer,"  Jacobian is never rebuilt\n");
424:     } else if (snes->lagjacobian > 1) {
425:       PetscViewerASCIIPrintf(viewer,"  Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);
426:     }
427:     SNESGetDM(snes,&dm);
428:     DMSNESGetJacobian(dm,&cJ,&ctx);
429:     if (cJ == SNESComputeJacobianDefault) {
430:       PetscViewerASCIIPrintf(viewer,"  Jacobian is built using finite differences one column at a time\n");
431:     } else if (cJ == SNESComputeJacobianDefaultColor) {
432:       PetscViewerASCIIPrintf(viewer,"  Jacobian is built using finite differences with coloring\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;
553:   KSP            ksp;
554:   PC             pc;
555:   PetscBool      match;
557:   MatNullSpace   nullsp;


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

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

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

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

591:     /* This version replaces the user provided Jacobian matrix with a
592:        matrix-free version but still employs the user-provided preconditioner matrix. */
593:     SNESSetJacobian(snes,J,0,0,0);
594:   } else {
595:     /* This version replaces both the user-provided Jacobian and the user-
596:      provided preconditioner Jacobian with the default matrix free version. */
597:     if ((snes->npcside== PC_LEFT) && snes->npc) {
598:       if (!snes->jacobian){SNESSetJacobian(snes,J,0,0,0);}
599:     } else {
600:       SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,0);
601:     }
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:   MatDestroy(&J);
612:   return(0);
613: }

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

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

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

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

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

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

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

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

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

704:    Collective

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

709:    Level: developer

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

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

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

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

777:    Collective on SNES

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

787:    Level: developer

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

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

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

821:    Collective on SNES

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

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

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

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

870:    Level: beginner

872: .keywords: SNES, nonlinear, set, options, database

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);}


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

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

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



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

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

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

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

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

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

1070:   for (i = 0; i < numberofsetfromoptions; i++) {
1071:     (*othersetfromoptions[i])(snes);
1072:   }

1074:   if (snes->ops->setfromoptions) {
1075:     (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
1076:   }

1078:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
1079:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)snes);
1080:   PetscOptionsEnd();

1082:   if (snes->linesearch) {
1083:     SNESGetLineSearch(snes, &snes->linesearch);
1084:     SNESLineSearchSetFromOptions(snes->linesearch);
1085:   }

1087:   if (snes->usesksp) {
1088:     if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
1089:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
1090:     KSPSetFromOptions(snes->ksp);
1091:   }

1093:   /* if someone has set the SNES NPC type, create it. */
1094:   SNESGetOptionsPrefix(snes, &optionsprefix);
1095:   PetscOptionsHasName(((PetscObject)snes)->options,optionsprefix, "-npc_snes_type", &pcset);
1096:   if (pcset && (!snes->npc)) {
1097:     SNESGetNPC(snes, &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: .keywords: SNES, nonlinear, set, options, database

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

1122:   if (snes->setfromoptionscalled) {SNESSetFromOptions(snes);}
1123:   return(0);
1124: }

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

1130:    Logically Collective on SNES

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

1137:    Level: intermediate

1139:    Notes:
1140:    This function is currently not available from Fortran.

1142: .keywords: SNES, nonlinear, set, application, context

1144: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
1145: @*/
1146: PetscErrorCode  SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
1147: {
1150:   snes->ops->usercompute = compute;
1151:   snes->ops->userdestroy = destroy;
1152:   return(0);
1153: }

1155: /*@
1156:    SNESSetApplicationContext - Sets the optional user-defined context for
1157:    the nonlinear solvers.

1159:    Logically Collective on SNES

1161:    Input Parameters:
1162: +  snes - the SNES context
1163: -  usrP - optional user context

1165:    Level: intermediate

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

1171: .keywords: SNES, nonlinear, set, application, context

1173: .seealso: SNESGetApplicationContext()
1174: @*/
1175: PetscErrorCode  SNESSetApplicationContext(SNES snes,void *usrP)
1176: {
1178:   KSP            ksp;

1182:   SNESGetKSP(snes,&ksp);
1183:   KSPSetApplicationContext(ksp,usrP);
1184:   snes->user = usrP;
1185:   return(0);
1186: }

1188: /*@
1189:    SNESGetApplicationContext - Gets the user-defined context for the
1190:    nonlinear solvers.

1192:    Not Collective

1194:    Input Parameter:
1195: .  snes - SNES context

1197:    Output Parameter:
1198: .  usrP - user context

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

1204:    Level: intermediate

1206: .keywords: SNES, nonlinear, get, application, context

1208: .seealso: SNESSetApplicationContext()
1209: @*/
1210: PetscErrorCode  SNESGetApplicationContext(SNES snes,void *usrP)
1211: {
1214:   *(void**)usrP = snes->user;
1215:   return(0);
1216: }

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

1222:    Collective on SNES

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

1229:    Options Database:
1230: + -snes_mf - use matrix free for both the mat and pmat operator
1231: - -snes_mf_operator - use matrix free only for the mat operator

1233:    Level: intermediate

1235: .keywords: SNES, nonlinear, get, iteration, number,

1237: .seealso:   SNESGetUseMatrixFree(), MatCreateSNESMF()
1238: @*/
1239: PetscErrorCode  SNESSetUseMatrixFree(SNES snes,PetscBool mf_operator,PetscBool mf)
1240: {
1245:   if (mf && !mf_operator) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"If using mf must also use mf_operator");
1246:   snes->mf          = mf;
1247:   snes->mf_operator = mf_operator;
1248:   return(0);
1249: }

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

1255:    Collective on SNES

1257:    Input Parameter:
1258: .  snes - SNES context

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

1264:    Options Database:
1265: + -snes_mf - use matrix free for both the mat and pmat operator
1266: - -snes_mf_operator - use matrix free only for the mat operator

1268:    Level: intermediate

1270: .keywords: SNES, nonlinear, get, iteration, number,

1272: .seealso:   SNESSetUseMatrixFree(), MatCreateSNESMF()
1273: @*/
1274: PetscErrorCode  SNESGetUseMatrixFree(SNES snes,PetscBool *mf_operator,PetscBool *mf)
1275: {
1278:   if (mf)          *mf          = snes->mf;
1279:   if (mf_operator) *mf_operator = snes->mf_operator;
1280:   return(0);
1281: }

1283: /*@
1284:    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1285:    at this time.

1287:    Not Collective

1289:    Input Parameter:
1290: .  snes - SNES context

1292:    Output Parameter:
1293: .  iter - iteration number

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

1298:    This is useful for using lagged Jacobians (where one does not recompute the
1299:    Jacobian at each SNES iteration). For example, the code
1300: .vb
1301:       SNESGetIterationNumber(snes,&it);
1302:       if (!(it % 2)) {
1303:         [compute Jacobian here]
1304:       }
1305: .ve
1306:    can be used in your ComputeJacobian() function to cause the Jacobian to be
1307:    recomputed every second SNES iteration.

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

1311:    Level: intermediate

1313: .keywords: SNES, nonlinear, get, iteration, number,

1315: .seealso:   SNESGetLinearSolveIterations()
1316: @*/
1317: PetscErrorCode  SNESGetIterationNumber(SNES snes,PetscInt *iter)
1318: {
1322:   *iter = snes->iter;
1323:   return(0);
1324: }

1326: /*@
1327:    SNESSetIterationNumber - Sets the current iteration number.

1329:    Not Collective

1331:    Input Parameter:
1332: .  snes - SNES context
1333: .  iter - iteration number

1335:    Level: developer

1337: .keywords: SNES, nonlinear, set, iteration, number,

1339: .seealso:   SNESGetLinearSolveIterations()
1340: @*/
1341: PetscErrorCode  SNESSetIterationNumber(SNES snes,PetscInt iter)
1342: {

1347:   PetscObjectSAWsTakeAccess((PetscObject)snes);
1348:   snes->iter = iter;
1349:   PetscObjectSAWsGrantAccess((PetscObject)snes);
1350:   return(0);
1351: }

1353: /*@
1354:    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1355:    attempted by the nonlinear solver.

1357:    Not Collective

1359:    Input Parameter:
1360: .  snes - SNES context

1362:    Output Parameter:
1363: .  nfails - number of unsuccessful steps attempted

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

1368:    Level: intermediate

1370: .keywords: SNES, nonlinear, get, number, unsuccessful, steps

1372: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1373:           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1374: @*/
1375: PetscErrorCode  SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1376: {
1380:   *nfails = snes->numFailures;
1381:   return(0);
1382: }

1384: /*@
1385:    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1386:    attempted by the nonlinear solver before it gives up.

1388:    Not Collective

1390:    Input Parameters:
1391: +  snes     - SNES context
1392: -  maxFails - maximum of unsuccessful steps

1394:    Level: intermediate

1396: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps

1398: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1399:           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1400: @*/
1401: PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1402: {
1405:   snes->maxFailures = maxFails;
1406:   return(0);
1407: }

1409: /*@
1410:    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1411:    attempted by the nonlinear solver before it gives up.

1413:    Not Collective

1415:    Input Parameter:
1416: .  snes     - SNES context

1418:    Output Parameter:
1419: .  maxFails - maximum of unsuccessful steps

1421:    Level: intermediate

1423: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps

1425: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1426:           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()

1428: @*/
1429: PetscErrorCode  SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1430: {
1434:   *maxFails = snes->maxFailures;
1435:   return(0);
1436: }

1438: /*@
1439:    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1440:      done by SNES.

1442:    Not Collective

1444:    Input Parameter:
1445: .  snes     - SNES context

1447:    Output Parameter:
1448: .  nfuncs - number of evaluations

1450:    Level: intermediate

1452:    Notes:
1453:     Reset every time SNESSolve is called unless SNESSetCountersReset() is used.

1455: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps

1457: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1458: @*/
1459: PetscErrorCode  SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1460: {
1464:   *nfuncs = snes->nfuncs;
1465:   return(0);
1466: }

1468: /*@
1469:    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1470:    linear solvers.

1472:    Not Collective

1474:    Input Parameter:
1475: .  snes - SNES context

1477:    Output Parameter:
1478: .  nfails - number of failed solves

1480:    Level: intermediate

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

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

1488: .keywords: SNES, nonlinear, get, number, unsuccessful, steps

1490: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1491: @*/
1492: PetscErrorCode  SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1493: {
1497:   *nfails = snes->numLinearSolveFailures;
1498:   return(0);
1499: }

1501: /*@
1502:    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1503:    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE

1505:    Logically Collective on SNES

1507:    Input Parameters:
1508: +  snes     - SNES context
1509: -  maxFails - maximum allowed linear solve failures

1511:    Level: intermediate

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

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

1519: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps

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

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

1536:    Not Collective

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

1541:    Output Parameter:
1542: .  maxFails - maximum of unsuccessful solves allowed

1544:    Level: intermediate

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

1549: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps

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

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

1566:    Not Collective

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

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

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

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

1580:    Level: intermediate

1582: .keywords: SNES, nonlinear, get, number, linear, iterations

1584: .seealso:  SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1585: @*/
1586: PetscErrorCode  SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1587: {
1591:   *lits = snes->linear_its;
1592:   return(0);
1593: }

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

1599:    Logically Collective on SNES

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

1605:    Notes:
1606:    This defaults to PETSC_TRUE

1608:    Level: developer

1610: .keywords: SNES, nonlinear, set, reset, number, linear, iterations

1612: .seealso:  SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC()
1613: @*/
1614: PetscErrorCode  SNESSetCountersReset(SNES snes,PetscBool reset)
1615: {
1619:   snes->counters_reset = reset;
1620:   return(0);
1621: }


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

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

1629:    Input Parameters:
1630: +  snes - the SNES context
1631: -  ksp - the KSP context

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

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

1640:    Level: developer

1642: .keywords: SNES, nonlinear, get, KSP, context

1644: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1645: @*/
1646: PetscErrorCode  SNESSetKSP(SNES snes,KSP ksp)
1647: {

1654:   PetscObjectReference((PetscObject)ksp);
1655:   if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1656:   snes->ksp = ksp;
1657:   return(0);
1658: }

1660: /* -----------------------------------------------------------*/
1661: /*@
1662:    SNESCreate - Creates a nonlinear solver context.

1664:    Collective on MPI_Comm

1666:    Input Parameters:
1667: .  comm - MPI communicator

1669:    Output Parameter:
1670: .  outsnes - the new SNES context

1672:    Options Database Keys:
1673: +   -snes_mf - Activates default matrix-free Jacobian-vector products,
1674:                and no preconditioning matrix
1675: .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
1676:                products, and a user-provided preconditioning matrix
1677:                as set by SNESSetJacobian()
1678: -   -snes_fd - Uses (slow!) finite differences to compute Jacobian

1680:    Level: beginner

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

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

1692: .keywords: SNES, nonlinear, create, context

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

1696: @*/
1697: PetscErrorCode  SNESCreate(MPI_Comm comm,SNES *outsnes)
1698: {
1700:   SNES           snes;
1701:   SNESKSPEW      *kctx;

1705:   *outsnes = NULL;
1706:   SNESInitializePackage();

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

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

1772:   snes->mf          = PETSC_FALSE;
1773:   snes->mf_operator = PETSC_FALSE;
1774:   snes->mf_version  = 1;

1776:   snes->numLinearSolveFailures = 0;
1777:   snes->maxLinearSolveFailures = 1;

1779:   snes->vizerotolerance = 1.e-8;
1780: #if defined(PETSC_USE_DEBUG)
1781:   snes->checkjacdomainerror = PETSC_TRUE;
1782: #else
1783:   snes->checkjacdomainerror = PETSC_FALSE;
1784: #endif

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

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

1792:   snes->kspconvctx  = (void*)kctx;
1793:   kctx->version     = 2;
1794:   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1795:                              this was too large for some test cases */
1796:   kctx->rtol_last   = 0.0;
1797:   kctx->rtol_max    = .9;
1798:   kctx->gamma       = 1.0;
1799:   kctx->alpha       = .5*(1.0 + PetscSqrtReal(5.0));
1800:   kctx->alpha2      = kctx->alpha;
1801:   kctx->threshold   = .1;
1802:   kctx->lresid_last = 0.0;
1803:   kctx->norm_last   = 0.0;

1805:   *outsnes = snes;
1806:   return(0);
1807: }

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

1812:      Synopsis:
1813:      #include "petscsnes.h"
1814:      PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);

1816:      Input Parameters:
1817: +     snes - the SNES context
1818: .     x    - state at which to evaluate residual
1819: -     ctx     - optional user-defined function context, passed in with SNESSetFunction()

1821:      Output Parameter:
1822: .     f  - vector to put residual (function value)

1824:    Level: intermediate

1826: .seealso:   SNESSetFunction(), SNESGetFunction()
1827: M*/

1829: /*@C
1830:    SNESSetFunction - Sets the function evaluation routine and function
1831:    vector for use by the SNES routines in solving systems of nonlinear
1832:    equations.

1834:    Logically Collective on SNES

1836:    Input Parameters:
1837: +  snes - the SNES context
1838: .  r - vector to store function value
1839: .  f - function evaluation routine; see SNESFunction for calling sequence details
1840: -  ctx - [optional] user-defined context for private data for the
1841:          function evaluation routine (may be NULL)

1843:    Notes:
1844:    The Newton-like methods typically solve linear systems of the form
1845: $      f'(x) x = -f(x),
1846:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

1848:    Level: beginner

1850: .keywords: SNES, nonlinear, set, function

1852: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1853: @*/
1854: PetscErrorCode  SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1855: {
1857:   DM             dm;

1861:   if (r) {
1864:     PetscObjectReference((PetscObject)r);
1865:     VecDestroy(&snes->vec_func);

1867:     snes->vec_func = r;
1868:   }
1869:   SNESGetDM(snes,&dm);
1870:   DMSNESSetFunction(dm,f,ctx);
1871:   return(0);
1872: }


1875: /*@C
1876:    SNESSetInitialFunction - Sets the function vector to be used as the
1877:    function norm at the initialization of the method.  In some
1878:    instances, the user has precomputed the function before calling
1879:    SNESSolve.  This function allows one to avoid a redundant call
1880:    to SNESComputeFunction in that case.

1882:    Logically Collective on SNES

1884:    Input Parameters:
1885: +  snes - the SNES context
1886: -  f - vector to store function value

1888:    Notes:
1889:    This should not be modified during the solution procedure.

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

1893:    Level: developer

1895: .keywords: SNES, nonlinear, set, function

1897: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1898: @*/
1899: PetscErrorCode  SNESSetInitialFunction(SNES snes, Vec f)
1900: {
1902:   Vec            vec_func;

1908:   if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1909:     snes->vec_func_init_set = PETSC_FALSE;
1910:     return(0);
1911:   }
1912:   SNESGetFunction(snes,&vec_func,NULL,NULL);
1913:   VecCopy(f, vec_func);

1915:   snes->vec_func_init_set = PETSC_TRUE;
1916:   return(0);
1917: }

1919: /*@
1920:    SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1921:    of the SNES method.

1923:    Logically Collective on SNES

1925:    Input Parameters:
1926: +  snes - the SNES context
1927: -  normschedule - the frequency of norm computation

1929:    Options Database Key:
1930: .  -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly>

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

1941:    Level: developer

1943: .keywords: SNES, nonlinear, set, function, norm, type

1945: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1946: @*/
1947: PetscErrorCode  SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1948: {
1951:   snes->normschedule = normschedule;
1952:   return(0);
1953: }


1956: /*@
1957:    SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1958:    of the SNES method.

1960:    Logically Collective on SNES

1962:    Input Parameters:
1963: +  snes - the SNES context
1964: -  normschedule - the type of the norm used

1966:    Level: advanced

1968: .keywords: SNES, nonlinear, set, function, norm, type

1970: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1971: @*/
1972: PetscErrorCode  SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1973: {
1976:   *normschedule = snes->normschedule;
1977:   return(0);
1978: }


1981: /*@
1982:   SNESSetFunctionNorm - Sets the last computed residual norm.

1984:   Logically Collective on SNES

1986:   Input Parameters:
1987: + snes - the SNES context

1989: - normschedule - the frequency of norm computation

1991:   Level: developer

1993: .keywords: SNES, nonlinear, set, function, norm, type
1994: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1995: @*/
1996: PetscErrorCode SNESSetFunctionNorm(SNES snes, PetscReal norm)
1997: {
2000:   snes->norm = norm;
2001:   return(0);
2002: }

2004: /*@
2005:   SNESGetFunctionNorm - Gets the last computed norm of the residual

2007:   Not Collective

2009:   Input Parameter:
2010: . snes - the SNES context

2012:   Output Parameter:
2013: . norm - the last computed residual norm

2015:   Level: developer

2017: .keywords: SNES, nonlinear, set, function, norm, type
2018: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2019: @*/
2020: PetscErrorCode SNESGetFunctionNorm(SNES snes, PetscReal *norm)
2021: {
2025:   *norm = snes->norm;
2026:   return(0);
2027: }

2029: /*@
2030:   SNESGetUpdateNorm - Gets the last computed norm of the Newton update

2032:   Not Collective

2034:   Input Parameter:
2035: . snes - the SNES context

2037:   Output Parameter:
2038: . ynorm - the last computed update norm

2040:   Level: developer

2042: .keywords: SNES, nonlinear, set, function, norm, type
2043: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm()
2044: @*/
2045: PetscErrorCode SNESGetUpdateNorm(SNES snes, PetscReal *ynorm)
2046: {
2050:   *ynorm = snes->ynorm;
2051:   return(0);
2052: }

2054: /*@
2055:   SNESGetSolutionNorm - Gets the last computed norm of the solution

2057:   Not Collective

2059:   Input Parameter:
2060: . snes - the SNES context

2062:   Output Parameter:
2063: . xnorm - the last computed solution norm

2065:   Level: developer

2067: .keywords: SNES, nonlinear, set, function, norm, type
2068: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm(), SNESGetUpdateNorm()
2069: @*/
2070: PetscErrorCode SNESGetSolutionNorm(SNES snes, PetscReal *xnorm)
2071: {
2075:   *xnorm = snes->xnorm;
2076:   return(0);
2077: }

2079: /*@C
2080:    SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
2081:    of the SNES method.

2083:    Logically Collective on SNES

2085:    Input Parameters:
2086: +  snes - the SNES context
2087: -  normschedule - the frequency of norm computation

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

2098:    Level: developer

2100: .keywords: SNES, nonlinear, set, function, norm, type

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


2113: /*@C
2114:    SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
2115:    of the SNES method.

2117:    Logically Collective on SNES

2119:    Input Parameters:
2120: +  snes - the SNES context
2121: -  normschedule - the type of the norm used

2123:    Level: advanced

2125: .keywords: SNES, nonlinear, set, function, norm, type

2127: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2128: @*/
2129: PetscErrorCode  SNESGetFunctionType(SNES snes, SNESFunctionType *type)
2130: {
2133:   *type = snes->functype;
2134:   return(0);
2135: }

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

2140:      Synopsis:
2141:      #include <petscsnes.h>
2142: $    SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);

2144: +  X   - solution vector
2145: .  B   - RHS vector
2146: -  ctx - optional user-defined Gauss-Seidel context

2148:    Level: intermediate

2150: .seealso:   SNESSetNGS(), SNESGetNGS()
2151: M*/

2153: /*@C
2154:    SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
2155:    use with composed nonlinear solvers.

2157:    Input Parameters:
2158: +  snes   - the SNES context
2159: .  f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
2160: -  ctx    - [optional] user-defined context for private data for the
2161:             smoother evaluation routine (may be NULL)

2163:    Notes:
2164:    The NGS routines are used by the composed nonlinear solver to generate
2165:     a problem appropriate update to the solution, particularly FAS.

2167:    Level: intermediate

2169: .keywords: SNES, nonlinear, set, Gauss-Seidel

2171: .seealso: SNESGetFunction(), SNESComputeNGS()
2172: @*/
2173: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
2174: {
2176:   DM             dm;

2180:   SNESGetDM(snes,&dm);
2181:   DMSNESSetNGS(dm,f,ctx);
2182:   return(0);
2183: }

2185: PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
2186: {
2188:   DM             dm;
2189:   DMSNES         sdm;

2192:   SNESGetDM(snes,&dm);
2193:   DMGetDMSNES(dm,&sdm);
2194:   if (!sdm->ops->computepfunction) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");
2195:   if (!sdm->ops->computepjacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard Jacobian.");
2196:   /*  A(x)*x - b(x) */
2197:   PetscStackPush("SNES Picard user function");
2198:   (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
2199:   PetscStackPop;
2200:   PetscStackPush("SNES Picard user Jacobian");
2201:   (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
2202:   PetscStackPop;
2203:   VecScale(f,-1.0);
2204:   MatMultAdd(snes->jacobian,x,f,f);
2205:   return(0);
2206: }

2208: PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
2209: {
2211:   /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
2212:   return(0);
2213: }

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

2218:    Logically Collective on SNES

2220:    Input Parameters:
2221: +  snes - the SNES context
2222: .  r - vector to store function value
2223: .  b - function evaluation routine
2224: .  Amat - matrix with which A(x) x - b(x) is to be computed
2225: .  Pmat - matrix from which preconditioner is computed (usually the same as Amat)
2226: .  J  - function to compute matrix value, see SNESJacobianFunction for details on its calling sequence
2227: -  ctx - [optional] user-defined context for private data for the
2228:          function evaluation routine (may be NULL)

2230:    Notes:
2231:     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
2232:     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.

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

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

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

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

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

2248:    Level: intermediate

2250: .keywords: SNES, nonlinear, set, function

2252: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
2253: @*/
2254: 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)
2255: {
2257:   DM             dm;

2261:   SNESGetDM(snes, &dm);
2262:   DMSNESSetPicard(dm,b,J,ctx);
2263:   SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
2264:   SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
2265:   return(0);
2266: }

2268: /*@C
2269:    SNESGetPicard - Returns the context for the Picard iteration

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

2273:    Input Parameter:
2274: .  snes - the SNES context

2276:    Output Parameter:
2277: +  r - the function (or NULL)
2278: .  f - the function (or NULL); see SNESFunction for calling sequence details
2279: .  Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
2280: .  Pmat  - the matrix from which the preconditioner will be constructed (or NULL)
2281: .  J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
2282: -  ctx - the function context (or NULL)

2284:    Level: advanced

2286: .keywords: SNES, nonlinear, get, function

2288: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
2289: @*/
2290: 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)
2291: {
2293:   DM             dm;

2297:   SNESGetFunction(snes,r,NULL,NULL);
2298:   SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
2299:   SNESGetDM(snes,&dm);
2300:   DMSNESGetPicard(dm,f,J,ctx);
2301:   return(0);
2302: }

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

2307:    Logically Collective on SNES

2309:    Input Parameters:
2310: +  snes - the SNES context
2311: .  func - function evaluation routine
2312: -  ctx - [optional] user-defined context for private data for the
2313:          function evaluation routine (may be NULL)

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

2318: .  f - function vector
2319: -  ctx - optional user-defined function context

2321:    Level: intermediate

2323: .keywords: SNES, nonlinear, set, function

2325: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
2326: @*/
2327: PetscErrorCode  SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
2328: {
2331:   if (func) snes->ops->computeinitialguess = func;
2332:   if (ctx)  snes->initialguessP            = ctx;
2333:   return(0);
2334: }

2336: /* --------------------------------------------------------------- */
2337: /*@C
2338:    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
2339:    it assumes a zero right hand side.

2341:    Logically Collective on SNES

2343:    Input Parameter:
2344: .  snes - the SNES context

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

2349:    Level: intermediate

2351: .keywords: SNES, nonlinear, get, function, right hand side

2353: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
2354: @*/
2355: PetscErrorCode  SNESGetRhs(SNES snes,Vec *rhs)
2356: {
2360:   *rhs = snes->vec_rhs;
2361:   return(0);
2362: }

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

2367:    Collective on SNES

2369:    Input Parameters:
2370: +  snes - the SNES context
2371: -  x - input vector

2373:    Output Parameter:
2374: .  y - function vector, as set by SNESSetFunction()

2376:    Notes:
2377:    SNESComputeFunction() is typically used within nonlinear solvers
2378:    implementations, so most users would not generally call this routine
2379:    themselves.

2381:    Level: developer

2383: .keywords: SNES, nonlinear, compute, function

2385: .seealso: SNESSetFunction(), SNESGetFunction()
2386: @*/
2387: PetscErrorCode  SNESComputeFunction(SNES snes,Vec x,Vec y)
2388: {
2390:   DM             dm;
2391:   DMSNES         sdm;

2399:   VecValidValues(x,2,PETSC_TRUE);

2401:   SNESGetDM(snes,&dm);
2402:   DMGetDMSNES(dm,&sdm);
2403:   if (sdm->ops->computefunction) {
2404:     if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2405:       PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2406:     }
2407:     VecLockReadPush(x);
2408:     PetscStackPush("SNES user function");
2409:     (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2410:     PetscStackPop;
2411:     VecLockReadPop(x);
2412:     if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2413:       PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2414:     }
2415:   } else if (snes->vec_rhs) {
2416:     MatMult(snes->jacobian, x, y);
2417:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2418:   if (snes->vec_rhs) {
2419:     VecAXPY(y,-1.0,snes->vec_rhs);
2420:   }
2421:   snes->nfuncs++;
2422:   /*
2423:      domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2424:      propagate the value to all processes
2425:   */
2426:   if (snes->domainerror) {
2427:     VecSetInf(y);
2428:   }
2429:   return(0);
2430: }

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

2435:    Collective on SNES

2437:    Input Parameters:
2438: +  snes - the SNES context
2439: .  x - input vector
2440: -  b - rhs vector

2442:    Output Parameter:
2443: .  x - new solution vector

2445:    Notes:
2446:    SNESComputeNGS() is typically used within composed nonlinear solver
2447:    implementations, so most users would not generally call this routine
2448:    themselves.

2450:    Level: developer

2452: .keywords: SNES, nonlinear, compute, function

2454: .seealso: SNESSetNGS(), SNESComputeFunction()
2455: @*/
2456: PetscErrorCode  SNESComputeNGS(SNES snes,Vec b,Vec x)
2457: {
2459:   DM             dm;
2460:   DMSNES         sdm;

2468:   if (b) {VecValidValues(b,2,PETSC_TRUE);}
2469:   PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2470:   SNESGetDM(snes,&dm);
2471:   DMGetDMSNES(dm,&sdm);
2472:   if (sdm->ops->computegs) {
2473:     if (b) {VecLockReadPush(b);}
2474:     PetscStackPush("SNES user NGS");
2475:     (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2476:     PetscStackPop;
2477:     if (b) {VecLockReadPop(b);}
2478:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2479:   PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2480:   return(0);
2481: }

2483: PetscErrorCode SNESTestJacobian(SNES snes)
2484: {
2485:   Mat               A,B,C,D,jacobian;
2486:   Vec               x = snes->vec_sol,f = snes->vec_func;
2487:   PetscErrorCode    ierr;
2488:   PetscReal         nrm,gnorm;
2489:   PetscReal         threshold = 1.e-5;
2490:   PetscInt          m,n,M,N;
2491:   void              *functx;
2492:   PetscBool         complete_print = PETSC_FALSE,threshold_print = PETSC_FALSE,test = PETSC_FALSE,flg;
2493:   PetscViewer       viewer,mviewer;
2494:   MPI_Comm          comm;
2495:   PetscInt          tabs;
2496:   static PetscBool  directionsprinted = PETSC_FALSE;
2497:   PetscViewerFormat format;

2500:   PetscObjectOptionsBegin((PetscObject)snes);
2501:   PetscOptionsName("-snes_test_jacobian","Compare hand-coded and finite difference Jacobians","None",&test);
2502:   PetscOptionsReal("-snes_test_jacobian", "Threshold for element difference between hand-coded and finite difference being meaningful", "None", threshold, &threshold,NULL);
2503:   PetscOptionsViewer("-snes_test_jacobian_view","View difference between hand-coded and finite difference Jacobians element entries","None",&mviewer,&format,&complete_print);
2504:   if (!complete_print) {
2505:     PetscOptionsViewer("-snes_test_jacobian_display","Display difference between hand-coded and finite difference Jacobians","None",&mviewer,&format,&complete_print);
2506:   }
2507:   /* for compatibility with PETSc 3.9 and older. */
2508:   PetscOptionsReal("-snes_test_jacobian_display_threshold", "Display difference between hand-coded and finite difference Jacobians which exceed input threshold", "None", threshold, &threshold, &threshold_print);
2509:   PetscOptionsEnd();
2510:   if (!test) return(0);

2512:   PetscObjectGetComm((PetscObject)snes,&comm);
2513:   PetscViewerASCIIGetStdout(comm,&viewer);
2514:   PetscViewerASCIIGetTab(viewer, &tabs);
2515:   PetscViewerASCIISetTab(viewer, ((PetscObject)snes)->tablevel);
2516:   PetscViewerASCIIPrintf(viewer,"  ---------- Testing Jacobian -------------\n");
2517:   if (!complete_print && !directionsprinted) {
2518:     PetscViewerASCIIPrintf(viewer,"  Run with -snes_test_jacobian_view and optionally -snes_test_jacobian <threshold> to show difference\n");
2519:     PetscViewerASCIIPrintf(viewer,"    of hand-coded and finite difference Jacobian entries greater than <threshold>.\n");
2520:   }
2521:   if (!directionsprinted) {
2522:     PetscViewerASCIIPrintf(viewer,"  Testing hand-coded Jacobian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");
2523:     PetscViewerASCIIPrintf(viewer,"    O(1.e-8), the hand-coded Jacobian is probably correct.\n");
2524:     directionsprinted = PETSC_TRUE;
2525:   }
2526:   if (complete_print) {
2527:     PetscViewerPushFormat(mviewer,format);
2528:   }

2530:   /* evaluate the function at this point because SNESComputeJacobianDefault() assumes that the function has been evaluated and put into snes->vec_func */
2531:   SNESComputeFunction(snes,x,f);

2533:   PetscObjectTypeCompare((PetscObject)snes->jacobian,MATMFFD,&flg);
2534:   if (!flg) jacobian = snes->jacobian;
2535:   else jacobian = snes->jacobian_pre;

2537:   while (jacobian) {
2538:     PetscObjectTypeCompareAny((PetscObject)jacobian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2539:     if (flg) {
2540:       A    = jacobian;
2541:       PetscObjectReference((PetscObject)A);
2542:     } else {
2543:       MatComputeOperator(jacobian,MATAIJ,&A);
2544:     }

2546:     MatCreate(PetscObjectComm((PetscObject)A),&B);
2547:     MatGetSize(A,&M,&N);
2548:     MatGetLocalSize(A,&m,&n);
2549:     MatSetSizes(B,m,n,M,N);
2550:     MatSetType(B,((PetscObject)A)->type_name);
2551:     MatSetUp(B);
2552:     MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);

2554:     SNESGetFunction(snes,NULL,NULL,&functx);
2555:     SNESComputeJacobianDefault(snes,x,B,B,functx);

2557:     MatDuplicate(B,MAT_COPY_VALUES,&D);
2558:     MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2559:     MatNorm(D,NORM_FROBENIUS,&nrm);
2560:     MatNorm(A,NORM_FROBENIUS,&gnorm);
2561:     MatDestroy(&D);
2562:     if (!gnorm) gnorm = 1; /* just in case */
2563:     PetscViewerASCIIPrintf(viewer,"  ||J - Jfd||_F/||J||_F = %g, ||J - Jfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);

2565:     if (complete_print) {
2566:       PetscViewerASCIIPrintf(viewer,"  Hand-coded Jacobian ----------\n");
2567:       MatView(jacobian,mviewer);
2568:       PetscViewerASCIIPrintf(viewer,"  Finite difference Jacobian ----------\n");
2569:       MatView(B,mviewer);
2570:     }

2572:     if (threshold_print || complete_print) {
2573:       PetscInt          Istart, Iend, *ccols, bncols, cncols, j, row;
2574:       PetscScalar       *cvals;
2575:       const PetscInt    *bcols;
2576:       const PetscScalar *bvals;

2578:       MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2579:       MatCreate(PetscObjectComm((PetscObject)A),&C);
2580:       MatSetSizes(C,m,n,M,N);
2581:       MatSetType(C,((PetscObject)A)->type_name);
2582:       MatSetUp(C);
2583:       MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
2584:       MatGetOwnershipRange(B,&Istart,&Iend);

2586:       for (row = Istart; row < Iend; row++) {
2587:         MatGetRow(B,row,&bncols,&bcols,&bvals);
2588:         PetscMalloc2(bncols,&ccols,bncols,&cvals);
2589:         for (j = 0, cncols = 0; j < bncols; j++) {
2590:           if (PetscAbsScalar(bvals[j]) > threshold) {
2591:             ccols[cncols] = bcols[j];
2592:             cvals[cncols] = bvals[j];
2593:             cncols += 1;
2594:           }
2595:         }
2596:         if (cncols) {
2597:           MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);
2598:         }
2599:         MatRestoreRow(B,row,&bncols,&bcols,&bvals);
2600:         PetscFree2(ccols,cvals);
2601:       }
2602:       MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2603:       MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2604:       PetscViewerASCIIPrintf(viewer,"  Hand-coded minus finite-difference Jacobian with tolerance %g ----------\n",(double)threshold);
2605:       MatView(C,complete_print ? mviewer : viewer);
2606:       MatDestroy(&C);
2607:     }
2608:     MatDestroy(&A);
2609:     MatDestroy(&B);

2611:     if (jacobian != snes->jacobian_pre) {
2612:       jacobian = snes->jacobian_pre;
2613:       PetscViewerASCIIPrintf(viewer,"  ---------- Testing Jacobian for preconditioner -------------\n");
2614:     }
2615:     else jacobian = NULL;
2616:   }
2617:   if (complete_print) {
2618:     PetscViewerPopFormat(mviewer);
2619:   }
2620:   if (mviewer) { PetscViewerDestroy(&mviewer); }
2621:   PetscViewerASCIISetTab(viewer,tabs);
2622:   return(0);
2623: }

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

2628:    Collective on SNES and Mat

2630:    Input Parameters:
2631: +  snes - the SNES context
2632: -  x - input vector

2634:    Output Parameters:
2635: +  A - Jacobian matrix
2636: -  B - optional preconditioning matrix

2638:   Options Database Keys:
2639: +    -snes_lag_preconditioner <lag>
2640: .    -snes_lag_jacobian <lag>
2641: .    -snes_test_jacobian - compare the user provided Jacobian with one compute via finite differences to check for errors
2642: .    -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
2643: .    -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
2644: .    -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2645: .    -snes_compare_explicit_draw  - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2646: .    -snes_compare_explicit_contour  - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2647: .    -snes_compare_operator  - Make the comparison options above use the operator instead of the preconditioning matrix
2648: .    -snes_compare_coloring - Compute the finite difference Jacobian using coloring and display norms of difference
2649: .    -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
2650: .    -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2651: .    -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2652: .    -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2653: .    -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
2654: -    -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences


2657:    Notes:
2658:    Most users should not need to explicitly call this routine, as it
2659:    is used internally within the nonlinear solvers.

2661:    Developer Notes:
2662:     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
2663:       for with the SNESType of test that has been removed.

2665:    Level: developer

2667: .keywords: SNES, compute, Jacobian, matrix

2669: .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2670: @*/
2671: PetscErrorCode  SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2672: {
2674:   PetscBool      flag;
2675:   DM             dm;
2676:   DMSNES         sdm;
2677:   KSP            ksp;

2683:   VecValidValues(X,2,PETSC_TRUE);
2684:   SNESGetDM(snes,&dm);
2685:   DMGetDMSNES(dm,&sdm);

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

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

2691:   if (snes->lagjacobian == -2) {
2692:     snes->lagjacobian = -1;

2694:     PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2695:   } else if (snes->lagjacobian == -1) {
2696:     PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2697:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2698:     if (flag) {
2699:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2700:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2701:     }
2702:     return(0);
2703:   } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2704:     PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2705:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2706:     if (flag) {
2707:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2708:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2709:     }
2710:     return(0);
2711:   }
2712:   if (snes->npc && snes->npcside== PC_LEFT) {
2713:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2714:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2715:       return(0);
2716:   }

2718:   PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2719:   VecLockReadPush(X);
2720:   PetscStackPush("SNES user Jacobian function");
2721:   (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2722:   PetscStackPop;
2723:   VecLockReadPop(X);
2724:   PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);

2726:   /* the next line ensures that snes->ksp exists */
2727:   SNESGetKSP(snes,&ksp);
2728:   if (snes->lagpreconditioner == -2) {
2729:     PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2730:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2731:     snes->lagpreconditioner = -1;
2732:   } else if (snes->lagpreconditioner == -1) {
2733:     PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2734:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2735:   } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2736:     PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2737:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2738:   } else {
2739:     PetscInfo(snes,"Rebuilding preconditioner\n");
2740:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2741:   }

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

2819:       MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2820:       MatColoringCreate(Bfd,&coloring);
2821:       MatColoringSetType(coloring,MATCOLORINGSL);
2822:       MatColoringSetFromOptions(coloring);
2823:       MatColoringApply(coloring,&iscoloring);
2824:       MatColoringDestroy(&coloring);
2825:       MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2826:       MatFDColoringSetFromOptions(matfdcoloring);
2827:       MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2828:       ISColoringDestroy(&iscoloring);

2830:       /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2831:       SNESGetFunction(snes,NULL,&func,&funcctx);
2832:       MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2833:       PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2834:       PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2835:       MatFDColoringSetFromOptions(matfdcoloring);
2836:       MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2837:       MatFDColoringDestroy(&matfdcoloring);

2839:       PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2840:       if (flag_draw || flag_contour) {
2841:         PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2842:         if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2843:       } else vdraw = NULL;
2844:       PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2845:       if (flag_display) {MatView(B,vstdout);}
2846:       if (vdraw) {MatView(B,vdraw);}
2847:       PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2848:       if (flag_display) {MatView(Bfd,vstdout);}
2849:       if (vdraw) {MatView(Bfd,vdraw);}
2850:       MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2851:       MatNorm(Bfd,NORM_1,&norm1);
2852:       MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2853:       MatNorm(Bfd,NORM_MAX,&normmax);
2854:       PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2855:       if (flag_display) {MatView(Bfd,vstdout);}
2856:       if (vdraw) {              /* Always use contour for the difference */
2857:         PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2858:         MatView(Bfd,vdraw);
2859:         PetscViewerPopFormat(vdraw);
2860:       }
2861:       if (flag_contour) {PetscViewerPopFormat(vdraw);}

2863:       if (flag_threshold) {
2864:         PetscInt bs,rstart,rend,i;
2865:         MatGetBlockSize(B,&bs);
2866:         MatGetOwnershipRange(B,&rstart,&rend);
2867:         for (i=rstart; i<rend; i++) {
2868:           const PetscScalar *ba,*ca;
2869:           const PetscInt    *bj,*cj;
2870:           PetscInt          bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
2871:           PetscReal         maxentry = 0,maxdiff = 0,maxrdiff = 0;
2872:           MatGetRow(B,i,&bn,&bj,&ba);
2873:           MatGetRow(Bfd,i,&cn,&cj,&ca);
2874:           if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
2875:           for (j=0; j<bn; j++) {
2876:             PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2877:             if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
2878:               maxentrycol = bj[j];
2879:               maxentry    = PetscRealPart(ba[j]);
2880:             }
2881:             if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
2882:               maxdiffcol = bj[j];
2883:               maxdiff    = PetscRealPart(ca[j]);
2884:             }
2885:             if (rdiff > maxrdiff) {
2886:               maxrdiffcol = bj[j];
2887:               maxrdiff    = rdiff;
2888:             }
2889:           }
2890:           if (maxrdiff > 1) {
2891:             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);
2892:             for (j=0; j<bn; j++) {
2893:               PetscReal rdiff;
2894:               rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2895:               if (rdiff > 1) {
2896:                 PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
2897:               }
2898:             }
2899:             PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
2900:           }
2901:           MatRestoreRow(B,i,&bn,&bj,&ba);
2902:           MatRestoreRow(Bfd,i,&cn,&cj,&ca);
2903:         }
2904:       }
2905:       PetscViewerDestroy(&vdraw);
2906:       MatDestroy(&Bfd);
2907:     }
2908:   }
2909:   return(0);
2910: }

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

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

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

2924:    Level: intermediate

2926: .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2927: M*/

2929: /*@C
2930:    SNESSetJacobian - Sets the function to compute Jacobian as well as the
2931:    location to store the matrix.

2933:    Logically Collective on SNES and Mat

2935:    Input Parameters:
2936: +  snes - the SNES context
2937: .  Amat - the matrix that defines the (approximate) Jacobian
2938: .  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2939: .  J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
2940: -  ctx - [optional] user-defined context for private data for the
2941:          Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)

2943:    Notes:
2944:    If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2945:    each matrix.

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

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

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

2956:    Level: beginner

2958: .keywords: SNES, nonlinear, set, Jacobian, matrix

2960: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J,
2961:           SNESSetPicard(), SNESJacobianFunction
2962: @*/
2963: PetscErrorCode  SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2964: {
2966:   DM             dm;

2974:   SNESGetDM(snes,&dm);
2975:   DMSNESSetJacobian(dm,J,ctx);
2976:   if (Amat) {
2977:     PetscObjectReference((PetscObject)Amat);
2978:     MatDestroy(&snes->jacobian);

2980:     snes->jacobian = Amat;
2981:   }
2982:   if (Pmat) {
2983:     PetscObjectReference((PetscObject)Pmat);
2984:     MatDestroy(&snes->jacobian_pre);

2986:     snes->jacobian_pre = Pmat;
2987:   }
2988:   return(0);
2989: }

2991: /*@C
2992:    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2993:    provided context for evaluating the Jacobian.

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

2997:    Input Parameter:
2998: .  snes - the nonlinear solver context

3000:    Output Parameters:
3001: +  Amat - location to stash (approximate) Jacobian matrix (or NULL)
3002: .  Pmat - location to stash matrix used to compute the preconditioner (or NULL)
3003: .  J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
3004: -  ctx - location to stash Jacobian ctx (or NULL)

3006:    Level: advanced

3008: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
3009: @*/
3010: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
3011: {
3013:   DM             dm;
3014:   DMSNES         sdm;

3018:   if (Amat) *Amat = snes->jacobian;
3019:   if (Pmat) *Pmat = snes->jacobian_pre;
3020:   SNESGetDM(snes,&dm);
3021:   DMGetDMSNES(dm,&sdm);
3022:   if (J) *J = sdm->ops->computejacobian;
3023:   if (ctx) *ctx = sdm->jacobianctx;
3024:   return(0);
3025: }

3027: /*@
3028:    SNESSetUp - Sets up the internal data structures for the later use
3029:    of a nonlinear solver.

3031:    Collective on SNES

3033:    Input Parameters:
3034: .  snes - the SNES context

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

3043:    Level: advanced

3045: .keywords: SNES, nonlinear, setup

3047: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
3048: @*/
3049: PetscErrorCode  SNESSetUp(SNES snes)
3050: {
3052:   DM             dm;
3053:   DMSNES         sdm;
3054:   SNESLineSearch linesearch, pclinesearch;
3055:   void           *lsprectx,*lspostctx;
3056:   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
3057:   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
3058:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
3059:   Vec            f,fpc;
3060:   void           *funcctx;
3061:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
3062:   void           *jacctx,*appctx;
3063:   Mat            j,jpre;

3067:   if (snes->setupcalled) return(0);
3068:   PetscLogEventBegin(SNES_Setup,snes,0,0,0);

3070:   if (!((PetscObject)snes)->type_name) {
3071:     SNESSetType(snes,SNESNEWTONLS);
3072:   }

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

3076:   SNESGetDM(snes,&dm);
3077:   DMGetDMSNES(dm,&sdm);
3078:   if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
3079:   if (!sdm->ops->computejacobian) {
3080:     DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
3081:   }
3082:   if (!snes->vec_func) {
3083:     DMCreateGlobalVector(dm,&snes->vec_func);
3084:   }

3086:   if (!snes->ksp) {
3087:     SNESGetKSP(snes, &snes->ksp);
3088:   }

3090:   if (snes->linesearch) {
3091:     SNESGetLineSearch(snes, &snes->linesearch);
3092:     SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);
3093:   }

3095:   if (snes->npc && (snes->npcside== PC_LEFT)) {
3096:     snes->mf          = PETSC_TRUE;
3097:     snes->mf_operator = PETSC_FALSE;
3098:   }

3100:   if (snes->npc) {
3101:     /* copy the DM over */
3102:     SNESGetDM(snes,&dm);
3103:     SNESSetDM(snes->npc,dm);

3105:     SNESGetFunction(snes,&f,&func,&funcctx);
3106:     VecDuplicate(f,&fpc);
3107:     SNESSetFunction(snes->npc,fpc,func,funcctx);
3108:     SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
3109:     SNESSetJacobian(snes->npc,j,jpre,jac,jacctx);
3110:     SNESGetApplicationContext(snes,&appctx);
3111:     SNESSetApplicationContext(snes->npc,appctx);
3112:     VecDestroy(&fpc);

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

3117:     /* default to 1 iteration */
3118:     SNESSetTolerances(snes->npc,0.0,0.0,0.0,1,snes->npc->max_funcs);
3119:     if (snes->npcside==PC_RIGHT) {
3120:       SNESSetNormSchedule(snes->npc,SNES_NORM_FINAL_ONLY);
3121:     } else {
3122:       SNESSetNormSchedule(snes->npc,SNES_NORM_NONE);
3123:     }
3124:     SNESSetFromOptions(snes->npc);

3126:     /* copy the line search context over */
3127:     if (snes->linesearch && snes->npc->linesearch) {
3128:       SNESGetLineSearch(snes,&linesearch);
3129:       SNESGetLineSearch(snes->npc,&pclinesearch);
3130:       SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
3131:       SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
3132:       SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
3133:       SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
3134:       PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
3135:     }
3136:   }
3137:   if (snes->mf) {
3138:     SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
3139:   }
3140:   if (snes->ops->usercompute && !snes->user) {
3141:     (*snes->ops->usercompute)(snes,(void**)&snes->user);
3142:   }

3144:   snes->jac_iter = 0;
3145:   snes->pre_iter = 0;

3147:   if (snes->ops->setup) {
3148:     (*snes->ops->setup)(snes);
3149:   }

3151:   if (snes->npc && (snes->npcside== PC_LEFT)) {
3152:     if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3153:       if (snes->linesearch){
3154:         SNESGetLineSearch(snes,&linesearch);
3155:         SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
3156:       }
3157:     }
3158:   }
3159:   PetscLogEventEnd(SNES_Setup,snes,0,0,0);
3160:   snes->setupcalled = PETSC_TRUE;
3161:   return(0);
3162: }

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

3167:    Collective on SNES

3169:    Input Parameter:
3170: .  snes - iterative context obtained from SNESCreate()

3172:    Level: intermediate

3174:    Notes:
3175:     Also calls the application context destroy routine set with SNESSetComputeApplicationContext()

3177: .keywords: SNES, destroy

3179: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
3180: @*/
3181: PetscErrorCode  SNESReset(SNES snes)
3182: {

3187:   if (snes->ops->userdestroy && snes->user) {
3188:     (*snes->ops->userdestroy)((void**)&snes->user);
3189:     snes->user = NULL;
3190:   }
3191:   if (snes->npc) {
3192:     SNESReset(snes->npc);
3193:   }

3195:   if (snes->ops->reset) {
3196:     (*snes->ops->reset)(snes);
3197:   }
3198:   if (snes->ksp) {
3199:     KSPReset(snes->ksp);
3200:   }

3202:   if (snes->linesearch) {
3203:     SNESLineSearchReset(snes->linesearch);
3204:   }

3206:   VecDestroy(&snes->vec_rhs);
3207:   VecDestroy(&snes->vec_sol);
3208:   VecDestroy(&snes->vec_sol_update);
3209:   VecDestroy(&snes->vec_func);
3210:   MatDestroy(&snes->jacobian);
3211:   MatDestroy(&snes->jacobian_pre);
3212:   VecDestroyVecs(snes->nwork,&snes->work);
3213:   VecDestroyVecs(snes->nvwork,&snes->vwork);

3215:   snes->alwayscomputesfinalresidual = PETSC_FALSE;

3217:   snes->nwork       = snes->nvwork = 0;
3218:   snes->setupcalled = PETSC_FALSE;
3219:   return(0);
3220: }

3222: /*@
3223:    SNESDestroy - Destroys the nonlinear solver context that was created
3224:    with SNESCreate().

3226:    Collective on SNES

3228:    Input Parameter:
3229: .  snes - the SNES context

3231:    Level: beginner

3233: .keywords: SNES, nonlinear, destroy

3235: .seealso: SNESCreate(), SNESSolve()
3236: @*/
3237: PetscErrorCode  SNESDestroy(SNES *snes)
3238: {

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

3246:   SNESReset((*snes));
3247:   SNESDestroy(&(*snes)->npc);

3249:   /* if memory was published with SAWs then destroy it */
3250:   PetscObjectSAWsViewOff((PetscObject)*snes);
3251:   if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}

3253:   if ((*snes)->dm) {DMCoarsenHookRemove((*snes)->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,*snes);}
3254:   DMDestroy(&(*snes)->dm);
3255:   KSPDestroy(&(*snes)->ksp);
3256:   SNESLineSearchDestroy(&(*snes)->linesearch);

3258:   PetscFree((*snes)->kspconvctx);
3259:   if ((*snes)->ops->convergeddestroy) {
3260:     (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
3261:   }
3262:   if ((*snes)->conv_malloc) {
3263:     PetscFree((*snes)->conv_hist);
3264:     PetscFree((*snes)->conv_hist_its);
3265:   }
3266:   SNESMonitorCancel((*snes));
3267:   PetscHeaderDestroy(snes);
3268:   return(0);
3269: }

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

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

3276:    Logically Collective on SNES

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

3283:    Options Database Keys:
3284: .    -snes_lag_preconditioner <lag>

3286:    Notes:
3287:    The default is 1
3288:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3289:    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use

3291:    Level: intermediate

3293: .keywords: SNES, nonlinear, set, convergence, tolerances

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

3297: @*/
3298: PetscErrorCode  SNESSetLagPreconditioner(SNES snes,PetscInt lag)
3299: {
3302:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3303:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3305:   snes->lagpreconditioner = lag;
3306:   return(0);
3307: }

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

3312:    Logically Collective on SNES

3314:    Input Parameters:
3315: +  snes - the SNES context
3316: -  steps - the number of refinements to do, defaults to 0

3318:    Options Database Keys:
3319: .    -snes_grid_sequence <steps>

3321:    Level: intermediate

3323:    Notes:
3324:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

3326: .keywords: SNES, nonlinear, set, convergence, tolerances

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

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

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

3343:    Logically Collective on SNES

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

3348:    Output Parameter:
3349: .  steps - the number of refinements to do, defaults to 0

3351:    Options Database Keys:
3352: .    -snes_grid_sequence <steps>

3354:    Level: intermediate

3356:    Notes:
3357:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

3359: .keywords: SNES, nonlinear, set, convergence, tolerances

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

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

3372: /*@
3373:    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt

3375:    Not Collective

3377:    Input Parameter:
3378: .  snes - the SNES context

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

3384:    Options Database Keys:
3385: .    -snes_lag_preconditioner <lag>

3387:    Notes:
3388:    The default is 1
3389:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

3391:    Level: intermediate

3393: .keywords: SNES, nonlinear, set, convergence, tolerances

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

3397: @*/
3398: PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
3399: {
3402:   *lag = snes->lagpreconditioner;
3403:   return(0);
3404: }

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

3410:    Logically Collective on SNES

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

3417:    Options Database Keys:
3418: .    -snes_lag_jacobian <lag>

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

3426:    Level: intermediate

3428: .keywords: SNES, nonlinear, set, convergence, tolerances

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

3432: @*/
3433: PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
3434: {
3437:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3438:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3440:   snes->lagjacobian = lag;
3441:   return(0);
3442: }

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

3447:    Not Collective

3449:    Input Parameter:
3450: .  snes - the SNES context

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

3456:    Options Database Keys:
3457: .    -snes_lag_jacobian <lag>

3459:    Notes:
3460:    The default is 1
3461:    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

3463:    Level: intermediate

3465: .keywords: SNES, nonlinear, set, convergence, tolerances

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

3469: @*/
3470: PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
3471: {
3474:   *lag = snes->lagjacobian;
3475:   return(0);
3476: }

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

3481:    Logically collective on SNES

3483:    Input Parameter:
3484: +  snes - the SNES context
3485: -   flg - jacobian lagging persists if true

3487:    Options Database Keys:
3488: .    -snes_lag_jacobian_persists <flg>

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

3495:    Level: developer

3497: .keywords: SNES, nonlinear, lag

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

3501: @*/
3502: PetscErrorCode  SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
3503: {
3507:   snes->lagjac_persist = flg;
3508:   return(0);
3509: }

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

3514:    Logically Collective on SNES

3516:    Input Parameter:
3517: +  snes - the SNES context
3518: -   flg - preconditioner lagging persists if true

3520:    Options Database Keys:
3521: .    -snes_lag_jacobian_persists <flg>

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

3528:    Level: developer

3530: .keywords: SNES, nonlinear, lag

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

3534: @*/
3535: PetscErrorCode  SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3536: {
3540:   snes->lagpre_persist = flg;
3541:   return(0);
3542: }

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

3547:    Logically Collective on SNES

3549:    Input Parameters:
3550: +  snes - the SNES context
3551: -  force - PETSC_TRUE require at least one iteration

3553:    Options Database Keys:
3554: .    -snes_force_iteration <force> - Sets forcing an iteration

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

3559:    Level: intermediate

3561: .keywords: SNES, nonlinear, set, convergence, tolerances

3563: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3564: @*/
3565: PetscErrorCode  SNESSetForceIteration(SNES snes,PetscBool force)
3566: {
3569:   snes->forceiteration = force;
3570:   return(0);
3571: }

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

3576:    Logically Collective on SNES

3578:    Input Parameters:
3579: .  snes - the SNES context

3581:    Output Parameter:
3582: .  force - PETSC_TRUE requires at least one iteration.

3584: .keywords: SNES, nonlinear, set, convergence, tolerances

3586:    Level: intermediate

3588: .seealso: SNESSetForceIteration(), SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3589: @*/
3590: PetscErrorCode  SNESGetForceIteration(SNES snes,PetscBool *force)
3591: {
3594:   *force = snes->forceiteration;
3595:   return(0);
3596: }

3598: /*@
3599:    SNESSetTolerances - Sets various parameters used in convergence tests.

3601:    Logically Collective on SNES

3603:    Input Parameters:
3604: +  snes - the SNES context
3605: .  abstol - absolute convergence tolerance
3606: .  rtol - relative convergence tolerance
3607: .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
3608: .  maxit - maximum number of iterations
3609: -  maxf - maximum number of function evaluations (-1 indicates no limit)

3611:    Options Database Keys:
3612: +    -snes_atol <abstol> - Sets abstol
3613: .    -snes_rtol <rtol> - Sets rtol
3614: .    -snes_stol <stol> - Sets stol
3615: .    -snes_max_it <maxit> - Sets maxit
3616: -    -snes_max_funcs <maxf> - Sets maxf

3618:    Notes:
3619:    The default maximum number of iterations is 50.
3620:    The default maximum number of function evaluations is 1000.

3622:    Level: intermediate

3624: .keywords: SNES, nonlinear, set, convergence, tolerances

3626: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance(), SNESSetForceIteration()
3627: @*/
3628: PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3629: {

3638:   if (abstol != PETSC_DEFAULT) {
3639:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3640:     snes->abstol = abstol;
3641:   }
3642:   if (rtol != PETSC_DEFAULT) {
3643:     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);
3644:     snes->rtol = rtol;
3645:   }
3646:   if (stol != PETSC_DEFAULT) {
3647:     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3648:     snes->stol = stol;
3649:   }
3650:   if (maxit != PETSC_DEFAULT) {
3651:     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3652:     snes->max_its = maxit;
3653:   }
3654:   if (maxf != PETSC_DEFAULT) {
3655:     if (maxf < -1) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be -1 or nonnegative",maxf);
3656:     snes->max_funcs = maxf;
3657:   }
3658:   snes->tolerancesset = PETSC_TRUE;
3659:   return(0);
3660: }

3662: /*@
3663:    SNESSetDivergenceTolerance - Sets the divergence tolerance used for the SNES divergence test.

3665:    Logically Collective on SNES

3667:    Input Parameters:
3668: +  snes - the SNES context
3669: -  divtol - the divergence tolerance. Use -1 to deactivate the test.

3671:    Options Database Keys:
3672: +    -snes_divergence_tolerance <divtol> - Sets divtol

3674:    Notes:
3675:    The default divergence tolerance is 1e4.

3677:    Level: intermediate

3679: .keywords: SNES, nonlinear, set, divergence, tolerance

3681: .seealso: SNESSetTolerances(), SNESGetDivergenceTolerance
3682: @*/
3683: PetscErrorCode  SNESSetDivergenceTolerance(SNES snes,PetscReal divtol)
3684: {

3689:   if (divtol != PETSC_DEFAULT) {
3690:     snes->divtol = divtol;
3691:   }
3692:   else {
3693:     snes->divtol = 1.0e4;
3694:   }
3695:   return(0);
3696: }

3698: /*@
3699:    SNESGetTolerances - Gets various parameters used in convergence tests.

3701:    Not Collective

3703:    Input Parameters:
3704: +  snes - the SNES context
3705: .  atol - absolute convergence tolerance
3706: .  rtol - relative convergence tolerance
3707: .  stol -  convergence tolerance in terms of the norm
3708:            of the change in the solution between steps
3709: .  maxit - maximum number of iterations
3710: -  maxf - maximum number of function evaluations

3712:    Notes:
3713:    The user can specify NULL for any parameter that is not needed.

3715:    Level: intermediate

3717: .keywords: SNES, nonlinear, get, convergence, tolerances

3719: .seealso: SNESSetTolerances()
3720: @*/
3721: PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3722: {
3725:   if (atol)  *atol  = snes->abstol;
3726:   if (rtol)  *rtol  = snes->rtol;
3727:   if (stol)  *stol  = snes->stol;
3728:   if (maxit) *maxit = snes->max_its;
3729:   if (maxf)  *maxf  = snes->max_funcs;
3730:   return(0);
3731: }

3733: /*@
3734:    SNESGetDivergenceTolerance - Gets divergence tolerance used in divergence test.

3736:    Not Collective

3738:    Input Parameters:
3739: +  snes - the SNES context
3740: -  divtol - divergence tolerance

3742:    Level: intermediate

3744: .keywords: SNES, nonlinear, get, divergence, tolerance

3746: .seealso: SNESSetDivergenceTolerance()
3747: @*/
3748: PetscErrorCode  SNESGetDivergenceTolerance(SNES snes,PetscReal *divtol)
3749: {
3752:   if (divtol) *divtol = snes->divtol;
3753:   return(0);
3754: }

3756: /*@
3757:    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.

3759:    Logically Collective on SNES

3761:    Input Parameters:
3762: +  snes - the SNES context
3763: -  tol - tolerance

3765:    Options Database Key:
3766: .  -snes_trtol <tol> - Sets tol

3768:    Level: intermediate

3770: .keywords: SNES, nonlinear, set, trust region, tolerance

3772: .seealso: SNESSetTolerances()
3773: @*/
3774: PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3775: {
3779:   snes->deltatol = tol;
3780:   return(0);
3781: }

3783: /*
3784:    Duplicate the lg monitors for SNES from KSP; for some reason with
3785:    dynamic libraries things don't work under Sun4 if we just use
3786:    macros instead of functions
3787: */
3788: PetscErrorCode  SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,void *ctx)
3789: {

3794:   KSPMonitorLGResidualNorm((KSP)snes,it,norm,ctx);
3795:   return(0);
3796: }

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

3803:   KSPMonitorLGResidualNormCreate(comm,host,label,x,y,m,n,lgctx);
3804:   return(0);
3805: }

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

3809: PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3810: {
3811:   PetscDrawLG      lg;
3812:   PetscErrorCode   ierr;
3813:   PetscReal        x,y,per;
3814:   PetscViewer      v = (PetscViewer)monctx;
3815:   static PetscReal prev; /* should be in the context */
3816:   PetscDraw        draw;

3820:   PetscViewerDrawGetDrawLG(v,0,&lg);
3821:   if (!n) {PetscDrawLGReset(lg);}
3822:   PetscDrawLGGetDraw(lg,&draw);
3823:   PetscDrawSetTitle(draw,"Residual norm");
3824:   x    = (PetscReal)n;
3825:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3826:   else y = -15.0;
3827:   PetscDrawLGAddPoint(lg,&x,&y);
3828:   if (n < 20 || !(n % 5) || snes->reason) {
3829:     PetscDrawLGDraw(lg);
3830:     PetscDrawLGSave(lg);
3831:   }

3833:   PetscViewerDrawGetDrawLG(v,1,&lg);
3834:   if (!n) {PetscDrawLGReset(lg);}
3835:   PetscDrawLGGetDraw(lg,&draw);
3836:   PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3837:    SNESMonitorRange_Private(snes,n,&per);
3838:   x    = (PetscReal)n;
3839:   y    = 100.0*per;
3840:   PetscDrawLGAddPoint(lg,&x,&y);
3841:   if (n < 20 || !(n % 5) || snes->reason) {
3842:     PetscDrawLGDraw(lg);
3843:     PetscDrawLGSave(lg);
3844:   }

3846:   PetscViewerDrawGetDrawLG(v,2,&lg);
3847:   if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3848:   PetscDrawLGGetDraw(lg,&draw);
3849:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3850:   x    = (PetscReal)n;
3851:   y    = (prev - rnorm)/prev;
3852:   PetscDrawLGAddPoint(lg,&x,&y);
3853:   if (n < 20 || !(n % 5) || snes->reason) {
3854:     PetscDrawLGDraw(lg);
3855:     PetscDrawLGSave(lg);
3856:   }

3858:   PetscViewerDrawGetDrawLG(v,3,&lg);
3859:   if (!n) {PetscDrawLGReset(lg);}
3860:   PetscDrawLGGetDraw(lg,&draw);
3861:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3862:   x    = (PetscReal)n;
3863:   y    = (prev - rnorm)/(prev*per);
3864:   if (n > 2) { /*skip initial crazy value */
3865:     PetscDrawLGAddPoint(lg,&x,&y);
3866:   }
3867:   if (n < 20 || !(n % 5) || snes->reason) {
3868:     PetscDrawLGDraw(lg);
3869:     PetscDrawLGSave(lg);
3870:   }
3871:   prev = rnorm;
3872:   return(0);
3873: }

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

3878:    Collective on SNES

3880:    Input Parameters:
3881: +  snes - nonlinear solver context obtained from SNESCreate()
3882: .  iter - iteration number
3883: -  rnorm - relative norm of the residual

3885:    Notes:
3886:    This routine is called by the SNES implementations.
3887:    It does not typically need to be called by the user.

3889:    Level: developer

3891: .seealso: SNESMonitorSet()
3892: @*/
3893: PetscErrorCode  SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3894: {
3896:   PetscInt       i,n = snes->numbermonitors;

3899:   VecLockReadPush(snes->vec_sol);
3900:   for (i=0; i<n; i++) {
3901:     (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3902:   }
3903:   VecLockReadPop(snes->vec_sol);
3904:   return(0);
3905: }

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

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

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

3916: +    snes - the SNES context
3917: .    its - iteration number
3918: .    norm - 2-norm function value (may be estimated)
3919: -    mctx - [optional] monitoring context

3921:    Level: advanced

3923: .seealso:   SNESMonitorSet(), SNESMonitorGet()
3924: M*/

3926: /*@C
3927:    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3928:    iteration of the nonlinear solver to display the iteration's
3929:    progress.

3931:    Logically Collective on SNES

3933:    Input Parameters:
3934: +  snes - the SNES context
3935: .  f - the monitor function, see SNESMonitorFunction for the calling sequence
3936: .  mctx - [optional] user-defined context for private data for the
3937:           monitor routine (use NULL if no context is desired)
3938: -  monitordestroy - [optional] routine that frees monitor context
3939:           (may be NULL)

3941:    Options Database Keys:
3942: +    -snes_monitor        - sets SNESMonitorDefault()
3943: .    -snes_monitor_lg_residualnorm    - sets line graph monitor,
3944:                             uses SNESMonitorLGCreate()
3945: -    -snes_monitor_cancel - cancels all monitors that have
3946:                             been hardwired into a code by
3947:                             calls to SNESMonitorSet(), but
3948:                             does not cancel those set via
3949:                             the options database.

3951:    Notes:
3952:    Several different monitoring routines may be set by calling
3953:    SNESMonitorSet() multiple times; all will be called in the
3954:    order in which they were set.

3956:    Fortran Notes:
3957:     Only a single monitor function can be set for each SNES object

3959:    Level: intermediate

3961: .keywords: SNES, nonlinear, set, monitor

3963: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3964: @*/
3965: PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3966: {
3967:   PetscInt       i;
3969:   PetscBool      identical;

3973:   for (i=0; i<snes->numbermonitors;i++) {
3974:     PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))snes->monitor[i],snes->monitorcontext[i],snes->monitordestroy[i],&identical);
3975:     if (identical) return(0);
3976:   }
3977:   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3978:   snes->monitor[snes->numbermonitors]          = f;
3979:   snes->monitordestroy[snes->numbermonitors]   = monitordestroy;
3980:   snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3981:   return(0);
3982: }

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

3987:    Logically Collective on SNES

3989:    Input Parameters:
3990: .  snes - the SNES context

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

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

4000:    Level: intermediate

4002: .keywords: SNES, nonlinear, set, monitor

4004: .seealso: SNESMonitorDefault(), SNESMonitorSet()
4005: @*/
4006: PetscErrorCode  SNESMonitorCancel(SNES snes)
4007: {
4009:   PetscInt       i;

4013:   for (i=0; i<snes->numbermonitors; i++) {
4014:     if (snes->monitordestroy[i]) {
4015:       (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
4016:     }
4017:   }
4018:   snes->numbermonitors = 0;
4019:   return(0);
4020: }

4022: /*MC
4023:     SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver

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

4029: +    snes - the SNES context
4030: .    it - current iteration (0 is the first and is before any Newton step)
4031: .    cctx - [optional] convergence context
4032: .    reason - reason for convergence/divergence
4033: .    xnorm - 2-norm of current iterate
4034: .    gnorm - 2-norm of current step
4035: -    f - 2-norm of function

4037:    Level: intermediate

4039: .seealso:   SNESSetConvergenceTest(), SNESGetConvergenceTest()
4040: M*/

4042: /*@C
4043:    SNESSetConvergenceTest - Sets the function that is to be used
4044:    to test for convergence of the nonlinear iterative solution.

4046:    Logically Collective on SNES

4048:    Input Parameters:
4049: +  snes - the SNES context
4050: .  SNESConvergenceTestFunction - routine to test for convergence
4051: .  cctx - [optional] context for private data for the convergence routine  (may be NULL)
4052: -  destroy - [optional] destructor for the context (may be NULL; PETSC_NULL_FUNCTION in Fortran)

4054:    Level: advanced

4056: .keywords: SNES, nonlinear, set, convergence, test

4058: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
4059: @*/
4060: PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
4061: {

4066:   if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
4067:   if (snes->ops->convergeddestroy) {
4068:     (*snes->ops->convergeddestroy)(snes->cnvP);
4069:   }
4070:   snes->ops->converged        = SNESConvergenceTestFunction;
4071:   snes->ops->convergeddestroy = destroy;
4072:   snes->cnvP                  = cctx;
4073:   return(0);
4074: }

4076: /*@
4077:    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.

4079:    Not Collective

4081:    Input Parameter:
4082: .  snes - the SNES context

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

4088:    Options Database:
4089: .   -snes_converged_reason - prints the reason to standard out

4091:    Level: intermediate

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

4096: .keywords: SNES, nonlinear, set, convergence, test

4098: .seealso: SNESSetConvergenceTest(), SNESSetConvergedReason(), SNESConvergedReason
4099: @*/
4100: PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
4101: {
4105:   *reason = snes->reason;
4106:   return(0);
4107: }

4109: /*@
4110:    SNESSetConvergedReason - Sets the reason the SNES iteration was stopped.

4112:    Not Collective

4114:    Input Parameters:
4115: +  snes - the SNES context
4116: -  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
4117:             manual pages for the individual convergence tests for complete lists

4119:    Level: intermediate

4121: .keywords: SNES, nonlinear, set, convergence, test
4122: .seealso: SNESGetConvergedReason(), SNESSetConvergenceTest(), SNESConvergedReason
4123: @*/
4124: PetscErrorCode SNESSetConvergedReason(SNES snes,SNESConvergedReason reason)
4125: {
4128:   snes->reason = reason;
4129:   return(0);
4130: }

4132: /*@
4133:    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.

4135:    Logically Collective on SNES

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

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

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

4153:    Level: intermediate

4155: .keywords: SNES, set, convergence, history

4157: .seealso: SNESGetConvergenceHistory()

4159: @*/
4160: PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
4161: {

4168:   if (!a) {
4169:     if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
4170:     PetscCalloc1(na,&a);
4171:     PetscCalloc1(na,&its);

4173:     snes->conv_malloc = PETSC_TRUE;
4174:   }
4175:   snes->conv_hist       = a;
4176:   snes->conv_hist_its   = its;
4177:   snes->conv_hist_max   = na;
4178:   snes->conv_hist_len   = 0;
4179:   snes->conv_hist_reset = reset;
4180:   return(0);
4181: }

4183: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4184: #include <engine.h>   /* MATLAB include file */
4185: #include <mex.h>      /* MATLAB include file */

4187: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
4188: {
4189:   mxArray   *mat;
4190:   PetscInt  i;
4191:   PetscReal *ar;

4194:   mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
4195:   ar  = (PetscReal*) mxGetData(mat);
4196:   for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
4197:   PetscFunctionReturn(mat);
4198: }
4199: #endif

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

4204:    Not Collective

4206:    Input Parameter:
4207: .  snes - iterative context obtained from SNESCreate()

4209:    Output Parameters:
4210: .  a   - array to hold history
4211: .  its - integer array holds the number of linear iterations (or
4212:          negative if not converged) for each solve.
4213: -  na  - size of a and its

4215:    Notes:
4216:     The calling sequence for this routine in Fortran is
4217: $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)

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

4223:    Level: intermediate

4225: .keywords: SNES, get, convergence, history

4227: .seealso: SNESSetConvergencHistory()

4229: @*/
4230: PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
4231: {
4234:   if (a)   *a   = snes->conv_hist;
4235:   if (its) *its = snes->conv_hist_its;
4236:   if (na)  *na  = snes->conv_hist_len;
4237:   return(0);
4238: }

4240: /*@C
4241:   SNESSetUpdate - Sets the general-purpose update function called
4242:   at the beginning of every iteration of the nonlinear solve. Specifically
4243:   it is called just before the Jacobian is "evaluated".

4245:   Logically Collective on SNES

4247:   Input Parameters:
4248: . snes - The nonlinear solver context
4249: . func - The function

4251:   Calling sequence of func:
4252: . func (SNES snes, PetscInt step);

4254: . step - The current step of the iteration

4256:   Level: advanced

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

4261: .keywords: SNES, update

4263: .seealso SNESSetJacobian(), SNESSolve()
4264: @*/
4265: PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
4266: {
4269:   snes->ops->update = func;
4270:   return(0);
4271: }

4273: /*
4274:    SNESScaleStep_Private - Scales a step so that its length is less than the
4275:    positive parameter delta.

4277:     Input Parameters:
4278: +   snes - the SNES context
4279: .   y - approximate solution of linear system
4280: .   fnorm - 2-norm of current function
4281: -   delta - trust region size

4283:     Output Parameters:
4284: +   gpnorm - predicted function norm at the new point, assuming local
4285:     linearization.  The value is zero if the step lies within the trust
4286:     region, and exceeds zero otherwise.
4287: -   ynorm - 2-norm of the step

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

4293: .keywords: SNES, nonlinear, scale, step
4294: */
4295: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
4296: {
4297:   PetscReal      nrm;
4298:   PetscScalar    cnorm;


4306:   VecNorm(y,NORM_2,&nrm);
4307:   if (nrm > *delta) {
4308:     nrm     = *delta/nrm;
4309:     *gpnorm = (1.0 - nrm)*(*fnorm);
4310:     cnorm   = nrm;
4311:     VecScale(y,cnorm);
4312:     *ynorm  = *delta;
4313:   } else {
4314:     *gpnorm = 0.0;
4315:     *ynorm  = nrm;
4316:   }
4317:   return(0);
4318: }

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

4323:    Collective on SNES

4325:    Parameter:
4326: +  snes - iterative context obtained from SNESCreate()
4327: -  viewer - the viewer to display the reason


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

4333:    Level: beginner

4335: .keywords: SNES, solve, linear system

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

4339: @*/
4340: PetscErrorCode  SNESReasonView(SNES snes,PetscViewer viewer)
4341: {
4342:   PetscViewerFormat format;
4343:   PetscBool         isAscii;
4344:   PetscErrorCode    ierr;

4347:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
4348:   if (isAscii) {
4349:     PetscViewerGetFormat(viewer, &format);
4350:     PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
4351:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
4352:       DM                dm;
4353:       Vec               u;
4354:       PetscDS           prob;
4355:       PetscInt          Nf, f;
4356:       PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
4357:       void            **exactCtx;
4358:       PetscReal         error;

4360:       SNESGetDM(snes, &dm);
4361:       SNESGetSolution(snes, &u);
4362:       DMGetDS(dm, &prob);
4363:       PetscDSGetNumFields(prob, &Nf);
4364:       PetscMalloc2(Nf, &exactSol, Nf, &exactCtx);
4365:       for (f = 0; f < Nf; ++f) {PetscDSGetExactSolution(prob, f, &exactSol[f], &exactCtx[f]);}
4366:       DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error);
4367:       PetscFree2(exactSol, exactCtx);
4368:       if (error < 1.0e-11) {PetscViewerASCIIPrintf(viewer, "L_2 Error: < 1.0e-11\n");}
4369:       else                 {PetscViewerASCIIPrintf(viewer, "L_2 Error: %g\n", error);}
4370:     }
4371:     if (snes->reason > 0) {
4372:       if (((PetscObject) snes)->prefix) {
4373:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4374:       } else {
4375:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4376:       }
4377:     } else {
4378:       if (((PetscObject) snes)->prefix) {
4379:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4380:       } else {
4381:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4382:       }
4383:     }
4384:     PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
4385:   }
4386:   return(0);
4387: }

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

4392:   Collective on SNES

4394:   Input Parameters:
4395: . snes   - the SNES object

4397:   Level: intermediate

4399: @*/
4400: PetscErrorCode SNESReasonViewFromOptions(SNES snes)
4401: {
4402:   PetscErrorCode    ierr;
4403:   PetscViewer       viewer;
4404:   PetscBool         flg;
4405:   static PetscBool  incall = PETSC_FALSE;
4406:   PetscViewerFormat format;

4409:   if (incall) return(0);
4410:   incall = PETSC_TRUE;
4411:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
4412:   if (flg) {
4413:     PetscViewerPushFormat(viewer,format);
4414:     SNESReasonView(snes,viewer);
4415:     PetscViewerPopFormat(viewer);
4416:     PetscViewerDestroy(&viewer);
4417:   }
4418:   incall = PETSC_FALSE;
4419:   return(0);
4420: }

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

4426:    Collective on SNES

4428:    Input Parameters:
4429: +  snes - the SNES context
4430: .  b - the constant part of the equation F(x) = b, or NULL to use zero.
4431: -  x - the solution vector.

4433:    Notes:
4434:    The user should initialize the vector,x, with the initial guess
4435:    for the nonlinear solve prior to calling SNESSolve.  In particular,
4436:    to employ an initial guess of zero, the user should explicitly set
4437:    this vector to zero by calling VecSet().

4439:    Level: beginner

4441: .keywords: SNES, nonlinear, solve

4443: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
4444: @*/
4445: PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
4446: {
4447:   PetscErrorCode    ierr;
4448:   PetscBool         flg;
4449:   PetscInt          grid;
4450:   Vec               xcreated = NULL;
4451:   DM                dm;


4460:   /* High level operations using the nonlinear solver */
4461:   {
4462:     PetscViewer       viewer;
4463:     PetscViewerFormat format;
4464:     PetscInt          num;
4465:     PetscBool         flg;
4466:     static PetscBool  incall = PETSC_FALSE;

4468:     if (!incall) {
4469:       /* Estimate the convergence rate of the discretization */
4470:       PetscOptionsGetViewer(PetscObjectComm((PetscObject) snes),((PetscObject)snes)->options, ((PetscObject) snes)->prefix, "-snes_convergence_estimate", &viewer, &format, &flg);
4471:       if (flg) {
4472:         PetscConvEst conv;
4473:         DM           dm;
4474:         PetscReal   *alpha; /* Convergence rate of the solution error for each field in the L_2 norm */
4475:         PetscInt     Nf;

4477:         incall = PETSC_TRUE;
4478:         SNESGetDM(snes, &dm);
4479:         DMGetNumFields(dm, &Nf);
4480:         PetscMalloc1(Nf, &alpha);
4481:         PetscConvEstCreate(PetscObjectComm((PetscObject) snes), &conv);
4482:         PetscConvEstSetSolver(conv, snes);
4483:         PetscConvEstSetFromOptions(conv);
4484:         PetscConvEstSetUp(conv);
4485:         PetscConvEstGetConvRate(conv, alpha);
4486:         PetscViewerPushFormat(viewer, format);
4487:         PetscConvEstRateView(conv, alpha, viewer);
4488:         PetscViewerPopFormat(viewer);
4489:         PetscViewerDestroy(&viewer);
4490:         PetscConvEstDestroy(&conv);
4491:         PetscFree(alpha);
4492:         incall = PETSC_FALSE;
4493:       }
4494:       /* Adaptively refine the initial grid */
4495:       num  = 1;
4496:       PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_initial", &num, &flg);
4497:       if (flg) {
4498:         DMAdaptor adaptor;

4500:         incall = PETSC_TRUE;
4501:         DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);
4502:         DMAdaptorSetSolver(adaptor, snes);
4503:         DMAdaptorSetSequenceLength(adaptor, num);
4504:         DMAdaptorSetFromOptions(adaptor);
4505:         DMAdaptorSetUp(adaptor);
4506:         DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_INITIAL, &dm, &x);
4507:         DMAdaptorDestroy(&adaptor);
4508:         incall = PETSC_FALSE;
4509:       }
4510:       /* Use grid sequencing to adapt */
4511:       num  = 0;
4512:       PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_sequence", &num, NULL);
4513:       if (num) {
4514:         DMAdaptor adaptor;

4516:         incall = PETSC_TRUE;
4517:         DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);
4518:         DMAdaptorSetSolver(adaptor, snes);
4519:         DMAdaptorSetSequenceLength(adaptor, num);
4520:         DMAdaptorSetFromOptions(adaptor);
4521:         DMAdaptorSetUp(adaptor);
4522:         DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_SEQUENTIAL, &dm, &x);
4523:         DMAdaptorDestroy(&adaptor);
4524:         incall = PETSC_FALSE;
4525:       }
4526:     }
4527:   }
4528:   if (!x) {
4529:     SNESGetDM(snes,&dm);
4530:     DMCreateGlobalVector(dm,&xcreated);
4531:     x    = xcreated;
4532:   }
4533:   SNESViewFromOptions(snes,NULL,"-snes_view_pre");

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

4538:     /* set solution vector */
4539:     if (!grid) {PetscObjectReference((PetscObject)x);}
4540:     VecDestroy(&snes->vec_sol);
4541:     snes->vec_sol = x;
4542:     SNESGetDM(snes,&dm);

4544:     /* set affine vector if provided */
4545:     if (b) { PetscObjectReference((PetscObject)b); }
4546:     VecDestroy(&snes->vec_rhs);
4547:     snes->vec_rhs = b;

4549:     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");
4550:     if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
4551:     if (snes->vec_rhs  == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
4552:     if (!snes->vec_sol_update /* && snes->vec_sol */) {
4553:       VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
4554:       PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
4555:     }
4556:     DMShellSetGlobalVector(dm,snes->vec_sol);
4557:     SNESSetUp(snes);

4559:     if (!grid) {
4560:       if (snes->ops->computeinitialguess) {
4561:         (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
4562:       }
4563:     }

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

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

4574:     if (snes->lagjac_persist) snes->jac_iter += snes->iter;
4575:     if (snes->lagpre_persist) snes->pre_iter += snes->iter;

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

4581:     if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
4582:     if (snes->reason < 0) break;
4583:     if (grid <  snes->gridsequence) {
4584:       DM  fine;
4585:       Vec xnew;
4586:       Mat interp;

4588:       DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
4589:       if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
4590:       DMCreateInterpolation(snes->dm,fine,&interp,NULL);
4591:       DMCreateGlobalVector(fine,&xnew);
4592:       MatInterpolate(interp,x,xnew);
4593:       DMInterpolate(snes->dm,interp,fine);
4594:       MatDestroy(&interp);
4595:       x    = xnew;

4597:       SNESReset(snes);
4598:       SNESSetDM(snes,fine);
4599:       SNESResetFromOptions(snes);
4600:       DMDestroy(&fine);
4601:       PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
4602:     }
4603:   }
4604:   SNESViewFromOptions(snes,NULL,"-snes_view");
4605:   VecViewFromOptions(snes->vec_sol,(PetscObject)snes,"-snes_view_solution");

4607:   VecDestroy(&xcreated);
4608:   PetscObjectSAWsBlock((PetscObject)snes);
4609:   return(0);
4610: }

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

4614: /*@C
4615:    SNESSetType - Sets the method for the nonlinear solver.

4617:    Collective on SNES

4619:    Input Parameters:
4620: +  snes - the SNES context
4621: -  type - a known method

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

4627:    Notes:
4628:    See "petsc/include/petscsnes.h" for available methods (for instance)
4629: +    SNESNEWTONLS - Newton's method with line search
4630:      (systems of nonlinear equations)
4631: .    SNESNEWTONTR - Newton's method with trust region
4632:      (systems of nonlinear equations)

4634:   Normally, it is best to use the SNESSetFromOptions() command and then
4635:   set the SNES solver type from the options database rather than by using
4636:   this routine.  Using the options database provides the user with
4637:   maximum flexibility in evaluating the many nonlinear solvers.
4638:   The SNESSetType() routine is provided for those situations where it
4639:   is necessary to set the nonlinear solver independently of the command
4640:   line or options database.  This might be the case, for example, when
4641:   the choice of solver changes during the execution of the program,
4642:   and the user's application is taking responsibility for choosing the
4643:   appropriate method.

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

4649:   Level: intermediate

4651: .keywords: SNES, set, type

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

4655: @*/
4656: PetscErrorCode  SNESSetType(SNES snes,SNESType type)
4657: {
4658:   PetscErrorCode ierr,(*r)(SNES);
4659:   PetscBool      match;


4665:   PetscObjectTypeCompare((PetscObject)snes,type,&match);
4666:   if (match) return(0);

4668:    PetscFunctionListFind(SNESList,type,&r);
4669:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
4670:   /* Destroy the previous private SNES context */
4671:   if (snes->ops->destroy) {
4672:     (*(snes)->ops->destroy)(snes);
4673:     snes->ops->destroy = NULL;
4674:   }
4675:   /* Reinitialize function pointers in SNESOps structure */
4676:   snes->ops->setup          = 0;
4677:   snes->ops->solve          = 0;
4678:   snes->ops->view           = 0;
4679:   snes->ops->setfromoptions = 0;
4680:   snes->ops->destroy        = 0;
4681:   SNESLineSearchDestroy(&snes->linesearch);
4682:   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
4683:   snes->setupcalled = PETSC_FALSE;

4685:   PetscObjectChangeTypeName((PetscObject)snes,type);
4686:   (*r)(snes);
4687:   return(0);
4688: }

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

4693:    Not Collective

4695:    Input Parameter:
4696: .  snes - nonlinear solver context

4698:    Output Parameter:
4699: .  type - SNES method (a character string)

4701:    Level: intermediate

4703: .keywords: SNES, nonlinear, get, type, name
4704: @*/
4705: PetscErrorCode  SNESGetType(SNES snes,SNESType *type)
4706: {
4710:   *type = ((PetscObject)snes)->type_name;
4711:   return(0);
4712: }

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

4717:   Logically Collective on SNES and Vec

4719:   Input Parameters:
4720: + snes - the SNES context obtained from SNESCreate()
4721: - u    - the solution vector

4723:   Level: beginner

4725: .keywords: SNES, set, solution
4726: @*/
4727: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4728: {
4729:   DM             dm;

4735:   PetscObjectReference((PetscObject) u);
4736:   VecDestroy(&snes->vec_sol);

4738:   snes->vec_sol = u;

4740:   SNESGetDM(snes, &dm);
4741:   DMShellSetGlobalVector(dm, u);
4742:   return(0);
4743: }

4745: /*@
4746:    SNESGetSolution - Returns the vector where the approximate solution is
4747:    stored. This is the fine grid solution when using SNESSetGridSequence().

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

4751:    Input Parameter:
4752: .  snes - the SNES context

4754:    Output Parameter:
4755: .  x - the solution

4757:    Level: intermediate

4759: .keywords: SNES, nonlinear, get, solution

4761: .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
4762: @*/
4763: PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
4764: {
4768:   *x = snes->vec_sol;
4769:   return(0);
4770: }

4772: /*@
4773:    SNESGetSolutionUpdate - Returns the vector where the solution update is
4774:    stored.

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

4778:    Input Parameter:
4779: .  snes - the SNES context

4781:    Output Parameter:
4782: .  x - the solution update

4784:    Level: advanced

4786: .keywords: SNES, nonlinear, get, solution, update

4788: .seealso: SNESGetSolution(), SNESGetFunction()
4789: @*/
4790: PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
4791: {
4795:   *x = snes->vec_sol_update;
4796:   return(0);
4797: }

4799: /*@C
4800:    SNESGetFunction - Returns the vector where the function is stored.

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

4804:    Input Parameter:
4805: .  snes - the SNES context

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

4812:    Level: advanced

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

4816: .keywords: SNES, nonlinear, get, function

4818: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4819: @*/
4820: PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4821: {
4823:   DM             dm;

4827:   if (r) {
4828:     if (!snes->vec_func) {
4829:       if (snes->vec_rhs) {
4830:         VecDuplicate(snes->vec_rhs,&snes->vec_func);
4831:       } else if (snes->vec_sol) {
4832:         VecDuplicate(snes->vec_sol,&snes->vec_func);
4833:       } else if (snes->dm) {
4834:         DMCreateGlobalVector(snes->dm,&snes->vec_func);
4835:       }
4836:     }
4837:     *r = snes->vec_func;
4838:   }
4839:   SNESGetDM(snes,&dm);
4840:   DMSNESGetFunction(dm,f,ctx);
4841:   return(0);
4842: }

4844: /*@C
4845:    SNESGetNGS - Returns the NGS function and context.

4847:    Input Parameter:
4848: .  snes - the SNES context

4850:    Output Parameter:
4851: +  f - the function (or NULL) see SNESNGSFunction for details
4852: -  ctx    - the function context (or NULL)

4854:    Level: advanced

4856: .keywords: SNES, nonlinear, get, function

4858: .seealso: SNESSetNGS(), SNESGetFunction()
4859: @*/

4861: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4862: {
4864:   DM             dm;

4868:   SNESGetDM(snes,&dm);
4869:   DMSNESGetNGS(dm,f,ctx);
4870:   return(0);
4871: }

4873: /*@C
4874:    SNESSetOptionsPrefix - Sets the prefix used for searching for all
4875:    SNES options in the database.

4877:    Logically Collective on SNES

4879:    Input Parameter:
4880: +  snes - the SNES context
4881: -  prefix - the prefix to prepend to all option names

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

4887:    Level: advanced

4889: .keywords: SNES, set, options, prefix, database

4891: .seealso: SNESSetFromOptions()
4892: @*/
4893: PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
4894: {

4899:   PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4900:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4901:   if (snes->linesearch) {
4902:     SNESGetLineSearch(snes,&snes->linesearch);
4903:     PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4904:   }
4905:   KSPSetOptionsPrefix(snes->ksp,prefix);
4906:   return(0);
4907: }

4909: /*@C
4910:    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4911:    SNES options in the database.

4913:    Logically Collective on SNES

4915:    Input Parameters:
4916: +  snes - the SNES context
4917: -  prefix - the prefix to prepend to all option names

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

4923:    Level: advanced

4925: .keywords: SNES, append, options, prefix, database

4927: .seealso: SNESGetOptionsPrefix()
4928: @*/
4929: PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4930: {

4935:   PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4936:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4937:   if (snes->linesearch) {
4938:     SNESGetLineSearch(snes,&snes->linesearch);
4939:     PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4940:   }
4941:   KSPAppendOptionsPrefix(snes->ksp,prefix);
4942:   return(0);
4943: }

4945: /*@C
4946:    SNESGetOptionsPrefix - Sets the prefix used for searching for all
4947:    SNES options in the database.

4949:    Not Collective

4951:    Input Parameter:
4952: .  snes - the SNES context

4954:    Output Parameter:
4955: .  prefix - pointer to the prefix string used

4957:    Notes:
4958:     On the fortran side, the user should pass in a string 'prefix' of
4959:    sufficient length to hold the prefix.

4961:    Level: advanced

4963: .keywords: SNES, get, options, prefix, database

4965: .seealso: SNESAppendOptionsPrefix()
4966: @*/
4967: PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4968: {

4973:   PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4974:   return(0);
4975: }


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

4981:    Not collective

4983:    Input Parameters:
4984: +  name_solver - name of a new user-defined solver
4985: -  routine_create - routine to create method context

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

4990:    Sample usage:
4991: .vb
4992:    SNESRegister("my_solver",MySolverCreate);
4993: .ve

4995:    Then, your solver can be chosen with the procedural interface via
4996: $     SNESSetType(snes,"my_solver")
4997:    or at runtime via the option
4998: $     -snes_type my_solver

5000:    Level: advanced

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

5004: .keywords: SNES, nonlinear, register

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

5008:   Level: advanced
5009: @*/
5010: PetscErrorCode  SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
5011: {

5015:   SNESInitializePackage();
5016:   PetscFunctionListAdd(&SNESList,sname,function);
5017:   return(0);
5018: }

5020: PetscErrorCode  SNESTestLocalMin(SNES snes)
5021: {
5023:   PetscInt       N,i,j;
5024:   Vec            u,uh,fh;
5025:   PetscScalar    value;
5026:   PetscReal      norm;

5029:   SNESGetSolution(snes,&u);
5030:   VecDuplicate(u,&uh);
5031:   VecDuplicate(u,&fh);

5033:   /* currently only works for sequential */
5034:   PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
5035:   VecGetSize(u,&N);
5036:   for (i=0; i<N; i++) {
5037:     VecCopy(u,uh);
5038:     PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
5039:     for (j=-10; j<11; j++) {
5040:       value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
5041:       VecSetValue(uh,i,value,ADD_VALUES);
5042:       SNESComputeFunction(snes,uh,fh);
5043:       VecNorm(fh,NORM_2,&norm);
5044:       PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);
5045:       value = -value;
5046:       VecSetValue(uh,i,value,ADD_VALUES);
5047:     }
5048:   }
5049:   VecDestroy(&uh);
5050:   VecDestroy(&fh);
5051:   return(0);
5052: }

5054: /*@
5055:    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
5056:    computing relative tolerance for linear solvers within an inexact
5057:    Newton method.

5059:    Logically Collective on SNES

5061:    Input Parameters:
5062: +  snes - SNES context
5063: -  flag - PETSC_TRUE or PETSC_FALSE

5065:     Options Database:
5066: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
5067: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
5068: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
5069: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
5070: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
5071: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
5072: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
5073: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

5075:    Notes:
5076:    Currently, the default is to use a constant relative tolerance for
5077:    the inner linear solvers.  Alternatively, one can use the
5078:    Eisenstat-Walker method, where the relative convergence tolerance
5079:    is reset at each Newton iteration according progress of the nonlinear
5080:    solver.

5082:    Level: advanced

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

5088: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton

5090: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
5091: @*/
5092: PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool flag)
5093: {
5097:   snes->ksp_ewconv = flag;
5098:   return(0);
5099: }

5101: /*@
5102:    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
5103:    for computing relative tolerance for linear solvers within an
5104:    inexact Newton method.

5106:    Not Collective

5108:    Input Parameter:
5109: .  snes - SNES context

5111:    Output Parameter:
5112: .  flag - PETSC_TRUE or PETSC_FALSE

5114:    Notes:
5115:    Currently, the default is to use a constant relative tolerance for
5116:    the inner linear solvers.  Alternatively, one can use the
5117:    Eisenstat-Walker method, where the relative convergence tolerance
5118:    is reset at each Newton iteration according progress of the nonlinear
5119:    solver.

5121:    Level: advanced

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

5127: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton

5129: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
5130: @*/
5131: PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
5132: {
5136:   *flag = snes->ksp_ewconv;
5137:   return(0);
5138: }

5140: /*@
5141:    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
5142:    convergence criteria for the linear solvers within an inexact
5143:    Newton method.

5145:    Logically Collective on SNES

5147:    Input Parameters:
5148: +    snes - SNES context
5149: .    version - version 1, 2 (default is 2) or 3
5150: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5151: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5152: .    gamma - multiplicative factor for version 2 rtol computation
5153:              (0 <= gamma2 <= 1)
5154: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
5155: .    alpha2 - power for safeguard
5156: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

5158:    Note:
5159:    Version 3 was contributed by Luis Chacon, June 2006.

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

5163:    Level: advanced

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

5170: .keywords: SNES, KSP, Eisenstat, Walker, set, parameters

5172: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
5173: @*/
5174: PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
5175: {
5176:   SNESKSPEW *kctx;

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

5190:   if (version != PETSC_DEFAULT)   kctx->version   = version;
5191:   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
5192:   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
5193:   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
5194:   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
5195:   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
5196:   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;

5198:   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);
5199:   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);
5200:   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);
5201:   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);
5202:   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);
5203:   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);
5204:   return(0);
5205: }

5207: /*@
5208:    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
5209:    convergence criteria for the linear solvers within an inexact
5210:    Newton method.

5212:    Not Collective

5214:    Input Parameters:
5215:      snes - SNES context

5217:    Output Parameters:
5218: +    version - version 1, 2 (default is 2) or 3
5219: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5220: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5221: .    gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
5222: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
5223: .    alpha2 - power for safeguard
5224: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

5226:    Level: advanced

5228: .keywords: SNES, KSP, Eisenstat, Walker, get, parameters

5230: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
5231: @*/
5232: PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
5233: {
5234:   SNESKSPEW *kctx;

5238:   kctx = (SNESKSPEW*)snes->kspconvctx;
5239:   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
5240:   if (version)   *version   = kctx->version;
5241:   if (rtol_0)    *rtol_0    = kctx->rtol_0;
5242:   if (rtol_max)  *rtol_max  = kctx->rtol_max;
5243:   if (gamma)     *gamma     = kctx->gamma;
5244:   if (alpha)     *alpha     = kctx->alpha;
5245:   if (alpha2)    *alpha2    = kctx->alpha2;
5246:   if (threshold) *threshold = kctx->threshold;
5247:   return(0);
5248: }

5250:  PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5251: {
5253:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
5254:   PetscReal      rtol  = PETSC_DEFAULT,stol;

5257:   if (!snes->ksp_ewconv) return(0);
5258:   if (!snes->iter) {
5259:     rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
5260:     VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
5261:   }
5262:   else {
5263:     if (kctx->version == 1) {
5264:       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
5265:       if (rtol < 0.0) rtol = -rtol;
5266:       stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
5267:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5268:     } else if (kctx->version == 2) {
5269:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5270:       stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
5271:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5272:     } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
5273:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5274:       /* safeguard: avoid sharp decrease of rtol */
5275:       stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
5276:       stol = PetscMax(rtol,stol);
5277:       rtol = PetscMin(kctx->rtol_0,stol);
5278:       /* safeguard: avoid oversolving */
5279:       stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
5280:       stol = PetscMax(rtol,stol);
5281:       rtol = PetscMin(kctx->rtol_0,stol);
5282:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
5283:   }
5284:   /* safeguard: avoid rtol greater than one */
5285:   rtol = PetscMin(rtol,kctx->rtol_max);
5286:   KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
5287:   PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
5288:   return(0);
5289: }

5291: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5292: {
5294:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
5295:   PCSide         pcside;
5296:   Vec            lres;

5299:   if (!snes->ksp_ewconv) return(0);
5300:   KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
5301:   kctx->norm_last = snes->norm;
5302:   if (kctx->version == 1) {
5303:     PC        pc;
5304:     PetscBool isNone;

5306:     KSPGetPC(ksp, &pc);
5307:     PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
5308:     KSPGetPCSide(ksp,&pcside);
5309:      if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
5310:       /* KSP residual is true linear residual */
5311:       KSPGetResidualNorm(ksp,&kctx->lresid_last);
5312:     } else {
5313:       /* KSP residual is preconditioned residual */
5314:       /* compute true linear residual norm */
5315:       VecDuplicate(b,&lres);
5316:       MatMult(snes->jacobian,x,lres);
5317:       VecAYPX(lres,-1.0,b);
5318:       VecNorm(lres,NORM_2,&kctx->lresid_last);
5319:       VecDestroy(&lres);
5320:     }
5321:   }
5322:   return(0);
5323: }

5325: /*@
5326:    SNESGetKSP - Returns the KSP context for a SNES solver.

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

5330:    Input Parameter:
5331: .  snes - the SNES context

5333:    Output Parameter:
5334: .  ksp - the KSP context

5336:    Notes:
5337:    The user can then directly manipulate the KSP context to set various
5338:    options, etc.  Likewise, the user can then extract and manipulate the
5339:    PC contexts as well.

5341:    Level: beginner

5343: .keywords: SNES, nonlinear, get, KSP, context

5345: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
5346: @*/
5347: PetscErrorCode  SNESGetKSP(SNES snes,KSP *ksp)
5348: {


5355:   if (!snes->ksp) {
5356:     PetscBool monitor = PETSC_FALSE;

5358:     KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
5359:     PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
5360:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);

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

5365:     PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
5366:     if (monitor) {
5367:       KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
5368:     }
5369:     monitor = PETSC_FALSE;
5370:     PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
5371:     if (monitor) {
5372:       PetscObject *objs;
5373:       KSPMonitorSNESLGResidualNormCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
5374:       objs[0] = (PetscObject) snes;
5375:       KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
5376:     }
5377:     PetscObjectSetOptions((PetscObject)snes->ksp,((PetscObject)snes)->options);
5378:   }
5379:   *ksp = snes->ksp;
5380:   return(0);
5381: }


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

5388:    Logically Collective on SNES

5390:    Input Parameters:
5391: +  snes - the nonlinear solver context
5392: -  dm - the dm, cannot be NULL

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

5399:    Level: intermediate

5401: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
5402: @*/
5403: PetscErrorCode  SNESSetDM(SNES snes,DM dm)
5404: {
5406:   KSP            ksp;
5407:   DMSNES         sdm;

5412:   PetscObjectReference((PetscObject)dm);
5413:   if (snes->dm) {               /* Move the DMSNES context over to the new DM unless the new DM already has one */
5414:     if (snes->dm->dmsnes && !dm->dmsnes) {
5415:       DMCopyDMSNES(snes->dm,dm);
5416:       DMGetDMSNES(snes->dm,&sdm);
5417:       if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
5418:     }
5419:     DMCoarsenHookRemove(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
5420:     DMDestroy(&snes->dm);
5421:   }
5422:   snes->dm     = dm;
5423:   snes->dmAuto = PETSC_FALSE;

5425:   SNESGetKSP(snes,&ksp);
5426:   KSPSetDM(ksp,dm);
5427:   KSPSetDMActive(ksp,PETSC_FALSE);
5428:   if (snes->npc) {
5429:     SNESSetDM(snes->npc, snes->dm);
5430:     SNESSetNPCSide(snes,snes->npcside);
5431:   }
5432:   return(0);
5433: }

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

5438:    Not Collective but DM obtained is parallel on SNES

5440:    Input Parameter:
5441: . snes - the preconditioner context

5443:    Output Parameter:
5444: .  dm - the dm

5446:    Level: intermediate

5448: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
5449: @*/
5450: PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
5451: {

5456:   if (!snes->dm) {
5457:     DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
5458:     snes->dmAuto = PETSC_TRUE;
5459:   }
5460:   *dm = snes->dm;
5461:   return(0);
5462: }

5464: /*@
5465:   SNESSetNPC - Sets the nonlinear preconditioner to be used.

5467:   Collective on SNES

5469:   Input Parameters:
5470: + snes - iterative context obtained from SNESCreate()
5471: - pc   - the preconditioner object

5473:   Notes:
5474:   Use SNESGetNPC() to retrieve the preconditioner context (for example,
5475:   to configure it using the API).

5477:   Level: developer

5479: .keywords: SNES, set, precondition
5480: .seealso: SNESGetNPC(), SNESHasNPC()
5481: @*/
5482: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
5483: {

5490:   PetscObjectReference((PetscObject) pc);
5491:   SNESDestroy(&snes->npc);
5492:   snes->npc = pc;
5493:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->npc);
5494:   return(0);
5495: }

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

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

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

5505:   Output Parameter:
5506: . pc - preconditioner context

5508:   Notes:
5509:     If a SNES was previously set with SNESSetNPC() then that SNES is returned.

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

5514:   Level: developer

5516: .keywords: SNES, get, preconditioner
5517: .seealso: SNESSetNPC(), SNESHasNPC(), SNES, SNESCreate()
5518: @*/
5519: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
5520: {
5522:   const char     *optionsprefix;

5527:   if (!snes->npc) {
5528:     SNESCreate(PetscObjectComm((PetscObject)snes),&snes->npc);
5529:     PetscObjectIncrementTabLevel((PetscObject)snes->npc,(PetscObject)snes,1);
5530:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->npc);
5531:     SNESGetOptionsPrefix(snes,&optionsprefix);
5532:     SNESSetOptionsPrefix(snes->npc,optionsprefix);
5533:     SNESAppendOptionsPrefix(snes->npc,"npc_");
5534:     SNESSetCountersReset(snes->npc,PETSC_FALSE);
5535:   }
5536:   *pc = snes->npc;
5537:   return(0);
5538: }

5540: /*@
5541:   SNESHasNPC - Returns whether a nonlinear preconditioner exists

5543:   Not Collective

5545:   Input Parameter:
5546: . snes - iterative context obtained from SNESCreate()

5548:   Output Parameter:
5549: . has_npc - whether the SNES has an NPC or not

5551:   Level: developer

5553: .keywords: SNES, has, preconditioner
5554: .seealso: SNESSetNPC(), SNESGetNPC()
5555: @*/
5556: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
5557: {
5560:   *has_npc = (PetscBool) (snes->npc ? PETSC_TRUE : PETSC_FALSE);
5561:   return(0);
5562: }

5564: /*@
5565:     SNESSetNPCSide - Sets the preconditioning side.

5567:     Logically Collective on SNES

5569:     Input Parameter:
5570: .   snes - iterative context obtained from SNESCreate()

5572:     Output Parameter:
5573: .   side - the preconditioning side, where side is one of
5574: .vb
5575:       PC_LEFT - left preconditioning
5576:       PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5577: .ve

5579:     Options Database Keys:
5580: .   -snes_pc_side <right,left>

5582:     Notes:
5583:     SNESNRICHARDSON and SNESNCG only support left preconditioning.

5585:     Level: intermediate

5587: .keywords: SNES, set, right, left, side, preconditioner, flag

5589: .seealso: SNESGetNPCSide(), KSPSetPCSide()
5590: @*/
5591: PetscErrorCode  SNESSetNPCSide(SNES snes,PCSide side)
5592: {
5596:   snes->npcside= side;
5597:   return(0);
5598: }

5600: /*@
5601:     SNESGetNPCSide - Gets the preconditioning side.

5603:     Not Collective

5605:     Input Parameter:
5606: .   snes - iterative context obtained from SNESCreate()

5608:     Output Parameter:
5609: .   side - the preconditioning side, where side is one of
5610: .vb
5611:       PC_LEFT - left preconditioning
5612:       PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5613: .ve

5615:     Level: intermediate

5617: .keywords: SNES, get, right, left, side, preconditioner, flag

5619: .seealso: SNESSetNPCSide(), KSPGetPCSide()
5620: @*/
5621: PetscErrorCode  SNESGetNPCSide(SNES snes,PCSide *side)
5622: {
5626:   *side = snes->npcside;
5627:   return(0);
5628: }

5630: /*@
5631:   SNESSetLineSearch - Sets the linesearch on the SNES instance.

5633:   Collective on SNES

5635:   Input Parameters:
5636: + snes - iterative context obtained from SNESCreate()
5637: - linesearch   - the linesearch object

5639:   Notes:
5640:   Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
5641:   to configure it using the API).

5643:   Level: developer

5645: .keywords: SNES, set, linesearch
5646: .seealso: SNESGetLineSearch()
5647: @*/
5648: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
5649: {

5656:   PetscObjectReference((PetscObject) linesearch);
5657:   SNESLineSearchDestroy(&snes->linesearch);

5659:   snes->linesearch = linesearch;

5661:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5662:   return(0);
5663: }

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

5669:   Not Collective

5671:   Input Parameter:
5672: . snes - iterative context obtained from SNESCreate()

5674:   Output Parameter:
5675: . linesearch - linesearch context

5677:   Level: beginner

5679: .keywords: SNES, get, linesearch
5680: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5681: @*/
5682: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5683: {
5685:   const char     *optionsprefix;

5690:   if (!snes->linesearch) {
5691:     SNESGetOptionsPrefix(snes, &optionsprefix);
5692:     SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5693:     SNESLineSearchSetSNES(snes->linesearch, snes);
5694:     SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5695:     PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5696:     PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5697:   }
5698:   *linesearch = snes->linesearch;
5699:   return(0);
5700: }

5702: #if defined(PETSC_HAVE_MATLAB_ENGINE)
5703: #include <mex.h>

5705: typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;

5707: /*
5708:    SNESComputeFunction_Matlab - Calls the function that has been set with SNESSetFunctionMatlab().

5710:    Collective on SNES

5712:    Input Parameters:
5713: +  snes - the SNES context
5714: -  x - input vector

5716:    Output Parameter:
5717: .  y - function vector, as set by SNESSetFunction()

5719:    Notes:
5720:    SNESComputeFunction() is typically used within nonlinear solvers
5721:    implementations, so most users would not generally call this routine
5722:    themselves.

5724:    Level: developer

5726: .keywords: SNES, nonlinear, compute, function

5728: .seealso: SNESSetFunction(), SNESGetFunction()
5729: */
5730: PetscErrorCode  SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
5731: {
5732:   PetscErrorCode    ierr;
5733:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5734:   int               nlhs  = 1,nrhs = 5;
5735:   mxArray           *plhs[1],*prhs[5];
5736:   long long int     lx = 0,ly = 0,ls = 0;


5745:   /* call Matlab function in ctx with arguments x and y */

5747:   PetscMemcpy(&ls,&snes,sizeof(snes));
5748:   PetscMemcpy(&lx,&x,sizeof(x));
5749:   PetscMemcpy(&ly,&y,sizeof(x));
5750:   prhs[0] = mxCreateDoubleScalar((double)ls);
5751:   prhs[1] = mxCreateDoubleScalar((double)lx);
5752:   prhs[2] = mxCreateDoubleScalar((double)ly);
5753:   prhs[3] = mxCreateString(sctx->funcname);
5754:   prhs[4] = sctx->ctx;
5755:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");
5756:   mxGetScalar(plhs[0]);
5757:   mxDestroyArray(prhs[0]);
5758:   mxDestroyArray(prhs[1]);
5759:   mxDestroyArray(prhs[2]);
5760:   mxDestroyArray(prhs[3]);
5761:   mxDestroyArray(plhs[0]);
5762:   return(0);
5763: }

5765: /*
5766:    SNESSetFunctionMatlab - Sets the function evaluation routine and function
5767:    vector for use by the SNES routines in solving systems of nonlinear
5768:    equations from MATLAB. Here the function is a string containing the name of a MATLAB function

5770:    Logically Collective on SNES

5772:    Input Parameters:
5773: +  snes - the SNES context
5774: .  r - vector to store function value
5775: -  f - function evaluation routine

5777:    Notes:
5778:    The Newton-like methods typically solve linear systems of the form
5779: $      f'(x) x = -f(x),
5780:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

5782:    Level: beginner

5784:    Developer Note:  This bleeds the allocated memory SNESMatlabContext *sctx;

5786: .keywords: SNES, nonlinear, set, function

5788: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5789: */
5790: PetscErrorCode  SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx)
5791: {
5792:   PetscErrorCode    ierr;
5793:   SNESMatlabContext *sctx;

5796:   /* currently sctx is memory bleed */
5797:   PetscNew(&sctx);
5798:   PetscStrallocpy(f,&sctx->funcname);
5799:   /*
5800:      This should work, but it doesn't
5801:   sctx->ctx = ctx;
5802:   mexMakeArrayPersistent(sctx->ctx);
5803:   */
5804:   sctx->ctx = mxDuplicateArray(ctx);
5805:   SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);
5806:   return(0);
5807: }

5809: /*
5810:    SNESComputeJacobian_Matlab - Calls the function that has been set with SNESSetJacobianMatlab().

5812:    Collective on SNES

5814:    Input Parameters:
5815: +  snes - the SNES context
5816: .  x - input vector
5817: .  A, B - the matrices
5818: -  ctx - user context

5820:    Level: developer

5822: .keywords: SNES, nonlinear, compute, function

5824: .seealso: SNESSetFunction(), SNESGetFunction()
5825: @*/
5826: PetscErrorCode  SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx)
5827: {
5828:   PetscErrorCode    ierr;
5829:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5830:   int               nlhs  = 2,nrhs = 6;
5831:   mxArray           *plhs[2],*prhs[6];
5832:   long long int     lx = 0,lA = 0,ls = 0, lB = 0;


5838:   /* call Matlab function in ctx with arguments x and y */

5840:   PetscMemcpy(&ls,&snes,sizeof(snes));
5841:   PetscMemcpy(&lx,&x,sizeof(x));
5842:   PetscMemcpy(&lA,A,sizeof(x));
5843:   PetscMemcpy(&lB,B,sizeof(x));
5844:   prhs[0] = mxCreateDoubleScalar((double)ls);
5845:   prhs[1] = mxCreateDoubleScalar((double)lx);
5846:   prhs[2] = mxCreateDoubleScalar((double)lA);
5847:   prhs[3] = mxCreateDoubleScalar((double)lB);
5848:   prhs[4] = mxCreateString(sctx->funcname);
5849:   prhs[5] = sctx->ctx;
5850:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");
5851:   mxGetScalar(plhs[0]);
5852:   mxDestroyArray(prhs[0]);
5853:   mxDestroyArray(prhs[1]);
5854:   mxDestroyArray(prhs[2]);
5855:   mxDestroyArray(prhs[3]);
5856:   mxDestroyArray(prhs[4]);
5857:   mxDestroyArray(plhs[0]);
5858:   mxDestroyArray(plhs[1]);
5859:   return(0);
5860: }

5862: /*
5863:    SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
5864:    vector for use by the SNES routines in solving systems of nonlinear
5865:    equations from MATLAB. Here the function is a string containing the name of a MATLAB function

5867:    Logically Collective on SNES

5869:    Input Parameters:
5870: +  snes - the SNES context
5871: .  A,B - Jacobian matrices
5872: .  J - function evaluation routine
5873: -  ctx - user context

5875:    Level: developer

5877:    Developer Note:  This bleeds the allocated memory SNESMatlabContext *sctx;

5879: .keywords: SNES, nonlinear, set, function

5881: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J
5882: */
5883: PetscErrorCode  SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx)
5884: {
5885:   PetscErrorCode    ierr;
5886:   SNESMatlabContext *sctx;

5889:   /* currently sctx is memory bleed */
5890:   PetscNew(&sctx);
5891:   PetscStrallocpy(J,&sctx->funcname);
5892:   /*
5893:      This should work, but it doesn't
5894:   sctx->ctx = ctx;
5895:   mexMakeArrayPersistent(sctx->ctx);
5896:   */
5897:   sctx->ctx = mxDuplicateArray(ctx);
5898:   SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);
5899:   return(0);
5900: }

5902: /*
5903:    SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().

5905:    Collective on SNES

5907: .seealso: SNESSetFunction(), SNESGetFunction()
5908: @*/
5909: PetscErrorCode  SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
5910: {
5911:   PetscErrorCode    ierr;
5912:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5913:   int               nlhs  = 1,nrhs = 6;
5914:   mxArray           *plhs[1],*prhs[6];
5915:   long long int     lx = 0,ls = 0;
5916:   Vec               x  = snes->vec_sol;


5921:   PetscMemcpy(&ls,&snes,sizeof(snes));
5922:   PetscMemcpy(&lx,&x,sizeof(x));
5923:   prhs[0] = mxCreateDoubleScalar((double)ls);
5924:   prhs[1] = mxCreateDoubleScalar((double)it);
5925:   prhs[2] = mxCreateDoubleScalar((double)fnorm);
5926:   prhs[3] = mxCreateDoubleScalar((double)lx);
5927:   prhs[4] = mxCreateString(sctx->funcname);
5928:   prhs[5] = sctx->ctx;
5929:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");
5930:   mxGetScalar(plhs[0]);
5931:   mxDestroyArray(prhs[0]);
5932:   mxDestroyArray(prhs[1]);
5933:   mxDestroyArray(prhs[2]);
5934:   mxDestroyArray(prhs[3]);
5935:   mxDestroyArray(prhs[4]);
5936:   mxDestroyArray(plhs[0]);
5937:   return(0);
5938: }

5940: /*
5941:    SNESMonitorSetMatlab - Sets the monitor function from MATLAB

5943:    Level: developer

5945:    Developer Note:  This bleeds the allocated memory SNESMatlabContext *sctx;

5947: .keywords: SNES, nonlinear, set, function

5949: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5950: */
5951: PetscErrorCode  SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx)
5952: {
5953:   PetscErrorCode    ierr;
5954:   SNESMatlabContext *sctx;

5957:   /* currently sctx is memory bleed */
5958:   PetscNew(&sctx);
5959:   PetscStrallocpy(f,&sctx->funcname);
5960:   /*
5961:      This should work, but it doesn't
5962:   sctx->ctx = ctx;
5963:   mexMakeArrayPersistent(sctx->ctx);
5964:   */
5965:   sctx->ctx = mxDuplicateArray(ctx);
5966:   SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);
5967:   return(0);
5968: }

5970: #endif