Actual source code: snes.c

petsc-master 2014-11-21
Report Typos and Errors
  2: #include <petsc-private/snesimpl.h>      /*I "petscsnes.h"  I*/
  3: #include <petscdmshell.h>

  5: PetscBool         SNESRegisterAllCalled = PETSC_FALSE;
  6: PetscFunctionList SNESList              = NULL;

  8: /* Logging support */
  9: PetscClassId  SNES_CLASSID, DMSNES_CLASSID;
 10: PetscLogEvent SNES_Solve, SNES_FunctionEval, SNES_JacobianEval, SNES_NGSEval, SNES_NGSFuncEval, SNES_NPCSolve;

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

 17:    Logically Collective on SNES

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

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

 26:    Level: intermediate

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

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

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

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

 50:    Not Collective

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

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

 58:    Level: intermediate

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

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

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

 79:    Logically Collective on SNES

 81:    Input Parameters:
 82: .  snes - the SNES context

 84:    Level: advanced

 86: .keywords: SNES, view

 88: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction
 89: @*/
 90: PetscErrorCode  SNESSetFunctionDomainError(SNES snes)
 91: {
 94:   snes->domainerror = PETSC_TRUE;
 95:   return(0);
 96: }

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

103:    Logically Collective on SNES

105:    Input Parameters:
106: .  snes - the SNES context

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

111:    Level: advanced

113: .keywords: SNES, view

115: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction()
116: @*/
117: PetscErrorCode  SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror)
118: {
122:   *domainerror = snes->domainerror;
123:   return(0);
124: }

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

131:   Collective on PetscViewer

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

138:    Level: intermediate

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

143:   Notes for advanced users:
144:   Most users should not need to know the details of the binary storage
145:   format, since SNESLoad() and TSView() completely hide these details.
146:   But for anyone who's interested, the standard binary matrix storage
147:   format is
148: .vb
149:      has not yet been determined
150: .ve

152: .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad()
153: @*/
154: PetscErrorCode  SNESLoad(SNES snes, PetscViewer viewer)
155: {
157:   PetscBool      isbinary;
158:   PetscInt       classid;
159:   char           type[256];
160:   KSP            ksp;
161:   DM             dm;
162:   DMSNES         dmsnes;

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

170:   PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);
171:   if (classid != SNES_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONG,"Not SNES next in file");
172:   PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);
173:   SNESSetType(snes, type);
174:   if (snes->ops->load) {
175:     (*snes->ops->load)(snes,viewer);
176:   }
177:   SNESGetDM(snes,&dm);
178:   DMGetDMSNES(dm,&dmsnes);
179:   DMSNESLoad(dmsnes,viewer);
180:   SNESGetKSP(snes,&ksp);
181:   KSPLoad(ksp,viewer);
182:   return(0);
183: }

185: #include <petscdraw.h>
186: #if defined(PETSC_HAVE_SAWS)
187: #include <petscviewersaws.h>
188: #endif
191: /*@C
192:    SNESView - Prints the SNES data structure.

194:    Collective on SNES

196:    Input Parameters:
197: +  SNES - the SNES context
198: -  viewer - visualization context

200:    Options Database Key:
201: .  -snes_view - Calls SNESView() at end of SNESSolve()

203:    Notes:
204:    The available visualization contexts include
205: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
206: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
207:          output where only the first processor opens
208:          the file.  All other processors send their
209:          data to the first processor to print.

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

214:    Level: beginner

216: .keywords: SNES, view

218: .seealso: PetscViewerASCIIOpen()
219: @*/
220: PetscErrorCode  SNESView(SNES snes,PetscViewer viewer)
221: {
222:   SNESKSPEW      *kctx;
224:   KSP            ksp;
225:   SNESLineSearch linesearch;
226:   PetscBool      iascii,isstring,isbinary,isdraw;
227:   DMSNES         dmsnes;
228: #if defined(PETSC_HAVE_SAWS)
229:   PetscBool      isams;
230: #endif

234:   if (!viewer) {
235:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
236:   }

240:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
241:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
242:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
243:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
244: #if defined(PETSC_HAVE_SAWS)
245:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&isams);
246: #endif
247:   if (iascii) {
248:     PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);
249:     if (!snes->setupcalled) {
250:       PetscViewerASCIIPrintf(viewer,"  SNES has not been set up so information may be incomplete\n");
251:     }
252:     if (snes->ops->view) {
253:       PetscViewerASCIIPushTab(viewer);
254:       (*snes->ops->view)(snes,viewer);
255:       PetscViewerASCIIPopTab(viewer);
256:     }
257:     PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);
258:     PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);
259:     PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);
260:     PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);
261:     if (snes->gridsequence) {
262:       PetscViewerASCIIPrintf(viewer,"  total number of grid sequence refinements=%D\n",snes->gridsequence);
263:     }
264:     if (snes->ksp_ewconv) {
265:       kctx = (SNESKSPEW*)snes->kspconvctx;
266:       if (kctx) {
267:         PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);
268:         PetscViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);
269:         PetscViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);
270:       }
271:     }
272:     if (snes->lagpreconditioner == -1) {
273:       PetscViewerASCIIPrintf(viewer,"  Preconditioned is never rebuilt\n");
274:     } else if (snes->lagpreconditioner > 1) {
275:       PetscViewerASCIIPrintf(viewer,"  Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);
276:     }
277:     if (snes->lagjacobian == -1) {
278:       PetscViewerASCIIPrintf(viewer,"  Jacobian is never rebuilt\n");
279:     } else if (snes->lagjacobian > 1) {
280:       PetscViewerASCIIPrintf(viewer,"  Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);
281:     }
282:   } else if (isstring) {
283:     const char *type;
284:     SNESGetType(snes,&type);
285:     PetscViewerStringSPrintf(viewer," %-3.3s",type);
286:   } else if (isbinary) {
287:     PetscInt    classid = SNES_FILE_CLASSID;
288:     MPI_Comm    comm;
289:     PetscMPIInt rank;
290:     char        type[256];

292:     PetscObjectGetComm((PetscObject)snes,&comm);
293:     MPI_Comm_rank(comm,&rank);
294:     if (!rank) {
295:       PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
296:       PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));
297:       PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR,PETSC_FALSE);
298:     }
299:     if (snes->ops->view) {
300:       (*snes->ops->view)(snes,viewer);
301:     }
302:   } else if (isdraw) {
303:     PetscDraw draw;
304:     char      str[36];
305:     PetscReal x,y,bottom,h;

307:     PetscViewerDrawGetDraw(viewer,0,&draw);
308:     PetscDrawGetCurrentPoint(draw,&x,&y);
309:     PetscStrcpy(str,"SNES: ");
310:     PetscStrcat(str,((PetscObject)snes)->type_name);
311:     PetscDrawBoxedString(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);
312:     bottom = y - h;
313:     PetscDrawPushCurrentPoint(draw,x,bottom);
314:     if (snes->ops->view) {
315:       (*snes->ops->view)(snes,viewer);
316:     }
317: #if defined(PETSC_HAVE_SAWS)
318:   } else if (isams) {
319:     PetscMPIInt rank;
320:     const char *name;

322:     PetscObjectGetName((PetscObject)snes,&name);
323:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
324:     if (!((PetscObject)snes)->amsmem && !rank) {
325:       char       dir[1024];

327:       PetscObjectViewSAWs((PetscObject)snes,viewer);
328:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
329:       PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT));
330:       if (!snes->conv_hist) {
331:         SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);
332:       }
333:       PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);
334:       PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE));
335:     }
336: #endif
337:   }
338:   if (snes->linesearch) {
339:     PetscViewerASCIIPushTab(viewer);
340:     SNESGetLineSearch(snes, &linesearch);
341:     SNESLineSearchView(linesearch, viewer);
342:     PetscViewerASCIIPopTab(viewer);
343:   }
344:   if (snes->pc && snes->usespc) {
345:     PetscViewerASCIIPushTab(viewer);
346:     SNESView(snes->pc, viewer);
347:     PetscViewerASCIIPopTab(viewer);
348:   }
349:   PetscViewerASCIIPushTab(viewer);
350:   DMGetDMSNES(snes->dm,&dmsnes);
351:   DMSNESView(dmsnes, viewer);
352:   PetscViewerASCIIPopTab(viewer);
353:   if (snes->usesksp) {
354:     SNESGetKSP(snes,&ksp);
355:     PetscViewerASCIIPushTab(viewer);
356:     KSPView(ksp,viewer);
357:     PetscViewerASCIIPopTab(viewer);
358:   }
359:   if (isdraw) {
360:     PetscDraw draw;
361:     PetscViewerDrawGetDraw(viewer,0,&draw);
362:     PetscDrawPopCurrentPoint(draw);
363:   }
364:   return(0);
365: }

367: /*
368:   We retain a list of functions that also take SNES command
369:   line options. These are called at the end SNESSetFromOptions()
370: */
371: #define MAXSETFROMOPTIONS 5
372: static PetscInt numberofsetfromoptions;
373: static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);

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

380:   Not Collective

382:   Input Parameter:
383: . snescheck - function that checks for options

385:   Level: developer

387: .seealso: SNESSetFromOptions()
388: @*/
389: PetscErrorCode  SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
390: {
392:   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
393:   othersetfromoptions[numberofsetfromoptions++] = snescheck;
394:   return(0);
395: }

397: extern PetscErrorCode  SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);

401: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
402: {
403:   Mat            J;
404:   KSP            ksp;
405:   PC             pc;
406:   PetscBool      match;


412:   if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
413:     Mat A = snes->jacobian, B = snes->jacobian_pre;
414:     MatCreateVecs(A ? A : B, NULL,&snes->vec_func);
415:   }

417:   if (version == 1) {
418:     MatCreateSNESMF(snes,&J);
419:     MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
420:     MatSetFromOptions(J);
421:   } else if (version == 2) {
422:     if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
423: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
424:     SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);
425: #else
426:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)");
427: #endif
428:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2");

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

433:     /* This version replaces the user provided Jacobian matrix with a
434:        matrix-free version but still employs the user-provided preconditioner matrix. */
435:     SNESSetJacobian(snes,J,0,0,0);
436:   } else {
437:     /* This version replaces both the user-provided Jacobian and the user-
438:      provided preconditioner Jacobian with the default matrix free version. */
439:     if ((snes->pcside == PC_LEFT) && snes->pc) {
440:       if (!snes->jacobian){SNESSetJacobian(snes,J,0,0,0);}
441:     } else {
442:       SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,0);
443:     }
444:     /* Force no preconditioner */
445:     SNESGetKSP(snes,&ksp);
446:     KSPGetPC(ksp,&pc);
447:     PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);
448:     if (!match) {
449:       PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");
450:       PCSetType(pc,PCNONE);
451:     }
452:   }
453:   MatDestroy(&J);
454:   return(0);
455: }

459: static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx)
460: {
461:   SNES           snes = (SNES)ctx;
463:   Vec            Xfine,Xfine_named = NULL,Xcoarse;

466:   if (PetscLogPrintInfo) {
467:     PetscInt finelevel,coarselevel,fineclevel,coarseclevel;
468:     DMGetRefineLevel(dmfine,&finelevel);
469:     DMGetCoarsenLevel(dmfine,&fineclevel);
470:     DMGetRefineLevel(dmcoarse,&coarselevel);
471:     DMGetCoarsenLevel(dmcoarse,&coarseclevel);
472:     PetscInfo4(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);
473:   }
474:   if (dmfine == snes->dm) Xfine = snes->vec_sol;
475:   else {
476:     DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
477:     Xfine = Xfine_named;
478:   }
479:   DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
480:   MatRestrict(Restrict,Xfine,Xcoarse);
481:   VecPointwiseMult(Xcoarse,Xcoarse,Rscale);
482:   DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
483:   if (Xfine_named) {DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);}
484:   return(0);
485: }

489: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
490: {

494:   DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
495:   return(0);
496: }

500: /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can
501:  * safely call SNESGetDM() in their residual evaluation routine. */
502: static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx)
503: {
504:   SNES           snes = (SNES)ctx;
506:   Mat            Asave = A,Bsave = B;
507:   Vec            X,Xnamed = NULL;
508:   DM             dmsave;
509:   void           *ctxsave;
510:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);

513:   dmsave = snes->dm;
514:   KSPGetDM(ksp,&snes->dm);
515:   if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */
516:   else {                                     /* We are on a coarser level, this vec was initialized using a DM restrict hook */
517:     DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
518:     X    = Xnamed;
519:     SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);
520:     /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */
521:     if (jac == SNESComputeJacobianDefaultColor) {
522:       SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,0);
523:     }
524:   }
525:   /* put the previous context back */

527:   SNESComputeJacobian(snes,X,A,B);
528:   if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
529:     SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
530:   }

532:   if (A != Asave || B != Bsave) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for changing matrices at this time");
533:   if (Xnamed) {
534:     DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
535:   }
536:   snes->dm = dmsave;
537:   return(0);
538: }

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

545:    Collective

547:    Input Arguments:
548: .  snes - snes to configure

550:    Level: developer

552: .seealso: SNESSetUp()
553: @*/
554: PetscErrorCode SNESSetUpMatrices(SNES snes)
555: {
557:   DM             dm;
558:   DMSNES         sdm;

561:   SNESGetDM(snes,&dm);
562:   DMGetDMSNES(dm,&sdm);
563:   if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"DMSNES not properly configured");
564:   else if (!snes->jacobian && snes->mf) {
565:     Mat  J;
566:     void *functx;
567:     MatCreateSNESMF(snes,&J);
568:     MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
569:     MatSetFromOptions(J);
570:     SNESGetFunction(snes,NULL,NULL,&functx);
571:     SNESSetJacobian(snes,J,J,0,0);
572:     MatDestroy(&J);
573:   } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
574:     Mat J,B;
575:     MatCreateSNESMF(snes,&J);
576:     MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
577:     MatSetFromOptions(J);
578:     DMCreateMatrix(snes->dm,&B);
579:     /* sdm->computejacobian was already set to reach here */
580:     SNESSetJacobian(snes,J,B,NULL,NULL);
581:     MatDestroy(&J);
582:     MatDestroy(&B);
583:   } else if (!snes->jacobian_pre) {
584:     Mat J,B;
585:     J    = snes->jacobian;
586:     DMCreateMatrix(snes->dm,&B);
587:     SNESSetJacobian(snes,J ? J : B,B,NULL,NULL);
588:     MatDestroy(&B);
589:   }
590:   {
591:     KSP ksp;
592:     SNESGetKSP(snes,&ksp);
593:     KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);
594:     DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
595:   }
596:   return(0);
597: }

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

604:    Collective on SNES

606:    Input Parameter:
607: .  snes - the SNES context

