Actual source code: snes.c

petsc-master 2015-04-27
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,NULL,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,NULL,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:     PetscDrawStringBoxed(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_field","Monitor norm of function (split into fields)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
782:   if (flg) {
783:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
784:     SNESMonitorSet(snes,SNESMonitorDefaultField,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
785:   }

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

790:   flg  = PETSC_FALSE;
791:   PetscOptionsBool("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,NULL);
792:   if (flg) {SNESMonitorSet(snes,SNESMonitorSolution,0,0);}
793:   flg  = PETSC_FALSE;
794:   PetscOptionsBool("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,NULL);
795:   if (flg) {SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);}
796:   flg  = PETSC_FALSE;
797:   PetscOptionsBool("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,NULL);
798:   if (flg) {SNESMonitorSet(snes,SNESMonitorResidual,0,0);}
799:   flg  = PETSC_FALSE;
800:   PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);
801:   if (flg) {
802:     PetscObject *objs;

804:     SNESMonitorLGCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&objs);
805:     SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorLGResidualNorm,objs,(PetscErrorCode (*)(void**))SNESMonitorLGDestroy);
806:   }
807:   flg  = PETSC_FALSE;
808:   PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
809:   if (flg) {
810:     PetscViewer ctx;

812:     PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
813:     SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
814:   }

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


821:   PetscOptionsString("-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
822:   if (flg) {
823:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
824:     SNESMonitorSet(snes,SNESMonitorFields,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
825:   }

827:   flg  = PETSC_FALSE;
828:   PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
829:   if (flg) {
830:     void *functx;
831:     SNESGetFunction(snes,NULL,NULL,&functx);
832:     SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
833:     PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
834:   }

836:   flg  = PETSC_FALSE;
837:   PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
838:   if (flg) {
839:     SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
840:   }

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

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

865:   flg  = PETSC_FALSE;
866:   SNESGetNPCSide(snes,&pcside);
867:   PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
868:   if (flg) {SNESSetNPCSide(snes,pcside);}

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

893:   for (i = 0; i < numberofsetfromoptions; i++) {
894:     (*othersetfromoptions[i])(snes);
895:   }

897:   if (snes->ops->setfromoptions) {
898:     (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
899:   }

901:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
902:   PetscObjectProcessOptionsHandlers((PetscObject)snes);
903:   PetscOptionsEnd();

905:   if (!snes->linesearch) {
906:     SNESGetLineSearch(snes, &snes->linesearch);
907:   }
908:   SNESLineSearchSetFromOptions(snes->linesearch);

910:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
911:   KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
912:   KSPSetFromOptions(snes->ksp);

914:   /* if someone has set the SNES NPC type, create it. */
915:   SNESGetOptionsPrefix(snes, &optionsprefix);
916:   PetscOptionsHasName(optionsprefix, "-npc_snes_type", &pcset);
917:   if (pcset && (!snes->pc)) {
918:     SNESGetNPC(snes, &snes->pc);
919:   }
920:   return(0);
921: }

925: /*@C
926:    SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
927:    the nonlinear solvers.

929:    Logically Collective on SNES

931:    Input Parameters:
932: +  snes - the SNES context
933: .  compute - function to compute the context
934: -  destroy - function to destroy the context

936:    Level: intermediate

938:    Notes:
939:    This function is currently not available from Fortran.

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

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

956: /*@
957:    SNESSetApplicationContext - Sets the optional user-defined context for
958:    the nonlinear solvers.

960:    Logically Collective on SNES

962:    Input Parameters:
963: +  snes - the SNES context
964: -  usrP - optional user context

966:    Level: intermediate

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

970: .seealso: SNESGetApplicationContext()
971: @*/
972: PetscErrorCode  SNESSetApplicationContext(SNES snes,void *usrP)
973: {
975:   KSP            ksp;

979:   SNESGetKSP(snes,&ksp);
980:   KSPSetApplicationContext(ksp,usrP);
981:   snes->user = usrP;
982:   return(0);
983: }

987: /*@
988:    SNESGetApplicationContext - Gets the user-defined context for the
989:    nonlinear solvers.

991:    Not Collective

993:    Input Parameter:
994: .  snes - SNES context

996:    Output Parameter:
997: .  usrP - user context

999:    Level: intermediate

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

1003: .seealso: SNESSetApplicationContext()
1004: @*/
1005: PetscErrorCode  SNESGetApplicationContext(SNES snes,void *usrP)
1006: {
1009:   *(void**)usrP = snes->user;
1010:   return(0);
1011: }

1015: /*@
1016:    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1017:    at this time.

1019:    Not Collective

1021:    Input Parameter:
1022: .  snes - SNES context

1024:    Output Parameter:
1025: .  iter - iteration number

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

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

1041:    Level: intermediate

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

1045: .seealso:   SNESGetLinearSolveIterations()
1046: @*/
1047: PetscErrorCode  SNESGetIterationNumber(SNES snes,PetscInt *iter)
1048: {
1052:   *iter = snes->iter;
1053:   return(0);
1054: }

1058: /*@
1059:    SNESSetIterationNumber - Sets the current iteration number.

1061:    Not Collective

1063:    Input Parameter:
1064: .  snes - SNES context
1065: .  iter - iteration number

1067:    Level: developer

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

1071: .seealso:   SNESGetLinearSolveIterations()
1072: @*/
1073: PetscErrorCode  SNESSetIterationNumber(SNES snes,PetscInt iter)
1074: {

1079:   PetscObjectSAWsTakeAccess((PetscObject)snes);
1080:   snes->iter = iter;
1081:   PetscObjectSAWsGrantAccess((PetscObject)snes);
1082:   return(0);
1083: }

1087: /*@
1088:    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1089:    attempted by the nonlinear solver.

1091:    Not Collective

1093:    Input Parameter:
1094: .  snes - SNES context

1096:    Output Parameter:
1097: .  nfails - number of unsuccessful steps attempted

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

1102:    Level: intermediate

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

1106: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1107:           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1108: @*/
1109: PetscErrorCode  SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1110: {
1114:   *nfails = snes->numFailures;
1115:   return(0);
1116: }

1120: /*@
1121:    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1122:    attempted by the nonlinear solver before it gives up.

1124:    Not Collective

1126:    Input Parameters:
1127: +  snes     - SNES context
1128: -  maxFails - maximum of unsuccessful steps

1130:    Level: intermediate

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

1134: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1135:           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1136: @*/
1137: PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1138: {
1141:   snes->maxFailures = maxFails;
1142:   return(0);
1143: }

1147: /*@
1148:    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1149:    attempted by the nonlinear solver before it gives up.

1151:    Not Collective

1153:    Input Parameter:
1154: .  snes     - SNES context

1156:    Output Parameter:
1157: .  maxFails - maximum of unsuccessful steps

1159:    Level: intermediate

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

1163: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1164:           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()

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

1178: /*@
1179:    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1180:      done by SNES.

1182:    Not Collective

1184:    Input Parameter:
1185: .  snes     - SNES context

1187:    Output Parameter:
1188: .  nfuncs - number of evaluations

1190:    Level: intermediate

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

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

1196: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1197: @*/
1198: PetscErrorCode  SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1199: {
1203:   *nfuncs = snes->nfuncs;
1204:   return(0);
1205: }

1209: /*@
1210:    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1211:    linear solvers.

1213:    Not Collective

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

1218:    Output Parameter:
1219: .  nfails - number of failed solves

1221:    Level: intermediate

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

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

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

1231: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1232: @*/
1233: PetscErrorCode  SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1234: {
1238:   *nfails = snes->numLinearSolveFailures;
1239:   return(0);
1240: }

1244: /*@
1245:    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1246:    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE

1248:    Logically Collective on SNES

1250:    Input Parameters:
1251: +  snes     - SNES context
1252: -  maxFails - maximum allowed linear solve failures

1254:    Level: intermediate

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

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

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

1263: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1264: @*/
1265: PetscErrorCode  SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1266: {
1270:   snes->maxLinearSolveFailures = maxFails;
1271:   return(0);
1272: }

1276: /*@
1277:    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1278:      are allowed before SNES terminates

1280:    Not Collective

1282:    Input Parameter:
1283: .  snes     - SNES context

1285:    Output Parameter:
1286: .  maxFails - maximum of unsuccessful solves allowed

1288:    Level: intermediate

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

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

1294: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1295: @*/
1296: PetscErrorCode  SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1297: {
1301:   *maxFails = snes->maxLinearSolveFailures;
1302:   return(0);
1303: }

1307: /*@
1308:    SNESGetLinearSolveIterations - Gets the total number of linear iterations
1309:    used by the nonlinear solver.

1311:    Not Collective

1313:    Input Parameter:
1314: .  snes - SNES context

1316:    Output Parameter:
1317: .  lits - number of linear iterations

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

1322:    Level: intermediate

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

1326: .seealso:  SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1327: @*/
1328: PetscErrorCode  SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1329: {
1333:   *lits = snes->linear_its;
1334:   return(0);
1335: }

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

1343:    Logically Collective on SNES

1345:    Input Parameter:
1346: +  snes - SNES context
1347: -  reset - whether to reset the counters or not

1349:    Notes:
1350:    This defaults to PETSC_TRUE

1352:    Level: developer

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

1356: .seealso:  SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC()
1357: @*/
1358: PetscErrorCode  SNESSetCountersReset(SNES snes,PetscBool reset)
1359: {
1363:   snes->counters_reset = reset;
1364:   return(0);
1365: }


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

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

1375:    Input Parameters:
1376: +  snes - the SNES context
1377: -  ksp - the KSP context

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

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

1386:    Level: developer

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

1390: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1391: @*/
1392: PetscErrorCode  SNESSetKSP(SNES snes,KSP ksp)
1393: {

1400:   PetscObjectReference((PetscObject)ksp);
1401:   if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1402:   snes->ksp = ksp;
1403:   return(0);
1404: }

1406: /* -----------------------------------------------------------*/
1409: /*@
1410:    SNESCreate - Creates a nonlinear solver context.

1412:    Collective on MPI_Comm

1414:    Input Parameters:
1415: .  comm - MPI communicator

1417:    Output Parameter:
1418: .  outsnes - the new SNES context

1420:    Options Database Keys:
1421: +   -snes_mf - Activates default matrix-free Jacobian-vector products,
1422:                and no preconditioning matrix
1423: .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
1424:                products, and a user-provided preconditioning matrix
1425:                as set by SNESSetJacobian()
1426: -   -snes_fd - Uses (slow!) finite differences to compute Jacobian

1428:    Level: beginner

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

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

1434: @*/
1435: PetscErrorCode  SNESCreate(MPI_Comm comm,SNES *outsnes)
1436: {
1438:   SNES           snes;
1439:   SNESKSPEW      *kctx;

1443:   *outsnes = NULL;
1444:   SNESInitializePackage();

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

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

1501:   snes->mf          = PETSC_FALSE;
1502:   snes->mf_operator = PETSC_FALSE;
1503:   snes->mf_version  = 1;

1505:   snes->numLinearSolveFailures = 0;
1506:   snes->maxLinearSolveFailures = 1;

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

1511:   snes->kspconvctx  = (void*)kctx;
1512:   kctx->version     = 2;
1513:   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1514:                              this was too large for some test cases */
1515:   kctx->rtol_last   = 0.0;
1516:   kctx->rtol_max    = .9;
1517:   kctx->gamma       = 1.0;
1518:   kctx->alpha       = .5*(1.0 + PetscSqrtReal(5.0));
1519:   kctx->alpha2      = kctx->alpha;
1520:   kctx->threshold   = .1;
1521:   kctx->lresid_last = 0.0;
1522:   kctx->norm_last   = 0.0;

1524:   *outsnes = snes;
1525:   return(0);
1526: }

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

1531:      Synopsis:
1532:      #include "petscsnes.h"
1533:      PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);

1535:      Input Parameters:
1536: +     snes - the SNES context
1537: .     x    - state at which to evaluate residual
1538: -     ctx     - optional user-defined function context, passed in with SNESSetFunction()

1540:      Output Parameter:
1541: .     f  - vector to put residual (function value)

1543:    Level: intermediate

1545: .seealso:   SNESSetFunction(), SNESGetFunction()
1546: M*/

1550: /*@C
1551:    SNESSetFunction - Sets the function evaluation routine and function
1552:    vector for use by the SNES routines in solving systems of nonlinear
1553:    equations.

1555:    Logically Collective on SNES

1557:    Input Parameters:
1558: +  snes - the SNES context
1559: .  r - vector to store function value
1560: .  f - function evaluation routine; see SNESFunction for calling sequence details
1561: -  ctx - [optional] user-defined context for private data for the
1562:          function evaluation routine (may be NULL)

1564:    Notes:
1565:    The Newton-like methods typically solve linear systems of the form
1566: $      f'(x) x = -f(x),
1567:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

1569:    Level: beginner

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

1573: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1574: @*/
1575: PetscErrorCode  SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1576: {
1578:   DM             dm;

1582:   if (r) {
1585:     PetscObjectReference((PetscObject)r);
1586:     VecDestroy(&snes->vec_func);

1588:     snes->vec_func = r;
1589:   }
1590:   SNESGetDM(snes,&dm);
1591:   DMSNESSetFunction(dm,f,ctx);
1592:   return(0);
1593: }


1598: /*@C
1599:    SNESSetInitialFunction - Sets the function vector to be used as the
1600:    function norm at the initialization of the method.  In some
1601:    instances, the user has precomputed the function before calling
1602:    SNESSolve.  This function allows one to avoid a redundant call
1603:    to SNESComputeFunction in that case.

1605:    Logically Collective on SNES

1607:    Input Parameters:
1608: +  snes - the SNES context
1609: -  f - vector to store function value

1611:    Notes:
1612:    This should not be modified during the solution procedure.

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

1616:    Level: developer

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

1620: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1621: @*/
1622: PetscErrorCode  SNESSetInitialFunction(SNES snes, Vec f)
1623: {
1625:   Vec            vec_func;

1631:   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1632:     snes->vec_func_init_set = PETSC_FALSE;
1633:     return(0);
1634:   }
1635:   SNESGetFunction(snes,&vec_func,NULL,NULL);
1636:   VecCopy(f, vec_func);

1638:   snes->vec_func_init_set = PETSC_TRUE;
1639:   return(0);
1640: }

1644: /*@
1645:    SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1646:    of the SNES method.

1648:    Logically Collective on SNES

1650:    Input Parameters:
1651: +  snes - the SNES context
1652: -  normschedule - the frequency of norm computation

1654:    Options Database Key:
1655: .  -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly>

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

1666:    Level: developer

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

1670: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1671: @*/
1672: PetscErrorCode  SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1673: {
1676:   snes->normschedule = normschedule;
1677:   return(0);
1678: }


1683: /*@
1684:    SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1685:    of the SNES method.

1687:    Logically Collective on SNES

1689:    Input Parameters:
1690: +  snes - the SNES context
1691: -  normschedule - the type of the norm used

1693:    Level: advanced

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

1697: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1698: @*/
1699: PetscErrorCode  SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1700: {
1703:   *normschedule = snes->normschedule;
1704:   return(0);
1705: }


1710: /*@C
1711:    SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
1712:    of the SNES method.

1714:    Logically Collective on SNES

1716:    Input Parameters:
1717: +  snes - the SNES context
1718: -  normschedule - the frequency of norm computation

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

1729:    Level: developer

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

1733: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1734: @*/
1735: PetscErrorCode  SNESSetFunctionType(SNES snes, SNESFunctionType type)
1736: {
1739:   snes->functype = type;
1740:   return(0);
1741: }


1746: /*@C
1747:    SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
1748:    of the SNES method.

1750:    Logically Collective on SNES

1752:    Input Parameters:
1753: +  snes - the SNES context
1754: -  normschedule - the type of the norm used

1756:    Level: advanced

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

1760: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1761: @*/
1762: PetscErrorCode  SNESGetFunctionType(SNES snes, SNESFunctionType *type)
1763: {
1766:   *type = snes->functype;
1767:   return(0);
1768: }

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

1773:      Synopsis:
1774:      #include <petscsnes.h>
1775: $    SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);

