Actual source code: snes.c

petsc-master 2015-03-04
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      issaws;
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,&issaws);
246: #endif
247:   if (iascii) {
248:     SNESNormSchedule normschedule;

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

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

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

326:     PetscObjectGetName((PetscObject)snes,&name);
327:     MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
328:     if (!((PetscObject)snes)->amsmem && !rank) {
329:       char       dir[1024];

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

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

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

384:   Not Collective

386:   Input Parameter:
387: . snescheck - function that checks for options

389:   Level: developer

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

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

405: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
406: {
407:   Mat            J;
408:   KSP            ksp;
409:   PC             pc;
410:   PetscBool      match;


416:   if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
417:     Mat A = snes->jacobian, B = snes->jacobian_pre;
418:     MatCreateVecs(A ? A : B, NULL,&snes->vec_func);
419:   }

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

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

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

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

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

497: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
498: {

502:   DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
503:   return(0);
504: }

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

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

535:   SNESComputeJacobian(snes,X,A,B);
536:   if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
537:     SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
538:   }

540:   if (A != Asave || B != Bsave) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for changing matrices at this time");
541:   if (Xnamed) {
542:     DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
543:   }
544:   snes->dm = dmsave;
545:   return(0);
546: }

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

553:    Collective

555:    Input Arguments:
556: .  snes - snes to configure

558:    Level: developer

560: .seealso: SNESSetUp()
561: @*/
562: PetscErrorCode SNESSetUpMatrices(SNES snes)
563: {
565:   DM             dm;
566:   DMSNES         sdm;

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

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

612:    Collective on SNES

614:    Input Parameter:
615: .  snes - the SNES context

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

646:     Options Database for Eisenstat-Walker method:
647: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
648: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
649: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
650: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
651: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
652: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
653: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
654: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

656:    Notes:
657:    To see all options, run your program with the -help option or consult
658:    Users-Manual: ch_snes

660:    Level: beginner

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

664: .seealso: SNESSetOptionsPrefix()
665: @*/
666: PetscErrorCode  SNESSetFromOptions(SNES snes)
667: {
668:   PetscBool      flg,pcset,persist,set;
669:   PetscInt       i,indx,lag,grids;
670:   const char     *deft        = SNESNEWTONLS;
671:   const char     *convtests[] = {"default","skip"};
672:   SNESKSPEW      *kctx        = NULL;
673:   char           type[256], monfilename[PETSC_MAX_PATH_LEN];
674:   PetscViewer    monviewer;
676:   PCSide         pcside;
677:   const char     *optionsprefix;

681:   SNESRegisterAll();
682:   PetscObjectOptionsBegin((PetscObject)snes);
683:   if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
684:   PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
685:   if (flg) {
686:     SNESSetType(snes,type);
687:   } else if (!((PetscObject)snes)->type_name) {
688:     SNESSetType(snes,deft);
689:   }
690:   PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);
691:   PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);

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

700:   PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
701:   if (flg) {
702:     SNESSetLagPreconditioner(snes,lag);
703:   }
704:   PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
705:   if (flg) {
706:     SNESSetLagPreconditionerPersists(snes,persist);
707:   }
708:   PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
709:   if (flg) {
710:     SNESSetLagJacobian(snes,lag);
711:   }
712:   PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
713:   if (flg) {
714:     SNESSetLagJacobianPersists(snes,persist);
715:   }

717:   PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
718:   if (flg) {
719:     SNESSetGridSequence(snes,grids);
720:   }

722:   PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);
723:   if (flg) {
724:     switch (indx) {
725:     case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
726:     case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL);    break;
727:     }
728:   }

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

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

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

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

740:   PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);
741:   PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);
742:   PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);
743:   PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);
744:   PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);
745:   PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);
746:   PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);

748:   flg  = PETSC_FALSE;
749:   PetscOptionsBool("-snes_check_jacobian","Check each Jacobian with a differenced one","SNESUpdateCheckJacobian",flg,&flg,&set);
750:   if (set && flg) {
751:     SNESSetUpdate(snes,SNESUpdateCheckJacobian);
752:   }

754:   flg  = PETSC_FALSE;
755:   PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);
756:   if (set && flg) {SNESMonitorCancel(snes);}

758:   PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
759:   if (flg) {
760:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
761:     SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
762:   }

764:   PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
765:   if (flg) {
766:     SNESMonitorSet(snes,SNESMonitorRange,0,0);
767:   }

769:   PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
770:   if (flg) {
771:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
772:     SNESMonitorSetRatio(snes,monviewer);
773:   }

775:   PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
776:   if (flg) {
777:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
778:     SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
779:   }

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

784:   flg  = PETSC_FALSE;
785:   PetscOptionsBool("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,NULL);
786:   if (flg) {SNESMonitorSet(snes,SNESMonitorSolution,0,0);}
787:   flg  = PETSC_FALSE;
788:   PetscOptionsBool("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,NULL);
789:   if (flg) {SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);}
790:   flg  = PETSC_FALSE;
791:   PetscOptionsBool("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,NULL);
792:   if (flg) {SNESMonitorSet(snes,SNESMonitorResidual,0,0);}
793:   flg  = PETSC_FALSE;
794:   PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);
795:   if (flg) {
796:     PetscObject *objs;

798:     SNESMonitorLGCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&objs);
799:     SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorLGResidualNorm,objs,(PetscErrorCode (*)(void**))SNESMonitorLGDestroy);
800:   }
801:   flg  = PETSC_FALSE;
802:   PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
803:   if (flg) {
804:     PetscViewer ctx;

806:     PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
807:     SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
808:   }

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


815:   PetscOptionsString("-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
816:   if (flg) {
817:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
818:     SNESMonitorSet(snes,SNESMonitorFields,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
819:   }

821:   flg  = PETSC_FALSE;
822:   PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
823:   if (flg) {
824:     void *functx;
825:     SNESGetFunction(snes,NULL,NULL,&functx);
826:     SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
827:     PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
828:   }

830:   flg  = PETSC_FALSE;
831:   PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
832:   if (flg) {
833:     SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
834:   }

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

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

859:   flg  = PETSC_FALSE;
860:   SNESGetNPCSide(snes,&pcside);
861:   PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
862:   if (flg) {SNESSetNPCSide(snes,pcside);}

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

887:   for (i = 0; i < numberofsetfromoptions; i++) {
888:     (*othersetfromoptions[i])(snes);
889:   }

891:   if (snes->ops->setfromoptions) {
892:     (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
893:   }

895:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
896:   PetscObjectProcessOptionsHandlers((PetscObject)snes);
897:   PetscOptionsEnd();

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

904:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
905:   KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
906:   KSPSetFromOptions(snes->ksp);

908:   /* if someone has set the SNES NPC type, create it. */
909:   SNESGetOptionsPrefix(snes, &optionsprefix);
910:   PetscOptionsHasName(optionsprefix, "-npc_snes_type", &pcset);
911:   if (pcset && (!snes->pc)) {
912:     SNESGetNPC(snes, &snes->pc);
913:   }
914:   return(0);
915: }

919: /*@C
920:    SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
921:    the nonlinear solvers.

923:    Logically Collective on SNES

925:    Input Parameters:
926: +  snes - the SNES context
927: .  compute - function to compute the context
928: -  destroy - function to destroy the context

930:    Level: intermediate

932:    Notes:
933:    This function is currently not available from Fortran.

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

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

950: /*@
951:    SNESSetApplicationContext - Sets the optional user-defined context for
952:    the nonlinear solvers.

954:    Logically Collective on SNES

956:    Input Parameters:
957: +  snes - the SNES context
958: -  usrP - optional user context

960:    Level: intermediate

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

964: .seealso: SNESGetApplicationContext()
965: @*/
966: PetscErrorCode  SNESSetApplicationContext(SNES snes,void *usrP)
967: {
969:   KSP            ksp;

973:   SNESGetKSP(snes,&ksp);
974:   KSPSetApplicationContext(ksp,usrP);
975:   snes->user = usrP;
976:   return(0);
977: }

981: /*@
982:    SNESGetApplicationContext - Gets the user-defined context for the
983:    nonlinear solvers.

985:    Not Collective

987:    Input Parameter:
988: .  snes - SNES context

990:    Output Parameter:
991: .  usrP - user context

993:    Level: intermediate

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

997: .seealso: SNESSetApplicationContext()
998: @*/
999: PetscErrorCode  SNESGetApplicationContext(SNES snes,void *usrP)
1000: {
1003:   *(void**)usrP = snes->user;
1004:   return(0);
1005: }

1009: /*@
1010:    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1011:    at this time.

1013:    Not Collective

1015:    Input Parameter:
1016: .  snes - SNES context

1018:    Output Parameter:
1019: .  iter - iteration number

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

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

1035:    Level: intermediate

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

1039: .seealso:   SNESGetLinearSolveIterations()
1040: @*/
1041: PetscErrorCode  SNESGetIterationNumber(SNES snes,PetscInt *iter)
1042: {
1046:   *iter = snes->iter;
1047:   return(0);
1048: }

1052: /*@
1053:    SNESSetIterationNumber - Sets the current iteration number.

1055:    Not Collective

1057:    Input Parameter:
1058: .  snes - SNES context
1059: .  iter - iteration number

1061:    Level: developer

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

1065: .seealso:   SNESGetLinearSolveIterations()
1066: @*/
1067: PetscErrorCode  SNESSetIterationNumber(SNES snes,PetscInt iter)
1068: {

1073:   PetscObjectSAWsTakeAccess((PetscObject)snes);
1074:   snes->iter = iter;
1075:   PetscObjectSAWsGrantAccess((PetscObject)snes);
1076:   return(0);
1077: }

1081: /*@
1082:    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1083:    attempted by the nonlinear solver.

1085:    Not Collective

1087:    Input Parameter:
1088: .  snes - SNES context

1090:    Output Parameter:
1091: .  nfails - number of unsuccessful steps attempted

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

1096:    Level: intermediate

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

1100: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1101:           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1102: @*/
1103: PetscErrorCode  SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1104: {
1108:   *nfails = snes->numFailures;
1109:   return(0);
1110: }

1114: /*@
1115:    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1116:    attempted by the nonlinear solver before it gives up.

1118:    Not Collective

1120:    Input Parameters:
1121: +  snes     - SNES context
1122: -  maxFails - maximum of unsuccessful steps

1124:    Level: intermediate

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

1128: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1129:           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1130: @*/
1131: PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1132: {
1135:   snes->maxFailures = maxFails;
1136:   return(0);
1137: }

1141: /*@
1142:    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1143:    attempted by the nonlinear solver before it gives up.

1145:    Not Collective

1147:    Input Parameter:
1148: .  snes     - SNES context

1150:    Output Parameter:
1151: .  maxFails - maximum of unsuccessful steps

1153:    Level: intermediate

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

1157: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1158:           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()

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

1172: /*@
1173:    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1174:      done by SNES.

1176:    Not Collective

1178:    Input Parameter:
1179: .  snes     - SNES context

1181:    Output Parameter:
1182: .  nfuncs - number of evaluations

1184:    Level: intermediate

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

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

1190: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1191: @*/
1192: PetscErrorCode  SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1193: {
1197:   *nfuncs = snes->nfuncs;
1198:   return(0);
1199: }

1203: /*@
1204:    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1205:    linear solvers.

1207:    Not Collective

1209:    Input Parameter:
1210: .  snes - SNES context

1212:    Output Parameter:
1213: .  nfails - number of failed solves

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

1218:    Level: intermediate

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

1222: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1223: @*/
1224: PetscErrorCode  SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1225: {
1229:   *nfails = snes->numLinearSolveFailures;
1230:   return(0);
1231: }

1235: /*@
1236:    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1237:    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE

1239:    Logically Collective on SNES

1241:    Input Parameters:
1242: +  snes     - SNES context
1243: -  maxFails - maximum allowed linear solve failures

1245:    Level: intermediate

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

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

1251: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1252: @*/
1253: PetscErrorCode  SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1254: {
1258:   snes->maxLinearSolveFailures = maxFails;
1259:   return(0);
1260: }

1264: /*@
1265:    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1266:      are allowed before SNES terminates

1268:    Not Collective

1270:    Input Parameter:
1271: .  snes     - SNES context

1273:    Output Parameter:
1274: .  maxFails - maximum of unsuccessful solves allowed

1276:    Level: intermediate

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

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

1282: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1283: @*/
1284: PetscErrorCode  SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1285: {
1289:   *maxFails = snes->maxLinearSolveFailures;
1290:   return(0);
1291: }

1295: /*@
1296:    SNESGetLinearSolveIterations - Gets the total number of linear iterations
1297:    used by the nonlinear solver.

1299:    Not Collective

1301:    Input Parameter:
1302: .  snes - SNES context

1304:    Output Parameter:
1305: .  lits - number of linear iterations

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

1310:    Level: intermediate

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

1314: .seealso:  SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1315: @*/
1316: PetscErrorCode  SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1317: {
1321:   *lits = snes->linear_its;
1322:   return(0);
1323: }

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

1331:    Logically Collective on SNES

1333:    Input Parameter:
1334: +  snes - SNES context
1335: -  reset - whether to reset the counters or not

1337:    Notes:
1338:    This is automatically called with FALSE

1340:    Level: developer

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

1344: .seealso:  SNESGetNumberFunctionEvals(), SNESGetNumberLinearSolveIterations(), SNESGetNPC()
1345: @*/
1346: PetscErrorCode  SNESSetCountersReset(SNES snes,PetscBool reset)
1347: {
1351:   snes->counters_reset = reset;
1352:   return(0);
1353: }


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

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

1363:    Input Parameters:
1364: +  snes - the SNES context
1365: -  ksp - the KSP context

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

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

1374:    Level: developer

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

1378: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1379: @*/
1380: PetscErrorCode  SNESSetKSP(SNES snes,KSP ksp)
1381: {

1388:   PetscObjectReference((PetscObject)ksp);
1389:   if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1390:   snes->ksp = ksp;
1391:   return(0);
1392: }

1394: /* -----------------------------------------------------------*/
1397: /*@
1398:    SNESCreate - Creates a nonlinear solver context.

1400:    Collective on MPI_Comm

1402:    Input Parameters:
1403: .  comm - MPI communicator

1405:    Output Parameter:
1406: .  outsnes - the new SNES context

1408:    Options Database Keys:
1409: +   -snes_mf - Activates default matrix-free Jacobian-vector products,
1410:                and no preconditioning matrix
1411: .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
1412:                products, and a user-provided preconditioning matrix
1413:                as set by SNESSetJacobian()
1414: -   -snes_fd - Uses (slow!) finite differences to compute Jacobian

1416:    Level: beginner

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

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

1422: @*/
1423: PetscErrorCode  SNESCreate(MPI_Comm comm,SNES *outsnes)
1424: {
1426:   SNES           snes;
1427:   SNESKSPEW      *kctx;

1431:   *outsnes = NULL;
1432:   SNESInitializePackage();

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

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

1489:   snes->mf          = PETSC_FALSE;
1490:   snes->mf_operator = PETSC_FALSE;
1491:   snes->mf_version  = 1;

1493:   snes->numLinearSolveFailures = 0;
1494:   snes->maxLinearSolveFailures = 1;

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

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

1512:   *outsnes = snes;
1513:   return(0);
1514: }

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

1519:      Synopsis:
1520:      #include "petscsnes.h"
1521:      PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);

1523:      Input Parameters:
1524: +     snes - the SNES context
1525: .     x    - state at which to evaluate residual
1526: -     ctx     - optional user-defined function context, passed in with SNESSetFunction()

1528:      Output Parameter:
1529: .     f  - vector to put residual (function value)

1531:    Level: intermediate

1533: .seealso:   SNESSetFunction(), SNESGetFunction()
1534: M*/

1538: /*@C
1539:    SNESSetFunction - Sets the function evaluation routine and function
1540:    vector for use by the SNES routines in solving systems of nonlinear
1541:    equations.

1543:    Logically Collective on SNES

1545:    Input Parameters:
1546: +  snes - the SNES context
1547: .  r - vector to store function value
1548: .  f - function evaluation routine; see SNESFunction for calling sequence details
1549: -  ctx - [optional] user-defined context for private data for the
1550:          function evaluation routine (may be NULL)

1552:    Notes:
1553:    The Newton-like methods typically solve linear systems of the form
1554: $      f'(x) x = -f(x),
1555:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

1557:    Level: beginner

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

1561: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1562: @*/
1563: PetscErrorCode  SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1564: {
1566:   DM             dm;

1570:   if (r) {
1573:     PetscObjectReference((PetscObject)r);
1574:     VecDestroy(&snes->vec_func);

1576:     snes->vec_func = r;
1577:   }
1578:   SNESGetDM(snes,&dm);
1579:   DMSNESSetFunction(dm,f,ctx);
1580:   return(0);
1581: }


1586: /*@C
1587:    SNESSetInitialFunction - Sets the function vector to be used as the
1588:    function norm at the initialization of the method.  In some
1589:    instances, the user has precomputed the function before calling
1590:    SNESSolve.  This function allows one to avoid a redundant call
1591:    to SNESComputeFunction in that case.

1593:    Logically Collective on SNES

1595:    Input Parameters:
1596: +  snes - the SNES context
1597: -  f - vector to store function value

1599:    Notes:
1600:    This should not be modified during the solution procedure.

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

1604:    Level: developer

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

1608: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1609: @*/
1610: PetscErrorCode  SNESSetInitialFunction(SNES snes, Vec f)
1611: {
1613:   Vec            vec_func;

1619:   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1620:     snes->vec_func_init_set = PETSC_FALSE;
1621:     return(0);
1622:   }
1623:   SNESGetFunction(snes,&vec_func,NULL,NULL);
1624:   VecCopy(f, vec_func);

1626:   snes->vec_func_init_set = PETSC_TRUE;
1627:   return(0);
1628: }

