Actual source code: snes.c

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

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

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

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

 17:    Logically Collective on SNES

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

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

 26:    Level: intermediate

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

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

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

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

 50:    Not Collective

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

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

 58:    Level: intermediate

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

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

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

 79:    Logically Collective on SNES

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

 84:    Level: advanced

 86: .keywords: SNES, view

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

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

103:    Logically Collective on SNES

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

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

111:    Level: advanced

113: .keywords: SNES, view

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

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

131:   Collective on PetscViewer

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

138:    Level: intermediate

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

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

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

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

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

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

194:    Collective on SNES

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

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

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

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

214:    Level: beginner

216: .keywords: SNES, view

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

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

240:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
241:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
242:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
243:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
244: #if defined(PETSC_HAVE_SAWS)
245:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
246: #endif
247:   if (iascii) {
248:     SNESNormSchedule normschedule;

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

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

311:     PetscViewerDrawGetDraw(viewer,0,&draw);
312:     PetscDrawGetCurrentPoint(draw,&x,&y);
313:     PetscStrcpy(str,"SNES: ");
314:     PetscStrcat(str,((PetscObject)snes)->type_name);
315:     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 is automatically called with FALSE

1352:    Level: developer

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

1356: .seealso:  SNESGetNumberFunctionEvals(), SNESGetNumberLinearSolveIterations(), 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,_p_SNES,struct _SNESOps,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:   VecValidValues(y,3,PETSC_FALSE);
2067:   return(0);
2068: }

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

2075:    Collective on SNES

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

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

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

2090:    Level: developer

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

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

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

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

2129:    Collective on SNES and Mat

2131:    Input Parameters:
2132: +  snes - the SNES context
2133: -  x - input vector

2135:    Output Parameters:
2136: +  A - Jacobian matrix
2137: -  B - optional preconditioning matrix

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


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

2159:    Level: developer

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

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

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

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

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

2185:   if (snes->lagjacobian == -2) {
2186:     snes->lagjacobian = -1;

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

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

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

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

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

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

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

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

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

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

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

2415:    Level: intermediate

2417: .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2418: M*/

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

2426:    Logically Collective on SNES and Mat

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

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

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

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

2446:    Level: beginner

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

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

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

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

2476:     snes->jacobian_pre = Pmat;
2477:   }
2478:   return(0);
2479: }

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

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

2489:    Input Parameter:
2490: .  snes - the nonlinear solver context

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

2498:    Level: advanced

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

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

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

2525:    Collective on SNES

2527:    Input Parameters:
2528: .  snes - the SNES context

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

2537:    Level: advanced

2539: .keywords: SNES, nonlinear, setup

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

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

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

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

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

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

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

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

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

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

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

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

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

2635:   snes->jac_iter = 0;
2636:   snes->pre_iter = 0;

2638:   if (snes->ops->setup) {
2639:     (*snes->ops->setup)(snes);
2640:   }

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

2649:   snes->setupcalled = PETSC_TRUE;
2650:   return(0);
2651: }

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

2658:    Collective on SNES

2660:    Input Parameter:
2661: .  snes - iterative context obtained from SNESCreate()

2663:    Level: intermediate

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

2667: .keywords: SNES, destroy

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

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

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

2692:   if (snes->linesearch) {
2693:     SNESLineSearchReset(snes->linesearch);
2694:   }

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

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

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

2716:    Collective on SNES

2718:    Input Parameter:
2719: .  snes - the SNES context

2721:    Level: beginner

2723: .keywords: SNES, nonlinear, destroy

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

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

2736:   SNESReset((*snes));
2737:   SNESDestroy(&(*snes)->pc);

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

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

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

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

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

2767:    Logically Collective on SNES

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

2774:    Options Database Keys:
2775: .    -snes_lag_preconditioner <lag>

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

2782:    Level: intermediate

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

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

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

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

2805:    Logically Collective on SNES

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

2811:    Options Database Keys:
2812: .    -snes_grid_sequence <steps>

2814:    Level: intermediate

2816:    Notes:
2817:    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.

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

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

2823: @*/
2824: PetscErrorCode  SNESSetGridSequence(SNES snes,PetscInt steps)
2825: {
2829:   snes->gridsequence = steps;
2830:   return(0);
2831: }

2835: /*@
2836:    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt

2838:    Not Collective

2840:    Input Parameter:
2841: .  snes - the SNES context

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

2847:    Options Database Keys:
2848: .    -snes_lag_preconditioner <lag>

2850:    Notes:
2851:    The default is 1
2852:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

2854:    Level: intermediate

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

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

2860: @*/
2861: PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
2862: {
2865:   *lag = snes->lagpreconditioner;
2866:   return(0);
2867: }

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

2875:    Logically Collective on SNES

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

2882:    Options Database Keys:
2883: .    -snes_lag_jacobian <lag>

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

2891:    Level: intermediate

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

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

2897: @*/
2898: PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
2899: {
2902:   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2903:   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2905:   snes->lagjacobian = lag;
2906:   return(0);
2907: }

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

2914:    Not Collective

2916:    Input Parameter:
2917: .  snes - the SNES context

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

2923:    Options Database Keys:
2924: .    -snes_lag_jacobian <lag>

2926:    Notes:
2927:    The default is 1
2928:    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1

2930:    Level: intermediate

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

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

2936: @*/
2937: PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
2938: {
2941:   *lag = snes->lagjacobian;
2942:   return(0);
2943: }

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

2950:    Logically collective on SNES

2952:    Input Parameter:
2953: +  snes - the SNES context
2954: -   flg - jacobian lagging persists if true

2956:    Options Database Keys:
2957: .    -snes_lag_jacobian_persists <flg>

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

2963:    Level: developer

2965: .keywords: SNES, nonlinear, lag

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

2969: @*/
2970: PetscErrorCode  SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
2971: {
2975:   snes->lagjac_persist = flg;
2976:   return(0);
2977: }

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

2984:    Logically Collective on SNES

2986:    Input Parameter:
2987: +  snes - the SNES context
2988: -   flg - preconditioner lagging persists if true

2990:    Options Database Keys:
2991: .    -snes_lag_jacobian_persists <flg>

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

2997:    Level: developer

2999: .keywords: SNES, nonlinear, lag

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

3003: @*/
3004: PetscErrorCode  SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3005: {
3009:   snes->lagpre_persist = flg;
3010:   return(0);
3011: }