609:    Options Database Keys:
610: +  -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
611: .  -snes_stol - convergence tolerance in terms of the norm
612:                 of the change in the solution between steps
613: .  -snes_atol <abstol> - absolute tolerance of residual norm
614: .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
615: .  -snes_max_it <max_it> - maximum number of iterations
616: .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
617: .  -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
618: .  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
619: .  -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
620: .  -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
621: .  -snes_trtol <trtol> - trust region tolerance
622: .  -snes_no_convergence_test - skip convergence test in nonlinear
623:                                solver; hence iterations will continue until max_it
624:                                or some other criterion is reached. Saves expense
625:                                of convergence test
626: .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
627:                                        filename given prints to stdout
628: .  -snes_monitor_solution - plots solution at each iteration
629: .  -snes_monitor_residual - plots residual (not its norm) at each iteration
630: .  -snes_monitor_solution_update - plots update to solution at each iteration
631: .  -snes_monitor_lg_residualnorm - plots residual norm at each iteration
632: .  -snes_monitor_lg_range - plots residual norm at each iteration
633: .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
634: .  -snes_fd_color - use finite differences with coloring to compute Jacobian
635: .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
636: -  -snes_converged_reason - print the reason for convergence/divergence after each solve

638:     Options Database for Eisenstat-Walker method:
639: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
640: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
641: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
642: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
643: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
644: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
645: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
646: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

648:    Notes:
649:    To see all options, run your program with the -help option or consult
650:    Users-Manual: ch_snes

652:    Level: beginner

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

656: .seealso: SNESSetOptionsPrefix()
657: @*/
658: PetscErrorCode  SNESSetFromOptions(SNES snes)
659: {
660:   PetscBool      flg,pcset,persist,set;
661:   PetscInt       i,indx,lag,grids;
662:   const char     *deft        = SNESNEWTONLS;
663:   const char     *convtests[] = {"default","skip"};
664:   SNESKSPEW      *kctx        = NULL;
665:   char           type[256], monfilename[PETSC_MAX_PATH_LEN];
666:   PetscViewer    monviewer;
668:   PCSide         pcside;
669:   const char     *optionsprefix;

673:   if (!SNESRegisterAllCalled) {SNESRegisterAll();}
674:   PetscObjectOptionsBegin((PetscObject)snes);
675:   if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
676:   PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
677:   if (flg) {
678:     SNESSetType(snes,type);
679:   } else if (!((PetscObject)snes)->type_name) {
680:     SNESSetType(snes,deft);
681:   }
682:   PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);
683:   PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);

685:   PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,NULL);
686:   PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,NULL);
687:   PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,NULL);
688:   PetscOptionsInt("-snes_max_fail","Maximum nonlinear step failures","SNESSetMaxNonlinearStepFailures",snes->maxFailures,&snes->maxFailures,NULL);
689:   PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,NULL);
690:   PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,NULL);

692:   PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
693:   if (flg) {
694:     SNESSetLagPreconditioner(snes,lag);
695:   }
696:   PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
697:   if (flg) {
698:     SNESSetLagPreconditionerPersists(snes,persist);
699:   }
700:   PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
701:   if (flg) {
702:     SNESSetLagJacobian(snes,lag);
703:   }
704:   PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
705:   if (flg) {
706:     SNESSetLagJacobianPersists(snes,persist);
707:   }

709:   PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
710:   if (flg) {
711:     SNESSetGridSequence(snes,grids);
712:   }

714:   PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);
715:   if (flg) {
716:     switch (indx) {
717:     case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
718:     case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL);    break;
719:     }
720:   }

722:   PetscOptionsBool("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",snes->printreason,&snes->printreason,NULL);

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

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

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

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

734:   PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);
735:   PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);
736:   PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);
737:   PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);
738:   PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);
739:   PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);
740:   PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);

742:   flg  = PETSC_FALSE;
743:   PetscOptionsBool("-snes_check_jacobian","Check each Jacobian with a differenced one","SNESUpdateCheckJacobian",flg,&flg,&set);
744:   if (set && flg) {
745:     SNESSetUpdate(snes,SNESUpdateCheckJacobian);
746:   }

748:   flg  = PETSC_FALSE;
749:   PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);
750:   if (set && flg) {SNESMonitorCancel(snes);}

752:   PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
753:   if (flg) {
754:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
755:     SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
756:   }

758:   PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
759:   if (flg) {
760:     SNESMonitorSet(snes,SNESMonitorRange,0,0);
761:   }

763:   PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
764:   if (flg) {
765:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
766:     SNESMonitorSetRatio(snes,monviewer);
767:   }

769:   PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
770:   if (flg) {
771:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
772:     SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
773:   }

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

778:   flg  = PETSC_FALSE;
779:   PetscOptionsBool("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,NULL);
780:   if (flg) {SNESMonitorSet(snes,SNESMonitorSolution,0,0);}
781:   flg  = PETSC_FALSE;
782:   PetscOptionsBool("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,NULL);
783:   if (flg) {SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);}
784:   flg  = PETSC_FALSE;
785:   PetscOptionsBool("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,NULL);
786:   if (flg) {SNESMonitorSet(snes,SNESMonitorResidual,0,0);}
787:   flg  = PETSC_FALSE;
788:   PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);
789:   if (flg) {
790:     PetscDrawLG ctx;

792:     SNESMonitorLGCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
793:     SNESMonitorSet(snes,SNESMonitorLGResidualNorm,ctx,(PetscErrorCode (*)(void**))SNESMonitorLGDestroy);
794:   }
795:   flg  = PETSC_FALSE;
796:   PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
797:   if (flg) {
798:     PetscViewer ctx;

800:     PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
801:     SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
802:   }

804:   flg  = PETSC_FALSE;
805:   PetscOptionsBool("-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",flg,&flg,NULL);
806:   if (flg) {SNESMonitorSet(snes,SNESMonitorJacUpdateSpectrum,0,0);}


809:   PetscOptionsString("-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
810:   if (flg) {
811:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
812:     SNESMonitorSet(snes,SNESMonitorFields,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
813:   }

815:   flg  = PETSC_FALSE;
816:   PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
817:   if (flg) {
818:     void *functx;
819:     SNESGetFunction(snes,NULL,NULL,&functx);
820:     SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
821:     PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
822:   }

824:   flg  = PETSC_FALSE;
825:   PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
826:   if (flg) {
827:     SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
828:   }

830:   flg  = PETSC_FALSE;
831:   PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);
832:   if (flg) {
833:     DM             dm;
834:     DMSNES         sdm;
835:     SNESGetDM(snes,&dm);
836:     DMGetDMSNES(dm,&sdm);
837:     sdm->jacobianctx = NULL;
838:     SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,0);
839:     PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
840:   }

842:   flg  = PETSC_FALSE;
843:   PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&snes->mf_operator,&flg);
844:   if (flg && snes->mf_operator) {
845:     snes->mf_operator = PETSC_TRUE;
846:     snes->mf          = PETSC_TRUE;
847:   }
848:   flg  = PETSC_FALSE;
849:   PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&snes->mf,&flg);
850:   if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE;
851:   PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,0);

853:   flg  = PETSC_FALSE;
854:   SNESGetNPCSide(snes,&pcside);
855:   PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
856:   if (flg) {SNESSetNPCSide(snes,pcside);}

858: #if defined(PETSC_HAVE_SAWS)
859:   /*
860:     Publish convergence information using SAWs
861:   */
862:   flg  = PETSC_FALSE;
863:   PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);
864:   if (flg) {
865:     void *ctx;
866:     SNESMonitorSAWsCreate(snes,&ctx);
867:     SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);
868:   }
869: #endif
870: #if defined(PETSC_HAVE_SAWS)
871:   {
872:   PetscBool set;
873:   flg  = PETSC_FALSE;
874:   PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);
875:   if (set) {
876:     PetscObjectSAWsSetBlock((PetscObject)snes,flg);
877:   }
878:   }
879: #endif

881:   for (i = 0; i < numberofsetfromoptions; i++) {
882:     (*othersetfromoptions[i])(snes);
883:   }

885:   if (snes->ops->setfromoptions) {
886:     (*snes->ops->setfromoptions)(snes);
887:   }

889:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
890:   PetscObjectProcessOptionsHandlers((PetscObject)snes);
891:   PetscOptionsEnd();

893:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
894:   KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
895:   KSPSetFromOptions(snes->ksp);

897:   if (!snes->linesearch) {
898:     SNESGetLineSearch(snes, &snes->linesearch);
899:   }
900:   SNESLineSearchSetFromOptions(snes->linesearch);

902:   /* if someone has set the SNES NPC type, create it. */
903:   SNESGetOptionsPrefix(snes, &optionsprefix);
904:   PetscOptionsHasName(optionsprefix, "-npc_snes_type", &pcset);
905:   if (pcset && (!snes->pc)) {
906:     SNESGetNPC(snes, &snes->pc);
907:   }
908:   return(0);
909: }

913: /*@C
914:    SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
915:    the nonlinear solvers.

917:    Logically Collective on SNES

919:    Input Parameters:
920: +  snes - the SNES context
921: .  compute - function to compute the context
922: -  destroy - function to destroy the context

924:    Level: intermediate

926:    Notes:
927:    This function is currently not available from Fortran.

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

931: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
932: @*/
933: PetscErrorCode  SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
934: {
937:   snes->ops->usercompute = compute;
938:   snes->ops->userdestroy = destroy;
939:   return(0);
940: }

944: /*@
945:    SNESSetApplicationContext - Sets the optional user-defined context for
946:    the nonlinear solvers.

948:    Logically Collective on SNES

950:    Input Parameters:
951: +  snes - the SNES context
952: -  usrP - optional user context

954:    Level: intermediate

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

958: .seealso: SNESGetApplicationContext()
959: @*/
960: PetscErrorCode  SNESSetApplicationContext(SNES snes,void *usrP)
961: {
963:   KSP            ksp;

967:   SNESGetKSP(snes,&ksp);
968:   KSPSetApplicationContext(ksp,usrP);
969:   snes->user = usrP;
970:   return(0);
971: }

975: /*@
976:    SNESGetApplicationContext - Gets the user-defined context for the
977:    nonlinear solvers.

979:    Not Collective

981:    Input Parameter:
982: .  snes - SNES context

984:    Output Parameter:
985: .  usrP - user context

987:    Level: intermediate

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

991: .seealso: SNESSetApplicationContext()
992: @*/
993: PetscErrorCode  SNESGetApplicationContext(SNES snes,void *usrP)
994: {
997:   *(void**)usrP = snes->user;
998:   return(0);
999: }

1003: /*@
1004:    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1005:    at this time.

1007:    Not Collective

1009:    Input Parameter:
1010: .  snes - SNES context

1012:    Output Parameter:
1013: .  iter - iteration number

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

1018:    This is useful for using lagged Jacobians (where one does not recompute the
1019:    Jacobian at each SNES iteration). For example, the code
1020: .vb
1021:       SNESGetIterationNumber(snes,&it);
1022:       if (!(it % 2)) {
1023:         [compute Jacobian here]
1024:       }
1025: .ve
1026:    can be used in your ComputeJacobian() function to cause the Jacobian to be
1027:    recomputed every second SNES iteration.

1029:    Level: intermediate

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

1033: .seealso:   SNESGetLinearSolveIterations()
1034: @*/
1035: PetscErrorCode  SNESGetIterationNumber(SNES snes,PetscInt *iter)
1036: {
1040:   *iter = snes->iter;
1041:   return(0);
1042: }

1046: /*@
1047:    SNESSetIterationNumber - Sets the current iteration number.

1049:    Not Collective

1051:    Input Parameter:
1052: .  snes - SNES context
1053: .  iter - iteration number

1055:    Level: developer

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

1059: .seealso:   SNESGetLinearSolveIterations()
1060: @*/
1061: PetscErrorCode  SNESSetIterationNumber(SNES snes,PetscInt iter)
1062: {

1067:   PetscObjectSAWsTakeAccess((PetscObject)snes);
1068:   snes->iter = iter;
1069:   PetscObjectSAWsGrantAccess((PetscObject)snes);
1070:   return(0);
1071: }

1075: /*@
1076:    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1077:    attempted by the nonlinear solver.

1079:    Not Collective

1081:    Input Parameter:
1082: .  snes - SNES context

1084:    Output Parameter:
1085: .  nfails - number of unsuccessful steps attempted

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

1090:    Level: intermediate

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

1094: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1095:           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1096: @*/
1097: PetscErrorCode  SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1098: {
1102:   *nfails = snes->numFailures;
1103:   return(0);
1104: }

1108: /*@
1109:    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1110:    attempted by the nonlinear solver before it gives up.

1112:    Not Collective

1114:    Input Parameters:
1115: +  snes     - SNES context
1116: -  maxFails - maximum of unsuccessful steps

1118:    Level: intermediate

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

1122: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1123:           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1124: @*/
1125: PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1126: {
1129:   snes->maxFailures = maxFails;
1130:   return(0);
1131: }

1135: /*@
1136:    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1137:    attempted by the nonlinear solver before it gives up.

1139:    Not Collective

1141:    Input Parameter:
1142: .  snes     - SNES context

1144:    Output Parameter:
1145: .  maxFails - maximum of unsuccessful steps

1147:    Level: intermediate

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

1151: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1152:           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()

1154: @*/
1155: PetscErrorCode  SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1156: {
1160:   *maxFails = snes->maxFailures;
1161:   return(0);
1162: }

1166: /*@
1167:    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1168:      done by SNES.

1170:    Not Collective

1172:    Input Parameter:
1173: .  snes     - SNES context

1175:    Output Parameter:
1176: .  nfuncs - number of evaluations

1178:    Level: intermediate

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

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

1184: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1185: @*/
1186: PetscErrorCode  SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1187: {
1191:   *nfuncs = snes->nfuncs;
1192:   return(0);
1193: }

1197: /*@
1198:    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1199:    linear solvers.

1201:    Not Collective

1203:    Input Parameter:
1204: .  snes - SNES context

1206:    Output Parameter:
1207: .  nfails - number of failed solves

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

1212:    Level: intermediate

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

1216: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1217: @*/
1218: PetscErrorCode  SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1219: {
1223:   *nfails = snes->numLinearSolveFailures;
1224:   return(0);
1225: }

1229: /*@
1230:    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1231:    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE

1233:    Logically Collective on SNES

1235:    Input Parameters:
1236: +  snes     - SNES context
1237: -  maxFails - maximum allowed linear solve failures

1239:    Level: intermediate

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

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

1245: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1246: @*/
1247: PetscErrorCode  SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1248: {
1252:   snes->maxLinearSolveFailures = maxFails;
1253:   return(0);
1254: }

1258: /*@
1259:    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1260:      are allowed before SNES terminates

1262:    Not Collective

1264:    Input Parameter:
1265: .  snes     - SNES context

1267:    Output Parameter:
1268: .  maxFails - maximum of unsuccessful solves allowed

1270:    Level: intermediate

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

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

1276: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1277: @*/
1278: PetscErrorCode  SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1279: {
1283:   *maxFails = snes->maxLinearSolveFailures;
1284:   return(0);
1285: }