1632: /*@
1633:    SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1634:    of the SNES method.

1636:    Logically Collective on SNES

1638:    Input Parameters:
1639: +  snes - the SNES context
1640: -  normschedule - the frequency of norm computation

1642:    Options Database Key:
1643: .  -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly>

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

1654:    Level: developer

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

1658: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1659: @*/
1660: PetscErrorCode  SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1661: {
1664:   snes->normschedule = normschedule;
1665:   return(0);
1666: }


1671: /*@
1672:    SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1673:    of the SNES method.

1675:    Logically Collective on SNES

1677:    Input Parameters:
1678: +  snes - the SNES context
1679: -  normschedule - the type of the norm used

1681:    Level: advanced

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

1685: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1686: @*/
1687: PetscErrorCode  SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1688: {
1691:   *normschedule = snes->normschedule;
1692:   return(0);
1693: }


1698: /*@C
1699:    SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
1700:    of the SNES method.

1702:    Logically Collective on SNES

1704:    Input Parameters:
1705: +  snes - the SNES context
1706: -  normschedule - the frequency of norm computation

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

1717:    Level: developer

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

1721: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1722: @*/
1723: PetscErrorCode  SNESSetFunctionType(SNES snes, SNESFunctionType type)
1724: {
1727:   snes->functype = type;
1728:   return(0);
1729: }


1734: /*@C
1735:    SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
1736:    of the SNES method.

1738:    Logically Collective on SNES

1740:    Input Parameters:
1741: +  snes - the SNES context
1742: -  normschedule - the type of the norm used

1744:    Level: advanced

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

1748: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1749: @*/
1750: PetscErrorCode  SNESGetFunctionType(SNES snes, SNESFunctionType *type)
1751: {
1754:   *type = snes->functype;
1755:   return(0);
1756: }

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

1761:      Synopsis:
1762:      #include <petscsnes.h>
1763: $    SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);

1765: +  X   - solution vector
1766: .  B   - RHS vector
1767: -  ctx - optional user-defined Gauss-Seidel context

1769:    Level: intermediate

1771: .seealso:   SNESSetNGS(), SNESGetNGS()
1772: M*/

1776: /*@C
1777:    SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
1778:    use with composed nonlinear solvers.

1780:    Input Parameters:
1781: +  snes   - the SNES context
1782: .  f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
1783: -  ctx    - [optional] user-defined context for private data for the
1784:             smoother evaluation routine (may be NULL)

1786:    Notes:
1787:    The NGS routines are used by the composed nonlinear solver to generate
1788:     a problem appropriate update to the solution, particularly FAS.

1790:    Level: intermediate

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

1794: .seealso: SNESGetFunction(), SNESComputeNGS()
1795: @*/
1796: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1797: {
1799:   DM             dm;

1803:   SNESGetDM(snes,&dm);
1804:   DMSNESSetNGS(dm,f,ctx);
1805:   return(0);
1806: }