1777: +  X   - solution vector
1778: .  B   - RHS vector
1779: -  ctx - optional user-defined Gauss-Seidel context

1781:    Level: intermediate

1783: .seealso:   SNESSetNGS(), SNESGetNGS()
1784: M*/

1788: /*@C
1789:    SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
1790:    use with composed nonlinear solvers.

1792:    Input Parameters:
1793: +  snes   - the SNES context
1794: .  f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
1795: -  ctx    - [optional] user-defined context for private data for the
1796:             smoother evaluation routine (may be NULL)

1798:    Notes:
1799:    The NGS routines are used by the composed nonlinear solver to generate
1800:     a problem appropriate update to the solution, particularly FAS.

1802:    Level: intermediate

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

1806: .seealso: SNESGetFunction(), SNESComputeNGS()
1807: @*/
1808: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1809: {
1811:   DM             dm;

1815:   SNESGetDM(snes,&dm);
1816:   DMSNESSetNGS(dm,f,ctx);
1817:   return(0);
1818: }

1822: PETSC_EXTERN PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
1823: {
1825:   DM             dm;
1826:   DMSNES         sdm;

1829:   SNESGetDM(snes,&dm);
1830:   DMGetDMSNES(dm,&sdm);
1831:   /*  A(x)*x - b(x) */
1832:   if (sdm->ops->computepfunction) {
1833:     (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
1834:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");

1836:   if (sdm->ops->computepjacobian) {
1837:     (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
1838:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard matrix.");
1839:   VecScale(f,-1.0);
1840:   MatMultAdd(snes->jacobian,x,f,f);
1841:   return(0);
1842: }

1846: PETSC_EXTERN PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
1847: {
1849:   /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
1850:   return(0);
1851: }

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

1858:    Logically Collective on SNES

1860:    Input Parameters:
1861: +  snes - the SNES context
1862: .  r - vector to store function value
1863: .  b - function evaluation routine
1864: .  Amat - matrix with which A(x) x - b(x) is to be computed
1865: .  Pmat - matrix from which preconditioner is computed (usually the same as Amat)
1866: .  J  - function to compute matrix value, see SNESJacobianFunction for details on its calling sequence
1867: -  ctx - [optional] user-defined context for private data for the
1868:          function evaluation routine (may be NULL)

1870:    Notes:
1871:     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
1872:     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.

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

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

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

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

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

1888:    Level: intermediate

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

1892: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
1893: @*/
1894: 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)
1895: {
1897:   DM             dm;

1901:   SNESGetDM(snes, &dm);
1902:   DMSNESSetPicard(dm,b,J,ctx);
1903:   SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
1904:   SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
1905:   return(0);
1906: }

1910: /*@C
1911:    SNESGetPicard - Returns the context for the Picard iteration

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

1915:    Input Parameter:
1916: .  snes - the SNES context

1918:    Output Parameter:
1919: +  r - the function (or NULL)
1920: .  f - the function (or NULL); see SNESFunction for calling sequence details
1921: .  Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
1922: .  Pmat  - the matrix from which the preconditioner will be constructed (or NULL)
1923: .  J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
1924: -  ctx - the function context (or NULL)

1926:    Level: advanced

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

1930: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
1931: @*/
1932: 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)
1933: {
1935:   DM             dm;

1939:   SNESGetFunction(snes,r,NULL,NULL);
1940:   SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
1941:   SNESGetDM(snes,&dm);
1942:   DMSNESGetPicard(dm,f,J,ctx);
1943:   return(0);
1944: }

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

1951:    Logically Collective on SNES

1953:    Input Parameters:
1954: +  snes - the SNES context
1955: .  func - function evaluation routine
1956: -  ctx - [optional] user-defined context for private data for the
1957:          function evaluation routine (may be NULL)

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

1962: .  f - function vector
1963: -  ctx - optional user-defined function context

1965:    Level: intermediate

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

1969: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
1970: @*/
1971: PetscErrorCode  SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
1972: {
1975:   if (func) snes->ops->computeinitialguess = func;
1976:   if (ctx)  snes->initialguessP            = ctx;
1977:   return(0);
1978: }

1980: /* --------------------------------------------------------------- */
1983: /*@C
1984:    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
1985:    it assumes a zero right hand side.

1987:    Logically Collective on SNES

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

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

1995:    Level: intermediate

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

1999: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
2000: @*/
2001: PetscErrorCode  SNESGetRhs(SNES snes,Vec *rhs)
2002: {
2006:   *rhs = snes->vec_rhs;
2007:   return(0);
2008: }

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

2015:    Collective on SNES

2017:    Input Parameters:
2018: +  snes - the SNES context
2019: -  x - input vector

2021:    Output Parameter:
2022: .  y - function vector, as set by SNESSetFunction()

2024:    Notes:
2025:    SNESComputeFunction() is typically used within nonlinear solvers
2026:    implementations, so most users would not generally call this routine
2027:    themselves.

2029:    Level: developer

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

2033: .seealso: SNESSetFunction(), SNESGetFunction()
2034: @*/
2035: PetscErrorCode  SNESComputeFunction(SNES snes,Vec x,Vec y)
2036: {
2038:   DM             dm;
2039:   DMSNES         sdm;

2047:   VecValidValues(x,2,PETSC_TRUE);

2049:   SNESGetDM(snes,&dm);
2050:   DMGetDMSNES(dm,&sdm);
2051:   if (sdm->ops->computefunction) {
2052:     PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2053:     VecLockPush(x);
2054:     PetscStackPush("SNES user function");
2055:     (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2056:     PetscStackPop;
2057:     VecLockPop(x);
2058:     PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2059:   } else if (snes->vec_rhs) {
2060:     MatMult(snes->jacobian, x, y);
2061:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2062:   if (snes->vec_rhs) {
2063:     VecAXPY(y,-1.0,snes->vec_rhs);
2064:   }
2065:   snes->nfuncs++;
2066:   return(0);
2067: }

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

2074:    Collective on SNES

2076:    Input Parameters:
2077: +  snes - the SNES context
2078: .  x - input vector
2079: -  b - rhs vector

2081:    Output Parameter:
2082: .  x - new solution vector

2084:    Notes:
2085:    SNESComputeNGS() is typically used within composed nonlinear solver
2086:    implementations, so most users would not generally call this routine
2087:    themselves.

2089:    Level: developer

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

2093: .seealso: SNESSetNGS(), SNESComputeFunction()
2094: @*/
2095: PetscErrorCode  SNESComputeNGS(SNES snes,Vec b,Vec x)
2096: {
2098:   DM             dm;
2099:   DMSNES         sdm;

2107:   if (b) {VecValidValues(b,2,PETSC_TRUE);}
2108:   PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2109:   SNESGetDM(snes,&dm);
2110:   DMGetDMSNES(dm,&sdm);
2111:   if (sdm->ops->computegs) {
2112:     if (b) {VecLockPush(b);}
2113:     PetscStackPush("SNES user NGS");
2114:     (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2115:     PetscStackPop;
2116:     if (b) {VecLockPop(b);}
2117:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2118:   PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2119:   return(0);
2120: }

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

2127:    Collective on SNES and Mat

2129:    Input Parameters:
2130: +  snes - the SNES context
2131: -  x - input vector

2133:    Output Parameters:
2134: +  A - Jacobian matrix
2135: -  B - optional preconditioning matrix

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


2153:    Notes:
2154:    Most users should not need to explicitly call this routine, as it
2155:    is used internally within the nonlinear solvers.

2157:    Level: developer

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

2161: .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2162: @*/
2163: PetscErrorCode  SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2164: {
2166:   PetscBool      flag;
2167:   DM             dm;
2168:   DMSNES         sdm;
2169:   KSP            ksp;

2175:   VecValidValues(X,2,PETSC_TRUE);
2176:   SNESGetDM(snes,&dm);
2177:   DMGetDMSNES(dm,&sdm);

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

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

2183:   if (snes->lagjacobian == -2) {
2184:     snes->lagjacobian = -1;

2186:     PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2187:   } else if (snes->lagjacobian == -1) {
2188:     PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2189:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2190:     if (flag) {
2191:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2192:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2193:     }
2194:     return(0);
2195:   } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2196:     PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2197:     PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2198:     if (flag) {
2199:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2200:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2201:     }
2202:     return(0);
2203:   }
2204:   if (snes->pc && snes->pcside == PC_LEFT) {
2205:       MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2206:       MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2207:       return(0);
2208:   }

2210:   PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2211:   VecLockPush(X);
2212:   PetscStackPush("SNES user Jacobian function");
2213:   (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2214:   PetscStackPop;
2215:   VecLockPop(X);
2216:   PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);

2218:   /* the next line ensures that snes->ksp exists */
2219:   SNESGetKSP(snes,&ksp);
2220:   if (snes->lagpreconditioner == -2) {
2221:     PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2222:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2223:     snes->lagpreconditioner = -1;
2224:   } else if (snes->lagpreconditioner == -1) {
2225:     PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2226:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2227:   } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2228:     PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2229:     KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2230:   } else {
2231:     PetscInfo(snes,"Rebuilding preconditioner\n");
2232:     KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2233:   }

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

2308:       MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2309:       MatColoringCreate(Bfd,&coloring);
2310:       MatColoringSetType(coloring,MATCOLORINGSL);
2311:       MatColoringSetFromOptions(coloring);
2312:       MatColoringApply(coloring,&iscoloring);
2313:       MatColoringDestroy(&coloring);
2314:       MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2315:       MatFDColoringSetFromOptions(matfdcoloring);
2316:       MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2317:       ISColoringDestroy(&iscoloring);

2319:       /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2320:       SNESGetFunction(snes,NULL,&func,&funcctx);
2321:       MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2322:       PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2323:       PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2324:       MatFDColoringSetFromOptions(matfdcoloring);
2325:       MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2326:       MatFDColoringDestroy(&matfdcoloring);

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

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

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

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

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

2413:    Level: intermediate

2415: .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2416: M*/

2420: /*@C
2421:    SNESSetJacobian - Sets the function to compute Jacobian as well as the
2422:    location to store the matrix.

2424:    Logically Collective on SNES and Mat

2426:    Input Parameters:
2427: +  snes - the SNES context
2428: .  Amat - the matrix that defines the (approximate) Jacobian
2429: .  Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2430: .  J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
2431: -  ctx - [optional] user-defined context for private data for the
2432:          Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)

2434:    Notes:
2435:    If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2436:    each matrix.

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

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

2444:    Level: beginner

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

2448: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J, 
2449:           SNESSetPicard(), SNESJacobianFunction
2450: @*/
2451: PetscErrorCode  SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2452: {
2454:   DM             dm;

2462:   SNESGetDM(snes,&dm);
2463:   DMSNESSetJacobian(dm,J,ctx);
2464:   if (Amat) {
2465:     PetscObjectReference((PetscObject)Amat);
2466:     MatDestroy(&snes->jacobian);

2468:     snes->jacobian = Amat;
2469:   }
2470:   if (Pmat) {
2471:     PetscObjectReference((PetscObject)Pmat);
2472:     MatDestroy(&snes->jacobian_pre);

2474:     snes->jacobian_pre = Pmat;
2475:   }
2476:   return(0);
2477: }

2481: /*@C
2482:    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2483:    provided context for evaluating the Jacobian.

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

2487:    Input Parameter:
2488: .  snes - the nonlinear solver context

2490:    Output Parameters:
2491: +  Amat - location to stash (approximate) Jacobian matrix (or NULL)
2492: .  Pmat - location to stash matrix used to compute the preconditioner (or NULL)
2493: .  J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
2494: -  ctx - location to stash Jacobian ctx (or NULL)

2496:    Level: advanced

2498: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
2499: @*/
2500: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2501: {
2503:   DM             dm;
2504:   DMSNES         sdm;

2508:   if (Amat) *Amat = snes->jacobian;
2509:   if (Pmat) *Pmat = snes->jacobian_pre;
2510:   SNESGetDM(snes,&dm);
2511:   DMGetDMSNES(dm,&sdm);
2512:   if (J) *J = sdm->ops->computejacobian;
2513:   if (ctx) *ctx = sdm->jacobianctx;
2514:   return(0);
2515: }

2519: /*@
2520:    SNESSetUp - Sets up the internal data structures for the later use
2521:    of a nonlinear solver.

2523:    Collective on SNES

2525:    Input Parameters:
2526: .  snes - the SNES context

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

2535:    Level: advanced

2537: .keywords: SNES, nonlinear, setup

2539: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
2540: @*/
2541: PetscErrorCode  SNESSetUp(SNES snes)
2542: {
2544:   DM             dm;
2545:   DMSNES         sdm;
2546:   SNESLineSearch linesearch, pclinesearch;
2547:   void           *lsprectx,*lspostctx;
2548:   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
2549:   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
2550:   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2551:   Vec            f,fpc;
2552:   void           *funcctx;
2553:   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2554:   void           *jacctx,*appctx;
2555:   Mat            j,jpre;

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

2561:   if (!((PetscObject)snes)->type_name) {
2562:     SNESSetType(snes,SNESNEWTONLS);
2563:   }

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

2567:   SNESGetDM(snes,&dm);
2568:   DMGetDMSNES(dm,&sdm);
2569:   if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
2570:   if (!sdm->ops->computejacobian) {
2571:     DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
2572:   }
2573:   if (!snes->vec_func) {
2574:     DMCreateGlobalVector(dm,&snes->vec_func);
2575:   }

2577:   if (!snes->ksp) {
2578:     SNESGetKSP(snes, &snes->ksp);
2579:   }

2581:   if (!snes->linesearch) {
2582:     SNESGetLineSearch(snes, &snes->linesearch);
2583:   }
2584:   SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);

2586:   if (snes->pc && (snes->pcside == PC_LEFT)) {
2587:     snes->mf          = PETSC_TRUE;
2588:     snes->mf_operator = PETSC_FALSE;
2589:   }

2591:   if (snes->pc) {
2592:     /* copy the DM over */
2593:     SNESGetDM(snes,&dm);
2594:     SNESSetDM(snes->pc,dm);

2596:     SNESGetFunction(snes,&f,&func,&funcctx);
2597:     VecDuplicate(f,&fpc);
2598:     SNESSetFunction(snes->pc,fpc,func,funcctx);
2599:     SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
2600:     SNESSetJacobian(snes->pc,j,jpre,jac,jacctx);
2601:     SNESGetApplicationContext(snes,&appctx);
2602:     SNESSetApplicationContext(snes->pc,appctx);
2603:     VecDestroy(&fpc);

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

2608:     /* default to 1 iteration */
2609:     SNESSetTolerances(snes->pc,0.0,0.0,0.0,1,snes->pc->max_funcs);
2610:     if (snes->pcside==PC_RIGHT) {
2611:       SNESSetNormSchedule(snes->pc,SNES_NORM_FINAL_ONLY);
2612:     } else {
2613:       SNESSetNormSchedule(snes->pc,SNES_NORM_NONE);
2614:     }
2615:     SNESSetFromOptions(snes->pc);

2617:     /* copy the line search context over */
2618:     SNESGetLineSearch(snes,&linesearch);
2619:     SNESGetLineSearch(snes->pc,&pclinesearch);
2620:     SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
2621:     SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
2622:     SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
2623:     SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
2624:     PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
2625:   }
2626:   if (snes->mf) {
2627:     SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
2628:   }
2629:   if (snes->ops->usercompute && !snes->user) {
2630:     (*snes->ops->usercompute)(snes,(void**)&snes->user);
2631:   }

2633:   snes->jac_iter = 0;
2634:   snes->pre_iter = 0;

2636:   if (snes->ops->setup) {
2637:     (*snes->ops->setup)(snes);
2638:   }

2640:   if (snes->pc && (snes->pcside == PC_LEFT)) {
2641:     if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2642:       SNESGetLineSearch(snes,&linesearch);
2643:       SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
2644:     }
2645:   }

2647:   snes->setupcalled = PETSC_TRUE;
2648:   return(0);
2649: }

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

2656:    Collective on SNES

2658:    Input Parameter:
2659: .  snes - iterative context obtained from SNESCreate()

2661:    Level: intermediate

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

2665: .keywords: SNES, destroy

2667: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
2668: @*/
2669: PetscErrorCode  SNESReset(SNES snes)
2670: {

2675:   if (snes->ops->userdestroy && snes->user) {
2676:     (*snes->ops->userdestroy)((void**)&snes->user);
2677:     snes->user = NULL;
2678:   }
2679:   if (snes->pc) {
2680:     SNESReset(snes->pc);
2681:   }

2683:   if (snes->ops->reset) {
2684:     (*snes->ops->reset)(snes);
2685:   }
2686:   if (snes->ksp) {
2687:     KSPReset(snes->ksp);
2688:   }

2690:   if (snes->linesearch) {
2691:     SNESLineSearchReset(snes->linesearch);
2692:   }

2694:   VecDestroy(&snes->vec_rhs);
2695:   VecDestroy(&snes->vec_sol);
2696:   VecDestroy(&snes->vec_sol_update);
2697:   VecDestroy(&snes->vec_func);
2698:   MatDestroy(&snes->jacobian);
2699:   MatDestroy(&snes->jacobian_pre);
2700:   VecDestroyVecs(snes->nwork,&snes->work);
2701:   VecDestroyVecs(snes->nvwork,&snes->vwork);

2703:   snes->nwork       = snes->nvwork = 0;
2704:   snes->setupcalled = PETSC_FALSE;
2705:   return(0);
2706: }

2710: /*@
2711:    SNESDestroy - Destroys the nonlinear solver context that was created
2712:    with SNESCreate().

2714:    Collective on SNES

2716:    Input Parameter:
2717: .  snes - the SNES context

2719:    Level: beginner

2721: .keywords: SNES, nonlinear, destroy

2723: .seealso: SNESCreate(), SNESSolve()
2724: @*/
2725: PetscErrorCode  SNESDestroy(SNES *snes)
2726: {

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

2734:   SNESReset((*snes));
2735:   SNESDestroy(&(*snes)->pc);

2737:   /* if memory was published with SAWs then destroy it */
2738:   PetscObjectSAWsViewOff((PetscObject)*snes);
2739:   if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}

2741:   DMDestroy(&(*snes)->dm);
2742:   KSPDestroy(&(*snes)->ksp);
2743:   SNESLineSearchDestroy(&(*snes)->linesearch);

2745:   PetscFree((*snes)->kspconvctx);
2746:   if ((*snes)->ops->convergeddestroy) {
2747:     (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
2748:   }
2749:   if ((*snes)->conv_malloc) {
2750:     PetscFree((*snes)->conv_hist);
2751:     PetscFree((*snes)->conv_hist_its);
2752:   }
2753:   SNESMonitorCancel((*snes));
2754:   PetscHeaderDestroy(snes);
2755:   return(0);
2756: }

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

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

2765:    Logically Collective on SNES

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

2772:    Options Database Keys:
2773: .    -snes_lag_preconditioner <lag>

2775:    Notes:
2776:    The default is 1
2777:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2778:    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use

2780:    Level: intermediate

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

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

2786: @*/
2787: PetscErrorCode  SNESSetLagPreconditioner(SNES snes,PetscInt lag)
2788: {
2791:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2792:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2794:   snes->lagpreconditioner = lag;
2795:   return(0);
2796: }

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

2803:    Logically Collective on SNES

2805:    Input Parameters:
2806: +  snes - the SNES context
2807: -  steps - the number of refinements to do, defaults to 0

2809:    Options Database Keys:
2810: .    -snes_grid_sequence <steps>

2812:    Level: intermediate

2814:    Notes:
2815:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

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

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

2821: @*/
2822: PetscErrorCode  SNESSetGridSequence(SNES snes,PetscInt steps)
2823: {
2827:   snes->gridsequence = steps;
2828:   return(0);
2829: }

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

2836:    Logically Collective on SNES

2838:    Input Parameter:
2839: .  snes - the SNES context

2841:    Output Parameter:
2842: .  steps - the number of refinements to do, defaults to 0

2844:    Options Database Keys:
2845: .    -snes_grid_sequence <steps>

2847:    Level: intermediate

2849:    Notes:
2850:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

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

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

2856: @*/
2857: PetscErrorCode  SNESGetGridSequence(SNES snes,PetscInt *steps)
2858: {
2861:   *steps = snes->gridsequence;
2862:   return(0);
2863: }

2867: /*@
2868:    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt

2870:    Not Collective

2872:    Input Parameter:
2873: .  snes - the SNES context

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

2879:    Options Database Keys:
2880: .    -snes_lag_preconditioner <lag>

2882:    Notes:
2883:    The default is 1
2884:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

2886:    Level: intermediate

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

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

2892: @*/
2893: PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
2894: {
2897:   *lag = snes->lagpreconditioner;
2898:   return(0);
2899: }

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

2907:    Logically Collective on SNES

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

2914:    Options Database Keys:
2915: .    -snes_lag_jacobian <lag>

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

2923:    Level: intermediate

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

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

2929: @*/
2930: PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
2931: {
2934:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2935:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2937:   snes->lagjacobian = lag;
2938:   return(0);
2939: }

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

2946:    Not Collective

2948:    Input Parameter:
2949: .  snes - the SNES context

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

2955:    Options Database Keys:
2956: .    -snes_lag_jacobian <lag>

2958:    Notes:
2959:    The default is 1
2960:    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

2962:    Level: intermediate

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

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

2968: @*/
2969: PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
2970: {
2973:   *lag = snes->lagjacobian;
2974:   return(0);
2975: }

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

2982:    Logically collective on SNES

2984:    Input Parameter:
2985: +  snes - the SNES context
2986: -   flg - jacobian lagging persists if true

2988:    Options Database Keys:
2989: .    -snes_lag_jacobian_persists <flg>

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

2995:    Level: developer

2997: .keywords: SNES, nonlinear, lag

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

3001: @*/
3002: PetscErrorCode  SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
3003: {
3007:   snes->lagjac_persist = flg;
3008:   return(0);
3009: }

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

3016:    Logically Collective on SNES

3018:    Input Parameter:
3019: +  snes - the SNES context
3020: -   flg - preconditioner lagging persists if true

3022:    Options Database Keys:
3023: .    -snes_lag_jacobian_persists <flg>

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

3029:    Level: developer

3031: .keywords: SNES, nonlinear, lag

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

3035: @*/
3036: PetscErrorCode  SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3037: {
3041:   snes->lagpre_persist = flg;
3042:   return(0);
3043: }

3047: /*@
3048:    SNESSetTolerances - Sets various parameters used in convergence tests.

3050:    Logically Collective on SNES

3052:    Input Parameters:
3053: +  snes - the SNES context
3054: .  abstol - absolute convergence tolerance
3055: .  rtol - relative convergence tolerance
3056: .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
3057: .  maxit - maximum number of iterations
3058: -  maxf - maximum number of function evaluations

3060:    Options Database Keys:
3061: +    -snes_atol <abstol> - Sets abstol
3062: .    -snes_rtol <rtol> - Sets rtol
3063: .    -snes_stol <stol> - Sets stol
3064: .    -snes_max_it <maxit> - Sets maxit
3065: -    -snes_max_funcs <maxf> - Sets maxf

3067:    Notes:
3068:    The default maximum number of iterations is 50.
3069:    The default maximum number of function evaluations is 1000.

3071:    Level: intermediate

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

3075: .seealso: SNESSetTrustRegionTolerance()
3076: @*/
3077: PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3078: {

3087:   if (abstol != PETSC_DEFAULT) {
3088:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3089:     snes->abstol = abstol;
3090:   }
3091:   if (rtol != PETSC_DEFAULT) {
3092:     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);
3093:     snes->rtol = rtol;
3094:   }
3095:   if (stol != PETSC_DEFAULT) {
3096:     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3097:     snes->stol = stol;
3098:   }
3099:   if (maxit != PETSC_DEFAULT) {
3100:     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3101:     snes->max_its = maxit;
3102:   }
3103:   if (maxf != PETSC_DEFAULT) {
3104:     if (maxf < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
3105:     snes->max_funcs = maxf;
3106:   }
3107:   snes->tolerancesset = PETSC_TRUE;
3108:   return(0);
3109: }

3113: /*@
3114:    SNESGetTolerances - Gets various parameters used in convergence tests.

3116:    Not Collective

3118:    Input Parameters:
3119: +  snes - the SNES context
3120: .  atol - absolute convergence tolerance
3121: .  rtol - relative convergence tolerance
3122: .  stol -  convergence tolerance in terms of the norm
3123:            of the change in the solution between steps
3124: .  maxit - maximum number of iterations
3125: -  maxf - maximum number of function evaluations

3127:    Notes:
3128:    The user can specify NULL for any parameter that is not needed.

3130:    Level: intermediate

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

3134: .seealso: SNESSetTolerances()
3135: @*/
3136: PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3137: {
3140:   if (atol)  *atol  = snes->abstol;
3141:   if (rtol)  *rtol  = snes->rtol;
3142:   if (stol)  *stol  = snes->stol;
3143:   if (maxit) *maxit = snes->max_its;
3144:   if (maxf)  *maxf  = snes->max_funcs;
3145:   return(0);
3146: }

3150: /*@
3151:    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.

3153:    Logically Collective on SNES

3155:    Input Parameters:
3156: +  snes - the SNES context
3157: -  tol - tolerance

3159:    Options Database Key:
3160: .  -snes_trtol <tol> - Sets tol

3162:    Level: intermediate

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

3166: .seealso: SNESSetTolerances()
3167: @*/
3168: PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3169: {
3173:   snes->deltatol = tol;
3174:   return(0);
3175: }

3177: /*
3178:    Duplicate the lg monitors for SNES from KSP; for some reason with
3179:    dynamic libraries things don't work under Sun4 if we just use
3180:    macros instead of functions
3181: */
3184: PetscErrorCode  SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,PetscObject *objs)
3185: {

3190:   KSPMonitorLGResidualNorm((KSP)snes,it,norm,objs);
3191:   return(0);
3192: }

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

3201:   KSPMonitorLGResidualNormCreate(host,label,x,y,m,n,draw);
3202:   return(0);
3203: }

3207: PetscErrorCode  SNESMonitorLGDestroy(PetscObject **objs)
3208: {

3212:   KSPMonitorLGResidualNormDestroy(objs);
3213:   return(0);
3214: }

3216: extern PetscErrorCode  SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3219: PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3220: {
3221:   PetscDrawLG      lg;
3222:   PetscErrorCode   ierr;
3223:   PetscReal        x,y,per;
3224:   PetscViewer      v = (PetscViewer)monctx;
3225:   static PetscReal prev; /* should be in the context */
3226:   PetscDraw        draw;

3229:   PetscViewerDrawGetDrawLG(v,0,&lg);
3230:   if (!n) {PetscDrawLGReset(lg);}
3231:   PetscDrawLGGetDraw(lg,&draw);
3232:   PetscDrawSetTitle(draw,"Residual norm");
3233:   x    = (PetscReal)n;
3234:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3235:   else y = -15.0;
3236:   PetscDrawLGAddPoint(lg,&x,&y);
3237:   if (n < 20 || !(n % 5)) {
3238:     PetscDrawLGDraw(lg);
3239:   }

3241:   PetscViewerDrawGetDrawLG(v,1,&lg);
3242:   if (!n) {PetscDrawLGReset(lg);}
3243:   PetscDrawLGGetDraw(lg,&draw);
3244:   PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3245:    SNESMonitorRange_Private(snes,n,&per);
3246:   x    = (PetscReal)n;
3247:   y    = 100.0*per;
3248:   PetscDrawLGAddPoint(lg,&x,&y);
3249:   if (n < 20 || !(n % 5)) {
3250:     PetscDrawLGDraw(lg);
3251:   }

3253:   PetscViewerDrawGetDrawLG(v,2,&lg);
3254:   if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3255:   PetscDrawLGGetDraw(lg,&draw);
3256:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3257:   x    = (PetscReal)n;
3258:   y    = (prev - rnorm)/prev;
3259:   PetscDrawLGAddPoint(lg,&x,&y);
3260:   if (n < 20 || !(n % 5)) {
3261:     PetscDrawLGDraw(lg);
3262:   }

3264:   PetscViewerDrawGetDrawLG(v,3,&lg);
3265:   if (!n) {PetscDrawLGReset(lg);}
3266:   PetscDrawLGGetDraw(lg,&draw);
3267:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3268:   x    = (PetscReal)n;
3269:   y    = (prev - rnorm)/(prev*per);
3270:   if (n > 2) { /*skip initial crazy value */
3271:     PetscDrawLGAddPoint(lg,&x,&y);
3272:   }
3273:   if (n < 20 || !(n % 5)) {
3274:     PetscDrawLGDraw(lg);
3275:   }
3276:   prev = rnorm;
3277:   return(0);
3278: }

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

3285:    Collective on SNES

3287:    Input Parameters:
3288: +  snes - nonlinear solver context obtained from SNESCreate()
3289: .  iter - iteration number
3290: -  rnorm - relative norm of the residual

3292:    Notes:
3293:    This routine is called by the SNES implementations.
3294:    It does not typically need to be called by the user.

3296:    Level: developer

3298: .seealso: SNESMonitorSet()
3299: @*/
3300: PetscErrorCode  SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3301: {
3303:   PetscInt       i,n = snes->numbermonitors;

3306:   VecLockPush(snes->vec_sol);
3307:   for (i=0; i<n; i++) {
3308:     (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3309:   }
3310:   VecLockPop(snes->vec_sol);
3311:   return(0);
3312: }

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

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

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

3323: +    snes - the SNES context
3324: .    its - iteration number
3325: .    norm - 2-norm function value (may be estimated)
3326: -    mctx - [optional] monitoring context

3328:    Level: advanced

3330: .seealso:   SNESMonitorSet(), SNESMonitorGet()
3331: M*/

3335: /*@C
3336:    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3337:    iteration of the nonlinear solver to display the iteration's
3338:    progress.

3340:    Logically Collective on SNES

3342:    Input Parameters:
3343: +  snes - the SNES context
3344: .  f - the monitor function, see SNESMonitorFunction for the calling sequence
3345: .  mctx - [optional] user-defined context for private data for the
3346:           monitor routine (use NULL if no context is desired)
3347: -  monitordestroy - [optional] routine that frees monitor context
3348:           (may be NULL)

3350:    Options Database Keys:
3351: +    -snes_monitor        - sets SNESMonitorDefault()
3352: .    -snes_monitor_lg_residualnorm    - sets line graph monitor,
3353:                             uses SNESMonitorLGCreate()
3354: -    -snes_monitor_cancel - cancels all monitors that have
3355:                             been hardwired into a code by
3356:                             calls to SNESMonitorSet(), but
3357:                             does not cancel those set via
3358:                             the options database.

3360:    Notes:
3361:    Several different monitoring routines may be set by calling
3362:    SNESMonitorSet() multiple times; all will be called in the
3363:    order in which they were set.

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

3367:    Level: intermediate

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

3371: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3372: @*/
3373: PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3374: {
3375:   PetscInt       i;

3380:   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3381:   for (i=0; i<snes->numbermonitors;i++) {
3382:     if (f == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) {
3383:       if (monitordestroy) {
3384:         (*monitordestroy)(&mctx);
3385:       }
3386:       return(0);
3387:     }
3388:   }
3389:   snes->monitor[snes->numbermonitors]          = f;
3390:   snes->monitordestroy[snes->numbermonitors]   = monitordestroy;
3391:   snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3392:   return(0);
3393: }

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

3400:    Logically Collective on SNES

3402:    Input Parameters:
3403: .  snes - the SNES context

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

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

3413:    Level: intermediate

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

3417: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3418: @*/
3419: PetscErrorCode  SNESMonitorCancel(SNES snes)
3420: {
3422:   PetscInt       i;

3426:   for (i=0; i<snes->numbermonitors; i++) {
3427:     if (snes->monitordestroy[i]) {
3428:       (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3429:     }
3430:   }
3431:   snes->numbermonitors = 0;
3432:   return(0);
3433: }

3435: /*MC
3436:     SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver

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

3442: +    snes - the SNES context
3443: .    it - current iteration (0 is the first and is before any Newton step)
3444: .    cctx - [optional] convergence context
3445: .    reason - reason for convergence/divergence
3446: .    xnorm - 2-norm of current iterate
3447: .    gnorm - 2-norm of current step
3448: -    f - 2-norm of function

3450:    Level: intermediate

3452: .seealso:   SNESSetConvergenceTest(), SNESGetConvergenceTest()
3453: M*/

3457: /*@C
3458:    SNESSetConvergenceTest - Sets the function that is to be used
3459:    to test for convergence of the nonlinear iterative solution.

3461:    Logically Collective on SNES

3463:    Input Parameters:
3464: +  snes - the SNES context
3465: .  SNESConvergenceTestFunction - routine to test for convergence
3466: .  cctx - [optional] context for private data for the convergence routine  (may be NULL)
3467: -  destroy - [optional] destructor for the context (may be NULL; NULL_FUNCTION in Fortran)

3469:    Level: advanced

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

3473: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3474: @*/
3475: PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3476: {

3481:   if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
3482:   if (snes->ops->convergeddestroy) {
3483:     (*snes->ops->convergeddestroy)(snes->cnvP);
3484:   }
3485:   snes->ops->converged        = SNESConvergenceTestFunction;
3486:   snes->ops->convergeddestroy = destroy;
3487:   snes->cnvP                  = cctx;
3488:   return(0);
3489: }

3493: /*@
3494:    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.

3496:    Not Collective

3498:    Input Parameter:
3499: .  snes - the SNES context

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

3505:    Level: intermediate

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

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

3511: .seealso: SNESSetConvergenceTest(), SNESConvergedReason
3512: @*/
3513: PetscErrorCode  SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
3514: {
3518:   *reason = snes->reason;
3519:   return(0);
3520: }

3524: /*@
3525:    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.

3527:    Logically Collective on SNES

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

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

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

3545:    Level: intermediate

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

3549: .seealso: SNESGetConvergenceHistory()

3551: @*/
3552: PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
3553: {

3560:   if (!a) {
3561:     if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
3562:     PetscCalloc1(na,&a);
3563:     PetscCalloc1(na,&its);

3565:     snes->conv_malloc = PETSC_TRUE;
3566:   }
3567:   snes->conv_hist       = a;
3568:   snes->conv_hist_its   = its;
3569:   snes->conv_hist_max   = na;
3570:   snes->conv_hist_len   = 0;
3571:   snes->conv_hist_reset = reset;
3572:   return(0);
3573: }

3575: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3576: #include <engine.h>   /* MATLAB include file */
3577: #include <mex.h>      /* MATLAB include file */

3581: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
3582: {
3583:   mxArray   *mat;
3584:   PetscInt  i;
3585:   PetscReal *ar;

3588:   mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
3589:   ar  = (PetscReal*) mxGetData(mat);
3590:   for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
3591:   PetscFunctionReturn(mat);
3592: }
3593: #endif

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

3600:    Not Collective

3602:    Input Parameter:
3603: .  snes - iterative context obtained from SNESCreate()

3605:    Output Parameters:
3606: .  a   - array to hold history
3607: .  its - integer array holds the number of linear iterations (or
3608:          negative if not converged) for each solve.
3609: -  na  - size of a and its

3611:    Notes:
3612:     The calling sequence for this routine in Fortran is
3613: $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)

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

3619:    Level: intermediate

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

3623: .seealso: SNESSetConvergencHistory()

3625: @*/
3626: PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
3627: {
3630:   if (a)   *a   = snes->conv_hist;
3631:   if (its) *its = snes->conv_hist_its;
3632:   if (na)  *na  = snes->conv_hist_len;
3633:   return(0);
3634: }

3638: /*@C
3639:   SNESSetUpdate - Sets the general-purpose update function called
3640:   at the beginning of every iteration of the nonlinear solve. Specifically
3641:   it is called just before the Jacobian is "evaluated".

3643:   Logically Collective on SNES

3645:   Input Parameters:
3646: . snes - The nonlinear solver context
3647: . func - The function

3649:   Calling sequence of func:
3650: . func (SNES snes, PetscInt step);

3652: . step - The current step of the iteration

3654:   Level: advanced

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

3659: .keywords: SNES, update

3661: .seealso SNESSetJacobian(), SNESSolve()
3662: @*/
3663: PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
3664: {
3667:   snes->ops->update = func;
3668:   return(0);
3669: }

3673: /*
3674:    SNESScaleStep_Private - Scales a step so that its length is less than the
3675:    positive parameter delta.

3677:     Input Parameters:
3678: +   snes - the SNES context
3679: .   y - approximate solution of linear system
3680: .   fnorm - 2-norm of current function
3681: -   delta - trust region size

3683:     Output Parameters:
3684: +   gpnorm - predicted function norm at the new point, assuming local
3685:     linearization.  The value is zero if the step lies within the trust
3686:     region, and exceeds zero otherwise.
3687: -   ynorm - 2-norm of the step

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

3693: .keywords: SNES, nonlinear, scale, step
3694: */
3695: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
3696: {
3697:   PetscReal      nrm;
3698:   PetscScalar    cnorm;


3706:   VecNorm(y,NORM_2,&nrm);
3707:   if (nrm > *delta) {
3708:     nrm     = *delta/nrm;
3709:     *gpnorm = (1.0 - nrm)*(*fnorm);
3710:     cnorm   = nrm;
3711:     VecScale(y,cnorm);
3712:     *ynorm  = *delta;
3713:   } else {
3714:     *gpnorm = 0.0;
3715:     *ynorm  = nrm;
3716:   }
3717:   return(0);
3718: }

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

3725:    Collective on SNES

3727:    Parameter:
3728: +  snes - iterative context obtained from SNESCreate()
3729: -  viewer - the viewer to display the reason


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

3735:    Level: beginner

3737: .keywords: SNES, solve, linear system

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

3741: @*/
3742: PetscErrorCode  SNESReasonView(SNES snes,PetscViewer viewer)
3743: {
3745:   PetscBool      isAscii;

3748:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
3749:   if (isAscii) {
3750:     PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
3751:     if (snes->reason > 0) {
3752:       if (((PetscObject) snes)->prefix) {
3753:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3754:       } else {
3755:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3756:       }
3757:     } else {
3758:       if (((PetscObject) snes)->prefix) {
3759:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3760:       } else {
3761:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3762:       }
3763:     }
3764:     PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
3765:   }
3766:   return(0);
3767: }

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

3774:   Collective on SNES

3776:   Input Parameters:
3777: . snes   - the SNES object

3779:   Level: intermediate

3781: @*/
3782: PetscErrorCode SNESReasonViewFromOptions(SNES snes)
3783: {
3784:   PetscErrorCode    ierr;
3785:   PetscViewer       viewer;
3786:   PetscBool         flg;
3787:   static PetscBool  incall = PETSC_FALSE;
3788:   PetscViewerFormat format;

3791:   if (incall) return(0);
3792:   incall = PETSC_TRUE;
3793:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
3794:   if (flg) {
3795:     PetscViewerPushFormat(viewer,format);
3796:     SNESReasonView(snes,viewer);
3797:     PetscViewerPopFormat(viewer);
3798:     PetscViewerDestroy(&viewer);
3799:   }
3800:   incall = PETSC_FALSE;
3801:   return(0);
3802: }

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

3810:    Collective on SNES

3812:    Input Parameters:
3813: +  snes - the SNES context
3814: .  b - the constant part of the equation F(x) = b, or NULL to use zero.
3815: -  x - the solution vector.

3817:    Notes:
3818:    The user should initialize the vector,x, with the initial guess
3819:    for the nonlinear solve prior to calling SNESSolve.  In particular,
3820:    to employ an initial guess of zero, the user should explicitly set
3821:    this vector to zero by calling VecSet().

3823:    Level: beginner

3825: .keywords: SNES, nonlinear, solve

3827: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
3828: @*/
3829: PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
3830: {
3831:   PetscErrorCode    ierr;
3832:   PetscBool         flg;
3833:   PetscInt          grid;
3834:   Vec               xcreated = NULL;
3835:   DM                dm;


3844:   if (!x) {
3845:     SNESGetDM(snes,&dm);
3846:     DMCreateGlobalVector(dm,&xcreated);
3847:     x    = xcreated;
3848:   }
3849:   SNESViewFromOptions(snes,NULL,"-snes_view_pre");

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

3854:     /* set solution vector */
3855:     if (!grid) {PetscObjectReference((PetscObject)x);}
3856:     VecDestroy(&snes->vec_sol);
3857:     snes->vec_sol = x;
3858:     SNESGetDM(snes,&dm);

3860:     /* set affine vector if provided */
3861:     if (b) { PetscObjectReference((PetscObject)b); }
3862:     VecDestroy(&snes->vec_rhs);
3863:     snes->vec_rhs = b;

3865:     if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
3866:     if (snes->vec_rhs  == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
3867:     if (!snes->vec_sol_update /* && snes->vec_sol */) {
3868:       VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
3869:       PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
3870:     }
3871:     DMShellSetGlobalVector(dm,snes->vec_sol);
3872:     SNESSetUp(snes);

3874:     if (!grid) {
3875:       if (snes->ops->computeinitialguess) {
3876:         (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
3877:       }
3878:     }

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

3883:     PetscLogEventBegin(SNES_Solve,snes,0,0,0);
3884:     (*snes->ops->solve)(snes);
3885:     PetscLogEventEnd(SNES_Solve,snes,0,0,0);
3886:     if (snes->domainerror) {
3887:       snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
3888:       snes->domainerror = PETSC_FALSE;
3889:     }
3890:     if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");

3892:     if (snes->lagjac_persist) snes->jac_iter += snes->iter;
3893:     if (snes->lagpre_persist) snes->pre_iter += snes->iter;

3895:     flg  = PETSC_FALSE;
3896:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,NULL);
3897:     if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
3898:     SNESReasonViewFromOptions(snes);

3900:     if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
3901:     if (grid <  snes->gridsequence) {
3902:       DM  fine;
3903:       Vec xnew;
3904:       Mat interp;

3906:       DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
3907:       if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
3908:       DMCreateInterpolation(snes->dm,fine,&interp,NULL);
3909:       DMCreateGlobalVector(fine,&xnew);
3910:       MatInterpolate(interp,x,xnew);
3911:       DMInterpolate(snes->dm,interp,fine);
3912:       MatDestroy(&interp);
3913:       x    = xnew;

3915:       SNESReset(snes);
3916:       SNESSetDM(snes,fine);
3917:       DMDestroy(&fine);
3918:       PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
3919:     }
3920:   }
3921:   SNESViewFromOptions(snes,NULL,"-snes_view");
3922:   VecViewFromOptions(snes->vec_sol,((PetscObject)snes)->prefix,"-snes_view_solution");

3924:   VecDestroy(&xcreated);
3925:   PetscObjectSAWsBlock((PetscObject)snes);
3926:   return(0);
3927: }

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

3933: /*@C
3934:    SNESSetType - Sets the method for the nonlinear solver.

3936:    Collective on SNES

3938:    Input Parameters:
3939: +  snes - the SNES context
3940: -  type - a known method

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

3946:    Notes:
3947:    See "petsc/include/petscsnes.h" for available methods (for instance)
3948: +    SNESNEWTONLS - Newton's method with line search
3949:      (systems of nonlinear equations)
3950: .    SNESNEWTONTR - Newton's method with trust region
3951:      (systems of nonlinear equations)

3953:   Normally, it is best to use the SNESSetFromOptions() command and then
3954:   set the SNES solver type from the options database rather than by using
3955:   this routine.  Using the options database provides the user with
3956:   maximum flexibility in evaluating the many nonlinear solvers.
3957:   The SNESSetType() routine is provided for those situations where it
3958:   is necessary to set the nonlinear solver independently of the command
3959:   line or options database.  This might be the case, for example, when
3960:   the choice of solver changes during the execution of the program,
3961:   and the user's application is taking responsibility for choosing the
3962:   appropriate method.

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

3967:   Level: intermediate

3969: .keywords: SNES, set, type

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

3973: @*/
3974: PetscErrorCode  SNESSetType(SNES snes,SNESType type)
3975: {
3976:   PetscErrorCode ierr,(*r)(SNES);
3977:   PetscBool      match;


3983:   PetscObjectTypeCompare((PetscObject)snes,type,&match);
3984:   if (match) return(0);

3986:    PetscFunctionListFind(SNESList,type,&r);
3987:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
3988:   /* Destroy the previous private SNES context */
3989:   if (snes->ops->destroy) {
3990:     (*(snes)->ops->destroy)(snes);
3991:     snes->ops->destroy = NULL;
3992:   }
3993:   /* Reinitialize function pointers in SNESOps structure */
3994:   snes->ops->setup          = 0;
3995:   snes->ops->solve          = 0;
3996:   snes->ops->view           = 0;
3997:   snes->ops->setfromoptions = 0;
3998:   snes->ops->destroy        = 0;
3999:   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
4000:   snes->setupcalled = PETSC_FALSE;

4002:   PetscObjectChangeTypeName((PetscObject)snes,type);
4003:   (*r)(snes);
4004:   return(0);
4005: }

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

4012:    Not Collective

4014:    Input Parameter:
4015: .  snes - nonlinear solver context

4017:    Output Parameter:
4018: .  type - SNES method (a character string)

4020:    Level: intermediate

4022: .keywords: SNES, nonlinear, get, type, name
4023: @*/
4024: PetscErrorCode  SNESGetType(SNES snes,SNESType *type)
4025: {
4029:   *type = ((PetscObject)snes)->type_name;
4030:   return(0);
4031: }

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

4038:   Logically Collective on SNES and Vec

4040:   Input Parameters:
4041: + snes - the SNES context obtained from SNESCreate()
4042: - u    - the solution vector

4044:   Level: beginner

4046: .keywords: SNES, set, solution
4047: @*/
4048: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4049: {
4050:   DM             dm;

4056:   PetscObjectReference((PetscObject) u);
4057:   VecDestroy(&snes->vec_sol);

4059:   snes->vec_sol = u;

4061:   SNESGetDM(snes, &dm);
4062:   DMShellSetGlobalVector(dm, u);
4063:   return(0);
4064: }

4068: /*@
4069:    SNESGetSolution - Returns the vector where the approximate solution is
4070:    stored. This is the fine grid solution when using SNESSetGridSequence().

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

4074:    Input Parameter:
4075: .  snes - the SNES context

4077:    Output Parameter:
4078: .  x - the solution

4080:    Level: intermediate

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

4084: .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
4085: @*/
4086: PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
4087: {
4091:   *x = snes->vec_sol;
4092:   return(0);
4093: }

4097: /*@
4098:    SNESGetSolutionUpdate - Returns the vector where the solution update is
4099:    stored.

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

4103:    Input Parameter:
4104: .  snes - the SNES context

4106:    Output Parameter:
4107: .  x - the solution update

4109:    Level: advanced

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

4113: .seealso: SNESGetSolution(), SNESGetFunction()
4114: @*/
4115: PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
4116: {
4120:   *x = snes->vec_sol_update;
4121:   return(0);
4122: }

4126: /*@C
4127:    SNESGetFunction - Returns the vector where the function is stored.

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

4131:    Input Parameter:
4132: .  snes - the SNES context

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

4139:    Level: advanced

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

4143: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4144: @*/
4145: PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4146: {
4148:   DM             dm;

4152:   if (r) {
4153:     if (!snes->vec_func) {
4154:       if (snes->vec_rhs) {
4155:         VecDuplicate(snes->vec_rhs,&snes->vec_func);
4156:       } else if (snes->vec_sol) {
4157:         VecDuplicate(snes->vec_sol,&snes->vec_func);
4158:       } else if (snes->dm) {
4159:         DMCreateGlobalVector(snes->dm,&snes->vec_func);
4160:       }
4161:     }
4162:     *r = snes->vec_func;
4163:   }
4164:   SNESGetDM(snes,&dm);
4165:   DMSNESGetFunction(dm,f,ctx);
4166:   return(0);
4167: }

4169: /*@C
4170:    SNESGetNGS - Returns the NGS function and context.

4172:    Input Parameter:
4173: .  snes - the SNES context

4175:    Output Parameter:
4176: +  f - the function (or NULL) see SNESNGSFunction for details
4177: -  ctx    - the function context (or NULL)

4179:    Level: advanced

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

4183: .seealso: SNESSetNGS(), SNESGetFunction()
4184: @*/

4188: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4189: {
4191:   DM             dm;

4195:   SNESGetDM(snes,&dm);
4196:   DMSNESGetNGS(dm,f,ctx);
4197:   return(0);
4198: }

4202: /*@C
4203:    SNESSetOptionsPrefix - Sets the prefix used for searching for all
4204:    SNES options in the database.

4206:    Logically Collective on SNES

4208:    Input Parameter:
4209: +  snes - the SNES context
4210: -  prefix - the prefix to prepend to all option names

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

4216:    Level: advanced

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

4220: .seealso: SNESSetFromOptions()
4221: @*/
4222: PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
4223: {

4228:   PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4229:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4230:   if (snes->linesearch) {
4231:     SNESGetLineSearch(snes,&snes->linesearch);
4232:     PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4233:   }
4234:   KSPSetOptionsPrefix(snes->ksp,prefix);
4235:   return(0);
4236: }

4240: /*@C
4241:    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4242:    SNES options in the database.

4244:    Logically Collective on SNES

4246:    Input Parameters:
4247: +  snes - the SNES context
4248: -  prefix - the prefix to prepend to all option names

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

4254:    Level: advanced

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

4258: .seealso: SNESGetOptionsPrefix()
4259: @*/
4260: PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4261: {

4266:   PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4267:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4268:   if (snes->linesearch) {
4269:     SNESGetLineSearch(snes,&snes->linesearch);
4270:     PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4271:   }
4272:   KSPAppendOptionsPrefix(snes->ksp,prefix);
4273:   return(0);
4274: }

4278: /*@C
4279:    SNESGetOptionsPrefix - Sets the prefix used for searching for all
4280:    SNES options in the database.

4282:    Not Collective

4284:    Input Parameter:
4285: .  snes - the SNES context

4287:    Output Parameter:
4288: .  prefix - pointer to the prefix string used

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

4293:    Level: advanced

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

4297: .seealso: SNESAppendOptionsPrefix()
4298: @*/
4299: PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4300: {

4305:   PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4306:   return(0);
4307: }


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

4315:    Not collective

4317:    Input Parameters:
4318: +  name_solver - name of a new user-defined solver
4319: -  routine_create - routine to create method context

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

4324:    Sample usage:
4325: .vb
4326:    SNESRegister("my_solver",MySolverCreate);
4327: .ve

4329:    Then, your solver can be chosen with the procedural interface via
4330: $     SNESSetType(snes,"my_solver")
4331:    or at runtime via the option
4332: $     -snes_type my_solver

4334:    Level: advanced

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

4338: .keywords: SNES, nonlinear, register

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

4342:   Level: advanced
4343: @*/
4344: PetscErrorCode  SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4345: {

4349:   PetscFunctionListAdd(&SNESList,sname,function);
4350:   return(0);
4351: }

4355: PetscErrorCode  SNESTestLocalMin(SNES snes)
4356: {
4358:   PetscInt       N,i,j;
4359:   Vec            u,uh,fh;
4360:   PetscScalar    value;
4361:   PetscReal      norm;

4364:   SNESGetSolution(snes,&u);
4365:   VecDuplicate(u,&uh);
4366:   VecDuplicate(u,&fh);

4368:   /* currently only works for sequential */
4369:   PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4370:   VecGetSize(u,&N);
4371:   for (i=0; i<N; i++) {
4372:     VecCopy(u,uh);
4373:     PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4374:     for (j=-10; j<11; j++) {
4375:       value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4376:       VecSetValue(uh,i,value,ADD_VALUES);
4377:       SNESComputeFunction(snes,uh,fh);
4378:       VecNorm(fh,NORM_2,&norm);
4379:       PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);
4380:       value = -value;
4381:       VecSetValue(uh,i,value,ADD_VALUES);
4382:     }
4383:   }
4384:   VecDestroy(&uh);
4385:   VecDestroy(&fh);
4386:   return(0);
4387: }

4391: /*@
4392:    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4393:    computing relative tolerance for linear solvers within an inexact
4394:    Newton method.

4396:    Logically Collective on SNES

4398:    Input Parameters:
4399: +  snes - SNES context
4400: -  flag - PETSC_TRUE or PETSC_FALSE

4402:     Options Database:
4403: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4404: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
4405: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4406: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4407: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
4408: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
4409: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4410: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

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

4419:    Level: advanced

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

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

4427: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4428: @*/
4429: PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool flag)
4430: {
4434:   snes->ksp_ewconv = flag;
4435:   return(0);
4436: }