3015: /*@
3016:    SNESSetTolerances - Sets various parameters used in convergence tests.

3018:    Logically Collective on SNES

3020:    Input Parameters:
3021: +  snes - the SNES context
3022: .  abstol - absolute convergence tolerance
3023: .  rtol - relative convergence tolerance
3024: .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
3025: .  maxit - maximum number of iterations
3026: -  maxf - maximum number of function evaluations

3028:    Options Database Keys:
3029: +    -snes_atol <abstol> - Sets abstol
3030: .    -snes_rtol <rtol> - Sets rtol
3031: .    -snes_stol <stol> - Sets stol
3032: .    -snes_max_it <maxit> - Sets maxit
3033: -    -snes_max_funcs <maxf> - Sets maxf

3035:    Notes:
3036:    The default maximum number of iterations is 50.
3037:    The default maximum number of function evaluations is 1000.

3039:    Level: intermediate

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

3043: .seealso: SNESSetTrustRegionTolerance()
3044: @*/
3045: PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3046: {

3055:   if (abstol != PETSC_DEFAULT) {
3056:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3057:     snes->abstol = abstol;
3058:   }
3059:   if (rtol != PETSC_DEFAULT) {
3060:     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);
3061:     snes->rtol = rtol;
3062:   }
3063:   if (stol != PETSC_DEFAULT) {
3064:     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3065:     snes->stol = stol;
3066:   }
3067:   if (maxit != PETSC_DEFAULT) {
3068:     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3069:     snes->max_its = maxit;
3070:   }
3071:   if (maxf != PETSC_DEFAULT) {
3072:     if (maxf < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
3073:     snes->max_funcs = maxf;
3074:   }
3075:   snes->tolerancesset = PETSC_TRUE;
3076:   return(0);
3077: }

3081: /*@
3082:    SNESGetTolerances - Gets various parameters used in convergence tests.

3084:    Not Collective

3086:    Input Parameters:
3087: +  snes - the SNES context
3088: .  atol - absolute convergence tolerance
3089: .  rtol - relative convergence tolerance
3090: .  stol -  convergence tolerance in terms of the norm
3091:            of the change in the solution between steps
3092: .  maxit - maximum number of iterations
3093: -  maxf - maximum number of function evaluations

3095:    Notes:
3096:    The user can specify NULL for any parameter that is not needed.

3098:    Level: intermediate

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

3102: .seealso: SNESSetTolerances()
3103: @*/
3104: PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3105: {
3108:   if (atol)  *atol  = snes->abstol;
3109:   if (rtol)  *rtol  = snes->rtol;
3110:   if (stol)  *stol  = snes->stol;
3111:   if (maxit) *maxit = snes->max_its;
3112:   if (maxf)  *maxf  = snes->max_funcs;
3113:   return(0);
3114: }

3118: /*@
3119:    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.

3121:    Logically Collective on SNES

3123:    Input Parameters:
3124: +  snes - the SNES context
3125: -  tol - tolerance

3127:    Options Database Key:
3128: .  -snes_trtol <tol> - Sets tol

3130:    Level: intermediate

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

3134: .seealso: SNESSetTolerances()
3135: @*/
3136: PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3137: {
3141:   snes->deltatol = tol;
3142:   return(0);
3143: }

3145: /*
3146:    Duplicate the lg monitors for SNES from KSP; for some reason with
3147:    dynamic libraries things don't work under Sun4 if we just use
3148:    macros instead of functions
3149: */
3152: PetscErrorCode  SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,PetscObject *objs)
3153: {

3158:   KSPMonitorLGResidualNorm((KSP)snes,it,norm,objs);
3159:   return(0);
3160: }

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

3169:   KSPMonitorLGResidualNormCreate(host,label,x,y,m,n,draw);
3170:   return(0);
3171: }

3175: PetscErrorCode  SNESMonitorLGDestroy(PetscObject **objs)
3176: {

3180:   KSPMonitorLGResidualNormDestroy(objs);
3181:   return(0);
3182: }

3184: extern PetscErrorCode  SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3187: PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3188: {
3189:   PetscDrawLG      lg;
3190:   PetscErrorCode   ierr;
3191:   PetscReal        x,y,per;
3192:   PetscViewer      v = (PetscViewer)monctx;
3193:   static PetscReal prev; /* should be in the context */
3194:   PetscDraw        draw;

3197:   PetscViewerDrawGetDrawLG(v,0,&lg);
3198:   if (!n) {PetscDrawLGReset(lg);}
3199:   PetscDrawLGGetDraw(lg,&draw);
3200:   PetscDrawSetTitle(draw,"Residual norm");
3201:   x    = (PetscReal)n;
3202:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3203:   else y = -15.0;
3204:   PetscDrawLGAddPoint(lg,&x,&y);
3205:   if (n < 20 || !(n % 5)) {
3206:     PetscDrawLGDraw(lg);
3207:   }

3209:   PetscViewerDrawGetDrawLG(v,1,&lg);
3210:   if (!n) {PetscDrawLGReset(lg);}
3211:   PetscDrawLGGetDraw(lg,&draw);
3212:   PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3213:    SNESMonitorRange_Private(snes,n,&per);
3214:   x    = (PetscReal)n;
3215:   y    = 100.0*per;
3216:   PetscDrawLGAddPoint(lg,&x,&y);
3217:   if (n < 20 || !(n % 5)) {
3218:     PetscDrawLGDraw(lg);
3219:   }

3221:   PetscViewerDrawGetDrawLG(v,2,&lg);
3222:   if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3223:   PetscDrawLGGetDraw(lg,&draw);
3224:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3225:   x    = (PetscReal)n;
3226:   y    = (prev - rnorm)/prev;
3227:   PetscDrawLGAddPoint(lg,&x,&y);
3228:   if (n < 20 || !(n % 5)) {
3229:     PetscDrawLGDraw(lg);
3230:   }

3232:   PetscViewerDrawGetDrawLG(v,3,&lg);
3233:   if (!n) {PetscDrawLGReset(lg);}
3234:   PetscDrawLGGetDraw(lg,&draw);
3235:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3236:   x    = (PetscReal)n;
3237:   y    = (prev - rnorm)/(prev*per);
3238:   if (n > 2) { /*skip initial crazy value */
3239:     PetscDrawLGAddPoint(lg,&x,&y);
3240:   }
3241:   if (n < 20 || !(n % 5)) {
3242:     PetscDrawLGDraw(lg);
3243:   }
3244:   prev = rnorm;
3245:   return(0);
3246: }

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

3253:    Collective on SNES

3255:    Input Parameters:
3256: +  snes - nonlinear solver context obtained from SNESCreate()
3257: .  iter - iteration number
3258: -  rnorm - relative norm of the residual

3260:    Notes:
3261:    This routine is called by the SNES implementations.
3262:    It does not typically need to be called by the user.

3264:    Level: developer