1810: PETSC_EXTERN PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
1811: {
1813:   DM             dm;
1814:   DMSNES         sdm;

1817:   SNESGetDM(snes,&dm);
1818:   DMGetDMSNES(dm,&sdm);
1819:   /*  A(x)*x - b(x) */
1820:   if (sdm->ops->computepfunction) {
1821:     (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
1822:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");

1824:   if (sdm->ops->computepjacobian) {
1825:     (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
1826:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard matrix.");
1827:   VecScale(f,-1.0);
1828:   MatMultAdd(snes->jacobian,x,f,f);
1829:   return(0);
1830: }

1834: PETSC_EXTERN PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
1835: {
1837:   /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
1838:   return(0);
1839: }

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

1846:    Logically Collective on SNES

1848:    Input Parameters:
1849: +  snes - the SNES context
1850: .  r - vector to store function value
1851: .  b - function evaluation routine
1852: .  Amat - matrix with which A(x) x - b(x) is to be computed
1853: .  Pmat - matrix from which preconditioner is computed (usually the same as Amat)
1854: .  J  - function to compute matrix value, see SNESJacobianFunction for details on its calling sequence
1855: -  ctx - [optional] user-defined context for private data for the
1856:          function evaluation routine (may be NULL)

1858:    Notes:
1859:     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
1860:     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.

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

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

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

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

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

1876:    Level: intermediate

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

1880: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
1881: @*/
1882: 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)
1883: {
1885:   DM             dm;

1889:   SNESGetDM(snes, &dm);
1890:   DMSNESSetPicard(dm,b,J,ctx);
1891:   SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
1892:   SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
1893:   return(0);
1894: }

1898: /*@C
1899:    SNESGetPicard - Returns the context for the Picard iteration

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

1903:    Input Parameter:
1904: .  snes - the SNES context

1906:    Output Parameter:
1907: +  r - the function (or NULL)
1908: .  f - the function (or NULL); see SNESFunction for calling sequence details
1909: .  Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
1910: .  Pmat  - the matrix from which the preconditioner will be constructed (or NULL)
1911: .  J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
1912: -  ctx - the function context (or NULL)

1914:    Level: advanced

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

1918: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
1919: @*/
1920: 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)
1921: {
1923:   DM             dm;

1927:   SNESGetFunction(snes,r,NULL,NULL);
1928:   SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
1929:   SNESGetDM(snes,&dm);
1930:   DMSNESGetPicard(dm,f,J,ctx);
1931:   return(0);
1932: }

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

1939:    Logically Collective on SNES

1941:    Input Parameters:
1942: +  snes - the SNES context
1943: .  func - function evaluation routine
1944: -  ctx - [optional] user-defined context for private data for the
1945:          function evaluation routine (may be NULL)

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

1950: .  f - function vector
1951: -  ctx - optional user-defined function context

1953:    Level: intermediate

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

1957: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
1958: @*/
1959: PetscErrorCode  SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
1960: {
1963:   if (func) snes->ops->computeinitialguess = func;
1964:   if (ctx)  snes->initialguessP            = ctx;
1965:   return(0);
1966: }

1968: /* --------------------------------------------------------------- */
1971: /*@C
1972:    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
1973:    it assumes a zero right hand side.

1975:    Logically Collective on SNES

1977:    Input Parameter:
1978: .  snes - the SNES context

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

1983:    Level: intermediate

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

1987: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
1988: @*/
1989: PetscErrorCode  SNESGetRhs(SNES snes,Vec *rhs)
1990: {
1994:   *rhs = snes->vec_rhs;
1995:   return(0);
1996: }

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

2003:    Collective on SNES

2005:    Input Parameters:
2006: +  snes - the SNES context
2007: -  x - input vector

2009:    Output Parameter:
2010: .  y - function vector, as set by SNESSetFunction()

2012:    Notes:
2013:    SNESComputeFunction() is typically used within nonlinear solvers
2014:    implementations, so most users would not generally call this routine
2015:    themselves.

2017:    Level: developer

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

2021: .seealso: SNESSetFunction(), SNESGetFunction()
2022: @*/
2023: PetscErrorCode  SNESComputeFunction(SNES snes,Vec x,Vec y)
2024: {
2026:   DM             dm;
2027:   DMSNES         sdm;

2035:   VecValidValues(x,2,PETSC_TRUE);

2037:   SNESGetDM(snes,&dm);
2038:   DMGetDMSNES(dm,&sdm);
2039:   if (sdm->ops->computefunction) {
2040:     PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2041:     VecLockPush(x);
2042:     PetscStackPush("SNES user function");
2043:     (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2044:     PetscStackPop;
2045:     VecLockPop(x);
2046:     PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2047:   } else if (snes->vec_rhs) {
2048:     MatMult(snes->jacobian, x, y);
2049:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2050:   if (snes->vec_rhs) {
2051:     VecAXPY(y,-1.0,snes->vec_rhs);
2052:   }
2053:   snes->nfuncs++;
2054:   VecValidValues(y,3,PETSC_FALSE);
2055:   return(0);
2056: }

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

2063:    Collective on SNES

2065:    Input Parameters:
2066: +  snes - the SNES context
2067: .  x - input vector
2068: -  b - rhs vector

2070:    Output Parameter:
2071: .  x - new solution vector

2073:    Notes:
2074:    SNESComputeNGS() is typically used within composed nonlinear solver
2075:    implementations, so most users would not generally call this routine
2076:    themselves.

2078:    Level: developer

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

2082: .seealso: SNESSetNGS(), SNESComputeFunction()
2083: @*/
2084: PetscErrorCode  SNESComputeNGS(SNES snes,Vec b,Vec x)
2085: {
2087:   DM             dm;
2088:   DMSNES         sdm;

2096:   if (b) {VecValidValues(b,2,PETSC_TRUE);}
2097:   PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2098:   SNESGetDM(snes,&dm);
2099:   DMGetDMSNES(dm,&sdm);
2100:   if (sdm->ops->computegs) {
2101:     if (b) {VecLockPush(b);}
2102:     PetscStackPush("SNES user NGS");
2103:     (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2104:     PetscStackPop;
2105:     if (b) {VecLockPop(b);}
2106:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2107:   PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2108:   VecValidValues(x,3,PETSC_FALSE);
2109:   return(0);
2110: }

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

2117:    Collective on SNES and Mat

2119:    Input Parameters:
2120: +  snes - the SNES context
2121: -  x - input vector

2123:    Output Parameters:
2124: +  A - Jacobian matrix
2125: -  B - optional preconditioning matrix

2127:   Options Database Keys:
2128: +    -snes_lag_preconditioner <lag>
2129: .    -snes_lag_jacobian <lag>
2130: .    -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2131: .    -snes_compare_explicit_draw  - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2132: .    -snes_compare_explicit_contour  - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2133: .    -snes_compare_operator  - Make the comparison options above use the operator instead of the preconditioning matrix
2134: .    -snes_compare_coloring - Compute the finite differece Jacobian using coloring and display norms of difference
2135: .    -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
2136: .    -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2137: .    -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2138: .    -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2139: .    -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
2140: -    -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences


2143:    Notes:
2144:    Most users should not need to explicitly call this routine, as it
2145:    is used internally within the nonlinear solvers.

2147:    Level: developer

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

2151: .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2152: @*/
2153: PetscErrorCode  SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2154: {
2156:   PetscBool      flag;
2157:   DM             dm;
2158:   DMSNES         sdm;
2159:   KSP            ksp;

2165:   VecValidValues(X,2,PETSC_TRUE);
2166:   SNESGetDM(snes,&dm);
2167:   DMGetDMSNES(dm,&sdm);

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

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

2173:   if (snes->lagjacobian == -2) {
2174:     snes->lagjacobian = -1;

2176:     PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2177:   } else if (snes->lagjacobian == -1) {
2178:     PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2179:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2180:     if (flag) {
2181:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2182:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2183:     }
2184:     return(0);
2185:   } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2186:     PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2187:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2188:     if (flag) {
2189:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2190:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2191:     }
2192:     return(0);
2193:   }
2194:   if (snes->pc && snes->pcside == PC_LEFT) {
2195:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2196:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2197:       return(0);
2198:   }

2200:   PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2201:   VecLockPush(X);
2202:   PetscStackPush("SNES user Jacobian function");
2203:   (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2204:   PetscStackPop;
2205:   VecLockPop(X);
2206:   PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);

2208:   /* the next line ensures that snes->ksp exists */
2209:   SNESGetKSP(snes,&ksp);
2210:   if (snes->lagpreconditioner == -2) {
2211:     PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2212:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2213:     snes->lagpreconditioner = -1;
2214:   } else if (snes->lagpreconditioner == -1) {
2215:     PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2216:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2217:   } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2218:     PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2219:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2220:   } else {
2221:     PetscInfo(snes,"Rebuilding preconditioner\n");
2222:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2223:   }

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

2298:       MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2299:       MatColoringCreate(Bfd,&coloring);
2300:       MatColoringSetType(coloring,MATCOLORINGSL);
2301:       MatColoringSetFromOptions(coloring);
2302:       MatColoringApply(coloring,&iscoloring);
2303:       MatColoringDestroy(&coloring);
2304:       MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2305:       MatFDColoringSetFromOptions(matfdcoloring);
2306:       MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2307:       ISColoringDestroy(&iscoloring);

2309:       /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2310:       SNESGetFunction(snes,NULL,&func,&funcctx);
2311:       MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2312:       PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2313:       PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2314:       MatFDColoringSetFromOptions(matfdcoloring);
2315:       MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2316:       MatFDColoringDestroy(&matfdcoloring);

2318:       PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2319:       if (flag_draw || flag_contour) {
2320:         PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2321:         if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2322:       } else vdraw = NULL;
2323:       PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2324:       if (flag_display) {MatView(B,vstdout);}
2325:       if (vdraw) {MatView(B,vdraw);}
2326:       PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2327:       if (flag_display) {MatView(Bfd,vstdout);}
2328:       if (vdraw) {MatView(Bfd,vdraw);}
2329:       MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2330:       MatNorm(Bfd,NORM_1,&norm1);
2331:       MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2332:       MatNorm(Bfd,NORM_MAX,&normmax);
2333:       PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2334:       if (flag_display) {MatView(Bfd,vstdout);}
2335:       if (vdraw) {              /* Always use contour for the difference */
2336:         PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2337:         MatView(Bfd,vdraw);
2338:         PetscViewerPopFormat(vdraw);
2339:       }
2340:       if (flag_contour) {PetscViewerPopFormat(vdraw);}

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

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

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

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

2403:    Level: intermediate

2405: .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2406: M*/

2410: /*@C
2411:    SNESSetJacobian - Sets the function to compute Jacobian as well as the
2412:    location to store the matrix.

2414:    Logically Collective on SNES and Mat

2416:    Input Parameters:
2417: +  snes - the SNES context
2418: .  Amat - the matrix that defines the (approximate) Jacobian
2419: .  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2420: .  J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
2421: -  ctx - [optional] user-defined context for private data for the
2422:          Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)

2424:    Notes:
2425:    If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2426:    each matrix.

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

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

2434:    Level: beginner

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

2438: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J, 
2439:           SNESSetPicard(), SNESJacobianFunction
2440: @*/
2441: PetscErrorCode  SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2442: {
2444:   DM             dm;

2452:   SNESGetDM(snes,&dm);
2453:   DMSNESSetJacobian(dm,J,ctx);
2454:   if (Amat) {
2455:     PetscObjectReference((PetscObject)Amat);
2456:     MatDestroy(&snes->jacobian);

2458:     snes->jacobian = Amat;
2459:   }
2460:   if (Pmat) {
2461:     PetscObjectReference((PetscObject)Pmat);
2462:     MatDestroy(&snes->jacobian_pre);

2464:     snes->jacobian_pre = Pmat;
2465:   }
2466:   return(0);
2467: }

2471: /*@C
2472:    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2473:    provided context for evaluating the Jacobian.

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

2477:    Input Parameter:
2478: .  snes - the nonlinear solver context

2480:    Output Parameters:
2481: +  Amat - location to stash (approximate) Jacobian matrix (or NULL)
2482: .  Pmat - location to stash matrix used to compute the preconditioner (or NULL)
2483: .  J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
2484: -  ctx - location to stash Jacobian ctx (or NULL)

2486:    Level: advanced

2488: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
2489: @*/
2490: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2491: {
2493:   DM             dm;
2494:   DMSNES         sdm;

2498:   if (Amat) *Amat = snes->jacobian;
2499:   if (Pmat) *Pmat = snes->jacobian_pre;
2500:   SNESGetDM(snes,&dm);
2501:   DMGetDMSNES(dm,&sdm);
2502:   if (J) *J = sdm->ops->computejacobian;
2503:   if (ctx) *ctx = sdm->jacobianctx;
2504:   return(0);
2505: }

2509: /*@
2510:    SNESSetUp - Sets up the internal data structures for the later use
2511:    of a nonlinear solver.

2513:    Collective on SNES

2515:    Input Parameters:
2516: .  snes - the SNES context

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

2525:    Level: advanced

2527: .keywords: SNES, nonlinear, setup

2529: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
2530: @*/
2531: PetscErrorCode  SNESSetUp(SNES snes)
2532: {
2534:   DM             dm;
2535:   DMSNES         sdm;
2536:   SNESLineSearch linesearch, pclinesearch;
2537:   void           *lsprectx,*lspostctx;
2538:   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
2539:   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
2540:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2541:   Vec            f,fpc;
2542:   void           *funcctx;
2543:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2544:   void           *jacctx,*appctx;
2545:   Mat            j,jpre;

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

2551:   if (!((PetscObject)snes)->type_name) {
2552:     SNESSetType(snes,SNESNEWTONLS);
2553:   }

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

2557:   SNESGetDM(snes,&dm);
2558:   DMGetDMSNES(dm,&sdm);
2559:   if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
2560:   if (!sdm->ops->computejacobian) {
2561:     DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
2562:   }
2563:   if (!snes->vec_func) {
2564:     DMCreateGlobalVector(dm,&snes->vec_func);
2565:   }

2567:   if (!snes->ksp) {
2568:     SNESGetKSP(snes, &snes->ksp);
2569:   }

2571:   if (!snes->linesearch) {
2572:     SNESGetLineSearch(snes, &snes->linesearch);
2573:   }
2574:   SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);

2576:   if (snes->pc && (snes->pcside == PC_LEFT)) {
2577:     snes->mf          = PETSC_TRUE;
2578:     snes->mf_operator = PETSC_FALSE;
2579:   }

2581:   if (snes->pc) {
2582:     /* copy the DM over */
2583:     SNESGetDM(snes,&dm);
2584:     SNESSetDM(snes->pc,dm);

2586:     SNESGetFunction(snes,&f,&func,&funcctx);
2587:     VecDuplicate(f,&fpc);
2588:     SNESSetFunction(snes->pc,fpc,func,funcctx);
2589:     SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
2590:     SNESSetJacobian(snes->pc,j,jpre,jac,jacctx);
2591:     SNESGetApplicationContext(snes,&appctx);
2592:     SNESSetApplicationContext(snes->pc,appctx);
2593:     VecDestroy(&fpc);

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

2598:     /* default to 1 iteration */
2599:     SNESSetTolerances(snes->pc,0.0,0.0,0.0,1,snes->pc->max_funcs);
2600:     if (snes->pcside==PC_RIGHT) {
2601:       SNESSetNormSchedule(snes->pc,SNES_NORM_FINAL_ONLY);
2602:     } else {
2603:       SNESSetNormSchedule(snes->pc,SNES_NORM_NONE);
2604:     }
2605:     SNESSetFromOptions(snes->pc);

2607:     /* copy the line search context over */
2608:     SNESGetLineSearch(snes,&linesearch);
2609:     SNESGetLineSearch(snes->pc,&pclinesearch);
2610:     SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
2611:     SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
2612:     SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
2613:     SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
2614:     PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
2615:   }
2616:   if (snes->mf) {
2617:     SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
2618:   }
2619:   if (snes->ops->usercompute && !snes->user) {
2620:     (*snes->ops->usercompute)(snes,(void**)&snes->user);
2621:   }

2623:   snes->jac_iter = 0;
2624:   snes->pre_iter = 0;

2626:   if (snes->ops->setup) {
2627:     (*snes->ops->setup)(snes);
2628:   }

2630:   if (snes->pc && (snes->pcside == PC_LEFT)) {
2631:     if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2632:       SNESGetLineSearch(snes,&linesearch);
2633:       SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
2634:     }
2635:   }

2637:   snes->setupcalled = PETSC_TRUE;
2638:   return(0);
2639: }

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

2646:    Collective on SNES

2648:    Input Parameter:
2649: .  snes - iterative context obtained from SNESCreate()

2651:    Level: intermediate

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

2655: .keywords: SNES, destroy

2657: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
2658: @*/
2659: PetscErrorCode  SNESReset(SNES snes)
2660: {

2665:   if (snes->ops->userdestroy && snes->user) {
2666:     (*snes->ops->userdestroy)((void**)&snes->user);
2667:     snes->user = NULL;
2668:   }
2669:   if (snes->pc) {
2670:     SNESReset(snes->pc);
2671:   }

2673:   if (snes->ops->reset) {
2674:     (*snes->ops->reset)(snes);
2675:   }
2676:   if (snes->ksp) {
2677:     KSPReset(snes->ksp);
2678:   }

2680:   if (snes->linesearch) {
2681:     SNESLineSearchReset(snes->linesearch);
2682:   }

2684:   VecDestroy(&snes->vec_rhs);
2685:   VecDestroy(&snes->vec_sol);
2686:   VecDestroy(&snes->vec_sol_update);
2687:   VecDestroy(&snes->vec_func);
2688:   MatDestroy(&snes->jacobian);
2689:   MatDestroy(&snes->jacobian_pre);
2690:   VecDestroyVecs(snes->nwork,&snes->work);
2691:   VecDestroyVecs(snes->nvwork,&snes->vwork);

2693:   snes->nwork       = snes->nvwork = 0;
2694:   snes->setupcalled = PETSC_FALSE;
2695:   return(0);
2696: }

2700: /*@
2701:    SNESDestroy - Destroys the nonlinear solver context that was created
2702:    with SNESCreate().

2704:    Collective on SNES

2706:    Input Parameter:
2707: .  snes - the SNES context

2709:    Level: beginner

2711: .keywords: SNES, nonlinear, destroy

2713: .seealso: SNESCreate(), SNESSolve()
2714: @*/
2715: PetscErrorCode  SNESDestroy(SNES *snes)
2716: {

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

2724:   SNESReset((*snes));
2725:   SNESDestroy(&(*snes)->pc);

2727:   /* if memory was published with SAWs then destroy it */
2728:   PetscObjectSAWsViewOff((PetscObject)*snes);
2729:   if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}

2731:   DMDestroy(&(*snes)->dm);
2732:   KSPDestroy(&(*snes)->ksp);
2733:   SNESLineSearchDestroy(&(*snes)->linesearch);

2735:   PetscFree((*snes)->kspconvctx);
2736:   if ((*snes)->ops->convergeddestroy) {
2737:     (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
2738:   }
2739:   if ((*snes)->conv_malloc) {
2740:     PetscFree((*snes)->conv_hist);
2741:     PetscFree((*snes)->conv_hist_its);
2742:   }
2743:   SNESMonitorCancel((*snes));
2744:   PetscHeaderDestroy(snes);
2745:   return(0);
2746: }

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

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

2755:    Logically Collective on SNES

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

2762:    Options Database Keys:
2763: .    -snes_lag_preconditioner <lag>

2765:    Notes:
2766:    The default is 1
2767:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2768:    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use

2770:    Level: intermediate

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

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

2776: @*/
2777: PetscErrorCode  SNESSetLagPreconditioner(SNES snes,PetscInt lag)
2778: {
2781:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2782:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2784:   snes->lagpreconditioner = lag;
2785:   return(0);
2786: }

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

2793:    Logically Collective on SNES

2795:    Input Parameters:
2796: +  snes - the SNES context
2797: -  steps - the number of refinements to do, defaults to 0

2799:    Options Database Keys:
2800: .    -snes_grid_sequence <steps>

2802:    Level: intermediate

2804:    Notes:
2805:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

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

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

2811: @*/
2812: PetscErrorCode  SNESSetGridSequence(SNES snes,PetscInt steps)
2813: {
2817:   snes->gridsequence = steps;
2818:   return(0);
2819: }

2823: /*@
2824:    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt

2826:    Not Collective

2828:    Input Parameter:
2829: .  snes - the SNES context

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

2835:    Options Database Keys:
2836: .    -snes_lag_preconditioner <lag>

2838:    Notes:
2839:    The default is 1
2840:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

2842:    Level: intermediate

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

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

2848: @*/
2849: PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
2850: {
2853:   *lag = snes->lagpreconditioner;
2854:   return(0);
2855: }

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

2863:    Logically Collective on SNES

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

2870:    Options Database Keys:
2871: .    -snes_lag_jacobian <lag>

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

2879:    Level: intermediate

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

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

2885: @*/
2886: PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
2887: {
2890:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2891:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2893:   snes->lagjacobian = lag;
2894:   return(0);
2895: }

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

2902:    Not Collective

2904:    Input Parameter:
2905: .  snes - the SNES context

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

2911:    Options Database Keys:
2912: .    -snes_lag_jacobian <lag>

2914:    Notes:
2915:    The default is 1
2916:    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

2918:    Level: intermediate

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

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

2924: @*/
2925: PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
2926: {
2929:   *lag = snes->lagjacobian;
2930:   return(0);
2931: }

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

2938:    Logically collective on SNES

2940:    Input Parameter:
2941: +  snes - the SNES context
2942: -   flg - jacobian lagging persists if true

2944:    Options Database Keys:
2945: .    -snes_lag_jacobian_persists <flg>

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

2951:    Level: developer

2953: .keywords: SNES, nonlinear, lag

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

2957: @*/
2958: PetscErrorCode  SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
2959: {
2963:   snes->lagjac_persist = flg;
2964:   return(0);
2965: }

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

2972:    Logically Collective on SNES

2974:    Input Parameter:
2975: +  snes - the SNES context
2976: -   flg - preconditioner lagging persists if true

2978:    Options Database Keys:
2979: .    -snes_lag_jacobian_persists <flg>

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

2985:    Level: developer

2987: .keywords: SNES, nonlinear, lag

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

2991: @*/
2992: PetscErrorCode  SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
2993: {
2997:   snes->lagpre_persist = flg;
2998:   return(0);
2999: }

3003: /*@
3004:    SNESSetTolerances - Sets various parameters used in convergence tests.

3006:    Logically Collective on SNES

3008:    Input Parameters:
3009: +  snes - the SNES context
3010: .  abstol - absolute convergence tolerance
3011: .  rtol - relative convergence tolerance
3012: .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
3013: .  maxit - maximum number of iterations
3014: -  maxf - maximum number of function evaluations

3016:    Options Database Keys:
3017: +    -snes_atol <abstol> - Sets abstol
3018: .    -snes_rtol <rtol> - Sets rtol
3019: .    -snes_stol <stol> - Sets stol
3020: .    -snes_max_it <maxit> - Sets maxit
3021: -    -snes_max_funcs <maxf> - Sets maxf

3023:    Notes:
3024:    The default maximum number of iterations is 50.
3025:    The default maximum number of function evaluations is 1000.

3027:    Level: intermediate

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

3031: .seealso: SNESSetTrustRegionTolerance()
3032: @*/
3033: PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3034: {

3043:   if (abstol != PETSC_DEFAULT) {
3044:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3045:     snes->abstol = abstol;
3046:   }
3047:   if (rtol != PETSC_DEFAULT) {
3048:     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);
3049:     snes->rtol = rtol;
3050:   }
3051:   if (stol != PETSC_DEFAULT) {
3052:     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3053:     snes->stol = stol;
3054:   }
3055:   if (maxit != PETSC_DEFAULT) {
3056:     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3057:     snes->max_its = maxit;
3058:   }
3059:   if (maxf != PETSC_DEFAULT) {
3060:     if (maxf < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
3061:     snes->max_funcs = maxf;
3062:   }
3063:   snes->tolerancesset = PETSC_TRUE;
3064:   return(0);
3065: }

3069: /*@
3070:    SNESGetTolerances - Gets various parameters used in convergence tests.

3072:    Not Collective

3074:    Input Parameters:
3075: +  snes - the SNES context
3076: .  atol - absolute convergence tolerance
3077: .  rtol - relative convergence tolerance
3078: .  stol -  convergence tolerance in terms of the norm
3079:            of the change in the solution between steps
3080: .  maxit - maximum number of iterations
3081: -  maxf - maximum number of function evaluations

3083:    Notes:
3084:    The user can specify NULL for any parameter that is not needed.

3086:    Level: intermediate

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

3090: .seealso: SNESSetTolerances()
3091: @*/
3092: PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3093: {
3096:   if (atol)  *atol  = snes->abstol;
3097:   if (rtol)  *rtol  = snes->rtol;
3098:   if (stol)  *stol  = snes->stol;
3099:   if (maxit) *maxit = snes->max_its;
3100:   if (maxf)  *maxf  = snes->max_funcs;
3101:   return(0);
3102: }

3106: /*@
3107:    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.

3109:    Logically Collective on SNES

3111:    Input Parameters:
3112: +  snes - the SNES context
3113: -  tol - tolerance

3115:    Options Database Key:
3116: .  -snes_trtol <tol> - Sets tol

3118:    Level: intermediate

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

3122: .seealso: SNESSetTolerances()
3123: @*/
3124: PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3125: {
3129:   snes->deltatol = tol;
3130:   return(0);
3131: }

3133: /*
3134:    Duplicate the lg monitors for SNES from KSP; for some reason with
3135:    dynamic libraries things don't work under Sun4 if we just use
3136:    macros instead of functions
3137: */
3140: PetscErrorCode  SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,PetscObject *objs)
3141: {

3146:   KSPMonitorLGResidualNorm((KSP)snes,it,norm,objs);
3147:   return(0);
3148: }

3152: PetscErrorCode  SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscObject **draw)
3153: {

3157:   KSPMonitorLGResidualNormCreate(host,label,x,y,m,n,draw);
3158:   return(0);
3159: }

3163: PetscErrorCode  SNESMonitorLGDestroy(PetscObject **objs)
3164: {

3168:   KSPMonitorLGResidualNormDestroy(objs);
3169:   return(0);
3170: }

3172: extern PetscErrorCode  SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3175: PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3176: {
3177:   PetscDrawLG      lg;
3178:   PetscErrorCode   ierr;
3179:   PetscReal        x,y,per;
3180:   PetscViewer      v = (PetscViewer)monctx;
3181:   static PetscReal prev; /* should be in the context */
3182:   PetscDraw        draw;

3185:   PetscViewerDrawGetDrawLG(v,0,&lg);
3186:   if (!n) {PetscDrawLGReset(lg);}
3187:   PetscDrawLGGetDraw(lg,&draw);
3188:   PetscDrawSetTitle(draw,"Residual norm");
3189:   x    = (PetscReal)n;
3190:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3191:   else y = -15.0;
3192:   PetscDrawLGAddPoint(lg,&x,&y);
3193:   if (n < 20 || !(n % 5)) {
3194:     PetscDrawLGDraw(lg);
3195:   }

3197:   PetscViewerDrawGetDrawLG(v,1,&lg);
3198:   if (!n) {PetscDrawLGReset(lg);}
3199:   PetscDrawLGGetDraw(lg,&draw);
3200:   PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3201:    SNESMonitorRange_Private(snes,n,&per);
3202:   x    = (PetscReal)n;
3203:   y    = 100.0*per;
3204:   PetscDrawLGAddPoint(lg,&x,&y);
3205:   if (n < 20 || !(n % 5)) {
3206:     PetscDrawLGDraw(lg);
3207:   }

3209:   PetscViewerDrawGetDrawLG(v,2,&lg);
3210:   if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3211:   PetscDrawLGGetDraw(lg,&draw);
3212:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3213:   x    = (PetscReal)n;
3214:   y    = (prev - rnorm)/prev;
3215:   PetscDrawLGAddPoint(lg,&x,&y);
3216:   if (n < 20 || !(n % 5)) {
3217:     PetscDrawLGDraw(lg);
3218:   }

3220:   PetscViewerDrawGetDrawLG(v,3,&lg);
3221:   if (!n) {PetscDrawLGReset(lg);}
3222:   PetscDrawLGGetDraw(lg,&draw);
3223:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3224:   x    = (PetscReal)n;
3225:   y    = (prev - rnorm)/(prev*per);
3226:   if (n > 2) { /*skip initial crazy value */
3227:     PetscDrawLGAddPoint(lg,&x,&y);
3228:   }
3229:   if (n < 20 || !(n % 5)) {
3230:     PetscDrawLGDraw(lg);
3231:   }
3232:   prev = rnorm;
3233:   return(0);
3234: }

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

3241:    Collective on SNES

3243:    Input Parameters:
3244: +  snes - nonlinear solver context obtained from SNESCreate()
3245: .  iter - iteration number
3246: -  rnorm - relative norm of the residual

3248:    Notes:
3249:    This routine is called by the SNES implementations.
3250:    It does not typically need to be called by the user.

3252:    Level: developer

3254: .seealso: SNESMonitorSet()
3255: @*/
3256: PetscErrorCode  SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3257: {
3259:   PetscInt       i,n = snes->numbermonitors;

3262:   VecLockPush(snes->vec_sol);
3263:   for (i=0; i<n; i++) {
3264:     (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3265:   }
3266:   VecLockPop(snes->vec_sol);
3267:   return(0);
3268: }

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

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

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

3279: +    snes - the SNES context
3280: .    its - iteration number
3281: .    norm - 2-norm function value (may be estimated)
3282: -    mctx - [optional] monitoring context

3284:    Level: advanced

3286: .seealso:   SNESMonitorSet(), SNESMonitorGet()
3287: M*/

3291: /*@C
3292:    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3293:    iteration of the nonlinear solver to display the iteration's
3294:    progress.

3296:    Logically Collective on SNES

3298:    Input Parameters:
3299: +  snes - the SNES context
3300: .  f - the monitor function, see SNESMonitorFunction for the calling sequence
3301: .  mctx - [optional] user-defined context for private data for the
3302:           monitor routine (use NULL if no context is desired)
3303: -  monitordestroy - [optional] routine that frees monitor context
3304:           (may be NULL)

3306:    Options Database Keys:
3307: +    -snes_monitor        - sets SNESMonitorDefault()
3308: .    -snes_monitor_lg_residualnorm    - sets line graph monitor,
3309:                             uses SNESMonitorLGCreate()
3310: -    -snes_monitor_cancel - cancels all monitors that have
3311:                             been hardwired into a code by
3312:                             calls to SNESMonitorSet(), but
3313:                             does not cancel those set via
3314:                             the options database.

3316:    Notes:
3317:    Several different monitoring routines may be set by calling
3318:    SNESMonitorSet() multiple times; all will be called in the
3319:    order in which they were set.

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

3323:    Level: intermediate

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

3327: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3328: @*/
3329: PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3330: {
3331:   PetscInt       i;

3336:   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3337:   for (i=0; i<snes->numbermonitors;i++) {
3338:     if (f == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) {
3339:       if (monitordestroy) {
3340:         (*monitordestroy)(&mctx);
3341:       }
3342:       return(0);
3343:     }
3344:   }
3345:   snes->monitor[snes->numbermonitors]          = f;
3346:   snes->monitordestroy[snes->numbermonitors]   = monitordestroy;
3347:   snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3348:   return(0);
3349: }

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

3356:    Logically Collective on SNES

3358:    Input Parameters:
3359: .  snes - the SNES context

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

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

3369:    Level: intermediate

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

3373: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3374: @*/
3375: PetscErrorCode  SNESMonitorCancel(SNES snes)
3376: {
3378:   PetscInt       i;

3382:   for (i=0; i<snes->numbermonitors; i++) {
3383:     if (snes->monitordestroy[i]) {
3384:       (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3385:     }
3386:   }
3387:   snes->numbermonitors = 0;
3388:   return(0);
3389: }

3391: /*MC
3392:     SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver

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

3398: +    snes - the SNES context
3399: .    it - current iteration (0 is the first and is before any Newton step)
3400: .    cctx - [optional] convergence context
3401: .    reason - reason for convergence/divergence
3402: .    xnorm - 2-norm of current iterate
3403: .    gnorm - 2-norm of current step
3404: -    f - 2-norm of function

3406:    Level: intermediate

3408: .seealso:   SNESSetConvergenceTest(), SNESGetConvergenceTest()
3409: M*/

3413: /*@C
3414:    SNESSetConvergenceTest - Sets the function that is to be used
3415:    to test for convergence of the nonlinear iterative solution.

3417:    Logically Collective on SNES

3419:    Input Parameters:
3420: +  snes - the SNES context
3421: .  SNESConvergenceTestFunction - routine to test for convergence
3422: .  cctx - [optional] context for private data for the convergence routine  (may be NULL)
3423: -  destroy - [optional] destructor for the context (may be NULL; NULL_FUNCTION in Fortran)

3425:    Level: advanced

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

3429: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3430: @*/
3431: PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3432: {

3437:   if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
3438:   if (snes->ops->convergeddestroy) {
3439:     (*snes->ops->convergeddestroy)(snes->cnvP);
3440:   }
3441:   snes->ops->converged        = SNESConvergenceTestFunction;
3442:   snes->ops->convergeddestroy = destroy;
3443:   snes->cnvP                  = cctx;
3444:   return(0);
3445: }

3449: /*@
3450:    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.

3452:    Not Collective

3454:    Input Parameter:
3455: .  snes - the SNES context

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

3461:    Level: intermediate

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

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

3467: .seealso: SNESSetConvergenceTest(), SNESConvergedReason
3468: @*/
3469: PetscErrorCode  SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
3470: {
3474:   *reason = snes->reason;
3475:   return(0);
3476: }

3480: /*@
3481:    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.

3483:    Logically Collective on SNES

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

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

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

3501:    Level: intermediate

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

3505: .seealso: SNESGetConvergenceHistory()

3507: @*/
3508: PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
3509: {

3516:   if (!a) {
3517:     if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
3518:     PetscCalloc1(na,&a);
3519:     PetscCalloc1(na,&its);

3521:     snes->conv_malloc = PETSC_TRUE;
3522:   }
3523:   snes->conv_hist       = a;
3524:   snes->conv_hist_its   = its;
3525:   snes->conv_hist_max   = na;
3526:   snes->conv_hist_len   = 0;
3527:   snes->conv_hist_reset = reset;
3528:   return(0);
3529: }

3531: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3532: #include <engine.h>   /* MATLAB include file */
3533: #include <mex.h>      /* MATLAB include file */

3537: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
3538: {
3539:   mxArray   *mat;
3540:   PetscInt  i;
3541:   PetscReal *ar;

3544:   mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
3545:   ar  = (PetscReal*) mxGetData(mat);
3546:   for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
3547:   PetscFunctionReturn(mat);
3548: }
3549: #endif

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

3556:    Not Collective

3558:    Input Parameter:
3559: .  snes - iterative context obtained from SNESCreate()

3561:    Output Parameters:
3562: .  a   - array to hold history
3563: .  its - integer array holds the number of linear iterations (or
3564:          negative if not converged) for each solve.
3565: -  na  - size of a and its

3567:    Notes:
3568:     The calling sequence for this routine in Fortran is
3569: $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)

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

3575:    Level: intermediate

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

3579: .seealso: SNESSetConvergencHistory()

3581: @*/
3582: PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
3583: {
3586:   if (a)   *a   = snes->conv_hist;
3587:   if (its) *its = snes->conv_hist_its;
3588:   if (na)  *na  = snes->conv_hist_len;
3589:   return(0);
3590: }

3594: /*@C
3595:   SNESSetUpdate - Sets the general-purpose update function called
3596:   at the beginning of every iteration of the nonlinear solve. Specifically
3597:   it is called just before the Jacobian is "evaluated".

3599:   Logically Collective on SNES

3601:   Input Parameters:
3602: . snes - The nonlinear solver context
3603: . func - The function

3605:   Calling sequence of func:
3606: . func (SNES snes, PetscInt step);

3608: . step - The current step of the iteration

3610:   Level: advanced

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

3615: .keywords: SNES, update

3617: .seealso SNESSetJacobian(), SNESSolve()
3618: @*/
3619: PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
3620: {
3623:   snes->ops->update = func;
3624:   return(0);
3625: }

3629: /*
3630:    SNESScaleStep_Private - Scales a step so that its length is less than the
3631:    positive parameter delta.

3633:     Input Parameters:
3634: +   snes - the SNES context
3635: .   y - approximate solution of linear system
3636: .   fnorm - 2-norm of current function
3637: -   delta - trust region size

3639:     Output Parameters:
3640: +   gpnorm - predicted function norm at the new point, assuming local
3641:     linearization.  The value is zero if the step lies within the trust
3642:     region, and exceeds zero otherwise.
3643: -   ynorm - 2-norm of the step

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

3649: .keywords: SNES, nonlinear, scale, step
3650: */
3651: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
3652: {
3653:   PetscReal      nrm;
3654:   PetscScalar    cnorm;


3662:   VecNorm(y,NORM_2,&nrm);
3663:   if (nrm > *delta) {
3664:     nrm     = *delta/nrm;
3665:     *gpnorm = (1.0 - nrm)*(*fnorm);
3666:     cnorm   = nrm;
3667:     VecScale(y,cnorm);
3668:     *ynorm  = *delta;
3669:   } else {
3670:     *gpnorm = 0.0;
3671:     *ynorm  = nrm;
3672:   }
3673:   return(0);
3674: }

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

3681:    Collective on SNES

3683:    Parameter:
3684: +  snes - iterative context obtained from SNESCreate()
3685: -  viewer - the viewer to display the reason


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

3691:    Level: beginner

3693: .keywords: SNES, solve, linear system

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

3697: @*/
3698: PetscErrorCode  SNESReasonView(SNES snes,PetscViewer viewer)
3699: {
3701:   PetscBool      isAscii;

3704:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
3705:   if (isAscii) {
3706:     PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
3707:     if (snes->reason > 0) {
3708:       if (((PetscObject) snes)->prefix) {
3709:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3710:       } else {
3711:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3712:       }
3713:     } else {
3714:       if (((PetscObject) snes)->prefix) {
3715:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3716:       } else {
3717:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3718:       }
3719:     }
3720:     PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
3721:   }
3722:   return(0);
3723: }

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

3730:   Collective on SNES

3732:   Input Parameters:
3733: . snes   - the SNES object

3735:   Level: intermediate

3737: @*/
3738: PetscErrorCode SNESReasonViewFromOptions(SNES snes)
3739: {
3740:   PetscErrorCode    ierr;
3741:   PetscViewer       viewer;
3742:   PetscBool         flg;
3743:   static PetscBool  incall = PETSC_FALSE;
3744:   PetscViewerFormat format;

3747:   if (incall) return(0);
3748:   incall = PETSC_TRUE;
3749:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
3750:   if (flg) {
3751:     PetscViewerPushFormat(viewer,format);
3752:     SNESReasonView(snes,viewer);
3753:     PetscViewerPopFormat(viewer);
3754:     PetscViewerDestroy(&viewer);
3755:   }
3756:   incall = PETSC_FALSE;
3757:   return(0);
3758: }

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

3766:    Collective on SNES

3768:    Input Parameters:
3769: +  snes - the SNES context
3770: .  b - the constant part of the equation F(x) = b, or NULL to use zero.
3771: -  x - the solution vector.

3773:    Notes:
3774:    The user should initialize the vector,x, with the initial guess
3775:    for the nonlinear solve prior to calling SNESSolve.  In particular,
3776:    to employ an initial guess of zero, the user should explicitly set
3777:    this vector to zero by calling VecSet().

3779:    Level: beginner

3781: .keywords: SNES, nonlinear, solve

3783: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
3784: @*/
3785: PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
3786: {
3787:   PetscErrorCode    ierr;
3788:   PetscBool         flg;
3789:   PetscInt          grid;
3790:   Vec               xcreated = NULL;
3791:   DM                dm;


3800:   if (!x) {
3801:     SNESGetDM(snes,&dm);
3802:     DMCreateGlobalVector(dm,&xcreated);
3803:     x    = xcreated;
3804:   }
3805:   SNESViewFromOptions(snes,NULL,"-snes_view_pre");

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

3810:     /* set solution vector */
3811:     if (!grid) {PetscObjectReference((PetscObject)x);}
3812:     VecDestroy(&snes->vec_sol);
3813:     snes->vec_sol = x;
3814:     SNESGetDM(snes,&dm);

3816:     /* set affine vector if provided */
3817:     if (b) { PetscObjectReference((PetscObject)b); }
3818:     VecDestroy(&snes->vec_rhs);
3819:     snes->vec_rhs = b;

3821:     if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
3822:     if (snes->vec_rhs  == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
3823:     if (!snes->vec_sol_update /* && snes->vec_sol */) {
3824:       VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
3825:       PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
3826:     }
3827:     DMShellSetGlobalVector(dm,snes->vec_sol);
3828:     SNESSetUp(snes);

3830:     if (!grid) {
3831:       if (snes->ops->computeinitialguess) {
3832:         (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
3833:       }
3834:     }

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

3839:     PetscLogEventBegin(SNES_Solve,snes,0,0,0);
3840:     (*snes->ops->solve)(snes);
3841:     PetscLogEventEnd(SNES_Solve,snes,0,0,0);
3842:     if (snes->domainerror) {
3843:       snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
3844:       snes->domainerror = PETSC_FALSE;
3845:     }
3846:     if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");

3848:     if (snes->lagjac_persist) snes->jac_iter += snes->iter;
3849:     if (snes->lagpre_persist) snes->pre_iter += snes->iter;

3851:     flg  = PETSC_FALSE;
3852:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,NULL);
3853:     if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
3854:     SNESReasonViewFromOptions(snes);

3856:     if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
3857:     if (grid <  snes->gridsequence) {
3858:       DM  fine;
3859:       Vec xnew;
3860:       Mat interp;

3862:       DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
3863:       if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
3864:       DMCreateInterpolation(snes->dm,fine,&interp,NULL);
3865:       DMCreateGlobalVector(fine,&xnew);
3866:       MatInterpolate(interp,x,xnew);
3867:       DMInterpolate(snes->dm,interp,fine);
3868:       MatDestroy(&interp);
3869:       x    = xnew;

3871:       SNESReset(snes);
3872:       SNESSetDM(snes,fine);
3873:       DMDestroy(&fine);
3874:       PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
3875:     }
3876:   }
3877:   SNESViewFromOptions(snes,NULL,"-snes_view");
3878:   VecViewFromOptions(snes->vec_sol,((PetscObject)snes)->prefix,"-snes_view_solution");

3880:   VecDestroy(&xcreated);
3881:   PetscObjectSAWsBlock((PetscObject)snes);
3882:   return(0);
3883: }

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

3889: /*@C
3890:    SNESSetType - Sets the method for the nonlinear solver.

3892:    Collective on SNES

3894:    Input Parameters:
3895: +  snes - the SNES context
3896: -  type - a known method

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

3902:    Notes:
3903:    See "petsc/include/petscsnes.h" for available methods (for instance)
3904: +    SNESNEWTONLS - Newton's method with line search
3905:      (systems of nonlinear equations)
3906: .    SNESNEWTONTR - Newton's method with trust region
3907:      (systems of nonlinear equations)

3909:   Normally, it is best to use the SNESSetFromOptions() command and then
3910:   set the SNES solver type from the options database rather than by using
3911:   this routine.  Using the options database provides the user with
3912:   maximum flexibility in evaluating the many nonlinear solvers.
3913:   The SNESSetType() routine is provided for those situations where it
3914:   is necessary to set the nonlinear solver independently of the command
3915:   line or options database.  This might be the case, for example, when
3916:   the choice of solver changes during the execution of the program,
3917:   and the user's application is taking responsibility for choosing the
3918:   appropriate method.

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

3923:   Level: intermediate

3925: .keywords: SNES, set, type

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

3929: @*/
3930: PetscErrorCode  SNESSetType(SNES snes,SNESType type)
3931: {
3932:   PetscErrorCode ierr,(*r)(SNES);
3933:   PetscBool      match;


3939:   PetscObjectTypeCompare((PetscObject)snes,type,&match);
3940:   if (match) return(0);

3942:    PetscFunctionListFind(SNESList,type,&r);
3943:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
3944:   /* Destroy the previous private SNES context */
3945:   if (snes->ops->destroy) {
3946:     (*(snes)->ops->destroy)(snes);
3947:     snes->ops->destroy = NULL;
3948:   }
3949:   /* Reinitialize function pointers in SNESOps structure */
3950:   snes->ops->setup          = 0;
3951:   snes->ops->solve          = 0;
3952:   snes->ops->view           = 0;
3953:   snes->ops->setfromoptions = 0;
3954:   snes->ops->destroy        = 0;
3955:   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
3956:   snes->setupcalled = PETSC_FALSE;

3958:   PetscObjectChangeTypeName((PetscObject)snes,type);
3959:   (*r)(snes);
3960:   return(0);
3961: }

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

3968:    Not Collective

3970:    Input Parameter:
3971: .  snes - nonlinear solver context

3973:    Output Parameter:
3974: .  type - SNES method (a character string)

3976:    Level: intermediate

3978: .keywords: SNES, nonlinear, get, type, name
3979: @*/
3980: PetscErrorCode  SNESGetType(SNES snes,SNESType *type)
3981: {
3985:   *type = ((PetscObject)snes)->type_name;
3986:   return(0);
3987: }

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

3994:   Logically Collective on SNES and Vec

3996:   Input Parameters:
3997: + snes - the SNES context obtained from SNESCreate()
3998: - u    - the solution vector

4000:   Level: beginner

4002: .keywords: SNES, set, solution
4003: @*/
4004: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4005: {
4006:   DM             dm;

4012:   PetscObjectReference((PetscObject) u);
4013:   VecDestroy(&snes->vec_sol);

4015:   snes->vec_sol = u;

4017:   SNESGetDM(snes, &dm);
4018:   DMShellSetGlobalVector(dm, u);
4019:   return(0);
4020: }

4024: /*@
4025:    SNESGetSolution - Returns the vector where the approximate solution is
4026:    stored. This is the fine grid solution when using SNESSetGridSequence().

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

4030:    Input Parameter:
4031: .  snes - the SNES context

4033:    Output Parameter:
4034: .  x - the solution

4036:    Level: intermediate

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

4040: .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
4041: @*/
4042: PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
4043: {
4047:   *x = snes->vec_sol;
4048:   return(0);
4049: }

4053: /*@
4054:    SNESGetSolutionUpdate - Returns the vector where the solution update is
4055:    stored.

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

4059:    Input Parameter:
4060: .  snes - the SNES context

4062:    Output Parameter:
4063: .  x - the solution update

4065:    Level: advanced

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

4069: .seealso: SNESGetSolution(), SNESGetFunction()
4070: @*/
4071: PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
4072: {
4076:   *x = snes->vec_sol_update;
4077:   return(0);
4078: }

4082: /*@C
4083:    SNESGetFunction - Returns the vector where the function is stored.

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

4087:    Input Parameter:
4088: .  snes - the SNES context

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

4095:    Level: advanced

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

4099: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4100: @*/
4101: PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4102: {
4104:   DM             dm;

4108:   if (r) {
4109:     if (!snes->vec_func) {
4110:       if (snes->vec_rhs) {
4111:         VecDuplicate(snes->vec_rhs,&snes->vec_func);
4112:       } else if (snes->vec_sol) {
4113:         VecDuplicate(snes->vec_sol,&snes->vec_func);
4114:       } else if (snes->dm) {
4115:         DMCreateGlobalVector(snes->dm,&snes->vec_func);
4116:       }
4117:     }
4118:     *r = snes->vec_func;
4119:   }
4120:   SNESGetDM(snes,&dm);
4121:   DMSNESGetFunction(dm,f,ctx);
4122:   return(0);
4123: }

4125: /*@C
4126:    SNESGetNGS - Returns the NGS function and context.

4128:    Input Parameter:
4129: .  snes - the SNES context

4131:    Output Parameter:
4132: +  f - the function (or NULL) see SNESNGSFunction for details
4133: -  ctx    - the function context (or NULL)

4135:    Level: advanced

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

4139: .seealso: SNESSetNGS(), SNESGetFunction()
4140: @*/

4144: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4145: {
4147:   DM             dm;

4151:   SNESGetDM(snes,&dm);
4152:   DMSNESGetNGS(dm,f,ctx);
4153:   return(0);
4154: }

4158: /*@C
4159:    SNESSetOptionsPrefix - Sets the prefix used for searching for all
4160:    SNES options in the database.

4162:    Logically Collective on SNES

4164:    Input Parameter:
4165: +  snes - the SNES context
4166: -  prefix - the prefix to prepend to all option names

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

4172:    Level: advanced

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

4176: .seealso: SNESSetFromOptions()
4177: @*/
4178: PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
4179: {

4184:   PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4185:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4186:   if (snes->linesearch) {
4187:     SNESGetLineSearch(snes,&snes->linesearch);
4188:     PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4189:   }
4190:   KSPSetOptionsPrefix(snes->ksp,prefix);
4191:   return(0);
4192: }

4196: /*@C
4197:    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4198:    SNES options in the database.

4200:    Logically Collective on SNES

4202:    Input Parameters:
4203: +  snes - the SNES context
4204: -  prefix - the prefix to prepend to all option names

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

4210:    Level: advanced

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

4214: .seealso: SNESGetOptionsPrefix()
4215: @*/
4216: PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4217: {

4222:   PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4223:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4224:   if (snes->linesearch) {
4225:     SNESGetLineSearch(snes,&snes->linesearch);
4226:     PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4227:   }
4228:   KSPAppendOptionsPrefix(snes->ksp,prefix);
4229:   return(0);
4230: }

4234: /*@C
4235:    SNESGetOptionsPrefix - Sets the prefix used for searching for all
4236:    SNES options in the database.

4238:    Not Collective

4240:    Input Parameter:
4241: .  snes - the SNES context

4243:    Output Parameter:
4244: .  prefix - pointer to the prefix string used

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

4249:    Level: advanced

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

4253: .seealso: SNESAppendOptionsPrefix()
4254: @*/
4255: PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4256: {

4261:   PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4262:   return(0);
4263: }


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

4271:    Not collective

4273:    Input Parameters:
4274: +  name_solver - name of a new user-defined solver
4275: -  routine_create - routine to create method context

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

4280:    Sample usage:
4281: .vb
4282:    SNESRegister("my_solver",MySolverCreate);
4283: .ve

4285:    Then, your solver can be chosen with the procedural interface via
4286: $     SNESSetType(snes,"my_solver")
4287:    or at runtime via the option
4288: $     -snes_type my_solver

4290:    Level: advanced

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

4294: .keywords: SNES, nonlinear, register

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

4298:   Level: advanced
4299: @*/
4300: PetscErrorCode  SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4301: {

4305:   PetscFunctionListAdd(&SNESList,sname,function);
4306:   return(0);
4307: }

4311: PetscErrorCode  SNESTestLocalMin(SNES snes)
4312: {
4314:   PetscInt       N,i,j;
4315:   Vec            u,uh,fh;
4316:   PetscScalar    value;
4317:   PetscReal      norm;

4320:   SNESGetSolution(snes,&u);
4321:   VecDuplicate(u,&uh);
4322:   VecDuplicate(u,&fh);

4324:   /* currently only works for sequential */
4325:   PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4326:   VecGetSize(u,&N);
4327:   for (i=0; i<N; i++) {
4328:     VecCopy(u,uh);
4329:     PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4330:     for (j=-10; j<11; j++) {
4331:       value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4332:       VecSetValue(uh,i,value,ADD_VALUES);
4333:       SNESComputeFunction(snes,uh,fh);
4334:       VecNorm(fh,NORM_2,&norm);
4335:       PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);
4336:       value = -value;
4337:       VecSetValue(uh,i,value,ADD_VALUES);
4338:     }
4339:   }
4340:   VecDestroy(&uh);
4341:   VecDestroy(&fh);
4342:   return(0);
4343: }

4347: /*@
4348:    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4349:    computing relative tolerance for linear solvers within an inexact
4350:    Newton method.

4352:    Logically Collective on SNES

4354:    Input Parameters:
4355: +  snes - SNES context
4356: -  flag - PETSC_TRUE or PETSC_FALSE

4358:     Options Database:
4359: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4360: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
4361: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4362: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4363: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
4364: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
4365: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4366: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

4368:    Notes:
4369:    Currently, the default is to use a constant relative tolerance for
4370:    the inner linear solvers.  Alternatively, one can use the
4371:    Eisenstat-Walker method, where the relative convergence tolerance
4372:    is reset at each Newton iteration according progress of the nonlinear
4373:    solver.

4375:    Level: advanced

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

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

4383: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4384: @*/
4385: PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool flag)
4386: {
4390:   snes->ksp_ewconv = flag;
4391:   return(0);
4392: }