1289: /*@
1290:    SNESGetLinearSolveIterations - Gets the total number of linear iterations
1291:    used by the nonlinear solver.

1293:    Not Collective

1295:    Input Parameter:
1296: .  snes - SNES context

1298:    Output Parameter:
1299: .  lits - number of linear iterations

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

1304:    Level: intermediate

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

1308: .seealso:  SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1309: @*/
1310: PetscErrorCode  SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1311: {
1315:   *lits = snes->linear_its;
1316:   return(0);
1317: }

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

1325:    Logically Collective on SNES

1327:    Input Parameter:
1328: +  snes - SNES context
1329: -  reset - whether to reset the counters or not

1331:    Notes:
1332:    This is automatically called with FALSE

1334:    Level: developer

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

1338: .seealso:  SNESGetNumberFunctionEvals(), SNESGetNumberLinearSolveIterations(), SNESGetNPC()
1339: @*/
1340: PetscErrorCode  SNESSetCountersReset(SNES snes,PetscBool reset)
1341: {
1345:   snes->counters_reset = reset;
1346:   return(0);
1347: }


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

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

1357:    Input Parameters:
1358: +  snes - the SNES context
1359: -  ksp - the KSP context

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

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

1368:    Level: developer

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

1372: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1373: @*/
1374: PetscErrorCode  SNESSetKSP(SNES snes,KSP ksp)
1375: {

1382:   PetscObjectReference((PetscObject)ksp);
1383:   if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1384:   snes->ksp = ksp;
1385:   return(0);
1386: }

1388: /* -----------------------------------------------------------*/
1391: /*@
1392:    SNESCreate - Creates a nonlinear solver context.

1394:    Collective on MPI_Comm

1396:    Input Parameters:
1397: .  comm - MPI communicator

1399:    Output Parameter:
1400: .  outsnes - the new SNES context

1402:    Options Database Keys:
1403: +   -snes_mf - Activates default matrix-free Jacobian-vector products,
1404:                and no preconditioning matrix
1405: .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
1406:                products, and a user-provided preconditioning matrix
1407:                as set by SNESSetJacobian()
1408: -   -snes_fd - Uses (slow!) finite differences to compute Jacobian

1410:    Level: beginner

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

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

1416: @*/
1417: PetscErrorCode  SNESCreate(MPI_Comm comm,SNES *outsnes)
1418: {
1420:   SNES           snes;
1421:   SNESKSPEW      *kctx;

1425:   *outsnes = NULL;
1426:   SNESInitializePackage();

1428:   PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_CLASSID,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);

1430:   snes->ops->converged    = SNESConvergedDefault;
1431:   snes->usesksp           = PETSC_TRUE;
1432:   snes->tolerancesset     = PETSC_FALSE;
1433:   snes->max_its           = 50;
1434:   snes->max_funcs         = 10000;
1435:   snes->norm              = 0.0;
1436:   snes->normschedule      = SNES_NORM_ALWAYS;
1437:   snes->functype          = SNES_FUNCTION_DEFAULT;
1438: #if defined(PETSC_USE_REAL_SINGLE)
1439:   snes->rtol              = 1.e-5;
1440: #else
1441:   snes->rtol              = 1.e-8;
1442: #endif
1443:   snes->ttol              = 0.0;
1444: #if defined(PETSC_USE_REAL_SINGLE)
1445:   snes->abstol            = 1.e-25;
1446: #else
1447:   snes->abstol            = 1.e-50;
1448: #endif
1449:   snes->stol              = 1.e-8;
1450: #if defined(PETSC_USE_REAL_SINGLE)
1451:   snes->deltatol          = 1.e-6;
1452: #else
1453:   snes->deltatol          = 1.e-12;
1454: #endif
1455:   snes->nfuncs            = 0;
1456:   snes->numFailures       = 0;
1457:   snes->maxFailures       = 1;
1458:   snes->linear_its        = 0;
1459:   snes->lagjacobian       = 1;
1460:   snes->jac_iter          = 0;
1461:   snes->lagjac_persist    = PETSC_FALSE;
1462:   snes->lagpreconditioner = 1;
1463:   snes->pre_iter          = 0;
1464:   snes->lagpre_persist    = PETSC_FALSE;
1465:   snes->numbermonitors    = 0;
1466:   snes->data              = 0;
1467:   snes->setupcalled       = PETSC_FALSE;
1468:   snes->ksp_ewconv        = PETSC_FALSE;
1469:   snes->nwork             = 0;
1470:   snes->work              = 0;
1471:   snes->nvwork            = 0;
1472:   snes->vwork             = 0;
1473:   snes->conv_hist_len     = 0;
1474:   snes->conv_hist_max     = 0;
1475:   snes->conv_hist         = NULL;
1476:   snes->conv_hist_its     = NULL;
1477:   snes->conv_hist_reset   = PETSC_TRUE;
1478:   snes->counters_reset    = PETSC_TRUE;
1479:   snes->vec_func_init_set = PETSC_FALSE;
1480:   snes->reason            = SNES_CONVERGED_ITERATING;
1481:   snes->pcside            = PC_RIGHT;

1483:   snes->mf          = PETSC_FALSE;
1484:   snes->mf_operator = PETSC_FALSE;
1485:   snes->mf_version  = 1;

1487:   snes->numLinearSolveFailures = 0;
1488:   snes->maxLinearSolveFailures = 1;

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

1493:   snes->kspconvctx  = (void*)kctx;
1494:   kctx->version     = 2;
1495:   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1496:                              this was too large for some test cases */
1497:   kctx->rtol_last   = 0.0;
1498:   kctx->rtol_max    = .9;
1499:   kctx->gamma       = 1.0;
1500:   kctx->alpha       = .5*(1.0 + PetscSqrtReal(5.0));
1501:   kctx->alpha2      = kctx->alpha;
1502:   kctx->threshold   = .1;
1503:   kctx->lresid_last = 0.0;
1504:   kctx->norm_last   = 0.0;

1506:   *outsnes = snes;
1507:   return(0);
1508: }

1510: /*MC
1511:     SNESFunction - functional form used to convey the nonlinear function to be solved by SNES

1513:      Synopsis:
1514:      #include <petscsnes.h>
1515:      SNESFunction(SNES snes,Vec x,Vec f,void *ctx);

1517:      Input Parameters:
1518: +     snes - the SNES context
1519: .     x    - state at which to evaluate residual
1520: -     ctx     - optional user-defined function context, passed in with SNESSetFunction()

1522:      Output Parameter:
1523: .     f  - vector to put residual (function value)

1525:    Level: intermediate

1527: .seealso:   SNESSetFunction(), SNESGetFunction()
1528: M*/

1532: /*@C
1533:    SNESSetFunction - Sets the function evaluation routine and function
1534:    vector for use by the SNES routines in solving systems of nonlinear
1535:    equations.

1537:    Logically Collective on SNES

1539:    Input Parameters:
1540: +  snes - the SNES context
1541: .  r - vector to store function value
1542: .  f - function evaluation routine; see SNESFunction for calling sequence details
1543: -  ctx - [optional] user-defined context for private data for the
1544:          function evaluation routine (may be NULL)

1546:    Notes:
1547:    The Newton-like methods typically solve linear systems of the form
1548: $      f'(x) x = -f(x),
1549:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

1551:    Level: beginner

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

1555: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1556: @*/
1557: PetscErrorCode  SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1558: {
1560:   DM             dm;

1564:   if (r) {
1567:     PetscObjectReference((PetscObject)r);
1568:     VecDestroy(&snes->vec_func);

1570:     snes->vec_func = r;
1571:   }
1572:   SNESGetDM(snes,&dm);
1573:   DMSNESSetFunction(dm,f,ctx);
1574:   return(0);
1575: }


1580: /*@C
1581:    SNESSetInitialFunction - Sets the function vector to be used as the
1582:    function norm at the initialization of the method.  In some
1583:    instances, the user has precomputed the function before calling
1584:    SNESSolve.  This function allows one to avoid a redundant call
1585:    to SNESComputeFunction in that case.

1587:    Logically Collective on SNES

1589:    Input Parameters:
1590: +  snes - the SNES context
1591: -  f - vector to store function value

1593:    Notes:
1594:    This should not be modified during the solution procedure.

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

1598:    Level: developer

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

1602: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1603: @*/
1604: PetscErrorCode  SNESSetInitialFunction(SNES snes, Vec f)
1605: {
1607:   Vec            vec_func;

1613:   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1614:     snes->vec_func_init_set = PETSC_FALSE;
1615:     return(0);
1616:   }
1617:   SNESGetFunction(snes,&vec_func,NULL,NULL);
1618:   VecCopy(f, vec_func);

1620:   snes->vec_func_init_set = PETSC_TRUE;
1621:   return(0);
1622: }

1626: /*@
1627:    SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1628:    of the SNES method.

1630:    Logically Collective on SNES

1632:    Input Parameters:
1633: +  snes - the SNES context
1634: -  normschedule - the frequency of norm computation

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

1645:    Level: developer

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

1649: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1650: @*/
1651: PetscErrorCode  SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1652: {
1655:   snes->normschedule = normschedule;
1656:   return(0);
1657: }


1662: /*@
1663:    SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1664:    of the SNES method.

1666:    Logically Collective on SNES

1668:    Input Parameters:
1669: +  snes - the SNES context
1670: -  normschedule - the type of the norm used

1672:    Level: advanced

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

1676: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1677: @*/
1678: PetscErrorCode  SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1679: {
1682:   *normschedule = snes->normschedule;
1683:   return(0);
1684: }


1689: /*@C
1690:    SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
1691:    of the SNES method.

1693:    Logically Collective on SNES

1695:    Input Parameters:
1696: +  snes - the SNES context
1697: -  normschedule - the frequency of norm computation

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

1708:    Level: developer

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

1712: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1713: @*/
1714: PetscErrorCode  SNESSetFunctionType(SNES snes, SNESFunctionType type)
1715: {
1718:   snes->functype = type;
1719:   return(0);
1720: }


1725: /*@C
1726:    SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
1727:    of the SNES method.

1729:    Logically Collective on SNES

1731:    Input Parameters:
1732: +  snes - the SNES context
1733: -  normschedule - the type of the norm used

1735:    Level: advanced

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

1739: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1740: @*/
1741: PetscErrorCode  SNESGetFunctionType(SNES snes, SNESFunctionType *type)
1742: {
1745:   *type = snes->functype;
1746:   return(0);
1747: }

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

1752:      Synopsis:
1753:      #include <petscsnes.h>
1754: $    SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);

1756: +  X   - solution vector
1757: .  B   - RHS vector
1758: -  ctx - optional user-defined Gauss-Seidel context

1760:    Level: intermediate

1762: .seealso:   SNESSetNGS(), SNESGetNGS()
1763: M*/

1767: /*@C
1768:    SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
1769:    use with composed nonlinear solvers.

1771:    Input Parameters:
1772: +  snes   - the SNES context
1773: .  f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
1774: -  ctx    - [optional] user-defined context for private data for the
1775:             smoother evaluation routine (may be NULL)

1777:    Notes:
1778:    The NGS routines are used by the composed nonlinear solver to generate
1779:     a problem appropriate update to the solution, particularly FAS.

1781:    Level: intermediate

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

1785: .seealso: SNESGetFunction(), SNESComputeNGS()
1786: @*/
1787: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1788: {
1790:   DM             dm;

1794:   SNESGetDM(snes,&dm);
1795:   DMSNESSetNGS(dm,f,ctx);
1796:   return(0);
1797: }