3266: .seealso: SNESMonitorSet()
3267: @*/
3268: PetscErrorCode  SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3269: {
3271:   PetscInt       i,n = snes->numbermonitors;

3274:   VecLockPush(snes->vec_sol);
3275:   for (i=0; i<n; i++) {
3276:     (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3277:   }
3278:   VecLockPop(snes->vec_sol);
3279:   return(0);
3280: }

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

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

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

3291: +    snes - the SNES context
3292: .    its - iteration number
3293: .    norm - 2-norm function value (may be estimated)
3294: -    mctx - [optional] monitoring context

3296:    Level: advanced

3298: .seealso:   SNESMonitorSet(), SNESMonitorGet()
3299: M*/

3303: /*@C
3304:    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3305:    iteration of the nonlinear solver to display the iteration's
3306:    progress.

3308:    Logically Collective on SNES

3310:    Input Parameters:
3311: +  snes - the SNES context
3312: .  f - the monitor function, see SNESMonitorFunction for the calling sequence
3313: .  mctx - [optional] user-defined context for private data for the
3314:           monitor routine (use NULL if no context is desired)
3315: -  monitordestroy - [optional] routine that frees monitor context
3316:           (may be NULL)

3318:    Options Database Keys:
3319: +    -snes_monitor        - sets SNESMonitorDefault()
3320: .    -snes_monitor_lg_residualnorm    - sets line graph monitor,
3321:                             uses SNESMonitorLGCreate()
3322: -    -snes_monitor_cancel - cancels all monitors that have
3323:                             been hardwired into a code by
3324:                             calls to SNESMonitorSet(), but
3325:                             does not cancel those set via
3326:                             the options database.

3328:    Notes:
3329:    Several different monitoring routines may be set by calling
3330:    SNESMonitorSet() multiple times; all will be called in the
3331:    order in which they were set.

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

3335:    Level: intermediate

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

3339: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3340: @*/
3341: PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3342: {
3343:   PetscInt       i;

3348:   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3349:   for (i=0; i<snes->numbermonitors;i++) {
3350:     if (f == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) {
3351:       if (monitordestroy) {
3352:         (*monitordestroy)(&mctx);
3353:       }
3354:       return(0);
3355:     }
3356:   }
3357:   snes->monitor[snes->numbermonitors]          = f;
3358:   snes->monitordestroy[snes->numbermonitors]   = monitordestroy;
3359:   snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3360:   return(0);
3361: }

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

3368:    Logically Collective on SNES

3370:    Input Parameters:
3371: .  snes - the SNES context

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

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

3381:    Level: intermediate

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

3385: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3386: @*/
3387: PetscErrorCode  SNESMonitorCancel(SNES snes)
3388: {
3390:   PetscInt       i;

3394:   for (i=0; i<snes->numbermonitors; i++) {
3395:     if (snes->monitordestroy[i]) {
3396:       (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3397:     }
3398:   }
3399:   snes->numbermonitors = 0;
3400:   return(0);
3401: }

3403: /*MC
3404:     SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver

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

3410: +    snes - the SNES context
3411: .    it - current iteration (0 is the first and is before any Newton step)
3412: .    cctx - [optional] convergence context
3413: .    reason - reason for convergence/divergence
3414: .    xnorm - 2-norm of current iterate
3415: .    gnorm - 2-norm of current step
3416: -    f - 2-norm of function

3418:    Level: intermediate

3420: .seealso:   SNESSetConvergenceTest(), SNESGetConvergenceTest()
3421: M*/

3425: /*@C
3426:    SNESSetConvergenceTest - Sets the function that is to be used
3427:    to test for convergence of the nonlinear iterative solution.

3429:    Logically Collective on SNES

3431:    Input Parameters:
3432: +  snes - the SNES context
3433: .  SNESConvergenceTestFunction - routine to test for convergence
3434: .  cctx - [optional] context for private data for the convergence routine  (may be NULL)
3435: -  destroy - [optional] destructor for the context (may be NULL; NULL_FUNCTION in Fortran)

3437:    Level: advanced

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

3441: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3442: @*/
3443: PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3444: {

3449:   if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
3450:   if (snes->ops->convergeddestroy) {
3451:     (*snes->ops->convergeddestroy)(snes->cnvP);
3452:   }
3453:   snes->ops->converged        = SNESConvergenceTestFunction;
3454:   snes->ops->convergeddestroy = destroy;
3455:   snes->cnvP                  = cctx;
3456:   return(0);
3457: }

3461: /*@
3462:    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.

3464:    Not Collective

3466:    Input Parameter:
3467: .  snes - the SNES context

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

3473:    Level: intermediate

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

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

3479: .seealso: SNESSetConvergenceTest(), SNESConvergedReason
3480: @*/
3481: PetscErrorCode  SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
3482: {
3486:   *reason = snes->reason;
3487:   return(0);
3488: }

3492: /*@
3493:    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.

3495:    Logically Collective on SNES

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

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

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

3513:    Level: intermediate

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

3517: .seealso: SNESGetConvergenceHistory()

3519: @*/
3520: PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
3521: {

3528:   if (!a) {
3529:     if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
3530:     PetscCalloc1(na,&a);
3531:     PetscCalloc1(na,&its);

3533:     snes->conv_malloc = PETSC_TRUE;
3534:   }
3535:   snes->conv_hist       = a;
3536:   snes->conv_hist_its   = its;
3537:   snes->conv_hist_max   = na;
3538:   snes->conv_hist_len   = 0;
3539:   snes->conv_hist_reset = reset;
3540:   return(0);
3541: }

3543: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3544: #include <engine.h>   /* MATLAB include file */
3545: #include <mex.h>      /* MATLAB include file */

3549: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
3550: {
3551:   mxArray   *mat;
3552:   PetscInt  i;
3553:   PetscReal *ar;

3556:   mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
3557:   ar  = (PetscReal*) mxGetData(mat);
3558:   for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
3559:   PetscFunctionReturn(mat);
3560: }
3561: #endif

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

3568:    Not Collective

3570:    Input Parameter:
3571: .  snes - iterative context obtained from SNESCreate()

3573:    Output Parameters:
3574: .  a   - array to hold history
3575: .  its - integer array holds the number of linear iterations (or
3576:          negative if not converged) for each solve.
3577: -  na  - size of a and its

3579:    Notes:
3580:     The calling sequence for this routine in Fortran is
3581: $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)

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

3587:    Level: intermediate

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

3591: .seealso: SNESSetConvergencHistory()

3593: @*/
3594: PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
3595: {
3598:   if (a)   *a   = snes->conv_hist;
3599:   if (its) *its = snes->conv_hist_its;
3600:   if (na)  *na  = snes->conv_hist_len;
3601:   return(0);
3602: }

3606: /*@C
3607:   SNESSetUpdate - Sets the general-purpose update function called
3608:   at the beginning of every iteration of the nonlinear solve. Specifically
3609:   it is called just before the Jacobian is "evaluated".

3611:   Logically Collective on SNES

3613:   Input Parameters:
3614: . snes - The nonlinear solver context
3615: . func - The function

3617:   Calling sequence of func:
3618: . func (SNES snes, PetscInt step);

3620: . step - The current step of the iteration

3622:   Level: advanced

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

3627: .keywords: SNES, update

3629: .seealso SNESSetJacobian(), SNESSolve()
3630: @*/
3631: PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
3632: {
3635:   snes->ops->update = func;
3636:   return(0);
3637: }

3641: /*
3642:    SNESScaleStep_Private - Scales a step so that its length is less than the
3643:    positive parameter delta.

3645:     Input Parameters:
3646: +   snes - the SNES context
3647: .   y - approximate solution of linear system
3648: .   fnorm - 2-norm of current function
3649: -   delta - trust region size

3651:     Output Parameters:
3652: +   gpnorm - predicted function norm at the new point, assuming local
3653:     linearization.  The value is zero if the step lies within the trust
3654:     region, and exceeds zero otherwise.
3655: -   ynorm - 2-norm of the step

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

3661: .keywords: SNES, nonlinear, scale, step
3662: */
3663: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
3664: {
3665:   PetscReal      nrm;
3666:   PetscScalar    cnorm;


3674:   VecNorm(y,NORM_2,&nrm);
3675:   if (nrm > *delta) {
3676:     nrm     = *delta/nrm;
3677:     *gpnorm = (1.0 - nrm)*(*fnorm);
3678:     cnorm   = nrm;
3679:     VecScale(y,cnorm);
3680:     *ynorm  = *delta;
3681:   } else {
3682:     *gpnorm = 0.0;
3683:     *ynorm  = nrm;
3684:   }
3685:   return(0);
3686: }

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

3693:    Collective on SNES

3695:    Parameter:
3696: +  snes - iterative context obtained from SNESCreate()
3697: -  viewer - the viewer to display the reason


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

3703:    Level: beginner

3705: .keywords: SNES, solve, linear system

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

3709: @*/
3710: PetscErrorCode  SNESReasonView(SNES snes,PetscViewer viewer)
3711: {
3713:   PetscBool      isAscii;

3716:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
3717:   if (isAscii) {
3718:     PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
3719:     if (snes->reason > 0) {
3720:       if (((PetscObject) snes)->prefix) {
3721:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3722:       } else {
3723:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3724:       }
3725:     } else {
3726:       if (((PetscObject) snes)->prefix) {
3727:         PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3728:       } else {
3729:         PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3730:       }
3731:     }
3732:     PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
3733:   }
3734:   return(0);
3735: }

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

3742:   Collective on SNES

3744:   Input Parameters:
3745: . snes   - the SNES object

3747:   Level: intermediate

3749: @*/
3750: PetscErrorCode SNESReasonViewFromOptions(SNES snes)
3751: {
3752:   PetscErrorCode    ierr;
3753:   PetscViewer       viewer;
3754:   PetscBool         flg;
3755:   static PetscBool  incall = PETSC_FALSE;
3756:   PetscViewerFormat format;

3759:   if (incall) return(0);
3760:   incall = PETSC_TRUE;
3761:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
3762:   if (flg) {
3763:     PetscViewerPushFormat(viewer,format);
3764:     SNESReasonView(snes,viewer);
3765:     PetscViewerPopFormat(viewer);
3766:     PetscViewerDestroy(&viewer);
3767:   }
3768:   incall = PETSC_FALSE;
3769:   return(0);
3770: }

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

3778:    Collective on SNES

3780:    Input Parameters:
3781: +  snes - the SNES context
3782: .  b - the constant part of the equation F(x) = b, or NULL to use zero.
3783: -  x - the solution vector.

3785:    Notes:
3786:    The user should initialize the vector,x, with the initial guess
3787:    for the nonlinear solve prior to calling SNESSolve.  In particular,
3788:    to employ an initial guess of zero, the user should explicitly set
3789:    this vector to zero by calling VecSet().

3791:    Level: beginner

3793: .keywords: SNES, nonlinear, solve

3795: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
3796: @*/
3797: PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
3798: {
3799:   PetscErrorCode    ierr;
3800:   PetscBool         flg;
3801:   PetscInt          grid;
3802:   Vec               xcreated = NULL;
3803:   DM                dm;


3812:   if (!x) {
3813:     SNESGetDM(snes,&dm);
3814:     DMCreateGlobalVector(dm,&xcreated);
3815:     x    = xcreated;
3816:   }
3817:   SNESViewFromOptions(snes,NULL,"-snes_view_pre");

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

3822:     /* set solution vector */
3823:     if (!grid) {PetscObjectReference((PetscObject)x);}
3824:     VecDestroy(&snes->vec_sol);
3825:     snes->vec_sol = x;
3826:     SNESGetDM(snes,&dm);

3828:     /* set affine vector if provided */
3829:     if (b) { PetscObjectReference((PetscObject)b); }
3830:     VecDestroy(&snes->vec_rhs);
3831:     snes->vec_rhs = b;

3833:     if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
3834:     if (snes->vec_rhs  == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
3835:     if (!snes->vec_sol_update /* && snes->vec_sol */) {
3836:       VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
3837:       PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
3838:     }
3839:     DMShellSetGlobalVector(dm,snes->vec_sol);
3840:     SNESSetUp(snes);

3842:     if (!grid) {
3843:       if (snes->ops->computeinitialguess) {
3844:         (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
3845:       }
3846:     }

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

3851:     PetscLogEventBegin(SNES_Solve,snes,0,0,0);
3852:     (*snes->ops->solve)(snes);
3853:     PetscLogEventEnd(SNES_Solve,snes,0,0,0);
3854:     if (snes->domainerror) {
3855:       snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
3856:       snes->domainerror = PETSC_FALSE;
3857:     }
3858:     if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");

3860:     if (snes->lagjac_persist) snes->jac_iter += snes->iter;
3861:     if (snes->lagpre_persist) snes->pre_iter += snes->iter;

3863:     flg  = PETSC_FALSE;
3864:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,NULL);
3865:     if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
3866:     SNESReasonViewFromOptions(snes);

3868:     if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
3869:     if (grid <  snes->gridsequence) {
3870:       DM  fine;
3871:       Vec xnew;
3872:       Mat interp;

3874:       DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
3875:       if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
3876:       DMCreateInterpolation(snes->dm,fine,&interp,NULL);
3877:       DMCreateGlobalVector(fine,&xnew);
3878:       MatInterpolate(interp,x,xnew);
3879:       DMInterpolate(snes->dm,interp,fine);
3880:       MatDestroy(&interp);
3881:       x    = xnew;

3883:       SNESReset(snes);
3884:       SNESSetDM(snes,fine);
3885:       DMDestroy(&fine);
3886:       PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
3887:     }
3888:   }
3889:   SNESViewFromOptions(snes,NULL,"-snes_view");
3890:   VecViewFromOptions(snes->vec_sol,((PetscObject)snes)->prefix,"-snes_view_solution");

3892:   VecDestroy(&xcreated);
3893:   PetscObjectSAWsBlock((PetscObject)snes);
3894:   return(0);
3895: }

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

3901: /*@C
3902:    SNESSetType - Sets the method for the nonlinear solver.

3904:    Collective on SNES

3906:    Input Parameters:
3907: +  snes - the SNES context
3908: -  type - a known method

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

3914:    Notes:
3915:    See "petsc/include/petscsnes.h" for available methods (for instance)
3916: +    SNESNEWTONLS - Newton's method with line search
3917:      (systems of nonlinear equations)
3918: .    SNESNEWTONTR - Newton's method with trust region
3919:      (systems of nonlinear equations)

3921:   Normally, it is best to use the SNESSetFromOptions() command and then
3922:   set the SNES solver type from the options database rather than by using
3923:   this routine.  Using the options database provides the user with
3924:   maximum flexibility in evaluating the many nonlinear solvers.
3925:   The SNESSetType() routine is provided for those situations where it
3926:   is necessary to set the nonlinear solver independently of the command
3927:   line or options database.  This might be the case, for example, when
3928:   the choice of solver changes during the execution of the program,
3929:   and the user's application is taking responsibility for choosing the
3930:   appropriate method.

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

3935:   Level: intermediate

3937: .keywords: SNES, set, type

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

3941: @*/
3942: PetscErrorCode  SNESSetType(SNES snes,SNESType type)
3943: {
3944:   PetscErrorCode ierr,(*r)(SNES);
3945:   PetscBool      match;


3951:   PetscObjectTypeCompare((PetscObject)snes,type,&match);
3952:   if (match) return(0);

3954:    PetscFunctionListFind(SNESList,type,&r);
3955:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
3956:   /* Destroy the previous private SNES context */
3957:   if (snes->ops->destroy) {
3958:     (*(snes)->ops->destroy)(snes);
3959:     snes->ops->destroy = NULL;
3960:   }
3961:   /* Reinitialize function pointers in SNESOps structure */
3962:   snes->ops->setup          = 0;
3963:   snes->ops->solve          = 0;
3964:   snes->ops->view           = 0;
3965:   snes->ops->setfromoptions = 0;
3966:   snes->ops->destroy        = 0;
3967:   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
3968:   snes->setupcalled = PETSC_FALSE;

3970:   PetscObjectChangeTypeName((PetscObject)snes,type);
3971:   (*r)(snes);
3972:   return(0);
3973: }

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

3980:    Not Collective

3982:    Input Parameter:
3983: .  snes - nonlinear solver context

3985:    Output Parameter:
3986: .  type - SNES method (a character string)

3988:    Level: intermediate

3990: .keywords: SNES, nonlinear, get, type, name
3991: @*/
3992: PetscErrorCode  SNESGetType(SNES snes,SNESType *type)
3993: {
3997:   *type = ((PetscObject)snes)->type_name;
3998:   return(0);
3999: }

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

4006:   Logically Collective on SNES and Vec

4008:   Input Parameters:
4009: + snes - the SNES context obtained from SNESCreate()
4010: - u    - the solution vector

4012:   Level: beginner

4014: .keywords: SNES, set, solution
4015: @*/
4016: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4017: {
4018:   DM             dm;

4024:   PetscObjectReference((PetscObject) u);
4025:   VecDestroy(&snes->vec_sol);

4027:   snes->vec_sol = u;

4029:   SNESGetDM(snes, &dm);
4030:   DMShellSetGlobalVector(dm, u);
4031:   return(0);
4032: }

4036: /*@
4037:    SNESGetSolution - Returns the vector where the approximate solution is
4038:    stored. This is the fine grid solution when using SNESSetGridSequence().

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

4042:    Input Parameter:
4043: .  snes - the SNES context

4045:    Output Parameter:
4046: .  x - the solution

4048:    Level: intermediate

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

4052: .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
4053: @*/
4054: PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
4055: {
4059:   *x = snes->vec_sol;
4060:   return(0);
4061: }

4065: /*@
4066:    SNESGetSolutionUpdate - Returns the vector where the solution update is
4067:    stored.

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

4071:    Input Parameter:
4072: .  snes - the SNES context

4074:    Output Parameter:
4075: .  x - the solution update

4077:    Level: advanced

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

4081: .seealso: SNESGetSolution(), SNESGetFunction()
4082: @*/
4083: PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
4084: {
4088:   *x = snes->vec_sol_update;
4089:   return(0);
4090: }

4094: /*@C
4095:    SNESGetFunction - Returns the vector where the function is stored.

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

4099:    Input Parameter:
4100: .  snes - the SNES context

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

4107:    Level: advanced

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

4111: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4112: @*/
4113: PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4114: {
4116:   DM             dm;

4120:   if (r) {
4121:     if (!snes->vec_func) {
4122:       if (snes->vec_rhs) {
4123:         VecDuplicate(snes->vec_rhs,&snes->vec_func);
4124:       } else if (snes->vec_sol) {
4125:         VecDuplicate(snes->vec_sol,&snes->vec_func);
4126:       } else if (snes->dm) {
4127:         DMCreateGlobalVector(snes->dm,&snes->vec_func);
4128:       }
4129:     }
4130:     *r = snes->vec_func;
4131:   }
4132:   SNESGetDM(snes,&dm);
4133:   DMSNESGetFunction(dm,f,ctx);
4134:   return(0);
4135: }

4137: /*@C
4138:    SNESGetNGS - Returns the NGS function and context.

4140:    Input Parameter:
4141: .  snes - the SNES context

4143:    Output Parameter:
4144: +  f - the function (or NULL) see SNESNGSFunction for details
4145: -  ctx    - the function context (or NULL)

4147:    Level: advanced

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

4151: .seealso: SNESSetNGS(), SNESGetFunction()
4152: @*/

4156: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4157: {
4159:   DM             dm;

4163:   SNESGetDM(snes,&dm);
4164:   DMSNESGetNGS(dm,f,ctx);
4165:   return(0);
4166: }

4170: /*@C
4171:    SNESSetOptionsPrefix - Sets the prefix used for searching for all
4172:    SNES options in the database.

4174:    Logically Collective on SNES

4176:    Input Parameter:
4177: +  snes - the SNES context
4178: -  prefix - the prefix to prepend to all option names

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

4184:    Level: advanced

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

4188: .seealso: SNESSetFromOptions()
4189: @*/
4190: PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
4191: {

4196:   PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4197:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4198:   if (snes->linesearch) {
4199:     SNESGetLineSearch(snes,&snes->linesearch);
4200:     PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4201:   }
4202:   KSPSetOptionsPrefix(snes->ksp,prefix);
4203:   return(0);
4204: }

4208: /*@C
4209:    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4210:    SNES options in the database.

4212:    Logically Collective on SNES

4214:    Input Parameters:
4215: +  snes - the SNES context
4216: -  prefix - the prefix to prepend to all option names

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

4222:    Level: advanced

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

4226: .seealso: SNESGetOptionsPrefix()
4227: @*/
4228: PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4229: {

4234:   PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4235:   if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4236:   if (snes->linesearch) {
4237:     SNESGetLineSearch(snes,&snes->linesearch);
4238:     PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4239:   }
4240:   KSPAppendOptionsPrefix(snes->ksp,prefix);
4241:   return(0);
4242: }

4246: /*@C
4247:    SNESGetOptionsPrefix - Sets the prefix used for searching for all
4248:    SNES options in the database.

4250:    Not Collective

4252:    Input Parameter:
4253: .  snes - the SNES context

4255:    Output Parameter:
4256: .  prefix - pointer to the prefix string used

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

4261:    Level: advanced

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

4265: .seealso: SNESAppendOptionsPrefix()
4266: @*/
4267: PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4268: {

4273:   PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4274:   return(0);
4275: }


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

4283:    Not collective

4285:    Input Parameters:
4286: +  name_solver - name of a new user-defined solver
4287: -  routine_create - routine to create method context

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

4292:    Sample usage:
4293: .vb
4294:    SNESRegister("my_solver",MySolverCreate);
4295: .ve

4297:    Then, your solver can be chosen with the procedural interface via
4298: $     SNESSetType(snes,"my_solver")
4299:    or at runtime via the option
4300: $     -snes_type my_solver

4302:    Level: advanced

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

4306: .keywords: SNES, nonlinear, register

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

4310:   Level: advanced
4311: @*/
4312: PetscErrorCode  SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4313: {

4317:   PetscFunctionListAdd(&SNESList,sname,function);
4318:   return(0);
4319: }

4323: PetscErrorCode  SNESTestLocalMin(SNES snes)
4324: {
4326:   PetscInt       N,i,j;
4327:   Vec            u,uh,fh;
4328:   PetscScalar    value;
4329:   PetscReal      norm;

4332:   SNESGetSolution(snes,&u);
4333:   VecDuplicate(u,&uh);
4334:   VecDuplicate(u,&fh);

4336:   /* currently only works for sequential */
4337:   PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4338:   VecGetSize(u,&N);
4339:   for (i=0; i<N; i++) {
4340:     VecCopy(u,uh);
4341:     PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4342:     for (j=-10; j<11; j++) {
4343:       value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4344:       VecSetValue(uh,i,value,ADD_VALUES);
4345:       SNESComputeFunction(snes,uh,fh);
4346:       VecNorm(fh,NORM_2,&norm);
4347:       PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);
4348:       value = -value;
4349:       VecSetValue(uh,i,value,ADD_VALUES);
4350:     }
4351:   }
4352:   VecDestroy(&uh);
4353:   VecDestroy(&fh);
4354:   return(0);
4355: }