4396: /*@
4397:    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
4398:    for computing relative tolerance for linear solvers within an
4399:    inexact Newton method.

4401:    Not Collective

4403:    Input Parameter:
4404: .  snes - SNES context

4406:    Output Parameter:
4407: .  flag - PETSC_TRUE or PETSC_FALSE

4409:    Notes:
4410:    Currently, the default is to use a constant relative tolerance for
4411:    the inner linear solvers.  Alternatively, one can use the
4412:    Eisenstat-Walker method, where the relative convergence tolerance
4413:    is reset at each Newton iteration according progress of the nonlinear
4414:    solver.

4416:    Level: advanced

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

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

4424: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4425: @*/
4426: PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
4427: {
4431:   *flag = snes->ksp_ewconv;
4432:   return(0);
4433: }

4437: /*@
4438:    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
4439:    convergence criteria for the linear solvers within an inexact
4440:    Newton method.

4442:    Logically Collective on SNES

4444:    Input Parameters:
4445: +    snes - SNES context
4446: .    version - version 1, 2 (default is 2) or 3
4447: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4448: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4449: .    gamma - multiplicative factor for version 2 rtol computation
4450:              (0 <= gamma2 <= 1)
4451: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
4452: .    alpha2 - power for safeguard
4453: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

4455:    Note:
4456:    Version 3 was contributed by Luis Chacon, June 2006.

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

4460:    Level: advanced

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

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

4469: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
4470: @*/
4471: PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
4472: {
4473:   SNESKSPEW *kctx;

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

4487:   if (version != PETSC_DEFAULT)   kctx->version   = version;
4488:   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
4489:   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
4490:   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
4491:   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
4492:   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
4493:   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;

4495:   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);
4496:   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);
4497:   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);
4498:   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);
4499:   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);
4500:   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);
4501:   return(0);
4502: }