1801: PETSC_EXTERN PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
1802: {
1804:   DM             dm;
1805:   DMSNES         sdm;

1808:   SNESGetDM(snes,&dm);
1809:   DMGetDMSNES(dm,&sdm);
1810:   /*  A(x)*x - b(x) */
1811:   if (sdm->ops->computepfunction) {
1812:     (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
1813:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");

1815:   if (sdm->ops->computepjacobian) {
1816:     (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
1817:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard matrix.");
1818:   VecScale(f,-1.0);
1819:   MatMultAdd(snes->jacobian,x,f,f);
1820:   return(0);
1821: }

1825: PETSC_EXTERN PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
1826: {
1828:   /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
1829:   return(0);
1830: }

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

1837:    Logically Collective on SNES

1839:    Input Parameters:
1840: +  snes - the SNES context
1841: .  r - vector to store function value
1842: .  b - function evaluation routine
1843: .  Amat - matrix with which A(x) x - b(x) is to be computed
1844: .  Pmat - matrix from which preconditioner is computed (usually the same as Amat)
1845: .  J  - function to compute matrix value
1846: -  ctx - [optional] user-defined context for private data for the
1847:          function evaluation routine (may be NULL)

1849:    Notes:
1850:     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
1851:     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.

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

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

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

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

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

1867:    Level: intermediate

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

1871: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard()
1872: @*/
1873: 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)
1874: {
1876:   DM             dm;

1880:   SNESGetDM(snes, &dm);
1881:   DMSNESSetPicard(dm,b,J,ctx);
1882:   SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
1883:   SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
1884:   return(0);
1885: }

1889: /*@C
1890:    SNESGetPicard - Returns the context for the Picard iteration

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

1894:    Input Parameter:
1895: .  snes - the SNES context

1897:    Output Parameter:
1898: +  r - the function (or NULL)
1899: .  f - the function (or NULL); see SNESFunction for calling sequence details
1900: .  Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
1901: .  Pmat  - the matrix from which the preconditioner will be constructed (or NULL)
1902: .  J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
1903: -  ctx - the function context (or NULL)

1905:    Level: advanced

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

1909: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
1910: @*/
1911: 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)
1912: {
1914:   DM             dm;

1918:   SNESGetFunction(snes,r,NULL,NULL);
1919:   SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
1920:   SNESGetDM(snes,&dm);
1921:   DMSNESGetPicard(dm,f,J,ctx);
1922:   return(0);
1923: }

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

1930:    Logically Collective on SNES

1932:    Input Parameters:
1933: +  snes - the SNES context
1934: .  func - function evaluation routine
1935: -  ctx - [optional] user-defined context for private data for the
1936:          function evaluation routine (may be NULL)

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

1941: .  f - function vector
1942: -  ctx - optional user-defined function context

1944:    Level: intermediate

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

1948: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
1949: @*/
1950: PetscErrorCode  SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
1951: {
1954:   if (func) snes->ops->computeinitialguess = func;
1955:   if (ctx)  snes->initialguessP            = ctx;
1956:   return(0);
1957: }

1959: /* --------------------------------------------------------------- */
1962: /*@C
1963:    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
1964:    it assumes a zero right hand side.

1966:    Logically Collective on SNES

1968:    Input Parameter:
1969: .  snes - the SNES context

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

1974:    Level: intermediate

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

1978: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
1979: @*/
1980: PetscErrorCode  SNESGetRhs(SNES snes,Vec *rhs)
1981: {
1985:   *rhs = snes->vec_rhs;
1986:   return(0);
1987: }

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

1994:    Collective on SNES

1996:    Input Parameters:
1997: +  snes - the SNES context
1998: -  x - input vector

2000:    Output Parameter:
2001: .  y - function vector, as set by SNESSetFunction()

2003:    Notes:
2004:    SNESComputeFunction() is typically used within nonlinear solvers
2005:    implementations, so most users would not generally call this routine
2006:    themselves.

2008:    Level: developer

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

2012: .seealso: SNESSetFunction(), SNESGetFunction()
2013: @*/
2014: PetscErrorCode  SNESComputeFunction(SNES snes,Vec x,Vec y)
2015: {
2017:   DM             dm;
2018:   DMSNES         sdm;

2026:   VecValidValues(x,2,PETSC_TRUE);

2028:   SNESGetDM(snes,&dm);
2029:   DMGetDMSNES(dm,&sdm);
2030:   if (sdm->ops->computefunction) {
2031:     PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2032:     PetscStackPush("SNES user function");
2033:     (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2034:     PetscStackPop;
2035:     PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2036:   } else if (snes->vec_rhs) {
2037:     MatMult(snes->jacobian, x, y);
2038:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2039:   if (snes->vec_rhs) {
2040:     VecAXPY(y,-1.0,snes->vec_rhs);
2041:   }
2042:   snes->nfuncs++;
2043:   VecValidValues(y,3,PETSC_FALSE);
2044:   return(0);
2045: }

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

2052:    Collective on SNES

2054:    Input Parameters:
2055: +  snes - the SNES context
2056: .  x - input vector
2057: -  b - rhs vector

2059:    Output Parameter:
2060: .  x - new solution vector

2062:    Notes:
2063:    SNESComputeNGS() is typically used within composed nonlinear solver
2064:    implementations, so most users would not generally call this routine
2065:    themselves.

2067:    Level: developer

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

2071: .seealso: SNESSetNGS(), SNESComputeFunction()
2072: @*/
2073: PetscErrorCode  SNESComputeNGS(SNES snes,Vec b,Vec x)
2074: {
2076:   DM             dm;
2077:   DMSNES         sdm;

2085:   if (b) {VecValidValues(b,2,PETSC_TRUE);}
2086:   PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2087:   SNESGetDM(snes,&dm);
2088:   DMGetDMSNES(dm,&sdm);
2089:   if (sdm->ops->computegs) {
2090:     PetscStackPush("SNES user NGS");
2091:     (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2092:     PetscStackPop;
2093:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2094:   PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2095:   VecValidValues(x,3,PETSC_FALSE);
2096:   return(0);
2097: }

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

2104:    Collective on SNES and Mat

2106:    Input Parameters:
2107: +  snes - the SNES context
2108: -  x - input vector

2110:    Output Parameters:
2111: +  A - Jacobian matrix
2112: -  B - optional preconditioning matrix

2114:   Options Database Keys:
2115: +    -snes_lag_preconditioner <lag>
2116: .    -snes_lag_jacobian <lag>
2117: .    -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2118: .    -snes_compare_explicit_draw  - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2119: .    -snes_compare_explicit_contour  - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2120: .    -snes_compare_operator  - Make the comparison options above use the operator instead of the preconditioning matrix
2121: .    -snes_compare_coloring - Compute the finite differece Jacobian using coloring and display norms of difference
2122: .    -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
2123: .    -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2124: .    -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2125: .    -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2126: .    -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
2127: -    -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences


2130:    Notes:
2131:    Most users should not need to explicitly call this routine, as it
2132:    is used internally within the nonlinear solvers.

2134:    Level: developer

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

2138: .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2139: @*/
2140: PetscErrorCode  SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2141: {
2143:   PetscBool      flag;
2144:   DM             dm;
2145:   DMSNES         sdm;
2146:   KSP            ksp;

2152:   VecValidValues(X,2,PETSC_TRUE);
2153:   SNESGetDM(snes,&dm);
2154:   DMGetDMSNES(dm,&sdm);

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

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

2160:   if (snes->lagjacobian == -2) {
2161:     snes->lagjacobian = -1;

2163:     PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2164:   } else if (snes->lagjacobian == -1) {
2165:     PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2166:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2167:     if (flag) {
2168:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2169:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2170:     }
2171:     return(0);
2172:   } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2173:     PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2174:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2175:     if (flag) {
2176:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2177:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2178:     }
2179:     return(0);
2180:   }
2181:   if (snes->pc && snes->pcside == PC_LEFT) {
2182:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2183:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2184:       return(0);
2185:   }

2187:   PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);

2189:   PetscStackPush("SNES user Jacobian function");
2190:   (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2191:   PetscStackPop;

2193:   PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);

2195:   /* the next line ensures that snes->ksp exists */
2196:   SNESGetKSP(snes,&ksp);
2197:   if (snes->lagpreconditioner == -2) {
2198:     PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2199:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2200:     snes->lagpreconditioner = -1;
2201:   } else if (snes->lagpreconditioner == -1) {
2202:     PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2203:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2204:   } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2205:     PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2206:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2207:   } else {
2208:     PetscInfo(snes,"Rebuilding preconditioner\n");
2209:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2210:   }

2212:   /* make sure user returned a correct Jacobian and preconditioner */
2215:   {
2216:     PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
2217:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit",&flag,NULL);
2218:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",&flag_draw,NULL);
2219:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",&flag_contour,NULL);
2220:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_operator",&flag_operator,NULL);
2221:     if (flag || flag_draw || flag_contour) {
2222:       Mat          Bexp_mine = NULL,Bexp,FDexp;
2223:       PetscViewer  vdraw,vstdout;
2224:       PetscBool    flg;
2225:       if (flag_operator) {
2226:         MatComputeExplicitOperator(A,&Bexp_mine);
2227:         Bexp = Bexp_mine;
2228:       } else {
2229:         /* See if the preconditioning matrix can be viewed and added directly */
2230:         PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2231:         if (flg) Bexp = B;
2232:         else {
2233:           /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
2234:           MatComputeExplicitOperator(B,&Bexp_mine);
2235:           Bexp = Bexp_mine;
2236:         }
2237:       }
2238:       MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);
2239:       SNESComputeJacobianDefault(snes,X,FDexp,FDexp,NULL);
2240:       PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2241:       if (flag_draw || flag_contour) {
2242:         PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2243:         if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2244:       } else vdraw = NULL;
2245:       PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator ? "Jacobian" : "preconditioning Jacobian");
2246:       if (flag) {MatView(Bexp,vstdout);}
2247:       if (vdraw) {MatView(Bexp,vdraw);}
2248:       PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");
2249:       if (flag) {MatView(FDexp,vstdout);}
2250:       if (vdraw) {MatView(FDexp,vdraw);}
2251:       MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);
2252:       PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");
2253:       if (flag) {MatView(FDexp,vstdout);}
2254:       if (vdraw) {              /* Always use contour for the difference */
2255:         PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2256:         MatView(FDexp,vdraw);
2257:         PetscViewerPopFormat(vdraw);
2258:       }
2259:       if (flag_contour) {PetscViewerPopFormat(vdraw);}
2260:       PetscViewerDestroy(&vdraw);
2261:       MatDestroy(&Bexp_mine);
2262:       MatDestroy(&FDexp);
2263:     }
2264:   }
2265:   {
2266:     PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
2267:     PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
2268:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring",&flag,NULL);
2269:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_display",&flag_display,NULL);
2270:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",&flag_draw,NULL);
2271:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",&flag_contour,NULL);
2272:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",&flag_threshold,NULL);
2273:     PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,NULL);
2274:     PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,NULL);
2275:     if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
2276:       Mat            Bfd;
2277:       PetscViewer    vdraw,vstdout;
2278:       MatColoring    coloring;
2279:       ISColoring     iscoloring;
2280:       MatFDColoring  matfdcoloring;
2281:       PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2282:       void           *funcctx;
2283:       PetscReal      norm1,norm2,normmax;

2285:       MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2286:       MatColoringCreate(Bfd,&coloring);
2287:       MatColoringSetType(coloring,MATCOLORINGSL);
2288:       MatColoringSetFromOptions(coloring);
2289:       MatColoringApply(coloring,&iscoloring);
2290:       MatColoringDestroy(&coloring);
2291:       MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2292:       MatFDColoringSetFromOptions(matfdcoloring);
2293:       MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2294:       ISColoringDestroy(&iscoloring);

2296:       /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2297:       SNESGetFunction(snes,NULL,&func,&funcctx);
2298:       MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2299:       PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2300:       PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2301:       MatFDColoringSetFromOptions(matfdcoloring);
2302:       MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2303:       MatFDColoringDestroy(&matfdcoloring);

2305:       PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2306:       if (flag_draw || flag_contour) {
2307:         PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2308:         if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2309:       } else vdraw = NULL;
2310:       PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2311:       if (flag_display) {MatView(B,vstdout);}
2312:       if (vdraw) {MatView(B,vdraw);}
2313:       PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2314:       if (flag_display) {MatView(Bfd,vstdout);}
2315:       if (vdraw) {MatView(Bfd,vdraw);}
2316:       MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2317:       MatNorm(Bfd,NORM_1,&norm1);
2318:       MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2319:       MatNorm(Bfd,NORM_MAX,&normmax);
2320:       PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2321:       if (flag_display) {MatView(Bfd,vstdout);}
2322:       if (vdraw) {              /* Always use contour for the difference */
2323:         PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2324:         MatView(Bfd,vdraw);
2325:         PetscViewerPopFormat(vdraw);
2326:       }
2327:       if (flag_contour) {PetscViewerPopFormat(vdraw);}

2329:       if (flag_threshold) {
2330:         PetscInt bs,rstart,rend,i;
2331:         MatGetBlockSize(B,&bs);
2332:         MatGetOwnershipRange(B,&rstart,&rend);
2333:         for (i=rstart; i<rend; i++) {
2334:           const PetscScalar *ba,*ca;
2335:           const PetscInt    *bj,*cj;
2336:           PetscInt          bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
2337:           PetscReal         maxentry = 0,maxdiff = 0,maxrdiff = 0;
2338:           MatGetRow(B,i,&bn,&bj,&ba);
2339:           MatGetRow(Bfd,i,&cn,&cj,&ca);
2340:           if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
2341:           for (j=0; j<bn; j++) {
2342:             PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2343:             if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
2344:               maxentrycol = bj[j];
2345:               maxentry    = PetscRealPart(ba[j]);
2346:             }
2347:             if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
2348:               maxdiffcol = bj[j];
2349:               maxdiff    = PetscRealPart(ca[j]);
2350:             }
2351:             if (rdiff > maxrdiff) {
2352:               maxrdiffcol = bj[j];
2353:               maxrdiff    = rdiff;
2354:             }
2355:           }
2356:           if (maxrdiff > 1) {
2357:             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);
2358:             for (j=0; j<bn; j++) {
2359:               PetscReal rdiff;
2360:               rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2361:               if (rdiff > 1) {
2362:                 PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
2363:               }
2364:             }
2365:             PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
2366:           }
2367:           MatRestoreRow(B,i,&bn,&bj,&ba);
2368:           MatRestoreRow(Bfd,i,&cn,&cj,&ca);
2369:         }
2370:       }
2371:       PetscViewerDestroy(&vdraw);
2372:       MatDestroy(&Bfd);
2373:     }
2374:   }
2375:   return(0);
2376: }

2378: /*MC
2379:     SNESJacobianFunction - function used to convey the nonlinear Jacobian of the function to be solved by SNES

2381:      Synopsis:
2382:      #include <petscsnes.h>
2383: $     SNESJacobianFunction(SNES snes,Vec x,Mat Amat,Mat Pmat,void *ctx);

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

2390:    Level: intermediate

2392: .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2393: M*/

2397: /*@C
2398:    SNESSetJacobian - Sets the function to compute Jacobian as well as the
2399:    location to store the matrix.

2401:    Logically Collective on SNES and Mat

2403:    Input Parameters:
2404: +  snes - the SNES context
2405: .  Amat - the matrix that defines the (approximate) Jacobian
2406: .  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2407: .  J - Jacobian evaluation routine (if NULL then SNES retains any previously set value)
2408: -  ctx - [optional] user-defined context for private data for the
2409:          Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)

2411:    Notes:
2412:    If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2413:    each matrix.

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

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

2421:    Level: beginner

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

2425: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J, SNESSetPicard()
2426: @*/
2427: PetscErrorCode  SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2428: {
2430:   DM             dm;

2438:   SNESGetDM(snes,&dm);
2439:   DMSNESSetJacobian(dm,J,ctx);
2440:   if (Amat) {
2441:     PetscObjectReference((PetscObject)Amat);
2442:     MatDestroy(&snes->jacobian);

2444:     snes->jacobian = Amat;
2445:   }
2446:   if (Pmat) {
2447:     PetscObjectReference((PetscObject)Pmat);
2448:     MatDestroy(&snes->jacobian_pre);

2450:     snes->jacobian_pre = Pmat;
2451:   }
2452:   return(0);
2453: }

2457: /*@C
2458:    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2459:    provided context for evaluating the Jacobian.

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

2463:    Input Parameter:
2464: .  snes - the nonlinear solver context

2466:    Output Parameters:
2467: +  Amat - location to stash (approximate) Jacobian matrix (or NULL)
2468: .  Pmat - location to stash matrix used to compute the preconditioner (or NULL)
2469: .  J - location to put Jacobian function (or NULL)
2470: -  ctx - location to stash Jacobian ctx (or NULL)

2472:    Level: advanced

2474: .seealso: SNESSetJacobian(), SNESComputeJacobian()
2475: @*/
2476: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2477: {
2479:   DM             dm;
2480:   DMSNES         sdm;

2484:   if (Amat) *Amat = snes->jacobian;
2485:   if (Pmat) *Pmat = snes->jacobian_pre;
2486:   SNESGetDM(snes,&dm);
2487:   DMGetDMSNES(dm,&sdm);
2488:   if (J) *J = sdm->ops->computejacobian;
2489:   if (ctx) *ctx = sdm->jacobianctx;
2490:   return(0);
2491: }

2495: /*@
2496:    SNESSetUp - Sets up the internal data structures for the later use
2497:    of a nonlinear solver.

2499:    Collective on SNES

2501:    Input Parameters:
2502: .  snes - the SNES context

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

2511:    Level: advanced

2513: .keywords: SNES, nonlinear, setup

2515: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
2516: @*/
2517: PetscErrorCode  SNESSetUp(SNES snes)
2518: {
2520:   DM             dm;
2521:   DMSNES         sdm;
2522:   SNESLineSearch linesearch, pclinesearch;
2523:   void           *lsprectx,*lspostctx;
2524:   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
2525:   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
2526:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2527:   Vec            f,fpc;
2528:   void           *funcctx;
2529:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2530:   void           *jacctx,*appctx;
2531:   Mat            j,jpre;

2535:   if (snes->setupcalled) return(0);

2537:   if (!((PetscObject)snes)->type_name) {
2538:     SNESSetType(snes,SNESNEWTONLS);
2539:   }

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

2543:   SNESGetDM(snes,&dm);
2544:   DMGetDMSNES(dm,&sdm);
2545:   if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
2546:   if (!sdm->ops->computejacobian) {
2547:     DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
2548:   }
2549:   if (!snes->vec_func) {
2550:     DMCreateGlobalVector(dm,&snes->vec_func);
2551:   }

2553:   if (!snes->ksp) {
2554:     SNESGetKSP(snes, &snes->ksp);
2555:   }

2557:   if (!snes->linesearch) {
2558:     SNESGetLineSearch(snes, &snes->linesearch);
2559:   }
2560:   SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);