4359: /*@
4360:    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4361:    computing relative tolerance for linear solvers within an inexact
4362:    Newton method.

4364:    Logically Collective on SNES

4366:    Input Parameters:
4367: +  snes - SNES context
4368: -  flag - PETSC_TRUE or PETSC_FALSE

4370:     Options Database:
4371: +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4372: .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
4373: .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4374: .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4375: .  -snes_ksp_ew_gamma <gamma> - Sets gamma
4376: .  -snes_ksp_ew_alpha <alpha> - Sets alpha
4377: .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4378: -  -snes_ksp_ew_threshold <threshold> - Sets threshold

4380:    Notes:
4381:    Currently, the default is to use a constant relative tolerance for
4382:    the inner linear solvers.  Alternatively, one can use the
4383:    Eisenstat-Walker method, where the relative convergence tolerance
4384:    is reset at each Newton iteration according progress of the nonlinear
4385:    solver.

4387:    Level: advanced

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

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

4395: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4396: @*/
4397: PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool flag)
4398: {
4402:   snes->ksp_ewconv = flag;
4403:   return(0);
4404: }

4408: /*@
4409:    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
4410:    for computing relative tolerance for linear solvers within an
4411:    inexact Newton method.

4413:    Not Collective

4415:    Input Parameter:
4416: .  snes - SNES context

4418:    Output Parameter:
4419: .  flag - PETSC_TRUE or PETSC_FALSE

4421:    Notes:
4422:    Currently, the default is to use a constant relative tolerance for
4423:    the inner linear solvers.  Alternatively, one can use the
4424:    Eisenstat-Walker method, where the relative convergence tolerance
4425:    is reset at each Newton iteration according progress of the nonlinear
4426:    solver.

4428:    Level: advanced

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

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

4436: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4437: @*/
4438: PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
4439: {
4443:   *flag = snes->ksp_ewconv;
4444:   return(0);
4445: }