4440: /*@
4441:    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
4442:    for computing relative tolerance for linear solvers within an
4443:    inexact Newton method.

4445:    Not Collective

4447:    Input Parameter:
4448: .  snes - SNES context

4450:    Output Parameter:
4451: .  flag - PETSC_TRUE or PETSC_FALSE

4453:    Notes:
4454:    Currently, the default is to use a constant relative tolerance for
4455:    the inner linear solvers.  Alternatively, one can use the
4456:    Eisenstat-Walker method, where the relative convergence tolerance
4457:    is reset at each Newton iteration according progress of the nonlinear
4458:    solver.

4460:    Level: advanced

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

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

4468: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4469: @*/
4470: PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
4471: {
4475:   *flag = snes->ksp_ewconv;
4476:   return(0);
4477: }

4481: /*@
4482:    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
4483:    convergence criteria for the linear solvers within an inexact
4484:    Newton method.

4486:    Logically Collective on SNES

4488:    Input Parameters:
4489: +    snes - SNES context
4490: .    version - version 1, 2 (default is 2) or 3
4491: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4492: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4493: .    gamma - multiplicative factor for version 2 rtol computation
4494:              (0 <= gamma2 <= 1)
4495: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
4496: .    alpha2 - power for safeguard
4497: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

4499:    Note:
4500:    Version 3 was contributed by Luis Chacon, June 2006.

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

4504:    Level: advanced

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

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

4513: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
4514: @*/
4515: PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
4516: {
4517:   SNESKSPEW *kctx;

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

4531:   if (version != PETSC_DEFAULT)   kctx->version   = version;
4532:   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
4533:   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
4534:   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
4535:   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
4536:   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
4537:   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;

4539:   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);
4540:   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);
4541:   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);
4542:   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);
4543:   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);
4544:   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);
4545:   return(0);
4546: }