2562:   if (snes->pc && (snes->pcside == PC_LEFT)) {
2563:     snes->mf          = PETSC_TRUE;
2564:     snes->mf_operator = PETSC_FALSE;
2565:   }

2567:   if (snes->pc) {
2568:     /* copy the DM over */
2569:     SNESGetDM(snes,&dm);
2570:     SNESSetDM(snes->pc,dm);

2572:     SNESGetFunction(snes,&f,&func,&funcctx);
2573:     VecDuplicate(f,&fpc);
2574:     SNESSetFunction(snes->pc,fpc,func,funcctx);
2575:     SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
2576:     SNESSetJacobian(snes->pc,j,jpre,jac,jacctx);
2577:     SNESGetApplicationContext(snes,&appctx);
2578:     SNESSetApplicationContext(snes->pc,appctx);
2579:     VecDestroy(&fpc);

2581:     /* copy the function pointers over */
2582:     PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)snes->pc);

2584:     /* default to 1 iteration */
2585:     SNESSetTolerances(snes->pc,0.0,0.0,0.0,1,snes->pc->max_funcs);
2586:     if (snes->pcside==PC_RIGHT) {
2587:       SNESSetNormSchedule(snes->pc,SNES_NORM_FINAL_ONLY);
2588:     } else {
2589:       SNESSetNormSchedule(snes->pc,SNES_NORM_NONE);
2590:     }
2591:     SNESSetFromOptions(snes->pc);

2593:     /* copy the line search context over */
2594:     SNESGetLineSearch(snes,&linesearch);
2595:     SNESGetLineSearch(snes->pc,&pclinesearch);
2596:     SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
2597:     SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
2598:     SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
2599:     SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
2600:     PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
2601:   }
2602:   if (snes->mf) {
2603:     SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
2604:   }
2605:   if (snes->ops->usercompute && !snes->user) {
2606:     (*snes->ops->usercompute)(snes,(void**)&snes->user);
2607:   }

2609:   snes->jac_iter = 0;
2610:   snes->pre_iter = 0;

2612:   if (snes->ops->setup) {
2613:     (*snes->ops->setup)(snes);
2614:   }

2616:   if (snes->pc && (snes->pcside == PC_LEFT)) {
2617:     if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2618:       SNESGetLineSearch(snes,&linesearch);
2619:       SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
2620:     }
2621:   }

2623:   snes->setupcalled = PETSC_TRUE;
2624:   return(0);
2625: }

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

2632:    Collective on SNES

2634:    Input Parameter:
2635: .  snes - iterative context obtained from SNESCreate()

2637:    Level: intermediate

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

2641: .keywords: SNES, destroy

2643: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
2644: @*/
2645: PetscErrorCode  SNESReset(SNES snes)
2646: {

2651:   if (snes->ops->userdestroy && snes->user) {
2652:     (*snes->ops->userdestroy)((void**)&snes->user);
2653:     snes->user = NULL;
2654:   }
2655:   if (snes->pc) {
2656:     SNESReset(snes->pc);
2657:   }

2659:   if (snes->ops->reset) {
2660:     (*snes->ops->reset)(snes);
2661:   }
2662:   if (snes->ksp) {
2663:     KSPReset(snes->ksp);
2664:   }

2666:   if (snes->linesearch) {
2667:     SNESLineSearchReset(snes->linesearch);
2668:   }

2670:   VecDestroy(&snes->vec_rhs);
2671:   VecDestroy(&snes->vec_sol);
2672:   VecDestroy(&snes->vec_sol_update);
2673:   VecDestroy(&snes->vec_func);
2674:   MatDestroy(&snes->jacobian);
2675:   MatDestroy(&snes->jacobian_pre);
2676:   VecDestroyVecs(snes->nwork,&snes->work);
2677:   VecDestroyVecs(snes->nvwork,&snes->vwork);

2679:   snes->nwork       = snes->nvwork = 0;
2680:   snes->setupcalled = PETSC_FALSE;
2681:   return(0);
2682: }

2686: /*@
2687:    SNESDestroy - Destroys the nonlinear solver context that was created
2688:    with SNESCreate().

2690:    Collective on SNES

2692:    Input Parameter:
2693: .  snes - the SNES context

2695:    Level: beginner

2697: .keywords: SNES, nonlinear, destroy

2699: .seealso: SNESCreate(), SNESSolve()
2700: @*/
2701: PetscErrorCode  SNESDestroy(SNES *snes)
2702: {

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

2710:   SNESReset((*snes));
2711:   SNESDestroy(&(*snes)->pc);

2713:   /* if memory was published with SAWs then destroy it */
2714:   PetscObjectSAWsViewOff((PetscObject)*snes);
2715:   if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}

2717:   DMDestroy(&(*snes)->dm);
2718:   KSPDestroy(&(*snes)->ksp);
2719:   SNESLineSearchDestroy(&(*snes)->linesearch);

2721:   PetscFree((*snes)->kspconvctx);
2722:   if ((*snes)->ops->convergeddestroy) {
2723:     (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
2724:   }
2725:   if ((*snes)->conv_malloc) {
2726:     PetscFree((*snes)->conv_hist);
2727:     PetscFree((*snes)->conv_hist_its);
2728:   }
2729:   SNESMonitorCancel((*snes));
2730:   PetscHeaderDestroy(snes);
2731:   return(0);
2732: }

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

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

2741:    Logically Collective on SNES

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

2748:    Options Database Keys:
2749: .    -snes_lag_preconditioner <lag>

2751:    Notes:
2752:    The default is 1
2753:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2754:    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use

2756:    Level: intermediate

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

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

2762: @*/
2763: PetscErrorCode  SNESSetLagPreconditioner(SNES snes,PetscInt lag)
2764: {
2767:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2768:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2770:   snes->lagpreconditioner = lag;
2771:   return(0);
2772: }

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

2779:    Logically Collective on SNES

2781:    Input Parameters:
2782: +  snes - the SNES context
2783: -  steps - the number of refinements to do, defaults to 0

2785:    Options Database Keys:
2786: .    -snes_grid_sequence <steps>

2788:    Level: intermediate

2790:    Notes:
2791:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

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

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

2797: @*/
2798: PetscErrorCode  SNESSetGridSequence(SNES snes,PetscInt steps)
2799: {
2803:   snes->gridsequence = steps;
2804:   return(0);
2805: }

2809: /*@
2810:    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt

2812:    Not Collective

2814:    Input Parameter:
2815: .  snes - the SNES context

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

2821:    Options Database Keys:
2822: .    -snes_lag_preconditioner <lag>

2824:    Notes:
2825:    The default is 1
2826:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

2828:    Level: intermediate

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

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

2834: @*/
2835: PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
2836: {
2839:   *lag = snes->lagpreconditioner;
2840:   return(0);
2841: }

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

2849:    Logically Collective on SNES

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

2856:    Options Database Keys:
2857: .    -snes_lag_jacobian <lag>

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

2865:    Level: intermediate

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

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

2871: @*/
2872: PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
2873: {
2876:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2877:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2879:   snes->lagjacobian = lag;
2880:   return(0);
2881: }

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

2888:    Not Collective

2890:    Input Parameter:
2891: .  snes - the SNES context

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

2897:    Options Database Keys:
2898: .    -snes_lag_jacobian <lag>

2900:    Notes:
2901:    The default is 1
2902:    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

2904:    Level: intermediate

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

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

2910: @*/
2911: PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
2912: {
2915:   *lag = snes->lagjacobian;
2916:   return(0);
2917: }

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

2924:    Logically collective on SNES

2926:    Input Parameter:
2927: +  snes - the SNES context
2928: -   flg - jacobian lagging persists if true

2930:    Options Database Keys:
2931: .    -snes_lag_jacobian_persists <flg>

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

2937:    Level: developer

2939: .keywords: SNES, nonlinear, lag

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

2943: @*/
2944: PetscErrorCode  SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
2945: {
2949:   snes->lagjac_persist = flg;
2950:   return(0);
2951: }

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

2958:    Logically Collective on SNES

2960:    Input Parameter:
2961: +  snes - the SNES context
2962: -   flg - preconditioner lagging persists if true

2964:    Options Database Keys:
2965: .    -snes_lag_jacobian_persists <flg>

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

2971:    Level: developer

2973: .keywords: SNES, nonlinear, lag

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

2977: @*/
2978: PetscErrorCode  SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
2979: {
2983:   snes->lagpre_persist = flg;
2984:   return(0);
2985: }

2989: /*@
2990:    SNESSetTolerances - Sets various parameters used in convergence tests.

2992:    Logically Collective on SNES

2994:    Input Parameters:
2995: +  snes - the SNES context
2996: .  abstol - absolute convergence tolerance
2997: .  rtol - relative convergence tolerance
2998: .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
2999: .  maxit - maximum number of iterations
3000: -  maxf - maximum number of function evaluations

3002:    Options Database Keys:
3003: +    -snes_atol <abstol> - Sets abstol
3004: .    -snes_rtol <rtol> - Sets rtol
3005: .    -snes_stol <stol> - Sets stol
3006: .    -snes_max_it <maxit> - Sets maxit
3007: -    -snes_max_funcs <maxf> - Sets maxf

3009:    Notes:
3010:    The default maximum number of iterations is 50.
3011:    The default maximum number of function evaluations is 1000.

3013:    Level: intermediate

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

3017: .seealso: SNESSetTrustRegionTolerance()
3018: @*/
3019: PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3020: {

3029:   if (abstol != PETSC_DEFAULT) {
3030:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3031:     snes->abstol = abstol;
3032:   }
3033:   if (rtol != PETSC_DEFAULT) {
3034:     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);
3035:     snes->rtol = rtol;
3036:   }
3037:   if (stol != PETSC_DEFAULT) {
3038:     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3039:     snes->stol = stol;
3040:   }
3041:   if (maxit != PETSC_DEFAULT) {
3042:     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3043:     snes->max_its = maxit;
3044:   }
3045:   if (maxf != PETSC_DEFAULT) {
3046:     if (maxf < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
3047:     snes->max_funcs = maxf;
3048:   }
3049:   snes->tolerancesset = PETSC_TRUE;
3050:   return(0);
3051: }

3055: /*@
3056:    SNESGetTolerances - Gets various parameters used in convergence tests.

3058:    Not Collective

3060:    Input Parameters:
3061: +  snes - the SNES context
3062: .  atol - absolute convergence tolerance
3063: .  rtol - relative convergence tolerance
3064: .  stol -  convergence tolerance in terms of the norm
3065:            of the change in the solution between steps
3066: .  maxit - maximum number of iterations
3067: -  maxf - maximum number of function evaluations

3069:    Notes:
3070:    The user can specify NULL for any parameter that is not needed.

3072:    Level: intermediate

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

3076: .seealso: SNESSetTolerances()
3077: @*/
3078: PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3079: {
3082:   if (atol)  *atol  = snes->abstol;
3083:   if (rtol)  *rtol  = snes->rtol;
3084:   if (stol)  *stol  = snes->stol;
3085:   if (maxit) *maxit = snes->max_its;
3086:   if (maxf)  *maxf  = snes->max_funcs;
3087:   return(0);
3088: }

3092: /*@
3093:    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.

3095:    Logically Collective on SNES

3097:    Input Parameters:
3098: +  snes - the SNES context
3099: -  tol - tolerance

3101:    Options Database Key:
3102: .  -snes_trtol <tol> - Sets tol

3104:    Level: intermediate

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

3108: .seealso: SNESSetTolerances()
3109: @*/
3110: PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3111: {
3115:   snes->deltatol = tol;
3116:   return(0);
3117: }

3119: /*
3120:    Duplicate the lg monitors for SNES from KSP; for some reason with
3121:    dynamic libraries things don't work under Sun4 if we just use
3122:    macros instead of functions
3123: */
3126: PetscErrorCode  SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,void *ctx)
3127: {

3132:   KSPMonitorLGResidualNorm((KSP)snes,it,norm,ctx);
3133:   return(0);
3134: }

3138: PetscErrorCode  SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
3139: {

3143:   KSPMonitorLGResidualNormCreate(host,label,x,y,m,n,draw);
3144:   return(0);
3145: }

3149: PetscErrorCode  SNESMonitorLGDestroy(PetscDrawLG *draw)
3150: {

3154:   KSPMonitorLGResidualNormDestroy(draw);
3155:   return(0);
3156: }

3158: extern PetscErrorCode  SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3161: PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3162: {
3163:   PetscDrawLG      lg;
3164:   PetscErrorCode   ierr;
3165:   PetscReal        x,y,per;
3166:   PetscViewer      v = (PetscViewer)monctx;
3167:   static PetscReal prev; /* should be in the context */
3168:   PetscDraw        draw;

3171:   PetscViewerDrawGetDrawLG(v,0,&lg);
3172:   if (!n) {PetscDrawLGReset(lg);}
3173:   PetscDrawLGGetDraw(lg,&draw);
3174:   PetscDrawSetTitle(draw,"Residual norm");
3175:   x    = (PetscReal)n;
3176:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3177:   else y = -15.0;
3178:   PetscDrawLGAddPoint(lg,&x,&y);
3179:   if (n < 20 || !(n % 5)) {
3180:     PetscDrawLGDraw(lg);
3181:   }

3183:   PetscViewerDrawGetDrawLG(v,1,&lg);
3184:   if (!n) {PetscDrawLGReset(lg);}
3185:   PetscDrawLGGetDraw(lg,&draw);
3186:   PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3187:    SNESMonitorRange_Private(snes,n,&per);
3188:   x    = (PetscReal)n;
3189:   y    = 100.0*per;
3190:   PetscDrawLGAddPoint(lg,&x,&y);
3191:   if (n < 20 || !(n % 5)) {
3192:     PetscDrawLGDraw(lg);
3193:   }

3195:   PetscViewerDrawGetDrawLG(v,2,&lg);
3196:   if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3197:   PetscDrawLGGetDraw(lg,&draw);
3198:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3199:   x    = (PetscReal)n;
3200:   y    = (prev - rnorm)/prev;
3201:   PetscDrawLGAddPoint(lg,&x,&y);
3202:   if (n < 20 || !(n % 5)) {
3203:     PetscDrawLGDraw(lg);
3204:   }

3206:   PetscViewerDrawGetDrawLG(v,3,&lg);
3207:   if (!n) {PetscDrawLGReset(lg);}
3208:   PetscDrawLGGetDraw(lg,&draw);
3209:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3210:   x    = (PetscReal)n;
3211:   y    = (prev - rnorm)/(prev*per);
3212:   if (n > 2) { /*skip initial crazy value */
3213:     PetscDrawLGAddPoint(lg,&x,&y);
3214:   }
3215:   if (n < 20 || !(n % 5)) {
3216:     PetscDrawLGDraw(lg);
3217:   }
3218:   prev = rnorm;
3219:   return(0);
3220: }

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

