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 <true,false> - cause an immediate error condition and stop the program if the solver does not converge
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: {
39: snes->errorifnotconverged = flg;
40: return 0;
41: }
43: /*@
44: SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge?
46: Not Collective
48: Input Parameter:
49: . snes - iterative context obtained from SNESCreate()
51: Output Parameter:
52: . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE
54: Level: intermediate
56: .seealso: SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIfNotConverged()
57: @*/
58: PetscErrorCode SNESGetErrorIfNotConverged(SNES snes,PetscBool *flag)
59: {
62: *flag = snes->errorifnotconverged;
63: return 0;
64: }
66: /*@
67: SNESSetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?
69: Logically Collective on SNES
71: Input Parameters:
72: + snes - the shell SNES
73: - flg - is the residual computed?
75: Level: advanced
77: .seealso: SNESGetAlwaysComputesFinalResidual()
78: @*/
79: PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES snes, PetscBool flg)
80: {
82: snes->alwayscomputesfinalresidual = flg;
83: return 0;
84: }
86: /*@
87: SNESGetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?
89: Logically Collective on SNES
91: Input Parameter:
92: . snes - the shell SNES
94: Output Parameter:
95: . flg - is the residual computed?
97: Level: advanced
99: .seealso: SNESSetAlwaysComputesFinalResidual()
100: @*/
101: PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES snes, PetscBool *flg)
102: {
104: *flg = snes->alwayscomputesfinalresidual;
105: return 0;
106: }
108: /*@
109: SNESSetFunctionDomainError - tells SNES that the input vector to your SNESFunction is not
110: in the functions domain. For example, negative pressure.
112: Logically Collective on SNES
114: Input Parameters:
115: . snes - the SNES context
117: Level: advanced
119: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction
120: @*/
121: PetscErrorCode SNESSetFunctionDomainError(SNES snes)
122: {
125: snes->domainerror = PETSC_TRUE;
126: return 0;
127: }
129: /*@
130: SNESSetJacobianDomainError - tells SNES that computeJacobian does not make sense any more. For example there is a negative element transformation.
132: Logically Collective on SNES
134: Input Parameters:
135: . snes - the SNES context
137: Level: advanced
139: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError()
140: @*/
141: PetscErrorCode SNESSetJacobianDomainError(SNES snes)
142: {
145: snes->jacobiandomainerror = PETSC_TRUE;
146: return 0;
147: }
149: /*@
150: SNESSetCheckJacobianDomainError - if or not to check jacobian domain error after each Jacobian evaluation. By default, we check Jacobian domain error
151: in the debug mode, and do not check it in the optimized mode.
153: Logically Collective on SNES
155: Input Parameters:
156: + snes - the SNES context
157: - flg - indicates if or not to check jacobian domain error after each Jacobian evaluation
159: Level: advanced
161: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESGetCheckJacobianDomainError()
162: @*/
163: PetscErrorCode SNESSetCheckJacobianDomainError(SNES snes, PetscBool flg)
164: {
166: snes->checkjacdomainerror = flg;
167: return 0;
168: }
170: /*@
171: SNESGetCheckJacobianDomainError - Get an indicator whether or not we are checking Jacobian domain errors after each Jacobian evaluation.
173: Logically Collective on SNES
175: Input Parameters:
176: . snes - the SNES context
178: Output Parameters:
179: . flg - PETSC_FALSE indicates that we don't check jacobian domain errors after each Jacobian evaluation
181: Level: advanced
183: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction(), SNESSetFunctionDomainError(), SNESSetCheckJacobianDomainError()
184: @*/
185: PetscErrorCode SNESGetCheckJacobianDomainError(SNES snes, PetscBool *flg)
186: {
189: *flg = snes->checkjacdomainerror;
190: return 0;
191: }
193: /*@
194: SNESGetFunctionDomainError - Gets the status of the domain error after a call to SNESComputeFunction;
196: Logically Collective on SNES
198: Input Parameters:
199: . snes - the SNES context
201: Output Parameters:
202: . domainerror - Set to PETSC_TRUE if there's a domain error; PETSC_FALSE otherwise.
204: Level: advanced
206: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction()
207: @*/
208: PetscErrorCode SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror)
209: {
212: *domainerror = snes->domainerror;
213: return 0;
214: }
216: /*@
217: SNESGetJacobianDomainError - Gets the status of the Jacobian domain error after a call to SNESComputeJacobian;
219: Logically Collective on SNES
221: Input Parameters:
222: . snes - the SNES context
224: Output Parameters:
225: . domainerror - Set to PETSC_TRUE if there's a jacobian domain error; PETSC_FALSE otherwise.
227: Level: advanced
229: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction(),SNESGetFunctionDomainError()
230: @*/
231: PetscErrorCode SNESGetJacobianDomainError(SNES snes, PetscBool *domainerror)
232: {
235: *domainerror = snes->jacobiandomainerror;
236: return 0;
237: }
239: /*@C
240: SNESLoad - Loads a SNES that has been stored in binary with SNESView().
242: Collective on PetscViewer
244: Input Parameters:
245: + newdm - the newly loaded SNES, this needs to have been created with SNESCreate() or
246: some related function before a call to SNESLoad().
247: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
249: Level: intermediate
251: Notes:
252: The type is determined by the data in the file, any type set into the SNES before this call is ignored.
254: Notes for advanced users:
255: Most users should not need to know the details of the binary storage
256: format, since SNESLoad() and TSView() completely hide these details.
257: But for anyone who's interested, the standard binary matrix storage
258: format is
259: .vb
260: has not yet been determined
261: .ve
263: .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad()
264: @*/
265: PetscErrorCode SNESLoad(SNES snes, PetscViewer viewer)
266: {
267: PetscBool isbinary;
268: PetscInt classid;
269: char type[256];
270: KSP ksp;
271: DM dm;
272: DMSNES dmsnes;
276: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
279: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
281: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
282: SNESSetType(snes, type);
283: if (snes->ops->load) {
284: (*snes->ops->load)(snes,viewer);
285: }
286: SNESGetDM(snes,&dm);
287: DMGetDMSNES(dm,&dmsnes);
288: DMSNESLoad(dmsnes,viewer);
289: SNESGetKSP(snes,&ksp);
290: KSPLoad(ksp,viewer);
291: return 0;
292: }
294: #include <petscdraw.h>
295: #if defined(PETSC_HAVE_SAWS)
296: #include <petscviewersaws.h>
297: #endif
299: /*@C
300: SNESViewFromOptions - View from Options
302: Collective on SNES
304: Input Parameters:
305: + A - the application ordering context
306: . obj - Optional object
307: - name - command line option
309: Level: intermediate
310: .seealso: SNES, SNESView, PetscObjectViewFromOptions(), SNESCreate()
311: @*/
312: PetscErrorCode SNESViewFromOptions(SNES A,PetscObject obj,const char name[])
313: {
315: PetscObjectViewFromOptions((PetscObject)A,obj,name);
316: return 0;
317: }
319: PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES,Vec,Mat,Mat,void*);
321: /*@C
322: SNESView - Prints the SNES data structure.
324: Collective on SNES
326: Input Parameters:
327: + SNES - the SNES context
328: - viewer - visualization context
330: Options Database Key:
331: . -snes_view - Calls SNESView() at end of SNESSolve()
333: Notes:
334: The available visualization contexts include
335: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
336: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
337: output where only the first processor opens
338: the file. All other processors send their
339: data to the first processor to print.
341: The available formats include
342: + PETSC_VIEWER_DEFAULT - standard output (default)
343: - PETSC_VIEWER_ASCII_INFO_DETAIL - more verbose output for SNESNASM
345: The user can open an alternative visualization context with
346: PetscViewerASCIIOpen() - output to a specified file.
348: In the debugger you can do "call SNESView(snes,0)" to display the SNES solver. (The same holds for any PETSc object viewer).
350: Level: beginner
352: .seealso: PetscViewerASCIIOpen()
353: @*/
354: PetscErrorCode SNESView(SNES snes,PetscViewer viewer)
355: {
356: SNESKSPEW *kctx;
357: KSP ksp;
358: SNESLineSearch linesearch;
359: PetscBool iascii,isstring,isbinary,isdraw;
360: DMSNES dmsnes;
361: #if defined(PETSC_HAVE_SAWS)
362: PetscBool issaws;
363: #endif
366: if (!viewer) {
367: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
368: }
372: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
373: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
374: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
375: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
376: #if defined(PETSC_HAVE_SAWS)
377: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
378: #endif
379: if (iascii) {
380: SNESNormSchedule normschedule;
381: DM dm;
382: PetscErrorCode (*cJ)(SNES,Vec,Mat,Mat,void*);
383: void *ctx;
384: const char *pre = "";
386: PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);
387: if (!snes->setupcalled) {
388: PetscViewerASCIIPrintf(viewer," SNES has not been set up so information may be incomplete\n");
389: }
390: if (snes->ops->view) {
391: PetscViewerASCIIPushTab(viewer);
392: (*snes->ops->view)(snes,viewer);
393: PetscViewerASCIIPopTab(viewer);
394: }
395: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);
396: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);
397: if (snes->usesksp) {
398: PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);
399: }
400: PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);
401: SNESGetNormSchedule(snes, &normschedule);
402: if (normschedule > 0) PetscViewerASCIIPrintf(viewer," norm schedule %s\n",SNESNormSchedules[normschedule]);
403: if (snes->gridsequence) {
404: PetscViewerASCIIPrintf(viewer," total number of grid sequence refinements=%D\n",snes->gridsequence);
405: }
406: if (snes->ksp_ewconv) {
407: kctx = (SNESKSPEW*)snes->kspconvctx;
408: if (kctx) {
409: PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);
410: PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);
411: PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);
412: }
413: }
414: if (snes->lagpreconditioner == -1) {
415: PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");
416: } else if (snes->lagpreconditioner > 1) {
417: PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);
418: }
419: if (snes->lagjacobian == -1) {
420: PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");
421: } else if (snes->lagjacobian > 1) {
422: PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);
423: }
424: SNESGetDM(snes,&dm);
425: DMSNESGetJacobian(dm,&cJ,&ctx);
426: if (snes->mf_operator) {
427: PetscViewerASCIIPrintf(viewer," Jacobian is applied matrix-free with differencing\n");
428: pre = "Preconditioning ";
429: }
430: if (cJ == SNESComputeJacobianDefault) {
431: PetscViewerASCIIPrintf(viewer," %sJacobian is built using finite differences one column at a time\n",pre);
432: } else if (cJ == SNESComputeJacobianDefaultColor) {
433: PetscViewerASCIIPrintf(viewer," %sJacobian is built using finite differences with coloring\n",pre);
434: /* it slightly breaks data encapsulation for access the DMDA information directly */
435: } else if (cJ == SNESComputeJacobian_DMDA) {
436: MatFDColoring fdcoloring;
437: PetscObjectQuery((PetscObject)dm,"DMDASNES_FDCOLORING",(PetscObject*)&fdcoloring);
438: if (fdcoloring) {
439: PetscViewerASCIIPrintf(viewer," %sJacobian is built using colored finite differences on a DMDA\n",pre);
440: } else {
441: PetscViewerASCIIPrintf(viewer," %sJacobian is built using a DMDA local Jacobian\n",pre);
442: }
443: } else if (snes->mf) {
444: PetscViewerASCIIPrintf(viewer," Jacobian is applied matrix-free with differencing, no explicit Jacobian\n");
445: }
446: } else if (isstring) {
447: const char *type;
448: SNESGetType(snes,&type);
449: PetscViewerStringSPrintf(viewer," SNESType: %-7.7s",type);
450: if (snes->ops->view) (*snes->ops->view)(snes,viewer);
451: } else if (isbinary) {
452: PetscInt classid = SNES_FILE_CLASSID;
453: MPI_Comm comm;
454: PetscMPIInt rank;
455: char type[256];
457: PetscObjectGetComm((PetscObject)snes,&comm);
458: MPI_Comm_rank(comm,&rank);
459: if (rank == 0) {
460: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT);
461: PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));
462: PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR);
463: }
464: if (snes->ops->view) {
465: (*snes->ops->view)(snes,viewer);
466: }
467: } else if (isdraw) {
468: PetscDraw draw;
469: char str[36];
470: PetscReal x,y,bottom,h;
472: PetscViewerDrawGetDraw(viewer,0,&draw);
473: PetscDrawGetCurrentPoint(draw,&x,&y);
474: PetscStrncpy(str,"SNES: ",sizeof(str));
475: PetscStrlcat(str,((PetscObject)snes)->type_name,sizeof(str));
476: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);
477: bottom = y - h;
478: PetscDrawPushCurrentPoint(draw,x,bottom);
479: if (snes->ops->view) {
480: (*snes->ops->view)(snes,viewer);
481: }
482: #if defined(PETSC_HAVE_SAWS)
483: } else if (issaws) {
484: PetscMPIInt rank;
485: const char *name;
487: PetscObjectGetName((PetscObject)snes,&name);
488: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
489: if (!((PetscObject)snes)->amsmem && rank == 0) {
490: char dir[1024];
492: PetscObjectViewSAWs((PetscObject)snes,viewer);
493: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
494: PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT));
495: if (!snes->conv_hist) {
496: SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);
497: }
498: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);
499: PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE));
500: }
501: #endif
502: }
503: if (snes->linesearch) {
504: SNESGetLineSearch(snes, &linesearch);
505: PetscViewerASCIIPushTab(viewer);
506: SNESLineSearchView(linesearch, viewer);
507: PetscViewerASCIIPopTab(viewer);
508: }
509: if (snes->npc && snes->usesnpc) {
510: PetscViewerASCIIPushTab(viewer);
511: SNESView(snes->npc, viewer);
512: PetscViewerASCIIPopTab(viewer);
513: }
514: PetscViewerASCIIPushTab(viewer);
515: DMGetDMSNES(snes->dm,&dmsnes);
516: DMSNESView(dmsnes, viewer);
517: PetscViewerASCIIPopTab(viewer);
518: if (snes->usesksp) {
519: SNESGetKSP(snes,&ksp);
520: PetscViewerASCIIPushTab(viewer);
521: KSPView(ksp,viewer);
522: PetscViewerASCIIPopTab(viewer);
523: }
524: if (isdraw) {
525: PetscDraw draw;
526: PetscViewerDrawGetDraw(viewer,0,&draw);
527: PetscDrawPopCurrentPoint(draw);
528: }
529: return 0;
530: }
532: /*
533: We retain a list of functions that also take SNES command
534: line options. These are called at the end SNESSetFromOptions()
535: */
536: #define MAXSETFROMOPTIONS 5
537: static PetscInt numberofsetfromoptions;
538: static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
540: /*@C
541: SNESAddOptionsChecker - Adds an additional function to check for SNES options.
543: Not Collective
545: Input Parameter:
546: . snescheck - function that checks for options
548: Level: developer
550: .seealso: SNESSetFromOptions()
551: @*/
552: PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
553: {
555: othersetfromoptions[numberofsetfromoptions++] = snescheck;
556: return 0;
557: }
559: PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
561: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
562: {
563: Mat J;
564: MatNullSpace nullsp;
568: if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
569: Mat A = snes->jacobian, B = snes->jacobian_pre;
570: MatCreateVecs(A ? A : B, NULL,&snes->vec_func);
571: }
573: if (version == 1) {
574: MatCreateSNESMF(snes,&J);
575: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
576: MatSetFromOptions(J);
577: /* TODO: the version 2 code should be merged into the MatCreateSNESMF() and MatCreateMFFD() infrastructure and then removed */
578: } else if (version == 2) {
580: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
581: SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);
582: #else
583: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator routines (version 2)");
584: #endif
585: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator routines, only version 1 and 2");
587: /* attach any user provided null space that was on Amat to the newly created matrix free matrix */
588: if (snes->jacobian) {
589: MatGetNullSpace(snes->jacobian,&nullsp);
590: if (nullsp) {
591: MatSetNullSpace(J,nullsp);
592: }
593: }
595: PetscInfo(snes,"Setting default matrix-free operator routines (version %D)\n", version);
596: if (hasOperator) {
598: /* This version replaces the user provided Jacobian matrix with a
599: matrix-free version but still employs the user-provided preconditioner matrix. */
600: SNESSetJacobian(snes,J,NULL,NULL,NULL);
601: } else {
602: /* This version replaces both the user-provided Jacobian and the user-
603: provided preconditioner Jacobian with the default matrix free version. */
604: if (snes->npcside == PC_LEFT && snes->npc) {
605: if (!snes->jacobian) SNESSetJacobian(snes,J,NULL,NULL,NULL);
606: } else {
607: KSP ksp;
608: PC pc;
609: PetscBool match;
611: SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,NULL);
612: /* Force no preconditioner */
613: SNESGetKSP(snes,&ksp);
614: KSPGetPC(ksp,&pc);
615: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);
616: if (!match) {
617: PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");
618: PCSetType(pc,PCNONE);
619: }
620: }
621: }
622: MatDestroy(&J);
623: return 0;
624: }
626: static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx)
627: {
628: SNES snes = (SNES)ctx;
629: Vec Xfine,Xfine_named = NULL,Xcoarse;
631: if (PetscLogPrintInfo) {
632: PetscInt finelevel,coarselevel,fineclevel,coarseclevel;
633: DMGetRefineLevel(dmfine,&finelevel);
634: DMGetCoarsenLevel(dmfine,&fineclevel);
635: DMGetRefineLevel(dmcoarse,&coarselevel);
636: DMGetCoarsenLevel(dmcoarse,&coarseclevel);
637: PetscInfo(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);
638: }
639: if (dmfine == snes->dm) Xfine = snes->vec_sol;
640: else {
641: DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
642: Xfine = Xfine_named;
643: }
644: DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
645: if (Inject) {
646: MatRestrict(Inject,Xfine,Xcoarse);
647: } else {
648: MatRestrict(Restrict,Xfine,Xcoarse);
649: VecPointwiseMult(Xcoarse,Xcoarse,Rscale);
650: }
651: DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
652: if (Xfine_named) DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
653: return 0;
654: }
656: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
657: {
658: DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
659: return 0;
660: }
662: /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can
663: * safely call SNESGetDM() in their residual evaluation routine. */
664: static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx)
665: {
666: SNES snes = (SNES)ctx;
667: Vec X,Xnamed = NULL;
668: DM dmsave;
669: void *ctxsave;
670: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL;
672: dmsave = snes->dm;
673: KSPGetDM(ksp,&snes->dm);
674: if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */
675: else { /* We are on a coarser level, this vec was initialized using a DM restrict hook */
676: DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
677: X = Xnamed;
678: SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);
679: /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */
680: if (jac == SNESComputeJacobianDefaultColor) {
681: SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,NULL);
682: }
683: }
684: /* Make sure KSP DM has the Jacobian computation routine */
685: {
686: DMSNES sdm;
688: DMGetDMSNES(snes->dm, &sdm);
689: if (!sdm->ops->computejacobian) {
690: DMCopyDMSNES(dmsave, snes->dm);
691: }
692: }
693: /* Compute the operators */
694: SNESComputeJacobian(snes,X,A,B);
695: /* Put the previous context back */
696: if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
697: SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
698: }
700: if (Xnamed) DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
701: snes->dm = dmsave;
702: return 0;
703: }
705: /*@
706: SNESSetUpMatrices - ensures that matrices are available for SNES, to be called by SNESSetUp_XXX()
708: Collective
710: Input Parameter:
711: . snes - snes to configure
713: Level: developer
715: .seealso: SNESSetUp()
716: @*/
717: PetscErrorCode SNESSetUpMatrices(SNES snes)
718: {
719: DM dm;
720: DMSNES sdm;
722: SNESGetDM(snes,&dm);
723: DMGetDMSNES(dm,&sdm);
724: if (!snes->jacobian && snes->mf) {
725: Mat J;
726: void *functx;
727: MatCreateSNESMF(snes,&J);
728: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
729: MatSetFromOptions(J);
730: SNESGetFunction(snes,NULL,NULL,&functx);
731: SNESSetJacobian(snes,J,J,NULL,NULL);
732: MatDestroy(&J);
733: } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
734: Mat J,B;
735: MatCreateSNESMF(snes,&J);
736: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
737: MatSetFromOptions(J);
738: DMCreateMatrix(snes->dm,&B);
739: /* sdm->computejacobian was already set to reach here */
740: SNESSetJacobian(snes,J,B,NULL,NULL);
741: MatDestroy(&J);
742: MatDestroy(&B);
743: } else if (!snes->jacobian_pre) {
744: PetscDS prob;
745: Mat J, B;
746: PetscBool hasPrec = PETSC_FALSE;
748: J = snes->jacobian;
749: DMGetDS(dm, &prob);
750: if (prob) PetscDSHasJacobianPreconditioner(prob, &hasPrec);
751: if (J) PetscObjectReference((PetscObject) J);
752: else if (hasPrec) DMCreateMatrix(snes->dm, &J);
753: DMCreateMatrix(snes->dm, &B);
754: SNESSetJacobian(snes, J ? J : B, B, NULL, NULL);
755: MatDestroy(&J);
756: MatDestroy(&B);
757: }
758: {
759: KSP ksp;
760: SNESGetKSP(snes,&ksp);
761: KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);
762: DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
763: }
764: return 0;
765: }
767: static PetscErrorCode SNESMonitorPauseFinal_Internal(SNES snes)
768: {
769: PetscInt i;
771: if (!snes->pauseFinal) return 0;
772: for (i = 0; i < snes->numbermonitors; ++i) {
773: PetscViewerAndFormat *vf = (PetscViewerAndFormat *) snes->monitorcontext[i];
774: PetscDraw draw;
775: PetscReal lpause;
777: if (!vf) continue;
778: if (vf->lg) {
780: if (((PetscObject) vf->lg)->classid != PETSC_DRAWLG_CLASSID) continue;
781: PetscDrawLGGetDraw(vf->lg, &draw);
782: PetscDrawGetPause(draw, &lpause);
783: PetscDrawSetPause(draw, -1.0);
784: PetscDrawPause(draw);
785: PetscDrawSetPause(draw, lpause);
786: } else {
787: PetscBool isdraw;
790: if (((PetscObject) vf->viewer)->classid != PETSC_VIEWER_CLASSID) continue;
791: PetscObjectTypeCompare((PetscObject) vf->viewer, PETSCVIEWERDRAW, &isdraw);
792: if (!isdraw) continue;
793: PetscViewerDrawGetDraw(vf->viewer, 0, &draw);
794: PetscDrawGetPause(draw, &lpause);
795: PetscDrawSetPause(draw, -1.0);
796: PetscDrawPause(draw);
797: PetscDrawSetPause(draw, lpause);
798: }
799: }
800: return 0;
801: }
803: /*@C
804: SNESMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
806: Collective on SNES
808: Input Parameters:
809: + snes - SNES object you wish to monitor
810: . name - the monitor type one is seeking
811: . help - message indicating what monitoring is done
812: . manual - manual page for the monitor
813: . monitor - the monitor function
814: - 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
816: Level: developer
818: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
819: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
820: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
821: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
822: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
823: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
824: PetscOptionsFList(), PetscOptionsEList()
825: @*/
826: PetscErrorCode SNESMonitorSetFromOptions(SNES snes,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNES,PetscViewerAndFormat*))
827: {
828: PetscViewer viewer;
829: PetscViewerFormat format;
830: PetscBool flg;
832: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,name,&viewer,&format,&flg);
833: if (flg) {
834: PetscViewerAndFormat *vf;
835: PetscViewerAndFormatCreate(viewer,format,&vf);
836: PetscObjectDereference((PetscObject)viewer);
837: if (monitorsetup) {
838: (*monitorsetup)(snes,vf);
839: }
840: SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
841: }
842: return 0;
843: }
845: /*@
846: SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
848: Collective on SNES
850: Input Parameter:
851: . snes - the SNES context
853: Options Database Keys:
854: + -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
855: . -snes_stol - convergence tolerance in terms of the norm
856: of the change in the solution between steps
857: . -snes_atol <abstol> - absolute tolerance of residual norm
858: . -snes_rtol <rtol> - relative decrease in tolerance norm from initial
859: . -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
860: . -snes_force_iteration <force> - force SNESSolve() to take at least one iteration
861: . -snes_max_it <max_it> - maximum number of iterations
862: . -snes_max_funcs <max_funcs> - maximum number of function evaluations
863: . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
864: . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
865: . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
866: . -snes_lag_preconditioner_persists <true,false> - retains the -snes_lag_preconditioner information across multiple SNESSolve()
867: . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
868: . -snes_lag_jacobian_persists <true,false> - retains the -snes_lag_jacobian information across multiple SNESSolve()
869: . -snes_trtol <trtol> - trust region tolerance
870: . -snes_convergence_test - <default,skip,correct_pressure> convergence test in nonlinear solver.
871: default SNESConvergedDefault(). skip SNESConvergedSkip() means continue iterating until max_it or some other criterion is reached, saving expense
872: of convergence test. correct_pressure SNESConvergedCorrectPressure() has special handling of a pressure null space.
873: . -snes_monitor [ascii][:filename][:viewer format] - prints residual norm at each iteration. if no filename given prints to stdout
874: . -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration
875: . -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration
876: . -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration
877: . -snes_monitor_lg_residualnorm - plots residual norm at each iteration
878: . -snes_monitor_lg_range - plots residual norm at each iteration
879: . -snes_monitor_pause_final - Pauses all monitor drawing after the solver ends
880: . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
881: . -snes_fd_color - use finite differences with coloring to compute Jacobian
882: . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
883: . -snes_converged_reason - print the reason for convergence/divergence after each solve
884: . -npc_snes_type <type> - the SNES type to use as a nonlinear preconditioner
885: . -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.
886: - -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.
888: Options Database for Eisenstat-Walker method:
889: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
890: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
891: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
892: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
893: . -snes_ksp_ew_gamma <gamma> - Sets gamma
894: . -snes_ksp_ew_alpha <alpha> - Sets alpha
895: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
896: - -snes_ksp_ew_threshold <threshold> - Sets threshold
898: Notes:
899: To see all options, run your program with the -help option or consult the users manual
901: Notes:
902: SNES supports three approaches for computing (approximate) Jacobians: user provided via SNESSetJacobian(), matrix free, and computing explicitly with
903: finite differences and coloring using MatFDColoring. It is also possible to use automatic differentiation and the MatFDColoring object.
905: Level: beginner
907: .seealso: SNESSetOptionsPrefix(), SNESResetFromOptions(), SNES, SNESCreate()
908: @*/
909: PetscErrorCode SNESSetFromOptions(SNES snes)
910: {
911: PetscBool flg,pcset,persist,set;
912: PetscInt i,indx,lag,grids;
913: const char *deft = SNESNEWTONLS;
914: const char *convtests[] = {"default","skip","correct_pressure"};
915: SNESKSPEW *kctx = NULL;
916: char type[256], monfilename[PETSC_MAX_PATH_LEN];
918: PCSide pcside;
919: const char *optionsprefix;
922: SNESRegisterAll();
923: PetscObjectOptionsBegin((PetscObject)snes);
924: if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
925: PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
926: if (flg) {
927: SNESSetType(snes,type);
928: } else if (!((PetscObject)snes)->type_name) {
929: SNESSetType(snes,deft);
930: }
931: PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);
932: PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);
934: PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,NULL);
935: PetscOptionsReal("-snes_divergence_tolerance","Stop if residual norm increases by this factor","SNESSetDivergenceTolerance",snes->divtol,&snes->divtol,NULL);
936: PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,NULL);
937: PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,NULL);
938: PetscOptionsInt("-snes_max_fail","Maximum nonlinear step failures","SNESSetMaxNonlinearStepFailures",snes->maxFailures,&snes->maxFailures,NULL);
939: PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,NULL);
940: PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,NULL);
941: PetscOptionsBool("-snes_force_iteration","Force SNESSolve() to take at least one iteration","SNESSetForceIteration",snes->forceiteration,&snes->forceiteration,NULL);
942: PetscOptionsBool("-snes_check_jacobian_domain_error","Check Jacobian domain error after Jacobian evaluation","SNESCheckJacobianDomainError",snes->checkjacdomainerror,&snes->checkjacdomainerror,NULL);
944: PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
945: if (flg) {
947: SNESSetLagPreconditioner(snes,lag);
948: }
949: PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple SNES solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
950: if (flg) {
951: SNESSetLagPreconditionerPersists(snes,persist);
952: }
953: PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
954: if (flg) {
956: SNESSetLagJacobian(snes,lag);
957: }
958: PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple SNES solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
959: if (flg) {
960: SNESSetLagJacobianPersists(snes,persist);
961: }
963: PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
964: if (flg) {
965: SNESSetGridSequence(snes,grids);
966: }
968: PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,sizeof(convtests)/sizeof(char*),"default",&indx,&flg);
969: if (flg) {
970: switch (indx) {
971: case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
972: case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL); break;
973: case 2: SNESSetConvergenceTest(snes,SNESConvergedCorrectPressure,NULL,NULL); break;
974: }
975: }
977: PetscOptionsEList("-snes_norm_schedule","SNES Norm schedule","SNESSetNormSchedule",SNESNormSchedules,5,"function",&indx,&flg);
978: if (flg) SNESSetNormSchedule(snes,(SNESNormSchedule)indx);
980: PetscOptionsEList("-snes_function_type","SNES Norm schedule","SNESSetFunctionType",SNESFunctionTypes,2,"unpreconditioned",&indx,&flg);
981: if (flg) SNESSetFunctionType(snes,(SNESFunctionType)indx);
983: kctx = (SNESKSPEW*)snes->kspconvctx;
985: PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,NULL);
987: PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);
988: PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);
989: PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);
990: PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);
991: PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);
992: PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);
993: PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);
995: flg = PETSC_FALSE;
996: PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);
997: if (set && flg) SNESMonitorCancel(snes);
999: SNESMonitorSetFromOptions(snes,"-snes_monitor","Monitor norm of function","SNESMonitorDefault",SNESMonitorDefault,SNESMonitorDefaultSetUp);
1000: SNESMonitorSetFromOptions(snes,"-snes_monitor_short","Monitor norm of function with fewer digits","SNESMonitorDefaultShort",SNESMonitorDefaultShort,NULL);
1001: SNESMonitorSetFromOptions(snes,"-snes_monitor_range","Monitor range of elements of function","SNESMonitorRange",SNESMonitorRange,NULL);
1003: SNESMonitorSetFromOptions(snes,"-snes_monitor_ratio","Monitor ratios of the norm of function for consecutive steps","SNESMonitorRatio",SNESMonitorRatio,SNESMonitorRatioSetUp);
1004: SNESMonitorSetFromOptions(snes,"-snes_monitor_field","Monitor norm of function (split into fields)","SNESMonitorDefaultField",SNESMonitorDefaultField,NULL);
1005: SNESMonitorSetFromOptions(snes,"-snes_monitor_solution","View solution at each iteration","SNESMonitorSolution",SNESMonitorSolution,NULL);
1006: SNESMonitorSetFromOptions(snes,"-snes_monitor_solution_update","View correction at each iteration","SNESMonitorSolutionUpdate",SNESMonitorSolutionUpdate,NULL);
1007: SNESMonitorSetFromOptions(snes,"-snes_monitor_residual","View residual at each iteration","SNESMonitorResidual",SNESMonitorResidual,NULL);
1008: SNESMonitorSetFromOptions(snes,"-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",SNESMonitorJacUpdateSpectrum,NULL);
1009: SNESMonitorSetFromOptions(snes,"-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet",SNESMonitorFields,NULL);
1010: PetscOptionsBool("-snes_monitor_pause_final", "Pauses all draw monitors at the final iterate", "SNESMonitorPauseFinal_Internal", PETSC_FALSE, &snes->pauseFinal, NULL);
1012: PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",NULL,monfilename,sizeof(monfilename),&flg);
1013: if (flg) PetscPythonMonitorSet((PetscObject)snes,monfilename);
1015: flg = PETSC_FALSE;
1016: PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
1017: if (flg) {
1018: PetscViewer ctx;
1020: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
1021: SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
1022: }
1024: flg = PETSC_FALSE;
1025: PetscOptionsBool("-snes_converged_reason_view_cancel","Remove all converged reason viewers","SNESConvergedReasonViewCancel",flg,&flg,&set);
1026: if (set && flg) SNESConvergedReasonViewCancel(snes);
1028: flg = PETSC_FALSE;
1029: PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
1030: if (flg) {
1031: void *functx;
1032: DM dm;
1033: DMSNES sdm;
1034: SNESGetDM(snes,&dm);
1035: DMGetDMSNES(dm,&sdm);
1036: sdm->jacobianctx = NULL;
1037: SNESGetFunction(snes,NULL,NULL,&functx);
1038: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
1039: PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
1040: }
1042: flg = PETSC_FALSE;
1043: PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
1044: if (flg) {
1045: SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
1046: }
1048: flg = PETSC_FALSE;
1049: PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);
1050: if (flg) {
1051: DM dm;
1052: DMSNES sdm;
1053: SNESGetDM(snes,&dm);
1054: DMGetDMSNES(dm,&sdm);
1055: sdm->jacobianctx = NULL;
1056: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,NULL);
1057: PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
1058: }
1060: flg = PETSC_FALSE;
1061: PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf_operator,&flg);
1062: if (flg && snes->mf_operator) {
1063: snes->mf_operator = PETSC_TRUE;
1064: snes->mf = PETSC_TRUE;
1065: }
1066: flg = PETSC_FALSE;
1067: PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf,&flg);
1068: if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE;
1069: PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,NULL);
1071: flg = PETSC_FALSE;
1072: SNESGetNPCSide(snes,&pcside);
1073: PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
1074: if (flg) SNESSetNPCSide(snes,pcside);
1076: #if defined(PETSC_HAVE_SAWS)
1077: /*
1078: Publish convergence information using SAWs
1079: */
1080: flg = PETSC_FALSE;
1081: PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);
1082: if (flg) {
1083: void *ctx;
1084: SNESMonitorSAWsCreate(snes,&ctx);
1085: SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);
1086: }
1087: #endif
1088: #if defined(PETSC_HAVE_SAWS)
1089: {
1090: PetscBool set;
1091: flg = PETSC_FALSE;
1092: PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);
1093: if (set) {
1094: PetscObjectSAWsSetBlock((PetscObject)snes,flg);
1095: }
1096: }
1097: #endif
1099: for (i = 0; i < numberofsetfromoptions; i++) {
1100: (*othersetfromoptions[i])(snes);
1101: }
1103: if (snes->ops->setfromoptions) {
1104: (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
1105: }
1107: /* process any options handlers added with PetscObjectAddOptionsHandler() */
1108: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)snes);
1109: PetscOptionsEnd();
1111: if (snes->linesearch) {
1112: SNESGetLineSearch(snes, &snes->linesearch);
1113: SNESLineSearchSetFromOptions(snes->linesearch);
1114: }
1116: if (snes->usesksp) {
1117: if (!snes->ksp) SNESGetKSP(snes,&snes->ksp);
1118: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
1119: KSPSetFromOptions(snes->ksp);
1120: }
1122: /* if user has set the SNES NPC type via options database, create it. */
1123: SNESGetOptionsPrefix(snes, &optionsprefix);
1124: PetscOptionsHasName(((PetscObject)snes)->options,optionsprefix, "-npc_snes_type", &pcset);
1125: if (pcset && (!snes->npc)) {
1126: SNESGetNPC(snes, &snes->npc);
1127: }
1128: if (snes->npc) {
1129: SNESSetFromOptions(snes->npc);
1130: }
1131: snes->setfromoptionscalled++;
1132: return 0;
1133: }
1135: /*@
1136: SNESResetFromOptions - Sets various SNES and KSP parameters from user options ONLY if the SNES was previously set from options
1138: Collective on SNES
1140: Input Parameter:
1141: . snes - the SNES context
1143: Level: beginner
1145: .seealso: SNESSetFromOptions(), SNESSetOptionsPrefix()
1146: @*/
1147: PetscErrorCode SNESResetFromOptions(SNES snes)
1148: {
1149: if (snes->setfromoptionscalled) SNESSetFromOptions(snes);
1150: return 0;
1151: }
1153: /*@C
1154: SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
1155: the nonlinear solvers.
1157: Logically Collective on SNES
1159: Input Parameters:
1160: + snes - the SNES context
1161: . compute - function to compute the context
1162: - destroy - function to destroy the context
1164: Level: intermediate
1166: Notes:
1167: This function is currently not available from Fortran.
1169: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
1170: @*/
1171: PetscErrorCode SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
1172: {
1174: snes->ops->usercompute = compute;
1175: snes->ops->userdestroy = destroy;
1176: return 0;
1177: }
1179: /*@
1180: SNESSetApplicationContext - Sets the optional user-defined context for
1181: the nonlinear solvers.
1183: Logically Collective on SNES
1185: Input Parameters:
1186: + snes - the SNES context
1187: - usrP - optional user context
1189: Level: intermediate
1191: Fortran Notes:
1192: To use this from Fortran you must write a Fortran interface definition for this
1193: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1195: .seealso: SNESGetApplicationContext()
1196: @*/
1197: PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP)
1198: {
1199: KSP ksp;
1202: SNESGetKSP(snes,&ksp);
1203: KSPSetApplicationContext(ksp,usrP);
1204: snes->user = usrP;
1205: return 0;
1206: }
1208: /*@
1209: SNESGetApplicationContext - Gets the user-defined context for the
1210: nonlinear solvers.
1212: Not Collective
1214: Input Parameter:
1215: . snes - SNES context
1217: Output Parameter:
1218: . usrP - user context
1220: Fortran Notes:
1221: To use this from Fortran you must write a Fortran interface definition for this
1222: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1224: Level: intermediate
1226: .seealso: SNESSetApplicationContext()
1227: @*/
1228: PetscErrorCode SNESGetApplicationContext(SNES snes,void *usrP)
1229: {
1231: *(void**)usrP = snes->user;
1232: return 0;
1233: }
1235: /*@
1236: SNESSetUseMatrixFree - indicates that SNES should use matrix free finite difference matrix vector products internally to apply the Jacobian.
1238: Collective on SNES
1240: Input Parameters:
1241: + snes - SNES context
1242: . mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used
1243: - mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1245: Options Database:
1246: + -snes_mf - use matrix free for both the mat and pmat operator
1247: . -snes_mf_operator - use matrix free only for the mat operator
1248: . -snes_fd_color - compute the Jacobian via coloring and finite differences.
1249: - -snes_fd - compute the Jacobian via finite differences (slow)
1251: Level: intermediate
1253: Notes:
1254: SNES supports three approaches for computing (approximate) Jacobians: user provided via SNESSetJacobian(), matrix free, and computing explicitly with
1255: finite differences and coloring using MatFDColoring. It is also possible to use automatic differentiation and the MatFDColoring object.
1257: .seealso: SNESGetUseMatrixFree(), MatCreateSNESMF(), SNESComputeJacobianDefaultColor()
1258: @*/
1259: PetscErrorCode SNESSetUseMatrixFree(SNES snes,PetscBool mf_operator,PetscBool mf)
1260: {
1264: snes->mf = mf_operator ? PETSC_TRUE : mf;
1265: snes->mf_operator = mf_operator;
1266: return 0;
1267: }
1269: /*@
1270: SNESGetUseMatrixFree - indicates if the SNES uses matrix free finite difference matrix vector products to apply the Jacobian.
1272: Collective on SNES
1274: Input Parameter:
1275: . snes - SNES context
1277: Output Parameters:
1278: + mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used
1279: - mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1281: Options Database:
1282: + -snes_mf - use matrix free for both the mat and pmat operator
1283: - -snes_mf_operator - use matrix free only for the mat operator
1285: Level: intermediate
1287: .seealso: SNESSetUseMatrixFree(), MatCreateSNESMF()
1288: @*/
1289: PetscErrorCode SNESGetUseMatrixFree(SNES snes,PetscBool *mf_operator,PetscBool *mf)
1290: {
1292: if (mf) *mf = snes->mf;
1293: if (mf_operator) *mf_operator = snes->mf_operator;
1294: return 0;
1295: }
1297: /*@
1298: SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1299: at this time.
1301: Not Collective
1303: Input Parameter:
1304: . snes - SNES context
1306: Output Parameter:
1307: . iter - iteration number
1309: Notes:
1310: For example, during the computation of iteration 2 this would return 1.
1312: This is useful for using lagged Jacobians (where one does not recompute the
1313: Jacobian at each SNES iteration). For example, the code
1314: .vb
1315: SNESGetIterationNumber(snes,&it);
1316: if (!(it % 2)) {
1317: [compute Jacobian here]
1318: }
1319: .ve
1320: can be used in your ComputeJacobian() function to cause the Jacobian to be
1321: recomputed every second SNES iteration.
1323: After the SNES solve is complete this will return the number of nonlinear iterations used.
1325: Level: intermediate
1327: .seealso: SNESGetLinearSolveIterations()
1328: @*/
1329: PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt *iter)
1330: {
1333: *iter = snes->iter;
1334: return 0;
1335: }
1337: /*@
1338: SNESSetIterationNumber - Sets the current iteration number.
1340: Not Collective
1342: Input Parameters:
1343: + snes - SNES context
1344: - iter - iteration number
1346: Level: developer
1348: .seealso: SNESGetLinearSolveIterations()
1349: @*/
1350: PetscErrorCode SNESSetIterationNumber(SNES snes,PetscInt iter)
1351: {
1353: PetscObjectSAWsTakeAccess((PetscObject)snes);
1354: snes->iter = iter;
1355: PetscObjectSAWsGrantAccess((PetscObject)snes);
1356: return 0;
1357: }
1359: /*@
1360: SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1361: attempted by the nonlinear solver.
1363: Not Collective
1365: Input Parameter:
1366: . snes - SNES context
1368: Output Parameter:
1369: . nfails - number of unsuccessful steps attempted
1371: Notes:
1372: This counter is reset to zero for each successive call to SNESSolve().
1374: Level: intermediate
1376: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1377: SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1378: @*/
1379: PetscErrorCode SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1380: {
1383: *nfails = snes->numFailures;
1384: return 0;
1385: }
1387: /*@
1388: SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1389: attempted by the nonlinear solver before it gives up.
1391: Not Collective
1393: Input Parameters:
1394: + snes - SNES context
1395: - maxFails - maximum of unsuccessful steps
1397: Level: intermediate
1399: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1400: SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1401: @*/
1402: PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1403: {
1405: snes->maxFailures = maxFails;
1406: return 0;
1407: }
1409: /*@
1410: SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1411: attempted by the nonlinear solver before it gives up.
1413: Not Collective
1415: Input Parameter:
1416: . snes - SNES context
1418: Output Parameter:
1419: . maxFails - maximum of unsuccessful steps
1421: Level: intermediate
1423: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1424: SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1426: @*/
1427: PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1428: {
1431: *maxFails = snes->maxFailures;
1432: return 0;
1433: }
1435: /*@
1436: SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1437: done by SNES.
1439: Not Collective
1441: Input Parameter:
1442: . snes - SNES context
1444: Output Parameter:
1445: . nfuncs - number of evaluations
1447: Level: intermediate
1449: Notes:
1450: Reset every time SNESSolve is called unless SNESSetCountersReset() is used.
1452: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1453: @*/
1454: PetscErrorCode SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1455: {
1458: *nfuncs = snes->nfuncs;
1459: return 0;
1460: }
1462: /*@
1463: SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1464: linear solvers.
1466: Not Collective
1468: Input Parameter:
1469: . snes - SNES context
1471: Output Parameter:
1472: . nfails - number of failed solves
1474: Level: intermediate
1476: Options Database Keys:
1477: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1479: Notes:
1480: This counter is reset to zero for each successive call to SNESSolve().
1482: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1483: @*/
1484: PetscErrorCode SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1485: {
1488: *nfails = snes->numLinearSolveFailures;
1489: return 0;
1490: }
1492: /*@
1493: SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1494: allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
1496: Logically Collective on SNES
1498: Input Parameters:
1499: + snes - SNES context
1500: - maxFails - maximum allowed linear solve failures
1502: Level: intermediate
1504: Options Database Keys:
1505: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1507: Notes:
1508: By default this is 0; that is SNES returns on the first failed linear solve
1510: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1511: @*/
1512: PetscErrorCode SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1513: {
1516: snes->maxLinearSolveFailures = maxFails;
1517: return 0;
1518: }
1520: /*@
1521: SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1522: are allowed before SNES terminates
1524: Not Collective
1526: Input Parameter:
1527: . snes - SNES context
1529: Output Parameter:
1530: . maxFails - maximum of unsuccessful solves allowed
1532: Level: intermediate
1534: Notes:
1535: By default this is 1; that is SNES returns on the first failed linear solve
1537: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1538: @*/
1539: PetscErrorCode SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1540: {
1543: *maxFails = snes->maxLinearSolveFailures;
1544: return 0;
1545: }
1547: /*@
1548: SNESGetLinearSolveIterations - Gets the total number of linear iterations
1549: used by the nonlinear solver.
1551: Not Collective
1553: Input Parameter:
1554: . snes - SNES context
1556: Output Parameter:
1557: . lits - number of linear iterations
1559: Notes:
1560: This counter is reset to zero for each successive call to SNESSolve() unless SNESSetCountersReset() is used.
1562: 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
1563: then call KSPGetIterationNumber() after the failed solve.
1565: Level: intermediate
1567: .seealso: SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1568: @*/
1569: PetscErrorCode SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1570: {
1573: *lits = snes->linear_its;
1574: return 0;
1575: }
1577: /*@
1578: SNESSetCountersReset - Sets whether or not the counters for linear iterations and function evaluations
1579: are reset every time SNESSolve() is called.
1581: Logically Collective on SNES
1583: Input Parameters:
1584: + snes - SNES context
1585: - reset - whether to reset the counters or not
1587: Notes:
1588: This defaults to PETSC_TRUE
1590: Level: developer
1592: .seealso: SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC()
1593: @*/
1594: PetscErrorCode SNESSetCountersReset(SNES snes,PetscBool reset)
1595: {
1598: snes->counters_reset = reset;
1599: return 0;
1600: }
1602: /*@
1603: SNESSetKSP - Sets a KSP context for the SNES object to use
1605: Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
1607: Input Parameters:
1608: + snes - the SNES context
1609: - ksp - the KSP context
1611: Notes:
1612: The SNES object already has its KSP object, you can obtain with SNESGetKSP()
1613: so this routine is rarely needed.
1615: The KSP object that is already in the SNES object has its reference count
1616: decreased by one.
1618: Level: developer
1620: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1621: @*/
1622: PetscErrorCode SNESSetKSP(SNES snes,KSP ksp)
1623: {
1627: PetscObjectReference((PetscObject)ksp);
1628: if (snes->ksp) PetscObjectDereference((PetscObject)snes->ksp);
1629: snes->ksp = ksp;
1630: return 0;
1631: }
1633: /* -----------------------------------------------------------*/
1634: /*@
1635: SNESCreate - Creates a nonlinear solver context.
1637: Collective
1639: Input Parameters:
1640: . comm - MPI communicator
1642: Output Parameter:
1643: . outsnes - the new SNES context
1645: Options Database Keys:
1646: + -snes_mf - Activates default matrix-free Jacobian-vector products,
1647: and no preconditioning matrix
1648: . -snes_mf_operator - Activates default matrix-free Jacobian-vector
1649: products, and a user-provided preconditioning matrix
1650: as set by SNESSetJacobian()
1651: - -snes_fd - Uses (slow!) finite differences to compute Jacobian
1653: Level: beginner
1655: Developer Notes:
1656: SNES always creates a KSP object even though many SNES methods do not use it. This is
1657: unfortunate and should be fixed at some point. The flag snes->usesksp indicates if the
1658: particular method does use KSP and regulates if the information about the KSP is printed
1659: in SNESView(). TSSetFromOptions() does call SNESSetFromOptions() which can lead to users being confused
1660: by help messages about meaningless SNES options.
1662: SNES always creates the snes->kspconvctx even though it is used by only one type. This should
1663: be fixed.
1665: .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner(), SNESSetLagJacobian()
1667: @*/
1668: PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes)
1669: {
1670: SNES snes;
1671: SNESKSPEW *kctx;
1674: *outsnes = NULL;
1675: SNESInitializePackage();
1677: PetscHeaderCreate(snes,SNES_CLASSID,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);
1679: snes->ops->converged = SNESConvergedDefault;
1680: snes->usesksp = PETSC_TRUE;
1681: snes->tolerancesset = PETSC_FALSE;
1682: snes->max_its = 50;
1683: snes->max_funcs = 10000;
1684: snes->norm = 0.0;
1685: snes->xnorm = 0.0;
1686: snes->ynorm = 0.0;
1687: snes->normschedule = SNES_NORM_ALWAYS;
1688: snes->functype = SNES_FUNCTION_DEFAULT;
1689: #if defined(PETSC_USE_REAL_SINGLE)
1690: snes->rtol = 1.e-5;
1691: #else
1692: snes->rtol = 1.e-8;
1693: #endif
1694: snes->ttol = 0.0;
1695: #if defined(PETSC_USE_REAL_SINGLE)
1696: snes->abstol = 1.e-25;
1697: #else
1698: snes->abstol = 1.e-50;
1699: #endif
1700: #if defined(PETSC_USE_REAL_SINGLE)
1701: snes->stol = 1.e-5;
1702: #else
1703: snes->stol = 1.e-8;
1704: #endif
1705: #if defined(PETSC_USE_REAL_SINGLE)
1706: snes->deltatol = 1.e-6;
1707: #else
1708: snes->deltatol = 1.e-12;
1709: #endif
1710: snes->divtol = 1.e4;
1711: snes->rnorm0 = 0;
1712: snes->nfuncs = 0;
1713: snes->numFailures = 0;
1714: snes->maxFailures = 1;
1715: snes->linear_its = 0;
1716: snes->lagjacobian = 1;
1717: snes->jac_iter = 0;
1718: snes->lagjac_persist = PETSC_FALSE;
1719: snes->lagpreconditioner = 1;
1720: snes->pre_iter = 0;
1721: snes->lagpre_persist = PETSC_FALSE;
1722: snes->numbermonitors = 0;
1723: snes->numberreasonviews = 0;
1724: snes->data = NULL;
1725: snes->setupcalled = PETSC_FALSE;
1726: snes->ksp_ewconv = PETSC_FALSE;
1727: snes->nwork = 0;
1728: snes->work = NULL;
1729: snes->nvwork = 0;
1730: snes->vwork = NULL;
1731: snes->conv_hist_len = 0;
1732: snes->conv_hist_max = 0;
1733: snes->conv_hist = NULL;
1734: snes->conv_hist_its = NULL;
1735: snes->conv_hist_reset = PETSC_TRUE;
1736: snes->counters_reset = PETSC_TRUE;
1737: snes->vec_func_init_set = PETSC_FALSE;
1738: snes->reason = SNES_CONVERGED_ITERATING;
1739: snes->npcside = PC_RIGHT;
1740: snes->setfromoptionscalled = 0;
1742: snes->mf = PETSC_FALSE;
1743: snes->mf_operator = PETSC_FALSE;
1744: snes->mf_version = 1;
1746: snes->numLinearSolveFailures = 0;
1747: snes->maxLinearSolveFailures = 1;
1749: snes->vizerotolerance = 1.e-8;
1750: snes->checkjacdomainerror = PetscDefined(USE_DEBUG) ? PETSC_TRUE : PETSC_FALSE;
1752: /* Set this to true if the implementation of SNESSolve_XXX does compute the residual at the final solution. */
1753: snes->alwayscomputesfinalresidual = PETSC_FALSE;
1755: /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
1756: PetscNewLog(snes,&kctx);
1758: snes->kspconvctx = (void*)kctx;
1759: kctx->version = 2;
1760: kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1761: this was too large for some test cases */
1762: kctx->rtol_last = 0.0;
1763: kctx->rtol_max = .9;
1764: kctx->gamma = 1.0;
1765: kctx->alpha = .5*(1.0 + PetscSqrtReal(5.0));
1766: kctx->alpha2 = kctx->alpha;
1767: kctx->threshold = .1;
1768: kctx->lresid_last = 0.0;
1769: kctx->norm_last = 0.0;
1771: *outsnes = snes;
1772: return 0;
1773: }
1775: /*MC
1776: SNESFunction - Functional form used to convey the nonlinear function to be solved by SNES
1778: Synopsis:
1779: #include "petscsnes.h"
1780: PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);
1782: Collective on snes
1784: Input Parameters:
1785: + snes - the SNES context
1786: . x - state at which to evaluate residual
1787: - ctx - optional user-defined function context, passed in with SNESSetFunction()
1789: Output Parameter:
1790: . f - vector to put residual (function value)
1792: Level: intermediate
1794: .seealso: SNESSetFunction(), SNESGetFunction()
1795: M*/
1797: /*@C
1798: SNESSetFunction - Sets the function evaluation routine and function
1799: vector for use by the SNES routines in solving systems of nonlinear
1800: equations.
1802: Logically Collective on SNES
1804: Input Parameters:
1805: + snes - the SNES context
1806: . r - vector to store function values, may be NULL
1807: . f - function evaluation routine; see SNESFunction for calling sequence details
1808: - ctx - [optional] user-defined context for private data for the
1809: function evaluation routine (may be NULL)
1811: Notes:
1812: The Newton-like methods typically solve linear systems of the form
1813: $ f'(x) x = -f(x),
1814: where f'(x) denotes the Jacobian matrix and f(x) is the function.
1816: Level: beginner
1818: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1819: @*/
1820: PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1821: {
1822: DM dm;
1825: if (r) {
1828: PetscObjectReference((PetscObject)r);
1829: VecDestroy(&snes->vec_func);
1830: snes->vec_func = r;
1831: }
1832: SNESGetDM(snes,&dm);
1833: DMSNESSetFunction(dm,f,ctx);
1834: if (f == SNESPicardComputeFunction) {
1835: DMSNESSetMFFunction(dm,SNESPicardComputeMFFunction,ctx);
1836: }
1837: return 0;
1838: }
1840: /*@C
1841: SNESSetInitialFunction - Sets the function vector to be used as the
1842: function norm at the initialization of the method. In some
1843: instances, the user has precomputed the function before calling
1844: SNESSolve. This function allows one to avoid a redundant call
1845: to SNESComputeFunction in that case.
1847: Logically Collective on SNES
1849: Input Parameters:
1850: + snes - the SNES context
1851: - f - vector to store function value
1853: Notes:
1854: This should not be modified during the solution procedure.
1856: This is used extensively in the SNESFAS hierarchy and in nonlinear preconditioning.
1858: Level: developer
1860: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1861: @*/
1862: PetscErrorCode SNESSetInitialFunction(SNES snes, Vec f)
1863: {
1864: Vec vec_func;
1869: if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1870: snes->vec_func_init_set = PETSC_FALSE;
1871: return 0;
1872: }
1873: SNESGetFunction(snes,&vec_func,NULL,NULL);
1874: VecCopy(f,vec_func);
1876: snes->vec_func_init_set = PETSC_TRUE;
1877: return 0;
1878: }
1880: /*@
1881: SNESSetNormSchedule - Sets the SNESNormSchedule used in convergence and monitoring
1882: of the SNES method.
1884: Logically Collective on SNES
1886: Input Parameters:
1887: + snes - the SNES context
1888: - normschedule - the frequency of norm computation
1890: Options Database Key:
1891: . -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - set the schedule
1893: Notes:
1894: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1895: of the nonlinear function and the taking of its norm at every iteration to
1896: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1897: (SNESNGS) and the like do not require the norm of the function to be computed, and therefore
1898: may either be monitored for convergence or not. As these are often used as nonlinear
1899: preconditioners, monitoring the norm of their error is not a useful enterprise within
1900: their solution.
1902: Level: developer
1904: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1905: @*/
1906: PetscErrorCode SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1907: {
1909: snes->normschedule = normschedule;
1910: return 0;
1911: }
1913: /*@
1914: SNESGetNormSchedule - Gets the SNESNormSchedule used in convergence and monitoring
1915: of the SNES method.
1917: Logically Collective on SNES
1919: Input Parameters:
1920: + snes - the SNES context
1921: - normschedule - the type of the norm used
1923: Level: advanced
1925: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1926: @*/
1927: PetscErrorCode SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1928: {
1930: *normschedule = snes->normschedule;
1931: return 0;
1932: }
1934: /*@
1935: SNESSetFunctionNorm - Sets the last computed residual norm.
1937: Logically Collective on SNES
1939: Input Parameters:
1940: + snes - the SNES context
1942: - normschedule - the frequency of norm computation
1944: Level: developer
1946: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1947: @*/
1948: PetscErrorCode SNESSetFunctionNorm(SNES snes, PetscReal norm)
1949: {
1951: snes->norm = norm;
1952: return 0;
1953: }
1955: /*@
1956: SNESGetFunctionNorm - Gets the last computed norm of the residual
1958: Not Collective
1960: Input Parameter:
1961: . snes - the SNES context
1963: Output Parameter:
1964: . norm - the last computed residual norm
1966: Level: developer
1968: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1969: @*/
1970: PetscErrorCode SNESGetFunctionNorm(SNES snes, PetscReal *norm)
1971: {
1974: *norm = snes->norm;
1975: return 0;
1976: }
1978: /*@
1979: SNESGetUpdateNorm - Gets the last computed norm of the Newton update
1981: Not Collective
1983: Input Parameter:
1984: . snes - the SNES context
1986: Output Parameter:
1987: . ynorm - the last computed update norm
1989: Level: developer
1991: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm()
1992: @*/
1993: PetscErrorCode SNESGetUpdateNorm(SNES snes, PetscReal *ynorm)
1994: {
1997: *ynorm = snes->ynorm;
1998: return 0;
1999: }
2001: /*@
2002: SNESGetSolutionNorm - Gets the last computed norm of the solution
2004: Not Collective
2006: Input Parameter:
2007: . snes - the SNES context
2009: Output Parameter:
2010: . xnorm - the last computed solution norm
2012: Level: developer
2014: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), SNESGetFunctionNorm(), SNESGetUpdateNorm()
2015: @*/
2016: PetscErrorCode SNESGetSolutionNorm(SNES snes, PetscReal *xnorm)
2017: {
2020: *xnorm = snes->xnorm;
2021: return 0;
2022: }
2024: /*@C
2025: SNESSetFunctionType - Sets the SNESNormSchedule used in convergence and monitoring
2026: of the SNES method.
2028: Logically Collective on SNES
2030: Input Parameters:
2031: + snes - the SNES context
2032: - normschedule - the frequency of norm computation
2034: Notes:
2035: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
2036: of the nonlinear function and the taking of its norm at every iteration to
2037: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
2038: (SNESNGS) and the like do not require the norm of the function to be computed, and therefore
2039: may either be monitored for convergence or not. As these are often used as nonlinear
2040: preconditioners, monitoring the norm of their error is not a useful enterprise within
2041: their solution.
2043: Level: developer
2045: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2046: @*/
2047: PetscErrorCode SNESSetFunctionType(SNES snes, SNESFunctionType type)
2048: {
2050: snes->functype = type;
2051: return 0;
2052: }
2054: /*@C
2055: SNESGetFunctionType - Gets the SNESNormSchedule used in convergence and monitoring
2056: of the SNES method.
2058: Logically Collective on SNES
2060: Input Parameters:
2061: + snes - the SNES context
2062: - normschedule - the type of the norm used
2064: Level: advanced
2066: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
2067: @*/
2068: PetscErrorCode SNESGetFunctionType(SNES snes, SNESFunctionType *type)
2069: {
2071: *type = snes->functype;
2072: return 0;
2073: }
2075: /*MC
2076: SNESNGSFunction - function used to convey a Gauss-Seidel sweep on the nonlinear function
2078: Synopsis:
2079: #include <petscsnes.h>
2080: $ SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);
2082: Collective on snes
2084: Input Parameters:
2085: + X - solution vector
2086: . B - RHS vector
2087: - ctx - optional user-defined Gauss-Seidel context
2089: Output Parameter:
2090: . X - solution vector
2092: Level: intermediate
2094: .seealso: SNESSetNGS(), SNESGetNGS()
2095: M*/
2097: /*@C
2098: SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
2099: use with composed nonlinear solvers.
2101: Input Parameters:
2102: + snes - the SNES context
2103: . f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
2104: - ctx - [optional] user-defined context for private data for the
2105: smoother evaluation routine (may be NULL)
2107: Notes:
2108: The NGS routines are used by the composed nonlinear solver to generate
2109: a problem appropriate update to the solution, particularly FAS.
2111: Level: intermediate
2113: .seealso: SNESGetFunction(), SNESComputeNGS()
2114: @*/
2115: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
2116: {
2117: DM dm;
2120: SNESGetDM(snes,&dm);
2121: DMSNESSetNGS(dm,f,ctx);
2122: return 0;
2123: }
2125: /*
2126: This is used for -snes_mf_operator; it uses a duplicate of snes->jacobian_pre because snes->jacobian_pre cannot be
2127: changed during the KSPSolve()
2128: */
2129: PetscErrorCode SNESPicardComputeMFFunction(SNES snes,Vec x,Vec f,void *ctx)
2130: {
2131: DM dm;
2132: DMSNES sdm;
2134: SNESGetDM(snes,&dm);
2135: DMGetDMSNES(dm,&sdm);
2137: /* A(x)*x - b(x) */
2138: if (sdm->ops->computepfunction) {
2139: PetscStackPush("SNES Picard user function");
2140: (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
2141: PetscStackPop;
2142: VecScale(f,-1.0);
2143: if (!snes->picard) {
2144: /* Cannot share nonzero pattern because of the possible use of SNESComputeJacobianDefault() */
2145: MatDuplicate(snes->jacobian_pre,MAT_DO_NOT_COPY_VALUES,&snes->picard);
2146: }
2147: PetscStackPush("SNES Picard user Jacobian");
2148: (*sdm->ops->computepjacobian)(snes,x,snes->picard,snes->picard,sdm->pctx);
2149: PetscStackPop;
2150: MatMultAdd(snes->picard,x,f,f);
2151: } else {
2152: PetscStackPush("SNES Picard user Jacobian");
2153: (*sdm->ops->computepjacobian)(snes,x,snes->picard,snes->picard,sdm->pctx);
2154: PetscStackPop;
2155: MatMult(snes->picard,x,f);
2156: }
2157: return 0;
2158: }
2160: PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
2161: {
2162: DM dm;
2163: DMSNES sdm;
2165: SNESGetDM(snes,&dm);
2166: DMGetDMSNES(dm,&sdm);
2168: /* A(x)*x - b(x) */
2169: if (sdm->ops->computepfunction) {
2170: PetscStackPush("SNES Picard user function");
2171: (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
2172: PetscStackPop;
2173: VecScale(f,-1.0);
2174: PetscStackPush("SNES Picard user Jacobian");
2175: (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
2176: PetscStackPop;
2177: MatMultAdd(snes->jacobian_pre,x,f,f);
2178: } else {
2179: PetscStackPush("SNES Picard user Jacobian");
2180: (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
2181: PetscStackPop;
2182: MatMult(snes->jacobian_pre,x,f);
2183: }
2184: return 0;
2185: }
2187: PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
2188: {
2189: /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
2190: /* must assembly if matrix-free to get the last SNES solution */
2191: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
2192: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
2193: return 0;
2194: }
2196: /*@C
2197: SNESSetPicard - Use SNES to solve the system A(x) x = bp(x) + b via a Picard type iteration (Picard linearization)
2199: Logically Collective on SNES
2201: Input Parameters:
2202: + snes - the SNES context
2203: . r - vector to store function values, may be NULL
2204: . bp - function evaluation routine, may be NULL
2205: . Amat - matrix with which A(x) x - bp(x) - b is to be computed
2206: . Pmat - matrix from which preconditioner is computed (usually the same as Amat)
2207: . J - function to compute matrix values, see SNESJacobianFunction() for details on its calling sequence
2208: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
2210: Notes:
2211: It is often better to provide the nonlinear function F() and some approximation to its Jacobian directly and use
2212: 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.
2214: One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both
2216: $ 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}
2217: $ 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.
2219: Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner.
2221: We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then
2222: the direct Picard iteration A(x^n) x^{n+1} = bp(x^n) + b
2224: 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
2225: 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
2226: different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).
2228: 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.
2230: When used with -snes_fd this will compute the true Jacobian (very slowly one column at at time) and thus represent Newton's method.
2232: 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
2233: 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
2234: 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.
2235: See the commment in src/snes/tutorials/ex15.c.
2237: Level: intermediate
2239: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
2240: @*/
2241: 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)
2242: {
2243: DM dm;
2246: SNESGetDM(snes, &dm);
2247: DMSNESSetPicard(dm,bp,J,ctx);
2248: DMSNESSetMFFunction(dm,SNESPicardComputeMFFunction,ctx);
2249: SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
2250: SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
2251: return 0;
2252: }
2254: /*@C
2255: SNESGetPicard - Returns the context for the Picard iteration
2257: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
2259: Input Parameter:
2260: . snes - the SNES context
2262: Output Parameters:
2263: + r - the function (or NULL)
2264: . f - the function (or NULL); see SNESFunction for calling sequence details
2265: . Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
2266: . Pmat - the matrix from which the preconditioner will be constructed (or NULL)
2267: . J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
2268: - ctx - the function context (or NULL)
2270: Level: advanced
2272: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
2273: @*/
2274: 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)
2275: {
2276: DM dm;
2279: SNESGetFunction(snes,r,NULL,NULL);
2280: SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
2281: SNESGetDM(snes,&dm);
2282: DMSNESGetPicard(dm,f,J,ctx);
2283: return 0;
2284: }
2286: /*@C
2287: SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem
2289: Logically Collective on SNES
2291: Input Parameters:
2292: + snes - the SNES context
2293: . func - function evaluation routine
2294: - ctx - [optional] user-defined context for private data for the
2295: function evaluation routine (may be NULL)
2297: Calling sequence of func:
2298: $ func (SNES snes,Vec x,void *ctx);
2300: . f - function vector
2301: - ctx - optional user-defined function context
2303: Level: intermediate
2305: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
2306: @*/
2307: PetscErrorCode SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
2308: {
2310: if (func) snes->ops->computeinitialguess = func;
2311: if (ctx) snes->initialguessP = ctx;
2312: return 0;
2313: }
2315: /* --------------------------------------------------------------- */
2316: /*@C
2317: SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
2318: it assumes a zero right hand side.
2320: Logically Collective on SNES
2322: Input Parameter:
2323: . snes - the SNES context
2325: Output Parameter:
2326: . rhs - the right hand side vector or NULL if the right hand side vector is null
2328: Level: intermediate
2330: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
2331: @*/
2332: PetscErrorCode SNESGetRhs(SNES snes,Vec *rhs)
2333: {
2336: *rhs = snes->vec_rhs;
2337: return 0;
2338: }
2340: /*@
2341: SNESComputeFunction - Calls the function that has been set with SNESSetFunction().
2343: Collective on SNES
2345: Input Parameters:
2346: + snes - the SNES context
2347: - x - input vector
2349: Output Parameter:
2350: . y - function vector, as set by SNESSetFunction()
2352: Notes:
2353: SNESComputeFunction() is typically used within nonlinear solvers
2354: implementations, so users would not generally call this routine themselves.
2356: Level: developer
2358: .seealso: SNESSetFunction(), SNESGetFunction(), SNESComputeMFFunction()
2359: @*/
2360: PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y)
2361: {
2362: DM dm;
2363: DMSNES sdm;
2370: VecValidValues(x,2,PETSC_TRUE);
2372: SNESGetDM(snes,&dm);
2373: DMGetDMSNES(dm,&sdm);
2374: if (sdm->ops->computefunction) {
2375: if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2376: PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2377: }
2378: VecLockReadPush(x);
2379: PetscStackPush("SNES user function");
2380: /* ensure domainerror is false prior to computefunction evaluation (may not have been reset) */
2381: snes->domainerror = PETSC_FALSE;
2382: (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2383: PetscStackPop;
2384: VecLockReadPop(x);
2385: if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2386: PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2387: }
2388: } else if (snes->vec_rhs) {
2389: MatMult(snes->jacobian, x, y);
2390: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2391: if (snes->vec_rhs) {
2392: VecAXPY(y,-1.0,snes->vec_rhs);
2393: }
2394: snes->nfuncs++;
2395: /*
2396: domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2397: propagate the value to all processes
2398: */
2399: if (snes->domainerror) {
2400: VecSetInf(y);
2401: }
2402: return 0;
2403: }
2405: /*@
2406: SNESComputeMFFunction - Calls the function that has been set with SNESSetMFFunction().
2408: Collective on SNES
2410: Input Parameters:
2411: + snes - the SNES context
2412: - x - input vector
2414: Output Parameter:
2415: . y - function vector, as set by SNESSetMFFunction()
2417: Notes:
2418: SNESComputeMFFunction() is used within the matrix vector products called by the matrix created with MatCreateSNESMF()
2419: so users would not generally call this routine themselves.
2421: Since this function is intended for use with finite differencing it does not subtract the right hand side vector provided with SNESSolve()
2422: while SNESComputeFunction() does. As such, this routine cannot be used with MatMFFDSetBase() with a provided F function value even if it applies the
2423: 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.
2425: Level: developer
2427: .seealso: SNESSetFunction(), SNESGetFunction(), SNESComputeFunction(), MatCreateSNESMF
2428: @*/
2429: PetscErrorCode SNESComputeMFFunction(SNES snes,Vec x,Vec y)
2430: {
2431: DM dm;
2432: DMSNES sdm;
2439: VecValidValues(x,2,PETSC_TRUE);
2441: SNESGetDM(snes,&dm);
2442: DMGetDMSNES(dm,&sdm);
2443: PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2444: VecLockReadPush(x);
2445: PetscStackPush("SNES user function");
2446: /* ensure domainerror is false prior to computefunction evaluation (may not have been reset) */
2447: snes->domainerror = PETSC_FALSE;
2448: (*sdm->ops->computemffunction)(snes,x,y,sdm->mffunctionctx);
2449: PetscStackPop;
2450: VecLockReadPop(x);
2451: PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2452: snes->nfuncs++;
2453: /*
2454: domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2455: propagate the value to all processes
2456: */
2457: if (snes->domainerror) {
2458: VecSetInf(y);
2459: }
2460: return 0;
2461: }
2463: /*@
2464: SNESComputeNGS - Calls the Gauss-Seidel function that has been set with SNESSetNGS().
2466: Collective on SNES
2468: Input Parameters:
2469: + snes - the SNES context
2470: . x - input vector
2471: - b - rhs vector
2473: Output Parameter:
2474: . x - new solution vector
2476: Notes:
2477: SNESComputeNGS() is typically used within composed nonlinear solver
2478: implementations, so most users would not generally call this routine
2479: themselves.
2481: Level: developer
2483: .seealso: SNESSetNGS(), SNESComputeFunction()
2484: @*/
2485: PetscErrorCode SNESComputeNGS(SNES snes,Vec b,Vec x)
2486: {
2487: DM dm;
2488: DMSNES sdm;
2495: if (b) VecValidValues(b,2,PETSC_TRUE);
2496: PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2497: SNESGetDM(snes,&dm);
2498: DMGetDMSNES(dm,&sdm);
2499: if (sdm->ops->computegs) {
2500: if (b) VecLockReadPush(b);
2501: PetscStackPush("SNES user NGS");
2502: (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2503: PetscStackPop;
2504: if (b) VecLockReadPop(b);
2505: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2506: PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2507: return 0;
2508: }
2510: PetscErrorCode SNESTestJacobian(SNES snes)
2511: {
2512: Mat A,B,C,D,jacobian;
2513: Vec x = snes->vec_sol,f = snes->vec_func;
2514: PetscErrorCode ierr;
2515: PetscReal nrm,gnorm;
2516: PetscReal threshold = 1.e-5;
2517: MatType mattype;
2518: PetscInt m,n,M,N;
2519: void *functx;
2520: PetscBool complete_print = PETSC_FALSE,threshold_print = PETSC_FALSE,test = PETSC_FALSE,flg,istranspose;
2521: PetscViewer viewer,mviewer;
2522: MPI_Comm comm;
2523: PetscInt tabs;
2524: static PetscBool directionsprinted = PETSC_FALSE;
2525: PetscViewerFormat format;
2527: PetscObjectOptionsBegin((PetscObject)snes);
2528: PetscOptionsName("-snes_test_jacobian","Compare hand-coded and finite difference Jacobians","None",&test);
2529: PetscOptionsReal("-snes_test_jacobian", "Threshold for element difference between hand-coded and finite difference being meaningful", "None", threshold, &threshold,NULL);
2530: PetscOptionsViewer("-snes_test_jacobian_view","View difference between hand-coded and finite difference Jacobians element entries","None",&mviewer,&format,&complete_print);
2531: if (!complete_print) {
2532: PetscOptionsDeprecated("-snes_test_jacobian_display","-snes_test_jacobian_view","3.13",NULL);
2533: PetscOptionsViewer("-snes_test_jacobian_display","Display difference between hand-coded and finite difference Jacobians","None",&mviewer,&format,&complete_print);
2534: }
2535: /* for compatibility with PETSc 3.9 and older. */
2536: PetscOptionsDeprecated("-snes_test_jacobian_display_threshold","-snes_test_jacobian","3.13","-snes_test_jacobian accepts an optional threshold (since v3.10)");
2537: PetscOptionsReal("-snes_test_jacobian_display_threshold", "Display difference between hand-coded and finite difference Jacobians which exceed input threshold", "None", threshold, &threshold, &threshold_print);
2538: PetscOptionsEnd();
2539: if (!test) return 0;
2541: PetscObjectGetComm((PetscObject)snes,&comm);
2542: PetscViewerASCIIGetStdout(comm,&viewer);
2543: PetscViewerASCIIGetTab(viewer, &tabs);
2544: PetscViewerASCIISetTab(viewer, ((PetscObject)snes)->tablevel);
2545: PetscViewerASCIIPrintf(viewer," ---------- Testing Jacobian -------------\n");
2546: if (!complete_print && !directionsprinted) {
2547: PetscViewerASCIIPrintf(viewer," Run with -snes_test_jacobian_view and optionally -snes_test_jacobian <threshold> to show difference\n");
2548: PetscViewerASCIIPrintf(viewer," of hand-coded and finite difference Jacobian entries greater than <threshold>.\n");
2549: }
2550: if (!directionsprinted) {
2551: PetscViewerASCIIPrintf(viewer," Testing hand-coded Jacobian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");
2552: PetscViewerASCIIPrintf(viewer," O(1.e-8), the hand-coded Jacobian is probably correct.\n");
2553: directionsprinted = PETSC_TRUE;
2554: }
2555: if (complete_print) {
2556: PetscViewerPushFormat(mviewer,format);
2557: }
2559: PetscObjectTypeCompare((PetscObject)snes->jacobian,MATMFFD,&flg);
2560: if (!flg) jacobian = snes->jacobian;
2561: else jacobian = snes->jacobian_pre;
2563: if (!x) {
2564: MatCreateVecs(jacobian, &x, NULL);
2565: } else {
2566: PetscObjectReference((PetscObject) x);
2567: }
2568: if (!f) {
2569: VecDuplicate(x, &f);
2570: } else {
2571: PetscObjectReference((PetscObject) f);
2572: }
2573: /* evaluate the function at this point because SNESComputeJacobianDefault() assumes that the function has been evaluated and put into snes->vec_func */
2574: SNESComputeFunction(snes,x,f);
2575: VecDestroy(&f);
2576: PetscObjectTypeCompare((PetscObject)snes,SNESKSPTRANSPOSEONLY,&istranspose);
2577: while (jacobian) {
2578: Mat JT = NULL, Jsave = NULL;
2580: if (istranspose) {
2581: MatCreateTranspose(jacobian,&JT);
2582: Jsave = jacobian;
2583: jacobian = JT;
2584: }
2585: PetscObjectBaseTypeCompareAny((PetscObject)jacobian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPISBAIJ,"");
2586: if (flg) {
2587: A = jacobian;
2588: PetscObjectReference((PetscObject)A);
2589: } else {
2590: MatComputeOperator(jacobian,MATAIJ,&A);
2591: }
2593: MatGetType(A,&mattype);
2594: MatGetSize(A,&M,&N);
2595: MatGetLocalSize(A,&m,&n);
2596: MatCreate(PetscObjectComm((PetscObject)A),&B);
2597: MatSetType(B,mattype);
2598: MatSetSizes(B,m,n,M,N);
2599: MatSetBlockSizesFromMats(B,A,A);
2600: MatSetUp(B);
2601: MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
2603: SNESGetFunction(snes,NULL,NULL,&functx);
2604: SNESComputeJacobianDefault(snes,x,B,B,functx);
2606: MatDuplicate(B,MAT_COPY_VALUES,&D);
2607: MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2608: MatNorm(D,NORM_FROBENIUS,&nrm);
2609: MatNorm(A,NORM_FROBENIUS,&gnorm);
2610: MatDestroy(&D);
2611: if (!gnorm) gnorm = 1; /* just in case */
2612: PetscViewerASCIIPrintf(viewer," ||J - Jfd||_F/||J||_F = %g, ||J - Jfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);
2614: if (complete_print) {
2615: PetscViewerASCIIPrintf(viewer," Hand-coded Jacobian ----------\n");
2616: MatView(A,mviewer);
2617: PetscViewerASCIIPrintf(viewer," Finite difference Jacobian ----------\n");
2618: MatView(B,mviewer);
2619: }
2621: if (threshold_print || complete_print) {
2622: PetscInt Istart, Iend, *ccols, bncols, cncols, j, row;
2623: PetscScalar *cvals;
2624: const PetscInt *bcols;
2625: const PetscScalar *bvals;
2627: MatCreate(PetscObjectComm((PetscObject)A),&C);
2628: MatSetType(C,mattype);
2629: MatSetSizes(C,m,n,M,N);
2630: MatSetBlockSizesFromMats(C,A,A);
2631: MatSetUp(C);
2632: MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
2634: MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2635: MatGetOwnershipRange(B,&Istart,&Iend);
2637: for (row = Istart; row < Iend; row++) {
2638: MatGetRow(B,row,&bncols,&bcols,&bvals);
2639: PetscMalloc2(bncols,&ccols,bncols,&cvals);
2640: for (j = 0, cncols = 0; j < bncols; j++) {
2641: if (PetscAbsScalar(bvals[j]) > threshold) {
2642: ccols[cncols] = bcols[j];
2643: cvals[cncols] = bvals[j];
2644: cncols += 1;
2645: }
2646: }
2647: if (cncols) {
2648: MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);
2649: }
2650: MatRestoreRow(B,row,&bncols,&bcols,&bvals);
2651: PetscFree2(ccols,cvals);
2652: }
2653: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2654: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2655: PetscViewerASCIIPrintf(viewer," Hand-coded minus finite-difference Jacobian with tolerance %g ----------\n",(double)threshold);
2656: MatView(C,complete_print ? mviewer : viewer);
2657: MatDestroy(&C);
2658: }
2659: MatDestroy(&A);
2660: MatDestroy(&B);
2661: MatDestroy(&JT);
2662: if (Jsave) jacobian = Jsave;
2663: if (jacobian != snes->jacobian_pre) {
2664: jacobian = snes->jacobian_pre;
2665: PetscViewerASCIIPrintf(viewer," ---------- Testing Jacobian for preconditioner -------------\n");
2666: }
2667: else jacobian = NULL;
2668: }
2669: VecDestroy(&x);
2670: if (complete_print) {
2671: PetscViewerPopFormat(mviewer);
2672: }
2673: if (mviewer) PetscViewerDestroy(&mviewer);
2674: PetscViewerASCIISetTab(viewer,tabs);
2675: return 0;
2676: }
2678: /*@
2679: SNESComputeJacobian - Computes the Jacobian matrix that has been set with SNESSetJacobian().
2681: Collective on SNES
2683: Input Parameters:
2684: + snes - the SNES context
2685: - x - input vector
2687: Output Parameters:
2688: + A - Jacobian matrix
2689: - B - optional preconditioning matrix
2691: Options Database Keys:
2692: + -snes_lag_preconditioner <lag> - how often to rebuild preconditioner
2693: . -snes_lag_jacobian <lag> - how often to rebuild Jacobian
2694: . -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.
2695: . -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
2696: . -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2697: . -snes_compare_explicit_draw - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2698: . -snes_compare_explicit_contour - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2699: . -snes_compare_operator - Make the comparison options above use the operator instead of the preconditioning matrix
2700: . -snes_compare_coloring - Compute the finite difference Jacobian using coloring and display norms of difference
2701: . -snes_compare_coloring_display - Compute the finite difference Jacobian using coloring and display verbose differences
2702: . -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2703: . -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2704: . -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2705: . -snes_compare_coloring_draw - Compute the finite difference Jacobian using coloring and draw differences
2706: - -snes_compare_coloring_draw_contour - Compute the finite difference Jacobian using coloring and show contours of matrices and differences
2708: Notes:
2709: Most users should not need to explicitly call this routine, as it
2710: is used internally within the nonlinear solvers.
2712: Developer Notes:
2713: 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
2714: for with the SNESType of test that has been removed.
2716: Level: developer
2718: .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2719: @*/
2720: PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2721: {
2722: PetscBool flag;
2723: DM dm;
2724: DMSNES sdm;
2725: KSP ksp;
2730: VecValidValues(X,2,PETSC_TRUE);
2731: SNESGetDM(snes,&dm);
2732: DMGetDMSNES(dm,&sdm);
2736: /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
2738: if (snes->lagjacobian == -2) {
2739: snes->lagjacobian = -1;
2741: PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2742: } else if (snes->lagjacobian == -1) {
2743: PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2744: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2745: if (flag) {
2746: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2747: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2748: }
2749: return 0;
2750: } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2751: PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2752: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2753: if (flag) {
2754: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2755: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2756: }
2757: return 0;
2758: }
2759: if (snes->npc && snes->npcside == PC_LEFT) {
2760: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2761: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2762: return 0;
2763: }
2765: PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2766: VecLockReadPush(X);
2767: PetscStackPush("SNES user Jacobian function");
2768: (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2769: PetscStackPop;
2770: VecLockReadPop(X);
2771: PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);
2773: /* attach latest linearization point to the preconditioning matrix */
2774: PetscObjectCompose((PetscObject)B,"__SNES_latest_X",(PetscObject)X);
2776: /* the next line ensures that snes->ksp exists */
2777: SNESGetKSP(snes,&ksp);
2778: if (snes->lagpreconditioner == -2) {
2779: PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2780: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2781: snes->lagpreconditioner = -1;
2782: } else if (snes->lagpreconditioner == -1) {
2783: PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2784: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2785: } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2786: PetscInfo(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2787: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2788: } else {
2789: PetscInfo(snes,"Rebuilding preconditioner\n");
2790: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2791: }
2793: SNESTestJacobian(snes);
2794: /* make sure user returned a correct Jacobian and preconditioner */
2797: {
2798: PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
2799: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit",NULL,NULL,&flag);
2800: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",NULL,NULL,&flag_draw);
2801: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",NULL,NULL,&flag_contour);
2802: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_compare_operator",NULL,NULL,&flag_operator);
2803: if (flag || flag_draw || flag_contour) {
2804: Mat Bexp_mine = NULL,Bexp,FDexp;
2805: PetscViewer vdraw,vstdout;
2806: PetscBool flg;
2807: if (flag_operator) {
2808: MatComputeOperator(A,MATAIJ,&Bexp_mine);
2809: Bexp = Bexp_mine;
2810: } else {
2811: /* See if the preconditioning matrix can be viewed and added directly */
2812: PetscObjectBaseTypeCompareAny((PetscObject)B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2813: if (flg) Bexp = B;
2814: else {
2815: /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
2816: MatComputeOperator(B,MATAIJ,&Bexp_mine);
2817: Bexp = Bexp_mine;
2818: }
2819: }
2820: MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);
2821: SNESComputeJacobianDefault(snes,X,FDexp,FDexp,NULL);
2822: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2823: if (flag_draw || flag_contour) {
2824: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2825: if (flag_contour) PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2826: } else vdraw = NULL;
2827: PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator ? "Jacobian" : "preconditioning Jacobian");
2828: if (flag) MatView(Bexp,vstdout);
2829: if (vdraw) MatView(Bexp,vdraw);
2830: PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");
2831: if (flag) MatView(FDexp,vstdout);
2832: if (vdraw) MatView(FDexp,vdraw);
2833: MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);
2834: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");
2835: if (flag) MatView(FDexp,vstdout);
2836: if (vdraw) { /* Always use contour for the difference */
2837: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2838: MatView(FDexp,vdraw);
2839: PetscViewerPopFormat(vdraw);
2840: }
2841: if (flag_contour) PetscViewerPopFormat(vdraw);
2842: PetscViewerDestroy(&vdraw);
2843: MatDestroy(&Bexp_mine);
2844: MatDestroy(&FDexp);
2845: }
2846: }
2847: {
2848: PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
2849: PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
2850: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring",NULL,NULL,&flag);
2851: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_display",NULL,NULL,&flag_display);
2852: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",NULL,NULL,&flag_draw);
2853: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",NULL,NULL,&flag_contour);
2854: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",NULL,NULL,&flag_threshold);
2855: if (flag_threshold) {
2856: PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,NULL);
2857: PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,NULL);
2858: }
2859: if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
2860: Mat Bfd;
2861: PetscViewer vdraw,vstdout;
2862: MatColoring coloring;
2863: ISColoring iscoloring;
2864: MatFDColoring matfdcoloring;
2865: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2866: void *funcctx;
2867: PetscReal norm1,norm2,normmax;
2869: MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2870: MatColoringCreate(Bfd,&coloring);
2871: MatColoringSetType(coloring,MATCOLORINGSL);
2872: MatColoringSetFromOptions(coloring);
2873: MatColoringApply(coloring,&iscoloring);
2874: MatColoringDestroy(&coloring);
2875: MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2876: MatFDColoringSetFromOptions(matfdcoloring);
2877: MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2878: ISColoringDestroy(&iscoloring);
2880: /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2881: SNESGetFunction(snes,NULL,&func,&funcctx);
2882: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2883: PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2884: PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2885: MatFDColoringSetFromOptions(matfdcoloring);
2886: MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2887: MatFDColoringDestroy(&matfdcoloring);
2889: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2890: if (flag_draw || flag_contour) {
2891: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2892: if (flag_contour) PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2893: } else vdraw = NULL;
2894: PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2895: if (flag_display) MatView(B,vstdout);
2896: if (vdraw) MatView(B,vdraw);
2897: PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2898: if (flag_display) MatView(Bfd,vstdout);
2899: if (vdraw) MatView(Bfd,vdraw);
2900: MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2901: MatNorm(Bfd,NORM_1,&norm1);
2902: MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2903: MatNorm(Bfd,NORM_MAX,&normmax);
2904: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2905: if (flag_display) MatView(Bfd,vstdout);
2906: if (vdraw) { /* Always use contour for the difference */
2907: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2908: MatView(Bfd,vdraw);
2909: PetscViewerPopFormat(vdraw);
2910: }
2911: if (flag_contour) PetscViewerPopFormat(vdraw);
2913: if (flag_threshold) {
2914: PetscInt bs,rstart,rend,i;
2915: MatGetBlockSize(B,&bs);
2916: MatGetOwnershipRange(B,&rstart,&rend);
2917: for (i=rstart; i<rend; i++) {
2918: const PetscScalar *ba,*ca;
2919: const PetscInt *bj,*cj;
2920: PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
2921: PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0;
2922: MatGetRow(B,i,&bn,&bj,&ba);
2923: MatGetRow(Bfd,i,&cn,&cj,&ca);
2925: for (j=0; j<bn; j++) {
2926: PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2927: if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
2928: maxentrycol = bj[j];
2929: maxentry = PetscRealPart(ba[j]);
2930: }
2931: if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
2932: maxdiffcol = bj[j];
2933: maxdiff = PetscRealPart(ca[j]);
2934: }
2935: if (rdiff > maxrdiff) {
2936: maxrdiffcol = bj[j];
2937: maxrdiff = rdiff;
2938: }
2939: }
2940: if (maxrdiff > 1) {
2941: 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);
2942: for (j=0; j<bn; j++) {
2943: PetscReal rdiff;
2944: rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2945: if (rdiff > 1) {
2946: PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
2947: }
2948: }
2949: PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
2950: }
2951: MatRestoreRow(B,i,&bn,&bj,&ba);
2952: MatRestoreRow(Bfd,i,&cn,&cj,&ca);
2953: }
2954: }
2955: PetscViewerDestroy(&vdraw);
2956: MatDestroy(&Bfd);
2957: }
2958: }
2959: return 0;
2960: }
2962: /*MC
2963: SNESJacobianFunction - Function used to convey the nonlinear Jacobian of the function to be solved by SNES
2965: Synopsis:
2966: #include "petscsnes.h"
2967: PetscErrorCode SNESJacobianFunction(SNES snes,Vec x,Mat Amat,Mat Pmat,void *ctx);
2969: Collective on snes
2971: Input Parameters:
2972: + x - input vector, the Jacobian is to be computed at this value
2973: - ctx - [optional] user-defined Jacobian context
2975: Output Parameters:
2976: + Amat - the matrix that defines the (approximate) Jacobian
2977: - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2979: Level: intermediate
2981: .seealso: SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2982: M*/
2984: /*@C
2985: SNESSetJacobian - Sets the function to compute Jacobian as well as the
2986: location to store the matrix.
2988: Logically Collective on SNES
2990: Input Parameters:
2991: + snes - the SNES context
2992: . Amat - the matrix that defines the (approximate) Jacobian
2993: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2994: . J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
2995: - ctx - [optional] user-defined context for private data for the
2996: Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)
2998: Notes:
2999: If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
3000: each matrix.
3002: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
3003: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
3005: If using SNESComputeJacobianDefaultColor() to assemble a Jacobian, the ctx argument
3006: must be a MatFDColoring.
3008: Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian. One common
3009: example is to use the "Picard linearization" which only differentiates through the highest order parts of each term.
3011: Level: beginner
3013: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J,
3014: SNESSetPicard(), SNESJacobianFunction
3015: @*/
3016: PetscErrorCode SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
3017: {
3018: DM dm;
3025: SNESGetDM(snes,&dm);
3026: DMSNESSetJacobian(dm,J,ctx);
3027: if (Amat) {
3028: PetscObjectReference((PetscObject)Amat);
3029: MatDestroy(&snes->jacobian);
3031: snes->jacobian = Amat;
3032: }
3033: if (Pmat) {
3034: PetscObjectReference((PetscObject)Pmat);
3035: MatDestroy(&snes->jacobian_pre);
3037: snes->jacobian_pre = Pmat;
3038: }
3039: return 0;
3040: }
3042: /*@C
3043: SNESGetJacobian - Returns the Jacobian matrix and optionally the user
3044: provided context for evaluating the Jacobian.
3046: Not Collective, but Mat object will be parallel if SNES object is
3048: Input Parameter:
3049: . snes - the nonlinear solver context
3051: Output Parameters:
3052: + Amat - location to stash (approximate) Jacobian matrix (or NULL)
3053: . Pmat - location to stash matrix used to compute the preconditioner (or NULL)
3054: . J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
3055: - ctx - location to stash Jacobian ctx (or NULL)
3057: Level: advanced
3059: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
3060: @*/
3061: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
3062: {
3063: DM dm;
3064: DMSNES sdm;
3067: if (Amat) *Amat = snes->jacobian;
3068: if (Pmat) *Pmat = snes->jacobian_pre;
3069: SNESGetDM(snes,&dm);
3070: DMGetDMSNES(dm,&sdm);
3071: if (J) *J = sdm->ops->computejacobian;
3072: if (ctx) *ctx = sdm->jacobianctx;
3073: return 0;
3074: }
3076: static PetscErrorCode SNESSetDefaultComputeJacobian(SNES snes)
3077: {
3078: DM dm;
3079: DMSNES sdm;
3081: SNESGetDM(snes,&dm);
3082: DMGetDMSNES(dm,&sdm);
3083: if (!sdm->ops->computejacobian && snes->jacobian_pre) {
3084: DM dm;
3085: PetscBool isdense,ismf;
3087: SNESGetDM(snes,&dm);
3088: PetscObjectTypeCompareAny((PetscObject)snes->jacobian_pre,&isdense,MATSEQDENSE,MATMPIDENSE,MATDENSE,NULL);
3089: PetscObjectTypeCompareAny((PetscObject)snes->jacobian_pre,&ismf,MATMFFD,MATSHELL,NULL);
3090: if (isdense) {
3091: DMSNESSetJacobian(dm,SNESComputeJacobianDefault,NULL);
3092: } else if (!ismf) {
3093: DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
3094: }
3095: }
3096: return 0;
3097: }
3099: /*@
3100: SNESSetUp - Sets up the internal data structures for the later use
3101: of a nonlinear solver.
3103: Collective on SNES
3105: Input Parameters:
3106: . snes - the SNES context
3108: Notes:
3109: For basic use of the SNES solvers the user need not explicitly call
3110: SNESSetUp(), since these actions will automatically occur during
3111: the call to SNESSolve(). However, if one wishes to control this
3112: phase separately, SNESSetUp() should be called after SNESCreate()
3113: and optional routines of the form SNESSetXXX(), but before SNESSolve().
3115: Level: advanced
3117: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
3118: @*/
3119: PetscErrorCode SNESSetUp(SNES snes)
3120: {
3121: DM dm;
3122: DMSNES sdm;
3123: SNESLineSearch linesearch, pclinesearch;
3124: void *lsprectx,*lspostctx;
3125: PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
3126: PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
3127: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
3128: Vec f,fpc;
3129: void *funcctx;
3130: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
3131: void *jacctx,*appctx;
3132: Mat j,jpre;
3135: if (snes->setupcalled) return 0;
3136: PetscLogEventBegin(SNES_Setup,snes,0,0,0);
3138: if (!((PetscObject)snes)->type_name) {
3139: SNESSetType(snes,SNESNEWTONLS);
3140: }
3142: SNESGetFunction(snes,&snes->vec_func,NULL,NULL);
3144: SNESGetDM(snes,&dm);
3145: DMGetDMSNES(dm,&sdm);
3147: SNESSetDefaultComputeJacobian(snes);
3149: if (!snes->vec_func) {
3150: DMCreateGlobalVector(dm,&snes->vec_func);
3151: }
3153: if (!snes->ksp) {
3154: SNESGetKSP(snes, &snes->ksp);
3155: }
3157: if (snes->linesearch) {
3158: SNESGetLineSearch(snes, &snes->linesearch);
3159: SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);
3160: }
3162: if (snes->npc && snes->npcside == PC_LEFT) {
3163: snes->mf = PETSC_TRUE;
3164: snes->mf_operator = PETSC_FALSE;
3165: }
3167: if (snes->npc) {
3168: /* copy the DM over */
3169: SNESGetDM(snes,&dm);
3170: SNESSetDM(snes->npc,dm);
3172: SNESGetFunction(snes,&f,&func,&funcctx);
3173: VecDuplicate(f,&fpc);
3174: SNESSetFunction(snes->npc,fpc,func,funcctx);
3175: SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
3176: SNESSetJacobian(snes->npc,j,jpre,jac,jacctx);
3177: SNESGetApplicationContext(snes,&appctx);
3178: SNESSetApplicationContext(snes->npc,appctx);
3179: VecDestroy(&fpc);
3181: /* copy the function pointers over */
3182: PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)snes->npc);
3184: /* default to 1 iteration */
3185: SNESSetTolerances(snes->npc,0.0,0.0,0.0,1,snes->npc->max_funcs);
3186: if (snes->npcside == PC_RIGHT) {
3187: SNESSetNormSchedule(snes->npc,SNES_NORM_FINAL_ONLY);
3188: } else {
3189: SNESSetNormSchedule(snes->npc,SNES_NORM_NONE);
3190: }
3191: SNESSetFromOptions(snes->npc);
3193: /* copy the line search context over */
3194: if (snes->linesearch && snes->npc->linesearch) {
3195: SNESGetLineSearch(snes,&linesearch);
3196: SNESGetLineSearch(snes->npc,&pclinesearch);
3197: SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
3198: SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
3199: SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
3200: SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
3201: PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
3202: }
3203: }
3204: if (snes->mf) {
3205: SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
3206: }
3207: if (snes->ops->usercompute && !snes->user) {
3208: (*snes->ops->usercompute)(snes,(void**)&snes->user);
3209: }
3211: snes->jac_iter = 0;
3212: snes->pre_iter = 0;
3214: if (snes->ops->setup) {
3215: (*snes->ops->setup)(snes);
3216: }
3218: SNESSetDefaultComputeJacobian(snes);
3220: if (snes->npc && snes->npcside == PC_LEFT) {
3221: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
3222: if (snes->linesearch) {
3223: SNESGetLineSearch(snes,&linesearch);
3224: SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
3225: }
3226: }
3227: }
3228: PetscLogEventEnd(SNES_Setup,snes,0,0,0);
3229: snes->setupcalled = PETSC_TRUE;
3230: return 0;
3231: }
3233: /*@
3234: SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats
3236: Collective on SNES
3238: Input Parameter:
3239: . snes - iterative context obtained from SNESCreate()
3241: Level: intermediate
3243: Notes:
3244: Also calls the application context destroy routine set with SNESSetComputeApplicationContext()
3246: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
3247: @*/
3248: PetscErrorCode SNESReset(SNES snes)
3249: {
3251: if (snes->ops->userdestroy && snes->user) {
3252: (*snes->ops->userdestroy)((void**)&snes->user);
3253: snes->user = NULL;
3254: }
3255: if (snes->npc) {
3256: SNESReset(snes->npc);
3257: }
3259: if (snes->ops->reset) {
3260: (*snes->ops->reset)(snes);
3261: }
3262: if (snes->ksp) {
3263: KSPReset(snes->ksp);
3264: }
3266: if (snes->linesearch) {
3267: SNESLineSearchReset(snes->linesearch);
3268: }
3270: VecDestroy(&snes->vec_rhs);
3271: VecDestroy(&snes->vec_sol);
3272: VecDestroy(&snes->vec_sol_update);
3273: VecDestroy(&snes->vec_func);
3274: MatDestroy(&snes->jacobian);
3275: MatDestroy(&snes->jacobian_pre);
3276: MatDestroy(&snes->picard);
3277: VecDestroyVecs(snes->nwork,&snes->work);
3278: VecDestroyVecs(snes->nvwork,&snes->vwork);
3280: snes->alwayscomputesfinalresidual = PETSC_FALSE;
3282: snes->nwork = snes->nvwork = 0;
3283: snes->setupcalled = PETSC_FALSE;
3284: return 0;
3285: }
3287: /*@
3288: SNESConvergedReasonViewCancel - Clears all the reasonview functions for a SNES object.
3290: Collective on SNES
3292: Input Parameter:
3293: . snes - iterative context obtained from SNESCreate()
3295: Level: intermediate
3297: .seealso: SNESCreate(), SNESDestroy(), SNESReset()
3298: @*/
3299: PetscErrorCode SNESConvergedReasonViewCancel(SNES snes)
3300: {
3301: PetscInt i;
3304: for (i=0; i<snes->numberreasonviews; i++) {
3305: if (snes->reasonviewdestroy[i]) {
3306: (*snes->reasonviewdestroy[i])(&snes->reasonviewcontext[i]);
3307: }
3308: }
3309: snes->numberreasonviews = 0;
3310: return 0;
3311: }
3313: /*@C
3314: SNESDestroy - Destroys the nonlinear solver context that was created
3315: with SNESCreate().
3317: Collective on SNES
3319: Input Parameter:
3320: . snes - the SNES context
3322: Level: beginner
3324: .seealso: SNESCreate(), SNESSolve()
3325: @*/
3326: PetscErrorCode SNESDestroy(SNES *snes)
3327: {
3328: if (!*snes) return 0;
3330: if (--((PetscObject)(*snes))->refct > 0) {*snes = NULL; return 0;}
3332: SNESReset((*snes));
3333: SNESDestroy(&(*snes)->npc);
3335: /* if memory was published with SAWs then destroy it */
3336: PetscObjectSAWsViewOff((PetscObject)*snes);
3337: if ((*snes)->ops->destroy) (*((*snes))->ops->destroy)((*snes));
3339: if ((*snes)->dm) DMCoarsenHookRemove((*snes)->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,*snes);
3340: DMDestroy(&(*snes)->dm);
3341: KSPDestroy(&(*snes)->ksp);
3342: SNESLineSearchDestroy(&(*snes)->linesearch);
3344: PetscFree((*snes)->kspconvctx);
3345: if ((*snes)->ops->convergeddestroy) {
3346: (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
3347: }
3348: if ((*snes)->conv_hist_alloc) {
3349: PetscFree2((*snes)->conv_hist,(*snes)->conv_hist_its);
3350: }
3351: SNESMonitorCancel((*snes));
3352: SNESConvergedReasonViewCancel((*snes));
3353: PetscHeaderDestroy(snes);
3354: return 0;
3355: }
3357: /* ----------- Routines to set solver parameters ---------- */
3359: /*@
3360: SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
3362: Logically Collective on SNES
3364: Input Parameters:
3365: + snes - the SNES context
3366: - lag - 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3367: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
3369: Options Database Keys:
3370: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3371: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3372: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3373: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag
3375: Notes:
3376: The default is 1
3377: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 or SNESSetLagPreconditionerPersists() was called
3379: SNESSetLagPreconditionerPersists() allows using the same uniform lagging (for example every second solve) across multiple solves.
3381: Level: intermediate
3383: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESSetLagPreconditionerPersists(),
3384: SNESSetLagJacobianPersists()
3386: @*/
3387: PetscErrorCode SNESSetLagPreconditioner(SNES snes,PetscInt lag)
3388: {
3393: snes->lagpreconditioner = lag;
3394: return 0;
3395: }
3397: /*@
3398: SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does
3400: Logically Collective on SNES
3402: Input Parameters:
3403: + snes - the SNES context
3404: - steps - the number of refinements to do, defaults to 0
3406: Options Database Keys:
3407: . -snes_grid_sequence <steps> - Use grid sequencing to generate initial guess
3409: Level: intermediate
3411: Notes:
3412: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
3414: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetGridSequence()
3416: @*/
3417: PetscErrorCode SNESSetGridSequence(SNES snes,PetscInt steps)
3418: {
3421: snes->gridsequence = steps;
3422: return 0;
3423: }
3425: /*@
3426: SNESGetGridSequence - gets the number of steps of grid sequencing that SNES does
3428: Logically Collective on SNES
3430: Input Parameter:
3431: . snes - the SNES context
3433: Output Parameter:
3434: . steps - the number of refinements to do, defaults to 0
3436: Options Database Keys:
3437: . -snes_grid_sequence <steps> - set number of refinements
3439: Level: intermediate
3441: Notes:
3442: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
3444: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESSetGridSequence()
3446: @*/
3447: PetscErrorCode SNESGetGridSequence(SNES snes,PetscInt *steps)
3448: {
3450: *steps = snes->gridsequence;
3451: return 0;
3452: }
3454: /*@
3455: SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
3457: Not Collective
3459: Input Parameter:
3460: . snes - the SNES context
3462: Output Parameter:
3463: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3464: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
3466: Options Database Keys:
3467: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3468: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3469: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3470: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag
3472: Notes:
3473: The default is 1
3474: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3476: Level: intermediate
3478: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner(), SNESSetLagJacobianPersists(), SNESSetLagPreconditionerPersists()
3480: @*/
3481: PetscErrorCode SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
3482: {
3484: *lag = snes->lagpreconditioner;
3485: return 0;
3486: }
3488: /*@
3489: SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
3490: often the preconditioner is rebuilt.
3492: Logically Collective on SNES
3494: Input Parameters:
3495: + snes - the SNES context
3496: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3497: the Jacobian is built etc. -2 means rebuild at next chance but then never again
3499: Options Database Keys:
3500: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3501: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3502: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3503: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag.
3505: Notes:
3506: The default is 1
3507: The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3508: 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
3509: at the next Newton step but never again (unless it is reset to another value)
3511: Level: intermediate
3513: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobianPersists(), SNESSetLagPreconditionerPersists()
3515: @*/
3516: PetscErrorCode SNESSetLagJacobian(SNES snes,PetscInt lag)
3517: {
3522: snes->lagjacobian = lag;
3523: return 0;
3524: }
3526: /*@
3527: SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
3529: Not Collective
3531: Input Parameter:
3532: . snes - the SNES context
3534: Output Parameter:
3535: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3536: the Jacobian is built etc.
3538: Notes:
3539: The default is 1
3540: The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 or SNESSetLagJacobianPersists() was called.
3542: Level: intermediate
3544: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner(), SNESSetLagJacobianPersists(), SNESSetLagPreconditionerPersists()
3546: @*/
3547: PetscErrorCode SNESGetLagJacobian(SNES snes,PetscInt *lag)
3548: {
3550: *lag = snes->lagjacobian;
3551: return 0;
3552: }
3554: /*@
3555: SNESSetLagJacobianPersists - Set whether or not the Jacobian lagging persists through multiple solves
3557: Logically collective on SNES
3559: Input Parameters:
3560: + snes - the SNES context
3561: - flg - jacobian lagging persists if true
3563: Options Database Keys:
3564: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3565: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3566: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3567: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag
3569: Notes:
3570: This is useful both for nonlinear preconditioning, where it's appropriate to have the Jacobian be stale by
3571: several solves, and for implicit time-stepping, where Jacobian lagging in the inner nonlinear solve over several
3572: timesteps may present huge efficiency gains.
3574: Level: developer
3576: .seealso: SNESSetLagPreconditionerPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC(), SNESSetLagJacobianPersists()
3578: @*/
3579: PetscErrorCode SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
3580: {
3583: snes->lagjac_persist = flg;
3584: return 0;
3585: }
3587: /*@
3588: SNESSetLagPreconditionerPersists - Set whether or not the preconditioner lagging persists through multiple nonlinear solves
3590: Logically Collective on SNES
3592: Input Parameters:
3593: + snes - the SNES context
3594: - flg - preconditioner lagging persists if true
3596: Options Database Keys:
3597: + -snes_lag_jacobian_persists <true,false> - sets the persistence
3598: . -snes_lag_jacobian <-2,1,2,...> - sets the lag
3599: . -snes_lag_preconditioner_persists <true,false> - sets the persistence
3600: - -snes_lag_preconditioner <-2,1,2,...> - sets the lag
3602: Notes:
3603: This is useful both for nonlinear preconditioning, where it's appropriate to have the preconditioner be stale
3604: by several solves, and for implicit time-stepping, where preconditioner lagging in the inner nonlinear solve over
3605: several timesteps may present huge efficiency gains.
3607: Level: developer
3609: .seealso: SNESSetLagJacobianPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC(), SNESSetLagPreconditioner()
3611: @*/
3612: PetscErrorCode SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3613: {
3616: snes->lagpre_persist = flg;
3617: return 0;
3618: }
3620: /*@
3621: SNESSetForceIteration - force SNESSolve() to take at least one iteration regardless of the initial residual norm
3623: Logically Collective on SNES
3625: Input Parameters:
3626: + snes - the SNES context
3627: - force - PETSC_TRUE require at least one iteration
3629: Options Database Keys:
3630: . -snes_force_iteration <force> - Sets forcing an iteration
3632: Notes:
3633: This is used sometimes with TS to prevent TS from detecting a false steady state solution
3635: Level: intermediate
3637: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3638: @*/
3639: PetscErrorCode SNESSetForceIteration(SNES snes,PetscBool force)
3640: {
3642: snes->forceiteration = force;
3643: return 0;
3644: }
3646: /*@
3647: SNESGetForceIteration - Whether or not to force SNESSolve() take at least one iteration regardless of the initial residual norm
3649: Logically Collective on SNES
3651: Input Parameters:
3652: . snes - the SNES context
3654: Output Parameter:
3655: . force - PETSC_TRUE requires at least one iteration.
3657: Level: intermediate
3659: .seealso: SNESSetForceIteration(), SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3660: @*/
3661: PetscErrorCode SNESGetForceIteration(SNES snes,PetscBool *force)
3662: {
3664: *force = snes->forceiteration;
3665: return 0;
3666: }
3668: /*@
3669: SNESSetTolerances - Sets various parameters used in convergence tests.
3671: Logically Collective on SNES
3673: Input Parameters:
3674: + snes - the SNES context
3675: . abstol - absolute convergence tolerance
3676: . rtol - relative convergence tolerance
3677: . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x ||
3678: . maxit - maximum number of iterations
3679: - maxf - maximum number of function evaluations (-1 indicates no limit)
3681: Options Database Keys:
3682: + -snes_atol <abstol> - Sets abstol
3683: . -snes_rtol <rtol> - Sets rtol
3684: . -snes_stol <stol> - Sets stol
3685: . -snes_max_it <maxit> - Sets maxit
3686: - -snes_max_funcs <maxf> - Sets maxf
3688: Notes:
3689: The default maximum number of iterations is 50.
3690: The default maximum number of function evaluations is 1000.
3692: Level: intermediate
3694: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance(), SNESSetForceIteration()
3695: @*/
3696: PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3697: {
3705: if (abstol != PETSC_DEFAULT) {
3707: snes->abstol = abstol;
3708: }
3709: if (rtol != PETSC_DEFAULT) {
3711: snes->rtol = rtol;
3712: }
3713: if (stol != PETSC_DEFAULT) {
3715: snes->stol = stol;
3716: }
3717: if (maxit != PETSC_DEFAULT) {
3719: snes->max_its = maxit;
3720: }
3721: if (maxf != PETSC_DEFAULT) {
3723: snes->max_funcs = maxf;
3724: }
3725: snes->tolerancesset = PETSC_TRUE;
3726: return 0;
3727: }
3729: /*@
3730: SNESSetDivergenceTolerance - Sets the divergence tolerance used for the SNES divergence test.
3732: Logically Collective on SNES
3734: Input Parameters:
3735: + snes - the SNES context
3736: - divtol - the divergence tolerance. Use -1 to deactivate the test.
3738: Options Database Keys:
3739: . -snes_divergence_tolerance <divtol> - Sets divtol
3741: Notes:
3742: The default divergence tolerance is 1e4.
3744: Level: intermediate
3746: .seealso: SNESSetTolerances(), SNESGetDivergenceTolerance
3747: @*/
3748: PetscErrorCode SNESSetDivergenceTolerance(SNES snes,PetscReal divtol)
3749: {
3753: if (divtol != PETSC_DEFAULT) {
3754: snes->divtol = divtol;
3755: }
3756: else {
3757: snes->divtol = 1.0e4;
3758: }
3759: return 0;
3760: }
3762: /*@
3763: SNESGetTolerances - Gets various parameters used in convergence tests.
3765: Not Collective
3767: Input Parameters:
3768: + snes - the SNES context
3769: . atol - absolute convergence tolerance
3770: . rtol - relative convergence tolerance
3771: . stol - convergence tolerance in terms of the norm
3772: of the change in the solution between steps
3773: . maxit - maximum number of iterations
3774: - maxf - maximum number of function evaluations
3776: Notes:
3777: The user can specify NULL for any parameter that is not needed.
3779: Level: intermediate
3781: .seealso: SNESSetTolerances()
3782: @*/
3783: PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3784: {
3786: if (atol) *atol = snes->abstol;
3787: if (rtol) *rtol = snes->rtol;
3788: if (stol) *stol = snes->stol;
3789: if (maxit) *maxit = snes->max_its;
3790: if (maxf) *maxf = snes->max_funcs;
3791: return 0;
3792: }
3794: /*@
3795: SNESGetDivergenceTolerance - Gets divergence tolerance used in divergence test.
3797: Not Collective
3799: Input Parameters:
3800: + snes - the SNES context
3801: - divtol - divergence tolerance
3803: Level: intermediate
3805: .seealso: SNESSetDivergenceTolerance()
3806: @*/
3807: PetscErrorCode SNESGetDivergenceTolerance(SNES snes,PetscReal *divtol)
3808: {
3810: if (divtol) *divtol = snes->divtol;
3811: return 0;
3812: }
3814: /*@
3815: SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
3817: Logically Collective on SNES
3819: Input Parameters:
3820: + snes - the SNES context
3821: - tol - tolerance
3823: Options Database Key:
3824: . -snes_trtol <tol> - Sets tol
3826: Level: intermediate
3828: .seealso: SNESSetTolerances()
3829: @*/
3830: PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3831: {
3834: snes->deltatol = tol;
3835: return 0;
3836: }
3838: PETSC_INTERN PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3840: PetscErrorCode SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3841: {
3842: PetscDrawLG lg;
3843: PetscReal x,y,per;
3844: PetscViewer v = (PetscViewer)monctx;
3845: static PetscReal prev; /* should be in the context */
3846: PetscDraw draw;
3849: PetscViewerDrawGetDrawLG(v,0,&lg);
3850: if (!n) PetscDrawLGReset(lg);
3851: PetscDrawLGGetDraw(lg,&draw);
3852: PetscDrawSetTitle(draw,"Residual norm");
3853: x = (PetscReal)n;
3854: if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3855: else y = -15.0;
3856: PetscDrawLGAddPoint(lg,&x,&y);
3857: if (n < 20 || !(n % 5) || snes->reason) {
3858: PetscDrawLGDraw(lg);
3859: PetscDrawLGSave(lg);
3860: }
3862: PetscViewerDrawGetDrawLG(v,1,&lg);
3863: if (!n) PetscDrawLGReset(lg);
3864: PetscDrawLGGetDraw(lg,&draw);
3865: PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3866: SNESMonitorRange_Private(snes,n,&per);
3867: x = (PetscReal)n;
3868: y = 100.0*per;
3869: PetscDrawLGAddPoint(lg,&x,&y);
3870: if (n < 20 || !(n % 5) || snes->reason) {
3871: PetscDrawLGDraw(lg);
3872: PetscDrawLGSave(lg);
3873: }
3875: PetscViewerDrawGetDrawLG(v,2,&lg);
3876: if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3877: PetscDrawLGGetDraw(lg,&draw);
3878: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3879: x = (PetscReal)n;
3880: y = (prev - rnorm)/prev;
3881: PetscDrawLGAddPoint(lg,&x,&y);
3882: if (n < 20 || !(n % 5) || snes->reason) {
3883: PetscDrawLGDraw(lg);
3884: PetscDrawLGSave(lg);
3885: }
3887: PetscViewerDrawGetDrawLG(v,3,&lg);
3888: if (!n) PetscDrawLGReset(lg);
3889: PetscDrawLGGetDraw(lg,&draw);
3890: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3891: x = (PetscReal)n;
3892: y = (prev - rnorm)/(prev*per);
3893: if (n > 2) { /*skip initial crazy value */
3894: PetscDrawLGAddPoint(lg,&x,&y);
3895: }
3896: if (n < 20 || !(n % 5) || snes->reason) {
3897: PetscDrawLGDraw(lg);
3898: PetscDrawLGSave(lg);
3899: }
3900: prev = rnorm;
3901: return 0;
3902: }
3904: /*@
3905: SNESMonitor - runs the user provided monitor routines, if they exist
3907: Collective on SNES
3909: Input Parameters:
3910: + snes - nonlinear solver context obtained from SNESCreate()
3911: . iter - iteration number
3912: - rnorm - relative norm of the residual
3914: Notes:
3915: This routine is called by the SNES implementations.
3916: It does not typically need to be called by the user.
3918: Level: developer
3920: .seealso: SNESMonitorSet()
3921: @*/
3922: PetscErrorCode SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3923: {
3924: PetscInt i,n = snes->numbermonitors;
3926: VecLockReadPush(snes->vec_sol);
3927: for (i=0; i<n; i++) {
3928: (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3929: }
3930: VecLockReadPop(snes->vec_sol);
3931: return 0;
3932: }
3934: /* ------------ Routines to set performance monitoring options ----------- */
3936: /*MC
3937: SNESMonitorFunction - functional form passed to SNESMonitorSet() to monitor convergence of nonlinear solver
3939: Synopsis:
3940: #include <petscsnes.h>
3941: $ PetscErrorCode SNESMonitorFunction(SNES snes,PetscInt its, PetscReal norm,void *mctx)
3943: Collective on snes
3945: Input Parameters:
3946: + snes - the SNES context
3947: . its - iteration number
3948: . norm - 2-norm function value (may be estimated)
3949: - mctx - [optional] monitoring context
3951: Level: advanced
3953: .seealso: SNESMonitorSet(), SNESMonitorGet()
3954: M*/
3956: /*@C
3957: SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3958: iteration of the nonlinear solver to display the iteration's
3959: progress.
3961: Logically Collective on SNES
3963: Input Parameters:
3964: + snes - the SNES context
3965: . f - the monitor function, see SNESMonitorFunction for the calling sequence
3966: . mctx - [optional] user-defined context for private data for the
3967: monitor routine (use NULL if no context is desired)
3968: - monitordestroy - [optional] routine that frees monitor context
3969: (may be NULL)
3971: Options Database Keys:
3972: + -snes_monitor - sets SNESMonitorDefault()
3973: . -snes_monitor draw::draw_lg - sets line graph monitor,
3974: - -snes_monitor_cancel - cancels all monitors that have
3975: been hardwired into a code by
3976: calls to SNESMonitorSet(), but
3977: does not cancel those set via
3978: the options database.
3980: Notes:
3981: Several different monitoring routines may be set by calling
3982: SNESMonitorSet() multiple times; all will be called in the
3983: order in which they were set.
3985: Fortran Notes:
3986: Only a single monitor function can be set for each SNES object
3988: Level: intermediate
3990: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3991: @*/
3992: PetscErrorCode SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3993: {
3994: PetscInt i;
3995: PetscBool identical;
3998: for (i=0; i<snes->numbermonitors;i++) {
3999: PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))snes->monitor[i],snes->monitorcontext[i],snes->monitordestroy[i],&identical);
4000: if (identical) return 0;
4001: }
4003: snes->monitor[snes->numbermonitors] = f;
4004: snes->monitordestroy[snes->numbermonitors] = monitordestroy;
4005: snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
4006: return 0;
4007: }
4009: /*@
4010: SNESMonitorCancel - Clears all the monitor functions for a SNES object.
4012: Logically Collective on SNES
4014: Input Parameters:
4015: . snes - the SNES context
4017: Options Database Key:
4018: . -snes_monitor_cancel - cancels all monitors that have been hardwired
4019: into a code by calls to SNESMonitorSet(), but does not cancel those
4020: set via the options database
4022: Notes:
4023: There is no way to clear one specific monitor from a SNES object.
4025: Level: intermediate
4027: .seealso: SNESMonitorDefault(), SNESMonitorSet()
4028: @*/
4029: PetscErrorCode SNESMonitorCancel(SNES snes)
4030: {
4031: PetscInt i;
4034: for (i=0; i<snes->numbermonitors; i++) {
4035: if (snes->monitordestroy[i]) {
4036: (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
4037: }
4038: }
4039: snes->numbermonitors = 0;
4040: return 0;
4041: }
4043: /*MC
4044: SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver
4046: Synopsis:
4047: #include <petscsnes.h>
4048: $ PetscErrorCode SNESConvergenceTest(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
4050: Collective on snes
4052: Input Parameters:
4053: + snes - the SNES context
4054: . it - current iteration (0 is the first and is before any Newton step)
4055: . xnorm - 2-norm of current iterate
4056: . gnorm - 2-norm of current step
4057: . f - 2-norm of function
4058: - cctx - [optional] convergence context
4060: Output Parameter:
4061: . reason - reason for convergence/divergence, only needs to be set when convergence or divergence is detected
4063: Level: intermediate
4065: .seealso: SNESSetConvergenceTest(), SNESGetConvergenceTest()
4066: M*/
4068: /*@C
4069: SNESSetConvergenceTest - Sets the function that is to be used
4070: to test for convergence of the nonlinear iterative solution.
4072: Logically Collective on SNES
4074: Input Parameters:
4075: + snes - the SNES context
4076: . SNESConvergenceTestFunction - routine to test for convergence
4077: . cctx - [optional] context for private data for the convergence routine (may be NULL)
4078: - destroy - [optional] destructor for the context (may be NULL; PETSC_NULL_FUNCTION in Fortran)
4080: Level: advanced
4082: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
4083: @*/
4084: PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
4085: {
4087: if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
4088: if (snes->ops->convergeddestroy) {
4089: (*snes->ops->convergeddestroy)(snes->cnvP);
4090: }
4091: snes->ops->converged = SNESConvergenceTestFunction;
4092: snes->ops->convergeddestroy = destroy;
4093: snes->cnvP = cctx;
4094: return 0;
4095: }
4097: /*@
4098: SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
4100: Not Collective
4102: Input Parameter:
4103: . snes - the SNES context
4105: Output Parameter:
4106: . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
4107: manual pages for the individual convergence tests for complete lists
4109: Options Database:
4110: . -snes_converged_reason - prints the reason to standard out
4112: Level: intermediate
4114: Notes:
4115: Should only be called after the call the SNESSolve() is complete, if it is called earlier it returns the value SNES__CONVERGED_ITERATING.
4117: .seealso: SNESSetConvergenceTest(), SNESSetConvergedReason(), SNESConvergedReason
4118: @*/
4119: PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
4120: {
4123: *reason = snes->reason;
4124: return 0;
4125: }
4127: /*@C
4128: SNESGetConvergedReasonString - Return a human readable string for snes converged reason
4130: Not Collective
4132: Input Parameter:
4133: . snes - the SNES context
4135: Output Parameter:
4136: . strreason - a human readable string that describes SNES converged reason
4138: Level: beginner
4140: .seealso: SNESGetConvergedReason()
4141: @*/
4142: PetscErrorCode SNESGetConvergedReasonString(SNES snes, const char** strreason)
4143: {
4146: *strreason = SNESConvergedReasons[snes->reason];
4147: return 0;
4148: }
4150: /*@
4151: SNESSetConvergedReason - Sets the reason the SNES iteration was stopped.
4153: Not Collective
4155: Input Parameters:
4156: + snes - the SNES context
4157: - reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
4158: manual pages for the individual convergence tests for complete lists
4160: Level: intermediate
4162: .seealso: SNESGetConvergedReason(), SNESSetConvergenceTest(), SNESConvergedReason
4163: @*/
4164: PetscErrorCode SNESSetConvergedReason(SNES snes,SNESConvergedReason reason)
4165: {
4167: snes->reason = reason;
4168: return 0;
4169: }
4171: /*@
4172: SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
4174: Logically Collective on SNES
4176: Input Parameters:
4177: + snes - iterative context obtained from SNESCreate()
4178: . a - array to hold history, this array will contain the function norms computed at each step
4179: . its - integer array holds the number of linear iterations for each solve.
4180: . na - size of a and its
4181: - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
4182: else it continues storing new values for new nonlinear solves after the old ones
4184: Notes:
4185: If 'a' and 'its' are NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
4186: default array of length 10000 is allocated.
4188: This routine is useful, e.g., when running a code for purposes
4189: of accurate performance monitoring, when no I/O should be done
4190: during the section of code that is being timed.
4192: Level: intermediate
4194: .seealso: SNESGetConvergenceHistory()
4196: @*/
4197: PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
4198: {
4202: if (!a) {
4203: if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
4204: PetscCalloc2(na,&a,na,&its);
4205: snes->conv_hist_alloc = PETSC_TRUE;
4206: }
4207: snes->conv_hist = a;
4208: snes->conv_hist_its = its;
4209: snes->conv_hist_max = (size_t)na;
4210: snes->conv_hist_len = 0;
4211: snes->conv_hist_reset = reset;
4212: return 0;
4213: }
4215: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4216: #include <engine.h> /* MATLAB include file */
4217: #include <mex.h> /* MATLAB include file */
4219: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
4220: {
4221: mxArray *mat;
4222: PetscInt i;
4223: PetscReal *ar;
4225: mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
4226: ar = (PetscReal*) mxGetData(mat);
4227: for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
4228: return mat;
4229: }
4230: #endif
4232: /*@C
4233: SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
4235: Not Collective
4237: Input Parameter:
4238: . snes - iterative context obtained from SNESCreate()
4240: Output Parameters:
4241: + a - array to hold history
4242: . its - integer array holds the number of linear iterations (or
4243: negative if not converged) for each solve.
4244: - na - size of a and its
4246: Notes:
4247: The calling sequence for this routine in Fortran is
4248: $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
4250: This routine is useful, e.g., when running a code for purposes
4251: of accurate performance monitoring, when no I/O should be done
4252: during the section of code that is being timed.
4254: Level: intermediate
4256: .seealso: SNESSetConvergenceHistory()
4258: @*/
4259: PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
4260: {
4262: if (a) *a = snes->conv_hist;
4263: if (its) *its = snes->conv_hist_its;
4264: if (na) *na = (PetscInt) snes->conv_hist_len;
4265: return 0;
4266: }
4268: /*@C
4269: SNESSetUpdate - Sets the general-purpose update function called
4270: at the beginning of every iteration of the nonlinear solve. Specifically
4271: it is called just before the Jacobian is "evaluated".
4273: Logically Collective on SNES
4275: Input Parameters:
4276: + snes - The nonlinear solver context
4277: - func - The function
4279: Calling sequence of func:
4280: $ func (SNES snes, PetscInt step);
4282: . step - The current step of the iteration
4284: Level: advanced
4286: Note:
4287: 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()
4288: This is not used by most users.
4290: 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.
4292: .seealso SNESSetJacobian(), SNESSolve(), SNESLineSearchSetPreCheck(), SNESLineSearchSetPostCheck(), SNESNewtonTRSetPreCheck(), SNESNewtonTRSetPostCheck(),
4293: SNESMonitorSet(), SNESSetDivergenceTest()
4294: @*/
4295: PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
4296: {
4298: snes->ops->update = func;
4299: return 0;
4300: }
4302: /*
4303: SNESScaleStep_Private - Scales a step so that its length is less than the
4304: positive parameter delta.
4306: Input Parameters:
4307: + snes - the SNES context
4308: . y - approximate solution of linear system
4309: . fnorm - 2-norm of current function
4310: - delta - trust region size
4312: Output Parameters:
4313: + gpnorm - predicted function norm at the new point, assuming local
4314: linearization. The value is zero if the step lies within the trust
4315: region, and exceeds zero otherwise.
4316: - ynorm - 2-norm of the step
4318: Note:
4319: For non-trust region methods such as SNESNEWTONLS, the parameter delta
4320: is set to be the maximum allowable step size.
4322: */
4323: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
4324: {
4325: PetscReal nrm;
4326: PetscScalar cnorm;
4332: VecNorm(y,NORM_2,&nrm);
4333: if (nrm > *delta) {
4334: nrm = *delta/nrm;
4335: *gpnorm = (1.0 - nrm)*(*fnorm);
4336: cnorm = nrm;
4337: VecScale(y,cnorm);
4338: *ynorm = *delta;
4339: } else {
4340: *gpnorm = 0.0;
4341: *ynorm = nrm;
4342: }
4343: return 0;
4344: }
4346: /*@C
4347: SNESConvergedReasonView - Displays the reason a SNES solve converged or diverged to a viewer
4349: Collective on SNES
4351: Parameter:
4352: + snes - iterative context obtained from SNESCreate()
4353: - viewer - the viewer to display the reason
4355: Options Database Keys:
4356: + -snes_converged_reason - print reason for converged or diverged, also prints number of iterations
4357: - -snes_converged_reason ::failed - only print reason and number of iterations when diverged
4359: Notes:
4360: To change the format of the output call PetscViewerPushFormat(viewer,format) before this call. Use PETSC_VIEWER_DEFAULT for the default,
4361: use PETSC_VIEWER_FAILED to only display a reason if it fails.
4363: Level: beginner
4365: .seealso: SNESCreate(), SNESSetUp(), SNESDestroy(), SNESSetTolerances(), SNESConvergedDefault(), SNESGetConvergedReason(), SNESConvergedReasonViewFromOptions(),
4366: PetscViewerPushFormat(), PetscViewerPopFormat()
4368: @*/
4369: PetscErrorCode SNESConvergedReasonView(SNES snes,PetscViewer viewer)
4370: {
4371: PetscViewerFormat format;
4372: PetscBool isAscii;
4374: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
4375: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
4376: if (isAscii) {
4377: PetscViewerGetFormat(viewer, &format);
4378: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
4379: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
4380: DM dm;
4381: Vec u;
4382: PetscDS prob;
4383: PetscInt Nf, f;
4384: PetscErrorCode (**exactSol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
4385: void **exactCtx;
4386: PetscReal error;
4388: SNESGetDM(snes, &dm);
4389: SNESGetSolution(snes, &u);
4390: DMGetDS(dm, &prob);
4391: PetscDSGetNumFields(prob, &Nf);
4392: PetscMalloc2(Nf, &exactSol, Nf, &exactCtx);
4393: for (f = 0; f < Nf; ++f) PetscDSGetExactSolution(prob, f, &exactSol[f], &exactCtx[f]);
4394: DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error);
4395: PetscFree2(exactSol, exactCtx);
4396: if (error < 1.0e-11) PetscViewerASCIIPrintf(viewer, "L_2 Error: < 1.0e-11\n");
4397: else PetscViewerASCIIPrintf(viewer, "L_2 Error: %g\n", error);
4398: }
4399: if (snes->reason > 0 && format != PETSC_VIEWER_FAILED) {
4400: if (((PetscObject) snes)->prefix) {
4401: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4402: } else {
4403: PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4404: }
4405: } else if (snes->reason <= 0) {
4406: if (((PetscObject) snes)->prefix) {
4407: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4408: } else {
4409: PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4410: }
4411: }
4412: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
4413: }
4414: return 0;
4415: }
4417: /*@C
4418: SNESConvergedReasonViewSet - Sets an ADDITIONAL function that is to be used at the
4419: end of the nonlinear solver to display the conver reason of the nonlinear solver.
4421: Logically Collective on SNES
4423: Input Parameters:
4424: + snes - the SNES context
4425: . f - the snes converged reason view function
4426: . vctx - [optional] user-defined context for private data for the
4427: snes converged reason view routine (use NULL if no context is desired)
4428: - reasonviewdestroy - [optional] routine that frees reasonview context
4429: (may be NULL)
4431: Options Database Keys:
4432: + -snes_converged_reason - sets a default SNESConvergedReasonView()
4433: - -snes_converged_reason_view_cancel - cancels all converged reason viewers that have
4434: been hardwired into a code by
4435: calls to SNESConvergedReasonViewSet(), but
4436: does not cancel those set via
4437: the options database.
4439: Notes:
4440: Several different converged reason view routines may be set by calling
4441: SNESConvergedReasonViewSet() multiple times; all will be called in the
4442: order in which they were set.
4444: Level: intermediate
4446: .seealso: SNESConvergedReasonView(), SNESConvergedReasonViewCancel()
4447: @*/
4448: PetscErrorCode SNESConvergedReasonViewSet(SNES snes,PetscErrorCode (*f)(SNES,void*),void *vctx,PetscErrorCode (*reasonviewdestroy)(void**))
4449: {
4450: PetscInt i;
4451: PetscBool identical;
4454: for (i=0; i<snes->numberreasonviews;i++) {
4455: PetscMonitorCompare((PetscErrorCode (*)(void))f,vctx,reasonviewdestroy,(PetscErrorCode (*)(void))snes->reasonview[i],snes->reasonviewcontext[i],snes->reasonviewdestroy[i],&identical);
4456: if (identical) return 0;
4457: }
4459: snes->reasonview[snes->numberreasonviews] = f;
4460: snes->reasonviewdestroy[snes->numberreasonviews] = reasonviewdestroy;
4461: snes->reasonviewcontext[snes->numberreasonviews++] = (void*)vctx;
4462: return 0;
4463: }
4465: /*@
4466: SNESConvergedReasonViewFromOptions - Processes command line options to determine if/how a SNESReason is to be viewed.
4467: All the user-provided convergedReasonView routines will be involved as well, if they exist.
4469: Collective on SNES
4471: Input Parameters:
4472: . snes - the SNES object
4474: Level: intermediate
4476: .seealso: SNESCreate(), SNESSetUp(), SNESDestroy(), SNESSetTolerances(), SNESConvergedDefault(), SNESGetConvergedReason(), SNESConvergedReasonView()
4478: @*/
4479: PetscErrorCode SNESConvergedReasonViewFromOptions(SNES snes)
4480: {
4481: PetscViewer viewer;
4482: PetscBool flg;
4483: static PetscBool incall = PETSC_FALSE;
4484: PetscViewerFormat format;
4485: PetscInt i;
4487: if (incall) return 0;
4488: incall = PETSC_TRUE;
4490: /* All user-provided viewers are called first, if they exist. */
4491: for (i=0; i<snes->numberreasonviews; i++) {
4492: (*snes->reasonview[i])(snes,snes->reasonviewcontext[i]);
4493: }
4495: /* Call PETSc default routine if users ask for it */
4496: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
4497: if (flg) {
4498: PetscViewerPushFormat(viewer,format);
4499: SNESConvergedReasonView(snes,viewer);
4500: PetscViewerPopFormat(viewer);
4501: PetscViewerDestroy(&viewer);
4502: }
4503: incall = PETSC_FALSE;
4504: return 0;
4505: }
4507: /*@
4508: SNESSolve - Solves a nonlinear system F(x) = b.
4509: Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
4511: Collective on SNES
4513: Input Parameters:
4514: + snes - the SNES context
4515: . b - the constant part of the equation F(x) = b, or NULL to use zero.
4516: - x - the solution vector.
4518: Notes:
4519: The user should initialize the vector,x, with the initial guess
4520: for the nonlinear solve prior to calling SNESSolve(). In particular,
4521: to employ an initial guess of zero, the user should explicitly set
4522: this vector to zero by calling VecSet().
4524: Level: beginner
4526: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution(),
4527: SNESNewtonTRSetPreCheck(), SNESNewtonTRGetPreCheck(), SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck(),
4528: SNESLineSearchSetPostCheck(), SNESLineSearchGetPostCheck(), SNESLineSearchSetPreCheck(), SNESLineSearchGetPreCheck()
4529: @*/
4530: PetscErrorCode SNESSolve(SNES snes,Vec b,Vec x)
4531: {
4532: PetscBool flg;
4533: PetscInt grid;
4534: Vec xcreated = NULL;
4535: DM dm;
4543: /* High level operations using the nonlinear solver */
4544: {
4545: PetscViewer viewer;
4546: PetscViewerFormat format;
4547: PetscInt num;
4548: PetscBool flg;
4549: static PetscBool incall = PETSC_FALSE;
4551: if (!incall) {
4552: /* Estimate the convergence rate of the discretization */
4553: PetscOptionsGetViewer(PetscObjectComm((PetscObject) snes),((PetscObject)snes)->options, ((PetscObject) snes)->prefix, "-snes_convergence_estimate", &viewer, &format, &flg);
4554: if (flg) {
4555: PetscConvEst conv;
4556: DM dm;
4557: PetscReal *alpha; /* Convergence rate of the solution error for each field in the L_2 norm */
4558: PetscInt Nf;
4560: incall = PETSC_TRUE;
4561: SNESGetDM(snes, &dm);
4562: DMGetNumFields(dm, &Nf);
4563: PetscCalloc1(Nf, &alpha);
4564: PetscConvEstCreate(PetscObjectComm((PetscObject) snes), &conv);
4565: PetscConvEstSetSolver(conv, (PetscObject) snes);
4566: PetscConvEstSetFromOptions(conv);
4567: PetscConvEstSetUp(conv);
4568: PetscConvEstGetConvRate(conv, alpha);
4569: PetscViewerPushFormat(viewer, format);
4570: PetscConvEstRateView(conv, alpha, viewer);
4571: PetscViewerPopFormat(viewer);
4572: PetscViewerDestroy(&viewer);
4573: PetscConvEstDestroy(&conv);
4574: PetscFree(alpha);
4575: incall = PETSC_FALSE;
4576: }
4577: /* Adaptively refine the initial grid */
4578: num = 1;
4579: PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_initial", &num, &flg);
4580: if (flg) {
4581: DMAdaptor adaptor;
4583: incall = PETSC_TRUE;
4584: DMAdaptorCreate(PetscObjectComm((PetscObject)snes), &adaptor);
4585: DMAdaptorSetSolver(adaptor, snes);
4586: DMAdaptorSetSequenceLength(adaptor, num);
4587: DMAdaptorSetFromOptions(adaptor);
4588: DMAdaptorSetUp(adaptor);
4589: DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_INITIAL, &dm, &x);
4590: DMAdaptorDestroy(&adaptor);
4591: incall = PETSC_FALSE;
4592: }
4593: /* Use grid sequencing to adapt */
4594: num = 0;
4595: PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_sequence", &num, NULL);
4596: if (num) {
4597: DMAdaptor adaptor;
4599: incall = PETSC_TRUE;
4600: DMAdaptorCreate(PetscObjectComm((PetscObject)snes), &adaptor);
4601: DMAdaptorSetSolver(adaptor, snes);
4602: DMAdaptorSetSequenceLength(adaptor, num);
4603: DMAdaptorSetFromOptions(adaptor);
4604: DMAdaptorSetUp(adaptor);
4605: DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_SEQUENTIAL, &dm, &x);
4606: DMAdaptorDestroy(&adaptor);
4607: incall = PETSC_FALSE;
4608: }
4609: }
4610: }
4611: if (!x) { x = snes->vec_sol; }
4612: if (!x) {
4613: SNESGetDM(snes,&dm);
4614: DMCreateGlobalVector(dm,&xcreated);
4615: x = xcreated;
4616: }
4617: SNESViewFromOptions(snes,NULL,"-snes_view_pre");
4619: for (grid=0; grid<snes->gridsequence; grid++) PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
4620: for (grid=0; grid<snes->gridsequence+1; grid++) {
4622: /* set solution vector */
4623: if (!grid) PetscObjectReference((PetscObject)x);
4624: VecDestroy(&snes->vec_sol);
4625: snes->vec_sol = x;
4626: SNESGetDM(snes,&dm);
4628: /* set affine vector if provided */
4629: if (b) PetscObjectReference((PetscObject)b);
4630: VecDestroy(&snes->vec_rhs);
4631: snes->vec_rhs = b;
4636: if (!snes->vec_sol_update /* && snes->vec_sol */) {
4637: VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
4638: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
4639: }
4640: DMShellSetGlobalVector(dm,snes->vec_sol);
4641: SNESSetUp(snes);
4643: if (!grid) {
4644: if (snes->ops->computeinitialguess) {
4645: (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
4646: }
4647: }
4649: if (snes->conv_hist_reset) snes->conv_hist_len = 0;
4650: if (snes->counters_reset) {snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;}
4652: PetscLogEventBegin(SNES_Solve,snes,0,0,0);
4653: (*snes->ops->solve)(snes);
4654: PetscLogEventEnd(SNES_Solve,snes,0,0,0);
4656: snes->domainerror = PETSC_FALSE; /* clear the flag if it has been set */
4658: if (snes->lagjac_persist) snes->jac_iter += snes->iter;
4659: if (snes->lagpre_persist) snes->pre_iter += snes->iter;
4661: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_local_min",NULL,NULL,&flg);
4662: if (flg && !PetscPreLoadingOn) SNESTestLocalMin(snes);
4663: /* Call converged reason views. This may involve user-provided viewers as well */
4664: SNESConvergedReasonViewFromOptions(snes);
4667: if (snes->reason < 0) break;
4668: if (grid < snes->gridsequence) {
4669: DM fine;
4670: Vec xnew;
4671: Mat interp;
4673: DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
4675: DMCreateInterpolation(snes->dm,fine,&interp,NULL);
4676: DMCreateGlobalVector(fine,&xnew);
4677: MatInterpolate(interp,x,xnew);
4678: DMInterpolate(snes->dm,interp,fine);
4679: MatDestroy(&interp);
4680: x = xnew;
4682: SNESReset(snes);
4683: SNESSetDM(snes,fine);
4684: SNESResetFromOptions(snes);
4685: DMDestroy(&fine);
4686: PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
4687: }
4688: }
4689: SNESViewFromOptions(snes,NULL,"-snes_view");
4690: VecViewFromOptions(snes->vec_sol,(PetscObject)snes,"-snes_view_solution");
4691: DMMonitor(snes->dm);
4692: SNESMonitorPauseFinal_Internal(snes);
4694: VecDestroy(&xcreated);
4695: PetscObjectSAWsBlock((PetscObject)snes);
4696: return 0;
4697: }
4699: /* --------- Internal routines for SNES Package --------- */
4701: /*@C
4702: SNESSetType - Sets the method for the nonlinear solver.
4704: Collective on SNES
4706: Input Parameters:
4707: + snes - the SNES context
4708: - type - a known method
4710: Options Database Key:
4711: . -snes_type <type> - Sets the method; use -help for a list
4712: of available methods (for instance, newtonls or newtontr)
4714: Notes:
4715: See "petsc/include/petscsnes.h" for available methods (for instance)
4716: + SNESNEWTONLS - Newton's method with line search
4717: (systems of nonlinear equations)
4718: - SNESNEWTONTR - Newton's method with trust region
4719: (systems of nonlinear equations)
4721: Normally, it is best to use the SNESSetFromOptions() command and then
4722: set the SNES solver type from the options database rather than by using
4723: this routine. Using the options database provides the user with
4724: maximum flexibility in evaluating the many nonlinear solvers.
4725: The SNESSetType() routine is provided for those situations where it
4726: is necessary to set the nonlinear solver independently of the command
4727: line or options database. This might be the case, for example, when
4728: the choice of solver changes during the execution of the program,
4729: and the user's application is taking responsibility for choosing the
4730: appropriate method.
4732: Developer Notes:
4733: SNESRegister() adds a constructor for a new SNESType to SNESList, SNESSetType() locates
4734: the constructor in that list and calls it to create the spexific object.
4736: Level: intermediate
4738: .seealso: SNESType, SNESCreate(), SNESDestroy(), SNESGetType(), SNESSetFromOptions()
4740: @*/
4741: PetscErrorCode SNESSetType(SNES snes,SNESType type)
4742: {
4743: PetscBool match;
4744: PetscErrorCode (*r)(SNES);
4749: PetscObjectTypeCompare((PetscObject)snes,type,&match);
4750: if (match) return 0;
4752: PetscFunctionListFind(SNESList,type,&r);
4754: /* Destroy the previous private SNES context */
4755: if (snes->ops->destroy) {
4756: (*(snes)->ops->destroy)(snes);
4757: snes->ops->destroy = NULL;
4758: }
4759: /* Reinitialize function pointers in SNESOps structure */
4760: snes->ops->setup = NULL;
4761: snes->ops->solve = NULL;
4762: snes->ops->view = NULL;
4763: snes->ops->setfromoptions = NULL;
4764: snes->ops->destroy = NULL;
4766: /* It may happen the user has customized the line search before calling SNESSetType */
4767: if (((PetscObject)snes)->type_name) SNESLineSearchDestroy(&snes->linesearch);
4769: /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
4770: snes->setupcalled = PETSC_FALSE;
4772: PetscObjectChangeTypeName((PetscObject)snes,type);
4773: (*r)(snes);
4774: return 0;
4775: }
4777: /*@C
4778: SNESGetType - Gets the SNES method type and name (as a string).
4780: Not Collective
4782: Input Parameter:
4783: . snes - nonlinear solver context
4785: Output Parameter:
4786: . type - SNES method (a character string)
4788: Level: intermediate
4790: @*/
4791: PetscErrorCode SNESGetType(SNES snes,SNESType *type)
4792: {
4795: *type = ((PetscObject)snes)->type_name;
4796: return 0;
4797: }
4799: /*@
4800: SNESSetSolution - Sets the solution vector for use by the SNES routines.
4802: Logically Collective on SNES
4804: Input Parameters:
4805: + snes - the SNES context obtained from SNESCreate()
4806: - u - the solution vector
4808: Level: beginner
4810: @*/
4811: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4812: {
4813: DM dm;
4817: PetscObjectReference((PetscObject) u);
4818: VecDestroy(&snes->vec_sol);
4820: snes->vec_sol = u;
4822: SNESGetDM(snes, &dm);
4823: DMShellSetGlobalVector(dm, u);
4824: return 0;
4825: }
4827: /*@
4828: SNESGetSolution - Returns the vector where the approximate solution is
4829: stored. This is the fine grid solution when using SNESSetGridSequence().
4831: Not Collective, but Vec is parallel if SNES is parallel
4833: Input Parameter:
4834: . snes - the SNES context
4836: Output Parameter:
4837: . x - the solution
4839: Level: intermediate
4841: .seealso: SNESGetSolutionUpdate(), SNESGetFunction()
4842: @*/
4843: PetscErrorCode SNESGetSolution(SNES snes,Vec *x)
4844: {
4847: *x = snes->vec_sol;
4848: return 0;
4849: }
4851: /*@
4852: SNESGetSolutionUpdate - Returns the vector where the solution update is
4853: stored.
4855: Not Collective, but Vec is parallel if SNES is parallel
4857: Input Parameter:
4858: . snes - the SNES context
4860: Output Parameter:
4861: . x - the solution update
4863: Level: advanced
4865: .seealso: SNESGetSolution(), SNESGetFunction()
4866: @*/
4867: PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x)
4868: {
4871: *x = snes->vec_sol_update;
4872: return 0;
4873: }
4875: /*@C
4876: SNESGetFunction - Returns the vector where the function is stored.
4878: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
4880: Input Parameter:
4881: . snes - the SNES context
4883: Output Parameters:
4884: + r - the vector that is used to store residuals (or NULL if you don't want it)
4885: . f - the function (or NULL if you don't want it); see SNESFunction for calling sequence details
4886: - ctx - the function context (or NULL if you don't want it)
4888: Level: advanced
4890: Notes: The vector r DOES NOT, in general contain the current value of the SNES nonlinear function
4892: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4893: @*/
4894: PetscErrorCode SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4895: {
4896: DM dm;
4899: if (r) {
4900: if (!snes->vec_func) {
4901: if (snes->vec_rhs) {
4902: VecDuplicate(snes->vec_rhs,&snes->vec_func);
4903: } else if (snes->vec_sol) {
4904: VecDuplicate(snes->vec_sol,&snes->vec_func);
4905: } else if (snes->dm) {
4906: DMCreateGlobalVector(snes->dm,&snes->vec_func);
4907: }
4908: }
4909: *r = snes->vec_func;
4910: }
4911: SNESGetDM(snes,&dm);
4912: DMSNESGetFunction(dm,f,ctx);
4913: return 0;
4914: }
4916: /*@C
4917: SNESGetNGS - Returns the NGS function and context.
4919: Input Parameter:
4920: . snes - the SNES context
4922: Output Parameters:
4923: + f - the function (or NULL) see SNESNGSFunction for details
4924: - ctx - the function context (or NULL)
4926: Level: advanced
4928: .seealso: SNESSetNGS(), SNESGetFunction()
4929: @*/
4931: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4932: {
4933: DM dm;
4936: SNESGetDM(snes,&dm);
4937: DMSNESGetNGS(dm,f,ctx);
4938: return 0;
4939: }
4941: /*@C
4942: SNESSetOptionsPrefix - Sets the prefix used for searching for all
4943: SNES options in the database.
4945: Logically Collective on SNES
4947: Input Parameters:
4948: + snes - the SNES context
4949: - prefix - the prefix to prepend to all option names
4951: Notes:
4952: A hyphen (-) must NOT be given at the beginning of the prefix name.
4953: The first character of all runtime options is AUTOMATICALLY the hyphen.
4955: Level: advanced
4957: .seealso: SNESSetFromOptions()
4958: @*/
4959: PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[])
4960: {
4962: PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4963: if (!snes->ksp) SNESGetKSP(snes,&snes->ksp);
4964: if (snes->linesearch) {
4965: SNESGetLineSearch(snes,&snes->linesearch);
4966: PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4967: }
4968: KSPSetOptionsPrefix(snes->ksp,prefix);
4969: return 0;
4970: }
4972: /*@C
4973: SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4974: SNES options in the database.
4976: Logically Collective on SNES
4978: Input Parameters:
4979: + snes - the SNES context
4980: - prefix - the prefix to prepend to all option names
4982: Notes:
4983: A hyphen (-) must NOT be given at the beginning of the prefix name.
4984: The first character of all runtime options is AUTOMATICALLY the hyphen.
4986: Level: advanced
4988: .seealso: SNESGetOptionsPrefix()
4989: @*/
4990: PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4991: {
4993: PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4994: if (!snes->ksp) SNESGetKSP(snes,&snes->ksp);
4995: if (snes->linesearch) {
4996: SNESGetLineSearch(snes,&snes->linesearch);
4997: PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4998: }
4999: KSPAppendOptionsPrefix(snes->ksp,prefix);
5000: return 0;
5001: }
5003: /*@C
5004: SNESGetOptionsPrefix - Sets the prefix used for searching for all
5005: SNES options in the database.
5007: Not Collective
5009: Input Parameter:
5010: . snes - the SNES context
5012: Output Parameter:
5013: . prefix - pointer to the prefix string used
5015: Notes:
5016: On the fortran side, the user should pass in a string 'prefix' of
5017: sufficient length to hold the prefix.
5019: Level: advanced
5021: .seealso: SNESAppendOptionsPrefix()
5022: @*/
5023: PetscErrorCode SNESGetOptionsPrefix(SNES snes,const char *prefix[])
5024: {
5026: PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
5027: return 0;
5028: }
5030: /*@C
5031: SNESRegister - Adds a method to the nonlinear solver package.
5033: Not collective
5035: Input Parameters:
5036: + name_solver - name of a new user-defined solver
5037: - routine_create - routine to create method context
5039: Notes:
5040: SNESRegister() may be called multiple times to add several user-defined solvers.
5042: Sample usage:
5043: .vb
5044: SNESRegister("my_solver",MySolverCreate);
5045: .ve
5047: Then, your solver can be chosen with the procedural interface via
5048: $ SNESSetType(snes,"my_solver")
5049: or at runtime via the option
5050: $ -snes_type my_solver
5052: Level: advanced
5054: Note: If your function is not being put into a shared library then use SNESRegister() instead
5056: .seealso: SNESRegisterAll(), SNESRegisterDestroy()
5058: Level: advanced
5059: @*/
5060: PetscErrorCode SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
5061: {
5062: SNESInitializePackage();
5063: PetscFunctionListAdd(&SNESList,sname,function);
5064: return 0;
5065: }
5067: PetscErrorCode SNESTestLocalMin(SNES snes)
5068: {
5069: PetscInt N,i,j;
5070: Vec u,uh,fh;
5071: PetscScalar value;
5072: PetscReal norm;
5074: SNESGetSolution(snes,&u);
5075: VecDuplicate(u,&uh);
5076: VecDuplicate(u,&fh);
5078: /* currently only works for sequential */
5079: PetscPrintf(PetscObjectComm((PetscObject)snes),"Testing FormFunction() for local min\n");
5080: VecGetSize(u,&N);
5081: for (i=0; i<N; i++) {
5082: VecCopy(u,uh);
5083: PetscPrintf(PetscObjectComm((PetscObject)snes),"i = %D\n",i);
5084: for (j=-10; j<11; j++) {
5085: value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
5086: VecSetValue(uh,i,value,ADD_VALUES);
5087: SNESComputeFunction(snes,uh,fh);
5088: VecNorm(fh,NORM_2,&norm);
5089: PetscPrintf(PetscObjectComm((PetscObject)snes)," j norm %D %18.16e\n",j,norm);
5090: value = -value;
5091: VecSetValue(uh,i,value,ADD_VALUES);
5092: }
5093: }
5094: VecDestroy(&uh);
5095: VecDestroy(&fh);
5096: return 0;
5097: }
5099: /*@
5100: SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
5101: computing relative tolerance for linear solvers within an inexact
5102: Newton method.
5104: Logically Collective on SNES
5106: Input Parameters:
5107: + snes - SNES context
5108: - flag - PETSC_TRUE or PETSC_FALSE
5110: Options Database:
5111: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
5112: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
5113: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
5114: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
5115: . -snes_ksp_ew_gamma <gamma> - Sets gamma
5116: . -snes_ksp_ew_alpha <alpha> - Sets alpha
5117: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
5118: - -snes_ksp_ew_threshold <threshold> - Sets threshold
5120: Notes:
5121: Currently, the default is to use a constant relative tolerance for
5122: the inner linear solvers. Alternatively, one can use the
5123: Eisenstat-Walker method, where the relative convergence tolerance
5124: is reset at each Newton iteration according progress of the nonlinear
5125: solver.
5127: Level: advanced
5129: Reference:
5130: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
5131: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
5133: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
5134: @*/
5135: PetscErrorCode SNESKSPSetUseEW(SNES snes,PetscBool flag)
5136: {
5139: snes->ksp_ewconv = flag;
5140: return 0;
5141: }
5143: /*@
5144: SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
5145: for computing relative tolerance for linear solvers within an
5146: inexact Newton method.
5148: Not Collective
5150: Input Parameter:
5151: . snes - SNES context
5153: Output Parameter:
5154: . flag - PETSC_TRUE or PETSC_FALSE
5156: Notes:
5157: Currently, the default is to use a constant relative tolerance for
5158: the inner linear solvers. Alternatively, one can use the
5159: Eisenstat-Walker method, where the relative convergence tolerance
5160: is reset at each Newton iteration according progress of the nonlinear
5161: solver.
5163: Level: advanced
5165: Reference:
5166: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
5167: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
5169: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
5170: @*/
5171: PetscErrorCode SNESKSPGetUseEW(SNES snes, PetscBool *flag)
5172: {
5175: *flag = snes->ksp_ewconv;
5176: return 0;
5177: }
5179: /*@
5180: SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
5181: convergence criteria for the linear solvers within an inexact
5182: Newton method.
5184: Logically Collective on SNES
5186: Input Parameters:
5187: + snes - SNES context
5188: . version - version 1, 2 (default is 2) or 3
5189: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5190: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5191: . gamma - multiplicative factor for version 2 rtol computation
5192: (0 <= gamma2 <= 1)
5193: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
5194: . alpha2 - power for safeguard
5195: - threshold - threshold for imposing safeguard (0 < threshold < 1)
5197: Note:
5198: Version 3 was contributed by Luis Chacon, June 2006.
5200: Use PETSC_DEFAULT to retain the default for any of the parameters.
5202: Level: advanced
5204: Reference:
5205: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
5206: inexact Newton method", Utah State University Math. Stat. Dept. Res.
5207: Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
5209: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
5210: @*/
5211: PetscErrorCode SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
5212: {
5213: SNESKSPEW *kctx;
5216: kctx = (SNESKSPEW*)snes->kspconvctx;
5226: if (version != PETSC_DEFAULT) kctx->version = version;
5227: if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0;
5228: if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max;
5229: if (gamma != PETSC_DEFAULT) kctx->gamma = gamma;
5230: if (alpha != PETSC_DEFAULT) kctx->alpha = alpha;
5231: if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2;
5232: if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
5240: return 0;
5241: }
5243: /*@
5244: SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
5245: convergence criteria for the linear solvers within an inexact
5246: Newton method.
5248: Not Collective
5250: Input Parameter:
5251: . snes - SNES context
5253: Output Parameters:
5254: + version - version 1, 2 (default is 2) or 3
5255: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
5256: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
5257: . gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
5258: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
5259: . alpha2 - power for safeguard
5260: - threshold - threshold for imposing safeguard (0 < threshold < 1)
5262: Level: advanced
5264: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
5265: @*/
5266: PetscErrorCode SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
5267: {
5268: SNESKSPEW *kctx;
5271: kctx = (SNESKSPEW*)snes->kspconvctx;
5273: if (version) *version = kctx->version;
5274: if (rtol_0) *rtol_0 = kctx->rtol_0;
5275: if (rtol_max) *rtol_max = kctx->rtol_max;
5276: if (gamma) *gamma = kctx->gamma;
5277: if (alpha) *alpha = kctx->alpha;
5278: if (alpha2) *alpha2 = kctx->alpha2;
5279: if (threshold) *threshold = kctx->threshold;
5280: return 0;
5281: }
5283: PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5284: {
5285: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
5286: PetscReal rtol = PETSC_DEFAULT,stol;
5288: if (!snes->ksp_ewconv) return 0;
5289: if (!snes->iter) {
5290: rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
5291: VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
5292: }
5293: else {
5294: if (kctx->version == 1) {
5295: rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
5296: if (rtol < 0.0) rtol = -rtol;
5297: stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
5298: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5299: } else if (kctx->version == 2) {
5300: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5301: stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
5302: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5303: } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
5304: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5305: /* safeguard: avoid sharp decrease of rtol */
5306: stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
5307: stol = PetscMax(rtol,stol);
5308: rtol = PetscMin(kctx->rtol_0,stol);
5309: /* safeguard: avoid oversolving */
5310: stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
5311: stol = PetscMax(rtol,stol);
5312: rtol = PetscMin(kctx->rtol_0,stol);
5313: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
5314: }
5315: /* safeguard: avoid rtol greater than one */
5316: rtol = PetscMin(rtol,kctx->rtol_max);
5317: KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
5318: PetscInfo(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
5319: return 0;
5320: }
5322: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5323: {
5324: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
5325: PCSide pcside;
5326: Vec lres;
5328: if (!snes->ksp_ewconv) return 0;
5329: KSPGetTolerances(ksp,&kctx->rtol_last,NULL,NULL,NULL);
5330: kctx->norm_last = snes->norm;
5331: if (kctx->version == 1) {
5332: PC pc;
5333: PetscBool isNone;
5335: KSPGetPC(ksp, &pc);
5336: PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
5337: KSPGetPCSide(ksp,&pcside);
5338: if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
5339: /* KSP residual is true linear residual */
5340: KSPGetResidualNorm(ksp,&kctx->lresid_last);
5341: } else {
5342: /* KSP residual is preconditioned residual */
5343: /* compute true linear residual norm */
5344: VecDuplicate(b,&lres);
5345: MatMult(snes->jacobian,x,lres);
5346: VecAYPX(lres,-1.0,b);
5347: VecNorm(lres,NORM_2,&kctx->lresid_last);
5348: VecDestroy(&lres);
5349: }
5350: }
5351: return 0;
5352: }
5354: /*@
5355: SNESGetKSP - Returns the KSP context for a SNES solver.
5357: Not Collective, but if SNES object is parallel, then KSP object is parallel
5359: Input Parameter:
5360: . snes - the SNES context
5362: Output Parameter:
5363: . ksp - the KSP context
5365: Notes:
5366: The user can then directly manipulate the KSP context to set various
5367: options, etc. Likewise, the user can then extract and manipulate the
5368: PC contexts as well.
5370: Level: beginner
5372: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
5373: @*/
5374: PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp)
5375: {
5379: if (!snes->ksp) {
5380: KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
5381: PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
5382: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);
5384: KSPSetPreSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPreSolve_SNESEW,snes);
5385: KSPSetPostSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPostSolve_SNESEW,snes);
5387: KSPMonitorSetFromOptions(snes->ksp, "-snes_monitor_ksp", "snes_preconditioned_residual", snes);
5388: PetscObjectSetOptions((PetscObject)snes->ksp,((PetscObject)snes)->options);
5389: }
5390: *ksp = snes->ksp;
5391: return 0;
5392: }
5394: #include <petsc/private/dmimpl.h>
5395: /*@
5396: SNESSetDM - Sets the DM that may be used by some nonlinear solvers or their underlying preconditioners
5398: Logically Collective on SNES
5400: Input Parameters:
5401: + snes - the nonlinear solver context
5402: - dm - the dm, cannot be NULL
5404: Notes:
5405: A DM can only be used for solving one problem at a time because information about the problem is stored on the DM,
5406: even when not using interfaces like DMSNESSetFunction(). Use DMClone() to get a distinct DM when solving different
5407: problems using the same function space.
5409: Level: intermediate
5411: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
5412: @*/
5413: PetscErrorCode SNESSetDM(SNES snes,DM dm)
5414: {
5415: KSP ksp;
5416: DMSNES sdm;
5420: PetscObjectReference((PetscObject)dm);
5421: if (snes->dm) { /* Move the DMSNES context over to the new DM unless the new DM already has one */
5422: if (snes->dm->dmsnes && !dm->dmsnes) {
5423: DMCopyDMSNES(snes->dm,dm);
5424: DMGetDMSNES(snes->dm,&sdm);
5425: if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
5426: }
5427: DMCoarsenHookRemove(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
5428: DMDestroy(&snes->dm);
5429: }
5430: snes->dm = dm;
5431: snes->dmAuto = PETSC_FALSE;
5433: SNESGetKSP(snes,&ksp);
5434: KSPSetDM(ksp,dm);
5435: KSPSetDMActive(ksp,PETSC_FALSE);
5436: if (snes->npc) {
5437: SNESSetDM(snes->npc,snes->dm);
5438: SNESSetNPCSide(snes,snes->npcside);
5439: }
5440: return 0;
5441: }
5443: /*@
5444: SNESGetDM - Gets the DM that may be used by some preconditioners
5446: Not Collective but DM obtained is parallel on SNES
5448: Input Parameter:
5449: . snes - the preconditioner context
5451: Output Parameter:
5452: . dm - the dm
5454: Level: intermediate
5456: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
5457: @*/
5458: PetscErrorCode SNESGetDM(SNES snes,DM *dm)
5459: {
5461: if (!snes->dm) {
5462: DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
5463: snes->dmAuto = PETSC_TRUE;
5464: }
5465: *dm = snes->dm;
5466: return 0;
5467: }
5469: /*@
5470: SNESSetNPC - Sets the nonlinear preconditioner to be used.
5472: Collective on SNES
5474: Input Parameters:
5475: + snes - iterative context obtained from SNESCreate()
5476: - pc - the preconditioner object
5478: Notes:
5479: Use SNESGetNPC() to retrieve the preconditioner context (for example,
5480: to configure it using the API).
5482: Level: developer
5484: .seealso: SNESGetNPC(), SNESHasNPC()
5485: @*/
5486: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
5487: {
5491: PetscObjectReference((PetscObject) pc);
5492: SNESDestroy(&snes->npc);
5493: snes->npc = pc;
5494: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->npc);
5495: return 0;
5496: }
5498: /*@
5499: SNESGetNPC - Creates a nonlinear preconditioning solver (SNES) to be used to precondition the nonlinear solver.
5501: Not Collective; but any changes to the obtained SNES object must be applied collectively
5503: Input Parameter:
5504: . snes - iterative context obtained from SNESCreate()
5506: Output Parameter:
5507: . pc - preconditioner context
5509: Options Database:
5510: . -npc_snes_type <type> - set the type of the SNES to use as the nonlinear preconditioner
5512: Notes:
5513: If a SNES was previously set with SNESSetNPC() then that SNES is returned, otherwise a new SNES object is created.
5515: The (preconditioner) SNES returned automatically inherits the same nonlinear function and Jacobian supplied to the original
5516: SNES during SNESSetUp()
5518: Level: developer
5520: .seealso: SNESSetNPC(), SNESHasNPC(), SNES, SNESCreate()
5521: @*/
5522: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
5523: {
5524: const char *optionsprefix;
5528: if (!snes->npc) {
5529: SNESCreate(PetscObjectComm((PetscObject)snes),&snes->npc);
5530: PetscObjectIncrementTabLevel((PetscObject)snes->npc,(PetscObject)snes,1);
5531: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->npc);
5532: SNESGetOptionsPrefix(snes,&optionsprefix);
5533: SNESSetOptionsPrefix(snes->npc,optionsprefix);
5534: SNESAppendOptionsPrefix(snes->npc,"npc_");
5535: SNESSetCountersReset(snes->npc,PETSC_FALSE);
5536: }
5537: *pc = snes->npc;
5538: return 0;
5539: }
5541: /*@
5542: SNESHasNPC - Returns whether a nonlinear preconditioner exists
5544: Not Collective
5546: Input Parameter:
5547: . snes - iterative context obtained from SNESCreate()
5549: Output Parameter:
5550: . has_npc - whether the SNES has an NPC or not
5552: Level: developer
5554: .seealso: SNESSetNPC(), SNESGetNPC()
5555: @*/
5556: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
5557: {
5559: *has_npc = (PetscBool) (snes->npc ? PETSC_TRUE : PETSC_FALSE);
5560: return 0;
5561: }
5563: /*@
5564: SNESSetNPCSide - Sets the preconditioning side.
5566: Logically Collective on SNES
5568: Input Parameter:
5569: . snes - iterative context obtained from SNESCreate()
5571: Output Parameter:
5572: . side - the preconditioning side, where side is one of
5573: .vb
5574: PC_LEFT - left preconditioning
5575: PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5576: .ve
5578: Options Database Keys:
5579: . -snes_npc_side <right,left> - nonlinear preconditioner side
5581: Notes:
5582: SNESNRICHARDSON and SNESNCG only support left preconditioning.
5584: Level: intermediate
5586: .seealso: SNESGetNPCSide(), KSPSetPCSide()
5587: @*/
5588: PetscErrorCode SNESSetNPCSide(SNES snes,PCSide side)
5589: {
5592: if (side == PC_SIDE_DEFAULT) side = PC_RIGHT;
5594: snes->npcside = side;
5595: return 0;
5596: }
5598: /*@
5599: SNESGetNPCSide - Gets the preconditioning side.
5601: Not Collective
5603: Input Parameter:
5604: . snes - iterative context obtained from SNESCreate()
5606: Output Parameter:
5607: . side - the preconditioning side, where side is one of
5608: .vb
5609: PC_LEFT - left preconditioning
5610: PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5611: .ve
5613: Level: intermediate
5615: .seealso: SNESSetNPCSide(), KSPGetPCSide()
5616: @*/
5617: PetscErrorCode SNESGetNPCSide(SNES snes,PCSide *side)
5618: {
5621: *side = snes->npcside;
5622: return 0;
5623: }
5625: /*@
5626: SNESSetLineSearch - Sets the linesearch on the SNES instance.
5628: Collective on SNES
5630: Input Parameters:
5631: + snes - iterative context obtained from SNESCreate()
5632: - linesearch - the linesearch object
5634: Notes:
5635: Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
5636: to configure it using the API).
5638: Level: developer
5640: .seealso: SNESGetLineSearch()
5641: @*/
5642: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
5643: {
5647: PetscObjectReference((PetscObject) linesearch);
5648: SNESLineSearchDestroy(&snes->linesearch);
5650: snes->linesearch = linesearch;
5652: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5653: return 0;
5654: }
5656: /*@
5657: SNESGetLineSearch - Returns a pointer to the line search context set with SNESSetLineSearch()
5658: or creates a default line search instance associated with the SNES and returns it.
5660: Not Collective
5662: Input Parameter:
5663: . snes - iterative context obtained from SNESCreate()
5665: Output Parameter:
5666: . linesearch - linesearch context
5668: Level: beginner
5670: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5671: @*/
5672: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5673: {
5674: const char *optionsprefix;
5678: if (!snes->linesearch) {
5679: SNESGetOptionsPrefix(snes, &optionsprefix);
5680: SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5681: SNESLineSearchSetSNES(snes->linesearch, snes);
5682: SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5683: PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5684: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5685: }
5686: *linesearch = snes->linesearch;
5687: return 0;
5688: }