4506: /*@
4507:    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
4508:    convergence criteria for the linear solvers within an inexact
4509:    Newton method.

4511:    Not Collective

4513:    Input Parameters:
4514:      snes - SNES context

4516:    Output Parameters:
4517: +    version - version 1, 2 (default is 2) or 3
4518: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4519: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4520: .    gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
4521: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
4522: .    alpha2 - power for safeguard
4523: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

4525:    Level: advanced

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

4529: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
4530: @*/
4531: PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
4532: {
4533:   SNESKSPEW *kctx;

4537:   kctx = (SNESKSPEW*)snes->kspconvctx;
4538:   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4539:   if (version)   *version   = kctx->version;
4540:   if (rtol_0)    *rtol_0    = kctx->rtol_0;
4541:   if (rtol_max)  *rtol_max  = kctx->rtol_max;
4542:   if (gamma)     *gamma     = kctx->gamma;
4543:   if (alpha)     *alpha     = kctx->alpha;
4544:   if (alpha2)    *alpha2    = kctx->alpha2;
4545:   if (threshold) *threshold = kctx->threshold;
4546:   return(0);
4547: }

4551:  PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4552: {
4554:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
4555:   PetscReal      rtol  = PETSC_DEFAULT,stol;

4558:   if (!snes->ksp_ewconv) return(0);
4559:   if (!snes->iter) {
4560:     rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
4561:     VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
4562:   }
4563:   else {
4564:     if (kctx->version == 1) {
4565:       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
4566:       if (rtol < 0.0) rtol = -rtol;
4567:       stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
4568:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4569:     } else if (kctx->version == 2) {
4570:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4571:       stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
4572:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4573:     } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
4574:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4575:       /* safeguard: avoid sharp decrease of rtol */
4576:       stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
4577:       stol = PetscMax(rtol,stol);
4578:       rtol = PetscMin(kctx->rtol_0,stol);
4579:       /* safeguard: avoid oversolving */
4580:       stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
4581:       stol = PetscMax(rtol,stol);
4582:       rtol = PetscMin(kctx->rtol_0,stol);
4583:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
4584:   }
4585:   /* safeguard: avoid rtol greater than one */
4586:   rtol = PetscMin(rtol,kctx->rtol_max);
4587:   KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
4588:   PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
4589:   return(0);
4590: }