4550: /*@
4551:    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
4552:    convergence criteria for the linear solvers within an inexact
4553:    Newton method.

4555:    Not Collective

4557:    Input Parameters:
4558:      snes - SNES context

4560:    Output Parameters:
4561: +    version - version 1, 2 (default is 2) or 3
4562: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4563: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4564: .    gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
4565: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
4566: .    alpha2 - power for safeguard
4567: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

4569:    Level: advanced

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

4573: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
4574: @*/
4575: PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
4576: {
4577:   SNESKSPEW *kctx;

4581:   kctx = (SNESKSPEW*)snes->kspconvctx;
4582:   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4583:   if (version)   *version   = kctx->version;
4584:   if (rtol_0)    *rtol_0    = kctx->rtol_0;
4585:   if (rtol_max)  *rtol_max  = kctx->rtol_max;
4586:   if (gamma)     *gamma     = kctx->gamma;
4587:   if (alpha)     *alpha     = kctx->alpha;
4588:   if (alpha2)    *alpha2    = kctx->alpha2;
4589:   if (threshold) *threshold = kctx->threshold;
4590:   return(0);
4591: }

4595:  PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4596: {
4598:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
4599:   PetscReal      rtol  = PETSC_DEFAULT,stol;

4602:   if (!snes->ksp_ewconv) return(0);
4603:   if (!snes->iter) {
4604:     rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
4605:     VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
4606:   }
4607:   else {
4608:     if (kctx->version == 1) {
4609:       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
4610:       if (rtol < 0.0) rtol = -rtol;
4611:       stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
4612:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4613:     } else if (kctx->version == 2) {
4614:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4615:       stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
4616:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4617:     } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
4618:       rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4619:       /* safeguard: avoid sharp decrease of rtol */
4620:       stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
4621:       stol = PetscMax(rtol,stol);
4622:       rtol = PetscMin(kctx->rtol_0,stol);
4623:       /* safeguard: avoid oversolving */
4624:       stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
4625:       stol = PetscMax(rtol,stol);
4626:       rtol = PetscMin(kctx->rtol_0,stol);
4627:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
4628:   }
4629:   /* safeguard: avoid rtol greater than one */
4630:   rtol = PetscMin(rtol,kctx->rtol_max);
4631:   KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
4632:   PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
4633:   return(0);
4634: }