4449: /*@
4450:    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
4451:    convergence criteria for the linear solvers within an inexact
4452:    Newton method.

4454:    Logically Collective on SNES

4456:    Input Parameters:
4457: +    snes - SNES context
4458: .    version - version 1, 2 (default is 2) or 3
4459: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4460: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4461: .    gamma - multiplicative factor for version 2 rtol computation
4462:              (0 <= gamma2 <= 1)
4463: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
4464: .    alpha2 - power for safeguard
4465: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

4467:    Note:
4468:    Version 3 was contributed by Luis Chacon, June 2006.

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

4472:    Level: advanced

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

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

4481: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
4482: @*/
4483: PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
4484: {
4485:   SNESKSPEW *kctx;

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

4499:   if (version != PETSC_DEFAULT)   kctx->version   = version;
4500:   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
4501:   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
4502:   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
4503:   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
4504:   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
4505:   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;

4507:   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);
4508:   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);
4509:   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);
4510:   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);
4511:   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);
4512:   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);
4513:   return(0);
4514: }

4518: /*@
4519:    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
4520:    convergence criteria for the linear solvers within an inexact
4521:    Newton method.

4523:    Not Collective

4525:    Input Parameters:
4526:      snes - SNES context

4528:    Output Parameters:
4529: +    version - version 1, 2 (default is 2) or 3
4530: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4531: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4532: .    gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
4533: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
4534: .    alpha2 - power for safeguard
4535: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

4537:    Level: advanced

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

4541: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
4542: @*/
4543: PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
4544: {
4545:   SNESKSPEW *kctx;

4549:   kctx = (SNESKSPEW*)snes->kspconvctx;
4550:   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4551:   if (version)   *version   = kctx->version;
4552:   if (rtol_0)    *rtol_0    = kctx->rtol_0;
4553:   if (rtol_max)  *rtol_max  = kctx->rtol_max;
4554:   if (gamma)     *gamma     = kctx->gamma;
4555:   if (alpha)     *alpha     = kctx->alpha;
4556:   if (alpha2)    *alpha2    = kctx->alpha2;
4557:   if (threshold) *threshold = kctx->threshold;
4558:   return(0);
4559: }

