Actual source code: snes.c
1: #include <petsc/private/snesimpl.h>
2: #include <petscdmshell.h>
3: #include <petscdraw.h>
4: #include <petscds.h>
5: #include <petscdmadaptor.h>
6: #include <petscconvest.h>
8: PetscBool SNESRegisterAllCalled = PETSC_FALSE;
9: PetscFunctionList SNESList = NULL;
11: /* Logging support */
12: PetscClassId SNES_CLASSID, DMSNES_CLASSID;
13: PetscLogEvent SNES_Solve, SNES_Setup, SNES_FunctionEval, SNES_JacobianEval, SNES_NGSEval, SNES_NGSFuncEval, SNES_NPCSolve, SNES_ObjectiveEval;
15: /*@
16: SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged.
18: Logically Collective on SNES
20: Input Parameters:
21: + snes - iterative context obtained from SNESCreate()
22: - flg - PETSC_TRUE indicates you want the error generated
24: Options database keys:
25: . -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)
27: Level: intermediate
29: Notes:
30: Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve()
31: to determine if it has converged.
33: .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIfNotConverged()
34: @*/
35: PetscErrorCode SNESSetErrorIfNotConverged(SNES snes,PetscBool flg)
36: {
40: snes->errorifnotconverged = flg;
41: return(0);
42: }
44: /*@
45: SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge?
47: Not Collective
49: Input Parameter:
50: . snes - iterative context obtained from SNESCreate()
52: Output Parameter:
53: . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE
55: Level: intermediate
57: .seealso: SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIfNotConverged()
58: @*/
59: PetscErrorCode SNESGetErrorIfNotConverged(SNES snes,PetscBool *flag)
60: {
64: *flag = snes->errorifnotconverged;
65: return(0);
66: }
68: /*@
69: SNESSetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?
71: Logically Collective on SNES
73: Input Parameters:
74: + snes - the shell SNES
75: - flg - is the residual computed?
77: Level: advanced
79: .seealso: SNESGetAlwaysComputesFinalResidual()
80: @*/
81: PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES snes, PetscBool flg)
82: {
85: snes->alwayscomputesfinalresidual = flg;
86: return(0);
87: }
89: /*@
90: SNESGetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?
92: Logically Collective on SNES
94: Input Parameter:
95: . snes - the shell SNES
97: Output Parameter:
98: . flg - is the residual computed?
100: Level: advanced
102: .seealso: SNESSetAlwaysComputesFinalResidual()
103: @*/
104: PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES snes, PetscBool *flg)
105: {
108: *flg = snes->alwayscomputesfinalresidual;
109: return(0);
110: }
112: /*@
113: SNESSetFunctionDomainError - tells SNES that the input vector to your SNESFunction is not
114: in the functions domain. For example, negative pressure.
116: Logically Collective on SNES
118: Input Parameters:
119: . snes - the SNES context
121: Level: advanced
123: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction
124: @*/
125: PetscErrorCode SNESSetFunctionDomainError(SNES snes)
126: {
129: if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates input vector is not in the function domain");
130: snes->domainerror = PETSC_TRUE;
131: return(0);
132: }
134: /*@
135: SNESSetJacobianDomainError - tells SNES that computeJacobian does not make sense any more. For example there is a negative element transformation.
137: Logically Collective on SNES
139: Input Parameters:
140: . snes - the SNES context
142: Level: advanced
144: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError()
145: @*/
146: PetscErrorCode SNESSetJacobianDomainError(SNES snes)
147: {
150: if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates computeJacobian does not make sense");
151: snes->jacobiandomainerror = PETSC_TRUE;
152: return(0);
153: }
155: /*@
156: SNESSetCheckJacobianDomainError - if or not to check jacobian domain error after each Jacobian evaluation. By default, we check Jacobian domain error
157: in the debug mode, and do not check it in the optimized mode.
159: Logically Collective on SNES
161: Input Parameters:
162: + snes - the SNES context
163: - flg - indicates if or not to check jacobian domain error after each Jacobian evaluation
165: Level: advanced
167: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESGetCheckJacobianDomainError()
168: @*/
169: PetscErrorCode SNESSetCheckJacobianDomainError(SNES snes, PetscBool flg)
170: {
173: snes->checkjacdomainerror = flg;
174: return(0);
175: }
177: /*@
178: SNESGetCheckJacobianDomainError - Get an indicator whether or not we are checking Jacobian domain errors after each Jacobian evaluation.
180: Logically Collective on SNES
182: Input Parameters:
183: . snes - the SNES context
185: Output Parameters:
186: . flg - PETSC_FALSE indicates that we don't check jacobian domain errors after each Jacobian evaluation
188: Level: advanced
190: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESSetCheckJacobianDomainError()
191: @*/
192: PetscErrorCode SNESGetCheckJacobianDomainError(SNES snes, PetscBool *flg)
193: {
197: *flg = snes->checkjacdomainerror;
198: return(0);
199: }
201: /*@
202: SNESGetFunctionDomainError - Gets the status of the domain error after a call to SNESComputeFunction;
204: Logically Collective on SNES
206: Input Parameters:
207: . snes - the SNES context
209: Output Parameters:
210: . domainerror - Set to PETSC_TRUE if there's a domain error; PETSC_FALSE otherwise.
212: Level: advanced
214: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction()
215: @*/
216: PetscErrorCode SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror)
217: {
221: *domainerror = snes->domainerror;
222: return(0);
223: }
225: /*@
226: SNESGetJacobianDomainError - Gets the status of the Jacobian domain error after a call to SNESComputeJacobian;
228: Logically Collective on SNES
230: Input Parameters:
231: . snes - the SNES context
233: Output Parameters:
234: . domainerror - Set to PETSC_TRUE if there's a jacobian domain error; PETSC_FALSE otherwise.
236: Level: advanced
238: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction(),SNESGetFunctionDomainError()
239: @*/
240: PetscErrorCode SNESGetJacobianDomainError(SNES snes, PetscBool *domainerror)
241: {
245: *domainerror = snes->jacobiandomainerror;
246: return(0);
247: }
249: /*@C
250: SNESLoad - Loads a SNES that has been stored in binary with SNESView().
252: Collective on PetscViewer
254: Input Parameters:
255: + newdm - the newly loaded SNES, this needs to have been created with SNESCreate() or
256: some related function before a call to SNESLoad().
257: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
259: Level: intermediate
261: Notes:
262: The type is determined by the data in the file, any type set into the SNES before this call is ignored.
264: Notes for advanced users:
265: Most users should not need to know the details of the binary storage
266: format, since SNESLoad() and TSView() completely hide these details.
267: But for anyone who's interested, the standard binary matrix storage
268: format is
269: .vb
270: has not yet been determined
271: .ve
273: .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad()
274: @*/
275: PetscErrorCode SNESLoad(SNES snes, PetscViewer viewer)
276: {
278: PetscBool isbinary;
279: PetscInt classid;
280: char type[256];
281: KSP ksp;
282: DM dm;
283: DMSNES dmsnes;
288: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
289: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
291: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
292: if (classid != SNES_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONG,"Not SNES next in file");
293: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
294: SNESSetType(snes, type);
295: if (snes->ops->load) {
296: (*snes->ops->load)(snes,viewer);
297: }
298: SNESGetDM(snes,&dm);
299: DMGetDMSNES(dm,&dmsnes);
300: DMSNESLoad(dmsnes,viewer);
301: SNESGetKSP(snes,&ksp);
302: KSPLoad(ksp,viewer);
303: return(0);
304: }
306: #include <petscdraw.h>
307: #if defined(PETSC_HAVE_SAWS)
308: #include <petscviewersaws.h>
309: #endif
311: /*@C
312: SNESViewFromOptions - View from Options
314: Collective on SNES
316: Input Parameters:
317: + A - the application ordering context
318: . obj - Optional object
319: - name - command line option
321: Level: intermediate
322: .seealso: SNES, SNESView, PetscObjectViewFromOptions(), SNESCreate()
323: @*/
324: PetscErrorCode SNESViewFromOptions(SNES A,PetscObject obj,const char name[])
325: {
330: PetscObjectViewFromOptions((PetscObject)A,obj,name);
331: return(0);
332: }
334: PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES,Vec,Mat,Mat,void*);
336: /*@C
337: SNESView - Prints the SNES data structure.
339: Collective on SNES
341: Input Parameters:
342: + SNES - the SNES context
343: - viewer - visualization context
345: Options Database Key:
346: . -snes_view - Calls SNESView() at end of SNESSolve()
348: Notes:
349: The available visualization contexts include
350: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
351: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
352: output where only the first processor opens
353: the file. All other processors send their
354: data to the first processor to print.
356: The available formats include
357: + PETSC_VIEWER_DEFAULT - standard output (default)
358: - PETSC_VIEWER_ASCII_INFO_DETAIL - more verbose output for SNESNASM
360: The user can open an alternative visualization context with
361: PetscViewerASCIIOpen() - output to a specified file.
363: In the debugger you can do "call SNESView(snes,0)" to display the SNES solver. (The same holds for any PETSc object viewer).
365: Level: beginner
367: .seealso: PetscViewerASCIIOpen()
368: @*/
369: PetscErrorCode SNESView(SNES snes,PetscViewer viewer)
370: {
371: SNESKSPEW *kctx;
373: KSP ksp;
374: SNESLineSearch linesearch;
375: PetscBool iascii,isstring,isbinary,isdraw;
376: DMSNES dmsnes;
377: #if defined(PETSC_HAVE_SAWS)
378: PetscBool issaws;
379: #endif
383: if (!viewer) {
384: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
385: }
389: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
390: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
391: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
392: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
393: #if defined(PETSC_HAVE_SAWS)
394: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
395: #endif
396: if (iascii) {
397: SNESNormSchedule normschedule;
398: DM dm;
399: PetscErrorCode (*cJ)(SNES,Vec,Mat,Mat,void*);
400: void *ctx;
401: const char *pre = "";
403: PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);
404: if (!snes->setupcalled) {
405: PetscViewerASCIIPrintf(viewer," SNES has not been set up so information may be incomplete\n");
406: }
407: if (snes->ops->view) {
408: PetscViewerASCIIPushTab(viewer);
409: (*snes->ops->view)(snes,viewer);
410: PetscViewerASCIIPopTab(viewer);
411: }
412: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);
413: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);
414: if (snes->usesksp) {
415: PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);
416: }
417: PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);
418: SNESGetNormSchedule(snes, &normschedule);
419: if (normschedule > 0) {PetscViewerASCIIPrintf(viewer," norm schedule %s\n",SNESNormSchedules[normschedule]);}
420: if (snes->gridsequence) {
421: PetscViewerASCIIPrintf(viewer," total number of grid sequence refinements=%D\n",snes->gridsequence);
422: }
423: if (snes->ksp_ewconv) {
424: kctx = (SNESKSPEW*)snes->kspconvctx;
425: if (kctx) {
426: PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);
427: PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);
428: PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);
429: }
430: }
431: if (snes->lagpreconditioner == -1) {
432: PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");
433: } else if (snes->lagpreconditioner > 1) {
434: PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);
435: }
436: if (snes->lagjacobian == -1) {
437: PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");
438: } else if (snes->lagjacobian > 1) {
439: PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);
440: }
441: SNESGetDM(snes,&dm);
442: DMSNESGetJacobian(dm,&cJ,&ctx);
443: if (snes->mf_operator) {
444: PetscViewerASCIIPrintf(viewer," Jacobian is applied matrix-free with differencing\n");
445: pre = "Preconditioning ";
446: }
447: if (cJ == SNESComputeJacobianDefault) {
448: PetscViewerASCIIPrintf(viewer," %sJacobian is built using finite differences one column at a time\n",pre);
449: } else if (cJ == SNESComputeJacobianDefaultColor) {
450: PetscViewerASCIIPrintf(viewer," %sJacobian is built using finite differences with coloring\n",pre);
451: /* it slightly breaks data encapsulation for access the DMDA information directly */
452: } else if (cJ == SNESComputeJacobian_DMDA) {
453: MatFDColoring fdcoloring;
454: PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);
455: if (fdcoloring) {
456: PetscViewerASCIIPrintf(viewer," %sJacobian is built using colored finite differences on a DMDA\n",pre);
457: } else {
458: PetscViewerASCIIPrintf(viewer," %sJacobian is built using a DMDA local Jacobian\n",pre);
459: }
460: } else if (snes->mf) {
461: PetscViewerASCIIPrintf(viewer," Jacobian is applied matrix-free with differencing, no explicit Jacobian\n");
462: }
463: } else if (isstring) {
464: const char *type;
465: SNESGetType(snes,&type);
466: PetscViewerStringSPrintf(viewer," SNESType: %-7.7s",type);
467: if (snes->ops->view) {(*snes->ops->view)(snes,viewer);}
468: } else if (isbinary) {
469: PetscInt classid = SNES_FILE_CLASSID;
470: MPI_Comm comm;
471: PetscMPIInt rank;
472: char type[256];
474: PetscObjectGetComm((PetscObject)snes,&comm);
475: MPI_Comm_rank(comm,&rank);
476: if (rank == 0) {
477: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT);
478: PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));
479: PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR);
480: }
481: if (snes->ops->view) {
482: (*snes->ops->view)(snes,viewer);
483: }
484: } else if (isdraw) {
485: PetscDraw draw;
486: char str[36];
487: PetscReal x,y,bottom,h;
489: PetscViewerDrawGetDraw(viewer,0,&draw);
490: PetscDrawGetCurrentPoint(draw,&x,&y);
491: PetscStrncpy(str,"SNES: ",sizeof(str));
492: PetscStrlcat(str,((PetscObject)snes)->type_name,sizeof(str));
493: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);
494: bottom = y - h;
495: PetscDrawPushCurrentPoint(draw,x,bottom);
496: if (snes->ops->view) {
497: (*snes->ops->view)(snes,viewer);
498: }
499: #if defined(PETSC_HAVE_SAWS)
500: } else if (issaws) {
501: PetscMPIInt rank;
502: const char *name;
504: PetscObjectGetName((PetscObject)snes,&name);
505: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
506: if (!((PetscObject)snes)->amsmem && rank == 0) {
507: char dir[1024];
509: PetscObjectViewSAWs((PetscObject)snes,viewer);
510: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
511: PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT));
512: if (!snes->conv_hist) {
513: SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);
514: }
515: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);
516: PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE));
517: }
518: #endif
519: }
520: if (snes->linesearch) {
521: SNESGetLineSearch(snes, &linesearch);
522: PetscViewerASCIIPushTab(viewer);
523: SNESLineSearchView(linesearch, viewer);
524: PetscViewerASCIIPopTab(viewer);
525: }
526: if (snes->npc && snes->usesnpc) {
527: PetscViewerASCIIPushTab(viewer);
528: SNESView(snes->npc, viewer);
529: PetscViewerASCIIPopTab(viewer);
530: }
531: PetscViewerASCIIPushTab(viewer);
532: DMGetDMSNES(snes->dm,&dmsnes);
533: DMSNESView(dmsnes, viewer);
534: PetscViewerASCIIPopTab(viewer);
535: if (snes->usesksp) {
536: SNESGetKSP(snes,&ksp);
537: PetscViewerASCIIPushTab(viewer);
538: KSPView(ksp,viewer);
539: PetscViewerASCIIPopTab(viewer);
540: }
541: if (isdraw) {
542: PetscDraw draw;
543: PetscViewerDrawGetDraw(viewer,0,&draw);
544: PetscDrawPopCurrentPoint(draw);
545: }
546: return(0);
547: }
549: /*
550: We retain a list of functions that also take SNES command
551: line options. These are called at the end SNESSetFromOptions()
552: */
553: #define MAXSETFROMOPTIONS 5
554: static PetscInt numberofsetfromoptions;
555: static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
557: /*@C
558: SNESAddOptionsChecker - Adds an additional function to check for SNES options.
560: Not Collective
562: Input Parameter:
563: . snescheck - function that checks for options
565: Level: developer
567: .seealso: SNESSetFromOptions()
568: @*/
569: PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
570: {
572: if (numberofsetfromoptions >= MAXSETFROMOPTIONS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
573: othersetfromoptions[numberofsetfromoptions++] = snescheck;
574: return(0);
575: }
577: PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
579: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
580: {
581: Mat J;
583: MatNullSpace nullsp;
588: if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
589: Mat A = snes->jacobian, B = snes->jacobian_pre;
590: MatCreateVecs(A ? A : B, NULL,&snes->vec_func);
591: }
593: if (version == 1) {
594: MatCreateSNESMF(snes,&J);
595: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
596: MatSetFromOptions(J);
597: } else if (version == 2) {
598: if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
599: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
600: SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);
601: #else
602: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator routines (version 2)");
603: #endif
604: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator routines, only version 1 and 2");
606: /* attach any user provided null space that was on Amat to the newly created matrix free matrix */
607: if (snes->jacobian) {
608: MatGetNullSpace(snes->jacobian,&nullsp);
609: if (nullsp) {
610: MatSetNullSpace(J,nullsp);
611: }
612: }
614: PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);
615: if (hasOperator) {
617: /* This version replaces the user provided Jacobian matrix with a
618: matrix-free version but still employs the user-provided preconditioner matrix. */
619: SNESSetJacobian(snes,J,NULL,NULL,NULL);
620: } else {
621: /* This version replaces both the user-provided Jacobian and the user-
622: provided preconditioner Jacobian with the default matrix free version. */
623: if ((snes->npcside== PC_LEFT) && snes->npc) {
624: if (!snes->jacobian) {SNESSetJacobian(snes,J,NULL,NULL,NULL);}
625: } else {
626: KSP ksp;
627: PC pc;
628: PetscBool match;
630: SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,NULL);
631: /* Force no preconditioner */
632: SNESGetKSP(snes,&ksp);
633: KSPGetPC(ksp,&pc);
634: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);
635: if (!match) {
636: PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");
637: PCSetType(pc,PCNONE);
638: }
639: }
640: }
641: MatDestroy(&J);
642: return(0);
643: }
645: static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx)
646: {
647: SNES snes = (SNES)ctx;
649: Vec Xfine,Xfine_named = NULL,Xcoarse;
652: if (PetscLogPrintInfo) {
653: PetscInt finelevel,coarselevel,fineclevel,coarseclevel;
654: DMGetRefineLevel(dmfine,&finelevel);
655: DMGetCoarsenLevel(dmfine,&fineclevel);
656: DMGetRefineLevel(dmcoarse,&coarselevel);
657: DMGetCoarsenLevel(dmcoarse,&coarseclevel);
658: PetscInfo4(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);
659: }
660: if (dmfine == snes->dm) Xfine = snes->vec_sol;
661: else {
662: DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
663: Xfine = Xfine_named;
664: }
665: DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
666: if (Inject) {
667: MatRestrict(Inject,Xfine,Xcoarse);
668: } else {
669: MatRestrict(Restrict,Xfine,Xcoarse);
670: VecPointwiseMult(Xcoarse,Xcoarse,Rscale);
671: }
672: DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
673: if (Xfine_named) {DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);}
674: return(0);
675: }
677: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
678: {
682: DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
683: return(0);
684: }
686: /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can
687: * safely call SNESGetDM() in their residual evaluation routine. */
688: static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx)
689: {
690: SNES snes = (SNES)ctx;
692: Vec X,Xnamed = NULL;
693: DM dmsave;
694: void *ctxsave;
695: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL;
698: dmsave = snes->dm;
699: KSPGetDM(ksp,&snes->dm);
700: if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */
701: else { /* We are on a coarser level, this vec was initialized using a DM restrict hook */
702: DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
703: X = Xnamed;
704: SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);
705: /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */
706: if (jac == SNESComputeJacobianDefaultColor) {
707: SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,NULL);
708: }
709: }
710: /* Make sure KSP DM has the Jacobian computation routine */
711: {
712: DMSNES sdm;
714: DMGetDMSNES(snes->dm, &sdm);
715: if (!sdm->ops->computejacobian) {
716: DMCopyDMSNES(dmsave, snes->dm);
717: }
718: }
719: /* Compute the operators */
720: SNESComputeJacobian(snes,X,A,B);
721: /* Put the previous context back */
722: if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
723: SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
724: }
726: if (Xnamed) {DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);}
727: snes->dm = dmsave;
728: return(0);
729: }
731: /*@
732: SNESSetUpMatrices - ensures that matrices are available for SNES, to be called by SNESSetUp_XXX()
734: Collective
736: Input Parameter:
737: . snes - snes to configure
739: Level: developer
741: .seealso: SNESSetUp()
742: @*/
743: PetscErrorCode SNESSetUpMatrices(SNES snes)
744: {
746: DM dm;
747: DMSNES sdm;
750: SNESGetDM(snes,&dm);
751: DMGetDMSNES(dm,&sdm);
752: if (!snes->jacobian && snes->mf) {
753: Mat J;
754: void *functx;
755: MatCreateSNESMF(snes,&J);
756: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
757: MatSetFromOptions(J);
758: SNESGetFunction(snes,NULL,NULL,&functx);
759: SNESSetJacobian(snes,J,J,NULL,NULL);
760: MatDestroy(&J);
761: } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
762: Mat J,B;
763: MatCreateSNESMF(snes,&J);
764: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
765: MatSetFromOptions(J);
766: DMCreateMatrix(snes->dm,&B);
767: /* sdm->computejacobian was already set to reach here */
768: SNESSetJacobian(snes,J,B,NULL,NULL);
769: MatDestroy(&J);
770: MatDestroy(&B);
771: } else if (!snes->jacobian_pre) {
772: PetscDS prob;
773: Mat J, B;
774: PetscBool hasPrec = PETSC_FALSE;
776: J = snes->jacobian;
777: DMGetDS(dm, &prob);
778: if (prob) {PetscDSHasJacobianPreconditioner(prob, &hasPrec);}
779: if (J) {PetscObjectReference((PetscObject) J);}
780: else if (hasPrec) {DMCreateMatrix(snes->dm, &J);}
781: DMCreateMatrix(snes->dm, &B);
782: SNESSetJacobian(snes, J ? J : B, B, NULL, NULL);
783: MatDestroy(&J);
784: MatDestroy(&B);
785: }
786: {
787: KSP ksp;
788: SNESGetKSP(snes,&ksp);
789: KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);
790: DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
791: }
792: return(0);
793: }
795: static PetscErrorCode SNESMonitorPauseFinal_Internal(SNES snes)
796: {
797: PetscInt i;
801: if (!snes->pauseFinal) return(0);
802: for (i = 0; i < snes->numbermonitors; ++i) {
803: PetscViewerAndFormat *vf = (PetscViewerAndFormat *) snes->monitorcontext[i];
804: PetscDraw draw;
805: PetscReal lpause;
807: if (!vf) continue;
808: if (vf->lg) {
810: if (((PetscObject) vf->lg)->classid != PETSC_DRAWLG_CLASSID) continue;
811: PetscDrawLGGetDraw(vf->lg, &draw);
812: PetscDrawGetPause(draw, &lpause);
813: PetscDrawSetPause(draw, -1.0);
814: PetscDrawPause(draw);
815: PetscDrawSetPause(draw, lpause);
816: } else {
817: PetscBool isdraw;
820: if (((PetscObject) vf->viewer)->classid != PETSC_VIEWER_CLASSID) continue;
821: PetscObjectTypeCompare((PetscObject) vf->viewer, PETSCVIEWERDRAW, &isdraw);
822: if (!isdraw) continue;
823: PetscViewerDrawGetDraw(vf->viewer, 0, &draw);
824: PetscDrawGetPause(draw, &lpause);
825: PetscDrawSetPause(draw, -1.0);
826: PetscDrawPause(draw);
827: PetscDrawSetPause(draw, lpause);
828: }
829: }
830: return(0);
831: }
833: /*@C
834: SNESMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
836: Collective on SNES
838: Input Parameters:
839: + snes - SNES object you wish to monitor
840: . name - the monitor type one is seeking
841: . help - message indicating what monitoring is done
842: . manual - manual page for the monitor
843: . monitor - the monitor function
844: - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the SNES or PetscViewer objects
846: Level: developer
848: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
849: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
850: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
851: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
852: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
853: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
854: PetscOptionsFList(), PetscOptionsEList()
855: @*/
856: PetscErrorCode SNESMonitorSetFromOptions(SNES snes,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNES,PetscViewerAndFormat*))
857: {
858: PetscErrorCode ierr;
859: PetscViewer viewer;
860: PetscViewerFormat format;
861: PetscBool flg;
864: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,name,&viewer,&format,&flg);
865: if (flg) {
866: PetscViewerAndFormat *vf;
867: PetscViewerAndFormatCreate(viewer,format,&vf);
868: PetscObjectDereference((PetscObject)viewer);
869: if (monitorsetup) {
870: (*monitorsetup)(snes,vf);
871: }
872: SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
873: }
874: return(0);
875: }
877: /*@
878: SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
880: Collective on SNES
882: Input Parameter:
883: . snes - the SNES context
885: Options Database Keys:
886: + -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
887: . -snes_stol - convergence tolerance in terms of the norm
888: of the change in the solution between steps
889: . -snes_atol <abstol> - absolute tolerance of residual norm
890: . -snes_rtol <rtol> - relative decrease in tolerance norm from initial
891: . -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
892: . -snes_force_iteration <force> - force SNESSolve() to take at least one iteration
893: . -snes_max_it <max_it> - maximum number of iterations
894: . -snes_max_funcs <max_funcs> - maximum number of function evaluations
895: . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
896: . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
897: . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
898: . -snes_lag_preconditioner_persists <true,false> - retains the -snes_lag_preconditioner information across multiple SNESSolve()
899: . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
900: . -snes_lag_jacobian_persists <true,false> - retains the -snes_lag_jacobian information across multiple SNESSolve()
901: . -snes_trtol <trtol> - trust region tolerance
902: . -snes_no_convergence_test - skip convergence test in nonlinear
903: solver; hence iterations will continue until max_it
904: or some other criterion is reached. Saves expense
905: of convergence test
906: . -snes_monitor [ascii][:filename][:viewer format] - prints residual norm at each iteration. if no filename given prints to stdout
907: . -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration
908: . -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration
909: . -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration
910: . -snes_monitor_lg_residualnorm - plots residual norm at each iteration
911: . -snes_monitor_lg_range - plots residual norm at each iteration
912: . -snes_monitor_pause_final - Pauses all monitor drawing after the solver ends
913: . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
914: . -snes_fd_color - use finite differences with coloring to compute Jacobian
915: . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
916: . -snes_converged_reason - print the reason for convergence/divergence after each solve
917: . -npc_snes_type <type> - the SNES type to use as a nonlinear preconditioner
918: . -snes_test_jacobian <optional threshold> - compare the user provided Jacobian with one computed via finite differences to check for errors. If a threshold is given, display only those entries whose difference is greater than the threshold.
919: - -snes_test_jacobian_view - display the user provided Jacobian, the finite difference Jacobian and the difference between them to help users detect the location of errors in the user provided Jacobian.
921: Options Database for Eisenstat-Walker method:
922: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
923: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
924: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
925: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
926: . -snes_ksp_ew_gamma <gamma> - Sets gamma
927: . -snes_ksp_ew_alpha <alpha> - Sets alpha
928: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
929: - -snes_ksp_ew_threshold <threshold> - Sets threshold
931: Notes:
932: To see all options, run your program with the -help option or consult the users manual
934: Notes:
935: SNES supports three approaches for computing (approximate) Jacobians: user provided via SNESSetJacobian(), matrix free, and computing explicitly with
936: finite differences and coloring using MatFDColoring. It is also possible to use automatic differentiation and the MatFDColoring object.
938: Level: beginner
940: .seealso: SNESSetOptionsPrefix(), SNESResetFromOptions(), SNES, SNESCreate()
941: @*/
942: PetscErrorCode SNESSetFromOptions(SNES snes)
943: {
944: PetscBool flg,pcset,persist,set;
945: PetscInt i,indx,lag,grids;
946: const char *deft = SNESNEWTONLS;
947: const char *convtests[] = {"default","skip","correct_pressure"};
948: SNESKSPEW *kctx = NULL;
949: char type[256], monfilename[PETSC_MAX_PATH_LEN];
951: PCSide pcside;
952: const char *optionsprefix;
956: SNESRegisterAll();
957: PetscObjectOptionsBegin((PetscObject)snes);
958: if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
959: PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
960: if (flg) {
961: SNESSetType(snes,type);
962: } else if (!((PetscObject)snes)->type_name) {
963: SNESSetType(snes,deft);
964: }
965: PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);
966: PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);
968: PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,NULL);
969: PetscOptionsReal("-snes_divergence_tolerance","Stop if residual norm increases by this factor","SNESSetDivergenceTolerance",snes->divtol,&snes->divtol,NULL);
970: PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,NULL);
971: PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,NULL);
972: PetscOptionsInt("-snes_max_fail","Maximum nonlinear step failures","SNESSetMaxNonlinearStepFailures",snes->maxFailures,&snes->maxFailures,NULL);
973: PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,NULL);
974: PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,NULL);
975: PetscOptionsBool("-snes_force_iteration","Force SNESSolve() to take at least one iteration","SNESSetForceIteration",snes->forceiteration,&snes->forceiteration,NULL);
976: PetscOptionsBool("-snes_check_jacobian_domain_error","Check Jacobian domain error after Jacobian evaluation","SNESCheckJacobianDomainError",snes->checkjacdomainerror,&snes->checkjacdomainerror,NULL);
978: PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
979: if (flg) {
980: if (lag == -1) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_USER,"Cannot set the lag to -1 from the command line since the preconditioner must be built as least once, perhaps you mean -2");
981: SNESSetLagPreconditioner(snes,lag);
982: }
983: PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple SNES solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
984: if (flg) {
985: SNESSetLagPreconditionerPersists(snes,persist);
986: }
987: PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
988: if (flg) {
989: if (lag == -1) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_USER,"Cannot set the lag to -1 from the command line since the Jacobian must be built as least once, perhaps you mean -2");
990: SNESSetLagJacobian(snes,lag);
991: }
992: PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple SNES solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
993: if (flg) {
994: SNESSetLagJacobianPersists(snes,persist);
995: }
997: PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
998: if (flg) {
999: SNESSetGridSequence(snes,grids);
1000: }
1002: PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,sizeof(convtests)/sizeof(char*),"default",&indx,&flg);
1003: if (flg) {
1004: switch (indx) {
1005: case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
1006: case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL); break;
1007: case 2: SNESSetConvergenceTest(snes,SNESConvergedCorrectPressure,NULL,NULL); break;
1008: }
1009: }
1011: PetscOptionsEList("-snes_norm_schedule","SNES Norm schedule","SNESSetNormSchedule",SNESNormSchedules,5,"function",&indx,&flg);
1012: if (flg) { SNESSetNormSchedule(snes,(SNESNormSchedule)indx); }
1014: PetscOptionsEList("-snes_function_type","SNES Norm schedule","SNESSetFunctionType",SNESFunctionTypes,2,"unpreconditioned",&indx,&flg);
1015: if (flg) { SNESSetFunctionType(snes,(SNESFunctionType)indx); }
1017: kctx = (SNESKSPEW*)snes->kspconvctx;
1019: PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,NULL);
1021: PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);
1022: PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);
1023: PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);
1024: PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);
1025: PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);
1026: PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);
1027: PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);
1029: flg = PETSC_FALSE;
1030: PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);
1031: if (set && flg) {SNESMonitorCancel(snes);}
1033: SNESMonitorSetFromOptions(snes,"-snes_monitor","Monitor norm of function","SNESMonitorDefault",SNESMonitorDefault,SNESMonitorDefaultSetUp);
1034: SNESMonitorSetFromOptions(snes,"-snes_monitor_short","Monitor norm of function with fewer digits","SNESMonitorDefaultShort",SNESMonitorDefaultShort,NULL);
1035: SNESMonitorSetFromOptions(snes,"-snes_monitor_range","Monitor range of elements of function","SNESMonitorRange",SNESMonitorRange,NULL);
1037: SNESMonitorSetFromOptions(snes,"-snes_monitor_ratio","Monitor ratios of the norm of function for consecutive steps","SNESMonitorRatio",SNESMonitorRatio,SNESMonitorRatioSetUp);
1038: SNESMonitorSetFromOptions(snes,"-snes_monitor_field","Monitor norm of function (split into fields)","SNESMonitorDefaultField",SNESMonitorDefaultField,NULL);
1039: SNESMonitorSetFromOptions(snes,"-snes_monitor_solution","View solution at each iteration","SNESMonitorSolution",SNESMonitorSolution,NULL);
1040: SNESMonitorSetFromOptions(snes,"-snes_monitor_solution_update","View correction at each iteration","SNESMonitorSolutionUpdate",SNESMonitorSolutionUpdate,NULL);
1041: SNESMonitorSetFromOptions(snes,"-snes_monitor_residual","View residual at each iteration","SNESMonitorResidual",SNESMonitorResidual,NULL);
1042: SNESMonitorSetFromOptions(snes,"-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",SNESMonitorJacUpdateSpectrum,NULL);
1043: SNESMonitorSetFromOptions(snes,"-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet",SNESMonitorFields,NULL);
1044: PetscOptionsBool("-snes_monitor_pause_final", "Pauses all draw monitors at the final iterate", "SNESMonitorPauseFinal_Internal", PETSC_FALSE, &snes->pauseFinal, NULL);
1046: PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",NULL,monfilename,sizeof(monfilename),&flg);
1047: if (flg) {PetscPythonMonitorSet((PetscObject)snes,monfilename);}
1049: flg = PETSC_FALSE;
1050: PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
1051: if (flg) {
1052: PetscViewer ctx;
1054: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
1055: SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
1056: }
1058: flg = PETSC_FALSE;
1059: PetscOptionsBool("-snes_converged_reason_view_cancel","Remove all converged reason viewers","SNESConvergedReasonViewCancel",flg,&flg,&set);
1060: if (set && flg) {SNESConvergedReasonViewCancel(snes);}
1062: flg = PETSC_FALSE;
1063: PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
1064: if (flg) {
1065: void *functx;
1066: DM dm;
1067: DMSNES sdm;
1068: SNESGetDM(snes,&dm);
1069: DMGetDMSNES(dm,&sdm);
1070: sdm->jacobianctx = NULL;
1071: SNESGetFunction(snes,NULL,NULL,&functx);
1072: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
1073: PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
1074: }
1076: flg = PETSC_FALSE;
1077: PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
1078: if (flg) {
1079: SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
1080: }
1082: flg = PETSC_FALSE;
1083: PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);
1084: if (flg) {
1085: DM dm;
1086: DMSNES sdm;
1087: SNESGetDM(snes,&dm);
1088: DMGetDMSNES(dm,&sdm);
1089: sdm->jacobianctx = NULL;
1090: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,NULL);
1091: PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
1092: }
1094: flg = PETSC_FALSE;
1095: PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf_operator,&flg);
1096: if (flg && snes->mf_operator) {
1097: snes->mf_operator = PETSC_TRUE;
1098: snes->mf = PETSC_TRUE;
1099: }
1100: flg = PETSC_FALSE;
1101: PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf,&flg);
1102: if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE;
1103: PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,NULL);
1105: flg = PETSC_FALSE;
1106: SNESGetNPCSide(snes,&pcside);
1107: PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
1108: if (flg) {SNESSetNPCSide(snes,pcside);}
1110: #if defined(PETSC_HAVE_SAWS)
1111: /*
1112: Publish convergence information using SAWs
1113: */
1114: flg = PETSC_FALSE;
1115: PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);
1116: if (flg) {
1117: void *ctx;
1118: SNESMonitorSAWsCreate(snes,&ctx);
1119: SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);
1120: }
1121: #endif
1122: #if defined(PETSC_HAVE_SAWS)
1123: {
1124: PetscBool set;
1125: flg = PETSC_FALSE;
1126: PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);
1127: if (set) {
1128: PetscObjectSAWsSetBlock((PetscObject)snes,flg);
1129: }
1130: }
1131: #endif
1133: for (i = 0; i < numberofsetfromoptions; i++) {
1134: (*othersetfromoptions[i])(snes);
1135: }
1137: if (snes->ops->setfromoptions) {
1138: (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
1139: }
1141: /* process any options handlers added with PetscObjectAddOptionsHandler() */
1142: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)snes);
1143: PetscOptionsEnd();
1145: if (snes->linesearch) {
1146: SNESGetLineSearch(snes, &snes->linesearch);
1147: SNESLineSearchSetFromOptions(snes->linesearch);
1148: }
1150: if (snes->usesksp) {
1151: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
1152: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
1153: KSPSetFromOptions(snes->ksp);
1154: }
1156: /* if user has set the SNES NPC type via options database, create it. */
1157: SNESGetOptionsPrefix(snes, &optionsprefix);
1158: PetscOptionsHasName(((PetscObject)snes)->options,optionsprefix, "-npc_snes_type", &pcset);
1159: if (pcset && (!snes->npc)) {
1160: SNESGetNPC(snes, &snes->npc);
1161: }
1162: if (snes->npc) {
1163: SNESSetFromOptions(snes->npc);
1164: }
1165: snes->setfromoptionscalled++;
1166: return(0);
1167: }
1169: /*@
1170: SNESResetFromOptions - Sets various SNES and KSP parameters from user options ONLY if the SNES was previously set from options
1172: Collective on SNES
1174: Input Parameter:
1175: . snes - the SNES context
1177: Level: beginner
1179: .seealso: SNESSetFromOptions(), SNESSetOptionsPrefix()
1180: @*/
1181: PetscErrorCode SNESResetFromOptions(SNES snes)
1182: {
1186: if (snes->setfromoptionscalled) {SNESSetFromOptions(snes);}
1187: return(0);
1188: }
1190: /*@C
1191: SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
1192: the nonlinear solvers.
1194: Logically Collective on SNES
1196: Input Parameters:
1197: + snes - the SNES context
1198: . compute - function to compute the context
1199: - destroy - function to destroy the context
1201: Level: intermediate
1203: Notes:
1204: This function is currently not available from Fortran.
1206: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
1207: @*/
1208: PetscErrorCode SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
1209: {
1212: snes->ops->usercompute = compute;
1213: snes->ops->userdestroy = destroy;
1214: return(0);
1215: }
1217: /*@
1218: SNESSetApplicationContext - Sets the optional user-defined context for
1219: the nonlinear solvers.
1221: Logically Collective on SNES
1223: Input Parameters:
1224: + snes - the SNES context
1225: - usrP - optional user context
1227: Level: intermediate
1229: Fortran Notes:
1230: To use this from Fortran you must write a Fortran interface definition for this
1231: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1233: .seealso: SNESGetApplicationContext()
1234: @*/
1235: PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP)
1236: {
1238: KSP ksp;
1242: SNESGetKSP(snes,&ksp);
1243: KSPSetApplicationContext(ksp,usrP);
1244: snes->user = usrP;
1245: return(0);
1246: }
1248: /*@
1249: SNESGetApplicationContext - Gets the user-defined context for the
1250: nonlinear solvers.
1252: Not Collective
1254: Input Parameter:
1255: . snes - SNES context
1257: Output Parameter:
1258: . usrP - user context
1260: Fortran Notes:
1261: To use this from Fortran you must write a Fortran interface definition for this
1262: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1264: Level: intermediate
1266: .seealso: SNESSetApplicationContext()
1267: @*/
1268: PetscErrorCode SNESGetApplicationContext(SNES snes,void *usrP)
1269: {
1272: *(void**)usrP = snes->user;
1273: return(0);
1274: }
1276: /*@
1277: SNESSetUseMatrixFree - indicates that SNES should use matrix free finite difference matrix vector products internally to apply the Jacobian.
1279: Collective on SNES
1281: Input Parameters:
1282: + snes - SNES context
1283: . mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used
1284: - mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1286: Options Database:
1287: + -snes_mf - use matrix free for both the mat and pmat operator
1288: . -snes_mf_operator - use matrix free only for the mat operator
1289: . -snes_fd_color - compute the Jacobian via coloring and finite differences.
1290: - -snes_fd - compute the Jacobian via finite differences (slow)
1292: Level: intermediate
1294: Notes:
1295: SNES supports three approaches for computing (approximate) Jacobians: user provided via SNESSetJacobian(), matrix free, and computing explicitly with
1296: finite differences and coloring using MatFDColoring. It is also possible to use automatic differentiation and the MatFDColoring object.
1298: .seealso: SNESGetUseMatrixFree(), MatCreateSNESMF(), SNESComputeJacobianDefaultColor()
1299: @*/
1300: PetscErrorCode SNESSetUseMatrixFree(SNES snes,PetscBool mf_operator,PetscBool mf)
1301: {
1306: snes->mf = mf_operator ? PETSC_TRUE : mf;
1307: snes->mf_operator = mf_operator;
1308: return(0);
1309: }
1311: /*@
1312: SNESGetUseMatrixFree - indicates if the SNES uses matrix free finite difference matrix vector products to apply the Jacobian.
1314: Collective on SNES
1316: Input Parameter:
1317: . snes - SNES context
1319: Output Parameters:
1320: + mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used
1321: - mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1323: Options Database:
1324: + -snes_mf - use matrix free for both the mat and pmat operator
1325: - -snes_mf_operator - use matrix free only for the mat operator
1327: Level: intermediate
1329: .seealso: SNESSetUseMatrixFree(), MatCreateSNESMF()
1330: @*/
1331: PetscErrorCode SNESGetUseMatrixFree(SNES snes,PetscBool *mf_operator,PetscBool *mf)
1332: {
1335: if (mf) *mf = snes->mf;
1336: if (mf_operator) *mf_operator = snes->mf_operator;
1337: return(0);
1338: }
1340: /*@
1341: SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1342: at this time.
1344: Not Collective
1346: Input Parameter:
1347: . snes - SNES context
1349: Output Parameter:
1350: . iter - iteration number
1352: Notes:
1353: For example, during the computation of iteration 2 this would return 1.
1355: This is useful for using lagged Jacobians (where one does not recompute the
1356: Jacobian at each SNES iteration). For example, the code
1357: .vb
1358: SNESGetIterationNumber(snes,&it);
1359: if (!(it % 2)) {
1360: [compute Jacobian here]
1361: }
1362: .ve
1363: can be used in your ComputeJacobian() function to cause the Jacobian to be
1364: recomputed every second SNES iteration.
1366: After the SNES solve is complete this will return the number of nonlinear iterations used.
1368: Level: intermediate
1370: .seealso: SNESGetLinearSolveIterations()
1371: @*/
1372: PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt *iter)
1373: {
1377: *iter = snes->iter;
1378: return(0);
1379: }
1381: /*@
1382: SNESSetIterationNumber - Sets the current iteration number.
1384: Not Collective
1386: Input Parameters:
1387: + snes - SNES context
1388: - iter - iteration number
1390: Level: developer
1392: .seealso: SNESGetLinearSolveIterations()
1393: @*/
1394: PetscErrorCode SNESSetIterationNumber(SNES snes,PetscInt iter)
1395: {
1400: PetscObjectSAWsTakeAccess((PetscObject)snes);
1401: snes->iter = iter;
1402: PetscObjectSAWsGrantAccess((PetscObject)snes);
1403: return(0);
1404: }
1406: /*@
1407: SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1408: attempted by the nonlinear solver.
1410: Not Collective
1412: Input Parameter:
1413: . snes - SNES context
1415: Output Parameter:
1416: . nfails - number of unsuccessful steps attempted
1418: Notes:
1419: This counter is reset to zero for each successive call to SNESSolve().
1421: Level: intermediate
1423: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1424: SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1425: @*/
1426: PetscErrorCode SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1427: {
1431: *nfails = snes->numFailures;
1432: return(0);
1433: }
1435: /*@
1436: SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1437: attempted by the nonlinear solver before it gives up.
1439: Not Collective
1441: Input Parameters:
1442: + snes - SNES context
1443: - maxFails - maximum of unsuccessful steps
1445: Level: intermediate
1447: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1448: SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1449: @*/
1450: PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1451: {
1454: snes->maxFailures = maxFails;
1455: return(0);
1456: }
1458: /*@
1459: SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1460: attempted by the nonlinear solver before it gives up.
1462: Not Collective
1464: Input Parameter:
1465: . snes - SNES context
1467: Output Parameter:
1468: . maxFails - maximum of unsuccessful steps
1470: Level: intermediate
1472: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1473: SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1475: @*/
1476: PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1477: {
1481: *maxFails = snes->maxFailures;
1482: return(0);
1483: }
1485: /*@
1486: SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1487: done by SNES.
1489: Not Collective
1491: Input Parameter:
1492: . snes - SNES context
1494: Output Parameter:
1495: . nfuncs - number of evaluations
1497: Level: intermediate
1499: Notes:
1500: Reset every time SNESSolve is called unless SNESSetCountersReset() is used.
1502: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1503: @*/
1504: PetscErrorCode SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1505: {
1509: *nfuncs = snes->nfuncs;
1510: return(0);
1511: }
1513: /*@
1514: SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1515: linear solvers.
1517: Not Collective
1519: Input Parameter:
1520: . snes - SNES context
1522: Output Parameter:
1523: . nfails - number of failed solves
1525: Level: intermediate
1527: Options Database Keys:
1528: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1530: Notes:
1531: This counter is reset to zero for each successive call to SNESSolve().
1533: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1534: @*/
1535: PetscErrorCode SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1536: {
1540: *nfails = snes->numLinearSolveFailures;
1541: return(0);
1542: }
1544: /*@
1545: SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1546: allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
1548: Logically Collective on SNES
1550: Input Parameters:
1551: + snes - SNES context
1552: - maxFails - maximum allowed linear solve failures
1554: Level: intermediate
1556: Options Database Keys:
1557: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1559: Notes:
1560: By default this is 0; that is SNES returns on the first failed linear solve
1562: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1563: @*/
1564: PetscErrorCode SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1565: {
1569: snes->maxLinearSolveFailures = maxFails;
1570: return(0);
1571: }
1573: /*@
1574: SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1575: are allowed before SNES terminates
1577: Not Collective
1579: Input Parameter:
1580: . snes - SNES context
1582: Output Parameter:
1583: . maxFails - maximum of unsuccessful solves allowed
1585: Level: intermediate
1587: Notes:
1588: By default this is 1; that is SNES returns on the first failed linear solve
1590: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1591: @*/
1592: PetscErrorCode SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1593: {
1597: *maxFails = snes->maxLinearSolveFailures;
1598: return(0);
1599: }
1601: /*@
1602: SNESGetLinearSolveIterations - Gets the total number of linear iterations
1603: used by the nonlinear solver.
1605: Not Collective
1607: Input Parameter:
1608: . snes - SNES context
1610: Output Parameter:
1611: . lits - number of linear iterations
1613: Notes:
1614: This counter is reset to zero for each successive call to SNESSolve() unless SNESSetCountersReset() is used.
1616: If the linear solver fails inside the SNESSolve() the iterations for that call to the linear solver are not included. If you wish to count them
1617: then call KSPGetIterationNumber() after the failed solve.
1619: Level: intermediate
1621: .seealso: SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1622: @*/
1623: PetscErrorCode SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1624: {
1628: *lits = snes->linear_its;
1629: return(0);
1630: }
1632: /*@
1633: SNESSetCountersReset - Sets whether or not the counters for linear iterations and function evaluations
1634: are reset every time SNESSolve() is called.
1636: Logically Collective on SNES
1638: Input Parameters:
1639: + snes - SNES context
1640: - reset - whether to reset the counters or not
1642: Notes:
1643: This defaults to PETSC_TRUE
1645: Level: developer
1647: .seealso: SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC()
1648: @*/
1649: PetscErrorCode SNESSetCountersReset(SNES snes,PetscBool reset)
1650: {
1654: snes->counters_reset = reset;
1655: return(0);
1656: }
1658: /*@
1659: SNESSetKSP - Sets a KSP context for the SNES object to use
1661: Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
1663: Input Parameters:
1664: + snes - the SNES context
1665: - ksp - the KSP context
1667: Notes:
1668: The SNES object already has its KSP object, you can obtain with SNESGetKSP()
1669: so this routine is rarely needed.
1671: The KSP object that is already in the SNES object has its reference count
1672: decreased by one.
1674: Level: developer
1676: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1677: @*/
1678: PetscErrorCode SNESSetKSP(SNES snes,KSP ksp)
1679: {
1686: PetscObjectReference((PetscObject)ksp);
1687: if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1688: snes->ksp = ksp;
1689: return(0);
1690: }
1692: /* -----------------------------------------------------------*/
1693: /*@
1694: SNESCreate - Creates a nonlinear solver context.
1696: Collective
1698: Input Parameters:
1699: . comm - MPI communicator
1701: Output Parameter:
1702: . outsnes - the new SNES context
1704: Options Database Keys:
1705: + -snes_mf - Activates default matrix-free Jacobian-vector products,
1706: and no preconditioning matrix
1707: . -snes_mf_operator - Activates default matrix-free Jacobian-vector
1708: products, and a user-provided preconditioning matrix
1709: as set by SNESSetJacobian()
1710: - -snes_fd - Uses (slow!) finite differences to compute Jacobian
1712: Level: beginner
1714: Developer Notes:
1715: SNES always creates a KSP object even though many SNES methods do not use it. This is
1716: unfortunate and should be fixed at some point. The flag snes->usesksp indicates if the
1717: particular method does use KSP and regulates if the information about the KSP is printed
1718: in SNESView(). TSSetFromOptions() does call SNESSetFromOptions() which can lead to users being confused
1719: by help messages about meaningless SNES options.
1721: SNES always creates the snes->kspconvctx even though it is used by only one type. This should
1722: be fixed.
1724: .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner(), SNESSetLagJacobian()
1726: @*/
1727: PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes)
1728: {
1730: SNES snes;
1731: SNESKSPEW *kctx;
1735: *outsnes = NULL;
1736: SNESInitializePackage();
1738: PetscHeaderCreate(snes,SNES_CLASSID,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);
1740: snes->ops->converged = SNESConvergedDefault;
1741: snes->usesksp = PETSC_TRUE;
1742: snes->tolerancesset = PETSC_FALSE;
1743: snes->max_its = 50;
1744: snes->max_funcs = 10000;
1745: snes->norm = 0.0;
1746: snes->xnorm = 0.0;
1747: snes->ynorm = 0.0;
1748: snes->normschedule = SNES_NORM_ALWAYS;
1749: snes->functype = SNES_FUNCTION_DEFAULT;
1750: #if defined(PETSC_USE_REAL_SINGLE)
1751: snes->rtol = 1.e-5;
1752: #else
1753: snes->rtol = 1.e-8;
1754: #endif
1755: snes->ttol = 0.0;
1756: #if defined(PETSC_USE_REAL_SINGLE)
1757: snes->abstol = 1.e-25;
1758: #else
1759: snes->abstol = 1.e-50;
1760: #endif
1761: #if defined(PETSC_USE_REAL_SINGLE)
1762: snes->stol = 1.e-5;
1763: #else
1764: snes->stol = 1.e-8;
1765: #endif
1766: #if defined(PETSC_USE_REAL_SINGLE)
1767: snes->deltatol = 1.e-6;
1768: #else
1769: snes->deltatol = 1.e-12;
1770: #endif
1771: snes->divtol = 1.e4;
1772: snes->rnorm0 = 0;
1773: snes->nfuncs = 0;
1774: snes->numFailures = 0;
1775: snes->maxFailures = 1;
1776: snes->linear_its = 0;
1777: snes->lagjacobian = 1;
1778: snes->jac_iter = 0;
1779: snes->lagjac_persist = PETSC_FALSE;
1780: snes->lagpreconditioner = 1;
1781: snes->pre_iter = 0;
1782: snes->lagpre_persist = PETSC_FALSE;
1783: snes->numbermonitors = 0;
1784: snes->numberreasonviews = 0;
1785: snes->data = NULL;
1786: snes->setupcalled = PETSC_FALSE;
1787: snes->ksp_ewconv = PETSC_FALSE;
1788: snes->nwork = 0;
1789: snes->work = NULL;
1790: snes->nvwork = 0;
1791: snes->vwork = NULL;
1792: snes->conv_hist_len = 0;
1793: snes->conv_hist_max = 0;
1794: snes->conv_hist = NULL;
1795: snes->conv_hist_its = NULL;
1796: snes->conv_hist_reset = PETSC_TRUE;
1797: snes->counters_reset = PETSC_TRUE;
1798: snes->vec_func_init_set = PETSC_FALSE;
1799: snes->reason = SNES_CONVERGED_ITERATING;
1800: snes->npcside = PC_RIGHT;
1801: snes->setfromoptionscalled = 0;
1803: snes->mf = PETSC_FALSE;
1804: snes->mf_operator = PETSC_FALSE;
1805: snes->mf_version = 1;
1807: snes->numLinearSolveFailures = 0;
1808: snes->maxLinearSolveFailures = 1;
1810: snes->vizerotolerance = 1.e-8;
1811: snes->checkjacdomainerror = PetscDefined(USE_DEBUG) ? PETSC_TRUE : PETSC_FALSE;
1813: /* Set this to true if the implementation of SNESSolve_XXX does compute the residual at the final solution. */
1814: snes->alwayscomputesfinalresidual = PETSC_FALSE;
1816: /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
1817: PetscNewLog(snes,&kctx);
1819: snes->kspconvctx = (void*)kctx;
1820: kctx->version = 2;
1821: kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1822: this was too large for some test cases */
1823: kctx->rtol_last = 0.0;
1824: kctx->rtol_max = .9;
1825: kctx->gamma = 1.0;
1826: kctx->alpha = .5*(1.0 + PetscSqrtReal(5.0));
1827: kctx->alpha2 = kctx->alpha;
1828: kctx->threshold = .1;
1829: kctx->lresid_last = 0.0;
1830: kctx->norm_last = 0.0;
1832: *outsnes = snes;
1833: return(0);
1834: }
1836: /*MC
1837: SNESFunction - Functional form used to convey the nonlinear function to be solved by SNES
1839: Synopsis:
1840: #include "petscsnes.h"
1841: PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);
1843: Collective on snes
1845: Input Parameters:
1846: + snes - the SNES context
1847: . x - state at which to evaluate residual
1848: - ctx - optional user-defined function context, passed in with SNESSetFunction()
1850: Output Parameter:
1851: . f - vector to put residual (function value)
1853: Level: intermediate
1855: .seealso: SNESSetFunction(), SNESGetFunction()
1856: M*/
1858: /*@C
1859: SNESSetFunction - Sets the function evaluation routine and function
1860: vector for use by the SNES routines in solving systems of nonlinear
1861: equations.
1863: Logically Collective on SNES
1865: Input Parameters:
1866: + snes - the SNES context
1867: . r - vector to store function values, may be NULL
1868: . f - function evaluation routine; see SNESFunction for calling sequence details
1869: - ctx - [optional] user-defined context for private data for the
1870: function evaluation routine (may be NULL)
1872: Notes:
1873: The Newton-like methods typically solve linear systems of the form
1874: $ f'(x) x = -f(x),
1875: where f'(x) denotes the Jacobian matrix and f(x) is the function.
1877: Level: beginner
1879: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1880: @*/
1881: PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1882: {
1884: DM dm;
1888: if (r) {
1891: PetscObjectReference((PetscObject)r);
1892: VecDestroy(&snes->vec_func);
1893: snes->vec_func = r;
1894: }
1895: SNESGetDM(snes,&dm);
1896: DMSNESSetFunction(dm,f,ctx);
1897: if (f == SNESPicardComputeFunction) {
1898: DMSNESSetMFFunction(dm,SNESPicardComputeMFFunction,ctx);
1899: }
1900: return(0);
1901: }
1903: /*@C
1904: SNESSetInitialFunction - Sets the function vector to be used as the
1905: function norm at the initialization of the method. In some
1906: instances, the user has precomputed the function before calling
1907: SNESSolve. This function allows one to avoid a redundant call
1908: to SNESComputeFunction in that case.
1910: Logically Collective on SNES
1912: Input Parameters:
1913: + snes - the SNES context
1914: - f - vector to store function value
1916: Notes:
1917: This should not be modified during the solution procedure.
1919: This is used extensively in the SNESFAS hierarchy and in nonlinear preconditioning.
1921: Level: developer
1923: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1924: @*/
1925: PetscErrorCode SNESSetInitialFunction(SNES snes, Vec f)
1926: {
1928: Vec vec_func;
1934: if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1935: snes->vec_func_init_set = PETSC_FALSE;
1936: return(0);
1937: }
1938: SNESGetFunction(snes,&vec_func,NULL,NULL);
1939: VecCopy(f, vec_func);
1941: snes->vec_func_init_set = PETSC_TRUE;
1942: return(0);
1943: }
1945: /*@
1946: SNESSetNormSchedule - Sets the SNESNormSchedule used in convergence and monitoring
1947: of the SNES method.
1949: Logically Collective on SNES
1951: Input Parameters:
1952: + snes - the SNES context
1953: - normschedule - the frequency of norm computation
1955: Options Database Key:
1956: . -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly>
1958: Notes:
1959: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1960: of the nonlinear function and the taking of its norm at every iteration to
1961: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1962: (SNESNGS) and the like do not require the norm of the function to be computed, and therefore
1963: may either be monitored for convergence or not. As these are often used as nonlinear
1964: preconditioners, monitoring the norm of their error is not a useful enterprise within
1965: their solution.
1967: Level: developer
1969: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1970: @*/
1971: PetscErrorCode SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1972: {
1975: snes->normschedule = normschedule;
1976: return(0);
1977: }
1979: /*@
1980: SNESGetNormSchedule - Gets the SNESNormSchedule used in convergence and monitoring
1981: of the SNES method.
1983: Logically Collective on SNES
1985: Input Parameters:
1986: + snes - the SNES context
1987: - normschedule - the type of the norm used
1989: Level: advanced
1991: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1992: @*/
1993: PetscErrorCode SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1994: {
1997: *normschedule = snes->normschedule;
1998: return(0);
1999: }
2001: /*@
2002: SNESSetFunctionNorm - Sets the last computed residual norm.
2004: Logically Collective on SNES
2006: Input Parameters:
2007: + snes - the SNES context
2009: - normschedule - the frequency of norm computation
2011: Level: developer
2013: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2014: @*/
2015: PetscErrorCode SNESSetFunctionNorm(SNES snes, PetscReal norm)
2016: {
2019: snes->norm = norm;
2020: return(0);
2021: }
2023: /*@
2024: SNESGetFunctionNorm - Gets the last computed norm of the residual
2026: Not Collective
2028: Input Parameter:
2029: . snes - the SNES context
2031: Output Parameter:
2032: . norm - the last computed residual norm
2034: Level: developer
2036: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2037: @*/
2038: PetscErrorCode SNESGetFunctionNorm(SNES snes, PetscReal *norm)
2039: {
2043: *norm = snes->norm;
2044: return(0);
2045: }
2047: /*@
2048: SNESGetUpdateNorm - Gets the last computed norm of the Newton update
2050: Not Collective
2052: Input Parameter:
2053: . snes - the SNES context
2055: Output Parameter:
2056: . ynorm - the last computed update norm
2058: Level: developer
2060: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm()
2061: @*/
2062: PetscErrorCode SNESGetUpdateNorm(SNES snes, PetscReal *ynorm)
2063: {
2067: *ynorm = snes->ynorm;
2068: return(0);
2069: }
2071: /*@
2072: SNESGetSolutionNorm - Gets the last computed norm of the solution
2074: Not Collective
2076: Input Parameter:
2077: . snes - the SNES context
2079: Output Parameter:
2080: . xnorm - the last computed solution norm
2082: Level: developer
2084: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm(), SNESGetUpdateNorm()
2085: @*/
2086: PetscErrorCode SNESGetSolutionNorm(SNES snes, PetscReal *xnorm)
2087: {
2091: *xnorm = snes->xnorm;
2092: return(0);
2093: }
2095: /*@C
2096: SNESSetFunctionType - Sets the SNESNormSchedule used in convergence and monitoring
2097: of the SNES method.
2099: Logically Collective on SNES
2101: Input Parameters:
2102: + snes - the SNES context
2103: - normschedule - the frequency of norm computation
2105: Notes:
2106: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
2107: of the nonlinear function and the taking of its norm at every iteration to
2108: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
2109: (SNESNGS) and the like do not require the norm of the function to be computed, and therefore
2110: may either be monitored for convergence or not. As these are often used as nonlinear
2111: preconditioners, monitoring the norm of their error is not a useful enterprise within
2112: their solution.
2114: Level: developer
2116: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2117: @*/
2118: PetscErrorCode SNESSetFunctionType(SNES snes, SNESFunctionType type)
2119: {
2122: snes->functype = type;
2123: return(0);
2124: }
2126: /*@C
2127: SNESGetFunctionType - Gets the SNESNormSchedule used in convergence and monitoring
2128: of the SNES method.
2130: Logically Collective on SNES
2132: Input Parameters:
2133: + snes - the SNES context
2134: - normschedule - the type of the norm used
2136: Level: advanced
2138: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2139: @*/
2140: PetscErrorCode SNESGetFunctionType(SNES snes, SNESFunctionType *type)
2141: {
2144: *type = snes->functype;
2145: return(0);
2146: }
2148: /*MC
2149: SNESNGSFunction - function used to convey a Gauss-Seidel sweep on the nonlinear function
2151: Synopsis:
2152: #include <petscsnes.h>
2153: $ SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);
2155: Collective on snes
2157: Input Parameters:
2158: + X - solution vector
2159: . B - RHS vector
2160: - ctx - optional user-defined Gauss-Seidel context
2162: Output Parameter:
2163: . X - solution vector
2165: Level: intermediate
2167: .seealso: SNESSetNGS(), SNESGetNGS()
2168: M*/
2170: /*@C
2171: SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
2172: use with composed nonlinear solvers.
2174: Input Parameters:
2175: + snes - the SNES context
2176: . f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
2177: - ctx - [optional] user-defined context for private data for the
2178: smoother evaluation routine (may be NULL)
2180: Notes:
2181: The NGS routines are used by the composed nonlinear solver to generate
2182: a problem appropriate update to the solution, particularly FAS.
2184: Level: intermediate
2186: .seealso: SNESGetFunction(), SNESComputeNGS()
2187: @*/
2188: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
2189: {
2191: DM dm;
2195: SNESGetDM(snes,&dm);
2196: DMSNESSetNGS(dm,f,ctx);
2197: return(0);
2198: }
2200: /*
2201: This is used for -snes_mf_operator; it uses a duplicate of snes->jacobian_pre because snes->jacobian_pre cannot be
2202: changed during the KSPSolve()
2203: */
2204: PetscErrorCode SNESPicardComputeMFFunction(SNES snes,Vec x,Vec f,void *ctx)
2205: {
2207: DM dm;
2208: DMSNES sdm;
2211: SNESGetDM(snes,&dm);
2212: DMGetDMSNES(dm,&sdm);
2213: if (!sdm->ops->computepjacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard Jacobian.");
2214: /* A(x)*x - b(x) */
2215: if (sdm->ops->computepfunction) {
2216: PetscStackPush("SNES Picard user function");
2217: (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
2218: PetscStackPop;
2219: VecScale(f,-1.0);
2220: if (!snes->picard) {
2221: /* Cannot share nonzero pattern because of the possible use of SNESComputeJacobianDefault() */
2222: MatDuplicate(snes->jacobian_pre,MAT_DO_NOT_COPY_VALUES,&snes->picard);
2223: }
2224: PetscStackPush("SNES Picard user Jacobian");
2225: (*sdm->ops->computepjacobian)(snes,x,snes->picard,snes->picard,sdm->pctx);
2226: PetscStackPop;
2227: MatMultAdd(snes->picard,x,f,f);
2228: } else {
2229: PetscStackPush("SNES Picard user Jacobian");
2230: (*sdm->ops->computepjacobian)(snes,x,snes->picard,snes->picard,sdm->pctx);
2231: PetscStackPop;
2232: MatMult(snes->picard,x,f);
2233: }
2234: return(0);
2235: }
2237: PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
2238: {
2240: DM dm;
2241: DMSNES sdm;
2244: SNESGetDM(snes,&dm);
2245: DMGetDMSNES(dm,&sdm);
2246: if (!sdm->ops->computepjacobian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard Jacobian.");
2247: /* A(x)*x - b(x) */
2248: if (sdm->ops->computepfunction) {
2249: PetscStackPush("SNES Picard user function");
2250: (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
2251: PetscStackPop;
2252: VecScale(f,-1.0);
2253: PetscStackPush("SNES Picard user Jacobian");
2254: (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
2255: PetscStackPop;
2256: MatMultAdd(snes->jacobian_pre,x,f,f);
2257: } else {
2258: PetscStackPush("SNES Picard user Jacobian");
2259: (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
2260: PetscStackPop;
2261: MatMult(snes->jacobian_pre,x,f);
2262: }
2263: return(0);
2264: }
2266: PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
2267: {
2270: /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
2271: /* must assembly if matrix-free to get the last SNES solution */
2272: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
2273: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
2274: return(0);
2275: }
2277: /*@C
2278: SNESSetPicard - Use SNES to solve the system A(x) x = bp(x) + b via a Picard type iteration (Picard linearization)
2280: Logically Collective on SNES
2282: Input Parameters:
2283: + snes - the SNES context
2284: . r - vector to store function values, may be NULL
2285: . bp - function evaluation routine, may be NULL
2286: . Amat - matrix with which A(x) x - bp(x) - b is to be computed
2287: . Pmat - matrix from which preconditioner is computed (usually the same as Amat)
2288: . J - function to compute matrix values, see SNESJacobianFunction() for details on its calling sequence
2289: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
2291: Notes:
2292: It is often better to provide the nonlinear function F() and some approximation to its Jacobian directly and use
2293: 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.
2295: One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both
2297: $ Solves the equation A(x) x = bp(x) - b via the defect correction algorithm A(x^{n}) (x^{n+1} - x^{n}) = bp(x^{n}) + b - A(x^{n})x^{n}
2298: $ Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = bp(x^{n}) + b iteration.
2300: Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner.
2302: We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then
2303: the direct Picard iteration A(x^n) x^{n+1} = bp(x^n) + b
2305: 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
2306: 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
2307: different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).
2309: When used with -snes_mf_operator this will run matrix-free Newton's method where the matrix-vector product is of the true Jacobian of A(x)x - bp(x) -b.
2311: When used with -snes_fd this will compute the true Jacobian (very slowly one column at at time) and thus represent Newton's method.
2313: When used with -snes_fd_coloring this will compute the Jacobian via coloring and thus represent a faster implementation of Newton's method. But the
2314: the nonzero structure of the Jacobian is, in general larger than that of the Picard matrix A so you must provide in A the needed nonzero structure for the correct
2315: coloring. When using DMDA this may mean creating the matrix A with DMCreateMatrix() using a wider stencil than strictly needed for A or with a DMDA_STENCIL_BOX.
2316: See the commment in src/snes/tutorials/ex15.c.
2318: Level: intermediate
2320: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
2321: @*/
2322: PetscErrorCode SNESSetPicard(SNES snes,Vec r,PetscErrorCode (*bp)(SNES,Vec,Vec,void*),Mat Amat, Mat Pmat, PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2323: {
2325: DM dm;
2329: SNESGetDM(snes, &dm);
2330: DMSNESSetPicard(dm,bp,J,ctx);
2331: DMSNESSetMFFunction(dm,SNESPicardComputeMFFunction,ctx);
2332: SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
2333: SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
2334: return(0);
2335: }
2337: /*@C
2338: SNESGetPicard - Returns the context for the Picard iteration
2340: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
2342: Input Parameter:
2343: . snes - the SNES context
2345: Output Parameters:
2346: + r - the function (or NULL)
2347: . f - the function (or NULL); see SNESFunction for calling sequence details
2348: . Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
2349: . Pmat - the matrix from which the preconditioner will be constructed (or NULL)
2350: . J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
2351: - ctx - the function context (or NULL)
2353: Level: advanced
2355: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
2356: @*/
2357: 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)
2358: {
2360: DM dm;
2364: SNESGetFunction(snes,r,NULL,NULL);
2365: SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
2366: SNESGetDM(snes,&dm);
2367: DMSNESGetPicard(dm,f,J,ctx);
2368: return(0);
2369: }
2371: /*@C
2372: SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem
2374: Logically Collective on SNES
2376: Input Parameters:
2377: + snes - the SNES context
2378: . func - function evaluation routine
2379: - ctx - [optional] user-defined context for private data for the
2380: function evaluation routine (may be NULL)
2382: Calling sequence of func:
2383: $ func (SNES snes,Vec x,void *ctx);
2385: . f - function vector
2386: - ctx - optional user-defined function context
2388: Level: intermediate
2390: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
2391: @*/
2392: PetscErrorCode SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
2393: {
2396: if (func) snes->ops->computeinitialguess = func;
2397: if (ctx) snes->initialguessP = ctx;
2398: return(0);
2399: }
2401: /* --------------------------------------------------------------- */
2402: /*@C
2403: SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
2404: it assumes a zero right hand side.
2406: Logically Collective on SNES
2408: Input Parameter:
2409: . snes - the SNES context
2411: Output Parameter:
2412: . rhs - the right hand side vector or NULL if the right hand side vector is null
2414: Level: intermediate
2416: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
2417: @*/
2418: PetscErrorCode SNESGetRhs(SNES snes,Vec *rhs)
2419: {
2423: *rhs = snes->vec_rhs;
2424: return(0);
2425: }
2427: /*@
2428: SNESComputeFunction - Calls the function that has been set with SNESSetFunction().
2430: Collective on SNES
2432: Input Parameters:
2433: + snes - the SNES context
2434: - x - input vector
2436: Output Parameter:
2437: . y - function vector, as set by SNESSetFunction()
2439: Notes:
2440: SNESComputeFunction() is typically used within nonlinear solvers
2441: implementations, so users would not generally call this routine themselves.
2443: Level: developer
2445: .seealso: SNESSetFunction(), SNESGetFunction(), SNESComputeMFFunction()
2446: @*/
2447: PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y)
2448: {
2450: DM dm;
2451: DMSNES sdm;
2459: VecValidValues(x,2,PETSC_TRUE);
2461: SNESGetDM(snes,&dm);
2462: DMGetDMSNES(dm,&sdm);
2463: if (sdm->ops->computefunction) {
2464: if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2465: PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2466: }
2467: VecLockReadPush(x);
2468: PetscStackPush("SNES user function");
2469: /* ensure domainerror is false prior to computefunction evaluation (may not have been reset) */
2470: snes->domainerror = PETSC_FALSE;
2471: (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2472: PetscStackPop;
2473: VecLockReadPop(x);
2474: if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2475: PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2476: }
2477: } else if (snes->vec_rhs) {
2478: MatMult(snes->jacobian, x, y);
2479: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2480: if (snes->vec_rhs) {
2481: VecAXPY(y,-1.0,snes->vec_rhs);
2482: }
2483: snes->nfuncs++;
2484: /*
2485: domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2486: propagate the value to all processes
2487: */
2488: if (snes->domainerror) {
2489: VecSetInf(y);
2490: }
2491: return(0);
2492: }
2494: /*@
2495: SNESComputeMFFunction - Calls the function that has been set with SNESSetMFFunction().
2497: Collective on SNES
2499: Input Parameters:
2500: + snes - the SNES context
2501: - x - input vector
2503: Output Parameter:
2504: . y - function vector, as set by SNESSetMFFunction()
2506: Notes:
2507: SNESComputeMFFunction() is used within the matrix vector products called by the matrix created with MatCreateSNESMF()
2508: so users would not generally call this routine themselves.
2510: Since this function is intended for use with finite differencing it does not subtract the right hand side vector provided with SNESSolve()
2511: while SNESComputeFunction() does. As such, this routine cannot be used with MatMFFDSetBase() with a provided F function value even if it applies the
2512: same function as SNESComputeFunction() if a SNESSolve() right hand side vector is use because the two functions difference would include this right hand side function.
2514: Level: developer
2516: .seealso: SNESSetFunction(), SNESGetFunction(), SNESComputeFunction(), MatCreateSNESMF
2517: @*/
2518: PetscErrorCode SNESComputeMFFunction(SNES snes,Vec x,Vec y)
2519: {
2521: DM dm;
2522: DMSNES sdm;
2530: VecValidValues(x,2,PETSC_TRUE);
2532: SNESGetDM(snes,&dm);
2533: DMGetDMSNES(dm,&sdm);
2534: PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2535: VecLockReadPush(x);
2536: PetscStackPush("SNES user function");
2537: /* ensure domainerror is false prior to computefunction evaluation (may not have been reset) */
2538: snes->domainerror = PETSC_FALSE;
2539: (*sdm->ops->computemffunction)(snes,x,y,sdm->mffunctionctx);
2540: PetscStackPop;
2541: VecLockReadPop(x);
2542: PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2543: snes->nfuncs++;
2544: /*
2545: domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2546: propagate the value to all processes
2547: */
2548: if (snes->domainerror) {
2549: VecSetInf(y);
2550: }
2551: return(0);
2552: }
2554: /*@
2555: SNESComputeNGS - Calls the Gauss-Seidel function that has been set with SNESSetNGS().
2557: Collective on SNES
2559: Input Parameters:
2560: + snes - the SNES context
2561: . x - input vector
2562: - b - rhs vector
2564: Output Parameter:
2565: . x - new solution vector
2567: Notes:
2568: SNESComputeNGS() is typically used within composed nonlinear solver
2569: implementations, so most users would not generally call this routine
2570: themselves.
2572: Level: developer
2574: .seealso: SNESSetNGS(), SNESComputeFunction()
2575: @*/
2576: PetscErrorCode SNESComputeNGS(SNES snes,Vec b,Vec x)
2577: {
2579: DM dm;
2580: DMSNES sdm;
2588: if (b) {VecValidValues(b,2,PETSC_TRUE);}
2589: PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2590: SNESGetDM(snes,&dm);
2591: DMGetDMSNES(dm,&sdm);
2592: if (sdm->ops->computegs) {
2593: if (b) {VecLockReadPush(b);}
2594: PetscStackPush("SNES user NGS");
2595: (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2596: PetscStackPop;
2597: if (b) {VecLockReadPop(b);}
2598: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2599: PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2600: return(0);
2601: }
2603: PetscErrorCode SNESTestJacobian(SNES snes)
2604: {
2605: Mat A,B,C,D,jacobian;
2606: Vec x = snes->vec_sol,f = snes->vec_func;
2607: PetscErrorCode ierr;
2608: PetscReal nrm,gnorm;
2609: PetscReal threshold = 1.e-5;
2610: MatType mattype;
2611: PetscInt m,n,M,N;
2612: void *functx;
2613: PetscBool complete_print = PETSC_FALSE,threshold_print = PETSC_FALSE,test = PETSC_FALSE,flg,istranspose;
2614: PetscViewer viewer,mviewer;
2615: MPI_Comm comm;
2616: PetscInt tabs;
2617: static PetscBool directionsprinted = PETSC_FALSE;
2618: PetscViewerFormat format;
2621: PetscObjectOptionsBegin((PetscObject)snes);
2622: PetscOptionsName("-snes_test_jacobian","Compare hand-coded and finite difference Jacobians","None",&test);
2623: PetscOptionsReal("-snes_test_jacobian", "Threshold for element difference between hand-coded and finite difference being meaningful", "None", threshold, &threshold,NULL);
2624: PetscOptionsViewer("-snes_test_jacobian_view","View difference between hand-coded and finite difference Jacobians element entries","None",&mviewer,&format,&complete_print);
2625: if (!complete_print) {
2626: PetscOptionsDeprecated("-snes_test_jacobian_display","-snes_test_jacobian_view","3.13",NULL);
2627: PetscOptionsViewer("-snes_test_jacobian_display","Display difference between hand-coded and finite difference Jacobians","None",&mviewer,&format,&complete_print);
2628: }
2629: /* for compatibility with PETSc 3.9 and older. */
2630: PetscOptionsDeprecated("-snes_test_jacobian_display_threshold","-snes_test_jacobian","3.13","-snes_test_jacobian accepts an optional threshold (since v3.10)");
2631: PetscOptionsReal("-snes_test_jacobian_display_threshold", "Display difference between hand-coded and finite difference Jacobians which exceed input threshold", "None", threshold, &threshold, &threshold_print);
2632: PetscOptionsEnd();
2633: if (!test) return(0);
2635: PetscObjectGetComm((PetscObject)snes,&comm);
2636: PetscViewerASCIIGetStdout(comm,&viewer);
2637: PetscViewerASCIIGetTab(viewer, &tabs);
2638: PetscViewerASCIISetTab(viewer, ((PetscObject)snes)->tablevel);
2639: PetscViewerASCIIPrintf(viewer," ---------- Testing Jacobian -------------\n");
2640: if (!complete_print && !directionsprinted) {
2641: PetscViewerASCIIPrintf(viewer," Run with -snes_test_jacobian_view and optionally -snes_test_jacobian <threshold> to show difference\n");
2642: PetscViewerASCIIPrintf(viewer," of hand-coded and finite difference Jacobian entries greater than <threshold>.\n");
2643: }
2644: if (!directionsprinted) {
2645: PetscViewerASCIIPrintf(viewer," Testing hand-coded Jacobian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");
2646: PetscViewerASCIIPrintf(viewer," O(1.e-8), the hand-coded Jacobian is probably correct.\n");
2647: directionsprinted = PETSC_TRUE;
2648: }
2649: if (complete_print) {
2650: PetscViewerPushFormat(mviewer,format);
2651: }
2653: PetscObjectTypeCompare((PetscObject)snes->jacobian,MATMFFD,&flg);
2654: if (!flg) jacobian = snes->jacobian;
2655: else jacobian = snes->jacobian_pre;
2657: if (!x) {
2658: MatCreateVecs(jacobian, &x, NULL);
2659: } else {
2660: PetscObjectReference((PetscObject) x);
2661: }
2662: if (!f) {
2663: VecDuplicate(x, &f);
2664: } else {
2665: PetscObjectReference((PetscObject) f);
2666: }
2667: /* evaluate the function at this point because SNESComputeJacobianDefault() assumes that the function has been evaluated and put into snes->vec_func */
2668: SNESComputeFunction(snes,x,f);
2669: VecDestroy(&f);
2670: PetscObjectTypeCompare((PetscObject)snes,SNESKSPTRANSPOSEONLY,&istranspose);
2671: while (jacobian) {
2672: Mat JT = NULL, Jsave = NULL;
2674: if (istranspose) {
2675: MatCreateTranspose(jacobian,&JT);
2676: Jsave = jacobian;
2677: jacobian = JT;
2678: }
2679: PetscObjectBaseTypeCompareAny((PetscObject)jacobian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPISBAIJ,"");
2680: if (flg) {
2681: A = jacobian;
2682: PetscObjectReference((PetscObject)A);
2683: } else {
2684: MatComputeOperator(jacobian,MATAIJ,&A);
2685: }
2687: MatGetType(A,&mattype);
2688: MatGetSize(A,&M,&N);
2689: MatGetLocalSize(A,&m,&n);
2690: MatCreate(PetscObjectComm((PetscObject)A),&B);
2691: MatSetType(B,mattype);
2692: MatSetSizes(B,m,n,M,N);
2693: MatSetBlockSizesFromMats(B,A,A);
2694: MatSetUp(B);
2695: MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
2697: SNESGetFunction(snes,NULL,NULL,&functx);
2698: SNESComputeJacobianDefault(snes,x,B,B,functx);
2700: MatDuplicate(B,MAT_COPY_VALUES,&D);
2701: MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2702: MatNorm(D,NORM_FROBENIUS,&nrm);
2703: MatNorm(A,NORM_FROBENIUS,&gnorm);
2704: MatDestroy(&D);
2705: if (!gnorm) gnorm = 1; /* just in case */
2706: PetscViewerASCIIPrintf(viewer," ||J - Jfd||_F/||J||_F = %g, ||J - Jfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);
2708: if (complete_print) {
2709: PetscViewerASCIIPrintf(viewer," Hand-coded Jacobian ----------\n");
2710: MatView(A,mviewer);
2711: PetscViewerASCIIPrintf(viewer," Finite difference Jacobian ----------\n");
2712: MatView(B,mviewer);
2713: }
2715: if (threshold_print || complete_print) {
2716: PetscInt Istart, Iend, *ccols, bncols, cncols, j, row;
2717: PetscScalar *cvals;
2718: const PetscInt *bcols;
2719: const PetscScalar *bvals;
2721: MatCreate(PetscObjectComm((PetscObject)A),&C);
2722: MatSetType(C,mattype);
2723: MatSetSizes(C,m,n,M,N);
2724: MatSetBlockSizesFromMats(C,A,A);
2725: MatSetUp(C);
2726: MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
2728: MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2729: MatGetOwnershipRange(B,&Istart,&Iend);
2731: for (row = Istart; row < Iend; row++) {
2732: MatGetRow(B,row,&bncols,&bcols,&bvals);
2733: PetscMalloc2(bncols,&ccols,bncols,&cvals);
2734: for (j = 0, cncols = 0; j < bncols; j++) {
2735: if (PetscAbsScalar(bvals[j]) > threshold) {
2736: ccols[cncols] = bcols[j];
2737: cvals[cncols] = bvals[j];
2738: cncols += 1;
2739: }
2740: }
2741: if (cncols) {
2742: MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);
2743: }
2744: MatRestoreRow(B,row,&bncols,&bcols,&bvals);
2745: PetscFree2(ccols,cvals);
2746: }
2747: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2748: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2749: PetscViewerASCIIPrintf(viewer," Hand-coded minus finite-difference Jacobian with tolerance %g ----------\n",(double)threshold);
2750: MatView(C,complete_print ? mviewer : viewer);
2751: MatDestroy(&C);
2752: }
2753: MatDestroy(&A);
2754: MatDestroy(&B);
2755: MatDestroy(&JT);
2756: if (Jsave) jacobian = Jsave;
2757: if (jacobian != snes->jacobian_pre) {
2758: jacobian = snes->jacobian_pre;
2759: PetscViewerASCIIPrintf(viewer," ---------- Testing Jacobian for preconditioner -------------\n");
2760: }
2761: else jacobian = NULL;
2762: }
2763: VecDestroy(&x);
2764: if (complete_print) {
2765: PetscViewerPopFormat(mviewer);
2766: }
2767: if (mviewer) { PetscViewerDestroy(&mviewer); }
2768: PetscViewerASCIISetTab(viewer,tabs);
2769: return(0);
2770: }
2772: /*@
2773: SNESComputeJacobian - Computes the Jacobian matrix that has been set with SNESSetJacobian().
2775: Collective on SNES
2777: Input Parameters:
2778: + snes - the SNES context
2779: - x - input vector
2781: Output Parameters:
2782: + A - Jacobian matrix
2783: - B - optional preconditioning matrix
2785: Options Database Keys:
2786: + -snes_lag_preconditioner <lag>
2787: . -snes_lag_jacobian <lag>
2788: . -snes_test_jacobian <optional threshold> - compare the user provided Jacobian with one compute via finite differences to check for errors. If a threshold is given, display only those entries whose difference is greater than the threshold.
2789: . -snes_test_jacobian_view - display the user provided Jacobian, the finite difference Jacobian and the difference between them to help users detect the location of errors in the user provided Jacobian
2790: . -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2791: . -snes_compare_explicit_draw - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2792: . -snes_compare_explicit_contour - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2793: . -snes_compare_operator - Make the comparison options above use the operator instead of the preconditioning matrix
2794: . -snes_compare_coloring - Compute the finite difference Jacobian using coloring and display norms of difference
2795: . -snes_compare_coloring_display - Compute the finite difference Jacobian using coloring and display verbose differences
2796: . -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2797: . -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2798: . -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2799: . -snes_compare_coloring_draw - Compute the finite difference Jacobian using coloring and draw differences
2800: - -snes_compare_coloring_draw_contour - Compute the finite difference Jacobian using coloring and show contours of matrices and differences
2802: Notes:
2803: Most users should not need to explicitly call this routine, as it
2804: is used internally within the nonlinear solvers.
2806: Developer Notes:
2807: This has duplicative ways of checking the accuracy of the user provided Jacobian (see the options above). This is for historical reasons, the routine SNESTestJacobian() use to used
2808: for with the SNESType of test that has been removed.
2810: Level: developer
2812: .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2813: @*/
2814: PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2815: {
2817: PetscBool flag;
2818: DM dm;
2819: DMSNES sdm;
2820: KSP ksp;
2826: VecValidValues(X,2,PETSC_TRUE);
2827: SNESGetDM(snes,&dm);
2828: DMGetDMSNES(dm,&sdm);
2830: if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_USER,"Must call SNESSetJacobian(), DMSNESSetJacobian(), DMDASNESSetJacobianLocal(), etc");
2832: /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
2834: if (snes->lagjacobian == -2) {
2835: snes->lagjacobian = -1;
2837: PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2838: } else if (snes->lagjacobian == -1) {
2839: PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2840: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2841: if (flag) {
2842: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2843: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2844: }
2845: return(0);
2846: } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2847: PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2848: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2849: if (flag) {
2850: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2851: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2852: }
2853: return(0);
2854: }
2855: if (snes->npc && snes->npcside== PC_LEFT) {
2856: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2857: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2858: return(0);
2859: }
2861: PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2862: VecLockReadPush(X);
2863: PetscStackPush("SNES user Jacobian function");
2864: (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2865: PetscStackPop;
2866: VecLockReadPop(X);
2867: PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);
2869: /* attach latest linearization point to the preconditioning matrix */
2870: PetscObjectCompose((PetscObject)B,"__SNES_latest_X",(PetscObject)X);
2872: /* the next line ensures that snes->ksp exists */
2873: SNESGetKSP(snes,&ksp);
2874: if (snes->lagpreconditioner == -2) {
2875: PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2876: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2877: snes->lagpreconditioner = -1;
2878: } else if (snes->lagpreconditioner == -1) {
2879: PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2880: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2881: } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2882: PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2883: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2884: } else {
2885: PetscInfo(snes,"Rebuilding preconditioner\n");
2886: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2887: }
2889: SNESTestJacobian(snes);
2890: /* make sure user returned a correct Jacobian and preconditioner */
2893: {
2894: PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
2895: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit",NULL,NULL,&flag);
2896: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",NULL,NULL,&flag_draw);
2897: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",NULL,NULL,&flag_contour);
2898: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_operator",NULL,NULL,&flag_operator);
2899: if (flag || flag_draw || flag_contour) {
2900: Mat Bexp_mine = NULL,Bexp,FDexp;
2901: PetscViewer vdraw,vstdout;
2902: PetscBool flg;
2903: if (flag_operator) {
2904: MatComputeOperator(A,MATAIJ,&Bexp_mine);
2905: Bexp = Bexp_mine;
2906: } else {
2907: /* See if the preconditioning matrix can be viewed and added directly */
2908: PetscObjectBaseTypeCompareAny((PetscObject)B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2909: if (flg) Bexp = B;
2910: else {
2911: /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
2912: MatComputeOperator(B,MATAIJ,&Bexp_mine);
2913: Bexp = Bexp_mine;
2914: }
2915: }
2916: MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);
2917: SNESComputeJacobianDefault(snes,X,FDexp,FDexp,NULL);
2918: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2919: if (flag_draw || flag_contour) {
2920: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2921: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2922: } else vdraw = NULL;
2923: PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator ? "Jacobian" : "preconditioning Jacobian");
2924: if (flag) {MatView(Bexp,vstdout);}
2925: if (vdraw) {MatView(Bexp,vdraw);}
2926: PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");
2927: if (flag) {MatView(FDexp,vstdout);}
2928: if (vdraw) {MatView(FDexp,vdraw);}
2929: MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);
2930: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");
2931: if (flag) {MatView(FDexp,vstdout);}
2932: if (vdraw) { /* Always use contour for the difference */
2933: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2934: MatView(FDexp,vdraw);
2935: PetscViewerPopFormat(vdraw);
2936: }
2937: if (flag_contour) {PetscViewerPopFormat(vdraw);}
2938: PetscViewerDestroy(&vdraw);
2939: MatDestroy(&Bexp_mine);
2940: MatDestroy(&FDexp);
2941: }
2942: }
2943: {
2944: PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
2945: PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
2946: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring",NULL,NULL,&flag);
2947: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_display",NULL,NULL,&flag_display);
2948: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",NULL,NULL,&flag_draw);
2949: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",NULL,NULL,&flag_contour);
2950: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",NULL,NULL,&flag_threshold);
2951: if (flag_threshold) {
2952: PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,NULL);
2953: PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,NULL);
2954: }
2955: if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
2956: Mat Bfd;
2957: PetscViewer vdraw,vstdout;
2958: MatColoring coloring;
2959: ISColoring iscoloring;
2960: MatFDColoring matfdcoloring;
2961: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2962: void *funcctx;
2963: PetscReal norm1,norm2,normmax;
2965: MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2966: MatColoringCreate(Bfd,&coloring);
2967: MatColoringSetType(coloring,MATCOLORINGSL);
2968: MatColoringSetFromOptions(coloring);
2969: MatColoringApply(coloring,&iscoloring);
2970: MatColoringDestroy(&coloring);
2971: MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2972: MatFDColoringSetFromOptions(matfdcoloring);
2973: MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2974: ISColoringDestroy(&iscoloring);
2976: /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2977: SNESGetFunction(snes,NULL,&func,&funcctx);
2978: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2979: PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2980: PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2981: MatFDColoringSetFromOptions(matfdcoloring);
2982: MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2983: MatFDColoringDestroy(&matfdcoloring);
2985: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2986: if (flag_draw || flag_contour) {
2987: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2988: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2989: } else vdraw = NULL;
2990: PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2991: if (flag_display) {MatView(B,vstdout);}
2992: if (vdraw) {MatView(B,vdraw);}
2993: PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2994: if (flag_display) {MatView(Bfd,vstdout);}
2995: if (vdraw) {MatView(Bfd,vdraw);}
2996: MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2997: MatNorm(Bfd,NORM_1,&norm1);
2998: MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2999: MatNorm(Bfd,NORM_MAX,&normmax);
3000: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
3001: if (flag_display) {MatView(Bfd,vstdout);}
3002: if (vdraw) { /* Always use contour for the difference */
3003: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
3004: MatView(Bfd,vdraw);
3005: PetscViewerPopFormat(vdraw);
3006: }
3007: if (flag_contour) {PetscViewerPopFormat(vdraw);}
3009: if (flag_threshold) {
3010: PetscInt bs,rstart,rend,i;
3011: MatGetBlockSize(B,&bs);
3012: MatGetOwnershipRange(B,&rstart,&rend);
3013: for (i=rstart; i<rend; i++) {
3014: const PetscScalar *ba,*ca;
3015: const PetscInt *bj,*cj;
3016: PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
3017: PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0;
3018: MatGetRow(B,i,&bn,&bj,&ba);
3019: MatGetRow(Bfd,i,&cn,&cj,&ca);
3020: if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
3021: for (j=0; j<bn; j++) {
3022: PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
3023: if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
3024: maxentrycol = bj[j];
3025: maxentry = PetscRealPart(ba[j]);
3026: }
3027: if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
3028: maxdiffcol = bj[j];
3029: maxdiff = PetscRealPart(ca[j]);
3030: }
3031: if (rdiff > maxrdiff) {
3032: maxrdiffcol = bj[j];
3033: maxrdiff = rdiff;
3034: }
3035: }
3036: if (maxrdiff > 1) {
3037: 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);
3038: for (j=0; j<bn; j++) {
3039: PetscReal rdiff;
3040: rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
3041: if (rdiff > 1) {
3042: PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
3043: }
3044: }
3045: PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
3046: }
3047: MatRestoreRow(B,i,&bn,&bj,&ba);
3048: MatRestoreRow(Bfd,i,&cn,&cj,&ca);
3049: }
3050: }
3051: PetscViewerDestroy(&vdraw);
3052: MatDestroy(&Bfd);
3053: }
3054: }
3055: return(0);
3056: }
3058: /*MC
3059: SNESJacobianFunction - Function used to convey the nonlinear Jacobian of the function to be solved by SNES
3061: Synopsis:
3062: #include "petscsnes.h"
3063: PetscErrorCode SNESJacobianFunction(SNES snes,Vec x,Mat Amat,Mat Pmat,void *ctx);
3065: Collective on snes
3067: Input Parameters:
3068: + x - input vector, the Jacobian is to be computed at this value
3069: - ctx - [optional] user-defined Jacobian context
3071: Output Parameters:
3072: + Amat - the matrix that defines the (approximate) Jacobian
3073: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
3075: Level: intermediate
3077: .seealso: SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
3078: M*/
3080: /*@C
3081: SNESSetJacobian - Sets the function to compute Jacobian as well as the
3082: location to store the matrix.
3084: Logically Collective on SNES
3086: Input Parameters:
3087: + snes - the SNES context
3088: . Amat - the matrix that defines the (approximate) Jacobian
3089: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
3090: . J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
3091: - ctx - [optional] user-defined context for private data for the
3092: Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)
3094: Notes:
3095: If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
3096: each matrix.
3098: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
3099: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
3101: If using SNESComputeJacobianDefaultColor() to assemble a Jacobian, the ctx argument
3102: must be a MatFDColoring.
3104: Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian. One common
3105: example is to use the "Picard linearization" which only differentiates through the highest order parts of each term.
3107: Level: beginner
3109: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J,
3110: SNESSetPicard(), SNESJacobianFunction
3111: @*/
3112: PetscErrorCode SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
3113: {
3115: DM dm;
3123: SNESGetDM(snes,&dm);
3124: DMSNESSetJacobian(dm,J,ctx);
3125: if (Amat) {
3126: PetscObjectReference((PetscObject)Amat);
3127: MatDestroy(&snes->jacobian);
3129: snes->jacobian = Amat;
3130: }
3131: if (Pmat) {
3132: PetscObjectReference((PetscObject)Pmat);
3133: MatDestroy(&snes->jacobian_pre);
3135: snes->jacobian_pre = Pmat;
3136: }
3137: return(0);
3138: }
3140: /*@C
3141: SNESGetJacobian - Returns the Jacobian matrix and optionally the user
3142: provided context for evaluating the Jacobian.
3144: Not Collective, but Mat object will be parallel if SNES object is
3146: Input Parameter:
3147: . snes - the nonlinear solver context
3149: Output Parameters:
3150: + Amat - location to stash (approximate) Jacobian matrix (or NULL)
3151: . Pmat - location to stash matrix used to compute the preconditioner (or NULL)
3152: . J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
3153: - ctx - location to stash Jacobian ctx (or NULL)
3155: Level: advanced
3157: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
3158: @*/
3159: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
3160: {
3162: DM dm;
3163: DMSNES sdm;
3167: if (Amat) *Amat = snes->jacobian;
3168: if (Pmat) *Pmat = snes->jacobian_pre;
3169: SNESGetDM(snes,&dm);
3170: DMGetDMSNES(dm,&sdm);
3171: if (J) *J = sdm->ops->computejacobian;
3172: if (ctx) *ctx = sdm->jacobianctx;
3173: return(0);
3174: }
3176: static PetscErrorCode SNESSetDefaultComputeJacobian(SNES snes)
3177: {
3179: DM dm;
3180: DMSNES sdm;
3183: SNESGetDM(snes,&dm);
3184: DMGetDMSNES(dm,&sdm);
3185: if (!sdm->ops->computejacobian && snes->jacobian_pre) {
3186: DM dm;
3187: PetscBool isdense,ismf;
3189: SNESGetDM(snes,&dm);
3190: PetscObjectTypeCompareAny((PetscObject)snes->jacobian_pre,&isdense,MATSEQDENSE,MATMPIDENSE,MATDENSE,NULL);
3191: PetscObjectTypeCompareAny((PetscObject)snes->jacobian_pre,&ismf,MATMFFD,MATSHELL,NULL);
3192: if (isdense) {
3193: DMSNESSetJacobian(dm,SNESComputeJacobianDefault,NULL);
3194: } else if (!ismf) {
3195: DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
3196: }
3197: }
3198: return(0);
3199: }
3201: /*@
3202: SNESSetUp - Sets up the internal data structures for the later use
3203: of a nonlinear solver.
3205: Collective on SNES
3207: Input Parameters:
3208: . snes - the SNES context
3210: Notes:
3211: For basic use of the SNES solvers the user need not explicitly call
3212: SNESSetUp(), since these actions will automatically occur during
3213: the call to SNESSolve(). However, if one wishes to control this
3214: phase separately, SNESSetUp() should be called after SNESCreate()
3215: and optional routines of the form SNESSetXXX(), but before SNESSolve().
3217: Level: advanced
3219: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
3220: @*/
3221: PetscErrorCode SNESSetUp(SNES snes)
3222: {
3224: DM dm;
3225: DMSNES sdm;
3226: SNESLineSearch linesearch, pclinesearch;
3227: void *lsprectx,*lspostctx;
3228: PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
3229: PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
3230: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
3231: Vec f,fpc;
3232: void *funcctx;
3233: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
3234: void *jacctx,*appctx;
3235: Mat j,jpre;
3239: if (snes->setupcalled) return(0);
3240: PetscLogEventBegin(SNES_Setup,snes,0,0,0);
3242: if (!((PetscObject)snes)->type_name) {
3243: SNESSetType(snes,SNESNEWTONLS);
3244: }
3246: SNESGetFunction(snes,&snes->vec_func,NULL,NULL);
3248: SNESGetDM(snes,&dm);
3249: DMGetDMSNES(dm,&sdm);
3250: if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
3251: SNESSetDefaultComputeJacobian(snes);
3253: if (!snes->vec_func) {
3254: DMCreateGlobalVector(dm,&snes->vec_func);
3255: }
3257: if (!snes->ksp) {
3258: SNESGetKSP(snes, &snes->ksp);
3259: }
3261: if (snes->linesearch) {
3262: SNESGetLineSearch(snes, &snes->linesearch);
3263: SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);
3264: }
3266: if (snes->npc && (snes->npcside== PC_LEFT)) {
3267: snes->mf = PETSC_TRUE;
3268: snes->mf_operator = PETSC_FALSE;
3269: }
3271: if (snes->npc) {
3272: /* copy the DM over */
3273: SNESGetDM(snes,&dm);
3274: SNESSetDM(snes->npc,dm);
3276: SNESGetFunction(snes,&f,&func,&funcctx);
3277: VecDuplicate(f,&fpc);
3278: SNESSetFunction(snes->npc,fpc,func,funcctx);
3279: SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
3280: SNESSetJacobian(snes->npc,j,jpre,jac,jacctx);
3281: SNESGetApplicationContext(snes,&appctx);
3282: SNESSetApplicationContext(snes->npc,appctx);
3283: VecDestroy(&fpc);
3285: /* copy the function pointers over */
3286: PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)snes->npc);
3288: /* default to 1 iteration */
3289: SNESSetTolerances(snes->npc,0.0,0.0,0.0,1,snes->npc->max_funcs);
3290: if (snes->npcside==PC_RIGHT) {
3291: SNESSetNormSchedule(snes->npc,SNES_NORM_FINAL_ONLY);
3292: } else {
3293: SNESSetNormSchedule(snes->npc,SNES_NORM_NONE);
3294: }
3295: SNESSetFromOptions(snes->npc);
3297: /* copy the line search context over */
3298: if (snes->linesearch && snes->npc->linesearch) {
3299: SNESGetLineSearch(snes,&linesearch);
3300: SNESGetLineSearch(snes->npc,&pclinesearch);
3301: SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
3302: SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
3303: SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
3304: SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
3305: PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
3306: }
3307: }
3308: if (snes->mf) {
3309: SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
3310: }
3311: if (snes->ops->usercompute && !snes->user) {
3312: (*snes->ops->usercompute)(snes,(void**)&snes->user);
3313: }
3315: snes->jac_iter = 0;
3316: snes->pre_iter = 0;
3318: if (snes->ops->setup) {
3319: (*snes->ops->setup)(snes);
3320: }
3322: SNESSetDefaultComputeJacobian(snes);
3324: if (snes->npc && (snes->npcside== PC_LEFT)) {
3325: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3326: if (snes->linesearch) {
3327: SNESGetLineSearch(snes,&linesearch);
3328: SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
3329: }
3330: }
3331: }
3332: PetscLogEventEnd(SNES_Setup,snes,0,0,0);
3333: snes->setupcalled = PETSC_TRUE;
3334: return(0);
3335: }
3337: /*@
3338: SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats
3340: Collective on SNES
3342: Input Parameter:
3343: . snes - iterative context obtained from SNESCreate()
3345: Level: intermediate
3347: Notes:
3348: Also calls the application context destroy routine set with SNESSetComputeApplicationContext()
3350: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
3351: @*/
3352: PetscErrorCode SNESReset(SNES snes)
3353: {
3358: if (snes->ops->userdestroy && snes->user) {
3359: (*snes->ops->userdestroy)((void**)&snes->user);
3360: snes->user = NULL;
3361: }
3362: if (snes->npc) {
3363: SNESReset(snes->npc);
3364: }
3366: if (snes->ops->reset) {
3367: (*snes->ops->reset)(snes);
3368: }
3369: if (snes->ksp) {
3370: KSPReset(snes->ksp);
3371: }
3373: if (snes->linesearch) {
3374: SNESLineSearchReset(snes->linesearch);
3375: }
3377: VecDestroy(&snes->vec_rhs);
3378: VecDestroy(&snes->vec_sol);
3379: VecDestroy(&snes->vec_sol_update);
3380: VecDestroy(&snes->vec_func);
3381: MatDestroy(&snes->jacobian);
3382: MatDestroy(&snes->jacobian_pre);
3383: MatDestroy(&snes->picard);
3384: VecDestroyVecs(snes->nwork,&snes->work);
3385: VecDestroyVecs(snes->nvwork,&snes->vwork);
3387: snes->alwayscomputesfinalresidual = PETSC_FALSE;
3389: snes->nwork = snes->nvwork = 0;
3390: snes->setupcalled = PETSC_FALSE;
3391: return(0);
3392: }
3394: /*@
3395: SNESConvergedReasonViewCancel - Clears all the reasonview functions for a SNES object.
3397: Collective on SNES
3399: Input Parameter:
3400: . snes - iterative context obtained from SNESCreate()
3402: Level: intermediate
3404: .seealso: SNESCreate(), SNESDestroy(), SNESReset()
3405: @*/
3406: PetscErrorCode SNESConvergedReasonViewCancel(SNES snes)
3407: {
3409: PetscInt i;
3413: for (i=0; i<snes->numberreasonviews; i++) {
3414: if (snes->reasonviewdestroy[i]) {
3415: (*snes->reasonviewdestroy[i])(&snes->reasonviewcontext[i]);
3416: }
3417: }
3418: snes->numberreasonviews = 0;
3419: return(0);
3420: }
3422: /*@C
3423: SNESDestroy - Destroys the nonlinear solver context that was created
3424: with SNESCreate().
3426: Collective on SNES
3428: Input Parameter:
3429: . snes - the SNES context
3431: Level: beginner
3433: .seealso: SNESCreate(), SNESSolve()
3434: @*/
3435: PetscErrorCode SNESDestroy(SNES *snes)
3436: {
3440: if (!*snes) return(0);
3442: if (--((PetscObject)(*snes))->refct > 0) {*snes = NULL; return(0);}
3444: SNESReset((*snes));
3445: SNESDestroy(&(*snes)->npc);
3447: /* if memory was published with SAWs then destroy it */
3448: PetscObjectSAWsViewOff((PetscObject)*snes);
3449: if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}
3451: if ((*snes)->dm) {DMCoarsenHookRemove((*snes)->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,*snes);}
3452: DMDestroy(&(*snes)->dm);
3453: KSPDestroy(&(*snes)->ksp);
3454: SNESLineSearchDestroy(&(*snes)->linesearch);
3456: PetscFree((*snes)->kspconvctx);
3457: if ((*snes)->ops->convergeddestroy) {
3458: (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
3459: }
3460: if ((*snes)->conv_hist_alloc) {
3461: PetscFree2((*snes)->conv_hist,(*snes)->conv_hist_its);
3462: }
3463: SNESMonitorCancel((*snes));
3464: SNESConvergedReasonViewCancel((*snes));
3465: PetscHeaderDestroy(snes);
3466: return(0);
3467: }
3469: /* ----------- Routines to set solver parameters ---------- */
3471: /*@
3472: SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
3474: Logically Collective on SNES
3476: Input Parameters:
3477: + snes - the SNES context
3478: - lag - 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3479: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
3481: Options Database Keys:
3482: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3483: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3484: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3485: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag
3487: Notes:
3488: The default is 1
3489: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 or SNESSetLagPreconditionerPersists() was called
3491: SNESSetLagPreconditionerPersists() allows using the same uniform lagging (for example every second solve) across multiple solves.
3493: Level: intermediate
3495: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESSetLagPreconditionerPersists(),
3496: SNESSetLagJacobianPersists()
3498: @*/
3499: PetscErrorCode SNESSetLagPreconditioner(SNES snes,PetscInt lag)
3500: {
3503: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3504: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3506: snes->lagpreconditioner = lag;
3507: return(0);
3508: }
3510: /*@
3511: SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does
3513: Logically Collective on SNES
3515: Input Parameters:
3516: + snes - the SNES context
3517: - steps - the number of refinements to do, defaults to 0
3519: Options Database Keys:
3520: . -snes_grid_sequence <steps>
3522: Level: intermediate
3524: Notes:
3525: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
3527: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetGridSequence()
3529: @*/
3530: PetscErrorCode SNESSetGridSequence(SNES snes,PetscInt steps)
3531: {
3535: snes->gridsequence = steps;
3536: return(0);
3537: }
3539: /*@
3540: SNESGetGridSequence - gets the number of steps of grid sequencing that SNES does
3542: Logically Collective on SNES
3544: Input Parameter:
3545: . snes - the SNES context
3547: Output Parameter:
3548: . steps - the number of refinements to do, defaults to 0
3550: Options Database Keys:
3551: . -snes_grid_sequence <steps>
3553: Level: intermediate
3555: Notes:
3556: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
3558: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESSetGridSequence()
3560: @*/
3561: PetscErrorCode SNESGetGridSequence(SNES snes,PetscInt *steps)
3562: {
3565: *steps = snes->gridsequence;
3566: return(0);
3567: }
3569: /*@
3570: SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
3572: Not Collective
3574: Input Parameter:
3575: . snes - the SNES context
3577: Output Parameter:
3578: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3579: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
3581: Options Database Keys:
3582: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3583: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3584: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3585: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag
3587: Notes:
3588: The default is 1
3589: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3591: Level: intermediate
3593: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner(), SNESSetLagJacobianPersists(), SNESSetLagPreconditionerPersists()
3595: @*/
3596: PetscErrorCode SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
3597: {
3600: *lag = snes->lagpreconditioner;
3601: return(0);
3602: }
3604: /*@
3605: SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
3606: often the preconditioner is rebuilt.
3608: Logically Collective on SNES
3610: Input Parameters:
3611: + snes - the SNES context
3612: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3613: the Jacobian is built etc. -2 means rebuild at next chance but then never again
3615: Options Database Keys:
3616: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3617: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3618: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3619: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag.
3621: Notes:
3622: The default is 1
3623: The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3624: 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
3625: at the next Newton step but never again (unless it is reset to another value)
3627: Level: intermediate
3629: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobianPersists(), SNESSetLagPreconditionerPersists()
3631: @*/
3632: PetscErrorCode SNESSetLagJacobian(SNES snes,PetscInt lag)
3633: {
3636: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3637: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3639: snes->lagjacobian = lag;
3640: return(0);
3641: }
3643: /*@
3644: SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
3646: Not Collective
3648: Input Parameter:
3649: . snes - the SNES context
3651: Output Parameter:
3652: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3653: the Jacobian is built etc.
3655: Notes:
3656: The default is 1
3657: The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 or SNESSetLagJacobianPersists() was called.
3659: Level: intermediate
3661: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner(), SNESSetLagJacobianPersists(), SNESSetLagPreconditionerPersists()
3663: @*/
3664: PetscErrorCode SNESGetLagJacobian(SNES snes,PetscInt *lag)
3665: {
3668: *lag = snes->lagjacobian;
3669: return(0);
3670: }
3672: /*@
3673: SNESSetLagJacobianPersists - Set whether or not the Jacobian lagging persists through multiple solves
3675: Logically collective on SNES
3677: Input Parameters:
3678: + snes - the SNES context
3679: - flg - jacobian lagging persists if true
3681: Options Database Keys:
3682: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3683: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3684: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3685: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag
3687: Notes:
3688: This is useful both for nonlinear preconditioning, where it's appropriate to have the Jacobian be stale by
3689: several solves, and for implicit time-stepping, where Jacobian lagging in the inner nonlinear solve over several
3690: timesteps may present huge efficiency gains.
3692: Level: developer
3694: .seealso: SNESSetLagPreconditionerPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC(), SNESSetLagJacobianPersists()
3696: @*/
3697: PetscErrorCode SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
3698: {
3702: snes->lagjac_persist = flg;
3703: return(0);
3704: }
3706: /*@
3707: SNESSetLagPreconditionerPersists - Set whether or not the preconditioner lagging persists through multiple nonlinear solves
3709: Logically Collective on SNES
3711: Input Parameters:
3712: + snes - the SNES context
3713: - flg - preconditioner lagging persists if true
3715: Options Database Keys:
3716: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3717: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3718: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3719: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag
3721: Notes:
3722: This is useful both for nonlinear preconditioning, where it's appropriate to have the preconditioner be stale
3723: by several solves, and for implicit time-stepping, where preconditioner lagging in the inner nonlinear solve over
3724: several timesteps may present huge efficiency gains.
3726: Level: developer
3728: .seealso: SNESSetLagJacobianPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC(), SNESSetLagPreconditioner()
3730: @*/
3731: PetscErrorCode SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3732: {
3736: snes->lagpre_persist = flg;
3737: return(0);
3738: }
3740: /*@
3741: SNESSetForceIteration - force SNESSolve() to take at least one iteration regardless of the initial residual norm
3743: Logically Collective on SNES
3745: Input Parameters:
3746: + snes - the SNES context
3747: - force - PETSC_TRUE require at least one iteration
3749: Options Database Keys:
3750: . -snes_force_iteration <force> - Sets forcing an iteration
3752: Notes:
3753: This is used sometimes with TS to prevent TS from detecting a false steady state solution
3755: Level: intermediate
3757: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3758: @*/
3759: PetscErrorCode SNESSetForceIteration(SNES snes,PetscBool force)
3760: {
3763: snes->forceiteration = force;
3764: return(0);
3765: }
3767: /*@
3768: SNESGetForceIteration - Whether or not to force SNESSolve() take at least one iteration regardless of the initial residual norm
3770: Logically Collective on SNES
3772: Input Parameters:
3773: . snes - the SNES context
3775: Output Parameter:
3776: . force - PETSC_TRUE requires at least one iteration.
3778: Level: intermediate
3780: .seealso: SNESSetForceIteration(), SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3781: @*/
3782: PetscErrorCode SNESGetForceIteration(SNES snes,PetscBool *force)
3783: {
3786: *force = snes->forceiteration;
3787: return(0);
3788: }
3790: /*@
3791: SNESSetTolerances - Sets various parameters used in convergence tests.
3793: Logically Collective on SNES
3795: Input Parameters:
3796: + snes - the SNES context
3797: . abstol - absolute convergence tolerance
3798: . rtol - relative convergence tolerance
3799: . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x ||
3800: . maxit - maximum number of iterations
3801: - maxf - maximum number of function evaluations (-1 indicates no limit)
3803: Options Database Keys:
3804: + -snes_atol <abstol> - Sets abstol
3805: . -snes_rtol <rtol> - Sets rtol
3806: . -snes_stol <stol> - Sets stol
3807: . -snes_max_it <maxit> - Sets maxit
3808: - -snes_max_funcs <maxf> - Sets maxf
3810: Notes:
3811: The default maximum number of iterations is 50.
3812: The default maximum number of function evaluations is 1000.
3814: Level: intermediate
3816: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance(), SNESSetForceIteration()
3817: @*/
3818: PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3819: {
3828: if (abstol != PETSC_DEFAULT) {
3829: if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3830: snes->abstol = abstol;
3831: }
3832: if (rtol != PETSC_DEFAULT) {
3833: 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);
3834: snes->rtol = rtol;
3835: }
3836: if (stol != PETSC_DEFAULT) {
3837: if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3838: snes->stol = stol;
3839: }
3840: if (maxit != PETSC_DEFAULT) {
3841: if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3842: snes->max_its = maxit;
3843: }
3844: if (maxf != PETSC_DEFAULT) {
3845: if (maxf < -1) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be -1 or nonnegative",maxf);
3846: snes->max_funcs = maxf;
3847: }
3848: snes->tolerancesset = PETSC_TRUE;
3849: return(0);
3850: }
3852: /*@
3853: SNESSetDivergenceTolerance - Sets the divergence tolerance used for the SNES divergence test.
3855: Logically Collective on SNES
3857: Input Parameters:
3858: + snes - the SNES context
3859: - divtol - the divergence tolerance. Use -1 to deactivate the test.
3861: Options Database Keys:
3862: . -snes_divergence_tolerance <divtol> - Sets divtol
3864: Notes:
3865: The default divergence tolerance is 1e4.
3867: Level: intermediate
3869: .seealso: SNESSetTolerances(), SNESGetDivergenceTolerance
3870: @*/
3871: PetscErrorCode SNESSetDivergenceTolerance(SNES snes,PetscReal divtol)
3872: {
3877: if (divtol != PETSC_DEFAULT) {
3878: snes->divtol = divtol;
3879: }
3880: else {
3881: snes->divtol = 1.0e4;
3882: }
3883: return(0);
3884: }
3886: /*@
3887: SNESGetTolerances - Gets various parameters used in convergence tests.
3889: Not Collective
3891: Input Parameters:
3892: + snes - the SNES context
3893: . atol - absolute convergence tolerance
3894: . rtol - relative convergence tolerance
3895: . stol - convergence tolerance in terms of the norm
3896: of the change in the solution between steps
3897: . maxit - maximum number of iterations
3898: - maxf - maximum number of function evaluations
3900: Notes:
3901: The user can specify NULL for any parameter that is not needed.
3903: Level: intermediate
3905: .seealso: SNESSetTolerances()
3906: @*/
3907: PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3908: {
3911: if (atol) *atol = snes->abstol;
3912: if (rtol) *rtol = snes->rtol;
3913: if (stol) *stol = snes->stol;
3914: if (maxit) *maxit = snes->max_its;
3915: if (maxf) *maxf = snes->max_funcs;
3916: return(0);
3917: }
3919: /*@
3920: SNESGetDivergenceTolerance - Gets divergence tolerance used in divergence test.
3922: Not Collective
3924: Input Parameters:
3925: + snes - the SNES context
3926: - divtol - divergence tolerance
3928: Level: intermediate
3930: .seealso: SNESSetDivergenceTolerance()
3931: @*/
3932: PetscErrorCode SNESGetDivergenceTolerance(SNES snes,PetscReal *divtol)
3933: {
3936: if (divtol) *divtol = snes->divtol;
3937: return(0);
3938: }
3940: /*@
3941: SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
3943: Logically Collective on SNES
3945: Input Parameters:
3946: + snes - the SNES context
3947: - tol - tolerance
3949: Options Database Key:
3950: . -snes_trtol <tol> - Sets tol
3952: Level: intermediate
3954: .seealso: SNESSetTolerances()
3955: @*/
3956: PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3957: {
3961: snes->deltatol = tol;
3962: return(0);
3963: }
3965: PETSC_INTERN PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3967: PetscErrorCode SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3968: {
3969: PetscDrawLG lg;
3970: PetscErrorCode ierr;
3971: PetscReal x,y,per;
3972: PetscViewer v = (PetscViewer)monctx;
3973: static PetscReal prev; /* should be in the context */
3974: PetscDraw draw;
3978: PetscViewerDrawGetDrawLG(v,0,&lg);
3979: if (!n) {PetscDrawLGReset(lg);}
3980: PetscDrawLGGetDraw(lg,&draw);
3981: PetscDrawSetTitle(draw,"Residual norm");
3982: x = (PetscReal)n;
3983: if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3984: else y = -15.0;
3985: PetscDrawLGAddPoint(lg,&x,&y);
3986: if (n < 20 || !(n % 5) || snes->reason) {
3987: PetscDrawLGDraw(lg);
3988: PetscDrawLGSave(lg);
3989: }
3991: PetscViewerDrawGetDrawLG(v,1,&lg);
3992: if (!n) {PetscDrawLGReset(lg);}
3993: PetscDrawLGGetDraw(lg,&draw);
3994: PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3995: SNESMonitorRange_Private(snes,n,&per);
3996: x = (PetscReal)n;
3997: y = 100.0*per;
3998: PetscDrawLGAddPoint(lg,&x,&y);
3999: if (n < 20 || !(n % 5) || snes->reason) {
4000: PetscDrawLGDraw(lg);
4001: PetscDrawLGSave(lg);
4002: }
4004: PetscViewerDrawGetDrawLG(v,2,&lg);
4005: if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
4006: PetscDrawLGGetDraw(lg,&draw);
4007: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
4008: x = (PetscReal)n;
4009: y = (prev - rnorm)/prev;
4010: PetscDrawLGAddPoint(lg,&x,&y);
4011: if (n < 20 || !(n % 5) || snes->reason) {
4012: PetscDrawLGDraw(lg);
4013: PetscDrawLGSave(lg);
4014: }
4016: PetscViewerDrawGetDrawLG(v,3,&lg);
4017: if (!n) {PetscDrawLGReset(lg);}
4018: PetscDrawLGGetDraw(lg,&draw);
4019: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
4020: x = (PetscReal)n;
4021: y = (prev - rnorm)/(prev*per);
4022: if (n > 2) { /*skip initial crazy value */
4023: PetscDrawLGAddPoint(lg,&x,&y);
4024: }
4025: if (n < 20 || !(n % 5) || snes->reason) {
4026: PetscDrawLGDraw(lg);
4027: PetscDrawLGSave(lg);
4028: }
4029: prev = rnorm;
4030: return(0);
4031: }
4033: /*@
4034: SNESMonitor - runs the user provided monitor routines, if they exist
4036: Collective on SNES
4038: Input Parameters:
4039: + snes - nonlinear solver context obtained from SNESCreate()
4040: . iter - iteration number
4041: - rnorm - relative norm of the residual
4043: Notes:
4044: This routine is called by the SNES implementations.
4045: It does not typically need to be called by the user.
4047: Level: developer
4049: .seealso: SNESMonitorSet()
4050: @*/
4051: PetscErrorCode SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
4052: {
4054: PetscInt i,n = snes->numbermonitors;
4057: VecLockReadPush(snes->vec_sol);
4058: for (i=0; i<n; i++) {
4059: (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
4060: }
4061: VecLockReadPop(snes->vec_sol);
4062: return(0);
4063: }
4065: /* ------------ Routines to set performance monitoring options ----------- */
4067: /*MC
4068: SNESMonitorFunction - functional form passed to SNESMonitorSet() to monitor convergence of nonlinear solver
4070: Synopsis:
4071: #include <petscsnes.h>
4072: $ PetscErrorCode SNESMonitorFunction(SNES snes,PetscInt its, PetscReal norm,void *mctx)
4074: Collective on snes
4076: Input Parameters:
4077: + snes - the SNES context
4078: . its - iteration number
4079: . norm - 2-norm function value (may be estimated)
4080: - mctx - [optional] monitoring context
4082: Level: advanced
4084: .seealso: SNESMonitorSet(), SNESMonitorGet()
4085: M*/
4087: /*@C
4088: SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
4089: iteration of the nonlinear solver to display the iteration's
4090: progress.
4092: Logically Collective on SNES
4094: Input Parameters:
4095: + snes - the SNES context
4096: . f - the monitor function, see SNESMonitorFunction for the calling sequence
4097: . mctx - [optional] user-defined context for private data for the
4098: monitor routine (use NULL if no context is desired)
4099: - monitordestroy - [optional] routine that frees monitor context
4100: (may be NULL)
4102: Options Database Keys:
4103: + -snes_monitor - sets SNESMonitorDefault()
4104: . -snes_monitor draw::draw_lg - sets line graph monitor,
4105: - -snes_monitor_cancel - cancels all monitors that have
4106: been hardwired into a code by
4107: calls to SNESMonitorSet(), but
4108: does not cancel those set via
4109: the options database.
4111: Notes:
4112: Several different monitoring routines may be set by calling
4113: SNESMonitorSet() multiple times; all will be called in the
4114: order in which they were set.
4116: Fortran Notes:
4117: Only a single monitor function can be set for each SNES object
4119: Level: intermediate
4121: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
4122: @*/
4123: PetscErrorCode SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
4124: {
4125: PetscInt i;
4127: PetscBool identical;
4131: for (i=0; i<snes->numbermonitors;i++) {
4132: PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))snes->monitor[i],snes->monitorcontext[i],snes->monitordestroy[i],&identical);
4133: if (identical) return(0);
4134: }
4135: if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
4136: snes->monitor[snes->numbermonitors] = f;
4137: snes->monitordestroy[snes->numbermonitors] = monitordestroy;
4138: snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
4139: return(0);
4140: }
4142: /*@
4143: SNESMonitorCancel - Clears all the monitor functions for a SNES object.
4145: Logically Collective on SNES
4147: Input Parameters:
4148: . snes - the SNES context
4150: Options Database Key:
4151: . -snes_monitor_cancel - cancels all monitors that have been hardwired
4152: into a code by calls to SNESMonitorSet(), but does not cancel those
4153: set via the options database
4155: Notes:
4156: There is no way to clear one specific monitor from a SNES object.
4158: Level: intermediate
4160: .seealso: SNESMonitorDefault(), SNESMonitorSet()
4161: @*/
4162: PetscErrorCode SNESMonitorCancel(SNES snes)
4163: {
4165: PetscInt i;
4169: for (i=0; i<snes->numbermonitors; i++) {
4170: if (snes->monitordestroy[i]) {
4171: (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
4172: }
4173: }
4174: snes->numbermonitors = 0;
4175: return(0);
4176: }
4178: /*MC
4179: SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver
4181: Synopsis:
4182: #include <petscsnes.h>
4183: $ PetscErrorCode SNESConvergenceTest(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
4185: Collective on snes
4187: Input Parameters:
4188: + snes - the SNES context
4189: . it - current iteration (0 is the first and is before any Newton step)
4190: . xnorm - 2-norm of current iterate
4191: . gnorm - 2-norm of current step
4192: . f - 2-norm of function
4193: - cctx - [optional] convergence context
4195: Output Parameter:
4196: . reason - reason for convergence/divergence, only needs to be set when convergence or divergence is detected
4198: Level: intermediate
4200: .seealso: SNESSetConvergenceTest(), SNESGetConvergenceTest()
4201: M*/
4203: /*@C
4204: SNESSetConvergenceTest - Sets the function that is to be used
4205: to test for convergence of the nonlinear iterative solution.
4207: Logically Collective on SNES
4209: Input Parameters:
4210: + snes - the SNES context
4211: . SNESConvergenceTestFunction - routine to test for convergence
4212: . cctx - [optional] context for private data for the convergence routine (may be NULL)
4213: - destroy - [optional] destructor for the context (may be NULL; PETSC_NULL_FUNCTION in Fortran)
4215: Level: advanced
4217: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
4218: @*/
4219: PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
4220: {
4225: if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
4226: if (snes->ops->convergeddestroy) {
4227: (*snes->ops->convergeddestroy)(snes->cnvP);
4228: }
4229: snes->ops->converged = SNESConvergenceTestFunction;
4230: snes->ops->convergeddestroy = destroy;
4231: snes->cnvP = cctx;
4232: return(0);
4233: }
4235: /*@
4236: SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
4238: Not Collective
4240: Input Parameter:
4241: . snes - the SNES context
4243: Output Parameter:
4244: . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
4245: manual pages for the individual convergence tests for complete lists
4247: Options Database:
4248: . -snes_converged_reason - prints the reason to standard out
4250: Level: intermediate
4252: Notes:
4253: Should only be called after the call the SNESSolve() is complete, if it is called earlier it returns the value SNES__CONVERGED_ITERATING.
4255: .seealso: SNESSetConvergenceTest(), SNESSetConvergedReason(), SNESConvergedReason
4256: @*/
4257: PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
4258: {
4262: *reason = snes->reason;
4263: return(0);
4264: }
4266: /*@C
4267: SNESGetConvergedReasonString - Return a human readable string for snes converged reason
4269: Not Collective
4271: Input Parameter:
4272: . snes - the SNES context
4274: Output Parameter:
4275: . strreason - a human readable string that describes SNES converged reason
4277: Level: beginner
4279: .seealso: SNESGetConvergedReason()
4280: @*/
4281: PetscErrorCode SNESGetConvergedReasonString(SNES snes, const char** strreason)
4282: {
4286: *strreason = SNESConvergedReasons[snes->reason];
4287: return(0);
4288: }
4290: /*@
4291: SNESSetConvergedReason - Sets the reason the SNES iteration was stopped.
4293: Not Collective
4295: Input Parameters:
4296: + snes - the SNES context
4297: - reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
4298: manual pages for the individual convergence tests for complete lists
4300: Level: intermediate
4302: .seealso: SNESGetConvergedReason(), SNESSetConvergenceTest(), SNESConvergedReason
4303: @*/
4304: PetscErrorCode SNESSetConvergedReason(SNES snes,SNESConvergedReason reason)
4305: {
4308: snes->reason = reason;
4309: return(0);
4310: }
4312: /*@
4313: SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
4315: Logically Collective on SNES
4317: Input Parameters:
4318: + snes - iterative context obtained from SNESCreate()
4319: . a - array to hold history, this array will contain the function norms computed at each step
4320: . its - integer array holds the number of linear iterations for each solve.
4321: . na - size of a and its
4322: - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
4323: else it continues storing new values for new nonlinear solves after the old ones
4325: Notes:
4326: If 'a' and 'its' are NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
4327: default array of length 10000 is allocated.
4329: This routine is useful, e.g., when running a code for purposes
4330: of accurate performance monitoring, when no I/O should be done
4331: during the section of code that is being timed.
4333: Level: intermediate
4335: .seealso: SNESGetConvergenceHistory()
4337: @*/
4338: PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
4339: {
4346: if (!a) {
4347: if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
4348: PetscCalloc2(na,&a,na,&its);
4349: snes->conv_hist_alloc = PETSC_TRUE;
4350: }
4351: snes->conv_hist = a;
4352: snes->conv_hist_its = its;
4353: snes->conv_hist_max = na;
4354: snes->conv_hist_len = 0;
4355: snes->conv_hist_reset = reset;
4356: return(0);
4357: }
4359: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4360: #include <engine.h> /* MATLAB include file */
4361: #include <mex.h> /* MATLAB include file */
4363: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
4364: {
4365: mxArray *mat;
4366: PetscInt i;
4367: PetscReal *ar;
4370: mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
4371: ar = (PetscReal*) mxGetData(mat);
4372: for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
4373: PetscFunctionReturn(mat);
4374: }
4375: #endif
4377: /*@C
4378: SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
4380: Not Collective
4382: Input Parameter:
4383: . snes - iterative context obtained from SNESCreate()
4385: Output Parameters:
4386: + a - array to hold history
4387: . its - integer array holds the number of linear iterations (or
4388: negative if not converged) for each solve.
4389: - na - size of a and its
4391: Notes:
4392: The calling sequence for this routine in Fortran is
4393: $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
4395: This routine is useful, e.g., when running a code for purposes
4396: of accurate performance monitoring, when no I/O should be done
4397: during the section of code that is being timed.
4399: Level: intermediate
4401: .seealso: SNESSetConvergenceHistory()
4403: @*/
4404: PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
4405: {
4408: if (a) *a = snes->conv_hist;
4409: if (its) *its = snes->conv_hist_its;
4410: if (na) *na = snes->conv_hist_len;
4411: return(0);
4412: }
4414: /*@C
4415: SNESSetUpdate - Sets the general-purpose update function called
4416: at the beginning of every iteration of the nonlinear solve. Specifically
4417: it is called just before the Jacobian is "evaluated".
4419: Logically Collective on SNES
4421: Input Parameters:
4422: + snes - The nonlinear solver context
4423: - func - The function
4425: Calling sequence of func:
4426: $ func (SNES snes, PetscInt step);
4428: . step - The current step of the iteration
4430: Level: advanced
4432: Note:
4433: 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()
4434: This is not used by most users.
4436: There are a varity of function hooks one many set that are called at different stages of the nonlinear solution process, see the functions listed below.
4438: .seealso SNESSetJacobian(), SNESSolve(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck(), SNESNewtonTRSetPreCheck(), SNESNewtonTRSetPostCheck(),
4439: SNESMonitorSet(), SNESSetDivergenceTest()
4440: @*/
4441: PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
4442: {
4445: snes->ops->update = func;
4446: return(0);
4447: }
4449: /*
4450: SNESScaleStep_Private - Scales a step so that its length is less than the
4451: positive parameter delta.
4453: Input Parameters:
4454: + snes - the SNES context
4455: . y - approximate solution of linear system
4456: . fnorm - 2-norm of current function
4457: - delta - trust region size
4459: Output Parameters:
4460: + gpnorm - predicted function norm at the new point, assuming local
4461: linearization. The value is zero if the step lies within the trust
4462: region, and exceeds zero otherwise.
4463: - ynorm - 2-norm of the step
4465: Note:
4466: For non-trust region methods such as SNESNEWTONLS, the parameter delta
4467: is set to be the maximum allowable step size.
4469: */
4470: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
4471: {
4472: PetscReal nrm;
4473: PetscScalar cnorm;
4481: VecNorm(y,NORM_2,&nrm);
4482: if (nrm > *delta) {
4483: nrm = *delta/nrm;
4484: *gpnorm = (1.0 - nrm)*(*fnorm);
4485: cnorm = nrm;
4486: VecScale(y,cnorm);
4487: *ynorm = *delta;
4488: } else {
4489: *gpnorm = 0.0;
4490: *ynorm = nrm;
4491: }
4492: return(0);
4493: }
4495: /*@C
4496: SNESConvergedReasonView - Displays the reason a SNES solve converged or diverged to a viewer
4498: Collective on SNES
4500: Parameter:
4501: + snes - iterative context obtained from SNESCreate()
4502: - viewer - the viewer to display the reason
4504: Options Database Keys:
4505: + -snes_converged_reason - print reason for converged or diverged, also prints number of iterations
4506: - -snes_converged_reason ::failed - only print reason and number of iterations when diverged
4508: Notes:
4509: To change the format of the output call PetscViewerPushFormat(viewer,format) before this call. Use PETSC_VIEWER_DEFAULT for the default,
4510: use PETSC_VIEWER_FAILED to only display a reason if it fails.
4512: Level: beginner
4514: .seealso: SNESCreate(), SNESSetUp(), SNESDestroy(), SNESSetTolerances(), SNESConvergedDefault(), SNESGetConvergedReason(), SNESConvergedReasonViewFromOptions(),
4515: PetscViewerPushFormat(), PetscViewerPopFormat()
4517: @*/
4518: PetscErrorCode SNESConvergedReasonView(SNES snes,PetscViewer viewer)
4519: {
4520: PetscViewerFormat format;
4521: PetscBool isAscii;
4522: PetscErrorCode ierr;
4525: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
4526: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
4527: if (isAscii) {
4528: PetscViewerGetFormat(viewer, &format);
4529: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
4530: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
4531: DM dm;
4532: Vec u;
4533: PetscDS prob;
4534: PetscInt Nf, f;
4535: PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
4536: void **exactCtx;
4537: PetscReal error;
4539: SNESGetDM(snes, &dm);
4540: SNESGetSolution(snes, &u);
4541: DMGetDS(dm, &prob);
4542: PetscDSGetNumFields(prob, &Nf);
4543: PetscMalloc2(Nf, &exactSol, Nf, &exactCtx);
4544: for (f = 0; f < Nf; ++f) {PetscDSGetExactSolution(prob, f, &exactSol[f], &exactCtx[f]);}
4545: DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error);
4546: PetscFree2(exactSol, exactCtx);
4547: if (error < 1.0e-11) {PetscViewerASCIIPrintf(viewer, "L_2 Error: < 1.0e-11\n");}
4548: else {PetscViewerASCIIPrintf(viewer, "L_2 Error: %g\n", error);}
4549: }
4550: if (snes->reason > 0 && format != PETSC_VIEWER_FAILED) {
4551: if (((PetscObject) snes)->prefix) {
4552: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4553: } else {
4554: PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4555: }
4556: } else if (snes->reason <= 0) {
4557: if (((PetscObject) snes)->prefix) {
4558: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4559: } else {
4560: PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4561: }
4562: }
4563: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
4564: }
4565: return(0);
4566: }
4568: /*@C
4569: SNESConvergedReasonViewSet - Sets an ADDITIONAL function that is to be used at the
4570: end of the nonlinear solver to display the conver reason of the nonlinear solver.
4572: Logically Collective on SNES
4574: Input Parameters:
4575: + snes - the SNES context
4576: . f - the snes converged reason view function
4577: . vctx - [optional] user-defined context for private data for the
4578: snes converged reason view routine (use NULL if no context is desired)
4579: - reasonviewdestroy - [optional] routine that frees reasonview context
4580: (may be NULL)
4582: Options Database Keys:
4583: + -snes_converged_reason - sets a default SNESConvergedReasonView()
4584: - -snes_converged_reason_view_cancel - cancels all converged reason viewers that have
4585: been hardwired into a code by
4586: calls to SNESConvergedReasonViewSet(), but
4587: does not cancel those set via
4588: the options database.
4590: Notes:
4591: Several different converged reason view routines may be set by calling
4592: SNESConvergedReasonViewSet() multiple times; all will be called in the
4593: order in which they were set.
4595: Level: intermediate
4597: .seealso: SNESConvergedReasonView(), SNESConvergedReasonViewCancel()
4598: @*/
4599: PetscErrorCode SNESConvergedReasonViewSet(SNES snes,PetscErrorCode (*f)(SNES,void*),void *vctx,PetscErrorCode (*reasonviewdestroy)(void**))
4600: {
4601: PetscInt i;
4603: PetscBool identical;
4607: for (i=0; i<snes->numberreasonviews;i++) {
4608: PetscMonitorCompare((PetscErrorCode (*)(void))f,vctx,reasonviewdestroy,(PetscErrorCode (*)(void))snes->reasonview[i],snes->reasonviewcontext[i],snes->reasonviewdestroy[i],&identical);
4609: if (identical) return(0);
4610: }
4611: if (snes->numberreasonviews >= MAXSNESREASONVIEWS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many SNES reasonview set");
4612: snes->reasonview[snes->numberreasonviews] = f;
4613: snes->reasonviewdestroy[snes->numberreasonviews] = reasonviewdestroy;
4614: snes->reasonviewcontext[snes->numberreasonviews++] = (void*)vctx;
4615: return(0);
4616: }
4618: /*@
4619: SNESConvergedReasonViewFromOptions - Processes command line options to determine if/how a SNESReason is to be viewed.
4620: All the user-provided convergedReasonView routines will be involved as well, if they exist.
4622: Collective on SNES
4624: Input Parameters:
4625: . snes - the SNES object
4627: Level: intermediate
4629: .seealso: SNESCreate(), SNESSetUp(), SNESDestroy(), SNESSetTolerances(), SNESConvergedDefault(), SNESGetConvergedReason(), SNESConvergedReasonView()
4631: @*/
4632: PetscErrorCode SNESConvergedReasonViewFromOptions(SNES snes)
4633: {
4634: PetscErrorCode ierr;
4635: PetscViewer viewer;
4636: PetscBool flg;
4637: static PetscBool incall = PETSC_FALSE;
4638: PetscViewerFormat format;
4639: PetscInt i;
4642: if (incall) return(0);
4643: incall = PETSC_TRUE;
4645: /* All user-provided viewers are called first, if they exist. */
4646: for (i=0; i<snes->numberreasonviews; i++) {
4647: (*snes->reasonview[i])(snes,snes->reasonviewcontext[i]);
4648: }
4650: /* Call PETSc default routine if users ask for it */
4651: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
4652: if (flg) {
4653: PetscViewerPushFormat(viewer,format);
4654: SNESConvergedReasonView(snes,viewer);
4655: PetscViewerPopFormat(viewer);
4656: PetscViewerDestroy(&viewer);
4657: }
4658: incall = PETSC_FALSE;
4659: return(0);
4660: }
4662: /*@
4663: SNESSolve - Solves a nonlinear system F(x) = b.
4664: Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
4666: Collective on SNES
4668: Input Parameters:
4669: + snes - the SNES context
4670: . b - the constant part of the equation F(x) = b, or NULL to use zero.
4671: - x - the solution vector.
4673: Notes:
4674: The user should initialize the vector,x, with the initial guess
4675: for the nonlinear solve prior to calling SNESSolve(). In particular,
4676: to employ an initial guess of zero, the user should explicitly set
4677: this vector to zero by calling VecSet().
4679: Level: beginner
4681: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution(),
4682: SNESNewtonTRSetPreCheck(), SNESNewtonTRGetPreCheck(), SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck(),
4683: SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck()
4684: @*/
4685: PetscErrorCode SNESSolve(SNES snes,Vec b,Vec x)
4686: {
4687: PetscErrorCode ierr;
4688: PetscBool flg;
4689: PetscInt grid;
4690: Vec xcreated = NULL;
4691: DM dm;
4700: /* High level operations using the nonlinear solver */
4701: {
4702: PetscViewer viewer;
4703: PetscViewerFormat format;
4704: PetscInt num;
4705: PetscBool flg;
4706: static PetscBool incall = PETSC_FALSE;
4708: if (!incall) {
4709: /* Estimate the convergence rate of the discretization */
4710: PetscOptionsGetViewer(PetscObjectComm((PetscObject) snes),((PetscObject)snes)->options, ((PetscObject) snes)->prefix, "-snes_convergence_estimate", &viewer, &format, &flg);
4711: if (flg) {
4712: PetscConvEst conv;
4713: DM dm;
4714: PetscReal *alpha; /* Convergence rate of the solution error for each field in the L_2 norm */
4715: PetscInt Nf;
4717: incall = PETSC_TRUE;
4718: SNESGetDM(snes, &dm);
4719: DMGetNumFields(dm, &Nf);
4720: PetscCalloc1(Nf, &alpha);
4721: PetscConvEstCreate(PetscObjectComm((PetscObject) snes), &conv);
4722: PetscConvEstSetSolver(conv, (PetscObject) snes);
4723: PetscConvEstSetFromOptions(conv);
4724: PetscConvEstSetUp(conv);
4725: PetscConvEstGetConvRate(conv, alpha);
4726: PetscViewerPushFormat(viewer, format);
4727: PetscConvEstRateView(conv, alpha, viewer);
4728: PetscViewerPopFormat(viewer);
4729: PetscViewerDestroy(&viewer);
4730: PetscConvEstDestroy(&conv);
4731: PetscFree(alpha);
4732: incall = PETSC_FALSE;
4733: }
4734: /* Adaptively refine the initial grid */
4735: num = 1;
4736: PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_initial", &num, &flg);
4737: if (flg) {
4738: DMAdaptor adaptor;
4740: incall = PETSC_TRUE;
4741: DMAdaptorCreate(PetscObjectComm((PetscObject)snes), &adaptor);
4742: DMAdaptorSetSolver(adaptor, snes);
4743: DMAdaptorSetSequenceLength(adaptor, num);
4744: DMAdaptorSetFromOptions(adaptor);
4745: DMAdaptorSetUp(adaptor);
4746: DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_INITIAL, &dm, &x);
4747: DMAdaptorDestroy(&adaptor);
4748: incall = PETSC_FALSE;
4749: }
4750: /* Use grid sequencing to adapt */
4751: num = 0;
4752: PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_sequence", &num, NULL);
4753: if (num) {
4754: DMAdaptor adaptor;
4756: incall = PETSC_TRUE;
4757: DMAdaptorCreate(PetscObjectComm((PetscObject)snes), &adaptor);
4758: DMAdaptorSetSolver(adaptor, snes);
4759: DMAdaptorSetSequenceLength(adaptor, num);
4760: DMAdaptorSetFromOptions(adaptor);
4761: DMAdaptorSetUp(adaptor);
4762: DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_SEQUENTIAL, &dm, &x);
4763: DMAdaptorDestroy(&adaptor);
4764: incall = PETSC_FALSE;
4765: }
4766: }
4767: }
4768: if (!x) {
4769: SNESGetDM(snes,&dm);
4770: DMCreateGlobalVector(dm,&xcreated);
4771: x = xcreated;
4772: }
4773: SNESViewFromOptions(snes,NULL,"-snes_view_pre");
4775: for (grid=0; grid<snes->gridsequence; grid++) {PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));}
4776: for (grid=0; grid<snes->gridsequence+1; grid++) {
4778: /* set solution vector */
4779: if (!grid) {PetscObjectReference((PetscObject)x);}
4780: VecDestroy(&snes->vec_sol);
4781: snes->vec_sol = x;
4782: SNESGetDM(snes,&dm);
4784: /* set affine vector if provided */
4785: if (b) { PetscObjectReference((PetscObject)b); }
4786: VecDestroy(&snes->vec_rhs);
4787: snes->vec_rhs = b;
4789: if (snes->vec_rhs && (snes->vec_func == snes->vec_rhs)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Right hand side vector cannot be function vector");
4790: if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
4791: if (snes->vec_rhs == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
4792: if (!snes->vec_sol_update /* && snes->vec_sol */) {
4793: VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
4794: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
4795: }
4796: DMShellSetGlobalVector(dm,snes->vec_sol);
4797: SNESSetUp(snes);
4799: if (!grid) {
4800: if (snes->ops->computeinitialguess) {
4801: (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
4802: }
4803: }
4805: if (snes->conv_hist_reset) snes->conv_hist_len = 0;
4806: if (snes->counters_reset) {snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;}
4808: PetscLogEventBegin(SNES_Solve,snes,0,0,0);
4809: (*snes->ops->solve)(snes);
4810: PetscLogEventEnd(SNES_Solve,snes,0,0,0);
4811: if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
4812: snes->domainerror = PETSC_FALSE; /* clear the flag if it has been set */
4814: if (snes->lagjac_persist) snes->jac_iter += snes->iter;
4815: if (snes->lagpre_persist) snes->pre_iter += snes->iter;
4817: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_local_min",NULL,NULL,&flg);
4818: if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
4819: /* Call converged reason views. This may involve user-provided viewers as well */
4820: SNESConvergedReasonViewFromOptions(snes);
4822: if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
4823: if (snes->reason < 0) break;
4824: if (grid < snes->gridsequence) {
4825: DM fine;
4826: Vec xnew;
4827: Mat interp;
4829: DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
4830: if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
4831: DMCreateInterpolation(snes->dm,fine,&interp,NULL);
4832: DMCreateGlobalVector(fine,&xnew);
4833: MatInterpolate(interp,x,xnew);
4834: DMInterpolate(snes->dm,interp,fine);
4835: MatDestroy(&interp);
4836: x = xnew;
4838: SNESReset(snes);
4839: SNESSetDM(snes,fine);
4840: SNESResetFromOptions(snes);
4841: DMDestroy(&fine);
4842: PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
4843: }
4844: }
4845: SNESViewFromOptions(snes,NULL,"-snes_view");
4846: VecViewFromOptions(snes->vec_sol,(PetscObject)snes,"-snes_view_solution");
4847: DMMonitor(snes->dm);
4848: SNESMonitorPauseFinal_Internal(snes);
4850: VecDestroy(&xcreated);
4851: PetscObjectSAWsBlock((PetscObject)snes);
4852: return(0);
4853: }
4855: /* --------- Internal routines for SNES Package --------- */
4857: /*@C
4858: SNESSetType - Sets the method for the nonlinear solver.
4860: Collective on SNES
4862: Input Parameters:
4863: + snes - the SNES context
4864: - type - a known method
4866: Options Database Key:
4867: . -snes_type <type> - Sets the method; use -help for a list
4868: of available methods (for instance, newtonls or newtontr)
4870: Notes:
4871: See "petsc/include/petscsnes.h" for available methods (for instance)
4872: + SNESNEWTONLS - Newton's method with line search
4873: (systems of nonlinear equations)
4874: - SNESNEWTONTR - Newton's method with trust region
4875: (systems of nonlinear equations)
4877: Normally, it is best to use the SNESSetFromOptions() command and then
4878: set the SNES solver type from the options database rather than by using
4879: this routine. Using the options database provides the user with
4880: maximum flexibility in evaluating the many nonlinear solvers.
4881: The SNESSetType() routine is provided for those situations where it
4882: is necessary to set the nonlinear solver independently of the command
4883: line or options database. This might be the case, for example, when
4884: the choice of solver changes during the execution of the program,
4885: and the user's application is taking responsibility for choosing the
4886: appropriate method.
4888: Developer Notes:
4889: SNESRegister() adds a constructor for a new SNESType to SNESList, SNESSetType() locates
4890: the constructor in that list and calls it to create the spexific object.
4892: Level: intermediate
4894: .seealso: SNESType, SNESCreate(), SNESDestroy(), SNESGetType(), SNESSetFromOptions()
4896: @*/
4897: PetscErrorCode SNESSetType(SNES snes,SNESType type)
4898: {
4899: PetscErrorCode ierr,(*r)(SNES);
4900: PetscBool match;
4906: PetscObjectTypeCompare((PetscObject)snes,type,&match);
4907: if (match) return(0);
4909: PetscFunctionListFind(SNESList,type,&r);
4910: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
4911: /* Destroy the previous private SNES context */
4912: if (snes->ops->destroy) {
4913: (*(snes)->ops->destroy)(snes);
4914: snes->ops->destroy = NULL;
4915: }
4916: /* Reinitialize function pointers in SNESOps structure */
4917: snes->ops->setup = NULL;
4918: snes->ops->solve = NULL;
4919: snes->ops->view = NULL;
4920: snes->ops->setfromoptions = NULL;
4921: snes->ops->destroy = NULL;
4923: /* It may happen the user has customized the line search before calling SNESSetType */
4924: if (((PetscObject)snes)->type_name) {
4925: SNESLineSearchDestroy(&snes->linesearch);
4926: }
4928: /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
4929: snes->setupcalled = PETSC_FALSE;
4931: PetscObjectChangeTypeName((PetscObject)snes,type);
4932: (*r)(snes);
4933: return(0);
4934: }
4936: /*@C
4937: SNESGetType - Gets the SNES method type and name (as a string).
4939: Not Collective
4941: Input Parameter:
4942: . snes - nonlinear solver context
4944: Output Parameter:
4945: . type - SNES method (a character string)
4947: Level: intermediate
4949: @*/
4950: PetscErrorCode SNESGetType(SNES snes,SNESType *type)
4951: {
4955: *type = ((PetscObject)snes)->type_name;
4956: return(0);
4957: }
4959: /*@
4960: SNESSetSolution - Sets the solution vector for use by the SNES routines.
4962: Logically Collective on SNES
4964: Input Parameters:
4965: + snes - the SNES context obtained from SNESCreate()
4966: - u - the solution vector
4968: Level: beginner
4970: @*/
4971: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4972: {
4973: DM dm;
4979: PetscObjectReference((PetscObject) u);
4980: VecDestroy(&snes->vec_sol);
4982: snes->vec_sol = u;
4984: SNESGetDM(snes, &dm);
4985: DMShellSetGlobalVector(dm, u);
4986: return(0);
4987: }
4989: /*@
4990: SNESGetSolution - Returns the vector where the approximate solution is
4991: stored. This is the fine grid solution when using SNESSetGridSequence().
4993: Not Collective, but Vec is parallel if SNES is parallel
4995: Input Parameter:
4996: . snes - the SNES context
4998: Output Parameter:
4999: . x - the solution
5001: Level: intermediate
5003: .seealso: SNESGetSolutionUpdate(), SNESGetFunction()
5004: @*/
5005: PetscErrorCode SNESGetSolution(SNES snes,Vec *x)
5006: {
5010: *x = snes->vec_sol;
5011: return(0);
5012: }
5014: /*@
5015: SNESGetSolutionUpdate - Returns the vector where the solution update is
5016: stored.
5018: Not Collective, but Vec is parallel if SNES is parallel
5020: Input Parameter:
5021: . snes - the SNES context
5023: Output Parameter:
5024: . x - the solution update
5026: Level: advanced
5028: .seealso: SNESGetSolution(), SNESGetFunction()
5029: @*/
5030: PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x)
5031: {
5035: *x = snes->vec_sol_update;
5036: return(0);
5037: }
5039: /*@C
5040: SNESGetFunction - Returns the vector where the function is stored.
5042: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
5044: Input Parameter:
5045: . snes - the SNES context
5047: Output Parameters:
5048: + r - the vector that is used to store residuals (or NULL if you don't want it)
5049: . f - the function (or NULL if you don't want it); see SNESFunction for calling sequence details
5050: - ctx - the function context (or NULL if you don't want it)
5052: Level: advanced
5054: Notes: The vector r DOES NOT, in general contain the current value of the SNES nonlinear function
5056: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
5057: @*/
5058: PetscErrorCode SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
5059: {
5061: DM dm;
5065: if (r) {
5066: if (!snes->vec_func) {
5067: if (snes->vec_rhs) {
5068: VecDuplicate(snes->vec_rhs,&snes->vec_func);
5069: } else if (snes->vec_sol) {
5070: VecDuplicate(snes->vec_sol,&snes->vec_func);
5071: } else if (snes->dm) {
5072: DMCreateGlobalVector(snes->dm,&snes->vec_func);
5073: }
5074: }
5075: *r = snes->vec_func;
5076: }
5077: SNESGetDM(snes,&dm);
5078: DMSNESGetFunction(dm,f,ctx);
5079: return(0);
5080: }
5082: /*@C
5083: SNESGetNGS - Returns the NGS function and context.
5085: Input Parameter:
5086: . snes - the SNES context
5088: Output Parameters:
5089: + f - the function (or NULL) see SNESNGSFunction for details
5090: - ctx - the function context (or NULL)
5092: Level: advanced
5094: .seealso: SNESSetNGS(), SNESGetFunction()
5095: @*/
5097: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
5098: {
5100: DM dm;
5104: SNESGetDM(snes,&dm);
5105: DMSNESGetNGS(dm,f,ctx);
5106: return(0);
5107: }
5109: /*@C
5110: SNESSetOptionsPrefix - Sets the prefix used for searching for all
5111: SNES options in the database.
5113: Logically Collective on SNES
5115: Input Parameters:
5116: + snes - the SNES context
5117: - prefix - the prefix to prepend to all option names
5119: Notes:
5120: A hyphen (-) must NOT be given at the beginning of the prefix name.
5121: The first character of all runtime options is AUTOMATICALLY the hyphen.
5123: Level: advanced
5125: .seealso: SNESSetFromOptions()
5126: @*/
5127: PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[])
5128: {
5133: PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
5134: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
5135: if (snes->linesearch) {
5136: SNESGetLineSearch(snes,&snes->linesearch);
5137: PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
5138: }
5139: KSPSetOptionsPrefix(snes->ksp,prefix);
5140: return(0);
5141: }
5143: /*@C
5144: SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
5145: SNES options in the database.
5147: Logically Collective on SNES
5149: Input Parameters:
5150: + snes - the SNES context
5151: - prefix - the prefix to prepend to all option names
5153: Notes:
5154: A hyphen (-) must NOT be given at the beginning of the prefix name.
5155: The first character of all runtime options is AUTOMATICALLY the hyphen.
5157: Level: advanced
5159: .seealso: SNESGetOptionsPrefix()
5160: @*/
5161: PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[])
5162: {
5167: PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
5168: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
5169: if (snes->linesearch) {
5170: SNESGetLineSearch(snes,&snes->linesearch);
5171: PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
5172: }
5173: KSPAppendOptionsPrefix(snes->ksp,prefix);
5174: return(0);
5175: }
5177: /*@C
5178: SNESGetOptionsPrefix - Sets the prefix used for searching for all
5179: SNES options in the database.
5181: Not Collective
5183: Input Parameter:
5184: . snes - the SNES context
5186: Output Parameter:
5187: . prefix - pointer to the prefix string used
5189: Notes:
5190: On the fortran side, the user should pass in a string 'prefix' of
5191: sufficient length to hold the prefix.
5193: Level: advanced
5195: .seealso: SNESAppendOptionsPrefix()
5196: @*/
5197: PetscErrorCode SNESGetOptionsPrefix(SNES snes,const char *prefix[])
5198: {
5203: PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
5204: return(0);
5205: }
5207: /*@C
5208: SNESRegister - Adds a method to the nonlinear solver package.
5210: Not collective
5212: Input Parameters:
5213: + name_solver - name of a new user-defined solver
5214: - routine_create - routine to create method context
5216: Notes:
5217: SNESRegister() may be called multiple times to add several user-defined solvers.
5219: Sample usage:
5220: .vb
5221: SNESRegister("my_solver",MySolverCreate);
5222: .ve
5224: Then, your solver can be chosen with the procedural interface via
5225: $ SNESSetType(snes,"my_solver")
5226: or at runtime via the option
5227: $ -snes_type my_solver
5229: Level: advanced
5231: Note: If your function is not being put into a shared library then use SNESRegister() instead
5233: .seealso: SNESRegisterAll(), SNESRegisterDestroy()
5235: Level: advanced
5236: @*/
5237: PetscErrorCode SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
5238: {
5242: SNESInitializePackage();
5243: PetscFunctionListAdd(&SNESList,sname,function);
5244: return(0);
5245: }
5247: PetscErrorCode SNESTestLocalMin(SNES snes)
5248: {
5250: PetscInt N,i,j;
5251: Vec u,uh,fh;
5252: PetscScalar value;
5253: PetscReal norm;
5256: SNESGetSolution(snes,&u);
5257: VecDuplicate(u,&uh);
5258: VecDuplicate(u,&fh);
5260: /* currently only works for sequential */
5261: PetscPrintf(PetscObjectComm((PetscObject)snes),"Testing FormFunction() for local min\n");
5262: VecGetSize(u,&N);
5263: for (i=0; i<N; i++) {
5264: VecCopy(u,uh);
5265: PetscPrintf(PetscObjectComm((PetscObject)snes),"i = %D\n",i);
5266: for (j=-10; j<11; j++) {
5267: value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
5268: VecSetValue(uh,i,value,ADD_VALUES);
5269: SNESComputeFunction(snes,uh,fh);
5270: VecNorm(fh,NORM_2,&norm);
5271: PetscPrintf(PetscObjectComm((PetscObject)snes)," j norm %D %18.16e\n",j,norm);
5272: value = -value;
5273: VecSetValue(uh,i,value,ADD_VALUES);
5274: }
5275: }
5276: VecDestroy(&uh);
5277: VecDestroy(&fh);
5278: return(0);
5279: }
5281: /*@
5282: SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
5283: computing relative tolerance for linear solvers within an inexact
5284: Newton method.
5286: Logically Collective on SNES
5288: Input Parameters:
5289: + snes - SNES context
5290: - flag - PETSC_TRUE or PETSC_FALSE
5292: Options Database:
5293: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
5294: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
5295: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
5296: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
5297: . -snes_ksp_ew_gamma <gamma> - Sets gamma
5298: . -snes_ksp_ew_alpha <alpha> - Sets alpha
5299: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
5300: - -snes_ksp_ew_threshold <threshold> - Sets threshold
5302: Notes:
5303: Currently, the default is to use a constant relative tolerance for
5304: the inner linear solvers. Alternatively, one can use the
5305: Eisenstat-Walker method, where the relative convergence tolerance
5306: is reset at each Newton iteration according progress of the nonlinear
5307: solver.
5309: Level: advanced
5311: Reference:
5312: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
5313: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
5315: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
5316: @*/
5317: PetscErrorCode SNESKSPSetUseEW(SNES snes,PetscBool flag)
5318: {
5322: snes->ksp_ewconv = flag;
5323: return(0);
5324: }
5326: /*@
5327: SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
5328: for computing relative tolerance for linear solvers within an
5329: inexact Newton method.
5331: Not Collective
5333: Input Parameter:
5334: . snes - SNES context
5336: Output Parameter:
5337: . flag - PETSC_TRUE or PETSC_FALSE
5339: Notes:
5340: Currently, the default is to use a constant relative tolerance for
5341: the inner linear solvers. Alternatively, one can use the
5342: Eisenstat-Walker method, where the relative convergence tolerance
5343: is reset at each Newton iteration according progress of the nonlinear
5344: solver.
5346: Level: advanced
5348: Reference:
5349: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
5350: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
5352: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
5353: @*/
5354: PetscErrorCode SNESKSPGetUseEW(SNES snes, PetscBool *flag)
5355: {
5359: *flag = snes->ksp_ewconv;
5360: return(0);
5361: }
5363: /*@
5364: SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
5365: convergence criteria for the linear solvers within an inexact
5366: Newton method.
5368: Logically Collective on SNES
5370: Input Parameters:
5371: + snes - SNES context
5372: . version - version 1, 2 (default is 2) or 3
5373: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5374: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5375: . gamma - multiplicative factor for version 2 rtol computation
5376: (0 <= gamma2 <= 1)
5377: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
5378: . alpha2 - power for safeguard
5379: - threshold - threshold for imposing safeguard (0 < threshold < 1)
5381: Note:
5382: Version 3 was contributed by Luis Chacon, June 2006.
5384: Use PETSC_DEFAULT to retain the default for any of the parameters.
5386: Level: advanced
5388: Reference:
5389: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
5390: inexact Newton method", Utah State University Math. Stat. Dept. Res.
5391: Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
5393: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
5394: @*/
5395: PetscErrorCode SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
5396: {
5397: SNESKSPEW *kctx;
5401: kctx = (SNESKSPEW*)snes->kspconvctx;
5402: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
5411: if (version != PETSC_DEFAULT) kctx->version = version;
5412: if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0;
5413: if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max;
5414: if (gamma != PETSC_DEFAULT) kctx->gamma = gamma;
5415: if (alpha != PETSC_DEFAULT) kctx->alpha = alpha;
5416: if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2;
5417: if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
5419: 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);
5420: 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);
5421: 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);
5422: 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);
5423: 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);
5424: 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);
5425: return(0);
5426: }
5428: /*@
5429: SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
5430: convergence criteria for the linear solvers within an inexact
5431: Newton method.
5433: Not Collective
5435: Input Parameter:
5436: . snes - SNES context
5438: Output Parameters:
5439: + version - version 1, 2 (default is 2) or 3
5440: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5441: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5442: . gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
5443: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
5444: . alpha2 - power for safeguard
5445: - threshold - threshold for imposing safeguard (0 < threshold < 1)
5447: Level: advanced
5449: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
5450: @*/
5451: PetscErrorCode SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
5452: {
5453: SNESKSPEW *kctx;
5457: kctx = (SNESKSPEW*)snes->kspconvctx;
5458: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
5459: if (version) *version = kctx->version;
5460: if (rtol_0) *rtol_0 = kctx->rtol_0;
5461: if (rtol_max) *rtol_max = kctx->rtol_max;
5462: if (gamma) *gamma = kctx->gamma;
5463: if (alpha) *alpha = kctx->alpha;
5464: if (alpha2) *alpha2 = kctx->alpha2;
5465: if (threshold) *threshold = kctx->threshold;
5466: return(0);
5467: }
5469: PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5470: {
5472: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
5473: PetscReal rtol = PETSC_DEFAULT,stol;
5476: if (!snes->ksp_ewconv) return(0);
5477: if (!snes->iter) {
5478: rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
5479: VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
5480: }
5481: else {
5482: if (kctx->version == 1) {
5483: rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
5484: if (rtol < 0.0) rtol = -rtol;
5485: stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
5486: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5487: } else if (kctx->version == 2) {
5488: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5489: stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
5490: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5491: } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
5492: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5493: /* safeguard: avoid sharp decrease of rtol */
5494: stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
5495: stol = PetscMax(rtol,stol);
5496: rtol = PetscMin(kctx->rtol_0,stol);
5497: /* safeguard: avoid oversolving */
5498: stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
5499: stol = PetscMax(rtol,stol);
5500: rtol = PetscMin(kctx->rtol_0,stol);
5501: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
5502: }
5503: /* safeguard: avoid rtol greater than one */
5504: rtol = PetscMin(rtol,kctx->rtol_max);
5505: KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
5506: PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
5507: return(0);
5508: }
5510: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5511: {
5513: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
5514: PCSide pcside;
5515: Vec lres;
5518: if (!snes->ksp_ewconv) return(0);
5519: KSPGetTolerances(ksp,&kctx->rtol_last,NULL,NULL,NULL);
5520: kctx->norm_last = snes->norm;
5521: if (kctx->version == 1) {
5522: PC pc;
5523: PetscBool isNone;
5525: KSPGetPC(ksp, &pc);
5526: PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
5527: KSPGetPCSide(ksp,&pcside);
5528: if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
5529: /* KSP residual is true linear residual */
5530: KSPGetResidualNorm(ksp,&kctx->lresid_last);
5531: } else {
5532: /* KSP residual is preconditioned residual */
5533: /* compute true linear residual norm */
5534: VecDuplicate(b,&lres);
5535: MatMult(snes->jacobian,x,lres);
5536: VecAYPX(lres,-1.0,b);
5537: VecNorm(lres,NORM_2,&kctx->lresid_last);
5538: VecDestroy(&lres);
5539: }
5540: }
5541: return(0);
5542: }
5544: /*@
5545: SNESGetKSP - Returns the KSP context for a SNES solver.
5547: Not Collective, but if SNES object is parallel, then KSP object is parallel
5549: Input Parameter:
5550: . snes - the SNES context
5552: Output Parameter:
5553: . ksp - the KSP context
5555: Notes:
5556: The user can then directly manipulate the KSP context to set various
5557: options, etc. Likewise, the user can then extract and manipulate the
5558: PC contexts as well.
5560: Level: beginner
5562: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
5563: @*/
5564: PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp)
5565: {
5572: if (!snes->ksp) {
5573: KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
5574: PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
5575: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);
5577: KSPSetPreSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPreSolve_SNESEW,snes);
5578: KSPSetPostSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPostSolve_SNESEW,snes);
5580: KSPMonitorSetFromOptions(snes->ksp, "-snes_monitor_ksp", "snes_preconditioned_residual", snes);
5581: PetscObjectSetOptions((PetscObject)snes->ksp,((PetscObject)snes)->options);
5582: }
5583: *ksp = snes->ksp;
5584: return(0);
5585: }
5587: #include <petsc/private/dmimpl.h>
5588: /*@
5589: SNESSetDM - Sets the DM that may be used by some nonlinear solvers or their underlying preconditioners
5591: Logically Collective on SNES
5593: Input Parameters:
5594: + snes - the nonlinear solver context
5595: - dm - the dm, cannot be NULL
5597: Notes:
5598: A DM can only be used for solving one problem at a time because information about the problem is stored on the DM,
5599: even when not using interfaces like DMSNESSetFunction(). Use DMClone() to get a distinct DM when solving different
5600: problems using the same function space.
5602: Level: intermediate
5604: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
5605: @*/
5606: PetscErrorCode SNESSetDM(SNES snes,DM dm)
5607: {
5609: KSP ksp;
5610: DMSNES sdm;
5615: PetscObjectReference((PetscObject)dm);
5616: if (snes->dm) { /* Move the DMSNES context over to the new DM unless the new DM already has one */
5617: if (snes->dm->dmsnes && !dm->dmsnes) {
5618: DMCopyDMSNES(snes->dm,dm);
5619: DMGetDMSNES(snes->dm,&sdm);
5620: if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
5621: }
5622: DMCoarsenHookRemove(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
5623: DMDestroy(&snes->dm);
5624: }
5625: snes->dm = dm;
5626: snes->dmAuto = PETSC_FALSE;
5628: SNESGetKSP(snes,&ksp);
5629: KSPSetDM(ksp,dm);
5630: KSPSetDMActive(ksp,PETSC_FALSE);
5631: if (snes->npc) {
5632: SNESSetDM(snes->npc, snes->dm);
5633: SNESSetNPCSide(snes,snes->npcside);
5634: }
5635: return(0);
5636: }
5638: /*@
5639: SNESGetDM - Gets the DM that may be used by some preconditioners
5641: Not Collective but DM obtained is parallel on SNES
5643: Input Parameter:
5644: . snes - the preconditioner context
5646: Output Parameter:
5647: . dm - the dm
5649: Level: intermediate
5651: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
5652: @*/
5653: PetscErrorCode SNESGetDM(SNES snes,DM *dm)
5654: {
5659: if (!snes->dm) {
5660: DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
5661: snes->dmAuto = PETSC_TRUE;
5662: }
5663: *dm = snes->dm;
5664: return(0);
5665: }
5667: /*@
5668: SNESSetNPC - Sets the nonlinear preconditioner to be used.
5670: Collective on SNES
5672: Input Parameters:
5673: + snes - iterative context obtained from SNESCreate()
5674: - pc - the preconditioner object
5676: Notes:
5677: Use SNESGetNPC() to retrieve the preconditioner context (for example,
5678: to configure it using the API).
5680: Level: developer
5682: .seealso: SNESGetNPC(), SNESHasNPC()
5683: @*/
5684: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
5685: {
5692: PetscObjectReference((PetscObject) pc);
5693: SNESDestroy(&snes->npc);
5694: snes->npc = pc;
5695: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->npc);
5696: return(0);
5697: }
5699: /*@
5700: SNESGetNPC - Creates a nonlinear preconditioning solver (SNES) to be used to precondition the nonlinear solver.
5702: Not Collective; but any changes to the obtained SNES object must be applied collectively
5704: Input Parameter:
5705: . snes - iterative context obtained from SNESCreate()
5707: Output Parameter:
5708: . pc - preconditioner context
5710: Options Database:
5711: . -npc_snes_type <type> - set the type of the SNES to use as the nonlinear preconditioner
5713: Notes:
5714: If a SNES was previously set with SNESSetNPC() then that SNES is returned, otherwise a new SNES object is created.
5716: The (preconditioner) SNES returned automatically inherits the same nonlinear function and Jacobian supplied to the original
5717: SNES during SNESSetUp()
5719: Level: developer
5721: .seealso: SNESSetNPC(), SNESHasNPC(), SNES, SNESCreate()
5722: @*/
5723: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
5724: {
5726: const char *optionsprefix;
5731: if (!snes->npc) {
5732: SNESCreate(PetscObjectComm((PetscObject)snes),&snes->npc);
5733: PetscObjectIncrementTabLevel((PetscObject)snes->npc,(PetscObject)snes,1);
5734: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->npc);
5735: SNESGetOptionsPrefix(snes,&optionsprefix);
5736: SNESSetOptionsPrefix(snes->npc,optionsprefix);
5737: SNESAppendOptionsPrefix(snes->npc,"npc_");
5738: SNESSetCountersReset(snes->npc,PETSC_FALSE);
5739: }
5740: *pc = snes->npc;
5741: return(0);
5742: }
5744: /*@
5745: SNESHasNPC - Returns whether a nonlinear preconditioner exists
5747: Not Collective
5749: Input Parameter:
5750: . snes - iterative context obtained from SNESCreate()
5752: Output Parameter:
5753: . has_npc - whether the SNES has an NPC or not
5755: Level: developer
5757: .seealso: SNESSetNPC(), SNESGetNPC()
5758: @*/
5759: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
5760: {
5763: *has_npc = (PetscBool) (snes->npc ? PETSC_TRUE : PETSC_FALSE);
5764: return(0);
5765: }
5767: /*@
5768: SNESSetNPCSide - Sets the preconditioning side.
5770: Logically Collective on SNES
5772: Input Parameter:
5773: . snes - iterative context obtained from SNESCreate()
5775: Output Parameter:
5776: . side - the preconditioning side, where side is one of
5777: .vb
5778: PC_LEFT - left preconditioning
5779: PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5780: .ve
5782: Options Database Keys:
5783: . -snes_npc_side <right,left>
5785: Notes:
5786: SNESNRICHARDSON and SNESNCG only support left preconditioning.
5788: Level: intermediate
5790: .seealso: SNESGetNPCSide(), KSPSetPCSide()
5791: @*/
5792: PetscErrorCode SNESSetNPCSide(SNES snes,PCSide side)
5793: {
5797: snes->npcside= side;
5798: return(0);
5799: }
5801: /*@
5802: SNESGetNPCSide - Gets the preconditioning side.
5804: Not Collective
5806: Input Parameter:
5807: . snes - iterative context obtained from SNESCreate()
5809: Output Parameter:
5810: . side - the preconditioning side, where side is one of
5811: .vb
5812: PC_LEFT - left preconditioning
5813: PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5814: .ve
5816: Level: intermediate
5818: .seealso: SNESSetNPCSide(), KSPGetPCSide()
5819: @*/
5820: PetscErrorCode SNESGetNPCSide(SNES snes,PCSide *side)
5821: {
5825: *side = snes->npcside;
5826: return(0);
5827: }
5829: /*@
5830: SNESSetLineSearch - Sets the linesearch on the SNES instance.
5832: Collective on SNES
5834: Input Parameters:
5835: + snes - iterative context obtained from SNESCreate()
5836: - linesearch - the linesearch object
5838: Notes:
5839: Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
5840: to configure it using the API).
5842: Level: developer
5844: .seealso: SNESGetLineSearch()
5845: @*/
5846: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
5847: {
5854: PetscObjectReference((PetscObject) linesearch);
5855: SNESLineSearchDestroy(&snes->linesearch);
5857: snes->linesearch = linesearch;
5859: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5860: return(0);
5861: }
5863: /*@
5864: SNESGetLineSearch - Returns a pointer to the line search context set with SNESSetLineSearch()
5865: or creates a default line search instance associated with the SNES and returns it.
5867: Not Collective
5869: Input Parameter:
5870: . snes - iterative context obtained from SNESCreate()
5872: Output Parameter:
5873: . linesearch - linesearch context
5875: Level: beginner
5877: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5878: @*/
5879: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5880: {
5882: const char *optionsprefix;
5887: if (!snes->linesearch) {
5888: SNESGetOptionsPrefix(snes, &optionsprefix);
5889: SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5890: SNESLineSearchSetSNES(snes->linesearch, snes);
5891: SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5892: PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5893: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5894: }
5895: *linesearch = snes->linesearch;
5896: return(0);
5897: }