Actual source code: snes.c
petsc-3.6.1 2015-08-06
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: if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates input vector is not in the function domain");
95: snes->domainerror = PETSC_TRUE;
96: return(0);
97: }
101: /*@
102: SNESGetFunctionDomainError - Gets the status of the domain error after a call to SNESComputeFunction;
104: Logically Collective on SNES
106: Input Parameters:
107: . snes - the SNES context
109: Output Parameters:
110: . domainerror - Set to PETSC_TRUE if there's a domain error; PETSC_FALSE otherwise.
112: Level: advanced
114: .keywords: SNES, view
116: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction()
117: @*/
118: PetscErrorCode SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror)
119: {
123: *domainerror = snes->domainerror;
124: return(0);
125: }
129: /*@C
130: SNESLoad - Loads a SNES that has been stored in binary with SNESView().
132: Collective on PetscViewer
134: Input Parameters:
135: + newdm - the newly loaded SNES, this needs to have been created with SNESCreate() or
136: some related function before a call to SNESLoad().
137: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
139: Level: intermediate
141: Notes:
142: The type is determined by the data in the file, any type set into the SNES before this call is ignored.
144: Notes for advanced users:
145: Most users should not need to know the details of the binary storage
146: format, since SNESLoad() and TSView() completely hide these details.
147: But for anyone who's interested, the standard binary matrix storage
148: format is
149: .vb
150: has not yet been determined
151: .ve
153: .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad()
154: @*/
155: PetscErrorCode SNESLoad(SNES snes, PetscViewer viewer)
156: {
158: PetscBool isbinary;
159: PetscInt classid;
160: char type[256];
161: KSP ksp;
162: DM dm;
163: DMSNES dmsnes;
168: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
169: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
171: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
172: if (classid != SNES_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONG,"Not SNES next in file");
173: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
174: SNESSetType(snes, type);
175: if (snes->ops->load) {
176: (*snes->ops->load)(snes,viewer);
177: }
178: SNESGetDM(snes,&dm);
179: DMGetDMSNES(dm,&dmsnes);
180: DMSNESLoad(dmsnes,viewer);
181: SNESGetKSP(snes,&ksp);
182: KSPLoad(ksp,viewer);
183: return(0);
184: }
186: #include <petscdraw.h>
187: #if defined(PETSC_HAVE_SAWS)
188: #include <petscviewersaws.h>
189: #endif
192: /*@C
193: SNESView - Prints the SNES data structure.
195: Collective on SNES
197: Input Parameters:
198: + SNES - the SNES context
199: - viewer - visualization context
201: Options Database Key:
202: . -snes_view - Calls SNESView() at end of SNESSolve()
204: Notes:
205: The available visualization contexts include
206: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
207: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
208: output where only the first processor opens
209: the file. All other processors send their
210: data to the first processor to print.
212: The user can open an alternative visualization context with
213: PetscViewerASCIIOpen() - output to a specified file.
215: Level: beginner
217: .keywords: SNES, view
219: .seealso: PetscViewerASCIIOpen()
220: @*/
221: PetscErrorCode SNESView(SNES snes,PetscViewer viewer)
222: {
223: SNESKSPEW *kctx;
225: KSP ksp;
226: SNESLineSearch linesearch;
227: PetscBool iascii,isstring,isbinary,isdraw;
228: DMSNES dmsnes;
229: #if defined(PETSC_HAVE_SAWS)
230: PetscBool issaws;
231: #endif
235: if (!viewer) {
236: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
237: }
241: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
242: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
243: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
244: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
245: #if defined(PETSC_HAVE_SAWS)
246: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
247: #endif
248: if (iascii) {
249: SNESNormSchedule normschedule;
251: PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);
252: if (!snes->setupcalled) {
253: PetscViewerASCIIPrintf(viewer," SNES has not been set up so information may be incomplete\n");
254: }
255: if (snes->ops->view) {
256: PetscViewerASCIIPushTab(viewer);
257: (*snes->ops->view)(snes,viewer);
258: PetscViewerASCIIPopTab(viewer);
259: }
260: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);
261: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);
262: PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);
263: PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);
264: SNESGetNormSchedule(snes, &normschedule);
265: if (normschedule > 0) {PetscViewerASCIIPrintf(viewer," norm schedule %s\n",SNESNormSchedules[normschedule]);}
266: if (snes->gridsequence) {
267: PetscViewerASCIIPrintf(viewer," total number of grid sequence refinements=%D\n",snes->gridsequence);
268: }
269: if (snes->ksp_ewconv) {
270: kctx = (SNESKSPEW*)snes->kspconvctx;
271: if (kctx) {
272: PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);
273: PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);
274: PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);
275: }
276: }
277: if (snes->lagpreconditioner == -1) {
278: PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");
279: } else if (snes->lagpreconditioner > 1) {
280: PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);
281: }
282: if (snes->lagjacobian == -1) {
283: PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");
284: } else if (snes->lagjacobian > 1) {
285: PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);
286: }
287: } else if (isstring) {
288: const char *type;
289: SNESGetType(snes,&type);
290: PetscViewerStringSPrintf(viewer," %-3.3s",type);
291: } else if (isbinary) {
292: PetscInt classid = SNES_FILE_CLASSID;
293: MPI_Comm comm;
294: PetscMPIInt rank;
295: char type[256];
297: PetscObjectGetComm((PetscObject)snes,&comm);
298: MPI_Comm_rank(comm,&rank);
299: if (!rank) {
300: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
301: PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));
302: PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR,PETSC_FALSE);
303: }
304: if (snes->ops->view) {
305: (*snes->ops->view)(snes,viewer);
306: }
307: } else if (isdraw) {
308: PetscDraw draw;
309: char str[36];
310: PetscReal x,y,bottom,h;
312: PetscViewerDrawGetDraw(viewer,0,&draw);
313: PetscDrawGetCurrentPoint(draw,&x,&y);
314: PetscStrcpy(str,"SNES: ");
315: PetscStrcat(str,((PetscObject)snes)->type_name);
316: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);
317: bottom = y - h;
318: PetscDrawPushCurrentPoint(draw,x,bottom);
319: if (snes->ops->view) {
320: (*snes->ops->view)(snes,viewer);
321: }
322: #if defined(PETSC_HAVE_SAWS)
323: } else if (issaws) {
324: PetscMPIInt rank;
325: const char *name;
327: PetscObjectGetName((PetscObject)snes,&name);
328: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
329: if (!((PetscObject)snes)->amsmem && !rank) {
330: char dir[1024];
332: PetscObjectViewSAWs((PetscObject)snes,viewer);
333: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
334: PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT));
335: if (!snes->conv_hist) {
336: SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);
337: }
338: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);
339: PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE));
340: }
341: #endif
342: }
343: if (snes->linesearch) {
344: PetscViewerASCIIPushTab(viewer);
345: SNESGetLineSearch(snes, &linesearch);
346: SNESLineSearchView(linesearch, viewer);
347: PetscViewerASCIIPopTab(viewer);
348: }
349: if (snes->pc && snes->usespc) {
350: PetscViewerASCIIPushTab(viewer);
351: SNESView(snes->pc, viewer);
352: PetscViewerASCIIPopTab(viewer);
353: }
354: PetscViewerASCIIPushTab(viewer);
355: DMGetDMSNES(snes->dm,&dmsnes);
356: DMSNESView(dmsnes, viewer);
357: PetscViewerASCIIPopTab(viewer);
358: if (snes->usesksp) {
359: SNESGetKSP(snes,&ksp);
360: PetscViewerASCIIPushTab(viewer);
361: KSPView(ksp,viewer);
362: PetscViewerASCIIPopTab(viewer);
363: }
364: if (isdraw) {
365: PetscDraw draw;
366: PetscViewerDrawGetDraw(viewer,0,&draw);
367: PetscDrawPopCurrentPoint(draw);
368: }
369: return(0);
370: }
372: /*
373: We retain a list of functions that also take SNES command
374: line options. These are called at the end SNESSetFromOptions()
375: */
376: #define MAXSETFROMOPTIONS 5
377: static PetscInt numberofsetfromoptions;
378: static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
382: /*@C
383: SNESAddOptionsChecker - Adds an additional function to check for SNES options.
385: Not Collective
387: Input Parameter:
388: . snescheck - function that checks for options
390: Level: developer
392: .seealso: SNESSetFromOptions()
393: @*/
394: PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
395: {
397: if (numberofsetfromoptions >= MAXSETFROMOPTIONS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
398: othersetfromoptions[numberofsetfromoptions++] = snescheck;
399: return(0);
400: }
402: extern PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
406: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
407: {
408: Mat J;
409: KSP ksp;
410: PC pc;
411: PetscBool match;
413: MatNullSpace nullsp;
418: if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
419: Mat A = snes->jacobian, B = snes->jacobian_pre;
420: MatCreateVecs(A ? A : B, NULL,&snes->vec_func);
421: }
423: if (version == 1) {
424: MatCreateSNESMF(snes,&J);
425: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
426: MatSetFromOptions(J);
427: } else if (version == 2) {
428: if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
429: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
430: SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);
431: #else
432: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)");
433: #endif
434: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2");
436: /* attach any user provided null space that was on Amat to the newly created matrix free matrix */
437: if (snes->jacobian) {
438: MatGetNullSpace(snes->jacobian,&nullsp);
439: if (nullsp) {
440: MatSetNullSpace(J,nullsp);
441: }
442: }
444: PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);
445: if (hasOperator) {
447: /* This version replaces the user provided Jacobian matrix with a
448: matrix-free version but still employs the user-provided preconditioner matrix. */
449: SNESSetJacobian(snes,J,0,0,0);
450: } else {
451: /* This version replaces both the user-provided Jacobian and the user-
452: provided preconditioner Jacobian with the default matrix free version. */
453: if ((snes->pcside == PC_LEFT) && snes->pc) {
454: if (!snes->jacobian){SNESSetJacobian(snes,J,0,0,0);}
455: } else {
456: SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,0);
457: }
458: /* Force no preconditioner */
459: SNESGetKSP(snes,&ksp);
460: KSPGetPC(ksp,&pc);
461: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);
462: if (!match) {
463: PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");
464: PCSetType(pc,PCNONE);
465: }
466: }
467: MatDestroy(&J);
468: return(0);
469: }
473: static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx)
474: {
475: SNES snes = (SNES)ctx;
477: Vec Xfine,Xfine_named = NULL,Xcoarse;
480: if (PetscLogPrintInfo) {
481: PetscInt finelevel,coarselevel,fineclevel,coarseclevel;
482: DMGetRefineLevel(dmfine,&finelevel);
483: DMGetCoarsenLevel(dmfine,&fineclevel);
484: DMGetRefineLevel(dmcoarse,&coarselevel);
485: DMGetCoarsenLevel(dmcoarse,&coarseclevel);
486: PetscInfo4(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);
487: }
488: if (dmfine == snes->dm) Xfine = snes->vec_sol;
489: else {
490: DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
491: Xfine = Xfine_named;
492: }
493: DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
494: if (Inject) {
495: MatRestrict(Inject,Xfine,Xcoarse);
496: } else {
497: MatRestrict(Restrict,Xfine,Xcoarse);
498: VecPointwiseMult(Xcoarse,Xcoarse,Rscale);
499: }
500: DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
501: if (Xfine_named) {DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);}
502: return(0);
503: }
507: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
508: {
512: DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
513: return(0);
514: }
518: /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can
519: * safely call SNESGetDM() in their residual evaluation routine. */
520: static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx)
521: {
522: SNES snes = (SNES)ctx;
524: Mat Asave = A,Bsave = B;
525: Vec X,Xnamed = NULL;
526: DM dmsave;
527: void *ctxsave;
528: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
531: dmsave = snes->dm;
532: KSPGetDM(ksp,&snes->dm);
533: if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */
534: else { /* We are on a coarser level, this vec was initialized using a DM restrict hook */
535: DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
536: X = Xnamed;
537: SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);
538: /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */
539: if (jac == SNESComputeJacobianDefaultColor) {
540: SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,0);
541: }
542: }
543: /* put the previous context back */
545: SNESComputeJacobian(snes,X,A,B);
546: if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
547: SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
548: }
550: if (A != Asave || B != Bsave) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for changing matrices at this time");
551: if (Xnamed) {
552: DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
553: }
554: snes->dm = dmsave;
555: return(0);
556: }
560: /*@
561: SNESSetUpMatrices - ensures that matrices are available for SNES, to be called by SNESSetUp_XXX()
563: Collective
565: Input Arguments:
566: . snes - snes to configure
568: Level: developer
570: .seealso: SNESSetUp()
571: @*/
572: PetscErrorCode SNESSetUpMatrices(SNES snes)
573: {
575: DM dm;
576: DMSNES sdm;
579: SNESGetDM(snes,&dm);
580: DMGetDMSNES(dm,&sdm);
581: if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"DMSNES not properly configured");
582: else if (!snes->jacobian && snes->mf) {
583: Mat J;
584: void *functx;
585: MatCreateSNESMF(snes,&J);
586: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
587: MatSetFromOptions(J);
588: SNESGetFunction(snes,NULL,NULL,&functx);
589: SNESSetJacobian(snes,J,J,0,0);
590: MatDestroy(&J);
591: } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
592: Mat J,B;
593: MatCreateSNESMF(snes,&J);
594: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
595: MatSetFromOptions(J);
596: DMCreateMatrix(snes->dm,&B);
597: /* sdm->computejacobian was already set to reach here */
598: SNESSetJacobian(snes,J,B,NULL,NULL);
599: MatDestroy(&J);
600: MatDestroy(&B);
601: } else if (!snes->jacobian_pre) {
602: Mat J,B;
603: J = snes->jacobian;
604: DMCreateMatrix(snes->dm,&B);
605: SNESSetJacobian(snes,J ? J : B,B,NULL,NULL);
606: MatDestroy(&B);
607: }
608: {
609: KSP ksp;
610: SNESGetKSP(snes,&ksp);
611: KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);
612: DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
613: }
614: return(0);
615: }
619: /*@
620: SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
622: Collective on SNES
624: Input Parameter:
625: . snes - the SNES context
627: Options Database Keys:
628: + -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
629: . -snes_stol - convergence tolerance in terms of the norm
630: of the change in the solution between steps
631: . -snes_atol <abstol> - absolute tolerance of residual norm
632: . -snes_rtol <rtol> - relative decrease in tolerance norm from initial
633: . -snes_max_it <max_it> - maximum number of iterations
634: . -snes_max_funcs <max_funcs> - maximum number of function evaluations
635: . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
636: . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
637: . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
638: . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
639: . -snes_trtol <trtol> - trust region tolerance
640: . -snes_no_convergence_test - skip convergence test in nonlinear
641: solver; hence iterations will continue until max_it
642: or some other criterion is reached. Saves expense
643: of convergence test
644: . -snes_monitor <optional filename> - prints residual norm at each iteration. if no
645: filename given prints to stdout
646: . -snes_monitor_solution - plots solution at each iteration
647: . -snes_monitor_residual - plots residual (not its norm) at each iteration
648: . -snes_monitor_solution_update - plots update to solution at each iteration
649: . -snes_monitor_lg_residualnorm - plots residual norm at each iteration
650: . -snes_monitor_lg_range - plots residual norm at each iteration
651: . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
652: . -snes_fd_color - use finite differences with coloring to compute Jacobian
653: . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
654: - -snes_converged_reason - print the reason for convergence/divergence after each solve
656: Options Database for Eisenstat-Walker method:
657: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
658: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
659: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
660: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
661: . -snes_ksp_ew_gamma <gamma> - Sets gamma
662: . -snes_ksp_ew_alpha <alpha> - Sets alpha
663: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
664: - -snes_ksp_ew_threshold <threshold> - Sets threshold
666: Notes:
667: To see all options, run your program with the -help option or consult
668: Users-Manual: ch_snes
670: Level: beginner
672: .keywords: SNES, nonlinear, set, options, database
674: .seealso: SNESSetOptionsPrefix()
675: @*/
676: PetscErrorCode SNESSetFromOptions(SNES snes)
677: {
678: PetscBool flg,pcset,persist,set;
679: PetscInt i,indx,lag,grids;
680: const char *deft = SNESNEWTONLS;
681: const char *convtests[] = {"default","skip"};
682: SNESKSPEW *kctx = NULL;
683: char type[256], monfilename[PETSC_MAX_PATH_LEN];
684: PetscViewer monviewer;
686: PCSide pcside;
687: const char *optionsprefix;
691: SNESRegisterAll();
692: PetscObjectOptionsBegin((PetscObject)snes);
693: if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
694: PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
695: if (flg) {
696: SNESSetType(snes,type);
697: } else if (!((PetscObject)snes)->type_name) {
698: SNESSetType(snes,deft);
699: }
700: PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);
701: PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);
703: PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,NULL);
704: PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,NULL);
705: PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,NULL);
706: PetscOptionsInt("-snes_max_fail","Maximum nonlinear step failures","SNESSetMaxNonlinearStepFailures",snes->maxFailures,&snes->maxFailures,NULL);
707: PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,NULL);
708: PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,NULL);
710: PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
711: if (flg) {
712: SNESSetLagPreconditioner(snes,lag);
713: }
714: PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
715: if (flg) {
716: SNESSetLagPreconditionerPersists(snes,persist);
717: }
718: PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
719: if (flg) {
720: SNESSetLagJacobian(snes,lag);
721: }
722: PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
723: if (flg) {
724: SNESSetLagJacobianPersists(snes,persist);
725: }
727: PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
728: if (flg) {
729: SNESSetGridSequence(snes,grids);
730: }
732: PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);
733: if (flg) {
734: switch (indx) {
735: case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
736: case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL); break;
737: }
738: }
740: PetscOptionsEList("-snes_norm_schedule","SNES Norm schedule","SNESSetNormSchedule",SNESNormSchedules,5,"function",&indx,&flg);
741: if (flg) { SNESSetNormSchedule(snes,(SNESNormSchedule)indx); }
743: PetscOptionsEList("-snes_function_type","SNES Norm schedule","SNESSetFunctionType",SNESFunctionTypes,2,"unpreconditioned",&indx,&flg);
744: if (flg) { SNESSetFunctionType(snes,(SNESFunctionType)indx); }
746: kctx = (SNESKSPEW*)snes->kspconvctx;
748: PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,NULL);
750: PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);
751: PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);
752: PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);
753: PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);
754: PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);
755: PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);
756: PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);
758: flg = PETSC_FALSE;
759: PetscOptionsBool("-snes_check_jacobian","Check each Jacobian with a differenced one","SNESUpdateCheckJacobian",flg,&flg,&set);
760: if (set && flg) {
761: SNESSetUpdate(snes,SNESUpdateCheckJacobian);
762: }
764: flg = PETSC_FALSE;
765: PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);
766: if (set && flg) {SNESMonitorCancel(snes);}
768: PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
769: if (flg) {
770: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
771: SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
772: }
774: PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
775: if (flg) {
776: SNESMonitorSet(snes,SNESMonitorRange,0,0);
777: }
779: PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
780: if (flg) {
781: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
782: SNESMonitorSetRatio(snes,monviewer);
783: }
785: PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
786: if (flg) {
787: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
788: SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
789: }
791: PetscOptionsString("-snes_monitor_field","Monitor norm of function (split into fields)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
792: if (flg) {
793: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
794: SNESMonitorSet(snes,SNESMonitorDefaultField,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
795: }
797: PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
798: if (flg) {PetscPythonMonitorSet((PetscObject)snes,monfilename);}
800: flg = PETSC_FALSE;
801: PetscOptionsBool("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,NULL);
802: if (flg) {SNESMonitorSet(snes,SNESMonitorSolution,0,0);}
803: flg = PETSC_FALSE;
804: PetscOptionsBool("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,NULL);
805: if (flg) {SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);}
806: flg = PETSC_FALSE;
807: PetscOptionsBool("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,NULL);
808: if (flg) {SNESMonitorSet(snes,SNESMonitorResidual,0,0);}
809: flg = PETSC_FALSE;
810: PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);
811: if (flg) {
812: PetscObject *objs;
814: SNESMonitorLGCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&objs);
815: SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorLGResidualNorm,objs,(PetscErrorCode (*)(void**))SNESMonitorLGDestroy);
816: }
817: flg = PETSC_FALSE;
818: PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
819: if (flg) {
820: PetscViewer ctx;
822: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&ctx);
823: SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
824: }
826: flg = PETSC_FALSE;
827: PetscOptionsBool("-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",flg,&flg,NULL);
828: if (flg) {SNESMonitorSet(snes,SNESMonitorJacUpdateSpectrum,0,0);}
831: PetscOptionsString("-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
832: if (flg) {
833: PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),monfilename,&monviewer);
834: SNESMonitorSet(snes,SNESMonitorFields,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
835: }
837: flg = PETSC_FALSE;
838: PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
839: if (flg) {
840: void *functx;
841: SNESGetFunction(snes,NULL,NULL,&functx);
842: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
843: PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
844: }
846: flg = PETSC_FALSE;
847: PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
848: if (flg) {
849: SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
850: }
852: flg = PETSC_FALSE;
853: PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);
854: if (flg) {
855: DM dm;
856: DMSNES sdm;
857: SNESGetDM(snes,&dm);
858: DMGetDMSNES(dm,&sdm);
859: sdm->jacobianctx = NULL;
860: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,0);
861: PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
862: }
864: flg = PETSC_FALSE;
865: PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&snes->mf_operator,&flg);
866: if (flg && snes->mf_operator) {
867: snes->mf_operator = PETSC_TRUE;
868: snes->mf = PETSC_TRUE;
869: }
870: flg = PETSC_FALSE;
871: PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&snes->mf,&flg);
872: if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE;
873: PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,0);
875: flg = PETSC_FALSE;
876: SNESGetNPCSide(snes,&pcside);
877: PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
878: if (flg) {SNESSetNPCSide(snes,pcside);}
880: #if defined(PETSC_HAVE_SAWS)
881: /*
882: Publish convergence information using SAWs
883: */
884: flg = PETSC_FALSE;
885: PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);
886: if (flg) {
887: void *ctx;
888: SNESMonitorSAWsCreate(snes,&ctx);
889: SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);
890: }
891: #endif
892: #if defined(PETSC_HAVE_SAWS)
893: {
894: PetscBool set;
895: flg = PETSC_FALSE;
896: PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);
897: if (set) {
898: PetscObjectSAWsSetBlock((PetscObject)snes,flg);
899: }
900: }
901: #endif
903: for (i = 0; i < numberofsetfromoptions; i++) {
904: (*othersetfromoptions[i])(snes);
905: }
907: if (snes->ops->setfromoptions) {
908: (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
909: }
911: /* process any options handlers added with PetscObjectAddOptionsHandler() */
912: PetscObjectProcessOptionsHandlers((PetscObject)snes);
913: PetscOptionsEnd();
915: if (!snes->linesearch) {
916: SNESGetLineSearch(snes, &snes->linesearch);
917: }
918: SNESLineSearchSetFromOptions(snes->linesearch);
920: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
921: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
922: KSPSetFromOptions(snes->ksp);
924: /* if someone has set the SNES NPC type, create it. */
925: SNESGetOptionsPrefix(snes, &optionsprefix);
926: PetscOptionsHasName(optionsprefix, "-npc_snes_type", &pcset);
927: if (pcset && (!snes->pc)) {
928: SNESGetNPC(snes, &snes->pc);
929: }
930: return(0);
931: }
935: /*@C
936: SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
937: the nonlinear solvers.
939: Logically Collective on SNES
941: Input Parameters:
942: + snes - the SNES context
943: . compute - function to compute the context
944: - destroy - function to destroy the context
946: Level: intermediate
948: Notes:
949: This function is currently not available from Fortran.
951: .keywords: SNES, nonlinear, set, application, context
953: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
954: @*/
955: PetscErrorCode SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
956: {
959: snes->ops->usercompute = compute;
960: snes->ops->userdestroy = destroy;
961: return(0);
962: }
966: /*@
967: SNESSetApplicationContext - Sets the optional user-defined context for
968: the nonlinear solvers.
970: Logically Collective on SNES
972: Input Parameters:
973: + snes - the SNES context
974: - usrP - optional user context
976: Level: intermediate
978: .keywords: SNES, nonlinear, set, application, context
980: .seealso: SNESGetApplicationContext()
981: @*/
982: PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP)
983: {
985: KSP ksp;
989: SNESGetKSP(snes,&ksp);
990: KSPSetApplicationContext(ksp,usrP);
991: snes->user = usrP;
992: return(0);
993: }
997: /*@
998: SNESGetApplicationContext - Gets the user-defined context for the
999: nonlinear solvers.
1001: Not Collective
1003: Input Parameter:
1004: . snes - SNES context
1006: Output Parameter:
1007: . usrP - user context
1009: Level: intermediate
1011: .keywords: SNES, nonlinear, get, application, context
1013: .seealso: SNESSetApplicationContext()
1014: @*/
1015: PetscErrorCode SNESGetApplicationContext(SNES snes,void *usrP)
1016: {
1019: *(void**)usrP = snes->user;
1020: return(0);
1021: }
1025: /*@
1026: SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1027: at this time.
1029: Not Collective
1031: Input Parameter:
1032: . snes - SNES context
1034: Output Parameter:
1035: . iter - iteration number
1037: Notes:
1038: For example, during the computation of iteration 2 this would return 1.
1040: This is useful for using lagged Jacobians (where one does not recompute the
1041: Jacobian at each SNES iteration). For example, the code
1042: .vb
1043: SNESGetIterationNumber(snes,&it);
1044: if (!(it % 2)) {
1045: [compute Jacobian here]
1046: }
1047: .ve
1048: can be used in your ComputeJacobian() function to cause the Jacobian to be
1049: recomputed every second SNES iteration.
1051: Level: intermediate
1053: .keywords: SNES, nonlinear, get, iteration, number,
1055: .seealso: SNESGetLinearSolveIterations()
1056: @*/
1057: PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt *iter)
1058: {
1062: *iter = snes->iter;
1063: return(0);
1064: }
1068: /*@
1069: SNESSetIterationNumber - Sets the current iteration number.
1071: Not Collective
1073: Input Parameter:
1074: . snes - SNES context
1075: . iter - iteration number
1077: Level: developer
1079: .keywords: SNES, nonlinear, set, iteration, number,
1081: .seealso: SNESGetLinearSolveIterations()
1082: @*/
1083: PetscErrorCode SNESSetIterationNumber(SNES snes,PetscInt iter)
1084: {
1089: PetscObjectSAWsTakeAccess((PetscObject)snes);
1090: snes->iter = iter;
1091: PetscObjectSAWsGrantAccess((PetscObject)snes);
1092: return(0);
1093: }
1097: /*@
1098: SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1099: attempted by the nonlinear solver.
1101: Not Collective
1103: Input Parameter:
1104: . snes - SNES context
1106: Output Parameter:
1107: . nfails - number of unsuccessful steps attempted
1109: Notes:
1110: This counter is reset to zero for each successive call to SNESSolve().
1112: Level: intermediate
1114: .keywords: SNES, nonlinear, get, number, unsuccessful, steps
1116: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1117: SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1118: @*/
1119: PetscErrorCode SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1120: {
1124: *nfails = snes->numFailures;
1125: return(0);
1126: }
1130: /*@
1131: SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1132: attempted by the nonlinear solver before it gives up.
1134: Not Collective
1136: Input Parameters:
1137: + snes - SNES context
1138: - maxFails - maximum of unsuccessful steps
1140: Level: intermediate
1142: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
1144: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1145: SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1146: @*/
1147: PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1148: {
1151: snes->maxFailures = maxFails;
1152: return(0);
1153: }
1157: /*@
1158: SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1159: attempted by the nonlinear solver before it gives up.
1161: Not Collective
1163: Input Parameter:
1164: . snes - SNES context
1166: Output Parameter:
1167: . maxFails - maximum of unsuccessful steps
1169: Level: intermediate
1171: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1173: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1174: SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1176: @*/
1177: PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1178: {
1182: *maxFails = snes->maxFailures;
1183: return(0);
1184: }
1188: /*@
1189: SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1190: done by SNES.
1192: Not Collective
1194: Input Parameter:
1195: . snes - SNES context
1197: Output Parameter:
1198: . nfuncs - number of evaluations
1200: Level: intermediate
1202: Notes: Reset every time SNESSolve is called unless SNESSetCountersReset() is used.
1204: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1206: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1207: @*/
1208: PetscErrorCode SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1209: {
1213: *nfuncs = snes->nfuncs;
1214: return(0);
1215: }
1219: /*@
1220: SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1221: linear solvers.
1223: Not Collective
1225: Input Parameter:
1226: . snes - SNES context
1228: Output Parameter:
1229: . nfails - number of failed solves
1231: Level: intermediate
1233: Options Database Keys:
1234: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1236: Notes:
1237: This counter is reset to zero for each successive call to SNESSolve().
1239: .keywords: SNES, nonlinear, get, number, unsuccessful, steps
1241: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1242: @*/
1243: PetscErrorCode SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1244: {
1248: *nfails = snes->numLinearSolveFailures;
1249: return(0);
1250: }
1254: /*@
1255: SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1256: allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
1258: Logically Collective on SNES
1260: Input Parameters:
1261: + snes - SNES context
1262: - maxFails - maximum allowed linear solve failures
1264: Level: intermediate
1266: Options Database Keys:
1267: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1269: Notes: By default this is 0; that is SNES returns on the first failed linear solve
1271: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
1273: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1274: @*/
1275: PetscErrorCode SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1276: {
1280: snes->maxLinearSolveFailures = maxFails;
1281: return(0);
1282: }
1286: /*@
1287: SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1288: are allowed before SNES terminates
1290: Not Collective
1292: Input Parameter:
1293: . snes - SNES context
1295: Output Parameter:
1296: . maxFails - maximum of unsuccessful solves allowed
1298: Level: intermediate
1300: Notes: By default this is 1; that is SNES returns on the first failed linear solve
1302: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1304: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1305: @*/
1306: PetscErrorCode SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1307: {
1311: *maxFails = snes->maxLinearSolveFailures;
1312: return(0);
1313: }
1317: /*@
1318: SNESGetLinearSolveIterations - Gets the total number of linear iterations
1319: used by the nonlinear solver.
1321: Not Collective
1323: Input Parameter:
1324: . snes - SNES context
1326: Output Parameter:
1327: . lits - number of linear iterations
1329: Notes:
1330: This counter is reset to zero for each successive call to SNESSolve() unless SNESSetCountersReset() is used.
1332: Level: intermediate
1334: .keywords: SNES, nonlinear, get, number, linear, iterations
1336: .seealso: SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1337: @*/
1338: PetscErrorCode SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1339: {
1343: *lits = snes->linear_its;
1344: return(0);
1345: }
1349: /*@
1350: SNESSetCountersReset - Sets whether or not the counters for linear iterations and function evaluations
1351: are reset every time SNESSolve() is called.
1353: Logically Collective on SNES
1355: Input Parameter:
1356: + snes - SNES context
1357: - reset - whether to reset the counters or not
1359: Notes:
1360: This defaults to PETSC_TRUE
1362: Level: developer
1364: .keywords: SNES, nonlinear, set, reset, number, linear, iterations
1366: .seealso: SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC()
1367: @*/
1368: PetscErrorCode SNESSetCountersReset(SNES snes,PetscBool reset)
1369: {
1373: snes->counters_reset = reset;
1374: return(0);
1375: }
1380: /*@
1381: SNESSetKSP - Sets a KSP context for the SNES object to use
1383: Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
1385: Input Parameters:
1386: + snes - the SNES context
1387: - ksp - the KSP context
1389: Notes:
1390: The SNES object already has its KSP object, you can obtain with SNESGetKSP()
1391: so this routine is rarely needed.
1393: The KSP object that is already in the SNES object has its reference count
1394: decreased by one.
1396: Level: developer
1398: .keywords: SNES, nonlinear, get, KSP, context
1400: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1401: @*/
1402: PetscErrorCode SNESSetKSP(SNES snes,KSP ksp)
1403: {
1410: PetscObjectReference((PetscObject)ksp);
1411: if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1412: snes->ksp = ksp;
1413: return(0);
1414: }
1416: /* -----------------------------------------------------------*/
1419: /*@
1420: SNESCreate - Creates a nonlinear solver context.
1422: Collective on MPI_Comm
1424: Input Parameters:
1425: . comm - MPI communicator
1427: Output Parameter:
1428: . outsnes - the new SNES context
1430: Options Database Keys:
1431: + -snes_mf - Activates default matrix-free Jacobian-vector products,
1432: and no preconditioning matrix
1433: . -snes_mf_operator - Activates default matrix-free Jacobian-vector
1434: products, and a user-provided preconditioning matrix
1435: as set by SNESSetJacobian()
1436: - -snes_fd - Uses (slow!) finite differences to compute Jacobian
1438: Level: beginner
1440: .keywords: SNES, nonlinear, create, context
1442: .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()
1444: @*/
1445: PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes)
1446: {
1448: SNES snes;
1449: SNESKSPEW *kctx;
1453: *outsnes = NULL;
1454: SNESInitializePackage();
1456: PetscHeaderCreate(snes,SNES_CLASSID,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);
1458: snes->ops->converged = SNESConvergedDefault;
1459: snes->usesksp = PETSC_TRUE;
1460: snes->tolerancesset = PETSC_FALSE;
1461: snes->max_its = 50;
1462: snes->max_funcs = 10000;
1463: snes->norm = 0.0;
1464: snes->normschedule = SNES_NORM_ALWAYS;
1465: snes->functype = SNES_FUNCTION_DEFAULT;
1466: #if defined(PETSC_USE_REAL_SINGLE)
1467: snes->rtol = 1.e-5;
1468: #else
1469: snes->rtol = 1.e-8;
1470: #endif
1471: snes->ttol = 0.0;
1472: #if defined(PETSC_USE_REAL_SINGLE)
1473: snes->abstol = 1.e-25;
1474: #else
1475: snes->abstol = 1.e-50;
1476: #endif
1477: snes->stol = 1.e-8;
1478: #if defined(PETSC_USE_REAL_SINGLE)
1479: snes->deltatol = 1.e-6;
1480: #else
1481: snes->deltatol = 1.e-12;
1482: #endif
1483: snes->nfuncs = 0;
1484: snes->numFailures = 0;
1485: snes->maxFailures = 1;
1486: snes->linear_its = 0;
1487: snes->lagjacobian = 1;
1488: snes->jac_iter = 0;
1489: snes->lagjac_persist = PETSC_FALSE;
1490: snes->lagpreconditioner = 1;
1491: snes->pre_iter = 0;
1492: snes->lagpre_persist = PETSC_FALSE;
1493: snes->numbermonitors = 0;
1494: snes->data = 0;
1495: snes->setupcalled = PETSC_FALSE;
1496: snes->ksp_ewconv = PETSC_FALSE;
1497: snes->nwork = 0;
1498: snes->work = 0;
1499: snes->nvwork = 0;
1500: snes->vwork = 0;
1501: snes->conv_hist_len = 0;
1502: snes->conv_hist_max = 0;
1503: snes->conv_hist = NULL;
1504: snes->conv_hist_its = NULL;
1505: snes->conv_hist_reset = PETSC_TRUE;
1506: snes->counters_reset = PETSC_TRUE;
1507: snes->vec_func_init_set = PETSC_FALSE;
1508: snes->reason = SNES_CONVERGED_ITERATING;
1509: snes->pcside = PC_RIGHT;
1511: snes->mf = PETSC_FALSE;
1512: snes->mf_operator = PETSC_FALSE;
1513: snes->mf_version = 1;
1515: snes->numLinearSolveFailures = 0;
1516: snes->maxLinearSolveFailures = 1;
1518: snes->vizerotolerance = 1.e-8;
1520: /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
1521: PetscNewLog(snes,&kctx);
1523: snes->kspconvctx = (void*)kctx;
1524: kctx->version = 2;
1525: kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1526: this was too large for some test cases */
1527: kctx->rtol_last = 0.0;
1528: kctx->rtol_max = .9;
1529: kctx->gamma = 1.0;
1530: kctx->alpha = .5*(1.0 + PetscSqrtReal(5.0));
1531: kctx->alpha2 = kctx->alpha;
1532: kctx->threshold = .1;
1533: kctx->lresid_last = 0.0;
1534: kctx->norm_last = 0.0;
1536: *outsnes = snes;
1537: return(0);
1538: }
1540: /*MC
1541: SNESFunction - Functional form used to convey the nonlinear function to be solved by SNES
1543: Synopsis:
1544: #include "petscsnes.h"
1545: PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);
1547: Input Parameters:
1548: + snes - the SNES context
1549: . x - state at which to evaluate residual
1550: - ctx - optional user-defined function context, passed in with SNESSetFunction()
1552: Output Parameter:
1553: . f - vector to put residual (function value)
1555: Level: intermediate
1557: .seealso: SNESSetFunction(), SNESGetFunction()
1558: M*/
1562: /*@C
1563: SNESSetFunction - Sets the function evaluation routine and function
1564: vector for use by the SNES routines in solving systems of nonlinear
1565: equations.
1567: Logically Collective on SNES
1569: Input Parameters:
1570: + snes - the SNES context
1571: . r - vector to store function value
1572: . f - function evaluation routine; see SNESFunction for calling sequence details
1573: - ctx - [optional] user-defined context for private data for the
1574: function evaluation routine (may be NULL)
1576: Notes:
1577: The Newton-like methods typically solve linear systems of the form
1578: $ f'(x) x = -f(x),
1579: where f'(x) denotes the Jacobian matrix and f(x) is the function.
1581: Level: beginner
1583: .keywords: SNES, nonlinear, set, function
1585: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1586: @*/
1587: PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1588: {
1590: DM dm;
1594: if (r) {
1597: PetscObjectReference((PetscObject)r);
1598: VecDestroy(&snes->vec_func);
1600: snes->vec_func = r;
1601: }
1602: SNESGetDM(snes,&dm);
1603: DMSNESSetFunction(dm,f,ctx);
1604: return(0);
1605: }
1610: /*@C
1611: SNESSetInitialFunction - Sets the function vector to be used as the
1612: function norm at the initialization of the method. In some
1613: instances, the user has precomputed the function before calling
1614: SNESSolve. This function allows one to avoid a redundant call
1615: to SNESComputeFunction in that case.
1617: Logically Collective on SNES
1619: Input Parameters:
1620: + snes - the SNES context
1621: - f - vector to store function value
1623: Notes:
1624: This should not be modified during the solution procedure.
1626: This is used extensively in the SNESFAS hierarchy and in nonlinear preconditioning.
1628: Level: developer
1630: .keywords: SNES, nonlinear, set, function
1632: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1633: @*/
1634: PetscErrorCode SNESSetInitialFunction(SNES snes, Vec f)
1635: {
1637: Vec vec_func;
1643: if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1644: snes->vec_func_init_set = PETSC_FALSE;
1645: return(0);
1646: }
1647: SNESGetFunction(snes,&vec_func,NULL,NULL);
1648: VecCopy(f, vec_func);
1650: snes->vec_func_init_set = PETSC_TRUE;
1651: return(0);
1652: }
1656: /*@
1657: SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1658: of the SNES method.
1660: Logically Collective on SNES
1662: Input Parameters:
1663: + snes - the SNES context
1664: - normschedule - the frequency of norm computation
1666: Options Database Key:
1667: . -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly>
1669: Notes:
1670: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1671: of the nonlinear function and the taking of its norm at every iteration to
1672: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1673: (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
1674: may either be monitored for convergence or not. As these are often used as nonlinear
1675: preconditioners, monitoring the norm of their error is not a useful enterprise within
1676: their solution.
1678: Level: developer
1680: .keywords: SNES, nonlinear, set, function, norm, type
1682: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1683: @*/
1684: PetscErrorCode SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1685: {
1688: snes->normschedule = normschedule;
1689: return(0);
1690: }
1695: /*@
1696: SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1697: of the SNES method.
1699: Logically Collective on SNES
1701: Input Parameters:
1702: + snes - the SNES context
1703: - normschedule - the type of the norm used
1705: Level: advanced
1707: .keywords: SNES, nonlinear, set, function, norm, type
1709: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1710: @*/
1711: PetscErrorCode SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1712: {
1715: *normschedule = snes->normschedule;
1716: return(0);
1717: }
1722: /*@C
1723: SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
1724: of the SNES method.
1726: Logically Collective on SNES
1728: Input Parameters:
1729: + snes - the SNES context
1730: - normschedule - the frequency of norm computation
1732: Notes:
1733: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1734: of the nonlinear function and the taking of its norm at every iteration to
1735: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1736: (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
1737: may either be monitored for convergence or not. As these are often used as nonlinear
1738: preconditioners, monitoring the norm of their error is not a useful enterprise within
1739: their solution.
1741: Level: developer
1743: .keywords: SNES, nonlinear, set, function, norm, type
1745: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1746: @*/
1747: PetscErrorCode SNESSetFunctionType(SNES snes, SNESFunctionType type)
1748: {
1751: snes->functype = type;
1752: return(0);
1753: }
1758: /*@C
1759: SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
1760: of the SNES method.
1762: Logically Collective on SNES
1764: Input Parameters:
1765: + snes - the SNES context
1766: - normschedule - the type of the norm used
1768: Level: advanced
1770: .keywords: SNES, nonlinear, set, function, norm, type
1772: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1773: @*/
1774: PetscErrorCode SNESGetFunctionType(SNES snes, SNESFunctionType *type)
1775: {
1778: *type = snes->functype;
1779: return(0);
1780: }
1782: /*MC
1783: SNESNGSFunction - function used to convey a Gauss-Seidel sweep on the nonlinear function
1785: Synopsis:
1786: #include <petscsnes.h>
1787: $ SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);
1789: + X - solution vector
1790: . B - RHS vector
1791: - ctx - optional user-defined Gauss-Seidel context
1793: Level: intermediate
1795: .seealso: SNESSetNGS(), SNESGetNGS()
1796: M*/
1800: /*@C
1801: SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
1802: use with composed nonlinear solvers.
1804: Input Parameters:
1805: + snes - the SNES context
1806: . f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
1807: - ctx - [optional] user-defined context for private data for the
1808: smoother evaluation routine (may be NULL)
1810: Notes:
1811: The NGS routines are used by the composed nonlinear solver to generate
1812: a problem appropriate update to the solution, particularly FAS.
1814: Level: intermediate
1816: .keywords: SNES, nonlinear, set, Gauss-Seidel
1818: .seealso: SNESGetFunction(), SNESComputeNGS()
1819: @*/
1820: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1821: {
1823: DM dm;
1827: SNESGetDM(snes,&dm);
1828: DMSNESSetNGS(dm,f,ctx);
1829: return(0);
1830: }
1834: PETSC_EXTERN PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
1835: {
1837: DM dm;
1838: DMSNES sdm;
1841: SNESGetDM(snes,&dm);
1842: DMGetDMSNES(dm,&sdm);
1843: /* A(x)*x - b(x) */
1844: if (sdm->ops->computepfunction) {
1845: (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
1846: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");
1848: if (sdm->ops->computepjacobian) {
1849: (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
1850: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard matrix.");
1851: VecScale(f,-1.0);
1852: MatMultAdd(snes->jacobian,x,f,f);
1853: return(0);
1854: }
1858: PETSC_EXTERN PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
1859: {
1861: /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
1862: return(0);
1863: }
1867: /*@C
1868: SNESSetPicard - Use SNES to solve the semilinear-system A(x) x = b(x) via a Picard type iteration (Picard linearization)
1870: Logically Collective on SNES
1872: Input Parameters:
1873: + snes - the SNES context
1874: . r - vector to store function value
1875: . b - function evaluation routine
1876: . Amat - matrix with which A(x) x - b(x) is to be computed
1877: . Pmat - matrix from which preconditioner is computed (usually the same as Amat)
1878: . J - function to compute matrix value, see SNESJacobianFunction for details on its calling sequence
1879: - ctx - [optional] user-defined context for private data for the
1880: function evaluation routine (may be NULL)
1882: Notes:
1883: 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
1884: 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.
1886: One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both
1888: $ 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}
1889: $ Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration.
1891: Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner.
1893: We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then
1894: the direct Picard iteration A(x^n) x^{n+1} = b(x^n)
1896: 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
1897: 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
1898: different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).
1900: Level: intermediate
1902: .keywords: SNES, nonlinear, set, function
1904: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
1905: @*/
1906: 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)
1907: {
1909: DM dm;
1913: SNESGetDM(snes, &dm);
1914: DMSNESSetPicard(dm,b,J,ctx);
1915: SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
1916: SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
1917: return(0);
1918: }
1922: /*@C
1923: SNESGetPicard - Returns the context for the Picard iteration
1925: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
1927: Input Parameter:
1928: . snes - the SNES context
1930: Output Parameter:
1931: + r - the function (or NULL)
1932: . f - the function (or NULL); see SNESFunction for calling sequence details
1933: . Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
1934: . Pmat - the matrix from which the preconditioner will be constructed (or NULL)
1935: . J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
1936: - ctx - the function context (or NULL)
1938: Level: advanced
1940: .keywords: SNES, nonlinear, get, function
1942: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
1943: @*/
1944: 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)
1945: {
1947: DM dm;
1951: SNESGetFunction(snes,r,NULL,NULL);
1952: SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
1953: SNESGetDM(snes,&dm);
1954: DMSNESGetPicard(dm,f,J,ctx);
1955: return(0);
1956: }
1960: /*@C
1961: SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem
1963: Logically Collective on SNES
1965: Input Parameters:
1966: + snes - the SNES context
1967: . func - function evaluation routine
1968: - ctx - [optional] user-defined context for private data for the
1969: function evaluation routine (may be NULL)
1971: Calling sequence of func:
1972: $ func (SNES snes,Vec x,void *ctx);
1974: . f - function vector
1975: - ctx - optional user-defined function context
1977: Level: intermediate
1979: .keywords: SNES, nonlinear, set, function
1981: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
1982: @*/
1983: PetscErrorCode SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
1984: {
1987: if (func) snes->ops->computeinitialguess = func;
1988: if (ctx) snes->initialguessP = ctx;
1989: return(0);
1990: }
1992: /* --------------------------------------------------------------- */
1995: /*@C
1996: SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
1997: it assumes a zero right hand side.
1999: Logically Collective on SNES
2001: Input Parameter:
2002: . snes - the SNES context
2004: Output Parameter:
2005: . rhs - the right hand side vector or NULL if the right hand side vector is null
2007: Level: intermediate
2009: .keywords: SNES, nonlinear, get, function, right hand side
2011: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
2012: @*/
2013: PetscErrorCode SNESGetRhs(SNES snes,Vec *rhs)
2014: {
2018: *rhs = snes->vec_rhs;
2019: return(0);
2020: }
2024: /*@
2025: SNESComputeFunction - Calls the function that has been set with SNESSetFunction().
2027: Collective on SNES
2029: Input Parameters:
2030: + snes - the SNES context
2031: - x - input vector
2033: Output Parameter:
2034: . y - function vector, as set by SNESSetFunction()
2036: Notes:
2037: SNESComputeFunction() is typically used within nonlinear solvers
2038: implementations, so most users would not generally call this routine
2039: themselves.
2041: Level: developer
2043: .keywords: SNES, nonlinear, compute, function
2045: .seealso: SNESSetFunction(), SNESGetFunction()
2046: @*/
2047: PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y)
2048: {
2050: DM dm;
2051: DMSNES sdm;
2059: VecValidValues(x,2,PETSC_TRUE);
2061: SNESGetDM(snes,&dm);
2062: DMGetDMSNES(dm,&sdm);
2063: if (sdm->ops->computefunction) {
2064: PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2065: VecLockPush(x);
2066: PetscStackPush("SNES user function");
2067: (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2068: PetscStackPop;
2069: VecLockPop(x);
2070: PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2071: } else if (snes->vec_rhs) {
2072: MatMult(snes->jacobian, x, y);
2073: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2074: if (snes->vec_rhs) {
2075: VecAXPY(y,-1.0,snes->vec_rhs);
2076: }
2077: snes->nfuncs++;
2078: /*
2079: domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2080: propagate the value to all processes
2081: */
2082: if (snes->domainerror) {
2083: VecSetInf(y);
2084: }
2085: return(0);
2086: }
2090: /*@
2091: SNESComputeNGS - Calls the Gauss-Seidel function that has been set with SNESSetNGS().
2093: Collective on SNES
2095: Input Parameters:
2096: + snes - the SNES context
2097: . x - input vector
2098: - b - rhs vector
2100: Output Parameter:
2101: . x - new solution vector
2103: Notes:
2104: SNESComputeNGS() is typically used within composed nonlinear solver
2105: implementations, so most users would not generally call this routine
2106: themselves.
2108: Level: developer
2110: .keywords: SNES, nonlinear, compute, function
2112: .seealso: SNESSetNGS(), SNESComputeFunction()
2113: @*/
2114: PetscErrorCode SNESComputeNGS(SNES snes,Vec b,Vec x)
2115: {
2117: DM dm;
2118: DMSNES sdm;
2126: if (b) {VecValidValues(b,2,PETSC_TRUE);}
2127: PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2128: SNESGetDM(snes,&dm);
2129: DMGetDMSNES(dm,&sdm);
2130: if (sdm->ops->computegs) {
2131: if (b) {VecLockPush(b);}
2132: PetscStackPush("SNES user NGS");
2133: (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2134: PetscStackPop;
2135: if (b) {VecLockPop(b);}
2136: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2137: PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2138: return(0);
2139: }
2143: /*@
2144: SNESComputeJacobian - Computes the Jacobian matrix that has been set with SNESSetJacobian().
2146: Collective on SNES and Mat
2148: Input Parameters:
2149: + snes - the SNES context
2150: - x - input vector
2152: Output Parameters:
2153: + A - Jacobian matrix
2154: - B - optional preconditioning matrix
2156: Options Database Keys:
2157: + -snes_lag_preconditioner <lag>
2158: . -snes_lag_jacobian <lag>
2159: . -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2160: . -snes_compare_explicit_draw - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2161: . -snes_compare_explicit_contour - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2162: . -snes_compare_operator - Make the comparison options above use the operator instead of the preconditioning matrix
2163: . -snes_compare_coloring - Compute the finite differece Jacobian using coloring and display norms of difference
2164: . -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
2165: . -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2166: . -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2167: . -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2168: . -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
2169: - -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences
2172: Notes:
2173: Most users should not need to explicitly call this routine, as it
2174: is used internally within the nonlinear solvers.
2176: Level: developer
2178: .keywords: SNES, compute, Jacobian, matrix
2180: .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2181: @*/
2182: PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2183: {
2185: PetscBool flag;
2186: DM dm;
2187: DMSNES sdm;
2188: KSP ksp;
2194: VecValidValues(X,2,PETSC_TRUE);
2195: SNESGetDM(snes,&dm);
2196: DMGetDMSNES(dm,&sdm);
2198: if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_USER,"Must call SNESSetJacobian(), DMSNESSetJacobian(), DMDASNESSetJacobianLocal(), etc");
2200: /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
2202: if (snes->lagjacobian == -2) {
2203: snes->lagjacobian = -1;
2205: PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2206: } else if (snes->lagjacobian == -1) {
2207: PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2208: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2209: if (flag) {
2210: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2211: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2212: }
2213: return(0);
2214: } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2215: PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2216: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2217: if (flag) {
2218: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2219: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2220: }
2221: return(0);
2222: }
2223: if (snes->pc && snes->pcside == PC_LEFT) {
2224: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2225: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2226: return(0);
2227: }
2229: PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2230: VecLockPush(X);
2231: PetscStackPush("SNES user Jacobian function");
2232: (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2233: PetscStackPop;
2234: VecLockPop(X);
2235: PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);
2237: /* the next line ensures that snes->ksp exists */
2238: SNESGetKSP(snes,&ksp);
2239: if (snes->lagpreconditioner == -2) {
2240: PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2241: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2242: snes->lagpreconditioner = -1;
2243: } else if (snes->lagpreconditioner == -1) {
2244: PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2245: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2246: } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2247: PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2248: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2249: } else {
2250: PetscInfo(snes,"Rebuilding preconditioner\n");
2251: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2252: }
2254: /* make sure user returned a correct Jacobian and preconditioner */
2257: {
2258: PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
2259: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit",&flag,NULL);
2260: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",&flag_draw,NULL);
2261: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",&flag_contour,NULL);
2262: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_operator",&flag_operator,NULL);
2263: if (flag || flag_draw || flag_contour) {
2264: Mat Bexp_mine = NULL,Bexp,FDexp;
2265: PetscViewer vdraw,vstdout;
2266: PetscBool flg;
2267: if (flag_operator) {
2268: MatComputeExplicitOperator(A,&Bexp_mine);
2269: Bexp = Bexp_mine;
2270: } else {
2271: /* See if the preconditioning matrix can be viewed and added directly */
2272: PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2273: if (flg) Bexp = B;
2274: else {
2275: /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
2276: MatComputeExplicitOperator(B,&Bexp_mine);
2277: Bexp = Bexp_mine;
2278: }
2279: }
2280: MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);
2281: SNESComputeJacobianDefault(snes,X,FDexp,FDexp,NULL);
2282: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2283: if (flag_draw || flag_contour) {
2284: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2285: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2286: } else vdraw = NULL;
2287: PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator ? "Jacobian" : "preconditioning Jacobian");
2288: if (flag) {MatView(Bexp,vstdout);}
2289: if (vdraw) {MatView(Bexp,vdraw);}
2290: PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");
2291: if (flag) {MatView(FDexp,vstdout);}
2292: if (vdraw) {MatView(FDexp,vdraw);}
2293: MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);
2294: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");
2295: if (flag) {MatView(FDexp,vstdout);}
2296: if (vdraw) { /* Always use contour for the difference */
2297: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2298: MatView(FDexp,vdraw);
2299: PetscViewerPopFormat(vdraw);
2300: }
2301: if (flag_contour) {PetscViewerPopFormat(vdraw);}
2302: PetscViewerDestroy(&vdraw);
2303: MatDestroy(&Bexp_mine);
2304: MatDestroy(&FDexp);
2305: }
2306: }
2307: {
2308: PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
2309: PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
2310: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring",&flag,NULL);
2311: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_display",&flag_display,NULL);
2312: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",&flag_draw,NULL);
2313: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",&flag_contour,NULL);
2314: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",&flag_threshold,NULL);
2315: PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,NULL);
2316: PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,NULL);
2317: if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
2318: Mat Bfd;
2319: PetscViewer vdraw,vstdout;
2320: MatColoring coloring;
2321: ISColoring iscoloring;
2322: MatFDColoring matfdcoloring;
2323: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2324: void *funcctx;
2325: PetscReal norm1,norm2,normmax;
2327: MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2328: MatColoringCreate(Bfd,&coloring);
2329: MatColoringSetType(coloring,MATCOLORINGSL);
2330: MatColoringSetFromOptions(coloring);
2331: MatColoringApply(coloring,&iscoloring);
2332: MatColoringDestroy(&coloring);
2333: MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2334: MatFDColoringSetFromOptions(matfdcoloring);
2335: MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2336: ISColoringDestroy(&iscoloring);
2338: /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2339: SNESGetFunction(snes,NULL,&func,&funcctx);
2340: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2341: PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2342: PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2343: MatFDColoringSetFromOptions(matfdcoloring);
2344: MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2345: MatFDColoringDestroy(&matfdcoloring);
2347: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2348: if (flag_draw || flag_contour) {
2349: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2350: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2351: } else vdraw = NULL;
2352: PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2353: if (flag_display) {MatView(B,vstdout);}
2354: if (vdraw) {MatView(B,vdraw);}
2355: PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2356: if (flag_display) {MatView(Bfd,vstdout);}
2357: if (vdraw) {MatView(Bfd,vdraw);}
2358: MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2359: MatNorm(Bfd,NORM_1,&norm1);
2360: MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2361: MatNorm(Bfd,NORM_MAX,&normmax);
2362: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2363: if (flag_display) {MatView(Bfd,vstdout);}
2364: if (vdraw) { /* Always use contour for the difference */
2365: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2366: MatView(Bfd,vdraw);
2367: PetscViewerPopFormat(vdraw);
2368: }
2369: if (flag_contour) {PetscViewerPopFormat(vdraw);}
2371: if (flag_threshold) {
2372: PetscInt bs,rstart,rend,i;
2373: MatGetBlockSize(B,&bs);
2374: MatGetOwnershipRange(B,&rstart,&rend);
2375: for (i=rstart; i<rend; i++) {
2376: const PetscScalar *ba,*ca;
2377: const PetscInt *bj,*cj;
2378: PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
2379: PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0;
2380: MatGetRow(B,i,&bn,&bj,&ba);
2381: MatGetRow(Bfd,i,&cn,&cj,&ca);
2382: if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
2383: for (j=0; j<bn; j++) {
2384: PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2385: if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
2386: maxentrycol = bj[j];
2387: maxentry = PetscRealPart(ba[j]);
2388: }
2389: if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
2390: maxdiffcol = bj[j];
2391: maxdiff = PetscRealPart(ca[j]);
2392: }
2393: if (rdiff > maxrdiff) {
2394: maxrdiffcol = bj[j];
2395: maxrdiff = rdiff;
2396: }
2397: }
2398: if (maxrdiff > 1) {
2399: 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);
2400: for (j=0; j<bn; j++) {
2401: PetscReal rdiff;
2402: rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2403: if (rdiff > 1) {
2404: PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
2405: }
2406: }
2407: PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
2408: }
2409: MatRestoreRow(B,i,&bn,&bj,&ba);
2410: MatRestoreRow(Bfd,i,&cn,&cj,&ca);
2411: }
2412: }
2413: PetscViewerDestroy(&vdraw);
2414: MatDestroy(&Bfd);
2415: }
2416: }
2417: return(0);
2418: }
2420: /*MC
2421: SNESJacobianFunction - Function used to convey the nonlinear Jacobian of the function to be solved by SNES
2423: Synopsis:
2424: #include "petscsnes.h"
2425: PetscErrorCode SNESJacobianFunction(SNES snes,Vec x,Mat Amat,Mat Pmat,void *ctx);
2427: + x - input vector
2428: . Amat - the matrix that defines the (approximate) Jacobian
2429: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2430: - ctx - [optional] user-defined Jacobian context
2432: Level: intermediate
2434: .seealso: SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2435: M*/
2439: /*@C
2440: SNESSetJacobian - Sets the function to compute Jacobian as well as the
2441: location to store the matrix.
2443: Logically Collective on SNES and Mat
2445: Input Parameters:
2446: + snes - the SNES context
2447: . Amat - the matrix that defines the (approximate) Jacobian
2448: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2449: . J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
2450: - ctx - [optional] user-defined context for private data for the
2451: Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)
2453: Notes:
2454: If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2455: each matrix.
2457: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
2458: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
2460: If using SNESComputeJacobianDefaultColor() to assemble a Jacobian, the ctx argument
2461: must be a MatFDColoring.
2463: Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian. One common
2464: example is to use the "Picard linearization" which only differentiates through the highest order parts of each term.
2466: Level: beginner
2468: .keywords: SNES, nonlinear, set, Jacobian, matrix
2470: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J,
2471: SNESSetPicard(), SNESJacobianFunction
2472: @*/
2473: PetscErrorCode SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2474: {
2476: DM dm;
2484: SNESGetDM(snes,&dm);
2485: DMSNESSetJacobian(dm,J,ctx);
2486: if (Amat) {
2487: PetscObjectReference((PetscObject)Amat);
2488: MatDestroy(&snes->jacobian);
2490: snes->jacobian = Amat;
2491: }
2492: if (Pmat) {
2493: PetscObjectReference((PetscObject)Pmat);
2494: MatDestroy(&snes->jacobian_pre);
2496: snes->jacobian_pre = Pmat;
2497: }
2498: return(0);
2499: }
2503: /*@C
2504: SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2505: provided context for evaluating the Jacobian.
2507: Not Collective, but Mat object will be parallel if SNES object is
2509: Input Parameter:
2510: . snes - the nonlinear solver context
2512: Output Parameters:
2513: + Amat - location to stash (approximate) Jacobian matrix (or NULL)
2514: . Pmat - location to stash matrix used to compute the preconditioner (or NULL)
2515: . J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
2516: - ctx - location to stash Jacobian ctx (or NULL)
2518: Level: advanced
2520: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
2521: @*/
2522: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2523: {
2525: DM dm;
2526: DMSNES sdm;
2530: if (Amat) *Amat = snes->jacobian;
2531: if (Pmat) *Pmat = snes->jacobian_pre;
2532: SNESGetDM(snes,&dm);
2533: DMGetDMSNES(dm,&sdm);
2534: if (J) *J = sdm->ops->computejacobian;
2535: if (ctx) *ctx = sdm->jacobianctx;
2536: return(0);
2537: }
2541: /*@
2542: SNESSetUp - Sets up the internal data structures for the later use
2543: of a nonlinear solver.
2545: Collective on SNES
2547: Input Parameters:
2548: . snes - the SNES context
2550: Notes:
2551: For basic use of the SNES solvers the user need not explicitly call
2552: SNESSetUp(), since these actions will automatically occur during
2553: the call to SNESSolve(). However, if one wishes to control this
2554: phase separately, SNESSetUp() should be called after SNESCreate()
2555: and optional routines of the form SNESSetXXX(), but before SNESSolve().
2557: Level: advanced
2559: .keywords: SNES, nonlinear, setup
2561: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
2562: @*/
2563: PetscErrorCode SNESSetUp(SNES snes)
2564: {
2566: DM dm;
2567: DMSNES sdm;
2568: SNESLineSearch linesearch, pclinesearch;
2569: void *lsprectx,*lspostctx;
2570: PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
2571: PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
2572: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2573: Vec f,fpc;
2574: void *funcctx;
2575: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2576: void *jacctx,*appctx;
2577: Mat j,jpre;
2581: if (snes->setupcalled) return(0);
2583: if (!((PetscObject)snes)->type_name) {
2584: SNESSetType(snes,SNESNEWTONLS);
2585: }
2587: SNESGetFunction(snes,&snes->vec_func,NULL,NULL);
2589: SNESGetDM(snes,&dm);
2590: DMGetDMSNES(dm,&sdm);
2591: if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
2592: if (!sdm->ops->computejacobian) {
2593: DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
2594: }
2595: if (!snes->vec_func) {
2596: DMCreateGlobalVector(dm,&snes->vec_func);
2597: }
2599: if (!snes->ksp) {
2600: SNESGetKSP(snes, &snes->ksp);
2601: }
2603: if (!snes->linesearch) {
2604: SNESGetLineSearch(snes, &snes->linesearch);
2605: }
2606: SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);
2608: if (snes->pc && (snes->pcside == PC_LEFT)) {
2609: snes->mf = PETSC_TRUE;
2610: snes->mf_operator = PETSC_FALSE;
2611: }
2613: if (snes->pc) {
2614: /* copy the DM over */
2615: SNESGetDM(snes,&dm);
2616: SNESSetDM(snes->pc,dm);
2618: SNESGetFunction(snes,&f,&func,&funcctx);
2619: VecDuplicate(f,&fpc);
2620: SNESSetFunction(snes->pc,fpc,func,funcctx);
2621: SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
2622: SNESSetJacobian(snes->pc,j,jpre,jac,jacctx);
2623: SNESGetApplicationContext(snes,&appctx);
2624: SNESSetApplicationContext(snes->pc,appctx);
2625: VecDestroy(&fpc);
2627: /* copy the function pointers over */
2628: PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)snes->pc);
2630: /* default to 1 iteration */
2631: SNESSetTolerances(snes->pc,0.0,0.0,0.0,1,snes->pc->max_funcs);
2632: if (snes->pcside==PC_RIGHT) {
2633: SNESSetNormSchedule(snes->pc,SNES_NORM_FINAL_ONLY);
2634: } else {
2635: SNESSetNormSchedule(snes->pc,SNES_NORM_NONE);
2636: }
2637: SNESSetFromOptions(snes->pc);
2639: /* copy the line search context over */
2640: SNESGetLineSearch(snes,&linesearch);
2641: SNESGetLineSearch(snes->pc,&pclinesearch);
2642: SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
2643: SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
2644: SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
2645: SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
2646: PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
2647: }
2648: if (snes->mf) {
2649: SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
2650: }
2651: if (snes->ops->usercompute && !snes->user) {
2652: (*snes->ops->usercompute)(snes,(void**)&snes->user);
2653: }
2655: snes->jac_iter = 0;
2656: snes->pre_iter = 0;
2658: if (snes->ops->setup) {
2659: (*snes->ops->setup)(snes);
2660: }
2662: if (snes->pc && (snes->pcside == PC_LEFT)) {
2663: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2664: SNESGetLineSearch(snes,&linesearch);
2665: SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
2666: }
2667: }
2669: snes->setupcalled = PETSC_TRUE;
2670: return(0);
2671: }
2675: /*@
2676: SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats
2678: Collective on SNES
2680: Input Parameter:
2681: . snes - iterative context obtained from SNESCreate()
2683: Level: intermediate
2685: Notes: Also calls the application context destroy routine set with SNESSetComputeApplicationContext()
2687: .keywords: SNES, destroy
2689: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
2690: @*/
2691: PetscErrorCode SNESReset(SNES snes)
2692: {
2697: if (snes->ops->userdestroy && snes->user) {
2698: (*snes->ops->userdestroy)((void**)&snes->user);
2699: snes->user = NULL;
2700: }
2701: if (snes->pc) {
2702: SNESReset(snes->pc);
2703: }
2705: if (snes->ops->reset) {
2706: (*snes->ops->reset)(snes);
2707: }
2708: if (snes->ksp) {
2709: KSPReset(snes->ksp);
2710: }
2712: if (snes->linesearch) {
2713: SNESLineSearchReset(snes->linesearch);
2714: }
2716: VecDestroy(&snes->vec_rhs);
2717: VecDestroy(&snes->vec_sol);
2718: VecDestroy(&snes->vec_sol_update);
2719: VecDestroy(&snes->vec_func);
2720: MatDestroy(&snes->jacobian);
2721: MatDestroy(&snes->jacobian_pre);
2722: VecDestroyVecs(snes->nwork,&snes->work);
2723: VecDestroyVecs(snes->nvwork,&snes->vwork);
2725: snes->nwork = snes->nvwork = 0;
2726: snes->setupcalled = PETSC_FALSE;
2727: return(0);
2728: }
2732: /*@
2733: SNESDestroy - Destroys the nonlinear solver context that was created
2734: with SNESCreate().
2736: Collective on SNES
2738: Input Parameter:
2739: . snes - the SNES context
2741: Level: beginner
2743: .keywords: SNES, nonlinear, destroy
2745: .seealso: SNESCreate(), SNESSolve()
2746: @*/
2747: PetscErrorCode SNESDestroy(SNES *snes)
2748: {
2752: if (!*snes) return(0);
2754: if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; return(0);}
2756: SNESReset((*snes));
2757: SNESDestroy(&(*snes)->pc);
2759: /* if memory was published with SAWs then destroy it */
2760: PetscObjectSAWsViewOff((PetscObject)*snes);
2761: if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}
2763: DMDestroy(&(*snes)->dm);
2764: KSPDestroy(&(*snes)->ksp);
2765: SNESLineSearchDestroy(&(*snes)->linesearch);
2767: PetscFree((*snes)->kspconvctx);
2768: if ((*snes)->ops->convergeddestroy) {
2769: (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
2770: }
2771: if ((*snes)->conv_malloc) {
2772: PetscFree((*snes)->conv_hist);
2773: PetscFree((*snes)->conv_hist_its);
2774: }
2775: SNESMonitorCancel((*snes));
2776: PetscHeaderDestroy(snes);
2777: return(0);
2778: }
2780: /* ----------- Routines to set solver parameters ---------- */
2784: /*@
2785: SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
2787: Logically Collective on SNES
2789: Input Parameters:
2790: + snes - the SNES context
2791: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2792: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
2794: Options Database Keys:
2795: . -snes_lag_preconditioner <lag>
2797: Notes:
2798: The default is 1
2799: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2800: If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
2802: Level: intermediate
2804: .keywords: SNES, nonlinear, set, convergence, tolerances
2806: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
2808: @*/
2809: PetscErrorCode SNESSetLagPreconditioner(SNES snes,PetscInt lag)
2810: {
2813: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2814: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2816: snes->lagpreconditioner = lag;
2817: return(0);
2818: }
2822: /*@
2823: SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does
2825: Logically Collective on SNES
2827: Input Parameters:
2828: + snes - the SNES context
2829: - steps - the number of refinements to do, defaults to 0
2831: Options Database Keys:
2832: . -snes_grid_sequence <steps>
2834: Level: intermediate
2836: Notes:
2837: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
2839: .keywords: SNES, nonlinear, set, convergence, tolerances
2841: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetGridSequence()
2843: @*/
2844: PetscErrorCode SNESSetGridSequence(SNES snes,PetscInt steps)
2845: {
2849: snes->gridsequence = steps;
2850: return(0);
2851: }
2855: /*@
2856: SNESGetGridSequence - gets the number of steps of grid sequencing that SNES does
2858: Logically Collective on SNES
2860: Input Parameter:
2861: . snes - the SNES context
2863: Output Parameter:
2864: . steps - the number of refinements to do, defaults to 0
2866: Options Database Keys:
2867: . -snes_grid_sequence <steps>
2869: Level: intermediate
2871: Notes:
2872: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
2874: .keywords: SNES, nonlinear, set, convergence, tolerances
2876: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESSetGridSequence()
2878: @*/
2879: PetscErrorCode SNESGetGridSequence(SNES snes,PetscInt *steps)
2880: {
2883: *steps = snes->gridsequence;
2884: return(0);
2885: }
2889: /*@
2890: SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
2892: Not Collective
2894: Input Parameter:
2895: . snes - the SNES context
2897: Output Parameter:
2898: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2899: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
2901: Options Database Keys:
2902: . -snes_lag_preconditioner <lag>
2904: Notes:
2905: The default is 1
2906: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2908: Level: intermediate
2910: .keywords: SNES, nonlinear, set, convergence, tolerances
2912: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
2914: @*/
2915: PetscErrorCode SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
2916: {
2919: *lag = snes->lagpreconditioner;
2920: return(0);
2921: }
2925: /*@
2926: SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
2927: often the preconditioner is rebuilt.
2929: Logically Collective on SNES
2931: Input Parameters:
2932: + snes - the SNES context
2933: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2934: the Jacobian is built etc. -2 means rebuild at next chance but then never again
2936: Options Database Keys:
2937: . -snes_lag_jacobian <lag>
2939: Notes:
2940: The default is 1
2941: The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2942: 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
2943: at the next Newton step but never again (unless it is reset to another value)
2945: Level: intermediate
2947: .keywords: SNES, nonlinear, set, convergence, tolerances
2949: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
2951: @*/
2952: PetscErrorCode SNESSetLagJacobian(SNES snes,PetscInt lag)
2953: {
2956: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2957: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2959: snes->lagjacobian = lag;
2960: return(0);
2961: }
2965: /*@
2966: SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
2968: Not Collective
2970: Input Parameter:
2971: . snes - the SNES context
2973: Output Parameter:
2974: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2975: the Jacobian is built etc.
2977: Options Database Keys:
2978: . -snes_lag_jacobian <lag>
2980: Notes:
2981: The default is 1
2982: The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2984: Level: intermediate
2986: .keywords: SNES, nonlinear, set, convergence, tolerances
2988: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
2990: @*/
2991: PetscErrorCode SNESGetLagJacobian(SNES snes,PetscInt *lag)
2992: {
2995: *lag = snes->lagjacobian;
2996: return(0);
2997: }
3001: /*@
3002: SNESSetLagJacobianPersists - Set whether or not the Jacobian lagging persists through multiple solves
3004: Logically collective on SNES
3006: Input Parameter:
3007: + snes - the SNES context
3008: - flg - jacobian lagging persists if true
3010: Options Database Keys:
3011: . -snes_lag_jacobian_persists <flg>
3013: Notes: This is useful both for nonlinear preconditioning, where it's appropriate to have the Jacobian be stale by
3014: several solves, and for implicit time-stepping, where Jacobian lagging in the inner nonlinear solve over several
3015: timesteps may present huge efficiency gains.
3017: Level: developer
3019: .keywords: SNES, nonlinear, lag
3021: .seealso: SNESSetLagPreconditionerPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()
3023: @*/
3024: PetscErrorCode SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
3025: {
3029: snes->lagjac_persist = flg;
3030: return(0);
3031: }
3035: /*@
3036: SNESSetLagPreconditionerPersists - Set whether or not the preconditioner lagging persists through multiple solves
3038: Logically Collective on SNES
3040: Input Parameter:
3041: + snes - the SNES context
3042: - flg - preconditioner lagging persists if true
3044: Options Database Keys:
3045: . -snes_lag_jacobian_persists <flg>
3047: Notes: This is useful both for nonlinear preconditioning, where it's appropriate to have the preconditioner be stale
3048: by several solves, and for implicit time-stepping, where preconditioner lagging in the inner nonlinear solve over
3049: several timesteps may present huge efficiency gains.
3051: Level: developer
3053: .keywords: SNES, nonlinear, lag
3055: .seealso: SNESSetLagJacobianPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()
3057: @*/
3058: PetscErrorCode SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3059: {
3063: snes->lagpre_persist = flg;
3064: return(0);
3065: }
3069: /*@
3070: SNESSetTolerances - Sets various parameters used in convergence tests.
3072: Logically Collective on SNES
3074: Input Parameters:
3075: + snes - the SNES context
3076: . abstol - absolute convergence tolerance
3077: . rtol - relative convergence tolerance
3078: . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x ||
3079: . maxit - maximum number of iterations
3080: - maxf - maximum number of function evaluations
3082: Options Database Keys:
3083: + -snes_atol <abstol> - Sets abstol
3084: . -snes_rtol <rtol> - Sets rtol
3085: . -snes_stol <stol> - Sets stol
3086: . -snes_max_it <maxit> - Sets maxit
3087: - -snes_max_funcs <maxf> - Sets maxf
3089: Notes:
3090: The default maximum number of iterations is 50.
3091: The default maximum number of function evaluations is 1000.
3093: Level: intermediate
3095: .keywords: SNES, nonlinear, set, convergence, tolerances
3097: .seealso: SNESSetTrustRegionTolerance()
3098: @*/
3099: PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3100: {
3109: if (abstol != PETSC_DEFAULT) {
3110: if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3111: snes->abstol = abstol;
3112: }
3113: if (rtol != PETSC_DEFAULT) {
3114: 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);
3115: snes->rtol = rtol;
3116: }
3117: if (stol != PETSC_DEFAULT) {
3118: if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3119: snes->stol = stol;
3120: }
3121: if (maxit != PETSC_DEFAULT) {
3122: if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3123: snes->max_its = maxit;
3124: }
3125: if (maxf != PETSC_DEFAULT) {
3126: if (maxf < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
3127: snes->max_funcs = maxf;
3128: }
3129: snes->tolerancesset = PETSC_TRUE;
3130: return(0);
3131: }
3135: /*@
3136: SNESGetTolerances - Gets various parameters used in convergence tests.
3138: Not Collective
3140: Input Parameters:
3141: + snes - the SNES context
3142: . atol - absolute convergence tolerance
3143: . rtol - relative convergence tolerance
3144: . stol - convergence tolerance in terms of the norm
3145: of the change in the solution between steps
3146: . maxit - maximum number of iterations
3147: - maxf - maximum number of function evaluations
3149: Notes:
3150: The user can specify NULL for any parameter that is not needed.
3152: Level: intermediate
3154: .keywords: SNES, nonlinear, get, convergence, tolerances
3156: .seealso: SNESSetTolerances()
3157: @*/
3158: PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3159: {
3162: if (atol) *atol = snes->abstol;
3163: if (rtol) *rtol = snes->rtol;
3164: if (stol) *stol = snes->stol;
3165: if (maxit) *maxit = snes->max_its;
3166: if (maxf) *maxf = snes->max_funcs;
3167: return(0);
3168: }
3172: /*@
3173: SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
3175: Logically Collective on SNES
3177: Input Parameters:
3178: + snes - the SNES context
3179: - tol - tolerance
3181: Options Database Key:
3182: . -snes_trtol <tol> - Sets tol
3184: Level: intermediate
3186: .keywords: SNES, nonlinear, set, trust region, tolerance
3188: .seealso: SNESSetTolerances()
3189: @*/
3190: PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3191: {
3195: snes->deltatol = tol;
3196: return(0);
3197: }
3199: /*
3200: Duplicate the lg monitors for SNES from KSP; for some reason with
3201: dynamic libraries things don't work under Sun4 if we just use
3202: macros instead of functions
3203: */
3206: PetscErrorCode SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,PetscObject *objs)
3207: {
3212: KSPMonitorLGResidualNorm((KSP)snes,it,norm,objs);
3213: return(0);
3214: }
3218: PetscErrorCode SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscObject **draw)
3219: {
3223: KSPMonitorLGResidualNormCreate(host,label,x,y,m,n,draw);
3224: return(0);
3225: }
3229: PetscErrorCode SNESMonitorLGDestroy(PetscObject **objs)
3230: {
3234: KSPMonitorLGResidualNormDestroy(objs);
3235: return(0);
3236: }
3238: extern PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3241: PetscErrorCode SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3242: {
3243: PetscDrawLG lg;
3244: PetscErrorCode ierr;
3245: PetscReal x,y,per;
3246: PetscViewer v = (PetscViewer)monctx;
3247: static PetscReal prev; /* should be in the context */
3248: PetscDraw draw;
3251: PetscViewerDrawGetDrawLG(v,0,&lg);
3252: if (!n) {PetscDrawLGReset(lg);}
3253: PetscDrawLGGetDraw(lg,&draw);
3254: PetscDrawSetTitle(draw,"Residual norm");
3255: x = (PetscReal)n;
3256: if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3257: else y = -15.0;
3258: PetscDrawLGAddPoint(lg,&x,&y);
3259: if (n < 20 || !(n % 5)) {
3260: PetscDrawLGDraw(lg);
3261: }
3263: PetscViewerDrawGetDrawLG(v,1,&lg);
3264: if (!n) {PetscDrawLGReset(lg);}
3265: PetscDrawLGGetDraw(lg,&draw);
3266: PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3267: SNESMonitorRange_Private(snes,n,&per);
3268: x = (PetscReal)n;
3269: y = 100.0*per;
3270: PetscDrawLGAddPoint(lg,&x,&y);
3271: if (n < 20 || !(n % 5)) {
3272: PetscDrawLGDraw(lg);
3273: }
3275: PetscViewerDrawGetDrawLG(v,2,&lg);
3276: if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3277: PetscDrawLGGetDraw(lg,&draw);
3278: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3279: x = (PetscReal)n;
3280: y = (prev - rnorm)/prev;
3281: PetscDrawLGAddPoint(lg,&x,&y);
3282: if (n < 20 || !(n % 5)) {
3283: PetscDrawLGDraw(lg);
3284: }
3286: PetscViewerDrawGetDrawLG(v,3,&lg);
3287: if (!n) {PetscDrawLGReset(lg);}
3288: PetscDrawLGGetDraw(lg,&draw);
3289: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3290: x = (PetscReal)n;
3291: y = (prev - rnorm)/(prev*per);
3292: if (n > 2) { /*skip initial crazy value */
3293: PetscDrawLGAddPoint(lg,&x,&y);
3294: }
3295: if (n < 20 || !(n % 5)) {
3296: PetscDrawLGDraw(lg);
3297: }
3298: prev = rnorm;
3299: return(0);
3300: }
3304: /*@
3305: SNESMonitor - runs the user provided monitor routines, if they exist
3307: Collective on SNES
3309: Input Parameters:
3310: + snes - nonlinear solver context obtained from SNESCreate()
3311: . iter - iteration number
3312: - rnorm - relative norm of the residual
3314: Notes:
3315: This routine is called by the SNES implementations.
3316: It does not typically need to be called by the user.
3318: Level: developer
3320: .seealso: SNESMonitorSet()
3321: @*/
3322: PetscErrorCode SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3323: {
3325: PetscInt i,n = snes->numbermonitors;
3328: VecLockPush(snes->vec_sol);
3329: for (i=0; i<n; i++) {
3330: (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3331: }
3332: VecLockPop(snes->vec_sol);
3333: return(0);
3334: }
3336: /* ------------ Routines to set performance monitoring options ----------- */
3338: /*MC
3339: SNESMonitorFunction - functional form passed to SNESMonitorSet() to monitor convergence of nonlinear solver
3341: Synopsis:
3342: #include <petscsnes.h>
3343: $ PetscErrorCode SNESMonitorFunction(SNES snes,PetscInt its, PetscReal norm,void *mctx)
3345: + snes - the SNES context
3346: . its - iteration number
3347: . norm - 2-norm function value (may be estimated)
3348: - mctx - [optional] monitoring context
3350: Level: advanced
3352: .seealso: SNESMonitorSet(), SNESMonitorGet()
3353: M*/
3357: /*@C
3358: SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3359: iteration of the nonlinear solver to display the iteration's
3360: progress.
3362: Logically Collective on SNES
3364: Input Parameters:
3365: + snes - the SNES context
3366: . f - the monitor function, see SNESMonitorFunction for the calling sequence
3367: . mctx - [optional] user-defined context for private data for the
3368: monitor routine (use NULL if no context is desired)
3369: - monitordestroy - [optional] routine that frees monitor context
3370: (may be NULL)
3372: Options Database Keys:
3373: + -snes_monitor - sets SNESMonitorDefault()
3374: . -snes_monitor_lg_residualnorm - sets line graph monitor,
3375: uses SNESMonitorLGCreate()
3376: - -snes_monitor_cancel - cancels all monitors that have
3377: been hardwired into a code by
3378: calls to SNESMonitorSet(), but
3379: does not cancel those set via
3380: the options database.
3382: Notes:
3383: Several different monitoring routines may be set by calling
3384: SNESMonitorSet() multiple times; all will be called in the
3385: order in which they were set.
3387: Fortran notes: Only a single monitor function can be set for each SNES object
3389: Level: intermediate
3391: .keywords: SNES, nonlinear, set, monitor
3393: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3394: @*/
3395: PetscErrorCode SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3396: {
3397: PetscInt i;
3402: if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3403: for (i=0; i<snes->numbermonitors;i++) {
3404: if (f == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) {
3405: if (monitordestroy) {
3406: (*monitordestroy)(&mctx);
3407: }
3408: return(0);
3409: }
3410: }
3411: snes->monitor[snes->numbermonitors] = f;
3412: snes->monitordestroy[snes->numbermonitors] = monitordestroy;
3413: snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3414: return(0);
3415: }
3419: /*@
3420: SNESMonitorCancel - Clears all the monitor functions for a SNES object.
3422: Logically Collective on SNES
3424: Input Parameters:
3425: . snes - the SNES context
3427: Options Database Key:
3428: . -snes_monitor_cancel - cancels all monitors that have been hardwired
3429: into a code by calls to SNESMonitorSet(), but does not cancel those
3430: set via the options database
3432: Notes:
3433: There is no way to clear one specific monitor from a SNES object.
3435: Level: intermediate
3437: .keywords: SNES, nonlinear, set, monitor
3439: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3440: @*/
3441: PetscErrorCode SNESMonitorCancel(SNES snes)
3442: {
3444: PetscInt i;
3448: for (i=0; i<snes->numbermonitors; i++) {
3449: if (snes->monitordestroy[i]) {
3450: (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3451: }
3452: }
3453: snes->numbermonitors = 0;
3454: return(0);
3455: }
3457: /*MC
3458: SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver
3460: Synopsis:
3461: #include <petscsnes.h>
3462: $ PetscErrorCode SNESConvergenceTest(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
3464: + snes - the SNES context
3465: . it - current iteration (0 is the first and is before any Newton step)
3466: . cctx - [optional] convergence context
3467: . reason - reason for convergence/divergence
3468: . xnorm - 2-norm of current iterate
3469: . gnorm - 2-norm of current step
3470: - f - 2-norm of function
3472: Level: intermediate
3474: .seealso: SNESSetConvergenceTest(), SNESGetConvergenceTest()
3475: M*/
3479: /*@C
3480: SNESSetConvergenceTest - Sets the function that is to be used
3481: to test for convergence of the nonlinear iterative solution.
3483: Logically Collective on SNES
3485: Input Parameters:
3486: + snes - the SNES context
3487: . SNESConvergenceTestFunction - routine to test for convergence
3488: . cctx - [optional] context for private data for the convergence routine (may be NULL)
3489: - destroy - [optional] destructor for the context (may be NULL; NULL_FUNCTION in Fortran)
3491: Level: advanced
3493: .keywords: SNES, nonlinear, set, convergence, test
3495: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3496: @*/
3497: PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3498: {
3503: if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
3504: if (snes->ops->convergeddestroy) {
3505: (*snes->ops->convergeddestroy)(snes->cnvP);
3506: }
3507: snes->ops->converged = SNESConvergenceTestFunction;
3508: snes->ops->convergeddestroy = destroy;
3509: snes->cnvP = cctx;
3510: return(0);
3511: }
3515: /*@
3516: SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
3518: Not Collective
3520: Input Parameter:
3521: . snes - the SNES context
3523: Output Parameter:
3524: . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
3525: manual pages for the individual convergence tests for complete lists
3527: Level: intermediate
3529: Notes: Can only be called after the call the SNESSolve() is complete.
3531: .keywords: SNES, nonlinear, set, convergence, test
3533: .seealso: SNESSetConvergenceTest(), SNESConvergedReason
3534: @*/
3535: PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
3536: {
3540: *reason = snes->reason;
3541: return(0);
3542: }
3546: /*@
3547: SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
3549: Logically Collective on SNES
3551: Input Parameters:
3552: + snes - iterative context obtained from SNESCreate()
3553: . a - array to hold history, this array will contain the function norms computed at each step
3554: . its - integer array holds the number of linear iterations for each solve.
3555: . na - size of a and its
3556: - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
3557: else it continues storing new values for new nonlinear solves after the old ones
3559: Notes:
3560: If 'a' and 'its' are NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
3561: default array of length 10000 is allocated.
3563: This routine is useful, e.g., when running a code for purposes
3564: of accurate performance monitoring, when no I/O should be done
3565: during the section of code that is being timed.
3567: Level: intermediate
3569: .keywords: SNES, set, convergence, history
3571: .seealso: SNESGetConvergenceHistory()
3573: @*/
3574: PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
3575: {
3582: if (!a) {
3583: if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
3584: PetscCalloc1(na,&a);
3585: PetscCalloc1(na,&its);
3587: snes->conv_malloc = PETSC_TRUE;
3588: }
3589: snes->conv_hist = a;
3590: snes->conv_hist_its = its;
3591: snes->conv_hist_max = na;
3592: snes->conv_hist_len = 0;
3593: snes->conv_hist_reset = reset;
3594: return(0);
3595: }
3597: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3598: #include <engine.h> /* MATLAB include file */
3599: #include <mex.h> /* MATLAB include file */
3603: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
3604: {
3605: mxArray *mat;
3606: PetscInt i;
3607: PetscReal *ar;
3610: mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
3611: ar = (PetscReal*) mxGetData(mat);
3612: for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
3613: PetscFunctionReturn(mat);
3614: }
3615: #endif
3619: /*@C
3620: SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
3622: Not Collective
3624: Input Parameter:
3625: . snes - iterative context obtained from SNESCreate()
3627: Output Parameters:
3628: . a - array to hold history
3629: . its - integer array holds the number of linear iterations (or
3630: negative if not converged) for each solve.
3631: - na - size of a and its
3633: Notes:
3634: The calling sequence for this routine in Fortran is
3635: $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
3637: This routine is useful, e.g., when running a code for purposes
3638: of accurate performance monitoring, when no I/O should be done
3639: during the section of code that is being timed.
3641: Level: intermediate
3643: .keywords: SNES, get, convergence, history
3645: .seealso: SNESSetConvergencHistory()
3647: @*/
3648: PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
3649: {
3652: if (a) *a = snes->conv_hist;
3653: if (its) *its = snes->conv_hist_its;
3654: if (na) *na = snes->conv_hist_len;
3655: return(0);
3656: }
3660: /*@C
3661: SNESSetUpdate - Sets the general-purpose update function called
3662: at the beginning of every iteration of the nonlinear solve. Specifically
3663: it is called just before the Jacobian is "evaluated".
3665: Logically Collective on SNES
3667: Input Parameters:
3668: . snes - The nonlinear solver context
3669: . func - The function
3671: Calling sequence of func:
3672: . func (SNES snes, PetscInt step);
3674: . step - The current step of the iteration
3676: Level: advanced
3678: 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()
3679: This is not used by most users.
3681: .keywords: SNES, update
3683: .seealso SNESSetJacobian(), SNESSolve()
3684: @*/
3685: PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
3686: {
3689: snes->ops->update = func;
3690: return(0);
3691: }
3695: /*
3696: SNESScaleStep_Private - Scales a step so that its length is less than the
3697: positive parameter delta.
3699: Input Parameters:
3700: + snes - the SNES context
3701: . y - approximate solution of linear system
3702: . fnorm - 2-norm of current function
3703: - delta - trust region size
3705: Output Parameters:
3706: + gpnorm - predicted function norm at the new point, assuming local
3707: linearization. The value is zero if the step lies within the trust
3708: region, and exceeds zero otherwise.
3709: - ynorm - 2-norm of the step
3711: Note:
3712: For non-trust region methods such as SNESNEWTONLS, the parameter delta
3713: is set to be the maximum allowable step size.
3715: .keywords: SNES, nonlinear, scale, step
3716: */
3717: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
3718: {
3719: PetscReal nrm;
3720: PetscScalar cnorm;
3728: VecNorm(y,NORM_2,&nrm);
3729: if (nrm > *delta) {
3730: nrm = *delta/nrm;
3731: *gpnorm = (1.0 - nrm)*(*fnorm);
3732: cnorm = nrm;
3733: VecScale(y,cnorm);
3734: *ynorm = *delta;
3735: } else {
3736: *gpnorm = 0.0;
3737: *ynorm = nrm;
3738: }
3739: return(0);
3740: }
3744: /*@
3745: SNESReasonView - Displays the reason a SNES solve converged or diverged to a viewer
3747: Collective on SNES
3749: Parameter:
3750: + snes - iterative context obtained from SNESCreate()
3751: - viewer - the viewer to display the reason
3754: Options Database Keys:
3755: . -snes_converged_reason - print reason for converged or diverged, also prints number of iterations
3757: Level: beginner
3759: .keywords: SNES, solve, linear system
3761: .seealso: SNESCreate(), SNESSetUp(), SNESDestroy(), SNESSetTolerances(), SNESConvergedDefault()
3763: @*/
3764: PetscErrorCode SNESReasonView(SNES snes,PetscViewer viewer)
3765: {
3767: PetscBool isAscii;
3770: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
3771: if (isAscii) {
3772: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
3773: if (snes->reason > 0) {
3774: if (((PetscObject) snes)->prefix) {
3775: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3776: } else {
3777: PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3778: }
3779: } else {
3780: if (((PetscObject) snes)->prefix) {
3781: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3782: } else {
3783: PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3784: }
3785: }
3786: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
3787: }
3788: return(0);
3789: }
3793: /*@C
3794: SNESReasonViewFromOptions - Processes command line options to determine if/how a SNESReason is to be viewed.
3796: Collective on SNES
3798: Input Parameters:
3799: . snes - the SNES object
3801: Level: intermediate
3803: @*/
3804: PetscErrorCode SNESReasonViewFromOptions(SNES snes)
3805: {
3806: PetscErrorCode ierr;
3807: PetscViewer viewer;
3808: PetscBool flg;
3809: static PetscBool incall = PETSC_FALSE;
3810: PetscViewerFormat format;
3813: if (incall) return(0);
3814: incall = PETSC_TRUE;
3815: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
3816: if (flg) {
3817: PetscViewerPushFormat(viewer,format);
3818: SNESReasonView(snes,viewer);
3819: PetscViewerPopFormat(viewer);
3820: PetscViewerDestroy(&viewer);
3821: }
3822: incall = PETSC_FALSE;
3823: return(0);
3824: }
3828: /*@C
3829: SNESSolve - Solves a nonlinear system F(x) = b.
3830: Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
3832: Collective on SNES
3834: Input Parameters:
3835: + snes - the SNES context
3836: . b - the constant part of the equation F(x) = b, or NULL to use zero.
3837: - x - the solution vector.
3839: Notes:
3840: The user should initialize the vector,x, with the initial guess
3841: for the nonlinear solve prior to calling SNESSolve. In particular,
3842: to employ an initial guess of zero, the user should explicitly set
3843: this vector to zero by calling VecSet().
3845: Level: beginner
3847: .keywords: SNES, nonlinear, solve
3849: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
3850: @*/
3851: PetscErrorCode SNESSolve(SNES snes,Vec b,Vec x)
3852: {
3853: PetscErrorCode ierr;
3854: PetscBool flg;
3855: PetscInt grid;
3856: Vec xcreated = NULL;
3857: DM dm;
3866: if (!x) {
3867: SNESGetDM(snes,&dm);
3868: DMCreateGlobalVector(dm,&xcreated);
3869: x = xcreated;
3870: }
3871: SNESViewFromOptions(snes,NULL,"-snes_view_pre");
3873: for (grid=0; grid<snes->gridsequence; grid++) {PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));}
3874: for (grid=0; grid<snes->gridsequence+1; grid++) {
3876: /* set solution vector */
3877: if (!grid) {PetscObjectReference((PetscObject)x);}
3878: VecDestroy(&snes->vec_sol);
3879: snes->vec_sol = x;
3880: SNESGetDM(snes,&dm);
3882: /* set affine vector if provided */
3883: if (b) { PetscObjectReference((PetscObject)b); }
3884: VecDestroy(&snes->vec_rhs);
3885: snes->vec_rhs = b;
3887: if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
3888: if (snes->vec_rhs == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
3889: if (!snes->vec_sol_update /* && snes->vec_sol */) {
3890: VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
3891: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
3892: }
3893: DMShellSetGlobalVector(dm,snes->vec_sol);
3894: SNESSetUp(snes);
3896: if (!grid) {
3897: if (snes->ops->computeinitialguess) {
3898: (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
3899: }
3900: }
3902: if (snes->conv_hist_reset) snes->conv_hist_len = 0;
3903: if (snes->counters_reset) {snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;}
3905: PetscLogEventBegin(SNES_Solve,snes,0,0,0);
3906: (*snes->ops->solve)(snes);
3907: PetscLogEventEnd(SNES_Solve,snes,0,0,0);
3908: if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
3909: snes->domainerror = PETSC_FALSE; /* clear the flag if it has been set */
3911: if (snes->lagjac_persist) snes->jac_iter += snes->iter;
3912: if (snes->lagpre_persist) snes->pre_iter += snes->iter;
3914: flg = PETSC_FALSE;
3915: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,NULL);
3916: if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
3917: SNESReasonViewFromOptions(snes);
3919: if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
3920: if (grid < snes->gridsequence) {
3921: DM fine;
3922: Vec xnew;
3923: Mat interp;
3925: DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
3926: if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
3927: DMCreateInterpolation(snes->dm,fine,&interp,NULL);
3928: DMCreateGlobalVector(fine,&xnew);
3929: MatInterpolate(interp,x,xnew);
3930: DMInterpolate(snes->dm,interp,fine);
3931: MatDestroy(&interp);
3932: x = xnew;
3934: SNESReset(snes);
3935: SNESSetDM(snes,fine);
3936: DMDestroy(&fine);
3937: PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
3938: }
3939: }
3940: SNESViewFromOptions(snes,NULL,"-snes_view");
3941: VecViewFromOptions(snes->vec_sol,(PetscObject)snes,"-snes_view_solution");
3943: VecDestroy(&xcreated);
3944: PetscObjectSAWsBlock((PetscObject)snes);
3945: return(0);
3946: }
3948: /* --------- Internal routines for SNES Package --------- */
3952: /*@C
3953: SNESSetType - Sets the method for the nonlinear solver.
3955: Collective on SNES
3957: Input Parameters:
3958: + snes - the SNES context
3959: - type - a known method
3961: Options Database Key:
3962: . -snes_type <type> - Sets the method; use -help for a list
3963: of available methods (for instance, newtonls or newtontr)
3965: Notes:
3966: See "petsc/include/petscsnes.h" for available methods (for instance)
3967: + SNESNEWTONLS - Newton's method with line search
3968: (systems of nonlinear equations)
3969: . SNESNEWTONTR - Newton's method with trust region
3970: (systems of nonlinear equations)
3972: Normally, it is best to use the SNESSetFromOptions() command and then
3973: set the SNES solver type from the options database rather than by using
3974: this routine. Using the options database provides the user with
3975: maximum flexibility in evaluating the many nonlinear solvers.
3976: The SNESSetType() routine is provided for those situations where it
3977: is necessary to set the nonlinear solver independently of the command
3978: line or options database. This might be the case, for example, when
3979: the choice of solver changes during the execution of the program,
3980: and the user's application is taking responsibility for choosing the
3981: appropriate method.
3983: Developer Notes: SNESRegister() adds a constructor for a new SNESType to SNESList, SNESSetType() locates
3984: the constructor in that list and calls it to create the spexific object.
3986: Level: intermediate
3988: .keywords: SNES, set, type
3990: .seealso: SNESType, SNESCreate(), SNESDestroy(), SNESGetType(), SNESSetFromOptions()
3992: @*/
3993: PetscErrorCode SNESSetType(SNES snes,SNESType type)
3994: {
3995: PetscErrorCode ierr,(*r)(SNES);
3996: PetscBool match;
4002: PetscObjectTypeCompare((PetscObject)snes,type,&match);
4003: if (match) return(0);
4005: PetscFunctionListFind(SNESList,type,&r);
4006: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
4007: /* Destroy the previous private SNES context */
4008: if (snes->ops->destroy) {
4009: (*(snes)->ops->destroy)(snes);
4010: snes->ops->destroy = NULL;
4011: }
4012: /* Reinitialize function pointers in SNESOps structure */
4013: snes->ops->setup = 0;
4014: snes->ops->solve = 0;
4015: snes->ops->view = 0;
4016: snes->ops->setfromoptions = 0;
4017: snes->ops->destroy = 0;
4018: /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
4019: snes->setupcalled = PETSC_FALSE;
4021: PetscObjectChangeTypeName((PetscObject)snes,type);
4022: (*r)(snes);
4023: return(0);
4024: }
4028: /*@C
4029: SNESGetType - Gets the SNES method type and name (as a string).
4031: Not Collective
4033: Input Parameter:
4034: . snes - nonlinear solver context
4036: Output Parameter:
4037: . type - SNES method (a character string)
4039: Level: intermediate
4041: .keywords: SNES, nonlinear, get, type, name
4042: @*/
4043: PetscErrorCode SNESGetType(SNES snes,SNESType *type)
4044: {
4048: *type = ((PetscObject)snes)->type_name;
4049: return(0);
4050: }
4054: /*@
4055: SNESSetSolution - Sets the solution vector for use by the SNES routines.
4057: Logically Collective on SNES and Vec
4059: Input Parameters:
4060: + snes - the SNES context obtained from SNESCreate()
4061: - u - the solution vector
4063: Level: beginner
4065: .keywords: SNES, set, solution
4066: @*/
4067: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4068: {
4069: DM dm;
4075: PetscObjectReference((PetscObject) u);
4076: VecDestroy(&snes->vec_sol);
4078: snes->vec_sol = u;
4080: SNESGetDM(snes, &dm);
4081: DMShellSetGlobalVector(dm, u);
4082: return(0);
4083: }
4087: /*@
4088: SNESGetSolution - Returns the vector where the approximate solution is
4089: stored. This is the fine grid solution when using SNESSetGridSequence().
4091: Not Collective, but Vec is parallel if SNES is parallel
4093: Input Parameter:
4094: . snes - the SNES context
4096: Output Parameter:
4097: . x - the solution
4099: Level: intermediate
4101: .keywords: SNES, nonlinear, get, solution
4103: .seealso: SNESGetSolutionUpdate(), SNESGetFunction()
4104: @*/
4105: PetscErrorCode SNESGetSolution(SNES snes,Vec *x)
4106: {
4110: *x = snes->vec_sol;
4111: return(0);
4112: }
4116: /*@
4117: SNESGetSolutionUpdate - Returns the vector where the solution update is
4118: stored.
4120: Not Collective, but Vec is parallel if SNES is parallel
4122: Input Parameter:
4123: . snes - the SNES context
4125: Output Parameter:
4126: . x - the solution update
4128: Level: advanced
4130: .keywords: SNES, nonlinear, get, solution, update
4132: .seealso: SNESGetSolution(), SNESGetFunction()
4133: @*/
4134: PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x)
4135: {
4139: *x = snes->vec_sol_update;
4140: return(0);
4141: }
4145: /*@C
4146: SNESGetFunction - Returns the vector where the function is stored.
4148: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
4150: Input Parameter:
4151: . snes - the SNES context
4153: Output Parameter:
4154: + r - the vector that is used to store residuals (or NULL if you don't want it)
4155: . f - the function (or NULL if you don't want it); see SNESFunction for calling sequence details
4156: - ctx - the function context (or NULL if you don't want it)
4158: Level: advanced
4160: .keywords: SNES, nonlinear, get, function
4162: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4163: @*/
4164: PetscErrorCode SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4165: {
4167: DM dm;
4171: if (r) {
4172: if (!snes->vec_func) {
4173: if (snes->vec_rhs) {
4174: VecDuplicate(snes->vec_rhs,&snes->vec_func);
4175: } else if (snes->vec_sol) {
4176: VecDuplicate(snes->vec_sol,&snes->vec_func);
4177: } else if (snes->dm) {
4178: DMCreateGlobalVector(snes->dm,&snes->vec_func);
4179: }
4180: }
4181: *r = snes->vec_func;
4182: }
4183: SNESGetDM(snes,&dm);
4184: DMSNESGetFunction(dm,f,ctx);
4185: return(0);
4186: }
4188: /*@C
4189: SNESGetNGS - Returns the NGS function and context.
4191: Input Parameter:
4192: . snes - the SNES context
4194: Output Parameter:
4195: + f - the function (or NULL) see SNESNGSFunction for details
4196: - ctx - the function context (or NULL)
4198: Level: advanced
4200: .keywords: SNES, nonlinear, get, function
4202: .seealso: SNESSetNGS(), SNESGetFunction()
4203: @*/
4207: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4208: {
4210: DM dm;
4214: SNESGetDM(snes,&dm);
4215: DMSNESGetNGS(dm,f,ctx);
4216: return(0);
4217: }
4221: /*@C
4222: SNESSetOptionsPrefix - Sets the prefix used for searching for all
4223: SNES options in the database.
4225: Logically Collective on SNES
4227: Input Parameter:
4228: + snes - the SNES context
4229: - prefix - the prefix to prepend to all option names
4231: Notes:
4232: A hyphen (-) must NOT be given at the beginning of the prefix name.
4233: The first character of all runtime options is AUTOMATICALLY the hyphen.
4235: Level: advanced
4237: .keywords: SNES, set, options, prefix, database
4239: .seealso: SNESSetFromOptions()
4240: @*/
4241: PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[])
4242: {
4247: PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4248: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4249: if (snes->linesearch) {
4250: SNESGetLineSearch(snes,&snes->linesearch);
4251: PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4252: }
4253: KSPSetOptionsPrefix(snes->ksp,prefix);
4254: return(0);
4255: }
4259: /*@C
4260: SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4261: SNES options in the database.
4263: Logically Collective on SNES
4265: Input Parameters:
4266: + snes - the SNES context
4267: - prefix - the prefix to prepend to all option names
4269: Notes:
4270: A hyphen (-) must NOT be given at the beginning of the prefix name.
4271: The first character of all runtime options is AUTOMATICALLY the hyphen.
4273: Level: advanced
4275: .keywords: SNES, append, options, prefix, database
4277: .seealso: SNESGetOptionsPrefix()
4278: @*/
4279: PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4280: {
4285: PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4286: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4287: if (snes->linesearch) {
4288: SNESGetLineSearch(snes,&snes->linesearch);
4289: PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4290: }
4291: KSPAppendOptionsPrefix(snes->ksp,prefix);
4292: return(0);
4293: }
4297: /*@C
4298: SNESGetOptionsPrefix - Sets the prefix used for searching for all
4299: SNES options in the database.
4301: Not Collective
4303: Input Parameter:
4304: . snes - the SNES context
4306: Output Parameter:
4307: . prefix - pointer to the prefix string used
4309: Notes: On the fortran side, the user should pass in a string 'prefix' of
4310: sufficient length to hold the prefix.
4312: Level: advanced
4314: .keywords: SNES, get, options, prefix, database
4316: .seealso: SNESAppendOptionsPrefix()
4317: @*/
4318: PetscErrorCode SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4319: {
4324: PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4325: return(0);
4326: }
4331: /*@C
4332: SNESRegister - Adds a method to the nonlinear solver package.
4334: Not collective
4336: Input Parameters:
4337: + name_solver - name of a new user-defined solver
4338: - routine_create - routine to create method context
4340: Notes:
4341: SNESRegister() may be called multiple times to add several user-defined solvers.
4343: Sample usage:
4344: .vb
4345: SNESRegister("my_solver",MySolverCreate);
4346: .ve
4348: Then, your solver can be chosen with the procedural interface via
4349: $ SNESSetType(snes,"my_solver")
4350: or at runtime via the option
4351: $ -snes_type my_solver
4353: Level: advanced
4355: Note: If your function is not being put into a shared library then use SNESRegister() instead
4357: .keywords: SNES, nonlinear, register
4359: .seealso: SNESRegisterAll(), SNESRegisterDestroy()
4361: Level: advanced
4362: @*/
4363: PetscErrorCode SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4364: {
4368: PetscFunctionListAdd(&SNESList,sname,function);
4369: return(0);
4370: }
4374: PetscErrorCode SNESTestLocalMin(SNES snes)
4375: {
4377: PetscInt N,i,j;
4378: Vec u,uh,fh;
4379: PetscScalar value;
4380: PetscReal norm;
4383: SNESGetSolution(snes,&u);
4384: VecDuplicate(u,&uh);
4385: VecDuplicate(u,&fh);
4387: /* currently only works for sequential */
4388: PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4389: VecGetSize(u,&N);
4390: for (i=0; i<N; i++) {
4391: VecCopy(u,uh);
4392: PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4393: for (j=-10; j<11; j++) {
4394: value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4395: VecSetValue(uh,i,value,ADD_VALUES);
4396: SNESComputeFunction(snes,uh,fh);
4397: VecNorm(fh,NORM_2,&norm);
4398: PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);
4399: value = -value;
4400: VecSetValue(uh,i,value,ADD_VALUES);
4401: }
4402: }
4403: VecDestroy(&uh);
4404: VecDestroy(&fh);
4405: return(0);
4406: }
4410: /*@
4411: SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4412: computing relative tolerance for linear solvers within an inexact
4413: Newton method.
4415: Logically Collective on SNES
4417: Input Parameters:
4418: + snes - SNES context
4419: - flag - PETSC_TRUE or PETSC_FALSE
4421: Options Database:
4422: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4423: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
4424: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4425: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4426: . -snes_ksp_ew_gamma <gamma> - Sets gamma
4427: . -snes_ksp_ew_alpha <alpha> - Sets alpha
4428: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4429: - -snes_ksp_ew_threshold <threshold> - Sets threshold
4431: Notes:
4432: Currently, the default is to use a constant relative tolerance for
4433: the inner linear solvers. Alternatively, one can use the
4434: Eisenstat-Walker method, where the relative convergence tolerance
4435: is reset at each Newton iteration according progress of the nonlinear
4436: solver.
4438: Level: advanced
4440: Reference:
4441: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4442: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
4444: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
4446: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4447: @*/
4448: PetscErrorCode SNESKSPSetUseEW(SNES snes,PetscBool flag)
4449: {
4453: snes->ksp_ewconv = flag;
4454: return(0);
4455: }
4459: /*@
4460: SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
4461: for computing relative tolerance for linear solvers within an
4462: inexact Newton method.
4464: Not Collective
4466: Input Parameter:
4467: . snes - SNES context
4469: Output Parameter:
4470: . flag - PETSC_TRUE or PETSC_FALSE
4472: Notes:
4473: Currently, the default is to use a constant relative tolerance for
4474: the inner linear solvers. Alternatively, one can use the
4475: Eisenstat-Walker method, where the relative convergence tolerance
4476: is reset at each Newton iteration according progress of the nonlinear
4477: solver.
4479: Level: advanced
4481: Reference:
4482: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4483: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
4485: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
4487: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4488: @*/
4489: PetscErrorCode SNESKSPGetUseEW(SNES snes, PetscBool *flag)
4490: {
4494: *flag = snes->ksp_ewconv;
4495: return(0);
4496: }
4500: /*@
4501: SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
4502: convergence criteria for the linear solvers within an inexact
4503: Newton method.
4505: Logically Collective on SNES
4507: Input Parameters:
4508: + snes - SNES context
4509: . version - version 1, 2 (default is 2) or 3
4510: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4511: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4512: . gamma - multiplicative factor for version 2 rtol computation
4513: (0 <= gamma2 <= 1)
4514: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
4515: . alpha2 - power for safeguard
4516: - threshold - threshold for imposing safeguard (0 < threshold < 1)
4518: Note:
4519: Version 3 was contributed by Luis Chacon, June 2006.
4521: Use PETSC_DEFAULT to retain the default for any of the parameters.
4523: Level: advanced
4525: Reference:
4526: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4527: inexact Newton method", Utah State University Math. Stat. Dept. Res.
4528: Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
4530: .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
4532: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
4533: @*/
4534: PetscErrorCode SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
4535: {
4536: SNESKSPEW *kctx;
4540: kctx = (SNESKSPEW*)snes->kspconvctx;
4541: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4550: if (version != PETSC_DEFAULT) kctx->version = version;
4551: if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0;
4552: if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max;
4553: if (gamma != PETSC_DEFAULT) kctx->gamma = gamma;
4554: if (alpha != PETSC_DEFAULT) kctx->alpha = alpha;
4555: if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2;
4556: if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
4558: 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);
4559: 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);
4560: 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);
4561: 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);
4562: 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);
4563: 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);
4564: return(0);
4565: }
4569: /*@
4570: SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
4571: convergence criteria for the linear solvers within an inexact
4572: Newton method.
4574: Not Collective
4576: Input Parameters:
4577: snes - SNES context
4579: Output Parameters:
4580: + version - version 1, 2 (default is 2) or 3
4581: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4582: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4583: . gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
4584: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
4585: . alpha2 - power for safeguard
4586: - threshold - threshold for imposing safeguard (0 < threshold < 1)
4588: Level: advanced
4590: .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
4592: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
4593: @*/
4594: PetscErrorCode SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
4595: {
4596: SNESKSPEW *kctx;
4600: kctx = (SNESKSPEW*)snes->kspconvctx;
4601: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4602: if (version) *version = kctx->version;
4603: if (rtol_0) *rtol_0 = kctx->rtol_0;
4604: if (rtol_max) *rtol_max = kctx->rtol_max;
4605: if (gamma) *gamma = kctx->gamma;
4606: if (alpha) *alpha = kctx->alpha;
4607: if (alpha2) *alpha2 = kctx->alpha2;
4608: if (threshold) *threshold = kctx->threshold;
4609: return(0);
4610: }
4614: PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4615: {
4617: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
4618: PetscReal rtol = PETSC_DEFAULT,stol;
4621: if (!snes->ksp_ewconv) return(0);
4622: if (!snes->iter) {
4623: rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
4624: VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
4625: }
4626: else {
4627: if (kctx->version == 1) {
4628: rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
4629: if (rtol < 0.0) rtol = -rtol;
4630: stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
4631: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4632: } else if (kctx->version == 2) {
4633: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4634: stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
4635: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4636: } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
4637: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4638: /* safeguard: avoid sharp decrease of rtol */
4639: stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
4640: stol = PetscMax(rtol,stol);
4641: rtol = PetscMin(kctx->rtol_0,stol);
4642: /* safeguard: avoid oversolving */
4643: stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
4644: stol = PetscMax(rtol,stol);
4645: rtol = PetscMin(kctx->rtol_0,stol);
4646: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
4647: }
4648: /* safeguard: avoid rtol greater than one */
4649: rtol = PetscMin(rtol,kctx->rtol_max);
4650: KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
4651: PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
4652: return(0);
4653: }
4657: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4658: {
4660: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
4661: PCSide pcside;
4662: Vec lres;
4665: if (!snes->ksp_ewconv) return(0);
4666: KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
4667: kctx->norm_last = snes->norm;
4668: if (kctx->version == 1) {
4669: PC pc;
4670: PetscBool isNone;
4672: KSPGetPC(ksp, &pc);
4673: PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
4674: KSPGetPCSide(ksp,&pcside);
4675: if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
4676: /* KSP residual is true linear residual */
4677: KSPGetResidualNorm(ksp,&kctx->lresid_last);
4678: } else {
4679: /* KSP residual is preconditioned residual */
4680: /* compute true linear residual norm */
4681: VecDuplicate(b,&lres);
4682: MatMult(snes->jacobian,x,lres);
4683: VecAYPX(lres,-1.0,b);
4684: VecNorm(lres,NORM_2,&kctx->lresid_last);
4685: VecDestroy(&lres);
4686: }
4687: }
4688: return(0);
4689: }
4693: /*@
4694: SNESGetKSP - Returns the KSP context for a SNES solver.
4696: Not Collective, but if SNES object is parallel, then KSP object is parallel
4698: Input Parameter:
4699: . snes - the SNES context
4701: Output Parameter:
4702: . ksp - the KSP context
4704: Notes:
4705: The user can then directly manipulate the KSP context to set various
4706: options, etc. Likewise, the user can then extract and manipulate the
4707: PC contexts as well.
4709: Level: beginner
4711: .keywords: SNES, nonlinear, get, KSP, context
4713: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
4714: @*/
4715: PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp)
4716: {
4723: if (!snes->ksp) {
4724: PetscBool monitor = PETSC_FALSE;
4726: KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
4727: PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
4728: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);
4730: KSPSetPreSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPreSolve_SNESEW,snes);
4731: KSPSetPostSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPostSolve_SNESEW,snes);
4733: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
4734: if (monitor) {
4735: KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
4736: }
4737: monitor = PETSC_FALSE;
4738: PetscOptionsGetBool(((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
4739: if (monitor) {
4740: PetscObject *objs;
4741: KSPMonitorSNESLGResidualNormCreate(0,0,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
4742: objs[0] = (PetscObject) snes;
4743: KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
4744: }
4745: }
4746: *ksp = snes->ksp;
4747: return(0);
4748: }
4751: #include <petsc/private/dmimpl.h>
4754: /*@
4755: SNESSetDM - Sets the DM that may be used by some preconditioners
4757: Logically Collective on SNES
4759: Input Parameters:
4760: + snes - the preconditioner context
4761: - dm - the dm
4763: Level: intermediate
4765: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
4766: @*/
4767: PetscErrorCode SNESSetDM(SNES snes,DM dm)
4768: {
4770: KSP ksp;
4771: DMSNES sdm;
4775: if (dm) {PetscObjectReference((PetscObject)dm);}
4776: if (snes->dm) { /* Move the DMSNES context over to the new DM unless the new DM already has one */
4777: if (snes->dm->dmsnes && snes->dmAuto && !dm->dmsnes) {
4778: DMCopyDMSNES(snes->dm,dm);
4779: DMGetDMSNES(snes->dm,&sdm);
4780: if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
4781: }
4782: DMDestroy(&snes->dm);
4783: }
4784: snes->dm = dm;
4785: snes->dmAuto = PETSC_FALSE;
4787: SNESGetKSP(snes,&ksp);
4788: KSPSetDM(ksp,dm);
4789: KSPSetDMActive(ksp,PETSC_FALSE);
4790: if (snes->pc) {
4791: SNESSetDM(snes->pc, snes->dm);
4792: SNESSetNPCSide(snes,snes->pcside);
4793: }
4794: return(0);
4795: }
4799: /*@
4800: SNESGetDM - Gets the DM that may be used by some preconditioners
4802: Not Collective but DM obtained is parallel on SNES
4804: Input Parameter:
4805: . snes - the preconditioner context
4807: Output Parameter:
4808: . dm - the dm
4810: Level: intermediate
4812: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
4813: @*/
4814: PetscErrorCode SNESGetDM(SNES snes,DM *dm)
4815: {
4820: if (!snes->dm) {
4821: DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
4822: snes->dmAuto = PETSC_TRUE;
4823: }
4824: *dm = snes->dm;
4825: return(0);
4826: }
4830: /*@
4831: SNESSetNPC - Sets the nonlinear preconditioner to be used.
4833: Collective on SNES
4835: Input Parameters:
4836: + snes - iterative context obtained from SNESCreate()
4837: - pc - the preconditioner object
4839: Notes:
4840: Use SNESGetNPC() to retrieve the preconditioner context (for example,
4841: to configure it using the API).
4843: Level: developer
4845: .keywords: SNES, set, precondition
4846: .seealso: SNESGetNPC(), SNESHasNPC()
4847: @*/
4848: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
4849: {
4856: PetscObjectReference((PetscObject) pc);
4857: SNESDestroy(&snes->pc);
4858: snes->pc = pc;
4859: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->pc);
4860: return(0);
4861: }
4865: /*@
4866: SNESGetNPC - Creates a nonlinear preconditioning solver (SNES) to be used to precondition the nonlinear solver.
4868: Not Collective
4870: Input Parameter:
4871: . snes - iterative context obtained from SNESCreate()
4873: Output Parameter:
4874: . pc - preconditioner context
4876: Notes: If a SNES was previously set with SNESSetNPC() then that SNES is returned.
4878: Level: developer
4880: .keywords: SNES, get, preconditioner
4881: .seealso: SNESSetNPC(), SNESHasNPC()
4882: @*/
4883: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
4884: {
4886: const char *optionsprefix;
4891: if (!snes->pc) {
4892: SNESCreate(PetscObjectComm((PetscObject)snes),&snes->pc);
4893: PetscObjectIncrementTabLevel((PetscObject)snes->pc,(PetscObject)snes,1);
4894: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->pc);
4895: SNESGetOptionsPrefix(snes,&optionsprefix);
4896: SNESSetOptionsPrefix(snes->pc,optionsprefix);
4897: SNESAppendOptionsPrefix(snes->pc,"npc_");
4898: SNESSetCountersReset(snes->pc,PETSC_FALSE);
4899: }
4900: *pc = snes->pc;
4901: return(0);
4902: }
4906: /*@
4907: SNESHasNPC - Returns whether a nonlinear preconditioner exists
4909: Not Collective
4911: Input Parameter:
4912: . snes - iterative context obtained from SNESCreate()
4914: Output Parameter:
4915: . has_npc - whether the SNES has an NPC or not
4917: Level: developer
4919: .keywords: SNES, has, preconditioner
4920: .seealso: SNESSetNPC(), SNESGetNPC()
4921: @*/
4922: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
4923: {
4926: *has_npc = (PetscBool) (snes->pc != NULL);
4927: return(0);
4928: }
4932: /*@
4933: SNESSetNPCSide - Sets the preconditioning side.
4935: Logically Collective on SNES
4937: Input Parameter:
4938: . snes - iterative context obtained from SNESCreate()
4940: Output Parameter:
4941: . side - the preconditioning side, where side is one of
4942: .vb
4943: PC_LEFT - left preconditioning (default)
4944: PC_RIGHT - right preconditioning
4945: .ve
4947: Options Database Keys:
4948: . -snes_pc_side <right,left>
4950: Level: intermediate
4952: .keywords: SNES, set, right, left, side, preconditioner, flag
4954: .seealso: SNESGetNPCSide(), KSPSetPCSide()
4955: @*/
4956: PetscErrorCode SNESSetNPCSide(SNES snes,PCSide side)
4957: {
4961: snes->pcside = side;
4962: return(0);
4963: }
4967: /*@
4968: SNESGetNPCSide - Gets the preconditioning side.
4970: Not Collective
4972: Input Parameter:
4973: . snes - iterative context obtained from SNESCreate()
4975: Output Parameter:
4976: . side - the preconditioning side, where side is one of
4977: .vb
4978: PC_LEFT - left preconditioning (default)
4979: PC_RIGHT - right preconditioning
4980: .ve
4982: Level: intermediate
4984: .keywords: SNES, get, right, left, side, preconditioner, flag
4986: .seealso: SNESSetNPCSide(), KSPGetPCSide()
4987: @*/
4988: PetscErrorCode SNESGetNPCSide(SNES snes,PCSide *side)
4989: {
4993: *side = snes->pcside;
4994: return(0);
4995: }
4999: /*@
5000: SNESSetLineSearch - Sets the linesearch on the SNES instance.
5002: Collective on SNES
5004: Input Parameters:
5005: + snes - iterative context obtained from SNESCreate()
5006: - linesearch - the linesearch object
5008: Notes:
5009: Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
5010: to configure it using the API).
5012: Level: developer
5014: .keywords: SNES, set, linesearch
5015: .seealso: SNESGetLineSearch()
5016: @*/
5017: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
5018: {
5025: PetscObjectReference((PetscObject) linesearch);
5026: SNESLineSearchDestroy(&snes->linesearch);
5028: snes->linesearch = linesearch;
5030: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5031: return(0);
5032: }
5036: /*@
5037: SNESGetLineSearch - Returns a pointer to the line search context set with SNESSetLineSearch()
5038: or creates a default line search instance associated with the SNES and returns it.
5040: Not Collective
5042: Input Parameter:
5043: . snes - iterative context obtained from SNESCreate()
5045: Output Parameter:
5046: . linesearch - linesearch context
5048: Level: beginner
5050: .keywords: SNES, get, linesearch
5051: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5052: @*/
5053: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5054: {
5056: const char *optionsprefix;
5061: if (!snes->linesearch) {
5062: SNESGetOptionsPrefix(snes, &optionsprefix);
5063: SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5064: SNESLineSearchSetSNES(snes->linesearch, snes);
5065: SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5066: PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5067: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5068: }
5069: *linesearch = snes->linesearch;
5070: return(0);
5071: }
5073: #if defined(PETSC_HAVE_MATLAB_ENGINE)
5074: #include <mex.h>
5076: typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
5080: /*
5081: SNESComputeFunction_Matlab - Calls the function that has been set with SNESSetFunctionMatlab().
5083: Collective on SNES
5085: Input Parameters:
5086: + snes - the SNES context
5087: - x - input vector
5089: Output Parameter:
5090: . y - function vector, as set by SNESSetFunction()
5092: Notes:
5093: SNESComputeFunction() is typically used within nonlinear solvers
5094: implementations, so most users would not generally call this routine
5095: themselves.
5097: Level: developer
5099: .keywords: SNES, nonlinear, compute, function
5101: .seealso: SNESSetFunction(), SNESGetFunction()
5102: */
5103: PetscErrorCode SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
5104: {
5105: PetscErrorCode ierr;
5106: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5107: int nlhs = 1,nrhs = 5;
5108: mxArray *plhs[1],*prhs[5];
5109: long long int lx = 0,ly = 0,ls = 0;
5118: /* call Matlab function in ctx with arguments x and y */
5120: PetscMemcpy(&ls,&snes,sizeof(snes));
5121: PetscMemcpy(&lx,&x,sizeof(x));
5122: PetscMemcpy(&ly,&y,sizeof(x));
5123: prhs[0] = mxCreateDoubleScalar((double)ls);
5124: prhs[1] = mxCreateDoubleScalar((double)lx);
5125: prhs[2] = mxCreateDoubleScalar((double)ly);
5126: prhs[3] = mxCreateString(sctx->funcname);
5127: prhs[4] = sctx->ctx;
5128: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");
5129: mxGetScalar(plhs[0]);
5130: mxDestroyArray(prhs[0]);
5131: mxDestroyArray(prhs[1]);
5132: mxDestroyArray(prhs[2]);
5133: mxDestroyArray(prhs[3]);
5134: mxDestroyArray(plhs[0]);
5135: return(0);
5136: }
5140: /*
5141: SNESSetFunctionMatlab - Sets the function evaluation routine and function
5142: vector for use by the SNES routines in solving systems of nonlinear
5143: equations from MATLAB. Here the function is a string containing the name of a MATLAB function
5145: Logically Collective on SNES
5147: Input Parameters:
5148: + snes - the SNES context
5149: . r - vector to store function value
5150: - f - function evaluation routine
5152: Notes:
5153: The Newton-like methods typically solve linear systems of the form
5154: $ f'(x) x = -f(x),
5155: where f'(x) denotes the Jacobian matrix and f(x) is the function.
5157: Level: beginner
5159: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5161: .keywords: SNES, nonlinear, set, function
5163: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5164: */
5165: PetscErrorCode SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx)
5166: {
5167: PetscErrorCode ierr;
5168: SNESMatlabContext *sctx;
5171: /* currently sctx is memory bleed */
5172: PetscNew(&sctx);
5173: PetscStrallocpy(f,&sctx->funcname);
5174: /*
5175: This should work, but it doesn't
5176: sctx->ctx = ctx;
5177: mexMakeArrayPersistent(sctx->ctx);
5178: */
5179: sctx->ctx = mxDuplicateArray(ctx);
5180: SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);
5181: return(0);
5182: }
5186: /*
5187: SNESComputeJacobian_Matlab - Calls the function that has been set with SNESSetJacobianMatlab().
5189: Collective on SNES
5191: Input Parameters:
5192: + snes - the SNES context
5193: . x - input vector
5194: . A, B - the matrices
5195: - ctx - user context
5197: Level: developer
5199: .keywords: SNES, nonlinear, compute, function
5201: .seealso: SNESSetFunction(), SNESGetFunction()
5202: @*/
5203: PetscErrorCode SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx)
5204: {
5205: PetscErrorCode ierr;
5206: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5207: int nlhs = 2,nrhs = 6;
5208: mxArray *plhs[2],*prhs[6];
5209: long long int lx = 0,lA = 0,ls = 0, lB = 0;
5215: /* call Matlab function in ctx with arguments x and y */
5217: PetscMemcpy(&ls,&snes,sizeof(snes));
5218: PetscMemcpy(&lx,&x,sizeof(x));
5219: PetscMemcpy(&lA,A,sizeof(x));
5220: PetscMemcpy(&lB,B,sizeof(x));
5221: prhs[0] = mxCreateDoubleScalar((double)ls);
5222: prhs[1] = mxCreateDoubleScalar((double)lx);
5223: prhs[2] = mxCreateDoubleScalar((double)lA);
5224: prhs[3] = mxCreateDoubleScalar((double)lB);
5225: prhs[4] = mxCreateString(sctx->funcname);
5226: prhs[5] = sctx->ctx;
5227: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");
5228: mxGetScalar(plhs[0]);
5229: mxDestroyArray(prhs[0]);
5230: mxDestroyArray(prhs[1]);
5231: mxDestroyArray(prhs[2]);
5232: mxDestroyArray(prhs[3]);
5233: mxDestroyArray(prhs[4]);
5234: mxDestroyArray(plhs[0]);
5235: mxDestroyArray(plhs[1]);
5236: return(0);
5237: }
5241: /*
5242: SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
5243: vector for use by the SNES routines in solving systems of nonlinear
5244: equations from MATLAB. Here the function is a string containing the name of a MATLAB function
5246: Logically Collective on SNES
5248: Input Parameters:
5249: + snes - the SNES context
5250: . A,B - Jacobian matrices
5251: . J - function evaluation routine
5252: - ctx - user context
5254: Level: developer
5256: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5258: .keywords: SNES, nonlinear, set, function
5260: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J
5261: */
5262: PetscErrorCode SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx)
5263: {
5264: PetscErrorCode ierr;
5265: SNESMatlabContext *sctx;
5268: /* currently sctx is memory bleed */
5269: PetscNew(&sctx);
5270: PetscStrallocpy(J,&sctx->funcname);
5271: /*
5272: This should work, but it doesn't
5273: sctx->ctx = ctx;
5274: mexMakeArrayPersistent(sctx->ctx);
5275: */
5276: sctx->ctx = mxDuplicateArray(ctx);
5277: SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);
5278: return(0);
5279: }
5283: /*
5284: SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().
5286: Collective on SNES
5288: .seealso: SNESSetFunction(), SNESGetFunction()
5289: @*/
5290: PetscErrorCode SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
5291: {
5292: PetscErrorCode ierr;
5293: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5294: int nlhs = 1,nrhs = 6;
5295: mxArray *plhs[1],*prhs[6];
5296: long long int lx = 0,ls = 0;
5297: Vec x = snes->vec_sol;
5302: PetscMemcpy(&ls,&snes,sizeof(snes));
5303: PetscMemcpy(&lx,&x,sizeof(x));
5304: prhs[0] = mxCreateDoubleScalar((double)ls);
5305: prhs[1] = mxCreateDoubleScalar((double)it);
5306: prhs[2] = mxCreateDoubleScalar((double)fnorm);
5307: prhs[3] = mxCreateDoubleScalar((double)lx);
5308: prhs[4] = mxCreateString(sctx->funcname);
5309: prhs[5] = sctx->ctx;
5310: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");
5311: mxGetScalar(plhs[0]);
5312: mxDestroyArray(prhs[0]);
5313: mxDestroyArray(prhs[1]);
5314: mxDestroyArray(prhs[2]);
5315: mxDestroyArray(prhs[3]);
5316: mxDestroyArray(prhs[4]);
5317: mxDestroyArray(plhs[0]);
5318: return(0);
5319: }
5323: /*
5324: SNESMonitorSetMatlab - Sets the monitor function from MATLAB
5326: Level: developer
5328: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5330: .keywords: SNES, nonlinear, set, function
5332: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5333: */
5334: PetscErrorCode SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx)
5335: {
5336: PetscErrorCode ierr;
5337: SNESMatlabContext *sctx;
5340: /* currently sctx is memory bleed */
5341: PetscNew(&sctx);
5342: PetscStrallocpy(f,&sctx->funcname);
5343: /*
5344: This should work, but it doesn't
5345: sctx->ctx = ctx;
5346: mexMakeArrayPersistent(sctx->ctx);
5347: */
5348: sctx->ctx = mxDuplicateArray(ctx);
5349: SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);
5350: return(0);
5351: }
5353: #endif