4563:  PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4564: {
4566:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
4567:   PetscReal      rtol  = PETSC_DEFAULT,stol;

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

4606: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4607: {
4609:   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
4610:   PCSide         pcside;
4611:   Vec            lres;

4614:   if (!snes->ksp_ewconv) return(0);
4615:   KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
4616:   kctx->norm_last = snes->norm;
4617:   if (kctx->version == 1) {
4618:     PC        pc;
4619:     PetscBool isNone;

4621:     KSPGetPC(ksp, &pc);
4622:     PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
4623:     KSPGetPCSide(ksp,&pcside);
4624:      if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
4625:       /* KSP residual is true linear residual */
4626:       KSPGetResidualNorm(ksp,&kctx->lresid_last);
4627:     } else {
4628:       /* KSP residual is preconditioned residual */
4629:       /* compute true linear residual norm */
4630:       VecDuplicate(b,&lres);
4631:       MatMult(snes->jacobian,x,lres);
4632:       VecAYPX(lres,-1.0,b);
4633:       VecNorm(lres,NORM_2,&kctx->lresid_last);
4634:       VecDestroy(&lres);
4635:     }
4636:   }
4637:   return(0);
4638: }

4642: /*@
4643:    SNESGetKSP - Returns the KSP context for a SNES solver.

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

4647:    Input Parameter:
4648: .  snes - the SNES context

4650:    Output Parameter:
4651: .  ksp - the KSP context

4653:    Notes:
4654:    The user can then directly manipulate the KSP context to set various
4655:    options, etc.  Likewise, the user can then extract and manipulate the
4656:    PC contexts as well.

4658:    Level: beginner

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

4662: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
4663: @*/
4664: PetscErrorCode  SNESGetKSP(SNES snes,KSP *ksp)
4665: {


4672:   if (!snes->ksp) {
4673:     PetscBool monitor = PETSC_FALSE;

4675:     KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
4676:     PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
4677:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);

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

4682:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
4683:     if (monitor) {
4684:       KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
4685:     }
4686:     monitor = PETSC_FALSE;
4687:     PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
4688:     if (monitor) {
4689:       PetscObject *objs;
4690:       KSPMonitorSNESLGResidualNormCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
4691:       objs[0] = (PetscObject) snes;
4692:       KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
4693:     }
4694:   }
4695:   *ksp = snes->ksp;
4696:   return(0);
4697: }


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

4706:    Logically Collective on SNES