4638: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4639: {
4641:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
4642:   PCSide         pcside;
4643:   Vec            lres;

4646:   if (!snes->ksp_ewconv) return(0);
4647:   KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
4648:   kctx->norm_last = snes->norm;
4649:   if (kctx->version == 1) {
4650:     PC        pc;
4651:     PetscBool isNone;

4653:     KSPGetPC(ksp, &pc);
4654:     PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
4655:     KSPGetPCSide(ksp,&pcside);
4656:      if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
4657:       /* KSP residual is true linear residual */
4658:       KSPGetResidualNorm(ksp,&kctx->lresid_last);
4659:     } else {
4660:       /* KSP residual is preconditioned residual */
4661:       /* compute true linear residual norm */
4662:       VecDuplicate(b,&lres);
4663:       MatMult(snes->jacobian,x,lres);
4664:       VecAYPX(lres,-1.0,b);
4665:       VecNorm(lres,NORM_2,&kctx->lresid_last);
4666:       VecDestroy(&lres);
4667:     }
4668:   }
4669:   return(0);
4670: }

4674: /*@
4675:    SNESGetKSP - Returns the KSP context for a SNES solver.

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

4679:    Input Parameter:
4680: .  snes - the SNES context

4682:    Output Parameter:
4683: .  ksp - the KSP context

4685:    Notes:
4686:    The user can then directly manipulate the KSP context to set various
4687:    options, etc.  Likewise, the user can then extract and manipulate the
4688:    PC contexts as well.

4690:    Level: beginner

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

4694: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
4695: @*/
4696: PetscErrorCode  SNESGetKSP(SNES snes,KSP *ksp)
4697: {


4704:   if (!snes->ksp) {
4705:     PetscBool monitor = PETSC_FALSE;

4707:     KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
4708:     PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
4709:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);

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

4714:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
4715:     if (monitor) {
4716:       KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
4717:     }
4718:     monitor = PETSC_FALSE;
4719:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
4720:     if (monitor) {
4721:       PetscObject *objs;
4722:       KSPMonitorSNESLGResidualNormCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
4723:       objs[0] = (PetscObject) snes;
4724:       KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
4725:     }
4726:   }
4727:   *ksp = snes->ksp;
4728:   return(0);
4729: }


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