4594: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4595: {
4597:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
4598:   PCSide         pcside;
4599:   Vec            lres;

4602:   if (!snes->ksp_ewconv) return(0);
4603:   KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
4604:   kctx->norm_last = snes->norm;
4605:   if (kctx->version == 1) {
4606:     PC        pc;
4607:     PetscBool isNone;

4609:     KSPGetPC(ksp, &pc);
4610:     PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
4611:     KSPGetPCSide(ksp,&pcside);
4612:      if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
4613:       /* KSP residual is true linear residual */
4614:       KSPGetResidualNorm(ksp,&kctx->lresid_last);
4615:     } else {
4616:       /* KSP residual is preconditioned residual */
4617:       /* compute true linear residual norm */
4618:       VecDuplicate(b,&lres);
4619:       MatMult(snes->jacobian,x,lres);
4620:       VecAYPX(lres,-1.0,b);
4621:       VecNorm(lres,NORM_2,&kctx->lresid_last);
4622:       VecDestroy(&lres);
4623:     }
4624:   }
4625:   return(0);
4626: }

4630: /*@
4631:    SNESGetKSP - Returns the KSP context for a SNES solver.

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

4635:    Input Parameter:
4636: .  snes - the SNES context

4638:    Output Parameter:
4639: .  ksp - the KSP context

4641:    Notes:
4642:    The user can then directly manipulate the KSP context to set various
4643:    options, etc.  Likewise, the user can then extract and manipulate the
4644:    PC contexts as well.

4646:    Level: beginner

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

4650: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
4651: @*/
4652: PetscErrorCode  SNESGetKSP(SNES snes,KSP *ksp)
4653: {


4660:   if (!snes->ksp) {
4661:     PetscBool monitor = PETSC_FALSE;

4663:     KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
4664:     PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
4665:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);

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

4670:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
4671:     if (monitor) {
4672:       KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
4673:     }
4674:     monitor = PETSC_FALSE;
4675:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
4676:     if (monitor) {
4677:       PetscObject *objs;
4678:       KSPMonitorSNESLGResidualNormCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
4679:       objs[0] = (PetscObject) snes;
4680:       KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
4681:     }
4682:   }
4683:   *ksp = snes->ksp;
4684:   return(0);
4685: }


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