4708:    Input Parameters:
4709: +  snes - the preconditioner context
4710: -  dm - the dm

4712:    Level: intermediate

4714: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
4715: @*/
4716: PetscErrorCode  SNESSetDM(SNES snes,DM dm)
4717: {
4719:   KSP            ksp;
4720:   DMSNES         sdm;

4724:   if (dm) {PetscObjectReference((PetscObject)dm);}
4725:   if (snes->dm) {               /* Move the DMSNES context over to the new DM unless the new DM already has one */
4726:     if (snes->dm->dmsnes && snes->dmAuto && !dm->dmsnes) {
4727:       DMCopyDMSNES(snes->dm,dm);
4728:       DMGetDMSNES(snes->dm,&sdm);
4729:       if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
4730:     }
4731:     DMDestroy(&snes->dm);
4732:   }
4733:   snes->dm     = dm;
4734:   snes->dmAuto = PETSC_FALSE;

4736:   SNESGetKSP(snes,&ksp);
4737:   KSPSetDM(ksp,dm);
4738:   KSPSetDMActive(ksp,PETSC_FALSE);
4739:   if (snes->pc) {
4740:     SNESSetDM(snes->pc, snes->dm);
4741:     SNESSetNPCSide(snes,snes->pcside);
4742:   }
4743:   return(0);
4744: }

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

4751:    Not Collective but DM obtained is parallel on SNES

4753:    Input Parameter:
4754: . snes - the preconditioner context

4756:    Output Parameter:
4757: .  dm - the dm

4759:    Level: intermediate

4761: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
4762: @*/
4763: PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
4764: {

4769:   if (!snes->dm) {
4770:     DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
4771:     snes->dmAuto = PETSC_TRUE;
4772:   }
4773:   *dm = snes->dm;
4774:   return(0);
4775: }

4779: /*@
4780:   SNESSetNPC - Sets the nonlinear preconditioner to be used.

4782:   Collective on SNES

4784:   Input Parameters:
4785: + snes - iterative context obtained from SNESCreate()
4786: - pc   - the preconditioner object

4788:   Notes:
4789:   Use SNESGetNPC() to retrieve the preconditioner context (for example,
4790:   to configure it using the API).

4792:   Level: developer

4794: .keywords: SNES, set, precondition
4795: .seealso: SNESGetNPC(), SNESHasNPC()
4796: @*/
4797: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
4798: {

4805:   PetscObjectReference((PetscObject) pc);
4806:   SNESDestroy(&snes->pc);
4807:   snes->pc = pc;
4808:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->pc);
4809:   return(0);
4810: }

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

4817:   Not Collective

4819:   Input Parameter:
4820: . snes - iterative context obtained from SNESCreate()

4822:   Output Parameter:
4823: . pc - preconditioner context

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

4827:   Level: developer

4829: .keywords: SNES, get, preconditioner
4830: .seealso: SNESSetNPC(), SNESHasNPC()
4831: @*/
4832: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
4833: {
4835:   const char     *optionsprefix;

4840:   if (!snes->pc) {
4841:     SNESCreate(PetscObjectComm((PetscObject)snes),&snes->pc);
4842:     PetscObjectIncrementTabLevel((PetscObject)snes->pc,(PetscObject)snes,1);
4843:     PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->pc);
4844:     SNESGetOptionsPrefix(snes,&optionsprefix);
4845:     SNESSetOptionsPrefix(snes->pc,optionsprefix);
4846:     SNESAppendOptionsPrefix(snes->pc,"npc_");
4847:     SNESSetCountersReset(snes->pc,PETSC_FALSE);
4848:   }
4849:   *pc = snes->pc;
4850:   return(0);
4851: }

4855: /*@
4856:   SNESHasNPC - Returns whether a nonlinear preconditioner exists

4858:   Not Collective

4860:   Input Parameter:
4861: . snes - iterative context obtained from SNESCreate()

4863:   Output Parameter:
4864: . has_npc - whether the SNES has an NPC or not

4866:   Level: developer

4868: .keywords: SNES, has, preconditioner
4869: .seealso: SNESSetNPC(), SNESGetNPC()
4870: @*/
4871: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
4872: {
4875:   *has_npc = (PetscBool) (snes->pc != NULL);
4876:   return(0);
4877: }

4881: /*@
4882:     SNESSetNPCSide - Sets the preconditioning side.

4884:     Logically Collective on SNES

4886:     Input Parameter:
4887: .   snes - iterative context obtained from SNESCreate()

4889:     Output Parameter:
4890: .   side - the preconditioning side, where side is one of
4891: .vb
4892:       PC_LEFT - left preconditioning (default)
4893:       PC_RIGHT - right preconditioning
4894: .ve

4896:     Options Database Keys:
4897: .   -snes_pc_side <right,left>

4899:     Level: intermediate

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

4903: .seealso: SNESGetNPCSide(), KSPSetPCSide()
4904: @*/
4905: PetscErrorCode  SNESSetNPCSide(SNES snes,PCSide side)
4906: {
4910:   snes->pcside = side;
4911:   return(0);
4912: }

4916: /*@
4917:     SNESGetNPCSide - Gets the preconditioning side.

4919:     Not Collective

4921:     Input Parameter:
4922: .   snes - iterative context obtained from SNESCreate()

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

4931:     Level: intermediate

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

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

4948: /*@
4949:   SNESSetLineSearch - Sets the linesearch on the SNES instance.

4951:   Collective on SNES

4953:   Input Parameters:
4954: + snes - iterative context obtained from SNESCreate()
4955: - linesearch   - the linesearch object

4957:   Notes:
4958:   Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
4959:   to configure it using the API).

4961:   Level: developer

4963: .keywords: SNES, set, linesearch
4964: .seealso: SNESGetLineSearch()
4965: @*/
4966: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
4967: {

4974:   PetscObjectReference((PetscObject) linesearch);
4975:   SNESLineSearchDestroy(&snes->linesearch);

4977:   snes->linesearch = linesearch;

4979:   PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
4980:   return(0);
4981: }

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

4989:   Not Collective

4991:   Input Parameter:
4992: . snes - iterative context obtained from SNESCreate()

4994:   Output Parameter:
4995: . linesearch - linesearch context

4997:   Level: beginner

4999: .keywords: SNES, get, linesearch
5000: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5001: @*/
5002: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5003: {
5005:   const char     *optionsprefix;

5010:   if (!snes->linesearch) {
5011:     SNESGetOptionsPrefix(snes, &optionsprefix);
5012:     SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5013:     SNESLineSearchSetSNES(snes->linesearch, snes);
5014:     SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5015:     PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5016:     PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5017:   }
5018:   *linesearch = snes->linesearch;
5019:   return(0);
5020: }

5022: #if defined(PETSC_HAVE_MATLAB_ENGINE)
5023: #include <mex.h>

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

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

5032:    Collective on SNES

5034:    Input Parameters:
5035: +  snes - the SNES context
5036: -  x - input vector

5038:    Output Parameter:
5039: .  y - function vector, as set by SNESSetFunction()

5041:    Notes:
5042:    SNESComputeFunction() is typically used within nonlinear solvers
5043:    implementations, so most users would not generally call this routine
5044:    themselves.

5046:    Level: developer

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