4738:    Logically Collective on SNES

4740:    Input Parameters:
4741: +  snes - the preconditioner context
4742: -  dm - the dm

4744:    Level: intermediate

4746: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
4747: @*/
4748: PetscErrorCode  SNESSetDM(SNES snes,DM dm)
4749: {
4751:   KSP            ksp;
4752:   DMSNES         sdm;

4756:   if (dm) {PetscObjectReference((PetscObject)dm);}
4757:   if (snes->dm) {               /* Move the DMSNES context over to the new DM unless the new DM already has one */
4758:     if (snes->dm->dmsnes && snes->dmAuto && !dm->dmsnes) {
4759:       DMCopyDMSNES(snes->dm,dm);
4760:       DMGetDMSNES(snes->dm,&sdm);
4761:       if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
4762:     }
4763:     DMDestroy(&snes->dm);
4764:   }
4765:   snes->dm     = dm;
4766:   snes->dmAuto = PETSC_FALSE;

4768:   SNESGetKSP(snes,&ksp);
4769:   KSPSetDM(ksp,dm);
4770:   KSPSetDMActive(ksp,PETSC_FALSE);
4771:   if (snes->pc) {
4772:     SNESSetDM(snes->pc, snes->dm);
4773:     SNESSetNPCSide(snes,snes->pcside);
4774:   }
4775:   return(0);
4776: }

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