3227:    Collective on SNES

3229:    Input Parameters:
3230: +  snes - nonlinear solver context obtained from SNESCreate()
3231: .  iter - iteration number
3232: -  rnorm - relative norm of the residual

3234:    Notes:
3235:    This routine is called by the SNES implementations.
3236:    It does not typically need to be called by the user.

3238:    Level: developer

3240: .seealso: SNESMonitorSet()
3241: @*/
3242: PetscErrorCode  SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3243: {
3245:   PetscInt       i,n = snes->numbermonitors;

3248:   for (i=0; i<n; i++) {
3249:     (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3250:   }
3251:   return(0);
3252: }

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

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

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

3263: +    snes - the SNES context
3264: .    its - iteration number
3265: .    norm - 2-norm function value (may be estimated)
3266: -    mctx - [optional] monitoring context

3268:    Level: advanced

3270: .seealso:   SNESMonitorSet(), SNESMonitorGet()
3271: M*/

3275: /*@C
3276:    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3277:    iteration of the nonlinear solver to display the iteration's
3278:    progress.

3280:    Logically Collective on SNES

3282:    Input Parameters:
3283: +  snes - the SNES context
3284: .  f - the monitor function, see SNESMonitorFunction for the calling sequence
3285: .  mctx - [optional] user-defined context for private data for the
3286:           monitor routine (use NULL if no context is desired)
3287: -  monitordestroy - [optional] routine that frees monitor context
3288:           (may be NULL)

3290:    Options Database Keys:
3291: +    -snes_monitor        - sets SNESMonitorDefault()
3292: .    -snes_monitor_lg_residualnorm    - sets line graph monitor,
3293:                             uses SNESMonitorLGCreate()
3294: -    -snes_monitor_cancel - cancels all monitors that have
3295:                             been hardwired into a code by
3296:                             calls to SNESMonitorSet(), but
3297:                             does not cancel those set via
3298:                             the options database.

3300:    Notes:
3301:    Several different monitoring routines may be set by calling
3302:    SNESMonitorSet() multiple times; all will be called in the
3303:    order in which they were set.

3305:    Fortran notes: Only a single monitor function can be set for each SNES object

3307:    Level: intermediate

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

3311: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3312: @*/
3313: PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3314: {
3315:   PetscInt       i;

3320:   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3321:   for (i=0; i<snes->numbermonitors;i++) {
3322:     if (f == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) {
3323:       if (monitordestroy) {
3324:         (*monitordestroy)(&mctx);
3325:       }
3326:       return(0);
3327:     }
3328:   }
3329:   snes->monitor[snes->numbermonitors]          = f;
3330:   snes->monitordestroy[snes->numbermonitors]   = monitordestroy;
3331:   snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3332:   return(0);
3333: }

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

3340:    Logically Collective on SNES

3342:    Input Parameters:
3343: .  snes - the SNES context

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

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

3353:    Level: intermediate

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

3357: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3358: @*/
3359: PetscErrorCode  SNESMonitorCancel(SNES snes)
3360: {
3362:   PetscInt       i;

3366:   for (i=0; i<snes->numbermonitors; i++) {
3367:     if (snes->monitordestroy[i]) {
3368:       (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3369:     }
3370:   }
3371:   snes->numbermonitors = 0;
3372:   return(0);
3373: }

3375: /*MC
3376:     SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver

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

3382: +    snes - the SNES context
3383: .    it - current iteration (0 is the first and is before any Newton step)
3384: .    cctx - [optional] convergence context
3385: .    reason - reason for convergence/divergence
3386: .    xnorm - 2-norm of current iterate
3387: .    gnorm - 2-norm of current step
3388: -    f - 2-norm of function

3390:    Level: intermediate

3392: .seealso:   SNESSetConvergenceTest(), SNESGetConvergenceTest()
3393: M*/

3397: /*@C
3398:    SNESSetConvergenceTest - Sets the function that is to be used
3399:    to test for convergence of the nonlinear iterative solution.

3401:    Logically Collective on SNES

3403:    Input Parameters:
3404: +  snes - the SNES context
3405: .  SNESConvergenceTestFunction - routine to test for convergence
3406: .  cctx - [optional] context for private data for the convergence routine  (may be NULL)
3407: -  destroy - [optional] destructor for the context (may be NULL; NULL_FUNCTION in Fortran)

3409:    Level: advanced

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

3413: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3414: @*/
3415: PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3416: {

3421:   if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
3422:   if (snes->ops->convergeddestroy) {
3423:     (*snes->ops->convergeddestroy)(snes->cnvP);
3424:   }
3425:   snes->ops->converged        = SNESConvergenceTestFunction;
3426:   snes->ops->convergeddestroy = destroy;
3427:   snes->cnvP                  = cctx;
3428:   return(0);
3429: }

3433: /*@
3434:    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.

3436:    Not Collective

3438:    Input Parameter:
3439: .  snes - the SNES context

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

3445:    Level: intermediate

3447:    Notes: Can only be called after the call the SNESSolve() is complete.

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

3451: .seealso: SNESSetConvergenceTest(), SNESConvergedReason
3452: @*/
3453: PetscErrorCode  SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
3454: {
3458:   *reason = snes->reason;
3459:   return(0);
3460: }

3464: /*@
3465:    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.

3467:    Logically Collective on SNES

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

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

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

3485:    Level: intermediate

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

3489: .seealso: SNESGetConvergenceHistory()

3491: @*/
3492: PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
3493: {

3500:   if (!a) {
3501:     if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
3502:     PetscCalloc1(na,&a);
3503:     PetscCalloc1(na,&its);

3505:     snes->conv_malloc = PETSC_TRUE;
3506:   }
3507:   snes->conv_hist       = a;
3508:   snes->conv_hist_its   = its;
3509:   snes->conv_hist_max   = na;
3510:   snes->conv_hist_len   = 0;
3511:   snes->conv_hist_reset = reset;
3512:   return(0);
3513: }

3515: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3516: #include <engine.h>   /* MATLAB include file */
3517: #include <mex.h>      /* MATLAB include file */

3521: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
3522: {
3523:   mxArray   *mat;
3524:   PetscInt  i;
3525:   PetscReal *ar;

3528:   mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
3529:   ar  = (PetscReal*) mxGetData(mat);
3530:   for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
3531:   PetscFunctionReturn(mat);
3532: }
3533: #endif

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

3540:    Not Collective

3542:    Input Parameter:
3543: .  snes - iterative context obtained from SNESCreate()

3545:    Output Parameters:
3546: .  a   - array to hold history
3547: .  its - integer array holds the number of linear iterations (or
3548:          negative if not converged) for each solve.
3549: -  na  - size of a and its

3551:    Notes:
3552:     The calling sequence for this routine in Fortran is
3553: $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)

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

3559:    Level: intermediate

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

3563: .seealso: SNESSetConvergencHistory()

3565: @*/
3566: PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
3567: {
3570:   if (a)   *a   = snes->conv_hist;
3571:   if (its) *its = snes->conv_hist_its;
3572:   if (na)  *na  = snes->conv_hist_len;
3573:   return(0);
3574: }

3578: /*@C
3579:   SNESSetUpdate - Sets the general-purpose update function called
3580:   at the beginning of every iteration of the nonlinear solve. Specifically
3581:   it is called just before the Jacobian is "evaluated".

3583:   Logically Collective on SNES

3585:   Input Parameters:
3586: . snes - The nonlinear solver context
3587: . func - The function

3589:   Calling sequence of func:
3590: . func (SNES snes, PetscInt step);

3592: . step - The current step of the iteration

3594:   Level: advanced

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

3599: .keywords: SNES, update

3601: .seealso SNESSetJacobian(), SNESSolve()
3602: @*/
3603: PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
3604: {
3607:   snes->ops->update = func;
3608:   return(0);
3609: }

3613: /*
3614:    SNESScaleStep_Private - Scales a step so that its length is less than the
3615:    positive parameter delta.

3617:     Input Parameters:
3618: +   snes - the SNES context
3619: .   y - approximate solution of linear system
3620: .   fnorm - 2-norm of current function
3621: -   delta - trust region size

3623:     Output Parameters:
3624: +   gpnorm - predicted function norm at the new point, assuming local
3625:     linearization.  The value is zero if the step lies within the trust
3626:     region, and exceeds zero otherwise.
3627: -   ynorm - 2-norm of the step

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

3633: .keywords: SNES, nonlinear, scale, step
3634: */
3635: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
3636: {
3637:   PetscReal      nrm;
3638:   PetscScalar    cnorm;


3646:   VecNorm(y,NORM_2,&nrm);
3647:   if (nrm > *delta) {
3648:     nrm     = *delta/nrm;
3649:     *gpnorm = (1.0 - nrm)*(*fnorm);
3650:     cnorm   = nrm;
3651:     VecScale(y,cnorm);
3652:     *ynorm  = *delta;
3653:   } else {
3654:     *gpnorm = 0.0;
3655:     *ynorm  = nrm;
3656:   }
3657:   return(0);
3658: }

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

3666:    Collective on SNES

3668:    Input Parameters:
3669: +  snes - the SNES context
3670: .  b - the constant part of the equation F(x) = b, or NULL to use zero.
3671: -  x - the solution vector.

3673:    Notes:
3674:    The user should initialize the vector,x, with the initial guess
3675:    for the nonlinear solve prior to calling SNESSolve.  In particular,
3676:    to employ an initial guess of zero, the user should explicitly set
3677:    this vector to zero by calling VecSet().

3679:    Level: beginner

3681: .keywords: SNES, nonlinear, solve

3683: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
3684: @*/
3685: PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
3686: {
3687:   PetscErrorCode    ierr;
3688:   PetscBool         flg;
3689:   PetscInt          grid;
3690:   Vec               xcreated = NULL;
3691:   DM                dm;


3700:   if (!x) {
3701:     SNESGetDM(snes,&dm);
3702:     DMCreateGlobalVector(dm,&xcreated);
3703:     x    = xcreated;
3704:   }
3705:   SNESViewFromOptions(snes,NULL,"-snes_view_pre");

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

3710:     /* set solution vector */
3711:     if (!grid) {PetscObjectReference((PetscObject)x);}
3712:     VecDestroy(&snes->vec_sol);
3713:     snes->vec_sol = x;
3714:     SNESGetDM(snes,&dm);

3716:     /* set affine vector if provided */
3717:     if (b) { PetscObjectReference((PetscObject)b); }
3718:     VecDestroy(&snes->vec_rhs);
3719:     snes->vec_rhs = b;

3721:     if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
3722:     if (snes->vec_rhs  == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
3723:     if (!snes->vec_sol_update /* && snes->vec_sol */) {
3724:       VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
3725:       PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
3726:     }
3727:     DMShellSetGlobalVector(dm,snes->vec_sol);
3728:     SNESSetUp(snes);

3730:     if (!grid) {
3731:       if (snes->ops->computeinitialguess) {
3732:         (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
3733:       }
3734:     }

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

3739:     PetscLogEventBegin(SNES_Solve,snes,0,0,0);
3740:     (*snes->ops->solve)(snes);
3741:     PetscLogEventEnd(SNES_Solve,snes,0,0,0);
3742:     if (snes->domainerror) {
3743:       snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
3744:       snes->domainerror = PETSC_FALSE;
3745:     }
3746:     if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");

3748:     if (snes->lagjac_persist) snes->jac_iter += snes->iter;
3749:     if (snes->lagpre_persist) snes->pre_iter += snes->iter;

3751:     flg  = PETSC_FALSE;
3752:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,NULL);
3753:     if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
3754:     if (snes->printreason) {
3755:       PetscViewerASCIIAddTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),((PetscObject)snes)->tablevel);
3756:       if (snes->reason > 0) {
3757:         if (((PetscObject) snes)->prefix) {
3758:           PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3759:         } else {
3760:           PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3761:         }
3762:       } else {
3763:         if (((PetscObject) snes)->prefix) {
3764:           PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3765:         } else {
3766:           PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3767:         }
3768:       }
3769:       PetscViewerASCIISubtractTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),((PetscObject)snes)->tablevel);
3770:     }

3772:     if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
3773:     if (grid <  snes->gridsequence) {
3774:       DM  fine;
3775:       Vec xnew;
3776:       Mat interp;

3778:       DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
3779:       if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
3780:       DMCreateInterpolation(snes->dm,fine,&interp,NULL);
3781:       DMCreateGlobalVector(fine,&xnew);
3782:       MatInterpolate(interp,x,xnew);
3783:       DMInterpolate(snes->dm,interp,fine);
3784:       MatDestroy(&interp);
3785:       x    = xnew;

3787:       SNESReset(snes);
3788:       SNESSetDM(snes,fine);
3789:       DMDestroy(&fine);
3790:       PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
3791:     }
3792:   }
3793:   SNESViewFromOptions(snes,NULL,"-snes_view");
3794:   VecViewFromOptions(snes->vec_sol,((PetscObject)snes)->prefix,"-snes_view_solution");

3796:   VecDestroy(&xcreated);
3797:   PetscObjectSAWsBlock((PetscObject)snes);
3798:   return(0);
3799: }

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

3805: /*@C
3806:    SNESSetType - Sets the method for the nonlinear solver.

3808:    Collective on SNES

3810:    Input Parameters:
3811: +  snes - the SNES context
3812: -  type - a known method

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

3818:    Notes:
3819:    See "petsc/include/petscsnes.h" for available methods (for instance)
3820: +    SNESNEWTONLS - Newton's method with line search
3821:      (systems of nonlinear equations)
3822: .    SNESNEWTONTR - Newton's method with trust region
3823:      (systems of nonlinear equations)

3825:   Normally, it is best to use the SNESSetFromOptions() command and then
3826:   set the SNES solver type from the options database rather than by using
3827:   this routine.  Using the options database provides the user with
3828:   maximum flexibility in evaluating the many nonlinear solvers.
3829:   The SNESSetType() routine is provided for those situations where it
3830:   is necessary to set the nonlinear solver independently of the command
3831:   line or options database.  This might be the case, for example, when
3832:   the choice of solver changes during the execution of the program,
3833:   and the user's application is taking responsibility for choosing the
3834:   appropriate method.

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

3839:   Level: intermediate

3841: .keywords: SNES, set, type

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

3845: @*/
3846: PetscErrorCode  SNESSetType(SNES snes,SNESType type)
3847: {
3848:   PetscErrorCode ierr,(*r)(SNES);
3849:   PetscBool      match;


3855:   PetscObjectTypeCompare((PetscObject)snes,type,&match);
3856:   if (match) return(0);

3858:    PetscFunctionListFind(SNESList,type,&r);
3859:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
3860:   /* Destroy the previous private SNES context */
3861:   if (snes->ops->destroy) {
3862:     (*(snes)->ops->destroy)(snes);
3863:     snes->ops->destroy = NULL;
3864:   }
3865:   /* Reinitialize function pointers in SNESOps structure */
3866:   snes->ops->setup          = 0;
3867:   snes->ops->solve          = 0;
3868:   snes->ops->view           = 0;
3869:   snes->ops->setfromoptions = 0;
3870:   snes->ops->destroy        = 0;
3871:   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
3872:   snes->setupcalled = PETSC_FALSE;

3874:   PetscObjectChangeTypeName((PetscObject)snes,type);
3875:   (*r)(snes);
3876:   return(0);
3877: }

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