4694:    Logically Collective on SNES

4696:    Input Parameters:
4697: +  snes - the preconditioner context
4698: -  dm - the dm

4700:    Level: intermediate

4702: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
4703: @*/
4704: PetscErrorCode  SNESSetDM(SNES snes,DM dm)
4705: {
4707:   KSP            ksp;
4708:   DMSNES         sdm;

4712:   if (dm) {PetscObjectReference((PetscObject)dm);}
4713:   if (snes->dm) {               /* Move the DMSNES context over to the new DM unless the new DM already has one */
4714:     if (snes->dm->dmsnes && snes->dmAuto && !dm->dmsnes) {
4715:       DMCopyDMSNES(snes->dm,dm);
4716:       DMGetDMSNES(snes->dm,&sdm);
4717:       if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
4718:     }
4719:     DMDestroy(&snes->dm);
4720:   }
4721:   snes->dm     = dm;
4722:   snes->dmAuto = PETSC_FALSE;

4724:   SNESGetKSP(snes,&ksp);
4725:   KSPSetDM(ksp,dm);
4726:   KSPSetDMActive(ksp,PETSC_FALSE);
4727:   if (snes->pc) {
4728:     SNESSetDM(snes->pc, snes->dm);
4729:     SNESSetNPCSide(snes,snes->pcside);
4730:   }
4731:   return(0);
4732: }

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

4739:    Not Collective but DM obtained is parallel on SNES

4741:    Input Parameter:
4742: . snes - the preconditioner context

4744:    Output Parameter:
4745: .  dm - the dm

4747:    Level: intermediate

4749: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
4750: @*/
4751: PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
4752: {

4757:   if (!snes->dm) {
4758:     DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
4759:     snes->dmAuto = PETSC_TRUE;
4760:   }
4761:   *dm = snes->dm;
4762:   return(0);
4763: }

4767: /*@
4768:   SNESSetNPC - Sets the nonlinear preconditioner to be used.

4770:   Collective on SNES

4772:   Input Parameters:
4773: + snes - iterative context obtained from SNESCreate()
4774: - pc   - the preconditioner object

4776:   Notes:
4777:   Use SNESGetNPC() to retrieve the preconditioner context (for example,
4778:   to configure it using the API).

4780:   Level: developer

4782: .keywords: SNES, set, precondition
4783: .seealso: SNESGetNPC(), SNESHasNPC()
4784: @*/
4785: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
4786: {

4793:   PetscObjectReference((PetscObject) pc);
4794:   SNESDestroy(&snes->pc);
4795:   snes->pc = pc;
4796:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->pc);
4797:   return(0);
4798: }

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

4805:   Not Collective

4807:   Input Parameter:
4808: . snes - iterative context obtained from SNESCreate()

4810:   Output Parameter:
4811: . pc - preconditioner context

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

4815:   Level: developer