4783:    Not Collective but DM obtained is parallel on SNES

4785:    Input Parameter:
4786: . snes - the preconditioner context

4788:    Output Parameter:
4789: .  dm - the dm

4791:    Level: intermediate

4793: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
4794: @*/
4795: PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
4796: {

4801:   if (!snes->dm) {
4802:     DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
4803:     snes->dmAuto = PETSC_TRUE;
4804:   }
4805:   *dm = snes->dm;
4806:   return(0);
4807: }

4811: /*@
4812:   SNESSetNPC - Sets the nonlinear preconditioner to be used.

4814:   Collective on SNES

4816:   Input Parameters:
4817: + snes - iterative context obtained from SNESCreate()
4818: - pc   - the preconditioner object

4820:   Notes:
4821:   Use SNESGetNPC() to retrieve the preconditioner context (for example,
4822:   to configure it using the API).

4824:   Level: developer

4826: .keywords: SNES, set, precondition
4827: .seealso: SNESGetNPC(), SNESHasNPC()
4828: @*/
4829: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
4830: {

4837:   PetscObjectReference((PetscObject) pc);
4838:   SNESDestroy(&snes->pc);
4839:   snes->pc = pc;
4840:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->pc);
4841:   return(0);
4842: }

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

4849:   Not Collective

4851:   Input Parameter:
4852: . snes - iterative context obtained from SNESCreate()

4854:   Output Parameter:
4855: . pc - preconditioner context

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

4859:   Level: developer

4861: .keywords: SNES, get, preconditioner
4862: .seealso: SNESSetNPC(), SNESHasNPC()
4863: @*/
4864: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
4865: {
4867:   const char     *optionsprefix;

4872:   if (!snes->pc) {
4873:     SNESCreate(PetscObjectComm((PetscObject)snes),&snes->pc);
4874:     PetscObjectIncrementTabLevel((PetscObject)snes->pc,(PetscObject)snes,1);
4875:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->pc);
4876:     SNESGetOptionsPrefix(snes,&optionsprefix);
4877:     SNESSetOptionsPrefix(snes->pc,optionsprefix);
4878:     SNESAppendOptionsPrefix(snes->pc,"npc_");
4879:     SNESSetCountersReset(snes->pc,PETSC_FALSE);
4880:   }
4881:   *pc = snes->pc;
4882:   return(0);
4883: }