3884:    Not Collective

3886:    Input Parameter:
3887: .  snes - nonlinear solver context

3889:    Output Parameter:
3890: .  type - SNES method (a character string)

3892:    Level: intermediate

3894: .keywords: SNES, nonlinear, get, type, name
3895: @*/
3896: PetscErrorCode  SNESGetType(SNES snes,SNESType *type)
3897: {
3901:   *type = ((PetscObject)snes)->type_name;
3902:   return(0);
3903: }

3907: /*@
3908:    SNESGetSolution - Returns the vector where the approximate solution is
3909:    stored. This is the fine grid solution when using SNESSetGridSequence().

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

3913:    Input Parameter:
3914: .  snes - the SNES context

3916:    Output Parameter:
3917: .  x - the solution

3919:    Level: intermediate

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

3923: .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
3924: @*/
3925: PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
3926: {
3930:   *x = snes->vec_sol;
3931:   return(0);
3932: }

3936: /*@
3937:    SNESGetSolutionUpdate - Returns the vector where the solution update is
3938:    stored.

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

3942:    Input Parameter:
3943: .  snes - the SNES context

3945:    Output Parameter:
3946: .  x - the solution update

3948:    Level: advanced

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

3952: .seealso: SNESGetSolution(), SNESGetFunction()
3953: @*/
3954: PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
3955: {
3959:   *x = snes->vec_sol_update;
3960:   return(0);
3961: }

3965: /*@C
3966:    SNESGetFunction - Returns the vector where the function is stored.

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

3970:    Input Parameter:
3971: .  snes - the SNES context

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

3978:    Level: advanced

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

3982: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
3983: @*/
3984: PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
3985: {
3987:   DM             dm;

3991:   if (r) {
3992:     if (!snes->vec_func) {
3993:       if (snes->vec_rhs) {
3994:         VecDuplicate(snes->vec_rhs,&snes->vec_func);
3995:       } else if (snes->vec_sol) {
3996:         VecDuplicate(snes->vec_sol,&snes->vec_func);
3997:       } else if (snes->dm) {
3998:         DMCreateGlobalVector(snes->dm,&snes->vec_func);
3999:       }
4000:     }
4001:     *r = snes->vec_func;
4002:   }
4003:   SNESGetDM(snes,&dm);
4004:   DMSNESGetFunction(dm,f,ctx);
4005:   return(0);
4006: }

4008: /*@C
4009:    SNESGetNGS - Returns the NGS function and context.

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

4014:    Output Parameter:
4015: +  f - the function (or NULL) see SNESNGSFunction for details
4016: -  ctx    - the function context (or NULL)

4018:    Level: advanced

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

4022: .seealso: SNESSetNGS(), SNESGetFunction()
4023: @*/

4027: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4028: {
4030:   DM             dm;

4034:   SNESGetDM(snes,&dm);
4035:   DMSNESGetNGS(dm,f,ctx);
4036:   return(0);
4037: }

4041: /*@C
4042:    SNESSetOptionsPrefix - Sets the prefix used for searching for all
4043:    SNES options in the database.

4045:    Logically Collective on SNES

4047:    Input Parameter:
4048: +  snes - the SNES context
4049: -  prefix - the prefix to prepend to all option names

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

4055:    Level: advanced

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

4059: .seealso: SNESSetFromOptions()
4060: @*/
4061: PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
4062: {

4067:   PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4068:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4069:   if (snes->linesearch) {
4070:     SNESGetLineSearch(snes,&snes->linesearch);
4071:     PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4072:   }
4073:   KSPSetOptionsPrefix(snes->ksp,prefix);
4074:   return(0);
4075: }

4079: /*@C
4080:    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4081:    SNES options in the database.

4083:    Logically Collective on SNES

4085:    Input Parameters:
4086: +  snes - the SNES context
4087: -  prefix - the prefix to prepend to all option names

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

4093:    Level: advanced

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

4097: .seealso: SNESGetOptionsPrefix()
4098: @*/
4099: PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4100: {

4105:   PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4106:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4107:   if (snes->linesearch) {
4108:     SNESGetLineSearch(snes,&snes->linesearch);
4109:     PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4110:   }
4111:   KSPAppendOptionsPrefix(snes->ksp,prefix);
4112:   return(0);
4113: }

4117: /*@C
4118:    SNESGetOptionsPrefix - Sets the prefix used for searching for all
4119:    SNES options in the database.

4121:    Not Collective

4123:    Input Parameter:
4124: .  snes - the SNES context

4126:    Output Parameter:
4127: .  prefix - pointer to the prefix string used

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

4132:    Level: advanced

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

4136: .seealso: SNESAppendOptionsPrefix()
4137: @*/
4138: PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4139: {

4144:   PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4145:   return(0);
4146: }


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

4154:    Not collective

4156:    Input Parameters:
4157: +  name_solver - name of a new user-defined solver
4158: -  routine_create - routine to create method context

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

4163:    Sample usage:
4164: .vb
4165:    SNESRegister("my_solver",MySolverCreate);
4166: .ve

4168:    Then, your solver can be chosen with the procedural interface via
4169: $     SNESSetType(snes,"my_solver")
4170:    or at runtime via the option
4171: $     -snes_type my_solver

4173:    Level: advanced

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

4177: .keywords: SNES, nonlinear, register

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

4181:   Level: advanced
4182: @*/
4183: PetscErrorCode  SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4184: {

4188:   PetscFunctionListAdd(&SNESList,sname,function);
4189:   return(0);
4190: }

4194: PetscErrorCode  SNESTestLocalMin(SNES snes)
4195: {
4197:   PetscInt       N,i,j;
4198:   Vec            u,uh,fh;
4199:   PetscScalar    value;
4200:   PetscReal      norm;

4203:   SNESGetSolution(snes,&u);
4204:   VecDuplicate(u,&uh);
4205:   VecDuplicate(u,&fh);

4207:   /* currently only works for sequential */
4208:   PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4209:   VecGetSize(u,&N);
4210:   for (i=0; i<N; i++) {
4211:     VecCopy(u,uh);
4212:     PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4213:     for (j=-10; j<11; j++) {
4214:       value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4215:       VecSetValue(uh,i,value,ADD_VALUES);
4216:       SNESComputeFunction(snes,uh,fh);
4217:       VecNorm(fh,NORM_2,&norm);
4218:       PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);
4219:       value = -value;
4220:       VecSetValue(uh,i,value,ADD_VALUES);
4221:     }
4222:   }
4223:   VecDestroy(&uh);
4224:   VecDestroy(&fh);
4225:   return(0);
4226: }

4230: /*@
4231:    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4232:    computing relative tolerance for linear solvers within an inexact
4233:    Newton method.

4235:    Logically Collective on SNES

4237:    Input Parameters:
4238: +  snes - SNES context
4239: -  flag - PETSC_TRUE or PETSC_FALSE

4241:     Options Database:
4242: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4243: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
4244: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4245: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4246: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
4247: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
4248: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4249: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

4251:    Notes:
4252:    Currently, the default is to use a constant relative tolerance for
4253:    the inner linear solvers.  Alternatively, one can use the
4254:    Eisenstat-Walker method, where the relative convergence tolerance
4255:    is reset at each Newton iteration according progress of the nonlinear
4256:    solver.

4258:    Level: advanced

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

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

4266: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4267: @*/
4268: PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool flag)
4269: {
4273:   snes->ksp_ewconv = flag;
4274:   return(0);
4275: }

4279: /*@
4280:    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
4281:    for computing relative tolerance for linear solvers within an
4282:    inexact Newton method.

4284:    Not Collective

4286:    Input Parameter:
4287: .  snes - SNES context

4289:    Output Parameter:
4290: .  flag - PETSC_TRUE or PETSC_FALSE

4292:    Notes:
4293:    Currently, the default is to use a constant relative tolerance for
4294:    the inner linear solvers.  Alternatively, one can use the
4295:    Eisenstat-Walker method, where the relative convergence tolerance
4296:    is reset at each Newton iteration according progress of the nonlinear
4297:    solver.

4299:    Level: advanced

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

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

4307: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4308: @*/
4309: PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
4310: {
4314:   *flag = snes->ksp_ewconv;
4315:   return(0);
4316: }

4320: /*@
4321:    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
4322:    convergence criteria for the linear solvers within an inexact
4323:    Newton method.

4325:    Logically Collective on SNES

4327:    Input Parameters:
4328: +    snes - SNES context
4329: .    version - version 1, 2 (default is 2) or 3
4330: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4331: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4332: .    gamma - multiplicative factor for version 2 rtol computation
4333:              (0 <= gamma2 <= 1)
4334: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
4335: .    alpha2 - power for safeguard
4336: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

4338:    Note:
4339:    Version 3 was contributed by Luis Chacon, June 2006.

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

4343:    Level: advanced

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

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

4352: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
4353: @*/
4354: PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
4355: {
4356:   SNESKSPEW *kctx;

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

4370:   if (version != PETSC_DEFAULT)   kctx->version   = version;
4371:   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
4372:   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
4373:   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
4374:   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
4375:   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
4376:   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;

4378:   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);
4379:   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);
4380:   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);
4381:   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);
4382:   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);
4383:   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);
4384:   return(0);
4385: }

4389: /*@
4390:    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
4391:    convergence criteria for the linear solvers within an inexact
4392:    Newton method.

4394:    Not Collective

4396:    Input Parameters:
4397:      snes - SNES context

4399:    Output Parameters:
4400: +    version - version 1, 2 (default is 2) or 3
4401: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4402: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4403: .    gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
4404: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
4405: .    alpha2 - power for safeguard
4406: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

4408:    Level: advanced

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

4412: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
4413: @*/
4414: PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
4415: {
4416:   SNESKSPEW *kctx;

4420:   kctx = (SNESKSPEW*)snes->kspconvctx;
4421:   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4422:   if (version)   *version   = kctx->version;
4423:   if (rtol_0)    *rtol_0    = kctx->rtol_0;
4424:   if (rtol_max)  *rtol_max  = kctx->rtol_max;
4425:   if (gamma)     *gamma     = kctx->gamma;
4426:   if (alpha)     *alpha     = kctx->alpha;
4427:   if (alpha2)    *alpha2    = kctx->alpha2;
4428:   if (threshold) *threshold = kctx->threshold;
4429:   return(0);
4430: }

4434:  PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4435: {
4437:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
4438:   PetscReal      rtol  = PETSC_DEFAULT,stol;

4441:   if (!snes->ksp_ewconv) return(0);
4442:   if (!snes->iter) {
4443:     rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
4444:     VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
4445:   }
4446:   else {
4447:     if (kctx->version == 1) {
4448:       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
4449:       if (rtol < 0.0) rtol = -rtol;
4450:       stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
4451:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4452:     } else if (kctx->version == 2) {
4453:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4454:       stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
4455:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4456:     } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
4457:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4458:       /* safeguard: avoid sharp decrease of rtol */
4459:       stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
4460:       stol = PetscMax(rtol,stol);
4461:       rtol = PetscMin(kctx->rtol_0,stol);
4462:       /* safeguard: avoid oversolving */
4463:       stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
4464:       stol = PetscMax(rtol,stol);
4465:       rtol = PetscMin(kctx->rtol_0,stol);
4466:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
4467:   }
4468:   /* safeguard: avoid rtol greater than one */
4469:   rtol = PetscMin(rtol,kctx->rtol_max);
4470:   KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
4471:   PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
4472:   return(0);
4473: }

4477: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4478: {
4480:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
4481:   PCSide         pcside;
4482:   Vec            lres;

4485:   if (!snes->ksp_ewconv) return(0);
4486:   KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
4487:   kctx->norm_last = snes->norm;
4488:   if (kctx->version == 1) {
4489:     PC        pc;
4490:     PetscBool isNone;

4492:     KSPGetPC(ksp, &pc);
4493:     PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
4494:     KSPGetPCSide(ksp,&pcside);
4495:      if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
4496:       /* KSP residual is true linear residual */
4497:       KSPGetResidualNorm(ksp,&kctx->lresid_last);
4498:     } else {
4499:       /* KSP residual is preconditioned residual */
4500:       /* compute true linear residual norm */
4501:       VecDuplicate(b,&lres);
4502:       MatMult(snes->jacobian,x,lres);
4503:       VecAYPX(lres,-1.0,b);
4504:       VecNorm(lres,NORM_2,&kctx->lresid_last);
4505:       VecDestroy(&lres);
4506:     }
4507:   }
4508:   return(0);
4509: }

4513: /*@
4514:    SNESGetKSP - Returns the KSP context for a SNES solver.

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

4518:    Input Parameter:
4519: .  snes - the SNES context

4521:    Output Parameter:
4522: .  ksp - the KSP context

4524:    Notes:
4525:    The user can then directly manipulate the KSP context to set various
4526:    options, etc.  Likewise, the user can then extract and manipulate the
4527:    PC contexts as well.

4529:    Level: beginner

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

4533: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
4534: @*/
4535: PetscErrorCode  SNESGetKSP(SNES snes,KSP *ksp)
4536: {


4543:   if (!snes->ksp) {
4544:     PetscBool monitor = PETSC_FALSE;

4546:     KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
4547:     PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
4548:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);

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

4553:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
4554:     if (monitor) {
4555:       KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
4556:     }
4557:     monitor = PETSC_FALSE;
4558:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
4559:     if (monitor) {
4560:       PetscObject *objs;
4561:       KSPMonitorSNESLGResidualNormCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
4562:       objs[0] = (PetscObject) snes;
4563:       KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
4564:     }
4565:   }
4566:   *ksp = snes->ksp;
4567:   return(0);
4568: }


