Actual source code: snes.c
petsc-3.5.4 2015-05-23
2: #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
3: #include <petscdmshell.h>
5: PetscBool SNESRegisterAllCalled = PETSC_FALSE;
6: PetscFunctionList SNESList = NULL;
8: /* Logging support */
9: PetscClassId SNES_CLASSID, DMSNES_CLASSID;
10: PetscLogEvent SNES_Solve, SNES_FunctionEval, SNES_JacobianEval, SNES_NGSEval, SNES_NGSFuncEval, SNES_NPCSolve;
14: /*@
15: SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged.
17: Logically Collective on SNES
19: Input Parameters:
20: + snes - iterative context obtained from SNESCreate()
21: - flg - PETSC_TRUE indicates you want the error generated
23: Options database keys:
24: . -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)
26: Level: intermediate
28: Notes:
29: Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve()
30: to determine if it has converged.
32: .keywords: SNES, set, initial guess, nonzero
34: .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
35: @*/
36: PetscErrorCode SNESSetErrorIfNotConverged(SNES snes,PetscBool flg)
37: {
41: snes->errorifnotconverged = flg;
42: return(0);
43: }
47: /*@
48: SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge?
50: Not Collective
52: Input Parameter:
53: . snes - iterative context obtained from SNESCreate()
55: Output Parameter:
56: . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE
58: Level: intermediate
60: .keywords: SNES, set, initial guess, nonzero
62: .seealso: SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
63: @*/
64: PetscErrorCode SNESGetErrorIfNotConverged(SNES snes,PetscBool *flag)
65: {
69: *flag = snes->errorifnotconverged;
70: return(0);
71: }
75: /*@
76: SNESSetFunctionDomainError - tells SNES that the input vector to your SNESFunction is not
77: in the functions domain. For example, negative pressure.
79: Logically Collective on SNES
81: Input Parameters:
82: . snes - the SNES context
84: Level: advanced
86: .keywords: SNES, view
88: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction
89: @*/
90: PetscErrorCode SNESSetFunctionDomainError(SNES snes)
91: {
94: snes->domainerror = PETSC_TRUE;
95: return(0);
96: }
100: /*@
101: SNESGetFunctionDomainError - Gets the status of the domain error after a call to SNESComputeFunction;
103: Logically Collective on SNES
105: Input Parameters:
106: . snes - the SNES context
108: Output Parameters:
109: . domainerror - Set to PETSC_TRUE if there's a domain error; PETSC_FALSE otherwise.
111: Level: advanced
113: .keywords: SNES, view
115: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction()
116: @*/
117: PetscErrorCode SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror)
118: {
122: *domainerror = snes->domainerror;
123: return(0);
124: }
128: /*@C
129: SNESLoad - Loads a SNES that has been stored in binary with SNESView().
131: Collective on PetscViewer
133: Input Parameters:
134: + newdm - the newly loaded SNES, this needs to have been created with SNESCreate() or
135: some related function before a call to SNESLoad().
136: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
138: Level: intermediate
140: Notes:
141: The type is determined by the data in the file, any type set into the SNES before this call is ignored.
143: Notes for advanced users:
144: Most users should not need to know the details of the binary storage
145: format, since SNESLoad() and TSView() completely hide these details.
146: But for anyone who's interested, the standard binary matrix storage
147: format is
148: .vb
149: has not yet been determined
150: .ve
152: .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad()
153: @*/
154: PetscErrorCode SNESLoad(SNES snes, PetscViewer viewer)
155: {
157: PetscBool isbinary;
158: PetscInt classid;
159: char type[256];
160: KSP ksp;
161: DM dm;
162: DMSNES dmsnes;
167: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
168: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
170: PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);
171: if (classid != SNES_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONG,"Not SNES next in file");
172: PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);
173: SNESSetType(snes, type);
174: if (snes->ops->load) {
175: (*snes->ops->load)(snes,viewer);
176: }
177: SNESGetDM(snes,&dm);
178: DMGetDMSNES(dm,&dmsnes);
179: DMSNESLoad(dmsnes,viewer);
180: SNESGetKSP(snes,&ksp);
181: KSPLoad(ksp,viewer);
182: return(0);
183: }
185: #include <petscdraw.h>
186: #if defined(PETSC_HAVE_SAWS)
187: #include <petscviewersaws.h>
188: #endif
191: /*@C
192: SNESView - Prints the SNES data structure.
194: Collective on SNES
196: Input Parameters:
197: + SNES - the SNES context
198: - viewer - visualization context
200: Options Database Key:
201: . -snes_view - Calls SNESView() at end of SNESSolve()
203: Notes:
204: The available visualization contexts include
205: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
206: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
207: output where only the first processor opens
208: the file. All other processors send their
209: data to the first processor to print.
211: The user can open an alternative visualization context with
212: PetscViewerASCIIOpen() - output to a specified file.
214: Level: beginner
216: .keywords: SNES, view
218: .seealso: PetscViewerASCIIOpen()
219: @*/
220: PetscErrorCode SNESView(SNES snes,PetscViewer viewer)
221: {
222: SNESKSPEW *kctx;
224: KSP ksp;
225: SNESLineSearch linesearch;
226: PetscBool iascii,isstring,isbinary,isdraw;
227: DMSNES dmsnes;
228: #if defined(PETSC_HAVE_SAWS)
229: PetscBool isams;
230: #endif
234: if (!viewer) {
235: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
236: }
240: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
241: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
242: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
243: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
244: #if defined(PETSC_HAVE_SAWS)
245: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&isams);
246: #endif
247: if (iascii) {
248: PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);
249: if (!snes->setupcalled) {
250: PetscViewerASCIIPrintf(viewer," SNES has not been set up so information may be incomplete\n");
251: }
252: if (snes->ops->view) {
253: PetscViewerASCIIPushTab(viewer);
254: (*snes->ops->view)(snes,viewer);
255: PetscViewerASCIIPopTab(viewer);
256: }
257: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);
258: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);
259: PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);
260: PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);
261: if (snes->gridsequence) {
262: PetscViewerASCIIPrintf(viewer," total number of grid sequence refinements=%D\n",snes->gridsequence);
263: }
264: if (snes->ksp_ewconv) {
265: kctx = (SNESKSPEW*)snes->kspconvctx;
266: if (kctx) {
267: PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);
268: PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);
269: PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);
270: }
271: }
272: if (snes->lagpreconditioner == -1) {
273: PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");
274: } else if (snes->lagpreconditioner > 1) {
275: PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);
276: }
277: if (snes->lagjacobian == -1) {
278: PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");
279: } else if (snes->lagjacobian > 1) {
280: PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);
281: }
282: } else if (isstring) {
283: const char *type;
284: SNESGetType(snes,&type);
285: PetscViewerStringSPrintf(viewer," %-3.3s",type);
286: } else if (isbinary) {
287: PetscInt classid = SNES_FILE_CLASSID;
288: MPI_Comm comm;
289: PetscMPIInt rank;
290: char type[256];
292: PetscObjectGetComm((PetscObject)snes,&comm);
293: MPI_Comm_rank(comm,&rank);
294: if (!rank) {
295: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
296: PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));
297: PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR,PETSC_FALSE);
298: }
299: if (snes->ops->view) {
300: (*snes->ops->view)(snes,viewer);
301: }
302: } else if (isdraw) {
303: PetscDraw draw;
304: char str[36];
305: PetscReal x,y,bottom,h;
307: PetscViewerDrawGetDraw(viewer,0,&draw);
308: PetscDrawGetCurrentPoint(draw,&x,&y);
309: PetscStrcpy(str,"SNES: ");
310: PetscStrcat(str,((PetscObject)snes)->type_name);
311: PetscDrawBoxedString(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);
312: bottom = y - h;
313: PetscDrawPushCurrentPoint(draw,x,bottom);
314: if (snes->ops->view) {
315: (*snes->ops->view)(snes,viewer);
316: }
317: #if defined(PETSC_HAVE_SAWS)
318: } else if (isams) {
319: PetscMPIInt rank;
320: const char *name;
322: PetscObjectGetName((PetscObject)snes,&name);
323: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
324: if (!((PetscObject)snes)->amsmem && !rank) {
325: char dir[1024];
327: PetscObjectViewSAWs((PetscObject)snes,viewer);
328: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
329: PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT));
330: if (!snes->conv_hist) {
331: SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);
332: }
333: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);
334: PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE));
335: }
336: #endif
337: }
338: if (snes->linesearch) {
339: PetscViewerASCIIPushTab(viewer);
340: SNESGetLineSearch(snes, &linesearch);
341: SNESLineSearchView(linesearch, viewer);
342: PetscViewerASCIIPopTab(viewer);
343: }
344: if (snes->pc && snes->usespc) {
345: PetscViewerASCIIPushTab(viewer);
346: SNESView(snes->pc, viewer);
347: PetscViewerASCIIPopTab(viewer);
348: }
349: PetscViewerASCIIPushTab(viewer);
350: DMGetDMSNES(snes->dm,&dmsnes);
351: DMSNESView(dmsnes, viewer);
352: PetscViewerASCIIPopTab(viewer);
353: if (snes->usesksp) {
354: SNESGetKSP(snes,&ksp);
355: PetscViewerASCIIPushTab(viewer);
356: KSPView(ksp,viewer);
357: PetscViewerASCIIPopTab(viewer);
358: }
359: if (isdraw) {
360: PetscDraw draw;
361: PetscViewerDrawGetDraw(viewer,0,&draw);
362: PetscDrawPopCurrentPoint(draw);
363: }
364: return(0);
365: }
367: /*
368: We retain a list of functions that also take SNES command
369: line options. These are called at the end SNESSetFromOptions()
370: */
371: #define MAXSETFROMOPTIONS 5
372: static PetscInt numberofsetfromoptions;
373: static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
377: /*@C
378: SNESAddOptionsChecker - Adds an additional function to check for SNES options.
380: Not Collective
382: Input Parameter:
383: . snescheck - function that checks for options
385: Level: developer
387: .seealso: SNESSetFromOptions()
388: @*/
389: PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
390: {
392: if (numberofsetfromoptions >= MAXSETFROMOPTIONS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
393: othersetfromoptions[numberofsetfromoptions++] = snescheck;
394: return(0);
395: }
397: extern PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
401: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
402: {
403: Mat J;
404: KSP ksp;
405: PC pc;
406: PetscBool match;
412: if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
413: Mat A = snes->jacobian, B = snes->jacobian_pre;
414: MatGetVecs(A ? A : B, NULL,&snes->vec_func);
415: }
417: if (version == 1) {
418: MatCreateSNESMF(snes,&J);
419: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
420: MatSetFromOptions(J);
421: } else if (version == 2) {
422: if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
423: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
424: SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);
425: #else
426: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)");
427: #endif
428: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2");
430: PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);
431: if (hasOperator) {
433: /* This version replaces the user provided Jacobian matrix with a
434: matrix-free version but still employs the user-provided preconditioner matrix. */
435: SNESSetJacobian(snes,J,0,0,0);
436: } else {
437: /* This version replaces both the user-provided Jacobian and the user-
438: provided preconditioner Jacobian with the default matrix free version. */
439: if ((snes->pcside == PC_LEFT) && snes->pc) {
440: if (!snes->jacobian){SNESSetJacobian(snes,J,0,0,0);}
441: } else {
442: SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,0);
443: }
444: /* Force no preconditioner */
445: SNESGetKSP(snes,&ksp);
446: KSPGetPC(ksp,&pc);
447: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);
448: if (!match) {
449: PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");
450: PCSetType(pc,PCNONE);
451: }
452: }
453: MatDestroy(&J);
454: return(0);
455: }
459: static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx)
460: {
461: SNES snes = (SNES)ctx;
463: Vec Xfine,Xfine_named = NULL,Xcoarse;
466: if (PetscLogPrintInfo) {
467: PetscInt finelevel,coarselevel,fineclevel,coarseclevel;
468: DMGetRefineLevel(dmfine,&finelevel);
469: DMGetCoarsenLevel(dmfine,&fineclevel);
470: DMGetRefineLevel(dmcoarse,&coarselevel);
471: DMGetCoarsenLevel(dmcoarse,&coarseclevel);
472: PetscInfo4(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);
473: }
474: if (dmfine == snes->dm) Xfine = snes->vec_sol;
475: else {
476: DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
477: Xfine = Xfine_named;
478: }
479: DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
480: MatRestrict(Restrict,Xfine,Xcoarse);
481: VecPointwiseMult(Xcoarse,Xcoarse,Rscale);
482: DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
483: if (Xfine_named) {DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);}
484: return(0);
485: }
489: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
490: {
494: DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
495: return(0);
496: }
500: /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can
501: * safely call SNESGetDM() in their residual evaluation routine. */
502: static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx)
503: {
504: SNES snes = (SNES)ctx;
506: Mat Asave = A,Bsave = B;
507: Vec X,Xnamed = NULL;
508: DM dmsave;
509: void *ctxsave;
510: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
513: dmsave = snes->dm;
514: KSPGetDM(ksp,&snes->dm);
515: if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */
516: else { /* We are on a coarser level, this vec was initialized using a DM restrict hook */
517: DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
518: X = Xnamed;
519: SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);
520: /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */
521: if (jac == SNESComputeJacobianDefaultColor) {
522: SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,0);
523: }
524: }
525: /* put the previous context back */
527: SNESComputeJacobian(snes,X,A,B);
528: if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
529: SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
530: }
532: if (A != Asave || B != Bsave) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for changing matrices at this time");
533: if (Xnamed) {
534: DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
535: }
536: snes->dm = dmsave;
537: return(0);
538: }
542: /*@
543: SNESSetUpMatrices - ensures that matrices are available for SNES, to be called by SNESSetUp_XXX()
545: Collective
547: Input Arguments:
548: . snes - snes to configure
550: Level: developer
552: .seealso: SNESSetUp()
553: @*/
554: PetscErrorCode SNESSetUpMatrices(SNES snes)
555: {
557: DM dm;
558: DMSNES sdm;
561: SNESGetDM(snes,&dm);
562: DMGetDMSNES(dm,&sdm);
563: if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"DMSNES not properly configured");
564: else if (!snes->jacobian && snes->mf) {
565: Mat J;
566: void *functx;
567: MatCreateSNESMF(snes,&J);
568: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
569: MatSetFromOptions(J);
570: SNESGetFunction(snes,NULL,NULL,&functx);
571: SNESSetJacobian(snes,J,J,0,0);
572: MatDestroy(&J);
573: } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
574: Mat J,B;
575: MatCreateSNESMF(snes,&J);
576: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
577: MatSetFromOptions(J);
578: DMCreateMatrix(snes->dm,&B);
579: /* sdm->computejacobian was already set to reach here */
580: SNESSetJacobian(snes,J,B,NULL,NULL);
581: MatDestroy(&J);
582: MatDestroy(&B);
583: } else if (!snes->jacobian_pre) {
584: Mat J,B;
585: J = snes->jacobian;
586: DMCreateMatrix(snes->dm,&B);
587: SNESSetJacobian(snes,J ? J : B,B,NULL,NULL);
588: MatDestroy(&B);
589: }
590: {
591: KSP ksp;
592: SNESGetKSP(snes,&ksp);
593: KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);
594: DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
595: }
596: return(0);
597: }
601: /*@
602: SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
604: Collective on SNES
606: Input Parameter:
607: . snes - the SNES context
609: Options Database Keys:
610: + -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
611: . -snes_stol - convergence tolerance in terms of the norm
612: of the change in the solution between steps
613: . -snes_atol <abstol> - absolute tolerance of residual norm
614: . -snes_rtol <rtol> - relative decrease in tolerance norm from initial
615: . -snes_max_it <max_it> - maximum number of iterations
616: . -snes_max_funcs <max_funcs> - maximum number of function evaluations
617: . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
618: . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
619: . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
620: . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
621: . -snes_trtol <trtol> - trust region tolerance
622: . -snes_no_convergence_test - skip convergence test in nonlinear
623: solver; hence iterations will continue until max_it
624: or some other criterion is reached. Saves expense
625: of convergence test
626: . -snes_monitor <optional filename> - prints residual norm at each iteration. if no
627: filename given prints to stdout
628: . -snes_monitor_solution - plots solution at each iteration
629: . -snes_monitor_residual - plots residual (not its norm) at each iteration
630: . -snes_monitor_solution_update - plots update to solution at each iteration
631: . -snes_monitor_lg_residualnorm - plots residual norm at each iteration
632: . -snes_monitor_lg_range - plots residual norm at each iteration
633: . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
634: . -snes_fd_color - use finite differences with coloring to compute Jacobian
635: . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
636: - -snes_converged_reason - print the reason for convergence/divergence after each solve
638: Options Database for Eisenstat-Walker method:
639: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
640: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
641: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
642: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
643: . -snes_ksp_ew_gamma <gamma> - Sets gamma
644: . -snes_ksp_ew_alpha <alpha> - Sets alpha
645: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
646: - -snes_ksp_ew_threshold <threshold> - Sets threshold
648: Notes:
649: To see all options, run your program with the -help option or consult
650: Users-Manual: ch_snes
652: Level: beginner
654: .keywords: SNES, nonlinear, set, options, database
656: .seealso: SNESSetOptionsPrefix()
657: @*/
658: PetscErrorCode SNESSetFromOptions(SNES snes)
659: {
660: PetscBool flg,pcset,persist;
661: PetscInt i,indx,lag,grids;
662: const char *deft = SNESNEWTONLS;
663: const char *convtests[] = {"default","skip"};
664: SNESKSPEW *kctx = NULL;
665: char type[256], monfilename[PETSC_MAX_PATH_LEN];
666: PetscViewer monviewer;
668: PCSide pcside;
669: const char *optionsprefix;
673: if (!SNESRegisterAllCalled) {SNESRegisterAll();}
674: PetscObjectOptionsBegin((PetscObject)snes);
675: if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
676: PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
677: if (flg) {
678: SNESSetType(snes,type);
679: } else if (!((PetscObject)snes)->type_name) {
680: SNESSetType(snes,deft);
681: }
682: PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,0);
683: PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,0);
685: PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,0);
686: PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,NULL);
687: PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,NULL);
688: PetscOptionsInt("-snes_max_fail","Maximum nonlinear step failures","SNESSetMaxNonlinearStepFailures",snes->maxFailures,&snes->maxFailures,NULL);
689: PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,NULL);
690: PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,NULL);
692: PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
693: if (flg) {
694: SNESSetLagPreconditioner(snes,lag);
695: }
696: PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
697: if (flg) {
698: SNESSetLagPreconditionerPersists(snes,persist);
699: }
700: PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
701: if (flg) {
702: SNESSetLagJacobian(snes,lag);
703: }
704: PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
705: if (flg) {
706: SNESSetLagJacobianPersists(snes,persist);
707: }
709: PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
710: if (flg) {
711: SNESSetGridSequence(snes,grids);
712: }
714: PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);
715: if (flg) {
716: switch (indx) {
717: case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
718: case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL); break;
719: }
720: }
722: PetscOptionsBool("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",snes->printreason,&snes->printreason,NULL);
724: PetscOptionsEList("-snes_norm_schedule","SNES Norm schedule","SNESSetNormSchedule",SNESNormSchedules,5,"function",&indx,&flg);
725: if (flg) { SNESSetNormSchedule(snes,(SNESNormSchedule)indx); }
727: PetscOptionsEList("-snes_function_type","SNES Norm schedule","SNESSetFunctionType",SNESFunctionTypes,2,"unpreconditioned",&indx,&flg);
728: if (flg) { SNESSetFunctionType(snes,(SNESFunctionType)indx); }
730: kctx = (SNESKSPEW*)snes->kspconvctx;
732: PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,NULL);
734: PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);
735: PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);
736: PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);
737: PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);
738: PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);
739: PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);
740: PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);
742: flg = PETSC_FALSE;
743: PetscOptionsBool("-snes_check_jacobian","Check each Jacobian with a differenced one","SNESUpdateCheckJacobian",flg,&flg,NULL);
744: if (flg) {
745: SNESSetUpdate(snes,SNESUpdateCheckJacobian);
746: }
748: flg = PETSC_FALSE;
749: PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,NULL);
750: if (flg) {SNESMonitorCancel(snes);}
752: PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
753: if (flg) {
754: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
755: SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
756: }
758: PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
759: if (flg) {
760: SNESMonitorSet(snes,SNESMonitorRange,0,0);
761: }
763: PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
764: if (flg) {
765: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
766: SNESMonitorSetRatio(snes,monviewer);
767: }
769: PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
770: if (flg) {
771: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
772: SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
773: }
775: PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
776: if (flg) {PetscPythonMonitorSet((PetscObject)snes,monfilename);}
778: flg = PETSC_FALSE;
779: PetscOptionsBool("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,NULL);
780: if (flg) {SNESMonitorSet(snes,SNESMonitorSolution,0,0);}
781: flg = PETSC_FALSE;
782: PetscOptionsBool("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,NULL);
783: if (flg) {SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);}
784: flg = PETSC_FALSE;
785: PetscOptionsBool("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,NULL);
786: if (flg) {SNESMonitorSet(snes,SNESMonitorResidual,0,0);}
787: flg = PETSC_FALSE;
788: PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);
789: if (flg) {
790: PetscDrawLG ctx;
792: SNESMonitorLGCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
793: SNESMonitorSet(snes,SNESMonitorLGResidualNorm,ctx,(PetscErrorCode (*)(void**))SNESMonitorLGDestroy);
794: }
795: flg = PETSC_FALSE;
796: PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
797: if (flg) {
798: PetscViewer ctx;
800: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
801: SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
802: }
804: flg = PETSC_FALSE;
805: PetscOptionsBool("-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",flg,&flg,NULL);
806: if (flg) {SNESMonitorSet(snes,SNESMonitorJacUpdateSpectrum,0,0);}
809: PetscOptionsString("-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
810: if (flg) {
811: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
812: SNESMonitorSet(snes,SNESMonitorFields,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
813: }
815: flg = PETSC_FALSE;
816: PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
817: if (flg) {
818: void *functx;
819: SNESGetFunction(snes,NULL,NULL,&functx);
820: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
821: PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
822: }
824: flg = PETSC_FALSE;
825: PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
826: if (flg) {
827: SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
828: }
830: flg = PETSC_FALSE;
831: PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);
832: if (flg) {
833: DM dm;
834: DMSNES sdm;
835: SNESGetDM(snes,&dm);
836: DMGetDMSNES(dm,&sdm);
837: sdm->jacobianctx = NULL;
838: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,0);
839: PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
840: }
842: flg = PETSC_FALSE;
843: PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&snes->mf_operator,&flg);
844: if (flg && snes->mf_operator) {
845: snes->mf_operator = PETSC_TRUE;
846: snes->mf = PETSC_TRUE;
847: }
848: flg = PETSC_FALSE;
849: PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&snes->mf,&flg);
850: if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE;
851: PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,0);
853: flg = PETSC_FALSE;
854: SNESGetNPCSide(snes,&pcside);
855: PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
856: if (flg) {SNESSetNPCSide(snes,pcside);}
858: #if defined(PETSC_HAVE_SAWS)
859: /*
860: Publish convergence information using SAWs
861: */
862: flg = PETSC_FALSE;
863: PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);
864: if (flg) {
865: void *ctx;
866: SNESMonitorSAWsCreate(snes,&ctx);
867: SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);
868: }
869: #endif
870: #if defined(PETSC_HAVE_SAWS)
871: {
872: PetscBool set;
873: flg = PETSC_FALSE;
874: PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);
875: if (set) {
876: PetscObjectSAWsSetBlock((PetscObject)snes,flg);
877: }
878: }
879: #endif
881: for (i = 0; i < numberofsetfromoptions; i++) {
882: (*othersetfromoptions[i])(snes);
883: }
885: if (snes->ops->setfromoptions) {
886: (*snes->ops->setfromoptions)(snes);
887: }
889: /* process any options handlers added with PetscObjectAddOptionsHandler() */
890: PetscObjectProcessOptionsHandlers((PetscObject)snes);
891: PetscOptionsEnd();
893: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
894: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
895: KSPSetFromOptions(snes->ksp);
897: if (!snes->linesearch) {
898: SNESGetLineSearch(snes, &snes->linesearch);
899: }
900: SNESLineSearchSetFromOptions(snes->linesearch);
902: /* if someone has set the SNES NPC type, create it. */
903: SNESGetOptionsPrefix(snes, &optionsprefix);
904: PetscOptionsHasName(optionsprefix, "-npc_snes_type", &pcset);
905: if (pcset && (!snes->pc)) {
906: SNESGetNPC(snes, &snes->pc);
907: }
908: return(0);
909: }
913: /*@C
914: SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
915: the nonlinear solvers.
917: Logically Collective on SNES
919: Input Parameters:
920: + snes - the SNES context
921: . compute - function to compute the context
922: - destroy - function to destroy the context
924: Level: intermediate
926: Notes:
927: This function is currently not available from Fortran.
929: .keywords: SNES, nonlinear, set, application, context
931: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
932: @*/
933: PetscErrorCode SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
934: {
937: snes->ops->usercompute = compute;
938: snes->ops->userdestroy = destroy;
939: return(0);
940: }
944: /*@
945: SNESSetApplicationContext - Sets the optional user-defined context for
946: the nonlinear solvers.
948: Logically Collective on SNES
950: Input Parameters:
951: + snes - the SNES context
952: - usrP - optional user context
954: Level: intermediate
956: .keywords: SNES, nonlinear, set, application, context
958: .seealso: SNESGetApplicationContext()
959: @*/
960: PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP)
961: {
963: KSP ksp;
967: SNESGetKSP(snes,&ksp);
968: KSPSetApplicationContext(ksp,usrP);
969: snes->user = usrP;
970: return(0);
971: }
975: /*@
976: SNESGetApplicationContext - Gets the user-defined context for the
977: nonlinear solvers.
979: Not Collective
981: Input Parameter:
982: . snes - SNES context
984: Output Parameter:
985: . usrP - user context
987: Level: intermediate
989: .keywords: SNES, nonlinear, get, application, context
991: .seealso: SNESSetApplicationContext()
992: @*/
993: PetscErrorCode SNESGetApplicationContext(SNES snes,void *usrP)
994: {
997: *(void**)usrP = snes->user;
998: return(0);
999: }
1003: /*@
1004: SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1005: at this time.
1007: Not Collective
1009: Input Parameter:
1010: . snes - SNES context
1012: Output Parameter:
1013: . iter - iteration number
1015: Notes:
1016: For example, during the computation of iteration 2 this would return 1.
1018: This is useful for using lagged Jacobians (where one does not recompute the
1019: Jacobian at each SNES iteration). For example, the code
1020: .vb
1021: SNESGetIterationNumber(snes,&it);
1022: if (!(it % 2)) {
1023: [compute Jacobian here]
1024: }
1025: .ve
1026: can be used in your ComputeJacobian() function to cause the Jacobian to be
1027: recomputed every second SNES iteration.
1029: Level: intermediate
1031: .keywords: SNES, nonlinear, get, iteration, number,
1033: .seealso: SNESGetLinearSolveIterations()
1034: @*/
1035: PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt *iter)
1036: {
1040: *iter = snes->iter;
1041: return(0);
1042: }
1046: /*@
1047: SNESSetIterationNumber - Sets the current iteration number.
1049: Not Collective
1051: Input Parameter:
1052: . snes - SNES context
1053: . iter - iteration number
1055: Level: developer
1057: .keywords: SNES, nonlinear, set, iteration, number,
1059: .seealso: SNESGetLinearSolveIterations()
1060: @*/
1061: PetscErrorCode SNESSetIterationNumber(SNES snes,PetscInt iter)
1062: {
1067: PetscObjectSAWsTakeAccess((PetscObject)snes);
1068: snes->iter = iter;
1069: PetscObjectSAWsGrantAccess((PetscObject)snes);
1070: return(0);
1071: }
1075: /*@
1076: SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1077: attempted by the nonlinear solver.
1079: Not Collective
1081: Input Parameter:
1082: . snes - SNES context
1084: Output Parameter:
1085: . nfails - number of unsuccessful steps attempted
1087: Notes:
1088: This counter is reset to zero for each successive call to SNESSolve().
1090: Level: intermediate
1092: .keywords: SNES, nonlinear, get, number, unsuccessful, steps
1094: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1095: SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1096: @*/
1097: PetscErrorCode SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1098: {
1102: *nfails = snes->numFailures;
1103: return(0);
1104: }
1108: /*@
1109: SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1110: attempted by the nonlinear solver before it gives up.
1112: Not Collective
1114: Input Parameters:
1115: + snes - SNES context
1116: - maxFails - maximum of unsuccessful steps
1118: Level: intermediate
1120: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
1122: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1123: SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1124: @*/
1125: PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1126: {
1129: snes->maxFailures = maxFails;
1130: return(0);
1131: }
1135: /*@
1136: SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1137: attempted by the nonlinear solver before it gives up.
1139: Not Collective
1141: Input Parameter:
1142: . snes - SNES context
1144: Output Parameter:
1145: . maxFails - maximum of unsuccessful steps
1147: Level: intermediate
1149: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1151: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1152: SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1154: @*/
1155: PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1156: {
1160: *maxFails = snes->maxFailures;
1161: return(0);
1162: }
1166: /*@
1167: SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1168: done by SNES.
1170: Not Collective
1172: Input Parameter:
1173: . snes - SNES context
1175: Output Parameter:
1176: . nfuncs - number of evaluations
1178: Level: intermediate
1180: Notes: Reset every time SNESSolve is called unless SNESSetCountersReset() is used.
1182: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1184: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1185: @*/
1186: PetscErrorCode SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1187: {
1191: *nfuncs = snes->nfuncs;
1192: return(0);
1193: }
1197: /*@
1198: SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1199: linear solvers.
1201: Not Collective
1203: Input Parameter:
1204: . snes - SNES context
1206: Output Parameter:
1207: . nfails - number of failed solves
1209: Notes:
1210: This counter is reset to zero for each successive call to SNESSolve().
1212: Level: intermediate
1214: .keywords: SNES, nonlinear, get, number, unsuccessful, steps
1216: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1217: @*/
1218: PetscErrorCode SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1219: {
1223: *nfails = snes->numLinearSolveFailures;
1224: return(0);
1225: }
1229: /*@
1230: SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1231: allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
1233: Logically Collective on SNES
1235: Input Parameters:
1236: + snes - SNES context
1237: - maxFails - maximum allowed linear solve failures
1239: Level: intermediate
1241: Notes: By default this is 0; that is SNES returns on the first failed linear solve
1243: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
1245: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1246: @*/
1247: PetscErrorCode SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1248: {
1252: snes->maxLinearSolveFailures = maxFails;
1253: return(0);
1254: }
1258: /*@
1259: SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1260: are allowed before SNES terminates
1262: Not Collective
1264: Input Parameter:
1265: . snes - SNES context
1267: Output Parameter:
1268: . maxFails - maximum of unsuccessful solves allowed
1270: Level: intermediate
1272: Notes: By default this is 1; that is SNES returns on the first failed linear solve
1274: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1276: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1277: @*/
1278: PetscErrorCode SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1279: {
1283: *maxFails = snes->maxLinearSolveFailures;
1284: return(0);
1285: }
1289: /*@
1290: SNESGetLinearSolveIterations - Gets the total number of linear iterations
1291: used by the nonlinear solver.
1293: Not Collective
1295: Input Parameter:
1296: . snes - SNES context
1298: Output Parameter:
1299: . lits - number of linear iterations
1301: Notes:
1302: This counter is reset to zero for each successive call to SNESSolve() unless SNESSetCountersReset() is used.
1304: Level: intermediate
1306: .keywords: SNES, nonlinear, get, number, linear, iterations
1308: .seealso: SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1309: @*/
1310: PetscErrorCode SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1311: {
1315: *lits = snes->linear_its;
1316: return(0);
1317: }
1321: /*@
1322: SNESSetCountersReset - Sets whether or not the counters for linear iterations and function evaluations
1323: are reset every time SNESSolve() is called.
1325: Logically Collective on SNES
1327: Input Parameter:
1328: + snes - SNES context
1329: - reset - whether to reset the counters or not
1331: Notes:
1332: This is automatically called with FALSE
1334: Level: developer
1336: .keywords: SNES, nonlinear, set, reset, number, linear, iterations
1338: .seealso: SNESGetNumberFunctionEvals(), SNESGetNumberLinearSolveIterations(), SNESGetNPC()
1339: @*/
1340: PetscErrorCode SNESSetCountersReset(SNES snes,PetscBool reset)
1341: {
1345: snes->counters_reset = reset;
1346: return(0);
1347: }
1352: /*@
1353: SNESSetKSP - Sets a KSP context for the SNES object to use
1355: Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
1357: Input Parameters:
1358: + snes - the SNES context
1359: - ksp - the KSP context
1361: Notes:
1362: The SNES object already has its KSP object, you can obtain with SNESGetKSP()
1363: so this routine is rarely needed.
1365: The KSP object that is already in the SNES object has its reference count
1366: decreased by one.
1368: Level: developer
1370: .keywords: SNES, nonlinear, get, KSP, context
1372: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1373: @*/
1374: PetscErrorCode SNESSetKSP(SNES snes,KSP ksp)
1375: {
1382: PetscObjectReference((PetscObject)ksp);
1383: if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1384: snes->ksp = ksp;
1385: return(0);
1386: }
1388: /* -----------------------------------------------------------*/
1391: /*@
1392: SNESCreate - Creates a nonlinear solver context.
1394: Collective on MPI_Comm
1396: Input Parameters:
1397: . comm - MPI communicator
1399: Output Parameter:
1400: . outsnes - the new SNES context
1402: Options Database Keys:
1403: + -snes_mf - Activates default matrix-free Jacobian-vector products,
1404: and no preconditioning matrix
1405: . -snes_mf_operator - Activates default matrix-free Jacobian-vector
1406: products, and a user-provided preconditioning matrix
1407: as set by SNESSetJacobian()
1408: - -snes_fd - Uses (slow!) finite differences to compute Jacobian
1410: Level: beginner
1412: .keywords: SNES, nonlinear, create, context
1414: .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()
1416: @*/
1417: PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes)
1418: {
1420: SNES snes;
1421: SNESKSPEW *kctx;
1425: *outsnes = NULL;
1426: SNESInitializePackage();
1428: PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_CLASSID,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);
1430: snes->ops->converged = SNESConvergedDefault;
1431: snes->usesksp = PETSC_TRUE;
1432: snes->tolerancesset = PETSC_FALSE;
1433: snes->max_its = 50;
1434: snes->max_funcs = 10000;
1435: snes->norm = 0.0;
1436: snes->normschedule = SNES_NORM_ALWAYS;
1437: snes->functype = SNES_FUNCTION_DEFAULT;
1438: #if defined(PETSC_USE_REAL_SINGLE)
1439: snes->rtol = 1.e-5;
1440: #else
1441: snes->rtol = 1.e-8;
1442: #endif
1443: snes->ttol = 0.0;
1444: #if defined(PETSC_USE_REAL_SINGLE)
1445: snes->abstol = 1.e-25;
1446: #else
1447: snes->abstol = 1.e-50;
1448: #endif
1449: snes->stol = 1.e-8;
1450: #if defined(PETSC_USE_REAL_SINGLE)
1451: snes->deltatol = 1.e-6;
1452: #else
1453: snes->deltatol = 1.e-12;
1454: #endif
1455: snes->nfuncs = 0;
1456: snes->numFailures = 0;
1457: snes->maxFailures = 1;
1458: snes->linear_its = 0;
1459: snes->lagjacobian = 1;
1460: snes->jac_iter = 0;
1461: snes->lagjac_persist = PETSC_FALSE;
1462: snes->lagpreconditioner = 1;
1463: snes->pre_iter = 0;
1464: snes->lagpre_persist = PETSC_FALSE;
1465: snes->numbermonitors = 0;
1466: snes->data = 0;
1467: snes->setupcalled = PETSC_FALSE;
1468: snes->ksp_ewconv = PETSC_FALSE;
1469: snes->nwork = 0;
1470: snes->work = 0;
1471: snes->nvwork = 0;
1472: snes->vwork = 0;
1473: snes->conv_hist_len = 0;
1474: snes->conv_hist_max = 0;
1475: snes->conv_hist = NULL;
1476: snes->conv_hist_its = NULL;
1477: snes->conv_hist_reset = PETSC_TRUE;
1478: snes->counters_reset = PETSC_TRUE;
1479: snes->vec_func_init_set = PETSC_FALSE;
1480: snes->reason = SNES_CONVERGED_ITERATING;
1481: snes->pcside = PC_RIGHT;
1483: snes->mf = PETSC_FALSE;
1484: snes->mf_operator = PETSC_FALSE;
1485: snes->mf_version = 1;
1487: snes->numLinearSolveFailures = 0;
1488: snes->maxLinearSolveFailures = 1;
1490: /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
1491: PetscNewLog(snes,&kctx);
1493: snes->kspconvctx = (void*)kctx;
1494: kctx->version = 2;
1495: kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1496: this was too large for some test cases */
1497: kctx->rtol_last = 0.0;
1498: kctx->rtol_max = .9;
1499: kctx->gamma = 1.0;
1500: kctx->alpha = .5*(1.0 + PetscSqrtReal(5.0));
1501: kctx->alpha2 = kctx->alpha;
1502: kctx->threshold = .1;
1503: kctx->lresid_last = 0.0;
1504: kctx->norm_last = 0.0;
1506: *outsnes = snes;
1507: return(0);
1508: }
1510: /*MC
1511: SNESFunction - functional form used to convey the nonlinear function to be solved by SNES
1513: Synopsis:
1514: #include <petscsnes.h>
1515: SNESFunction(SNES snes,Vec x,Vec f,void *ctx);
1517: Input Parameters:
1518: + snes - the SNES context
1519: . x - state at which to evaluate residual
1520: - ctx - optional user-defined function context, passed in with SNESSetFunction()
1522: Output Parameter:
1523: . f - vector to put residual (function value)
1525: Level: intermediate
1527: .seealso: SNESSetFunction(), SNESGetFunction()
1528: M*/
1532: /*@C
1533: SNESSetFunction - Sets the function evaluation routine and function
1534: vector for use by the SNES routines in solving systems of nonlinear
1535: equations.
1537: Logically Collective on SNES
1539: Input Parameters:
1540: + snes - the SNES context
1541: . r - vector to store function value
1542: . f - function evaluation routine; see SNESFunction for calling sequence details
1543: - ctx - [optional] user-defined context for private data for the
1544: function evaluation routine (may be NULL)
1546: Notes:
1547: The Newton-like methods typically solve linear systems of the form
1548: $ f'(x) x = -f(x),
1549: where f'(x) denotes the Jacobian matrix and f(x) is the function.
1551: Level: beginner
1553: .keywords: SNES, nonlinear, set, function
1555: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1556: @*/
1557: PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1558: {
1560: DM dm;
1564: if (r) {
1567: PetscObjectReference((PetscObject)r);
1568: VecDestroy(&snes->vec_func);
1570: snes->vec_func = r;
1571: }
1572: SNESGetDM(snes,&dm);
1573: DMSNESSetFunction(dm,f,ctx);
1574: return(0);
1575: }
1580: /*@C
1581: SNESSetInitialFunction - Sets the function vector to be used as the
1582: function norm at the initialization of the method. In some
1583: instances, the user has precomputed the function before calling
1584: SNESSolve. This function allows one to avoid a redundant call
1585: to SNESComputeFunction in that case.
1587: Logically Collective on SNES
1589: Input Parameters:
1590: + snes - the SNES context
1591: - f - vector to store function value
1593: Notes:
1594: This should not be modified during the solution procedure.
1596: This is used extensively in the SNESFAS hierarchy and in nonlinear preconditioning.
1598: Level: developer
1600: .keywords: SNES, nonlinear, set, function
1602: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1603: @*/
1604: PetscErrorCode SNESSetInitialFunction(SNES snes, Vec f)
1605: {
1607: Vec vec_func;
1613: if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1614: snes->vec_func_init_set = PETSC_FALSE;
1615: return(0);
1616: }
1617: SNESGetFunction(snes,&vec_func,NULL,NULL);
1618: VecCopy(f, vec_func);
1620: snes->vec_func_init_set = PETSC_TRUE;
1621: return(0);
1622: }
1626: /*@
1627: SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1628: of the SNES method.
1630: Logically Collective on SNES
1632: Input Parameters:
1633: + snes - the SNES context
1634: - normschedule - the frequency of norm computation
1636: Notes:
1637: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1638: of the nonlinear function and the taking of its norm at every iteration to
1639: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1640: (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
1641: may either be monitored for convergence or not. As these are often used as nonlinear
1642: preconditioners, monitoring the norm of their error is not a useful enterprise within
1643: their solution.
1645: Level: developer
1647: .keywords: SNES, nonlinear, set, function, norm, type
1649: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1650: @*/
1651: PetscErrorCode SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1652: {
1655: snes->normschedule = normschedule;
1656: return(0);
1657: }
1662: /*@
1663: SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1664: of the SNES method.
1666: Logically Collective on SNES
1668: Input Parameters:
1669: + snes - the SNES context
1670: - normschedule - the type of the norm used
1672: Level: advanced
1674: .keywords: SNES, nonlinear, set, function, norm, type
1676: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1677: @*/
1678: PetscErrorCode SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1679: {
1682: *normschedule = snes->normschedule;
1683: return(0);
1684: }
1689: /*@C
1690: SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
1691: of the SNES method.
1693: Logically Collective on SNES
1695: Input Parameters:
1696: + snes - the SNES context
1697: - normschedule - the frequency of norm computation
1699: Notes:
1700: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1701: of the nonlinear function and the taking of its norm at every iteration to
1702: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1703: (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
1704: may either be monitored for convergence or not. As these are often used as nonlinear
1705: preconditioners, monitoring the norm of their error is not a useful enterprise within
1706: their solution.
1708: Level: developer
1710: .keywords: SNES, nonlinear, set, function, norm, type
1712: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1713: @*/
1714: PetscErrorCode SNESSetFunctionType(SNES snes, SNESFunctionType type)
1715: {
1718: snes->functype = type;
1719: return(0);
1720: }
1725: /*@C
1726: SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
1727: of the SNES method.
1729: Logically Collective on SNES
1731: Input Parameters:
1732: + snes - the SNES context
1733: - normschedule - the type of the norm used
1735: Level: advanced
1737: .keywords: SNES, nonlinear, set, function, norm, type
1739: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1740: @*/
1741: PetscErrorCode SNESGetFunctionType(SNES snes, SNESFunctionType *type)
1742: {
1745: *type = snes->functype;
1746: return(0);
1747: }
1749: /*MC
1750: SNESNGSFunction - function used to convey a Gauss-Seidel sweep on the nonlinear function
1752: Synopsis:
1753: #include <petscsnes.h>
1754: $ SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);
1756: + X - solution vector
1757: . B - RHS vector
1758: - ctx - optional user-defined Gauss-Seidel context
1760: Level: intermediate
1762: .seealso: SNESSetNGS(), SNESGetNGS()
1763: M*/
1767: /*@C
1768: SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
1769: use with composed nonlinear solvers.
1771: Input Parameters:
1772: + snes - the SNES context
1773: . f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
1774: - ctx - [optional] user-defined context for private data for the
1775: smoother evaluation routine (may be NULL)
1777: Notes:
1778: The NGS routines are used by the composed nonlinear solver to generate
1779: a problem appropriate update to the solution, particularly FAS.
1781: Level: intermediate
1783: .keywords: SNES, nonlinear, set, Gauss-Seidel
1785: .seealso: SNESGetFunction(), SNESComputeNGS()
1786: @*/
1787: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1788: {
1790: DM dm;
1794: SNESGetDM(snes,&dm);
1795: DMSNESSetNGS(dm,f,ctx);
1796: return(0);
1797: }
1801: PETSC_EXTERN PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
1802: {
1804: DM dm;
1805: DMSNES sdm;
1808: SNESGetDM(snes,&dm);
1809: DMGetDMSNES(dm,&sdm);
1810: /* A(x)*x - b(x) */
1811: if (sdm->ops->computepfunction) {
1812: (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
1813: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");
1815: if (sdm->ops->computepjacobian) {
1816: (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
1817: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard matrix.");
1818: VecScale(f,-1.0);
1819: MatMultAdd(snes->jacobian,x,f,f);
1820: return(0);
1821: }
1825: PETSC_EXTERN PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
1826: {
1828: /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
1829: return(0);
1830: }
1834: /*@C
1835: SNESSetPicard - Use SNES to solve the semilinear-system A(x) x = b(x) via a Picard type iteration (Picard linearization)
1837: Logically Collective on SNES
1839: Input Parameters:
1840: + snes - the SNES context
1841: . r - vector to store function value
1842: . b - function evaluation routine
1843: . Amat - matrix with which A(x) x - b(x) is to be computed
1844: . Pmat - matrix from which preconditioner is computed (usually the same as Amat)
1845: . J - function to compute matrix value
1846: - ctx - [optional] user-defined context for private data for the
1847: function evaluation routine (may be NULL)
1849: Notes:
1850: We do not recomemend using this routine. It is far better to provide the nonlinear function F() and some approximation to the Jacobian and use
1851: 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.
1853: One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both
1855: $ Solves the equation A(x) x = b(x) via the defect correction algorithm A(x^{n}) (x^{n+1} - x^{n}) = b(x^{n}) - A(x^{n})x^{n}
1856: $ Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration.
1858: Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner.
1860: We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then
1861: the direct Picard iteration A(x^n) x^{n+1} = b(x^n)
1863: 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
1864: 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
1865: different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).
1867: Level: intermediate
1869: .keywords: SNES, nonlinear, set, function
1871: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard()
1872: @*/
1873: PetscErrorCode SNESSetPicard(SNES snes,Vec r,PetscErrorCode (*b)(SNES,Vec,Vec,void*),Mat Amat, Mat Pmat, PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
1874: {
1876: DM dm;
1880: SNESGetDM(snes, &dm);
1881: DMSNESSetPicard(dm,b,J,ctx);
1882: SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
1883: SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
1884: return(0);
1885: }
1889: /*@C
1890: SNESGetPicard - Returns the context for the Picard iteration
1892: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
1894: Input Parameter:
1895: . snes - the SNES context
1897: Output Parameter:
1898: + r - the function (or NULL)
1899: . f - the function (or NULL); see SNESFunction for calling sequence details
1900: . Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
1901: . Pmat - the matrix from which the preconditioner will be constructed (or NULL)
1902: . J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
1903: - ctx - the function context (or NULL)
1905: Level: advanced
1907: .keywords: SNES, nonlinear, get, function
1909: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
1910: @*/
1911: 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)
1912: {
1914: DM dm;
1918: SNESGetFunction(snes,r,NULL,NULL);
1919: SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
1920: SNESGetDM(snes,&dm);
1921: DMSNESGetPicard(dm,f,J,ctx);
1922: return(0);
1923: }
1927: /*@C
1928: SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem
1930: Logically Collective on SNES
1932: Input Parameters:
1933: + snes - the SNES context
1934: . func - function evaluation routine
1935: - ctx - [optional] user-defined context for private data for the
1936: function evaluation routine (may be NULL)
1938: Calling sequence of func:
1939: $ func (SNES snes,Vec x,void *ctx);
1941: . f - function vector
1942: - ctx - optional user-defined function context
1944: Level: intermediate
1946: .keywords: SNES, nonlinear, set, function
1948: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
1949: @*/
1950: PetscErrorCode SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
1951: {
1954: if (func) snes->ops->computeinitialguess = func;
1955: if (ctx) snes->initialguessP = ctx;
1956: return(0);
1957: }
1959: /* --------------------------------------------------------------- */
1962: /*@C
1963: SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
1964: it assumes a zero right hand side.
1966: Logically Collective on SNES
1968: Input Parameter:
1969: . snes - the SNES context
1971: Output Parameter:
1972: . rhs - the right hand side vector or NULL if the right hand side vector is null
1974: Level: intermediate
1976: .keywords: SNES, nonlinear, get, function, right hand side
1978: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
1979: @*/
1980: PetscErrorCode SNESGetRhs(SNES snes,Vec *rhs)
1981: {
1985: *rhs = snes->vec_rhs;
1986: return(0);
1987: }
1991: /*@
1992: SNESComputeFunction - Calls the function that has been set with SNESSetFunction().
1994: Collective on SNES
1996: Input Parameters:
1997: + snes - the SNES context
1998: - x - input vector
2000: Output Parameter:
2001: . y - function vector, as set by SNESSetFunction()
2003: Notes:
2004: SNESComputeFunction() is typically used within nonlinear solvers
2005: implementations, so most users would not generally call this routine
2006: themselves.
2008: Level: developer
2010: .keywords: SNES, nonlinear, compute, function
2012: .seealso: SNESSetFunction(), SNESGetFunction()
2013: @*/
2014: PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y)
2015: {
2017: DM dm;
2018: DMSNES sdm;
2026: VecValidValues(x,2,PETSC_TRUE);
2028: SNESGetDM(snes,&dm);
2029: DMGetDMSNES(dm,&sdm);
2030: if (sdm->ops->computefunction) {
2031: PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2032: PetscStackPush("SNES user function");
2033: (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2034: PetscStackPop;
2035: PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2036: } else if (snes->vec_rhs) {
2037: MatMult(snes->jacobian, x, y);
2038: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2039: if (snes->vec_rhs) {
2040: VecAXPY(y,-1.0,snes->vec_rhs);
2041: }
2042: snes->nfuncs++;
2043: VecValidValues(y,3,PETSC_FALSE);
2044: return(0);
2045: }
2049: /*@
2050: SNESComputeNGS - Calls the Gauss-Seidel function that has been set with SNESSetNGS().
2052: Collective on SNES
2054: Input Parameters:
2055: + snes - the SNES context
2056: . x - input vector
2057: - b - rhs vector
2059: Output Parameter:
2060: . x - new solution vector
2062: Notes:
2063: SNESComputeNGS() is typically used within composed nonlinear solver
2064: implementations, so most users would not generally call this routine
2065: themselves.
2067: Level: developer
2069: .keywords: SNES, nonlinear, compute, function
2071: .seealso: SNESSetNGS(), SNESComputeFunction()
2072: @*/
2073: PetscErrorCode SNESComputeNGS(SNES snes,Vec b,Vec x)
2074: {
2076: DM dm;
2077: DMSNES sdm;
2085: if (b) {VecValidValues(b,2,PETSC_TRUE);}
2086: PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2087: SNESGetDM(snes,&dm);
2088: DMGetDMSNES(dm,&sdm);
2089: if (sdm->ops->computegs) {
2090: PetscStackPush("SNES user NGS");
2091: (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2092: PetscStackPop;
2093: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2094: PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2095: VecValidValues(x,3,PETSC_FALSE);
2096: return(0);
2097: }
2101: /*@
2102: SNESComputeJacobian - Computes the Jacobian matrix that has been set with SNESSetJacobian().
2104: Collective on SNES and Mat
2106: Input Parameters:
2107: + snes - the SNES context
2108: - x - input vector
2110: Output Parameters:
2111: + A - Jacobian matrix
2112: - B - optional preconditioning matrix
2114: Options Database Keys:
2115: + -snes_lag_preconditioner <lag>
2116: . -snes_lag_jacobian <lag>
2117: . -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2118: . -snes_compare_explicit_draw - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2119: . -snes_compare_explicit_contour - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2120: . -snes_compare_operator - Make the comparison options above use the operator instead of the preconditioning matrix
2121: . -snes_compare_coloring - Compute the finite differece Jacobian using coloring and display norms of difference
2122: . -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
2123: . -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2124: . -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2125: . -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2126: . -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
2127: - -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences
2130: Notes:
2131: Most users should not need to explicitly call this routine, as it
2132: is used internally within the nonlinear solvers.
2134: See KSPSetOperators() for important information about setting the
2135: flag parameter.
2137: Level: developer
2139: .keywords: SNES, compute, Jacobian, matrix
2141: .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2142: @*/
2143: PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2144: {
2146: PetscBool flag;
2147: DM dm;
2148: DMSNES sdm;
2149: KSP ksp;
2155: VecValidValues(X,2,PETSC_TRUE);
2156: SNESGetDM(snes,&dm);
2157: DMGetDMSNES(dm,&sdm);
2159: if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_USER,"Must call SNESSetJacobian(), DMSNESSetJacobian(), DMDASNESSetJacobianLocal(), etc");
2161: /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
2163: if (snes->lagjacobian == -2) {
2164: snes->lagjacobian = -1;
2166: PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2167: } else if (snes->lagjacobian == -1) {
2168: PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2169: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2170: if (flag) {
2171: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2172: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2173: }
2174: return(0);
2175: } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2176: PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2177: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2178: if (flag) {
2179: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2180: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2181: }
2182: return(0);
2183: }
2184: if (snes->pc && snes->pcside == PC_LEFT) {
2185: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2186: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2187: return(0);
2188: }
2190: PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2192: PetscStackPush("SNES user Jacobian function");
2193: (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2194: PetscStackPop;
2196: PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);
2198: /* the next line ensures that snes->ksp exists */
2199: SNESGetKSP(snes,&ksp);
2200: if (snes->lagpreconditioner == -2) {
2201: PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2202: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2203: snes->lagpreconditioner = -1;
2204: } else if (snes->lagpreconditioner == -1) {
2205: PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2206: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2207: } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2208: PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2209: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2210: } else {
2211: PetscInfo(snes,"Rebuilding preconditioner\n");
2212: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2213: }
2215: /* make sure user returned a correct Jacobian and preconditioner */
2218: {
2219: PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
2220: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit",&flag,NULL);
2221: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",&flag_draw,NULL);
2222: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",&flag_contour,NULL);
2223: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_operator",&flag_operator,NULL);
2224: if (flag || flag_draw || flag_contour) {
2225: Mat Bexp_mine = NULL,Bexp,FDexp;
2226: PetscViewer vdraw,vstdout;
2227: PetscBool flg;
2228: if (flag_operator) {
2229: MatComputeExplicitOperator(A,&Bexp_mine);
2230: Bexp = Bexp_mine;
2231: } else {
2232: /* See if the preconditioning matrix can be viewed and added directly */
2233: PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2234: if (flg) Bexp = B;
2235: else {
2236: /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
2237: MatComputeExplicitOperator(B,&Bexp_mine);
2238: Bexp = Bexp_mine;
2239: }
2240: }
2241: MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);
2242: SNESComputeJacobianDefault(snes,X,FDexp,FDexp,NULL);
2243: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2244: if (flag_draw || flag_contour) {
2245: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2246: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2247: } else vdraw = NULL;
2248: PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator ? "Jacobian" : "preconditioning Jacobian");
2249: if (flag) {MatView(Bexp,vstdout);}
2250: if (vdraw) {MatView(Bexp,vdraw);}
2251: PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");
2252: if (flag) {MatView(FDexp,vstdout);}
2253: if (vdraw) {MatView(FDexp,vdraw);}
2254: MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);
2255: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");
2256: if (flag) {MatView(FDexp,vstdout);}
2257: if (vdraw) { /* Always use contour for the difference */
2258: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2259: MatView(FDexp,vdraw);
2260: PetscViewerPopFormat(vdraw);
2261: }
2262: if (flag_contour) {PetscViewerPopFormat(vdraw);}
2263: PetscViewerDestroy(&vdraw);
2264: MatDestroy(&Bexp_mine);
2265: MatDestroy(&FDexp);
2266: }
2267: }
2268: {
2269: PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
2270: PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
2271: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring",&flag,NULL);
2272: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_display",&flag_display,NULL);
2273: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",&flag_draw,NULL);
2274: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",&flag_contour,NULL);
2275: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",&flag_threshold,NULL);
2276: PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,NULL);
2277: PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,NULL);
2278: if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
2279: Mat Bfd;
2280: PetscViewer vdraw,vstdout;
2281: MatColoring coloring;
2282: ISColoring iscoloring;
2283: MatFDColoring matfdcoloring;
2284: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2285: void *funcctx;
2286: PetscReal norm1,norm2,normmax;
2288: MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2289: MatColoringCreate(Bfd,&coloring);
2290: MatColoringSetType(coloring,MATCOLORINGSL);
2291: MatColoringSetFromOptions(coloring);
2292: MatColoringApply(coloring,&iscoloring);
2293: MatColoringDestroy(&coloring);
2294: MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2295: MatFDColoringSetFromOptions(matfdcoloring);
2296: MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2297: ISColoringDestroy(&iscoloring);
2299: /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2300: SNESGetFunction(snes,NULL,&func,&funcctx);
2301: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2302: PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2303: PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2304: MatFDColoringSetFromOptions(matfdcoloring);
2305: MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2306: MatFDColoringDestroy(&matfdcoloring);
2308: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2309: if (flag_draw || flag_contour) {
2310: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2311: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2312: } else vdraw = NULL;
2313: PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2314: if (flag_display) {MatView(B,vstdout);}
2315: if (vdraw) {MatView(B,vdraw);}
2316: PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2317: if (flag_display) {MatView(Bfd,vstdout);}
2318: if (vdraw) {MatView(Bfd,vdraw);}
2319: MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2320: MatNorm(Bfd,NORM_1,&norm1);
2321: MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2322: MatNorm(Bfd,NORM_MAX,&normmax);
2323: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2324: if (flag_display) {MatView(Bfd,vstdout);}
2325: if (vdraw) { /* Always use contour for the difference */
2326: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2327: MatView(Bfd,vdraw);
2328: PetscViewerPopFormat(vdraw);
2329: }
2330: if (flag_contour) {PetscViewerPopFormat(vdraw);}
2332: if (flag_threshold) {
2333: PetscInt bs,rstart,rend,i;
2334: MatGetBlockSize(B,&bs);
2335: MatGetOwnershipRange(B,&rstart,&rend);
2336: for (i=rstart; i<rend; i++) {
2337: const PetscScalar *ba,*ca;
2338: const PetscInt *bj,*cj;
2339: PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
2340: PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0;
2341: MatGetRow(B,i,&bn,&bj,&ba);
2342: MatGetRow(Bfd,i,&cn,&cj,&ca);
2343: if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
2344: for (j=0; j<bn; j++) {
2345: PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2346: if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
2347: maxentrycol = bj[j];
2348: maxentry = PetscRealPart(ba[j]);
2349: }
2350: if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
2351: maxdiffcol = bj[j];
2352: maxdiff = PetscRealPart(ca[j]);
2353: }
2354: if (rdiff > maxrdiff) {
2355: maxrdiffcol = bj[j];
2356: maxrdiff = rdiff;
2357: }
2358: }
2359: if (maxrdiff > 1) {
2360: 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);
2361: for (j=0; j<bn; j++) {
2362: PetscReal rdiff;
2363: rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2364: if (rdiff > 1) {
2365: PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
2366: }
2367: }
2368: PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
2369: }
2370: MatRestoreRow(B,i,&bn,&bj,&ba);
2371: MatRestoreRow(Bfd,i,&cn,&cj,&ca);
2372: }
2373: }
2374: PetscViewerDestroy(&vdraw);
2375: MatDestroy(&Bfd);
2376: }
2377: }
2378: return(0);
2379: }
2381: /*MC
2382: SNESJacobianFunction - function used to convey the nonlinear Jacobian of the function to be solved by SNES
2384: Synopsis:
2385: #include <petscsnes.h>
2386: $ SNESJacobianFunction(SNES snes,Vec x,Mat Amat,Mat Pmat,void *ctx);
2388: + x - input vector
2389: . Amat - the matrix that defines the (approximate) Jacobian
2390: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2391: - ctx - [optional] user-defined Jacobian context
2393: Level: intermediate
2395: .seealso: SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2396: M*/
2400: /*@C
2401: SNESSetJacobian - Sets the function to compute Jacobian as well as the
2402: location to store the matrix.
2404: Logically Collective on SNES and Mat
2406: Input Parameters:
2407: + snes - the SNES context
2408: . Amat - the matrix that defines the (approximate) Jacobian
2409: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2410: . J - Jacobian evaluation routine (if NULL then SNES retains any previously set value)
2411: - ctx - [optional] user-defined context for private data for the
2412: Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)
2414: Notes:
2415: If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2416: each matrix.
2418: If using SNESComputeJacobianDefaultColor() to assemble a Jacobian, the ctx argument
2419: must be a MatFDColoring.
2421: Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian. One common
2422: example is to use the "Picard linearization" which only differentiates through the highest order parts of each term.
2424: Level: beginner
2426: .keywords: SNES, nonlinear, set, Jacobian, matrix
2428: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J, SNESSetPicard()
2429: @*/
2430: PetscErrorCode SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2431: {
2433: DM dm;
2441: SNESGetDM(snes,&dm);
2442: DMSNESSetJacobian(dm,J,ctx);
2443: if (Amat) {
2444: PetscObjectReference((PetscObject)Amat);
2445: MatDestroy(&snes->jacobian);
2447: snes->jacobian = Amat;
2448: }
2449: if (Pmat) {
2450: PetscObjectReference((PetscObject)Pmat);
2451: MatDestroy(&snes->jacobian_pre);
2453: snes->jacobian_pre = Pmat;
2454: }
2455: return(0);
2456: }
2460: /*@C
2461: SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2462: provided context for evaluating the Jacobian.
2464: Not Collective, but Mat object will be parallel if SNES object is
2466: Input Parameter:
2467: . snes - the nonlinear solver context
2469: Output Parameters:
2470: + Amat - location to stash (approximate) Jacobian matrix (or NULL)
2471: . Pmat - location to stash matrix used to compute the preconditioner (or NULL)
2472: . J - location to put Jacobian function (or NULL)
2473: - ctx - location to stash Jacobian ctx (or NULL)
2475: Level: advanced
2477: .seealso: SNESSetJacobian(), SNESComputeJacobian()
2478: @*/
2479: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2480: {
2482: DM dm;
2483: DMSNES sdm;
2487: if (Amat) *Amat = snes->jacobian;
2488: if (Pmat) *Pmat = snes->jacobian_pre;
2489: SNESGetDM(snes,&dm);
2490: DMGetDMSNES(dm,&sdm);
2491: if (J) *J = sdm->ops->computejacobian;
2492: if (ctx) *ctx = sdm->jacobianctx;
2493: return(0);
2494: }
2498: /*@
2499: SNESSetUp - Sets up the internal data structures for the later use
2500: of a nonlinear solver.
2502: Collective on SNES
2504: Input Parameters:
2505: . snes - the SNES context
2507: Notes:
2508: For basic use of the SNES solvers the user need not explicitly call
2509: SNESSetUp(), since these actions will automatically occur during
2510: the call to SNESSolve(). However, if one wishes to control this
2511: phase separately, SNESSetUp() should be called after SNESCreate()
2512: and optional routines of the form SNESSetXXX(), but before SNESSolve().
2514: Level: advanced
2516: .keywords: SNES, nonlinear, setup
2518: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
2519: @*/
2520: PetscErrorCode SNESSetUp(SNES snes)
2521: {
2523: DM dm;
2524: DMSNES sdm;
2525: SNESLineSearch linesearch, pclinesearch;
2526: void *lsprectx,*lspostctx;
2527: PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
2528: PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
2529: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2530: Vec f,fpc;
2531: void *funcctx;
2532: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2533: void *jacctx,*appctx;
2534: Mat j,jpre;
2538: if (snes->setupcalled) return(0);
2540: if (!((PetscObject)snes)->type_name) {
2541: SNESSetType(snes,SNESNEWTONLS);
2542: }
2544: SNESGetFunction(snes,&snes->vec_func,NULL,NULL);
2546: SNESGetDM(snes,&dm);
2547: DMGetDMSNES(dm,&sdm);
2548: if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
2549: if (!sdm->ops->computejacobian) {
2550: DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
2551: }
2552: if (!snes->vec_func) {
2553: DMCreateGlobalVector(dm,&snes->vec_func);
2554: }
2556: if (!snes->ksp) {
2557: SNESGetKSP(snes, &snes->ksp);
2558: }
2560: if (!snes->linesearch) {
2561: SNESGetLineSearch(snes, &snes->linesearch);
2562: }
2563: SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);
2565: if (snes->pc && (snes->pcside == PC_LEFT)) {
2566: snes->mf = PETSC_TRUE;
2567: snes->mf_operator = PETSC_FALSE;
2568: }
2570: if (snes->pc) {
2571: /* copy the DM over */
2572: SNESGetDM(snes,&dm);
2573: SNESSetDM(snes->pc,dm);
2575: SNESGetFunction(snes,&f,&func,&funcctx);
2576: VecDuplicate(f,&fpc);
2577: SNESSetFunction(snes->pc,fpc,func,funcctx);
2578: SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
2579: SNESSetJacobian(snes->pc,j,jpre,jac,jacctx);
2580: SNESGetApplicationContext(snes,&appctx);
2581: SNESSetApplicationContext(snes->pc,appctx);
2582: VecDestroy(&fpc);
2584: /* copy the function pointers over */
2585: PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)snes->pc);
2587: /* default to 1 iteration */
2588: SNESSetTolerances(snes->pc,0.0,0.0,0.0,1,snes->pc->max_funcs);
2589: if (snes->pcside==PC_RIGHT) {
2590: SNESSetNormSchedule(snes->pc,SNES_NORM_FINAL_ONLY);
2591: } else {
2592: SNESSetNormSchedule(snes->pc,SNES_NORM_NONE);
2593: }
2594: SNESSetFromOptions(snes->pc);
2596: /* copy the line search context over */
2597: SNESGetLineSearch(snes,&linesearch);
2598: SNESGetLineSearch(snes->pc,&pclinesearch);
2599: SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
2600: SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
2601: SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
2602: SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
2603: PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
2604: }
2605: if (snes->mf) {
2606: SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
2607: }
2608: if (snes->ops->usercompute && !snes->user) {
2609: (*snes->ops->usercompute)(snes,(void**)&snes->user);
2610: }
2612: snes->jac_iter = 0;
2613: snes->pre_iter = 0;
2615: if (snes->ops->setup) {
2616: (*snes->ops->setup)(snes);
2617: }
2619: if (snes->pc && (snes->pcside == PC_LEFT)) {
2620: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2621: SNESGetLineSearch(snes,&linesearch);
2622: SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
2623: }
2624: }
2626: snes->setupcalled = PETSC_TRUE;
2627: return(0);
2628: }
2632: /*@
2633: SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats
2635: Collective on SNES
2637: Input Parameter:
2638: . snes - iterative context obtained from SNESCreate()
2640: Level: intermediate
2642: Notes: Also calls the application context destroy routine set with SNESSetComputeApplicationContext()
2644: .keywords: SNES, destroy
2646: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
2647: @*/
2648: PetscErrorCode SNESReset(SNES snes)
2649: {
2654: if (snes->ops->userdestroy && snes->user) {
2655: (*snes->ops->userdestroy)((void**)&snes->user);
2656: snes->user = NULL;
2657: }
2658: if (snes->pc) {
2659: SNESReset(snes->pc);
2660: }
2662: if (snes->ops->reset) {
2663: (*snes->ops->reset)(snes);
2664: }
2665: if (snes->ksp) {
2666: KSPReset(snes->ksp);
2667: }
2669: if (snes->linesearch) {
2670: SNESLineSearchReset(snes->linesearch);
2671: }
2673: VecDestroy(&snes->vec_rhs);
2674: VecDestroy(&snes->vec_sol);
2675: VecDestroy(&snes->vec_sol_update);
2676: VecDestroy(&snes->vec_func);
2677: MatDestroy(&snes->jacobian);
2678: MatDestroy(&snes->jacobian_pre);
2679: VecDestroyVecs(snes->nwork,&snes->work);
2680: VecDestroyVecs(snes->nvwork,&snes->vwork);
2682: snes->nwork = snes->nvwork = 0;
2683: snes->setupcalled = PETSC_FALSE;
2684: return(0);
2685: }
2689: /*@
2690: SNESDestroy - Destroys the nonlinear solver context that was created
2691: with SNESCreate().
2693: Collective on SNES
2695: Input Parameter:
2696: . snes - the SNES context
2698: Level: beginner
2700: .keywords: SNES, nonlinear, destroy
2702: .seealso: SNESCreate(), SNESSolve()
2703: @*/
2704: PetscErrorCode SNESDestroy(SNES *snes)
2705: {
2709: if (!*snes) return(0);
2711: if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; return(0);}
2713: SNESReset((*snes));
2714: SNESDestroy(&(*snes)->pc);
2716: /* if memory was published with SAWs then destroy it */
2717: PetscObjectSAWsViewOff((PetscObject)*snes);
2718: if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}
2720: DMDestroy(&(*snes)->dm);
2721: KSPDestroy(&(*snes)->ksp);
2722: SNESLineSearchDestroy(&(*snes)->linesearch);
2724: PetscFree((*snes)->kspconvctx);
2725: if ((*snes)->ops->convergeddestroy) {
2726: (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
2727: }
2728: if ((*snes)->conv_malloc) {
2729: PetscFree((*snes)->conv_hist);
2730: PetscFree((*snes)->conv_hist_its);
2731: }
2732: SNESMonitorCancel((*snes));
2733: PetscHeaderDestroy(snes);
2734: return(0);
2735: }
2737: /* ----------- Routines to set solver parameters ---------- */
2741: /*@
2742: SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
2744: Logically Collective on SNES
2746: Input Parameters:
2747: + snes - the SNES context
2748: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2749: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
2751: Options Database Keys:
2752: . -snes_lag_preconditioner <lag>
2754: Notes:
2755: The default is 1
2756: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2757: If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
2759: Level: intermediate
2761: .keywords: SNES, nonlinear, set, convergence, tolerances
2763: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
2765: @*/
2766: PetscErrorCode SNESSetLagPreconditioner(SNES snes,PetscInt lag)
2767: {
2770: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2771: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2773: snes->lagpreconditioner = lag;
2774: return(0);
2775: }
2779: /*@
2780: SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does
2782: Logically Collective on SNES
2784: Input Parameters:
2785: + snes - the SNES context
2786: - steps - the number of refinements to do, defaults to 0
2788: Options Database Keys:
2789: . -snes_grid_sequence <steps>
2791: Level: intermediate
2793: Notes:
2794: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
2796: .keywords: SNES, nonlinear, set, convergence, tolerances
2798: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
2800: @*/
2801: PetscErrorCode SNESSetGridSequence(SNES snes,PetscInt steps)
2802: {
2806: snes->gridsequence = steps;
2807: return(0);
2808: }
2812: /*@
2813: SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
2815: Not Collective
2817: Input Parameter:
2818: . snes - the SNES context
2820: Output Parameter:
2821: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2822: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
2824: Options Database Keys:
2825: . -snes_lag_preconditioner <lag>
2827: Notes:
2828: The default is 1
2829: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2831: Level: intermediate
2833: .keywords: SNES, nonlinear, set, convergence, tolerances
2835: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
2837: @*/
2838: PetscErrorCode SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
2839: {
2842: *lag = snes->lagpreconditioner;
2843: return(0);
2844: }
2848: /*@
2849: SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
2850: often the preconditioner is rebuilt.
2852: Logically Collective on SNES
2854: Input Parameters:
2855: + snes - the SNES context
2856: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2857: the Jacobian is built etc. -2 means rebuild at next chance but then never again
2859: Options Database Keys:
2860: . -snes_lag_jacobian <lag>
2862: Notes:
2863: The default is 1
2864: The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2865: 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
2866: at the next Newton step but never again (unless it is reset to another value)
2868: Level: intermediate
2870: .keywords: SNES, nonlinear, set, convergence, tolerances
2872: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
2874: @*/
2875: PetscErrorCode SNESSetLagJacobian(SNES snes,PetscInt lag)
2876: {
2879: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2880: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2882: snes->lagjacobian = lag;
2883: return(0);
2884: }
2888: /*@
2889: SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
2891: Not Collective
2893: Input Parameter:
2894: . snes - the SNES context
2896: Output Parameter:
2897: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2898: the Jacobian is built etc.
2900: Options Database Keys:
2901: . -snes_lag_jacobian <lag>
2903: Notes:
2904: The default is 1
2905: The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2907: Level: intermediate
2909: .keywords: SNES, nonlinear, set, convergence, tolerances
2911: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
2913: @*/
2914: PetscErrorCode SNESGetLagJacobian(SNES snes,PetscInt *lag)
2915: {
2918: *lag = snes->lagjacobian;
2919: return(0);
2920: }
2924: /*@
2925: SNESSetLagJacobianPersists - Set whether or not the Jacobian lagging persists through multiple solves
2927: Logically collective on SNES
2929: Input Parameter:
2930: + snes - the SNES context
2931: - flg - jacobian lagging persists if true
2933: Options Database Keys:
2934: . -snes_lag_jacobian_persists <flg>
2936: Notes: This is useful both for nonlinear preconditioning, where it's appropriate to have the Jacobian be stale by
2937: several solves, and for implicit time-stepping, where Jacobian lagging in the inner nonlinear solve over several
2938: timesteps may present huge efficiency gains.
2940: Level: developer
2942: .keywords: SNES, nonlinear, lag
2944: .seealso: SNESSetLagPreconditionerPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()
2946: @*/
2947: PetscErrorCode SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
2948: {
2952: snes->lagjac_persist = flg;
2953: return(0);
2954: }
2958: /*@
2959: SNESSetLagPreconditionerPersists - Set whether or not the preconditioner lagging persists through multiple solves
2961: Logically Collective on SNES
2963: Input Parameter:
2964: + snes - the SNES context
2965: - flg - preconditioner lagging persists if true
2967: Options Database Keys:
2968: . -snes_lag_jacobian_persists <flg>
2970: Notes: This is useful both for nonlinear preconditioning, where it's appropriate to have the preconditioner be stale
2971: by several solves, and for implicit time-stepping, where preconditioner lagging in the inner nonlinear solve over
2972: several timesteps may present huge efficiency gains.
2974: Level: developer
2976: .keywords: SNES, nonlinear, lag
2978: .seealso: SNESSetLagJacobianPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()
2980: @*/
2981: PetscErrorCode SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
2982: {
2986: snes->lagpre_persist = flg;
2987: return(0);
2988: }
2992: /*@
2993: SNESSetTolerances - Sets various parameters used in convergence tests.
2995: Logically Collective on SNES
2997: Input Parameters:
2998: + snes - the SNES context
2999: . abstol - absolute convergence tolerance
3000: . rtol - relative convergence tolerance
3001: . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x ||
3002: . maxit - maximum number of iterations
3003: - maxf - maximum number of function evaluations
3005: Options Database Keys:
3006: + -snes_atol <abstol> - Sets abstol
3007: . -snes_rtol <rtol> - Sets rtol
3008: . -snes_stol <stol> - Sets stol
3009: . -snes_max_it <maxit> - Sets maxit
3010: - -snes_max_funcs <maxf> - Sets maxf
3012: Notes:
3013: The default maximum number of iterations is 50.
3014: The default maximum number of function evaluations is 1000.
3016: Level: intermediate
3018: .keywords: SNES, nonlinear, set, convergence, tolerances
3020: .seealso: SNESSetTrustRegionTolerance()
3021: @*/
3022: PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3023: {
3032: if (abstol != PETSC_DEFAULT) {
3033: if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3034: snes->abstol = abstol;
3035: }
3036: if (rtol != PETSC_DEFAULT) {
3037: if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %g must be non-negative and less than 1.0",(double)rtol);
3038: snes->rtol = rtol;
3039: }
3040: if (stol != PETSC_DEFAULT) {
3041: if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3042: snes->stol = stol;
3043: }
3044: if (maxit != PETSC_DEFAULT) {
3045: if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3046: snes->max_its = maxit;
3047: }
3048: if (maxf != PETSC_DEFAULT) {
3049: if (maxf < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
3050: snes->max_funcs = maxf;
3051: }
3052: snes->tolerancesset = PETSC_TRUE;
3053: return(0);
3054: }
3058: /*@
3059: SNESGetTolerances - Gets various parameters used in convergence tests.
3061: Not Collective
3063: Input Parameters:
3064: + snes - the SNES context
3065: . atol - absolute convergence tolerance
3066: . rtol - relative convergence tolerance
3067: . stol - convergence tolerance in terms of the norm
3068: of the change in the solution between steps
3069: . maxit - maximum number of iterations
3070: - maxf - maximum number of function evaluations
3072: Notes:
3073: The user can specify NULL for any parameter that is not needed.
3075: Level: intermediate
3077: .keywords: SNES, nonlinear, get, convergence, tolerances
3079: .seealso: SNESSetTolerances()
3080: @*/
3081: PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3082: {
3085: if (atol) *atol = snes->abstol;
3086: if (rtol) *rtol = snes->rtol;
3087: if (stol) *stol = snes->stol;
3088: if (maxit) *maxit = snes->max_its;
3089: if (maxf) *maxf = snes->max_funcs;
3090: return(0);
3091: }
3095: /*@
3096: SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
3098: Logically Collective on SNES
3100: Input Parameters:
3101: + snes - the SNES context
3102: - tol - tolerance
3104: Options Database Key:
3105: . -snes_trtol <tol> - Sets tol
3107: Level: intermediate
3109: .keywords: SNES, nonlinear, set, trust region, tolerance
3111: .seealso: SNESSetTolerances()
3112: @*/
3113: PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3114: {
3118: snes->deltatol = tol;
3119: return(0);
3120: }
3122: /*
3123: Duplicate the lg monitors for SNES from KSP; for some reason with
3124: dynamic libraries things don't work under Sun4 if we just use
3125: macros instead of functions
3126: */
3129: PetscErrorCode SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,void *ctx)
3130: {
3135: KSPMonitorLGResidualNorm((KSP)snes,it,norm,ctx);
3136: return(0);
3137: }
3141: PetscErrorCode SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
3142: {
3146: KSPMonitorLGResidualNormCreate(host,label,x,y,m,n,draw);
3147: return(0);
3148: }
3152: PetscErrorCode SNESMonitorLGDestroy(PetscDrawLG *draw)
3153: {
3157: KSPMonitorLGResidualNormDestroy(draw);
3158: return(0);
3159: }
3161: extern PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3164: PetscErrorCode SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3165: {
3166: PetscDrawLG lg;
3167: PetscErrorCode ierr;
3168: PetscReal x,y,per;
3169: PetscViewer v = (PetscViewer)monctx;
3170: static PetscReal prev; /* should be in the context */
3171: PetscDraw draw;
3174: PetscViewerDrawGetDrawLG(v,0,&lg);
3175: if (!n) {PetscDrawLGReset(lg);}
3176: PetscDrawLGGetDraw(lg,&draw);
3177: PetscDrawSetTitle(draw,"Residual norm");
3178: x = (PetscReal)n;
3179: if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3180: else y = -15.0;
3181: PetscDrawLGAddPoint(lg,&x,&y);
3182: if (n < 20 || !(n % 5)) {
3183: PetscDrawLGDraw(lg);
3184: }
3186: PetscViewerDrawGetDrawLG(v,1,&lg);
3187: if (!n) {PetscDrawLGReset(lg);}
3188: PetscDrawLGGetDraw(lg,&draw);
3189: PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3190: SNESMonitorRange_Private(snes,n,&per);
3191: x = (PetscReal)n;
3192: y = 100.0*per;
3193: PetscDrawLGAddPoint(lg,&x,&y);
3194: if (n < 20 || !(n % 5)) {
3195: PetscDrawLGDraw(lg);
3196: }
3198: PetscViewerDrawGetDrawLG(v,2,&lg);
3199: if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3200: PetscDrawLGGetDraw(lg,&draw);
3201: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3202: x = (PetscReal)n;
3203: y = (prev - rnorm)/prev;
3204: PetscDrawLGAddPoint(lg,&x,&y);
3205: if (n < 20 || !(n % 5)) {
3206: PetscDrawLGDraw(lg);
3207: }
3209: PetscViewerDrawGetDrawLG(v,3,&lg);
3210: if (!n) {PetscDrawLGReset(lg);}
3211: PetscDrawLGGetDraw(lg,&draw);
3212: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3213: x = (PetscReal)n;
3214: y = (prev - rnorm)/(prev*per);
3215: if (n > 2) { /*skip initial crazy value */
3216: PetscDrawLGAddPoint(lg,&x,&y);
3217: }
3218: if (n < 20 || !(n % 5)) {
3219: PetscDrawLGDraw(lg);
3220: }
3221: prev = rnorm;
3222: return(0);
3223: }
3227: /*@
3228: SNESMonitor - runs the user provided monitor routines, if they exist
3230: Collective on SNES
3232: Input Parameters:
3233: + snes - nonlinear solver context obtained from SNESCreate()
3234: . iter - iteration number
3235: - rnorm - relative norm of the residual
3237: Notes:
3238: This routine is called by the SNES implementations.
3239: It does not typically need to be called by the user.
3241: Level: developer
3243: .seealso: SNESMonitorSet()
3244: @*/
3245: PetscErrorCode SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3246: {
3248: PetscInt i,n = snes->numbermonitors;
3251: for (i=0; i<n; i++) {
3252: (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3253: }
3254: return(0);
3255: }
3257: /* ------------ Routines to set performance monitoring options ----------- */
3259: /*MC
3260: SNESMonitorFunction - functional form passed to SNESMonitorSet() to monitor convergence of nonlinear solver
3262: Synopsis:
3263: #include <petscsnes.h>
3264: $ PetscErrorCode SNESMonitorFunction(SNES snes,PetscInt its, PetscReal norm,void *mctx)
3266: + snes - the SNES context
3267: . its - iteration number
3268: . norm - 2-norm function value (may be estimated)
3269: - mctx - [optional] monitoring context
3271: Level: advanced
3273: .seealso: SNESMonitorSet(), SNESMonitorGet()
3274: M*/
3278: /*@C
3279: SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3280: iteration of the nonlinear solver to display the iteration's
3281: progress.
3283: Logically Collective on SNES
3285: Input Parameters:
3286: + snes - the SNES context
3287: . f - the monitor function, see SNESMonitorFunction for the calling sequence
3288: . mctx - [optional] user-defined context for private data for the
3289: monitor routine (use NULL if no context is desired)
3290: - monitordestroy - [optional] routine that frees monitor context
3291: (may be NULL)
3293: Options Database Keys:
3294: + -snes_monitor - sets SNESMonitorDefault()
3295: . -snes_monitor_lg_residualnorm - sets line graph monitor,
3296: uses SNESMonitorLGCreate()
3297: - -snes_monitor_cancel - cancels all monitors that have
3298: been hardwired into a code by
3299: calls to SNESMonitorSet(), but
3300: does not cancel those set via
3301: the options database.
3303: Notes:
3304: Several different monitoring routines may be set by calling
3305: SNESMonitorSet() multiple times; all will be called in the
3306: order in which they were set.
3308: Fortran notes: Only a single monitor function can be set for each SNES object
3310: Level: intermediate
3312: .keywords: SNES, nonlinear, set, monitor
3314: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3315: @*/
3316: PetscErrorCode SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3317: {
3318: PetscInt i;
3323: if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3324: for (i=0; i<snes->numbermonitors;i++) {
3325: if (f == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) {
3326: if (monitordestroy) {
3327: (*monitordestroy)(&mctx);
3328: }
3329: return(0);
3330: }
3331: }
3332: snes->monitor[snes->numbermonitors] = f;
3333: snes->monitordestroy[snes->numbermonitors] = monitordestroy;
3334: snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3335: return(0);
3336: }
3340: /*@
3341: SNESMonitorCancel - Clears all the monitor functions for a SNES object.
3343: Logically Collective on SNES
3345: Input Parameters:
3346: . snes - the SNES context
3348: Options Database Key:
3349: . -snes_monitor_cancel - cancels all monitors that have been hardwired
3350: into a code by calls to SNESMonitorSet(), but does not cancel those
3351: set via the options database
3353: Notes:
3354: There is no way to clear one specific monitor from a SNES object.
3356: Level: intermediate
3358: .keywords: SNES, nonlinear, set, monitor
3360: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3361: @*/
3362: PetscErrorCode SNESMonitorCancel(SNES snes)
3363: {
3365: PetscInt i;
3369: for (i=0; i<snes->numbermonitors; i++) {
3370: if (snes->monitordestroy[i]) {
3371: (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3372: }
3373: }
3374: snes->numbermonitors = 0;
3375: return(0);
3376: }
3378: /*MC
3379: SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver
3381: Synopsis:
3382: #include <petscsnes.h>
3383: $ PetscErrorCode SNESConvergenceTest(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
3385: + snes - the SNES context
3386: . it - current iteration (0 is the first and is before any Newton step)
3387: . cctx - [optional] convergence context
3388: . reason - reason for convergence/divergence
3389: . xnorm - 2-norm of current iterate
3390: . gnorm - 2-norm of current step
3391: - f - 2-norm of function
3393: Level: intermediate
3395: .seealso: SNESSetConvergenceTest(), SNESGetConvergenceTest()
3396: M*/
3400: /*@C
3401: SNESSetConvergenceTest - Sets the function that is to be used
3402: to test for convergence of the nonlinear iterative solution.
3404: Logically Collective on SNES
3406: Input Parameters:
3407: + snes - the SNES context
3408: . SNESConvergenceTestFunction - routine to test for convergence
3409: . cctx - [optional] context for private data for the convergence routine (may be NULL)
3410: - destroy - [optional] destructor for the context (may be NULL; NULL_FUNCTION in Fortran)
3412: Level: advanced
3414: .keywords: SNES, nonlinear, set, convergence, test
3416: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3417: @*/
3418: PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3419: {
3424: if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
3425: if (snes->ops->convergeddestroy) {
3426: (*snes->ops->convergeddestroy)(snes->cnvP);
3427: }
3428: snes->ops->converged = SNESConvergenceTestFunction;
3429: snes->ops->convergeddestroy = destroy;
3430: snes->cnvP = cctx;
3431: return(0);
3432: }
3436: /*@
3437: SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
3439: Not Collective
3441: Input Parameter:
3442: . snes - the SNES context
3444: Output Parameter:
3445: . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
3446: manual pages for the individual convergence tests for complete lists
3448: Level: intermediate
3450: Notes: Can only be called after the call the SNESSolve() is complete.
3452: .keywords: SNES, nonlinear, set, convergence, test
3454: .seealso: SNESSetConvergenceTest(), SNESConvergedReason
3455: @*/
3456: PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
3457: {
3461: *reason = snes->reason;
3462: return(0);
3463: }
3467: /*@
3468: SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
3470: Logically Collective on SNES
3472: Input Parameters:
3473: + snes - iterative context obtained from SNESCreate()
3474: . a - array to hold history, this array will contain the function norms computed at each step
3475: . its - integer array holds the number of linear iterations for each solve.
3476: . na - size of a and its
3477: - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
3478: else it continues storing new values for new nonlinear solves after the old ones
3480: Notes:
3481: If 'a' and 'its' are NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
3482: default array of length 10000 is allocated.
3484: This routine is useful, e.g., when running a code for purposes
3485: of accurate performance monitoring, when no I/O should be done
3486: during the section of code that is being timed.
3488: Level: intermediate
3490: .keywords: SNES, set, convergence, history
3492: .seealso: SNESGetConvergenceHistory()
3494: @*/
3495: PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
3496: {
3503: if (!a) {
3504: if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
3505: PetscCalloc1(na,&a);
3506: PetscCalloc1(na,&its);
3508: snes->conv_malloc = PETSC_TRUE;
3509: }
3510: snes->conv_hist = a;
3511: snes->conv_hist_its = its;
3512: snes->conv_hist_max = na;
3513: snes->conv_hist_len = 0;
3514: snes->conv_hist_reset = reset;
3515: return(0);
3516: }
3518: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3519: #include <engine.h> /* MATLAB include file */
3520: #include <mex.h> /* MATLAB include file */
3524: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
3525: {
3526: mxArray *mat;
3527: PetscInt i;
3528: PetscReal *ar;
3531: mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
3532: ar = (PetscReal*) mxGetData(mat);
3533: for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
3534: PetscFunctionReturn(mat);
3535: }
3536: #endif
3540: /*@C
3541: SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
3543: Not Collective
3545: Input Parameter:
3546: . snes - iterative context obtained from SNESCreate()
3548: Output Parameters:
3549: . a - array to hold history
3550: . its - integer array holds the number of linear iterations (or
3551: negative if not converged) for each solve.
3552: - na - size of a and its
3554: Notes:
3555: The calling sequence for this routine in Fortran is
3556: $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
3558: This routine is useful, e.g., when running a code for purposes
3559: of accurate performance monitoring, when no I/O should be done
3560: during the section of code that is being timed.
3562: Level: intermediate
3564: .keywords: SNES, get, convergence, history
3566: .seealso: SNESSetConvergencHistory()
3568: @*/
3569: PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
3570: {
3573: if (a) *a = snes->conv_hist;
3574: if (its) *its = snes->conv_hist_its;
3575: if (na) *na = snes->conv_hist_len;
3576: return(0);
3577: }
3581: /*@C
3582: SNESSetUpdate - Sets the general-purpose update function called
3583: at the beginning of every iteration of the nonlinear solve. Specifically
3584: it is called just before the Jacobian is "evaluated".
3586: Logically Collective on SNES
3588: Input Parameters:
3589: . snes - The nonlinear solver context
3590: . func - The function
3592: Calling sequence of func:
3593: . func (SNES snes, PetscInt step);
3595: . step - The current step of the iteration
3597: Level: advanced
3599: Note: This is NOT what one uses to update the ghost points before a function evaluation, that should be done at the beginning of your FormFunction()
3600: This is not used by most users.
3602: .keywords: SNES, update
3604: .seealso SNESSetJacobian(), SNESSolve()
3605: @*/
3606: PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
3607: {
3610: snes->ops->update = func;
3611: return(0);
3612: }
3616: /*
3617: SNESScaleStep_Private - Scales a step so that its length is less than the
3618: positive parameter delta.
3620: Input Parameters:
3621: + snes - the SNES context
3622: . y - approximate solution of linear system
3623: . fnorm - 2-norm of current function
3624: - delta - trust region size
3626: Output Parameters:
3627: + gpnorm - predicted function norm at the new point, assuming local
3628: linearization. The value is zero if the step lies within the trust
3629: region, and exceeds zero otherwise.
3630: - ynorm - 2-norm of the step
3632: Note:
3633: For non-trust region methods such as SNESNEWTONLS, the parameter delta
3634: is set to be the maximum allowable step size.
3636: .keywords: SNES, nonlinear, scale, step
3637: */
3638: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
3639: {
3640: PetscReal nrm;
3641: PetscScalar cnorm;
3649: VecNorm(y,NORM_2,&nrm);
3650: if (nrm > *delta) {
3651: nrm = *delta/nrm;
3652: *gpnorm = (1.0 - nrm)*(*fnorm);
3653: cnorm = nrm;
3654: VecScale(y,cnorm);
3655: *ynorm = *delta;
3656: } else {
3657: *gpnorm = 0.0;
3658: *ynorm = nrm;
3659: }
3660: return(0);
3661: }
3665: /*@C
3666: SNESSolve - Solves a nonlinear system F(x) = b.
3667: Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
3669: Collective on SNES
3671: Input Parameters:
3672: + snes - the SNES context
3673: . b - the constant part of the equation F(x) = b, or NULL to use zero.
3674: - x - the solution vector.
3676: Notes:
3677: The user should initialize the vector,x, with the initial guess
3678: for the nonlinear solve prior to calling SNESSolve. In particular,
3679: to employ an initial guess of zero, the user should explicitly set
3680: this vector to zero by calling VecSet().
3682: Level: beginner
3684: .keywords: SNES, nonlinear, solve
3686: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
3687: @*/
3688: PetscErrorCode SNESSolve(SNES snes,Vec b,Vec x)
3689: {
3690: PetscErrorCode ierr;
3691: PetscBool flg;
3692: PetscInt grid;
3693: Vec xcreated = NULL;
3694: DM dm;
3703: if (!x) {
3704: SNESGetDM(snes,&dm);
3705: DMCreateGlobalVector(dm,&xcreated);
3706: x = xcreated;
3707: }
3708: SNESViewFromOptions(snes,NULL,"-snes_view_pre");
3710: for (grid=0; grid<snes->gridsequence; grid++) {PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));}
3711: for (grid=0; grid<snes->gridsequence+1; grid++) {
3713: /* set solution vector */
3714: if (!grid) {PetscObjectReference((PetscObject)x);}
3715: VecDestroy(&snes->vec_sol);
3716: snes->vec_sol = x;
3717: SNESGetDM(snes,&dm);
3719: /* set affine vector if provided */
3720: if (b) { PetscObjectReference((PetscObject)b); }
3721: VecDestroy(&snes->vec_rhs);
3722: snes->vec_rhs = b;
3724: if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
3725: if (snes->vec_rhs == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
3726: if (!snes->vec_sol_update /* && snes->vec_sol */) {
3727: VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
3728: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
3729: }
3730: DMShellSetGlobalVector(dm,snes->vec_sol);
3731: SNESSetUp(snes);
3733: if (!grid) {
3734: if (snes->ops->computeinitialguess) {
3735: (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
3736: }
3737: }
3739: if (snes->conv_hist_reset) snes->conv_hist_len = 0;
3740: if (snes->counters_reset) {snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;}
3742: PetscLogEventBegin(SNES_Solve,snes,0,0,0);
3743: (*snes->ops->solve)(snes);
3744: PetscLogEventEnd(SNES_Solve,snes,0,0,0);
3745: if (snes->domainerror) {
3746: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
3747: snes->domainerror = PETSC_FALSE;
3748: }
3749: if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
3751: if (snes->lagjac_persist) snes->jac_iter += snes->iter;
3752: if (snes->lagpre_persist) snes->pre_iter += snes->iter;
3754: flg = PETSC_FALSE;
3755: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,NULL);
3756: if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
3757: if (snes->printreason) {
3758: PetscViewerASCIIAddTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),((PetscObject)snes)->tablevel);
3759: if (snes->reason > 0) {
3760: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3761: } else {
3762: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3763: }
3764: PetscViewerASCIISubtractTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)),((PetscObject)snes)->tablevel);
3765: }
3767: if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
3768: if (grid < snes->gridsequence) {
3769: DM fine;
3770: Vec xnew;
3771: Mat interp;
3773: DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
3774: if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
3775: DMCreateInterpolation(snes->dm,fine,&interp,NULL);
3776: DMCreateGlobalVector(fine,&xnew);
3777: MatInterpolate(interp,x,xnew);
3778: DMInterpolate(snes->dm,interp,fine);
3779: MatDestroy(&interp);
3780: x = xnew;
3782: SNESReset(snes);
3783: SNESSetDM(snes,fine);
3784: DMDestroy(&fine);
3785: PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
3786: }
3787: }
3788: SNESViewFromOptions(snes,NULL,"-snes_view");
3789: VecViewFromOptions(snes->vec_sol,((PetscObject)snes)->prefix,"-snes_view_solution");
3791: VecDestroy(&xcreated);
3792: PetscObjectSAWsBlock((PetscObject)snes);
3793: return(0);
3794: }
3796: /* --------- Internal routines for SNES Package --------- */
3800: /*@C
3801: SNESSetType - Sets the method for the nonlinear solver.
3803: Collective on SNES
3805: Input Parameters:
3806: + snes - the SNES context
3807: - type - a known method
3809: Options Database Key:
3810: . -snes_type <type> - Sets the method; use -help for a list
3811: of available methods (for instance, newtonls or newtontr)
3813: Notes:
3814: See "petsc/include/petscsnes.h" for available methods (for instance)
3815: + SNESNEWTONLS - Newton's method with line search
3816: (systems of nonlinear equations)
3817: . SNESNEWTONTR - Newton's method with trust region
3818: (systems of nonlinear equations)
3820: Normally, it is best to use the SNESSetFromOptions() command and then
3821: set the SNES solver type from the options database rather than by using
3822: this routine. Using the options database provides the user with
3823: maximum flexibility in evaluating the many nonlinear solvers.
3824: The SNESSetType() routine is provided for those situations where it
3825: is necessary to set the nonlinear solver independently of the command
3826: line or options database. This might be the case, for example, when
3827: the choice of solver changes during the execution of the program,
3828: and the user's application is taking responsibility for choosing the
3829: appropriate method.
3831: Developer Notes: SNESRegister() adds a constructor for a new SNESType to SNESList, SNESSetType() locates
3832: the constructor in that list and calls it to create the spexific object.
3834: Level: intermediate
3836: .keywords: SNES, set, type
3838: .seealso: SNESType, SNESCreate(), SNESDestroy(), SNESGetType(), SNESSetFromOptions()
3840: @*/
3841: PetscErrorCode SNESSetType(SNES snes,SNESType type)
3842: {
3843: PetscErrorCode ierr,(*r)(SNES);
3844: PetscBool match;
3850: PetscObjectTypeCompare((PetscObject)snes,type,&match);
3851: if (match) return(0);
3853: PetscFunctionListFind(SNESList,type,&r);
3854: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
3855: /* Destroy the previous private SNES context */
3856: if (snes->ops->destroy) {
3857: (*(snes)->ops->destroy)(snes);
3858: snes->ops->destroy = NULL;
3859: }
3860: /* Reinitialize function pointers in SNESOps structure */
3861: snes->ops->setup = 0;
3862: snes->ops->solve = 0;
3863: snes->ops->view = 0;
3864: snes->ops->setfromoptions = 0;
3865: snes->ops->destroy = 0;
3866: /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
3867: snes->setupcalled = PETSC_FALSE;
3869: PetscObjectChangeTypeName((PetscObject)snes,type);
3870: (*r)(snes);
3871: return(0);
3872: }
3876: /*@C
3877: SNESGetType - Gets the SNES method type and name (as a string).
3879: Not Collective
3881: Input Parameter:
3882: . snes - nonlinear solver context
3884: Output Parameter:
3885: . type - SNES method (a character string)
3887: Level: intermediate
3889: .keywords: SNES, nonlinear, get, type, name
3890: @*/
3891: PetscErrorCode SNESGetType(SNES snes,SNESType *type)
3892: {
3896: *type = ((PetscObject)snes)->type_name;
3897: return(0);
3898: }
3902: /*@
3903: SNESGetSolution - Returns the vector where the approximate solution is
3904: stored. This is the fine grid solution when using SNESSetGridSequence().
3906: Not Collective, but Vec is parallel if SNES is parallel
3908: Input Parameter:
3909: . snes - the SNES context
3911: Output Parameter:
3912: . x - the solution
3914: Level: intermediate
3916: .keywords: SNES, nonlinear, get, solution
3918: .seealso: SNESGetSolutionUpdate(), SNESGetFunction()
3919: @*/
3920: PetscErrorCode SNESGetSolution(SNES snes,Vec *x)
3921: {
3925: *x = snes->vec_sol;
3926: return(0);
3927: }
3931: /*@
3932: SNESGetSolutionUpdate - Returns the vector where the solution update is
3933: stored.
3935: Not Collective, but Vec is parallel if SNES is parallel
3937: Input Parameter:
3938: . snes - the SNES context
3940: Output Parameter:
3941: . x - the solution update
3943: Level: advanced
3945: .keywords: SNES, nonlinear, get, solution, update
3947: .seealso: SNESGetSolution(), SNESGetFunction()
3948: @*/
3949: PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x)
3950: {
3954: *x = snes->vec_sol_update;
3955: return(0);
3956: }
3960: /*@C
3961: SNESGetFunction - Returns the vector where the function is stored.
3963: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
3965: Input Parameter:
3966: . snes - the SNES context
3968: Output Parameter:
3969: + r - the vector that is used to store residuals (or NULL if you don't want it)
3970: . f - the function (or NULL if you don't want it); see SNESFunction for calling sequence details
3971: - ctx - the function context (or NULL if you don't want it)
3973: Level: advanced
3975: .keywords: SNES, nonlinear, get, function
3977: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
3978: @*/
3979: PetscErrorCode SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
3980: {
3982: DM dm;
3986: if (r) {
3987: if (!snes->vec_func) {
3988: if (snes->vec_rhs) {
3989: VecDuplicate(snes->vec_rhs,&snes->vec_func);
3990: } else if (snes->vec_sol) {
3991: VecDuplicate(snes->vec_sol,&snes->vec_func);
3992: } else if (snes->dm) {
3993: DMCreateGlobalVector(snes->dm,&snes->vec_func);
3994: }
3995: }
3996: *r = snes->vec_func;
3997: }
3998: SNESGetDM(snes,&dm);
3999: DMSNESGetFunction(dm,f,ctx);
4000: return(0);
4001: }
4003: /*@C
4004: SNESGetNGS - Returns the NGS function and context.
4006: Input Parameter:
4007: . snes - the SNES context
4009: Output Parameter:
4010: + f - the function (or NULL) see SNESNGSFunction for details
4011: - ctx - the function context (or NULL)
4013: Level: advanced
4015: .keywords: SNES, nonlinear, get, function
4017: .seealso: SNESSetNGS(), SNESGetFunction()
4018: @*/
4022: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4023: {
4025: DM dm;
4029: SNESGetDM(snes,&dm);
4030: DMSNESGetNGS(dm,f,ctx);
4031: return(0);
4032: }
4036: /*@C
4037: SNESSetOptionsPrefix - Sets the prefix used for searching for all
4038: SNES options in the database.
4040: Logically Collective on SNES
4042: Input Parameter:
4043: + snes - the SNES context
4044: - prefix - the prefix to prepend to all option names
4046: Notes:
4047: A hyphen (-) must NOT be given at the beginning of the prefix name.
4048: The first character of all runtime options is AUTOMATICALLY the hyphen.
4050: Level: advanced
4052: .keywords: SNES, set, options, prefix, database
4054: .seealso: SNESSetFromOptions()
4055: @*/
4056: PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[])
4057: {
4062: PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4063: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4064: if (snes->linesearch) {
4065: SNESGetLineSearch(snes,&snes->linesearch);
4066: PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4067: }
4068: KSPSetOptionsPrefix(snes->ksp,prefix);
4069: return(0);
4070: }
4074: /*@C
4075: SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4076: SNES options in the database.
4078: Logically Collective on SNES
4080: Input Parameters:
4081: + snes - the SNES context
4082: - prefix - the prefix to prepend to all option names
4084: Notes:
4085: A hyphen (-) must NOT be given at the beginning of the prefix name.
4086: The first character of all runtime options is AUTOMATICALLY the hyphen.
4088: Level: advanced
4090: .keywords: SNES, append, options, prefix, database
4092: .seealso: SNESGetOptionsPrefix()
4093: @*/
4094: PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4095: {
4100: PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4101: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4102: if (snes->linesearch) {
4103: SNESGetLineSearch(snes,&snes->linesearch);
4104: PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4105: }
4106: KSPAppendOptionsPrefix(snes->ksp,prefix);
4107: return(0);
4108: }
4112: /*@C
4113: SNESGetOptionsPrefix - Sets the prefix used for searching for all
4114: SNES options in the database.
4116: Not Collective
4118: Input Parameter:
4119: . snes - the SNES context
4121: Output Parameter:
4122: . prefix - pointer to the prefix string used
4124: Notes: On the fortran side, the user should pass in a string 'prefix' of
4125: sufficient length to hold the prefix.
4127: Level: advanced
4129: .keywords: SNES, get, options, prefix, database
4131: .seealso: SNESAppendOptionsPrefix()
4132: @*/
4133: PetscErrorCode SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4134: {
4139: PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4140: return(0);
4141: }
4146: /*@C
4147: SNESRegister - Adds a method to the nonlinear solver package.
4149: Not collective
4151: Input Parameters:
4152: + name_solver - name of a new user-defined solver
4153: - routine_create - routine to create method context
4155: Notes:
4156: SNESRegister() may be called multiple times to add several user-defined solvers.
4158: Sample usage:
4159: .vb
4160: SNESRegister("my_solver",MySolverCreate);
4161: .ve
4163: Then, your solver can be chosen with the procedural interface via
4164: $ SNESSetType(snes,"my_solver")
4165: or at runtime via the option
4166: $ -snes_type my_solver
4168: Level: advanced
4170: Note: If your function is not being put into a shared library then use SNESRegister() instead
4172: .keywords: SNES, nonlinear, register
4174: .seealso: SNESRegisterAll(), SNESRegisterDestroy()
4176: Level: advanced
4177: @*/
4178: PetscErrorCode SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4179: {
4183: PetscFunctionListAdd(&SNESList,sname,function);
4184: return(0);
4185: }
4189: PetscErrorCode SNESTestLocalMin(SNES snes)
4190: {
4192: PetscInt N,i,j;
4193: Vec u,uh,fh;
4194: PetscScalar value;
4195: PetscReal norm;
4198: SNESGetSolution(snes,&u);
4199: VecDuplicate(u,&uh);
4200: VecDuplicate(u,&fh);
4202: /* currently only works for sequential */
4203: PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4204: VecGetSize(u,&N);
4205: for (i=0; i<N; i++) {
4206: VecCopy(u,uh);
4207: PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4208: for (j=-10; j<11; j++) {
4209: value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4210: VecSetValue(uh,i,value,ADD_VALUES);
4211: SNESComputeFunction(snes,uh,fh);
4212: VecNorm(fh,NORM_2,&norm);
4213: PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);
4214: value = -value;
4215: VecSetValue(uh,i,value,ADD_VALUES);
4216: }
4217: }
4218: VecDestroy(&uh);
4219: VecDestroy(&fh);
4220: return(0);
4221: }
4225: /*@
4226: SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4227: computing relative tolerance for linear solvers within an inexact
4228: Newton method.
4230: Logically Collective on SNES
4232: Input Parameters:
4233: + snes - SNES context
4234: - flag - PETSC_TRUE or PETSC_FALSE
4236: Options Database:
4237: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4238: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
4239: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4240: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4241: . -snes_ksp_ew_gamma <gamma> - Sets gamma
4242: . -snes_ksp_ew_alpha <alpha> - Sets alpha
4243: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4244: - -snes_ksp_ew_threshold <threshold> - Sets threshold
4246: Notes:
4247: Currently, the default is to use a constant relative tolerance for
4248: the inner linear solvers. Alternatively, one can use the
4249: Eisenstat-Walker method, where the relative convergence tolerance
4250: is reset at each Newton iteration according progress of the nonlinear
4251: solver.
4253: Level: advanced
4255: Reference:
4256: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4257: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
4259: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
4261: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4262: @*/
4263: PetscErrorCode SNESKSPSetUseEW(SNES snes,PetscBool flag)
4264: {
4268: snes->ksp_ewconv = flag;
4269: return(0);
4270: }
4274: /*@
4275: SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
4276: for computing relative tolerance for linear solvers within an
4277: inexact Newton method.
4279: Not Collective
4281: Input Parameter:
4282: . snes - SNES context
4284: Output Parameter:
4285: . flag - PETSC_TRUE or PETSC_FALSE
4287: Notes:
4288: Currently, the default is to use a constant relative tolerance for
4289: the inner linear solvers. Alternatively, one can use the
4290: Eisenstat-Walker method, where the relative convergence tolerance
4291: is reset at each Newton iteration according progress of the nonlinear
4292: solver.
4294: Level: advanced
4296: Reference:
4297: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4298: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
4300: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
4302: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4303: @*/
4304: PetscErrorCode SNESKSPGetUseEW(SNES snes, PetscBool *flag)
4305: {
4309: *flag = snes->ksp_ewconv;
4310: return(0);
4311: }
4315: /*@
4316: SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
4317: convergence criteria for the linear solvers within an inexact
4318: Newton method.
4320: Logically Collective on SNES
4322: Input Parameters:
4323: + snes - SNES context
4324: . version - version 1, 2 (default is 2) or 3
4325: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4326: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4327: . gamma - multiplicative factor for version 2 rtol computation
4328: (0 <= gamma2 <= 1)
4329: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
4330: . alpha2 - power for safeguard
4331: - threshold - threshold for imposing safeguard (0 < threshold < 1)
4333: Note:
4334: Version 3 was contributed by Luis Chacon, June 2006.
4336: Use PETSC_DEFAULT to retain the default for any of the parameters.
4338: Level: advanced
4340: Reference:
4341: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4342: inexact Newton method", Utah State University Math. Stat. Dept. Res.
4343: Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
4345: .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
4347: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
4348: @*/
4349: PetscErrorCode SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
4350: {
4351: SNESKSPEW *kctx;
4355: kctx = (SNESKSPEW*)snes->kspconvctx;
4356: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4365: if (version != PETSC_DEFAULT) kctx->version = version;
4366: if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0;
4367: if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max;
4368: if (gamma != PETSC_DEFAULT) kctx->gamma = gamma;
4369: if (alpha != PETSC_DEFAULT) kctx->alpha = alpha;
4370: if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2;
4371: if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
4373: if (kctx->version < 1 || kctx->version > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
4374: if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %g",(double)kctx->rtol_0);
4375: if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%g) < 1.0\n",(double)kctx->rtol_max);
4376: if (kctx->gamma < 0.0 || kctx->gamma > 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%g) <= 1.0\n",(double)kctx->gamma);
4377: if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%g) <= 2.0\n",(double)kctx->alpha);
4378: if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%g) < 1.0\n",(double)kctx->threshold);
4379: return(0);
4380: }
4384: /*@
4385: SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
4386: convergence criteria for the linear solvers within an inexact
4387: Newton method.
4389: Not Collective
4391: Input Parameters:
4392: snes - SNES context
4394: Output Parameters:
4395: + version - version 1, 2 (default is 2) or 3
4396: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4397: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4398: . gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
4399: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
4400: . alpha2 - power for safeguard
4401: - threshold - threshold for imposing safeguard (0 < threshold < 1)
4403: Level: advanced
4405: .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
4407: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
4408: @*/
4409: PetscErrorCode SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
4410: {
4411: SNESKSPEW *kctx;
4415: kctx = (SNESKSPEW*)snes->kspconvctx;
4416: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4417: if (version) *version = kctx->version;
4418: if (rtol_0) *rtol_0 = kctx->rtol_0;
4419: if (rtol_max) *rtol_max = kctx->rtol_max;
4420: if (gamma) *gamma = kctx->gamma;
4421: if (alpha) *alpha = kctx->alpha;
4422: if (alpha2) *alpha2 = kctx->alpha2;
4423: if (threshold) *threshold = kctx->threshold;
4424: return(0);
4425: }
4429: PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4430: {
4432: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
4433: PetscReal rtol = PETSC_DEFAULT,stol;
4436: if (!snes->ksp_ewconv) return(0);
4437: if (!snes->iter) {
4438: rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
4439: VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
4440: }
4441: else {
4442: if (kctx->version == 1) {
4443: rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
4444: if (rtol < 0.0) rtol = -rtol;
4445: stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
4446: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4447: } else if (kctx->version == 2) {
4448: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4449: stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
4450: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4451: } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
4452: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4453: /* safeguard: avoid sharp decrease of rtol */
4454: stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
4455: stol = PetscMax(rtol,stol);
4456: rtol = PetscMin(kctx->rtol_0,stol);
4457: /* safeguard: avoid oversolving */
4458: stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
4459: stol = PetscMax(rtol,stol);
4460: rtol = PetscMin(kctx->rtol_0,stol);
4461: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
4462: }
4463: /* safeguard: avoid rtol greater than one */
4464: rtol = PetscMin(rtol,kctx->rtol_max);
4465: KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
4466: PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
4467: return(0);
4468: }
4472: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4473: {
4475: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
4476: PCSide pcside;
4477: Vec lres;
4480: if (!snes->ksp_ewconv) return(0);
4481: KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
4482: kctx->norm_last = snes->norm;
4483: if (kctx->version == 1) {
4484: KSPGetPCSide(ksp,&pcside);
4485: if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
4486: /* KSP residual is true linear residual */
4487: KSPGetResidualNorm(ksp,&kctx->lresid_last);
4488: } else {
4489: /* KSP residual is preconditioned residual */
4490: /* compute true linear residual norm */
4491: VecDuplicate(b,&lres);
4492: MatMult(snes->jacobian,x,lres);
4493: VecAYPX(lres,-1.0,b);
4494: VecNorm(lres,NORM_2,&kctx->lresid_last);
4495: VecDestroy(&lres);
4496: }
4497: }
4498: return(0);
4499: }
4503: /*@
4504: SNESGetKSP - Returns the KSP context for a SNES solver.
4506: Not Collective, but if SNES object is parallel, then KSP object is parallel
4508: Input Parameter:
4509: . snes - the SNES context
4511: Output Parameter:
4512: . ksp - the KSP context
4514: Notes:
4515: The user can then directly manipulate the KSP context to set various
4516: options, etc. Likewise, the user can then extract and manipulate the
4517: PC contexts as well.
4519: Level: beginner
4521: .keywords: SNES, nonlinear, get, KSP, context
4523: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
4524: @*/
4525: PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp)
4526: {
4533: if (!snes->ksp) {
4534: KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
4535: PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
4536: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);
4538: KSPSetPreSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPreSolve_SNESEW,snes);
4539: KSPSetPostSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPostSolve_SNESEW,snes);
4540: }
4541: *ksp = snes->ksp;
4542: return(0);
4543: }
4546: #include <petsc-private/dmimpl.h>
4549: /*@
4550: SNESSetDM - Sets the DM that may be used by some preconditioners
4552: Logically Collective on SNES
4554: Input Parameters:
4555: + snes - the preconditioner context
4556: - dm - the dm
4558: Level: intermediate
4560: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
4561: @*/
4562: PetscErrorCode SNESSetDM(SNES snes,DM dm)
4563: {
4565: KSP ksp;
4566: DMSNES sdm;
4570: if (dm) {PetscObjectReference((PetscObject)dm);}
4571: if (snes->dm) { /* Move the DMSNES context over to the new DM unless the new DM already has one */
4572: if (snes->dm->dmsnes && snes->dmAuto && !dm->dmsnes) {
4573: DMCopyDMSNES(snes->dm,dm);
4574: DMGetDMSNES(snes->dm,&sdm);
4575: if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
4576: }
4577: DMDestroy(&snes->dm);
4578: }
4579: snes->dm = dm;
4580: snes->dmAuto = PETSC_FALSE;
4582: SNESGetKSP(snes,&ksp);
4583: KSPSetDM(ksp,dm);
4584: KSPSetDMActive(ksp,PETSC_FALSE);
4585: if (snes->pc) {
4586: SNESSetDM(snes->pc, snes->dm);
4587: SNESSetNPCSide(snes,snes->pcside);
4588: }
4589: return(0);
4590: }
4594: /*@
4595: SNESGetDM - Gets the DM that may be used by some preconditioners
4597: Not Collective but DM obtained is parallel on SNES
4599: Input Parameter:
4600: . snes - the preconditioner context
4602: Output Parameter:
4603: . dm - the dm
4605: Level: intermediate
4607: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
4608: @*/
4609: PetscErrorCode SNESGetDM(SNES snes,DM *dm)
4610: {
4615: if (!snes->dm) {
4616: DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
4617: snes->dmAuto = PETSC_TRUE;
4618: }
4619: *dm = snes->dm;
4620: return(0);
4621: }
4625: /*@
4626: SNESSetNPC - Sets the nonlinear preconditioner to be used.
4628: Collective on SNES
4630: Input Parameters:
4631: + snes - iterative context obtained from SNESCreate()
4632: - pc - the preconditioner object
4634: Notes:
4635: Use SNESGetNPC() to retrieve the preconditioner context (for example,
4636: to configure it using the API).
4638: Level: developer
4640: .keywords: SNES, set, precondition
4641: .seealso: SNESGetNPC()
4642: @*/
4643: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
4644: {
4651: PetscObjectReference((PetscObject) pc);
4652: SNESDestroy(&snes->pc);
4653: snes->pc = pc;
4654: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->pc);
4655: return(0);
4656: }
4660: /*@
4661: SNESGetNPC - Creates a nonlinear preconditioning solver (SNES) to be used to precondition the nonlinear solver.
4663: Not Collective
4665: Input Parameter:
4666: . snes - iterative context obtained from SNESCreate()
4668: Output Parameter:
4669: . pc - preconditioner context
4671: Notes: If a SNES was previously set with SNESSetNPC() then that SNES is returned.
4673: Level: developer
4675: .keywords: SNES, get, preconditioner
4676: .seealso: SNESSetNPC()
4677: @*/
4678: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
4679: {
4681: const char *optionsprefix;
4686: if (!snes->pc) {
4687: SNESCreate(PetscObjectComm((PetscObject)snes),&snes->pc);
4688: PetscObjectIncrementTabLevel((PetscObject)snes->pc,(PetscObject)snes,1);
4689: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->pc);
4690: SNESGetOptionsPrefix(snes,&optionsprefix);
4691: SNESSetOptionsPrefix(snes->pc,optionsprefix);
4692: SNESAppendOptionsPrefix(snes->pc,"npc_");
4693: SNESSetCountersReset(snes->pc,PETSC_FALSE);
4694: }
4695: *pc = snes->pc;
4696: return(0);
4697: }
4701: /*@
4702: SNESSetNPCSide - Sets the preconditioning side.
4704: Logically Collective on SNES
4706: Input Parameter:
4707: . snes - iterative context obtained from SNESCreate()
4709: Output Parameter:
4710: . side - the preconditioning side, where side is one of
4711: .vb
4712: PC_LEFT - left preconditioning (default)
4713: PC_RIGHT - right preconditioning
4714: .ve
4716: Options Database Keys:
4717: . -snes_pc_side <right,left>
4719: Level: intermediate
4721: .keywords: SNES, set, right, left, side, preconditioner, flag
4723: .seealso: SNESGetNPCSide(), KSPSetPCSide()
4724: @*/
4725: PetscErrorCode SNESSetNPCSide(SNES snes,PCSide side)
4726: {
4730: snes->pcside = side;
4731: return(0);
4732: }
4736: /*@
4737: SNESGetNPCSide - Gets the preconditioning side.
4739: Not Collective
4741: Input Parameter:
4742: . snes - iterative context obtained from SNESCreate()
4744: Output Parameter:
4745: . side - the preconditioning side, where side is one of
4746: .vb
4747: PC_LEFT - left preconditioning (default)
4748: PC_RIGHT - right preconditioning
4749: .ve
4751: Level: intermediate
4753: .keywords: SNES, get, right, left, side, preconditioner, flag
4755: .seealso: SNESSetNPCSide(), KSPGetPCSide()
4756: @*/
4757: PetscErrorCode SNESGetNPCSide(SNES snes,PCSide *side)
4758: {
4762: *side = snes->pcside;
4763: return(0);
4764: }
4768: /*@
4769: SNESSetLineSearch - Sets the linesearch on the SNES instance.
4771: Collective on SNES
4773: Input Parameters:
4774: + snes - iterative context obtained from SNESCreate()
4775: - linesearch - the linesearch object
4777: Notes:
4778: Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
4779: to configure it using the API).
4781: Level: developer
4783: .keywords: SNES, set, linesearch
4784: .seealso: SNESGetLineSearch()
4785: @*/
4786: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
4787: {
4794: PetscObjectReference((PetscObject) linesearch);
4795: SNESLineSearchDestroy(&snes->linesearch);
4797: snes->linesearch = linesearch;
4799: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
4800: return(0);
4801: }
4805: /*@
4806: SNESGetLineSearch - Returns a pointer to the line search context set with SNESSetLineSearch()
4807: or creates a default line search instance associated with the SNES and returns it.
4809: Not Collective
4811: Input Parameter:
4812: . snes - iterative context obtained from SNESCreate()
4814: Output Parameter:
4815: . linesearch - linesearch context
4817: Level: beginner
4819: .keywords: SNES, get, linesearch
4820: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
4821: @*/
4822: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
4823: {
4825: const char *optionsprefix;
4830: if (!snes->linesearch) {
4831: SNESGetOptionsPrefix(snes, &optionsprefix);
4832: SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
4833: SNESLineSearchSetSNES(snes->linesearch, snes);
4834: SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
4835: PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
4836: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
4837: }
4838: *linesearch = snes->linesearch;
4839: return(0);
4840: }
4842: #if defined(PETSC_HAVE_MATLAB_ENGINE)
4843: #include <mex.h>
4845: typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
4849: /*
4850: SNESComputeFunction_Matlab - Calls the function that has been set with SNESSetFunctionMatlab().
4852: Collective on SNES
4854: Input Parameters:
4855: + snes - the SNES context
4856: - x - input vector
4858: Output Parameter:
4859: . y - function vector, as set by SNESSetFunction()
4861: Notes:
4862: SNESComputeFunction() is typically used within nonlinear solvers
4863: implementations, so most users would not generally call this routine
4864: themselves.
4866: Level: developer
4868: .keywords: SNES, nonlinear, compute, function
4870: .seealso: SNESSetFunction(), SNESGetFunction()
4871: */
4872: PetscErrorCode SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
4873: {
4874: PetscErrorCode ierr;
4875: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
4876: int nlhs = 1,nrhs = 5;
4877: mxArray *plhs[1],*prhs[5];
4878: long long int lx = 0,ly = 0,ls = 0;
4887: /* call Matlab function in ctx with arguments x and y */
4889: PetscMemcpy(&ls,&snes,sizeof(snes));
4890: PetscMemcpy(&lx,&x,sizeof(x));
4891: PetscMemcpy(&ly,&y,sizeof(x));
4892: prhs[0] = mxCreateDoubleScalar((double)ls);
4893: prhs[1] = mxCreateDoubleScalar((double)lx);
4894: prhs[2] = mxCreateDoubleScalar((double)ly);
4895: prhs[3] = mxCreateString(sctx->funcname);
4896: prhs[4] = sctx->ctx;
4897: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");
4898: mxGetScalar(plhs[0]);
4899: mxDestroyArray(prhs[0]);
4900: mxDestroyArray(prhs[1]);
4901: mxDestroyArray(prhs[2]);
4902: mxDestroyArray(prhs[3]);
4903: mxDestroyArray(plhs[0]);
4904: return(0);
4905: }
4909: /*
4910: SNESSetFunctionMatlab - Sets the function evaluation routine and function
4911: vector for use by the SNES routines in solving systems of nonlinear
4912: equations from MATLAB. Here the function is a string containing the name of a MATLAB function
4914: Logically Collective on SNES
4916: Input Parameters:
4917: + snes - the SNES context
4918: . r - vector to store function value
4919: - f - function evaluation routine
4921: Notes:
4922: The Newton-like methods typically solve linear systems of the form
4923: $ f'(x) x = -f(x),
4924: where f'(x) denotes the Jacobian matrix and f(x) is the function.
4926: Level: beginner
4928: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
4930: .keywords: SNES, nonlinear, set, function
4932: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
4933: */
4934: PetscErrorCode SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx)
4935: {
4936: PetscErrorCode ierr;
4937: SNESMatlabContext *sctx;
4940: /* currently sctx is memory bleed */
4941: PetscMalloc(sizeof(SNESMatlabContext),&sctx);
4942: PetscStrallocpy(f,&sctx->funcname);
4943: /*
4944: This should work, but it doesn't
4945: sctx->ctx = ctx;
4946: mexMakeArrayPersistent(sctx->ctx);
4947: */
4948: sctx->ctx = mxDuplicateArray(ctx);
4949: SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);
4950: return(0);
4951: }
4955: /*
4956: SNESComputeJacobian_Matlab - Calls the function that has been set with SNESSetJacobianMatlab().
4958: Collective on SNES
4960: Input Parameters:
4961: + snes - the SNES context
4962: . x - input vector
4963: . A, B - the matrices
4964: - ctx - user context
4966: Level: developer
4968: .keywords: SNES, nonlinear, compute, function
4970: .seealso: SNESSetFunction(), SNESGetFunction()
4971: @*/
4972: PetscErrorCode SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx)
4973: {
4974: PetscErrorCode ierr;
4975: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
4976: int nlhs = 2,nrhs = 6;
4977: mxArray *plhs[2],*prhs[6];
4978: long long int lx = 0,lA = 0,ls = 0, lB = 0;
4984: /* call Matlab function in ctx with arguments x and y */
4986: PetscMemcpy(&ls,&snes,sizeof(snes));
4987: PetscMemcpy(&lx,&x,sizeof(x));
4988: PetscMemcpy(&lA,A,sizeof(x));
4989: PetscMemcpy(&lB,B,sizeof(x));
4990: prhs[0] = mxCreateDoubleScalar((double)ls);
4991: prhs[1] = mxCreateDoubleScalar((double)lx);
4992: prhs[2] = mxCreateDoubleScalar((double)lA);
4993: prhs[3] = mxCreateDoubleScalar((double)lB);
4994: prhs[4] = mxCreateString(sctx->funcname);
4995: prhs[5] = sctx->ctx;
4996: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");
4997: mxGetScalar(plhs[0]);
4998: mxDestroyArray(prhs[0]);
4999: mxDestroyArray(prhs[1]);
5000: mxDestroyArray(prhs[2]);
5001: mxDestroyArray(prhs[3]);
5002: mxDestroyArray(prhs[4]);
5003: mxDestroyArray(plhs[0]);
5004: mxDestroyArray(plhs[1]);
5005: return(0);
5006: }
5010: /*
5011: SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
5012: vector for use by the SNES routines in solving systems of nonlinear
5013: equations from MATLAB. Here the function is a string containing the name of a MATLAB function
5015: Logically Collective on SNES
5017: Input Parameters:
5018: + snes - the SNES context
5019: . A,B - Jacobian matrices
5020: . J - function evaluation routine
5021: - ctx - user context
5023: Level: developer
5025: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5027: .keywords: SNES, nonlinear, set, function
5029: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J
5030: */
5031: PetscErrorCode SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx)
5032: {
5033: PetscErrorCode ierr;
5034: SNESMatlabContext *sctx;
5037: /* currently sctx is memory bleed */
5038: PetscMalloc(sizeof(SNESMatlabContext),&sctx);
5039: PetscStrallocpy(J,&sctx->funcname);
5040: /*
5041: This should work, but it doesn't
5042: sctx->ctx = ctx;
5043: mexMakeArrayPersistent(sctx->ctx);
5044: */
5045: sctx->ctx = mxDuplicateArray(ctx);
5046: SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);
5047: return(0);
5048: }
5052: /*
5053: SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().
5055: Collective on SNES
5057: .seealso: SNESSetFunction(), SNESGetFunction()
5058: @*/
5059: PetscErrorCode SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
5060: {
5061: PetscErrorCode ierr;
5062: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5063: int nlhs = 1,nrhs = 6;
5064: mxArray *plhs[1],*prhs[6];
5065: long long int lx = 0,ls = 0;
5066: Vec x = snes->vec_sol;
5071: PetscMemcpy(&ls,&snes,sizeof(snes));
5072: PetscMemcpy(&lx,&x,sizeof(x));
5073: prhs[0] = mxCreateDoubleScalar((double)ls);
5074: prhs[1] = mxCreateDoubleScalar((double)it);
5075: prhs[2] = mxCreateDoubleScalar((double)fnorm);
5076: prhs[3] = mxCreateDoubleScalar((double)lx);
5077: prhs[4] = mxCreateString(sctx->funcname);
5078: prhs[5] = sctx->ctx;
5079: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");
5080: mxGetScalar(plhs[0]);
5081: mxDestroyArray(prhs[0]);
5082: mxDestroyArray(prhs[1]);
5083: mxDestroyArray(prhs[2]);
5084: mxDestroyArray(prhs[3]);
5085: mxDestroyArray(prhs[4]);
5086: mxDestroyArray(plhs[0]);
5087: return(0);
5088: }
5092: /*
5093: SNESMonitorSetMatlab - Sets the monitor function from MATLAB
5095: Level: developer
5097: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5099: .keywords: SNES, nonlinear, set, function
5101: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5102: */
5103: PetscErrorCode SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx)
5104: {
5105: PetscErrorCode ierr;
5106: SNESMatlabContext *sctx;
5109: /* currently sctx is memory bleed */
5110: PetscMalloc(sizeof(SNESMatlabContext),&sctx);
5111: PetscStrallocpy(f,&sctx->funcname);
5112: /*
5113: This should work, but it doesn't
5114: sctx->ctx = ctx;
5115: mexMakeArrayPersistent(sctx->ctx);
5116: */
5117: sctx->ctx = mxDuplicateArray(ctx);
5118: SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);
5119: return(0);
5120: }
5122: #endif