4887: /*@
4888:   SNESHasNPC - Returns whether a nonlinear preconditioner exists

4890:   Not Collective

4892:   Input Parameter:
4893: . snes - iterative context obtained from SNESCreate()

4895:   Output Parameter:
4896: . has_npc - whether the SNES has an NPC or not

4898:   Level: developer

4900: .keywords: SNES, has, preconditioner
4901: .seealso: SNESSetNPC(), SNESGetNPC()
4902: @*/
4903: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
4904: {
4907:   *has_npc = (PetscBool) (snes->pc != NULL);
4908:   return(0);
4909: }

4913: /*@
4914:     SNESSetNPCSide - Sets the preconditioning side.

4916:     Logically Collective on SNES

4918:     Input Parameter:
4919: .   snes - iterative context obtained from SNESCreate()

4921:     Output Parameter:
4922: .   side - the preconditioning side, where side is one of
4923: .vb
4924:       PC_LEFT - left preconditioning (default)
4925:       PC_RIGHT - right preconditioning
4926: .ve

4928:     Options Database Keys:
4929: .   -snes_pc_side <right,left>

4931:     Level: intermediate

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

4935: .seealso: SNESGetNPCSide(), KSPSetPCSide()
4936: @*/
4937: PetscErrorCode  SNESSetNPCSide(SNES snes,PCSide side)
4938: {
4942:   snes->pcside = side;
4943:   return(0);
4944: }

4948: /*@
4949:     SNESGetNPCSide - Gets the preconditioning side.

4951:     Not Collective

4953:     Input Parameter:
4954: .   snes - iterative context obtained from SNESCreate()

4956:     Output Parameter:
4957: .   side - the preconditioning side, where side is one of
4958: .vb
4959:       PC_LEFT - left preconditioning (default)
4960:       PC_RIGHT - right preconditioning
4961: .ve

4963:     Level: intermediate

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

4967: .seealso: SNESSetNPCSide(), KSPGetPCSide()
4968: @*/
4969: PetscErrorCode  SNESGetNPCSide(SNES snes,PCSide *side)
4970: {
4974:   *side = snes->pcside;
4975:   return(0);
4976: }

4980: /*@
4981:   SNESSetLineSearch - Sets the linesearch on the SNES instance.

4983:   Collective on SNES

4985:   Input Parameters:
4986: + snes - iterative context obtained from SNESCreate()
4987: - linesearch   - the linesearch object

4989:   Notes:
4990:   Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
4991:   to configure it using the API).

4993:   Level: developer

4995: .keywords: SNES, set, linesearch
4996: .seealso: SNESGetLineSearch()
4997: @*/
4998: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
4999: {

5006:   PetscObjectReference((PetscObject) linesearch);
5007:   SNESLineSearchDestroy(&snes->linesearch);

5009:   snes->linesearch = linesearch;

5011:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5012:   return(0);
5013: }

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

5021:   Not Collective

5023:   Input Parameter:
5024: . snes - iterative context obtained from SNESCreate()

5026:   Output Parameter:
5027: . linesearch - linesearch context

5029:   Level: beginner

5031: .keywords: SNES, get, linesearch
5032: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5033: @*/
5034: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5035: {
5037:   const char     *optionsprefix;

5042:   if (!snes->linesearch) {
5043:     SNESGetOptionsPrefix(snes, &optionsprefix);
5044:     SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5045:     SNESLineSearchSetSNES(snes->linesearch, snes);
5046:     SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5047:     PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5048:     PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5049:   }
5050:   *linesearch = snes->linesearch;
5051:   return(0);
5052: }

5054: #if defined(PETSC_HAVE_MATLAB_ENGINE)
5055: #include <mex.h>

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

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

5064:    Collective on SNES

5066:    Input Parameters:
5067: +  snes - the SNES context
5068: -  x - input vector

5070:    Output Parameter:
5071: .  y - function vector, as set by SNESSetFunction()

5073:    Notes:
5074:    SNESComputeFunction() is typically used within nonlinear solvers
5075:    implementations, so most users would not generally call this routine
5076:    themselves.

5078:    Level: developer

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

5082: .seealso: SNESSetFunction(), SNESGetFunction()
5083: */
5084: PetscErrorCode  SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
5085: {
5086:   PetscErrorCode    ierr;
5087:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5088:   int               nlhs  = 1,nrhs = 5;
5089:   mxArray           *plhs[1],*prhs[5];
5090:   long long int     lx = 0,ly = 0,ls = 0;


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

5101:   PetscMemcpy(&ls,&snes,sizeof(snes));
5102:   PetscMemcpy(&lx,&x,sizeof(x));
5103:   PetscMemcpy(&ly,&y,sizeof(x));
5104:   prhs[0] = mxCreateDoubleScalar((double)ls);
5105:   prhs[1] = mxCreateDoubleScalar((double)lx);
5106:   prhs[2] = mxCreateDoubleScalar((double)ly);
5107:   prhs[3] = mxCreateString(sctx->funcname);
5108:   prhs[4] = sctx->ctx;
5109:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");
5110:   mxGetScalar(plhs[0]);
5111:   mxDestroyArray(prhs[0]);
5112:   mxDestroyArray(prhs[1]);
5113:   mxDestroyArray(prhs[2]);
5114:   mxDestroyArray(prhs[3]);
5115:   mxDestroyArray(plhs[0]);
5116:   return(0);
5117: }

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

5126:    Logically Collective on SNES

5128:    Input Parameters:
5129: +  snes - the SNES context
5130: .  r - vector to store function value
5131: -  f - function evaluation routine

5133:    Notes:
5134:    The Newton-like methods typically solve linear systems of the form
5135: $      f'(x) x = -f(x),
5136:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

5138:    Level: beginner

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

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

5144: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5145: */
5146: PetscErrorCode  SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx)
5147: {
5148:   PetscErrorCode    ierr;
5149:   SNESMatlabContext *sctx;

5152:   /* currently sctx is memory bleed */
5153:   PetscNew(&sctx);
5154:   PetscStrallocpy(f,&sctx->funcname);
5155:   /*
5156:      This should work, but it doesn't
5157:   sctx->ctx = ctx;
5158:   mexMakeArrayPersistent(sctx->ctx);
5159:   */
5160:   sctx->ctx = mxDuplicateArray(ctx);
5161:   SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);
5162:   return(0);
5163: }

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

5170:    Collective on SNES

5172:    Input Parameters:
5173: +  snes - the SNES context
5174: .  x - input vector
5175: .  A, B - the matrices
5176: -  ctx - user context

5178:    Level: developer

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

5182: .seealso: SNESSetFunction(), SNESGetFunction()
5183: @*/
5184: PetscErrorCode  SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx)
5185: {
5186:   PetscErrorCode    ierr;
5187:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5188:   int               nlhs  = 2,nrhs = 6;
5189:   mxArray           *plhs[2],*prhs[6];
5190:   long long int     lx = 0,lA = 0,ls = 0, lB = 0;


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

5198:   PetscMemcpy(&ls,&snes,sizeof(snes));
5199:   PetscMemcpy(&lx,&x,sizeof(x));
5200:   PetscMemcpy(&lA,A,sizeof(x));
5201:   PetscMemcpy(&lB,B,sizeof(x));
5202:   prhs[0] = mxCreateDoubleScalar((double)ls);
5203:   prhs[1] = mxCreateDoubleScalar((double)lx);
5204:   prhs[2] = mxCreateDoubleScalar((double)lA);
5205:   prhs[3] = mxCreateDoubleScalar((double)lB);
5206:   prhs[4] = mxCreateString(sctx->funcname);
5207:   prhs[5] = sctx->ctx;
5208:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");
5209:   mxGetScalar(plhs[0]);
5210:   mxDestroyArray(prhs[0]);
5211:   mxDestroyArray(prhs[1]);
5212:   mxDestroyArray(prhs[2]);
5213:   mxDestroyArray(prhs[3]);
5214:   mxDestroyArray(prhs[4]);
5215:   mxDestroyArray(plhs[0]);
5216:   mxDestroyArray(plhs[1]);
5217:   return(0);
5218: }

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

5227:    Logically Collective on SNES

5229:    Input Parameters:
5230: +  snes - the SNES context
5231: .  A,B - Jacobian matrices
5232: .  J - function evaluation routine
5233: -  ctx - user context

5235:    Level: developer

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

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

5241: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J
5242: */
5243: PetscErrorCode  SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx)
5244: {
5245:   PetscErrorCode    ierr;
5246:   SNESMatlabContext *sctx;

5249:   /* currently sctx is memory bleed */
5250:   PetscNew(&sctx);
5251:   PetscStrallocpy(J,&sctx->funcname);
5252:   /*
5253:      This should work, but it doesn't
5254:   sctx->ctx = ctx;
5255:   mexMakeArrayPersistent(sctx->ctx);
5256:   */
5257:   sctx->ctx = mxDuplicateArray(ctx);
5258:   SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);
5259:   return(0);
5260: }

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

5267:    Collective on SNES

5269: .seealso: SNESSetFunction(), SNESGetFunction()
5270: @*/
5271: PetscErrorCode  SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
5272: {
5273:   PetscErrorCode    ierr;
5274:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5275:   int               nlhs  = 1,nrhs = 6;
5276:   mxArray           *plhs[1],*prhs[6];
5277:   long long int     lx = 0,ls = 0;
5278:   Vec               x  = snes->vec_sol;


5283:   PetscMemcpy(&ls,&snes,sizeof(snes));
5284:   PetscMemcpy(&lx,&x,sizeof(x));
5285:   prhs[0] = mxCreateDoubleScalar((double)ls);
5286:   prhs[1] = mxCreateDoubleScalar((double)it);
5287:   prhs[2] = mxCreateDoubleScalar((double)fnorm);
5288:   prhs[3] = mxCreateDoubleScalar((double)lx);
5289:   prhs[4] = mxCreateString(sctx->funcname);
5290:   prhs[5] = sctx->ctx;
5291:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");
5292:   mxGetScalar(plhs[0]);
5293:   mxDestroyArray(prhs[0]);
5294:   mxDestroyArray(prhs[1]);
5295:   mxDestroyArray(prhs[2]);
5296:   mxDestroyArray(prhs[3]);
5297:   mxDestroyArray(prhs[4]);
5298:   mxDestroyArray(plhs[0]);
5299:   return(0);
5300: }

5304: /*
5305:    SNESMonitorSetMatlab - Sets the monitor function from MATLAB

5307:    Level: developer

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

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

5313: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5314: */
5315: PetscErrorCode  SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx)
5316: {
5317:   PetscErrorCode    ierr;
5318:   SNESMatlabContext *sctx;

5321:   /* currently sctx is memory bleed */
5322:   PetscNew(&sctx);
5323:   PetscStrallocpy(f,&sctx->funcname);
5324:   /*
5325:      This should work, but it doesn't
5326:   sctx->ctx = ctx;
5327:   mexMakeArrayPersistent(sctx->ctx);
5328:   */
5329:   sctx->ctx = mxDuplicateArray(ctx);
5330:   SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);
5331:   return(0);
5332: }

5334: #endif