4571: #include <petsc-private/dmimpl.h>
4574: /*@
4575:    SNESSetDM - Sets the DM that may be used by some preconditioners

4577:    Logically Collective on SNES

4579:    Input Parameters:
4580: +  snes - the preconditioner context
4581: -  dm - the dm

4583:    Level: intermediate

4585: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
4586: @*/
4587: PetscErrorCode  SNESSetDM(SNES snes,DM dm)
4588: {
4590:   KSP            ksp;
4591:   DMSNES         sdm;

4595:   if (dm) {PetscObjectReference((PetscObject)dm);}
4596:   if (snes->dm) {               /* Move the DMSNES context over to the new DM unless the new DM already has one */
4597:     if (snes->dm->dmsnes && snes->dmAuto && !dm->dmsnes) {
4598:       DMCopyDMSNES(snes->dm,dm);
4599:       DMGetDMSNES(snes->dm,&sdm);
4600:       if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
4601:     }
4602:     DMDestroy(&snes->dm);
4603:   }
4604:   snes->dm     = dm;
4605:   snes->dmAuto = PETSC_FALSE;

4607:   SNESGetKSP(snes,&ksp);
4608:   KSPSetDM(ksp,dm);
4609:   KSPSetDMActive(ksp,PETSC_FALSE);
4610:   if (snes->pc) {
4611:     SNESSetDM(snes->pc, snes->dm);
4612:     SNESSetNPCSide(snes,snes->pcside);
4613:   }
4614:   return(0);
4615: }

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

4622:    Not Collective but DM obtained is parallel on SNES

4624:    Input Parameter:
4625: . snes - the preconditioner context

4627:    Output Parameter:
4628: .  dm - the dm

4630:    Level: intermediate

4632: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
4633: @*/
4634: PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
4635: {

4640:   if (!snes->dm) {
4641:     DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
4642:     snes->dmAuto = PETSC_TRUE;
4643:   }
4644:   *dm = snes->dm;
4645:   return(0);
4646: }

4650: /*@
4651:   SNESSetNPC - Sets the nonlinear preconditioner to be used.

4653:   Collective on SNES

4655:   Input Parameters:
4656: + snes - iterative context obtained from SNESCreate()
4657: - pc   - the preconditioner object

4659:   Notes:
4660:   Use SNESGetNPC() to retrieve the preconditioner context (for example,
4661:   to configure it using the API).

4663:   Level: developer

4665: .keywords: SNES, set, precondition
4666: .seealso: SNESGetNPC()
4667: @*/
4668: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
4669: {

4676:   PetscObjectReference((PetscObject) pc);
4677:   SNESDestroy(&snes->pc);
4678:   snes->pc = pc;
4679:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->pc);
4680:   return(0);
4681: }

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

4688:   Not Collective

4690:   Input Parameter:
4691: . snes - iterative context obtained from SNESCreate()

4693:   Output Parameter:
4694: . pc - preconditioner context

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

4698:   Level: developer

4700: .keywords: SNES, get, preconditioner
4701: .seealso: SNESSetNPC()
4702: @*/
4703: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
4704: {
4706:   const char     *optionsprefix;

4711:   if (!snes->pc) {
4712:     SNESCreate(PetscObjectComm((PetscObject)snes),&snes->pc);
4713:     PetscObjectIncrementTabLevel((PetscObject)snes->pc,(PetscObject)snes,1);
4714:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->pc);
4715:     SNESGetOptionsPrefix(snes,&optionsprefix);
4716:     SNESSetOptionsPrefix(snes->pc,optionsprefix);
4717:     SNESAppendOptionsPrefix(snes->pc,"npc_");
4718:     SNESSetCountersReset(snes->pc,PETSC_FALSE);
4719:   }
4720:   *pc = snes->pc;
4721:   return(0);
4722: }

4726: /*@
4727:     SNESSetNPCSide - Sets the preconditioning side.

4729:     Logically Collective on SNES

4731:     Input Parameter:
4732: .   snes - iterative context obtained from SNESCreate()

4734:     Output Parameter:
4735: .   side - the preconditioning side, where side is one of
4736: .vb
4737:       PC_LEFT - left preconditioning (default)
4738:       PC_RIGHT - right preconditioning
4739: .ve

4741:     Options Database Keys:
4742: .   -snes_pc_side <right,left>

4744:     Level: intermediate

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

4748: .seealso: SNESGetNPCSide(), KSPSetPCSide()
4749: @*/
4750: PetscErrorCode  SNESSetNPCSide(SNES snes,PCSide side)
4751: {
4755:   snes->pcside = side;
4756:   return(0);
4757: }

4761: /*@
4762:     SNESGetNPCSide - Gets the preconditioning side.

4764:     Not Collective

4766:     Input Parameter:
4767: .   snes - iterative context obtained from SNESCreate()

4769:     Output Parameter:
4770: .   side - the preconditioning side, where side is one of
4771: .vb
4772:       PC_LEFT - left preconditioning (default)
4773:       PC_RIGHT - right preconditioning
4774: .ve

4776:     Level: intermediate

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

4780: .seealso: SNESSetNPCSide(), KSPGetPCSide()
4781: @*/
4782: PetscErrorCode  SNESGetNPCSide(SNES snes,PCSide *side)
4783: {
4787:   *side = snes->pcside;
4788:   return(0);
4789: }

4793: /*@
4794:   SNESSetLineSearch - Sets the linesearch on the SNES instance.

4796:   Collective on SNES

4798:   Input Parameters:
4799: + snes - iterative context obtained from SNESCreate()
4800: - linesearch   - the linesearch object

4802:   Notes:
4803:   Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
4804:   to configure it using the API).

4806:   Level: developer

4808: .keywords: SNES, set, linesearch
4809: .seealso: SNESGetLineSearch()
4810: @*/
4811: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
4812: {

4819:   PetscObjectReference((PetscObject) linesearch);
4820:   SNESLineSearchDestroy(&snes->linesearch);

4822:   snes->linesearch = linesearch;

4824:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
4825:   return(0);
4826: }

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

4834:   Not Collective

4836:   Input Parameter:
4837: . snes - iterative context obtained from SNESCreate()

4839:   Output Parameter:
4840: . linesearch - linesearch context

4842:   Level: beginner

4844: .keywords: SNES, get, linesearch
4845: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
4846: @*/
4847: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
4848: {
4850:   const char     *optionsprefix;

4855:   if (!snes->linesearch) {
4856:     SNESGetOptionsPrefix(snes, &optionsprefix);
4857:     SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
4858:     SNESLineSearchSetSNES(snes->linesearch, snes);
4859:     SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
4860:     PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
4861:     PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
4862:   }
4863:   *linesearch = snes->linesearch;
4864:   return(0);
4865: }

4867: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4868: #include <mex.h>

4870: typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;

4874: /*
4875:    SNESComputeFunction_Matlab - Calls the function that has been set with SNESSetFunctionMatlab().

4877:    Collective on SNES

4879:    Input Parameters:
4880: +  snes - the SNES context
4881: -  x - input vector

4883:    Output Parameter:
4884: .  y - function vector, as set by SNESSetFunction()

4886:    Notes:
4887:    SNESComputeFunction() is typically used within nonlinear solvers
4888:    implementations, so most users would not generally call this routine
4889:    themselves.

4891:    Level: developer

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

4895: .seealso: SNESSetFunction(), SNESGetFunction()
4896: */
4897: PetscErrorCode  SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
4898: {
4899:   PetscErrorCode    ierr;
4900:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
4901:   int               nlhs  = 1,nrhs = 5;
4902:   mxArray           *plhs[1],*prhs[5];
4903:   long long int     lx = 0,ly = 0,ls = 0;


4912:   /* call Matlab function in ctx with arguments x and y */

4914:   PetscMemcpy(&ls,&snes,sizeof(snes));
4915:   PetscMemcpy(&lx,&x,sizeof(x));
4916:   PetscMemcpy(&ly,&y,sizeof(x));
4917:   prhs[0] = mxCreateDoubleScalar((double)ls);
4918:   prhs[1] = mxCreateDoubleScalar((double)lx);
4919:   prhs[2] = mxCreateDoubleScalar((double)ly);
4920:   prhs[3] = mxCreateString(sctx->funcname);
4921:   prhs[4] = sctx->ctx;
4922:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");
4923:   mxGetScalar(plhs[0]);
4924:   mxDestroyArray(prhs[0]);
4925:   mxDestroyArray(prhs[1]);
4926:   mxDestroyArray(prhs[2]);
4927:   mxDestroyArray(prhs[3]);
4928:   mxDestroyArray(plhs[0]);
4929:   return(0);
4930: }

4934: /*
4935:    SNESSetFunctionMatlab - Sets the function evaluation routine and function
4936:    vector for use by the SNES routines in solving systems of nonlinear
4937:    equations from MATLAB. Here the function is a string containing the name of a MATLAB function

4939:    Logically Collective on SNES

4941:    Input Parameters:
4942: +  snes - the SNES context
4943: .  r - vector to store function value
4944: -  f - function evaluation routine

4946:    Notes:
4947:    The Newton-like methods typically solve linear systems of the form
4948: $      f'(x) x = -f(x),
4949:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

4951:    Level: beginner

4953:    Developer Note:  This bleeds the allocated memory SNESMatlabContext *sctx;

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

4957: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
4958: */
4959: PetscErrorCode  SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx)
4960: {
4961:   PetscErrorCode    ierr;
4962:   SNESMatlabContext *sctx;

4965:   /* currently sctx is memory bleed */
4966:   PetscMalloc(sizeof(SNESMatlabContext),&sctx);
4967:   PetscStrallocpy(f,&sctx->funcname);
4968:   /*
4969:      This should work, but it doesn't
4970:   sctx->ctx = ctx;
4971:   mexMakeArrayPersistent(sctx->ctx);
4972:   */
4973:   sctx->ctx = mxDuplicateArray(ctx);
4974:   SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);
4975:   return(0);
4976: }

4980: /*
4981:    SNESComputeJacobian_Matlab - Calls the function that has been set with SNESSetJacobianMatlab().

4983:    Collective on SNES

4985:    Input Parameters:
4986: +  snes - the SNES context
4987: .  x - input vector
4988: .  A, B - the matrices
4989: -  ctx - user context

4991:    Level: developer

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

4995: .seealso: SNESSetFunction(), SNESGetFunction()
4996: @*/
4997: PetscErrorCode  SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx)
4998: {
4999:   PetscErrorCode    ierr;
5000:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5001:   int               nlhs  = 2,nrhs = 6;
5002:   mxArray           *plhs[2],*prhs[6];
5003:   long long int     lx = 0,lA = 0,ls = 0, lB = 0;


5009:   /* call Matlab function in ctx with arguments x and y */

5011:   PetscMemcpy(&ls,&snes,sizeof(snes));
5012:   PetscMemcpy(&lx,&x,sizeof(x));
5013:   PetscMemcpy(&lA,A,sizeof(x));
5014:   PetscMemcpy(&lB,B,sizeof(x));
5015:   prhs[0] = mxCreateDoubleScalar((double)ls);
5016:   prhs[1] = mxCreateDoubleScalar((double)lx);
5017:   prhs[2] = mxCreateDoubleScalar((double)lA);
5018:   prhs[3] = mxCreateDoubleScalar((double)lB);
5019:   prhs[4] = mxCreateString(sctx->funcname);
5020:   prhs[5] = sctx->ctx;
5021:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");
5022:   mxGetScalar(plhs[0]);
5023:   mxDestroyArray(prhs[0]);
5024:   mxDestroyArray(prhs[1]);
5025:   mxDestroyArray(prhs[2]);
5026:   mxDestroyArray(prhs[3]);
5027:   mxDestroyArray(prhs[4]);
5028:   mxDestroyArray(plhs[0]);
5029:   mxDestroyArray(plhs[1]);
5030:   return(0);
5031: }

5035: /*
5036:    SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
5037:    vector for use by the SNES routines in solving systems of nonlinear
5038:    equations from MATLAB. Here the function is a string containing the name of a MATLAB function

5040:    Logically Collective on SNES

5042:    Input Parameters:
5043: +  snes - the SNES context
5044: .  A,B - Jacobian matrices
5045: .  J - function evaluation routine
5046: -  ctx - user context

5048:    Level: developer

5050:    Developer Note:  This bleeds the allocated memory SNESMatlabContext *sctx;

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

5054: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J
5055: */
5056: PetscErrorCode  SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx)
5057: {
5058:   PetscErrorCode    ierr;
5059:   SNESMatlabContext *sctx;

5062:   /* currently sctx is memory bleed */
5063:   PetscMalloc(sizeof(SNESMatlabContext),&sctx);
5064:   PetscStrallocpy(J,&sctx->funcname);
5065:   /*
5066:      This should work, but it doesn't
5067:   sctx->ctx = ctx;
5068:   mexMakeArrayPersistent(sctx->ctx);
5069:   */
5070:   sctx->ctx = mxDuplicateArray(ctx);
5071:   SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);
5072:   return(0);
5073: }

5077: /*
5078:    SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().

5080:    Collective on SNES

5082: .seealso: SNESSetFunction(), SNESGetFunction()
5083: @*/
5084: PetscErrorCode  SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
5085: {
5086:   PetscErrorCode    ierr;
5087:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5088:   int               nlhs  = 1,nrhs = 6;
5089:   mxArray           *plhs[1],*prhs[6];
5090:   long long int     lx = 0,ls = 0;
5091:   Vec               x  = snes->vec_sol;


5096:   PetscMemcpy(&ls,&snes,sizeof(snes));
5097:   PetscMemcpy(&lx,&x,sizeof(x));
5098:   prhs[0] = mxCreateDoubleScalar((double)ls);
5099:   prhs[1] = mxCreateDoubleScalar((double)it);
5100:   prhs[2] = mxCreateDoubleScalar((double)fnorm);
5101:   prhs[3] = mxCreateDoubleScalar((double)lx);
5102:   prhs[4] = mxCreateString(sctx->funcname);
5103:   prhs[5] = sctx->ctx;
5104:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");
5105:   mxGetScalar(plhs[0]);
5106:   mxDestroyArray(prhs[0]);
5107:   mxDestroyArray(prhs[1]);
5108:   mxDestroyArray(prhs[2]);
5109:   mxDestroyArray(prhs[3]);
5110:   mxDestroyArray(prhs[4]);
5111:   mxDestroyArray(plhs[0]);
5112:   return(0);
5113: }

5117: /*
5118:    SNESMonitorSetMatlab - Sets the monitor function from MATLAB

5120:    Level: developer

5122:    Developer Note:  This bleeds the allocated memory SNESMatlabContext *sctx;

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

5126: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5127: */
5128: PetscErrorCode  SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx)
5129: {
5130:   PetscErrorCode    ierr;
5131:   SNESMatlabContext *sctx;

5134:   /* currently sctx is memory bleed */
5135:   PetscMalloc(sizeof(SNESMatlabContext),&sctx);
5136:   PetscStrallocpy(f,&sctx->funcname);
5137:   /*
5138:      This should work, but it doesn't
5139:   sctx->ctx = ctx;
5140:   mexMakeArrayPersistent(sctx->ctx);
5141:   */
5142:   sctx->ctx = mxDuplicateArray(ctx);
5143:   SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);
5144:   return(0);
5145: }

5147: #endif