4817: .keywords: SNES, get, preconditioner
4818: .seealso: SNESSetNPC(), SNESHasNPC()
4819: @*/
4820: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
4821: {
4823:   const char     *optionsprefix;

4828:   if (!snes->pc) {
4829:     SNESCreate(PetscObjectComm((PetscObject)snes),&snes->pc);
4830:     PetscObjectIncrementTabLevel((PetscObject)snes->pc,(PetscObject)snes,1);
4831:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->pc);
4832:     SNESGetOptionsPrefix(snes,&optionsprefix);
4833:     SNESSetOptionsPrefix(snes->pc,optionsprefix);
4834:     SNESAppendOptionsPrefix(snes->pc,"npc_");
4835:     SNESSetCountersReset(snes->pc,PETSC_FALSE);
4836:   }
4837:   *pc = snes->pc;
4838:   return(0);
4839: }

4843: /*@
4844:   SNESHasNPC - Returns whether a nonlinear preconditioner exists

4846:   Not Collective

4848:   Input Parameter:
4849: . snes - iterative context obtained from SNESCreate()

4851:   Output Parameter:
4852: . has_npc - whether the SNES has an NPC or not

4854:   Level: developer

4856: .keywords: SNES, has, preconditioner
4857: .seealso: SNESSetNPC(), SNESGetNPC()
4858: @*/
4859: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
4860: {
4863:   *has_npc = (PetscBool) (snes->pc != NULL);
4864:   return(0);
4865: }

4869: /*@
4870:     SNESSetNPCSide - Sets the preconditioning side.

4872:     Logically Collective on SNES

4874:     Input Parameter:
4875: .   snes - iterative context obtained from SNESCreate()

4877:     Output Parameter:
4878: .   side - the preconditioning side, where side is one of
4879: .vb
4880:       PC_LEFT - left preconditioning (default)
4881:       PC_RIGHT - right preconditioning
4882: .ve

4884:     Options Database Keys:
4885: .   -snes_pc_side <right,left>

4887:     Level: intermediate

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

4891: .seealso: SNESGetNPCSide(), KSPSetPCSide()
4892: @*/
4893: PetscErrorCode  SNESSetNPCSide(SNES snes,PCSide side)
4894: {
4898:   snes->pcside = side;
4899:   return(0);
4900: }

4904: /*@
4905:     SNESGetNPCSide - Gets the preconditioning side.

4907:     Not Collective

4909:     Input Parameter:
4910: .   snes - iterative context obtained from SNESCreate()

4912:     Output Parameter:
4913: .   side - the preconditioning side, where side is one of
4914: .vb
4915:       PC_LEFT - left preconditioning (default)
4916:       PC_RIGHT - right preconditioning
4917: .ve

4919:     Level: intermediate

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

4923: .seealso: SNESSetNPCSide(), KSPGetPCSide()
4924: @*/
4925: PetscErrorCode  SNESGetNPCSide(SNES snes,PCSide *side)
4926: {
4930:   *side = snes->pcside;
4931:   return(0);
4932: }

4936: /*@
4937:   SNESSetLineSearch - Sets the linesearch on the SNES instance.

4939:   Collective on SNES

4941:   Input Parameters:
4942: + snes - iterative context obtained from SNESCreate()
4943: - linesearch   - the linesearch object

4945:   Notes:
4946:   Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
4947:   to configure it using the API).

4949:   Level: developer

4951: .keywords: SNES, set, linesearch
4952: .seealso: SNESGetLineSearch()
4953: @*/
4954: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
4955: {

4962:   PetscObjectReference((PetscObject) linesearch);
4963:   SNESLineSearchDestroy(&snes->linesearch);

4965:   snes->linesearch = linesearch;

4967:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
4968:   return(0);
4969: }

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

4977:   Not Collective

4979:   Input Parameter:
4980: . snes - iterative context obtained from SNESCreate()

4982:   Output Parameter:
4983: . linesearch - linesearch context

4985:   Level: beginner

4987: .keywords: SNES, get, linesearch
4988: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
4989: @*/
4990: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
4991: {
4993:   const char     *optionsprefix;

4998:   if (!snes->linesearch) {
4999:     SNESGetOptionsPrefix(snes, &optionsprefix);
5000:     SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5001:     SNESLineSearchSetSNES(snes->linesearch, snes);
5002:     SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5003:     PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5004:     PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5005:   }
5006:   *linesearch = snes->linesearch;
5007:   return(0);
5008: }

5010: #if defined(PETSC_HAVE_MATLAB_ENGINE)
5011: #include <mex.h>

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

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

5020:    Collective on SNES

5022:    Input Parameters:
5023: +  snes - the SNES context
5024: -  x - input vector

5026:    Output Parameter:
5027: .  y - function vector, as set by SNESSetFunction()

5029:    Notes:
5030:    SNESComputeFunction() is typically used within nonlinear solvers
5031:    implementations, so most users would not generally call this routine
5032:    themselves.

5034:    Level: developer

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

5038: .seealso: SNESSetFunction(), SNESGetFunction()
5039: */
5040: PetscErrorCode  SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
5041: {
5042:   PetscErrorCode    ierr;
5043:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5044:   int               nlhs  = 1,nrhs = 5;
5045:   mxArray           *plhs[1],*prhs[5];
5046:   long long int     lx = 0,ly = 0,ls = 0;


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

5057:   PetscMemcpy(&ls,&snes,sizeof(snes));
5058:   PetscMemcpy(&lx,&x,sizeof(x));
5059:   PetscMemcpy(&ly,&y,sizeof(x));
5060:   prhs[0] = mxCreateDoubleScalar((double)ls);
5061:   prhs[1] = mxCreateDoubleScalar((double)lx);
5062:   prhs[2] = mxCreateDoubleScalar((double)ly);
5063:   prhs[3] = mxCreateString(sctx->funcname);
5064:   prhs[4] = sctx->ctx;
5065:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");
5066:   mxGetScalar(plhs[0]);
5067:   mxDestroyArray(prhs[0]);
5068:   mxDestroyArray(prhs[1]);
5069:   mxDestroyArray(prhs[2]);
5070:   mxDestroyArray(prhs[3]);
5071:   mxDestroyArray(plhs[0]);
5072:   return(0);
5073: }

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

5082:    Logically Collective on SNES

5084:    Input Parameters:
5085: +  snes - the SNES context
5086: .  r - vector to store function value
5087: -  f - function evaluation routine

5089:    Notes:
5090:    The Newton-like methods typically solve linear systems of the form
5091: $      f'(x) x = -f(x),
5092:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

5094:    Level: beginner

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

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

5100: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5101: */
5102: PetscErrorCode  SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx)
5103: {
5104:   PetscErrorCode    ierr;
5105:   SNESMatlabContext *sctx;

5108:   /* currently sctx is memory bleed */
5109:   PetscNew(&sctx);
5110:   PetscStrallocpy(f,&sctx->funcname);
5111:   /*
5112:      This should work, but it doesn't
5113:   sctx->ctx = ctx;
5114:   mexMakeArrayPersistent(sctx->ctx);
5115:   */
5116:   sctx->ctx = mxDuplicateArray(ctx);
5117:   SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);
5118:   return(0);
5119: }

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

5126:    Collective on SNES

5128:    Input Parameters:
5129: +  snes - the SNES context
5130: .  x - input vector
5131: .  A, B - the matrices
5132: -  ctx - user context

5134:    Level: developer

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

5138: .seealso: SNESSetFunction(), SNESGetFunction()
5139: @*/
5140: PetscErrorCode  SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx)
5141: {
5142:   PetscErrorCode    ierr;
5143:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5144:   int               nlhs  = 2,nrhs = 6;
5145:   mxArray           *plhs[2],*prhs[6];
5146:   long long int     lx = 0,lA = 0,ls = 0, lB = 0;


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

5154:   PetscMemcpy(&ls,&snes,sizeof(snes));
5155:   PetscMemcpy(&lx,&x,sizeof(x));
5156:   PetscMemcpy(&lA,A,sizeof(x));
5157:   PetscMemcpy(&lB,B,sizeof(x));
5158:   prhs[0] = mxCreateDoubleScalar((double)ls);
5159:   prhs[1] = mxCreateDoubleScalar((double)lx);
5160:   prhs[2] = mxCreateDoubleScalar((double)lA);
5161:   prhs[3] = mxCreateDoubleScalar((double)lB);
5162:   prhs[4] = mxCreateString(sctx->funcname);
5163:   prhs[5] = sctx->ctx;
5164:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");
5165:   mxGetScalar(plhs[0]);
5166:   mxDestroyArray(prhs[0]);
5167:   mxDestroyArray(prhs[1]);
5168:   mxDestroyArray(prhs[2]);
5169:   mxDestroyArray(prhs[3]);
5170:   mxDestroyArray(prhs[4]);
5171:   mxDestroyArray(plhs[0]);
5172:   mxDestroyArray(plhs[1]);
5173:   return(0);
5174: }

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

5183:    Logically Collective on SNES

5185:    Input Parameters:
5186: +  snes - the SNES context
5187: .  A,B - Jacobian matrices
5188: .  J - function evaluation routine
5189: -  ctx - user context

5191:    Level: developer

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

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

5197: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J
5198: */
5199: PetscErrorCode  SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx)
5200: {
5201:   PetscErrorCode    ierr;
5202:   SNESMatlabContext *sctx;

5205:   /* currently sctx is memory bleed */
5206:   PetscNew(&sctx);
5207:   PetscStrallocpy(J,&sctx->funcname);
5208:   /*
5209:      This should work, but it doesn't
5210:   sctx->ctx = ctx;
5211:   mexMakeArrayPersistent(sctx->ctx);
5212:   */
5213:   sctx->ctx = mxDuplicateArray(ctx);
5214:   SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);
5215:   return(0);
5216: }

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

5223:    Collective on SNES

5225: .seealso: SNESSetFunction(), SNESGetFunction()
5226: @*/
5227: PetscErrorCode  SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
5228: {
5229:   PetscErrorCode    ierr;
5230:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5231:   int               nlhs  = 1,nrhs = 6;
5232:   mxArray           *plhs[1],*prhs[6];
5233:   long long int     lx = 0,ls = 0;
5234:   Vec               x  = snes->vec_sol;


5239:   PetscMemcpy(&ls,&snes,sizeof(snes));
5240:   PetscMemcpy(&lx,&x,sizeof(x));
5241:   prhs[0] = mxCreateDoubleScalar((double)ls);
5242:   prhs[1] = mxCreateDoubleScalar((double)it);
5243:   prhs[2] = mxCreateDoubleScalar((double)fnorm);
5244:   prhs[3] = mxCreateDoubleScalar((double)lx);
5245:   prhs[4] = mxCreateString(sctx->funcname);
5246:   prhs[5] = sctx->ctx;
5247:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");
5248:   mxGetScalar(plhs[0]);
5249:   mxDestroyArray(prhs[0]);
5250:   mxDestroyArray(prhs[1]);
5251:   mxDestroyArray(prhs[2]);
5252:   mxDestroyArray(prhs[3]);
5253:   mxDestroyArray(prhs[4]);
5254:   mxDestroyArray(plhs[0]);
5255:   return(0);
5256: }

5260: /*
5261:    SNESMonitorSetMatlab - Sets the monitor function from MATLAB

5263:    Level: developer

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

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

5269: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5270: */
5271: PetscErrorCode  SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx)
5272: {
5273:   PetscErrorCode    ierr;
5274:   SNESMatlabContext *sctx;

5277:   /* currently sctx is memory bleed */
5278:   PetscNew(&sctx);
5279:   PetscStrallocpy(f,&sctx->funcname);
5280:   /*
5281:      This should work, but it doesn't
5282:   sctx->ctx = ctx;
5283:   mexMakeArrayPersistent(sctx->ctx);
5284:   */
5285:   sctx->ctx = mxDuplicateArray(ctx);
5286:   SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);
5287:   return(0);
5288: }

5290: #endif