5050: .seealso: SNESSetFunction(), SNESGetFunction()
5051: */
5052: PetscErrorCode  SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
5053: {
5054:   PetscErrorCode    ierr;
5055:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5056:   int               nlhs  = 1,nrhs = 5;
5057:   mxArray           *plhs[1],*prhs[5];
5058:   long long int     lx = 0,ly = 0,ls = 0;


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

5069:   PetscMemcpy(&ls,&snes,sizeof(snes));
5070:   PetscMemcpy(&lx,&x,sizeof(x));
5071:   PetscMemcpy(&ly,&y,sizeof(x));
5072:   prhs[0] = mxCreateDoubleScalar((double)ls);
5073:   prhs[1] = mxCreateDoubleScalar((double)lx);
5074:   prhs[2] = mxCreateDoubleScalar((double)ly);
5075:   prhs[3] = mxCreateString(sctx->funcname);
5076:   prhs[4] = sctx->ctx;
5077:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");
5078:   mxGetScalar(plhs[0]);
5079:   mxDestroyArray(prhs[0]);
5080:   mxDestroyArray(prhs[1]);
5081:   mxDestroyArray(prhs[2]);
5082:   mxDestroyArray(prhs[3]);
5083:   mxDestroyArray(plhs[0]);
5084:   return(0);
5085: }

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

5094:    Logically Collective on SNES

5096:    Input Parameters:
5097: +  snes - the SNES context
5098: .  r - vector to store function value
5099: -  f - function evaluation routine

5101:    Notes:
5102:    The Newton-like methods typically solve linear systems of the form
5103: $      f'(x) x = -f(x),
5104:    where f'(x) denotes the Jacobian matrix and f(x) is the function.

5106:    Level: beginner

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

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

5112: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5113: */
5114: PetscErrorCode  SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx)
5115: {
5116:   PetscErrorCode    ierr;
5117:   SNESMatlabContext *sctx;

5120:   /* currently sctx is memory bleed */
5121:   PetscNew(&sctx);
5122:   PetscStrallocpy(f,&sctx->funcname);
5123:   /*
5124:      This should work, but it doesn't
5125:   sctx->ctx = ctx;
5126:   mexMakeArrayPersistent(sctx->ctx);
5127:   */
5128:   sctx->ctx = mxDuplicateArray(ctx);
5129:   SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);
5130:   return(0);
5131: }

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

5138:    Collective on SNES

5140:    Input Parameters:
5141: +  snes - the SNES context
5142: .  x - input vector
5143: .  A, B - the matrices
5144: -  ctx - user context

5146:    Level: developer

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

5150: .seealso: SNESSetFunction(), SNESGetFunction()
5151: @*/
5152: PetscErrorCode  SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx)
5153: {
5154:   PetscErrorCode    ierr;
5155:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5156:   int               nlhs  = 2,nrhs = 6;
5157:   mxArray           *plhs[2],*prhs[6];
5158:   long long int     lx = 0,lA = 0,ls = 0, lB = 0;


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

5166:   PetscMemcpy(&ls,&snes,sizeof(snes));
5167:   PetscMemcpy(&lx,&x,sizeof(x));
5168:   PetscMemcpy(&lA,A,sizeof(x));
5169:   PetscMemcpy(&lB,B,sizeof(x));
5170:   prhs[0] = mxCreateDoubleScalar((double)ls);
5171:   prhs[1] = mxCreateDoubleScalar((double)lx);
5172:   prhs[2] = mxCreateDoubleScalar((double)lA);
5173:   prhs[3] = mxCreateDoubleScalar((double)lB);
5174:   prhs[4] = mxCreateString(sctx->funcname);
5175:   prhs[5] = sctx->ctx;
5176:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");
5177:   mxGetScalar(plhs[0]);
5178:   mxDestroyArray(prhs[0]);
5179:   mxDestroyArray(prhs[1]);
5180:   mxDestroyArray(prhs[2]);
5181:   mxDestroyArray(prhs[3]);
5182:   mxDestroyArray(prhs[4]);
5183:   mxDestroyArray(plhs[0]);
5184:   mxDestroyArray(plhs[1]);
5185:   return(0);
5186: }

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

5195:    Logically Collective on SNES

5197:    Input Parameters:
5198: +  snes - the SNES context
5199: .  A,B - Jacobian matrices
5200: .  J - function evaluation routine
5201: -  ctx - user context

5203:    Level: developer

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

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

5209: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J
5210: */
5211: PetscErrorCode  SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx)
5212: {
5213:   PetscErrorCode    ierr;
5214:   SNESMatlabContext *sctx;

5217:   /* currently sctx is memory bleed */
5218:   PetscNew(&sctx);
5219:   PetscStrallocpy(J,&sctx->funcname);
5220:   /*
5221:      This should work, but it doesn't
5222:   sctx->ctx = ctx;
5223:   mexMakeArrayPersistent(sctx->ctx);
5224:   */
5225:   sctx->ctx = mxDuplicateArray(ctx);
5226:   SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);
5227:   return(0);
5228: }

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

5235:    Collective on SNES

5237: .seealso: SNESSetFunction(), SNESGetFunction()
5238: @*/
5239: PetscErrorCode  SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
5240: {
5241:   PetscErrorCode    ierr;
5242:   SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5243:   int               nlhs  = 1,nrhs = 6;
5244:   mxArray           *plhs[1],*prhs[6];
5245:   long long int     lx = 0,ls = 0;
5246:   Vec               x  = snes->vec_sol;


5251:   PetscMemcpy(&ls,&snes,sizeof(snes));
5252:   PetscMemcpy(&lx,&x,sizeof(x));
5253:   prhs[0] = mxCreateDoubleScalar((double)ls);
5254:   prhs[1] = mxCreateDoubleScalar((double)it);
5255:   prhs[2] = mxCreateDoubleScalar((double)fnorm);
5256:   prhs[3] = mxCreateDoubleScalar((double)lx);
5257:   prhs[4] = mxCreateString(sctx->funcname);
5258:   prhs[5] = sctx->ctx;
5259:   mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");
5260:   mxGetScalar(plhs[0]);
5261:   mxDestroyArray(prhs[0]);
5262:   mxDestroyArray(prhs[1]);
5263:   mxDestroyArray(prhs[2]);
5264:   mxDestroyArray(prhs[3]);
5265:   mxDestroyArray(prhs[4]);
5266:   mxDestroyArray(plhs[0]);
5267:   return(0);
5268: }

5272: /*
5273:    SNESMonitorSetMatlab - Sets the monitor function from MATLAB

5275:    Level: developer

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

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

5281: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5282: */
5283: PetscErrorCode  SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx)
5284: {
5285:   PetscErrorCode    ierr;
5286:   SNESMatlabContext *sctx;

5289:   /* currently sctx is memory bleed */
5290:   PetscNew(&sctx);
5291:   PetscStrallocpy(f,&sctx->funcname);
5292:   /*
5293:      This should work, but it doesn't
5294:   sctx->ctx = ctx;
5295:   mexMakeArrayPersistent(sctx->ctx);
5296:   */
5297:   sctx->ctx = mxDuplicateArray(ctx);
5298:   SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);
5299:   return(0);
5300: }

5302: #endif