Actual source code: snes.c
petsc-3.9.4 2018-09-11
2: #include <petsc/private/snesimpl.h>
3: #include <petscdmshell.h>
4: #include <petscdraw.h>
5: #include <petscds.h>
6: #include <petscdmadaptor.h>
7: #include <petscconvest.h>
9: PetscBool SNESRegisterAllCalled = PETSC_FALSE;
10: PetscFunctionList SNESList = NULL;
12: /* Logging support */
13: PetscClassId SNES_CLASSID, DMSNES_CLASSID;
14: PetscLogEvent SNES_Solve, SNES_FunctionEval, SNES_JacobianEval, SNES_NGSEval, SNES_NGSFuncEval, SNES_NPCSolve, SNES_ObjectiveEval;
16: /*@
17: SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged.
19: Logically Collective on SNES
21: Input Parameters:
22: + snes - iterative context obtained from SNESCreate()
23: - flg - PETSC_TRUE indicates you want the error generated
25: Options database keys:
26: . -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)
28: Level: intermediate
30: Notes:
31: Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve()
32: to determine if it has converged.
34: .keywords: SNES, set, initial guess, nonzero
36: .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
37: @*/
38: PetscErrorCode SNESSetErrorIfNotConverged(SNES snes,PetscBool flg)
39: {
43: snes->errorifnotconverged = flg;
44: return(0);
45: }
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: }
73: /*@
74: SNESSetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?
76: Logically Collective on SNES
78: Input Parameters:
79: + snes - the shell SNES
80: - flg - is the residual computed?
82: Level: advanced
84: .seealso: SNESGetAlwaysComputesFinalResidual()
85: @*/
86: PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES snes, PetscBool flg)
87: {
90: snes->alwayscomputesfinalresidual = flg;
91: return(0);
92: }
94: /*@
95: SNESGetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?
97: Logically Collective on SNES
99: Input Parameter:
100: . snes - the shell SNES
102: Output Parameter:
103: . flg - is the residual computed?
105: Level: advanced
107: .seealso: SNESSetAlwaysComputesFinalResidual()
108: @*/
109: PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES snes, PetscBool *flg)
110: {
113: *flg = snes->alwayscomputesfinalresidual;
114: return(0);
115: }
117: /*@
118: SNESSetFunctionDomainError - tells SNES that the input vector to your SNESFunction is not
119: in the functions domain. For example, negative pressure.
121: Logically Collective on SNES
123: Input Parameters:
124: . snes - the SNES context
126: Level: advanced
128: .keywords: SNES, view
130: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction
131: @*/
132: PetscErrorCode SNESSetFunctionDomainError(SNES snes)
133: {
136: if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates input vector is not in the function domain");
137: snes->domainerror = PETSC_TRUE;
138: return(0);
139: }
141: /*@
142: SNESGetFunctionDomainError - Gets the status of the domain error after a call to SNESComputeFunction;
144: Logically Collective on SNES
146: Input Parameters:
147: . snes - the SNES context
149: Output Parameters:
150: . domainerror - Set to PETSC_TRUE if there's a domain error; PETSC_FALSE otherwise.
152: Level: advanced
154: .keywords: SNES, view
156: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction()
157: @*/
158: PetscErrorCode SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror)
159: {
163: *domainerror = snes->domainerror;
164: return(0);
165: }
167: /*@C
168: SNESLoad - Loads a SNES that has been stored in binary with SNESView().
170: Collective on PetscViewer
172: Input Parameters:
173: + newdm - the newly loaded SNES, this needs to have been created with SNESCreate() or
174: some related function before a call to SNESLoad().
175: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
177: Level: intermediate
179: Notes:
180: The type is determined by the data in the file, any type set into the SNES before this call is ignored.
182: Notes for advanced users:
183: Most users should not need to know the details of the binary storage
184: format, since SNESLoad() and TSView() completely hide these details.
185: But for anyone who's interested, the standard binary matrix storage
186: format is
187: .vb
188: has not yet been determined
189: .ve
191: .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad()
192: @*/
193: PetscErrorCode SNESLoad(SNES snes, PetscViewer viewer)
194: {
196: PetscBool isbinary;
197: PetscInt classid;
198: char type[256];
199: KSP ksp;
200: DM dm;
201: DMSNES dmsnes;
206: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
207: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
209: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
210: if (classid != SNES_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONG,"Not SNES next in file");
211: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
212: SNESSetType(snes, type);
213: if (snes->ops->load) {
214: (*snes->ops->load)(snes,viewer);
215: }
216: SNESGetDM(snes,&dm);
217: DMGetDMSNES(dm,&dmsnes);
218: DMSNESLoad(dmsnes,viewer);
219: SNESGetKSP(snes,&ksp);
220: KSPLoad(ksp,viewer);
221: return(0);
222: }
224: #include <petscdraw.h>
225: #if defined(PETSC_HAVE_SAWS)
226: #include <petscviewersaws.h>
227: #endif
229: /*@C
230: SNESView - Prints the SNES data structure.
232: Collective on SNES
234: Input Parameters:
235: + SNES - the SNES context
236: - viewer - visualization context
238: Options Database Key:
239: . -snes_view - Calls SNESView() at end of SNESSolve()
241: Notes:
242: The available visualization contexts include
243: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
244: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
245: output where only the first processor opens
246: the file. All other processors send their
247: data to the first processor to print.
249: The user can open an alternative visualization context with
250: PetscViewerASCIIOpen() - output to a specified file.
252: Level: beginner
254: .keywords: SNES, view
256: .seealso: PetscViewerASCIIOpen()
257: @*/
258: PetscErrorCode SNESView(SNES snes,PetscViewer viewer)
259: {
260: SNESKSPEW *kctx;
262: KSP ksp;
263: SNESLineSearch linesearch;
264: PetscBool iascii,isstring,isbinary,isdraw;
265: DMSNES dmsnes;
266: #if defined(PETSC_HAVE_SAWS)
267: PetscBool issaws;
268: #endif
272: if (!viewer) {
273: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
274: }
278: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
279: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
280: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
281: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
282: #if defined(PETSC_HAVE_SAWS)
283: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
284: #endif
285: if (iascii) {
286: SNESNormSchedule normschedule;
287: DM dm;
288: PetscErrorCode (*cJ)(SNES,Vec,Mat,Mat,void*);
289: void *ctx;
291: PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);
292: if (!snes->setupcalled) {
293: PetscViewerASCIIPrintf(viewer," SNES has not been set up so information may be incomplete\n");
294: }
295: if (snes->ops->view) {
296: PetscViewerASCIIPushTab(viewer);
297: (*snes->ops->view)(snes,viewer);
298: PetscViewerASCIIPopTab(viewer);
299: }
300: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);
301: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);
302: if (snes->usesksp) {
303: PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);
304: }
305: PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);
306: SNESGetNormSchedule(snes, &normschedule);
307: if (normschedule > 0) {PetscViewerASCIIPrintf(viewer," norm schedule %s\n",SNESNormSchedules[normschedule]);}
308: if (snes->gridsequence) {
309: PetscViewerASCIIPrintf(viewer," total number of grid sequence refinements=%D\n",snes->gridsequence);
310: }
311: if (snes->ksp_ewconv) {
312: kctx = (SNESKSPEW*)snes->kspconvctx;
313: if (kctx) {
314: PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);
315: PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);
316: PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);
317: }
318: }
319: if (snes->lagpreconditioner == -1) {
320: PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");
321: } else if (snes->lagpreconditioner > 1) {
322: PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);
323: }
324: if (snes->lagjacobian == -1) {
325: PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");
326: } else if (snes->lagjacobian > 1) {
327: PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);
328: }
329: SNESGetDM(snes,&dm);
330: DMSNESGetJacobian(dm,&cJ,&ctx);
331: if (cJ == SNESComputeJacobianDefault) {
332: PetscViewerASCIIPrintf(viewer," Jacobian is built using finite differences one column at a time\n");
333: } else if (cJ == SNESComputeJacobianDefaultColor) {
334: PetscViewerASCIIPrintf(viewer," Jacobian is built using finite differences with coloring\n");
335: }
336: } else if (isstring) {
337: const char *type;
338: SNESGetType(snes,&type);
339: PetscViewerStringSPrintf(viewer," %-3.3s",type);
340: } else if (isbinary) {
341: PetscInt classid = SNES_FILE_CLASSID;
342: MPI_Comm comm;
343: PetscMPIInt rank;
344: char type[256];
346: PetscObjectGetComm((PetscObject)snes,&comm);
347: MPI_Comm_rank(comm,&rank);
348: if (!rank) {
349: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
350: PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));
351: PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR,PETSC_FALSE);
352: }
353: if (snes->ops->view) {
354: (*snes->ops->view)(snes,viewer);
355: }
356: } else if (isdraw) {
357: PetscDraw draw;
358: char str[36];
359: PetscReal x,y,bottom,h;
361: PetscViewerDrawGetDraw(viewer,0,&draw);
362: PetscDrawGetCurrentPoint(draw,&x,&y);
363: PetscStrncpy(str,"SNES: ",sizeof(str));
364: PetscStrlcat(str,((PetscObject)snes)->type_name,sizeof(str));
365: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);
366: bottom = y - h;
367: PetscDrawPushCurrentPoint(draw,x,bottom);
368: if (snes->ops->view) {
369: (*snes->ops->view)(snes,viewer);
370: }
371: #if defined(PETSC_HAVE_SAWS)
372: } else if (issaws) {
373: PetscMPIInt rank;
374: const char *name;
376: PetscObjectGetName((PetscObject)snes,&name);
377: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
378: if (!((PetscObject)snes)->amsmem && !rank) {
379: char dir[1024];
381: PetscObjectViewSAWs((PetscObject)snes,viewer);
382: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
383: PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT));
384: if (!snes->conv_hist) {
385: SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);
386: }
387: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);
388: PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE));
389: }
390: #endif
391: }
392: if (snes->linesearch) {
393: SNESGetLineSearch(snes, &linesearch);
394: PetscViewerASCIIPushTab(viewer);
395: SNESLineSearchView(linesearch, viewer);
396: PetscViewerASCIIPopTab(viewer);
397: }
398: if (snes->npc && snes->usesnpc) {
399: PetscViewerASCIIPushTab(viewer);
400: SNESView(snes->npc, viewer);
401: PetscViewerASCIIPopTab(viewer);
402: }
403: PetscViewerASCIIPushTab(viewer);
404: DMGetDMSNES(snes->dm,&dmsnes);
405: DMSNESView(dmsnes, viewer);
406: PetscViewerASCIIPopTab(viewer);
407: if (snes->usesksp) {
408: SNESGetKSP(snes,&ksp);
409: PetscViewerASCIIPushTab(viewer);
410: KSPView(ksp,viewer);
411: PetscViewerASCIIPopTab(viewer);
412: }
413: if (isdraw) {
414: PetscDraw draw;
415: PetscViewerDrawGetDraw(viewer,0,&draw);
416: PetscDrawPopCurrentPoint(draw);
417: }
418: return(0);
419: }
421: /*
422: We retain a list of functions that also take SNES command
423: line options. These are called at the end SNESSetFromOptions()
424: */
425: #define MAXSETFROMOPTIONS 5
426: static PetscInt numberofsetfromoptions;
427: static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
429: /*@C
430: SNESAddOptionsChecker - Adds an additional function to check for SNES options.
432: Not Collective
434: Input Parameter:
435: . snescheck - function that checks for options
437: Level: developer
439: .seealso: SNESSetFromOptions()
440: @*/
441: PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
442: {
444: if (numberofsetfromoptions >= MAXSETFROMOPTIONS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
445: othersetfromoptions[numberofsetfromoptions++] = snescheck;
446: return(0);
447: }
449: extern PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
451: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
452: {
453: Mat J;
454: KSP ksp;
455: PC pc;
456: PetscBool match;
458: MatNullSpace nullsp;
463: if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
464: Mat A = snes->jacobian, B = snes->jacobian_pre;
465: MatCreateVecs(A ? A : B, NULL,&snes->vec_func);
466: }
468: if (version == 1) {
469: MatCreateSNESMF(snes,&J);
470: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
471: MatSetFromOptions(J);
472: } else if (version == 2) {
473: if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
474: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
475: SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);
476: #else
477: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)");
478: #endif
479: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2");
481: /* attach any user provided null space that was on Amat to the newly created matrix free matrix */
482: if (snes->jacobian) {
483: MatGetNullSpace(snes->jacobian,&nullsp);
484: if (nullsp) {
485: MatSetNullSpace(J,nullsp);
486: }
487: }
489: PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);
490: if (hasOperator) {
492: /* This version replaces the user provided Jacobian matrix with a
493: matrix-free version but still employs the user-provided preconditioner matrix. */
494: SNESSetJacobian(snes,J,0,0,0);
495: } else {
496: /* This version replaces both the user-provided Jacobian and the user-
497: provided preconditioner Jacobian with the default matrix free version. */
498: if ((snes->npcside== PC_LEFT) && snes->npc) {
499: if (!snes->jacobian){SNESSetJacobian(snes,J,0,0,0);}
500: } else {
501: SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,0);
502: }
503: /* Force no preconditioner */
504: SNESGetKSP(snes,&ksp);
505: KSPGetPC(ksp,&pc);
506: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);
507: if (!match) {
508: PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");
509: PCSetType(pc,PCNONE);
510: }
511: }
512: MatDestroy(&J);
513: return(0);
514: }
516: static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx)
517: {
518: SNES snes = (SNES)ctx;
520: Vec Xfine,Xfine_named = NULL,Xcoarse;
523: if (PetscLogPrintInfo) {
524: PetscInt finelevel,coarselevel,fineclevel,coarseclevel;
525: DMGetRefineLevel(dmfine,&finelevel);
526: DMGetCoarsenLevel(dmfine,&fineclevel);
527: DMGetRefineLevel(dmcoarse,&coarselevel);
528: DMGetCoarsenLevel(dmcoarse,&coarseclevel);
529: PetscInfo4(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);
530: }
531: if (dmfine == snes->dm) Xfine = snes->vec_sol;
532: else {
533: DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
534: Xfine = Xfine_named;
535: }
536: DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
537: if (Inject) {
538: MatRestrict(Inject,Xfine,Xcoarse);
539: } else {
540: MatRestrict(Restrict,Xfine,Xcoarse);
541: VecPointwiseMult(Xcoarse,Xcoarse,Rscale);
542: }
543: DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
544: if (Xfine_named) {DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);}
545: return(0);
546: }
548: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
549: {
553: DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
554: return(0);
555: }
557: /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can
558: * safely call SNESGetDM() in their residual evaluation routine. */
559: static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx)
560: {
561: SNES snes = (SNES)ctx;
563: Mat Asave = A,Bsave = B;
564: Vec X,Xnamed = NULL;
565: DM dmsave;
566: void *ctxsave;
567: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL;
570: dmsave = snes->dm;
571: KSPGetDM(ksp,&snes->dm);
572: if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */
573: else { /* We are on a coarser level, this vec was initialized using a DM restrict hook */
574: DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
575: X = Xnamed;
576: SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);
577: /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */
578: if (jac == SNESComputeJacobianDefaultColor) {
579: SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,0);
580: }
581: }
582: /* put the previous context back */
584: SNESComputeJacobian(snes,X,A,B);
585: if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
586: SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
587: }
589: if (A != Asave || B != Bsave) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for changing matrices at this time");
590: if (Xnamed) {
591: DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
592: }
593: snes->dm = dmsave;
594: return(0);
595: }
597: /*@
598: SNESSetUpMatrices - ensures that matrices are available for SNES, to be called by SNESSetUp_XXX()
600: Collective
602: Input Arguments:
603: . snes - snes to configure
605: Level: developer
607: .seealso: SNESSetUp()
608: @*/
609: PetscErrorCode SNESSetUpMatrices(SNES snes)
610: {
612: DM dm;
613: DMSNES sdm;
616: SNESGetDM(snes,&dm);
617: DMGetDMSNES(dm,&sdm);
618: if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"DMSNES not properly configured");
619: else if (!snes->jacobian && snes->mf) {
620: Mat J;
621: void *functx;
622: MatCreateSNESMF(snes,&J);
623: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
624: MatSetFromOptions(J);
625: SNESGetFunction(snes,NULL,NULL,&functx);
626: SNESSetJacobian(snes,J,J,0,0);
627: MatDestroy(&J);
628: } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
629: Mat J,B;
630: MatCreateSNESMF(snes,&J);
631: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
632: MatSetFromOptions(J);
633: DMCreateMatrix(snes->dm,&B);
634: /* sdm->computejacobian was already set to reach here */
635: SNESSetJacobian(snes,J,B,NULL,NULL);
636: MatDestroy(&J);
637: MatDestroy(&B);
638: } else if (!snes->jacobian_pre) {
639: PetscDS prob;
640: Mat J, B;
641: PetscBool hasPrec = PETSC_FALSE;
643: J = snes->jacobian;
644: DMGetDS(dm, &prob);
645: if (prob) {PetscDSHasJacobianPreconditioner(prob, &hasPrec);}
646: if (J) {PetscObjectReference((PetscObject) J);}
647: else if (hasPrec) {DMCreateMatrix(snes->dm, &J);}
648: DMCreateMatrix(snes->dm, &B);
649: SNESSetJacobian(snes, J ? J : B, B, NULL, NULL);
650: MatDestroy(&J);
651: MatDestroy(&B);
652: }
653: {
654: KSP ksp;
655: SNESGetKSP(snes,&ksp);
656: KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);
657: DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
658: }
659: return(0);
660: }
662: /*@C
663: SNESMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
665: Collective on SNES
667: Input Parameters:
668: + snes - SNES object you wish to monitor
669: . name - the monitor type one is seeking
670: . help - message indicating what monitoring is done
671: . manual - manual page for the monitor
672: . monitor - the monitor function
673: - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the SNES or PetscViewer objects
675: Level: developer
677: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
678: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
679: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
680: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
681: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
682: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
683: PetscOptionsFList(), PetscOptionsEList()
684: @*/
685: PetscErrorCode SNESMonitorSetFromOptions(SNES snes,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNES,PetscViewerAndFormat*))
686: {
687: PetscErrorCode ierr;
688: PetscViewer viewer;
689: PetscViewerFormat format;
690: PetscBool flg;
693: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,name,&viewer,&format,&flg);
694: if (flg) {
695: PetscViewerAndFormat *vf;
696: PetscViewerAndFormatCreate(viewer,format,&vf);
697: PetscObjectDereference((PetscObject)viewer);
698: if (monitorsetup) {
699: (*monitorsetup)(snes,vf);
700: }
701: SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
702: }
703: return(0);
704: }
706: /*@
707: SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
709: Collective on SNES
711: Input Parameter:
712: . snes - the SNES context
714: Options Database Keys:
715: + -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
716: . -snes_stol - convergence tolerance in terms of the norm
717: of the change in the solution between steps
718: . -snes_atol <abstol> - absolute tolerance of residual norm
719: . -snes_rtol <rtol> - relative decrease in tolerance norm from initial
720: . -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
721: . -snes_force_iteration <force> - force SNESSolve() to take at least one iteration
722: . -snes_max_it <max_it> - maximum number of iterations
723: . -snes_max_funcs <max_funcs> - maximum number of function evaluations
724: . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
725: . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
726: . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
727: . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
728: . -snes_trtol <trtol> - trust region tolerance
729: . -snes_no_convergence_test - skip convergence test in nonlinear
730: solver; hence iterations will continue until max_it
731: or some other criterion is reached. Saves expense
732: of convergence test
733: . -snes_monitor [ascii][:filename][:viewer format] - prints residual norm at each iteration. if no filename given prints to stdout
734: . -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration
735: . -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration
736: . -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration
737: . -snes_monitor_lg_residualnorm - plots residual norm at each iteration
738: . -snes_monitor_lg_range - plots residual norm at each iteration
739: . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
740: . -snes_fd_color - use finite differences with coloring to compute Jacobian
741: . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
742: - -snes_converged_reason - print the reason for convergence/divergence after each solve
744: Options Database for Eisenstat-Walker method:
745: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
746: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
747: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
748: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
749: . -snes_ksp_ew_gamma <gamma> - Sets gamma
750: . -snes_ksp_ew_alpha <alpha> - Sets alpha
751: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
752: - -snes_ksp_ew_threshold <threshold> - Sets threshold
754: Notes:
755: To see all options, run your program with the -help option or consult
756: Users-Manual: ch_snes
758: Level: beginner
760: .keywords: SNES, nonlinear, set, options, database
762: .seealso: SNESSetOptionsPrefix()
763: @*/
764: PetscErrorCode SNESSetFromOptions(SNES snes)
765: {
766: PetscBool flg,pcset,persist,set;
767: PetscInt i,indx,lag,grids;
768: const char *deft = SNESNEWTONLS;
769: const char *convtests[] = {"default","skip"};
770: SNESKSPEW *kctx = NULL;
771: char type[256], monfilename[PETSC_MAX_PATH_LEN];
773: PCSide pcside;
774: const char *optionsprefix;
778: SNESRegisterAll();
779: PetscObjectOptionsBegin((PetscObject)snes);
780: if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
781: PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
782: if (flg) {
783: SNESSetType(snes,type);
784: } else if (!((PetscObject)snes)->type_name) {
785: SNESSetType(snes,deft);
786: }
787: PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);
788: PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);
790: PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,NULL);
791: PetscOptionsReal("-snes_divergence_tolerance","Stop if residual norm increases by this factor","SNESSetDivergenceTolerance",snes->divtol,&snes->divtol,NULL);
792: PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,NULL);
793: PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,NULL);
794: PetscOptionsInt("-snes_max_fail","Maximum nonlinear step failures","SNESSetMaxNonlinearStepFailures",snes->maxFailures,&snes->maxFailures,NULL);
795: PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,NULL);
796: PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,NULL);
797: PetscOptionsBool("-snes_force_iteration","Force SNESSolve() to take at least one iteration","SNESSetForceIteration",snes->forceiteration,&snes->forceiteration,NULL);
799: PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
800: if (flg) {
801: SNESSetLagPreconditioner(snes,lag);
802: }
803: PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
804: if (flg) {
805: SNESSetLagPreconditionerPersists(snes,persist);
806: }
807: PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
808: if (flg) {
809: SNESSetLagJacobian(snes,lag);
810: }
811: PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
812: if (flg) {
813: SNESSetLagJacobianPersists(snes,persist);
814: }
816: PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
817: if (flg) {
818: SNESSetGridSequence(snes,grids);
819: }
821: PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);
822: if (flg) {
823: switch (indx) {
824: case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
825: case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL); break;
826: }
827: }
829: PetscOptionsEList("-snes_norm_schedule","SNES Norm schedule","SNESSetNormSchedule",SNESNormSchedules,5,"function",&indx,&flg);
830: if (flg) { SNESSetNormSchedule(snes,(SNESNormSchedule)indx); }
832: PetscOptionsEList("-snes_function_type","SNES Norm schedule","SNESSetFunctionType",SNESFunctionTypes,2,"unpreconditioned",&indx,&flg);
833: if (flg) { SNESSetFunctionType(snes,(SNESFunctionType)indx); }
835: kctx = (SNESKSPEW*)snes->kspconvctx;
837: PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,NULL);
839: PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);
840: PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);
841: PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);
842: PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);
843: PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);
844: PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);
845: PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);
847: flg = PETSC_FALSE;
848: PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);
849: if (set && flg) {SNESMonitorCancel(snes);}
851: SNESMonitorSetFromOptions(snes,"-snes_monitor","Monitor norm of function","SNESMonitorDefault",SNESMonitorDefault,NULL);
852: SNESMonitorSetFromOptions(snes,"-snes_monitor_short","Monitor norm of function with fewer digits","SNESMonitorDefaultShort",SNESMonitorDefaultShort,NULL);
853: SNESMonitorSetFromOptions(snes,"-snes_monitor_range","Monitor range of elements of function","SNESMonitorRange",SNESMonitorRange,NULL);
855: SNESMonitorSetFromOptions(snes,"-snes_monitor_ratio","Monitor ratios of the norm of function for consecutive steps","SNESMonitorRatio",SNESMonitorRatio,SNESMonitorRatioSetUp);
856: SNESMonitorSetFromOptions(snes,"-snes_monitor_field","Monitor norm of function (split into fields)","SNESMonitorDefaultField",SNESMonitorDefaultField,NULL);
857: SNESMonitorSetFromOptions(snes,"-snes_monitor_solution","View solution at each iteration","SNESMonitorSolution",SNESMonitorSolution,NULL);
858: SNESMonitorSetFromOptions(snes,"-snes_monitor_solution_update","View correction at each iteration","SNESMonitorSolutionUpdate",SNESMonitorSolutionUpdate,NULL);
859: SNESMonitorSetFromOptions(snes,"-snes_monitor_residual","View residual at each iteration","SNESMonitorResidual",SNESMonitorResidual,NULL);
860: SNESMonitorSetFromOptions(snes,"-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",SNESMonitorJacUpdateSpectrum,NULL);
861: SNESMonitorSetFromOptions(snes,"-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet",SNESMonitorFields,NULL);
863: PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
864: if (flg) {PetscPythonMonitorSet((PetscObject)snes,monfilename);}
867: flg = PETSC_FALSE;
868: PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);
869: if (flg) {
870: PetscDrawLG ctx;
872: SNESMonitorLGCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
873: SNESMonitorSet(snes,SNESMonitorLGResidualNorm,ctx,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
874: }
875: flg = PETSC_FALSE;
876: PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
877: if (flg) {
878: PetscViewer ctx;
880: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
881: SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
882: }
886: flg = PETSC_FALSE;
887: PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
888: if (flg) {
889: void *functx;
890: DM dm;
891: DMSNES sdm;
892: SNESGetDM(snes,&dm);
893: DMGetDMSNES(dm,&sdm);
894: sdm->jacobianctx = NULL;
895: SNESGetFunction(snes,NULL,NULL,&functx);
896: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
897: PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
898: }
900: flg = PETSC_FALSE;
901: PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
902: if (flg) {
903: SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
904: }
906: flg = PETSC_FALSE;
907: PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);
908: if (flg) {
909: DM dm;
910: DMSNES sdm;
911: SNESGetDM(snes,&dm);
912: DMGetDMSNES(dm,&sdm);
913: sdm->jacobianctx = NULL;
914: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,0);
915: PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
916: }
918: flg = PETSC_FALSE;
919: PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf_operator,&flg);
920: if (flg && snes->mf_operator) {
921: snes->mf_operator = PETSC_TRUE;
922: snes->mf = PETSC_TRUE;
923: }
924: flg = PETSC_FALSE;
925: PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf,&flg);
926: if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE;
927: PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,0);
929: flg = PETSC_FALSE;
930: SNESGetNPCSide(snes,&pcside);
931: PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
932: if (flg) {SNESSetNPCSide(snes,pcside);}
934: #if defined(PETSC_HAVE_SAWS)
935: /*
936: Publish convergence information using SAWs
937: */
938: flg = PETSC_FALSE;
939: PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);
940: if (flg) {
941: void *ctx;
942: SNESMonitorSAWsCreate(snes,&ctx);
943: SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);
944: }
945: #endif
946: #if defined(PETSC_HAVE_SAWS)
947: {
948: PetscBool set;
949: flg = PETSC_FALSE;
950: PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);
951: if (set) {
952: PetscObjectSAWsSetBlock((PetscObject)snes,flg);
953: }
954: }
955: #endif
957: for (i = 0; i < numberofsetfromoptions; i++) {
958: (*othersetfromoptions[i])(snes);
959: }
961: if (snes->ops->setfromoptions) {
962: (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
963: }
965: /* process any options handlers added with PetscObjectAddOptionsHandler() */
966: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)snes);
967: PetscOptionsEnd();
969: if (!snes->linesearch) {
970: SNESGetLineSearch(snes, &snes->linesearch);
971: }
972: SNESLineSearchSetFromOptions(snes->linesearch);
974: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
975: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
976: KSPSetFromOptions(snes->ksp);
978: /* if someone has set the SNES NPC type, create it. */
979: SNESGetOptionsPrefix(snes, &optionsprefix);
980: PetscOptionsHasName(((PetscObject)snes)->options,optionsprefix, "-npc_snes_type", &pcset);
981: if (pcset && (!snes->npc)) {
982: SNESGetNPC(snes, &snes->npc);
983: }
984: return(0);
985: }
987: /*@C
988: SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
989: the nonlinear solvers.
991: Logically Collective on SNES
993: Input Parameters:
994: + snes - the SNES context
995: . compute - function to compute the context
996: - destroy - function to destroy the context
998: Level: intermediate
1000: Notes:
1001: This function is currently not available from Fortran.
1003: .keywords: SNES, nonlinear, set, application, context
1005: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
1006: @*/
1007: PetscErrorCode SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
1008: {
1011: snes->ops->usercompute = compute;
1012: snes->ops->userdestroy = destroy;
1013: return(0);
1014: }
1016: /*@
1017: SNESSetApplicationContext - Sets the optional user-defined context for
1018: the nonlinear solvers.
1020: Logically Collective on SNES
1022: Input Parameters:
1023: + snes - the SNES context
1024: - usrP - optional user context
1026: Level: intermediate
1028: Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
1029: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1031: .keywords: SNES, nonlinear, set, application, context
1033: .seealso: SNESGetApplicationContext()
1034: @*/
1035: PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP)
1036: {
1038: KSP ksp;
1042: SNESGetKSP(snes,&ksp);
1043: KSPSetApplicationContext(ksp,usrP);
1044: snes->user = usrP;
1045: return(0);
1046: }
1048: /*@
1049: SNESGetApplicationContext - Gets the user-defined context for the
1050: nonlinear solvers.
1052: Not Collective
1054: Input Parameter:
1055: . snes - SNES context
1057: Output Parameter:
1058: . usrP - user context
1060: Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
1061: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1063: Level: intermediate
1065: .keywords: SNES, nonlinear, get, application, context
1067: .seealso: SNESSetApplicationContext()
1068: @*/
1069: PetscErrorCode SNESGetApplicationContext(SNES snes,void *usrP)
1070: {
1073: *(void**)usrP = snes->user;
1074: return(0);
1075: }
1077: /*@
1078: SNESSetUseMatrixFree - indicates that SNES should use matrix free finite difference matrix vector products internally to apply
1079: the Jacobian.
1081: Collective on SNES
1083: Input Parameters:
1084: + snes - SNES context
1085: . mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1086: - mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used
1088: Options Database:
1089: + -snes_mf - use matrix free for both the mat and pmat operator
1090: - -snes_mf_operator - use matrix free only for the mat operator
1092: Level: intermediate
1094: .keywords: SNES, nonlinear, get, iteration, number,
1096: .seealso: SNESGetUseMatrixFree(), MatCreateSNESMF()
1097: @*/
1098: PetscErrorCode SNESSetUseMatrixFree(SNES snes,PetscBool mf_operator,PetscBool mf)
1099: {
1104: if (mf && !mf_operator) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"If using mf must also use mf_operator");
1105: snes->mf = mf;
1106: snes->mf_operator = mf_operator;
1107: return(0);
1108: }
1110: /*@
1111: SNESGetUseMatrixFree - indicates if the SNES uses matrix free finite difference matrix vector products to apply
1112: the Jacobian.
1114: Collective on SNES
1116: Input Parameter:
1117: . snes - SNES context
1119: Output Parameters:
1120: + mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1121: - mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used
1123: Options Database:
1124: + -snes_mf - use matrix free for both the mat and pmat operator
1125: - -snes_mf_operator - use matrix free only for the mat operator
1127: Level: intermediate
1129: .keywords: SNES, nonlinear, get, iteration, number,
1131: .seealso: SNESSetUseMatrixFree(), MatCreateSNESMF()
1132: @*/
1133: PetscErrorCode SNESGetUseMatrixFree(SNES snes,PetscBool *mf_operator,PetscBool *mf)
1134: {
1137: if (mf) *mf = snes->mf;
1138: if (mf_operator) *mf_operator = snes->mf_operator;
1139: return(0);
1140: }
1142: /*@
1143: SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1144: at this time.
1146: Not Collective
1148: Input Parameter:
1149: . snes - SNES context
1151: Output Parameter:
1152: . iter - iteration number
1154: Notes:
1155: For example, during the computation of iteration 2 this would return 1.
1157: This is useful for using lagged Jacobians (where one does not recompute the
1158: Jacobian at each SNES iteration). For example, the code
1159: .vb
1160: SNESGetIterationNumber(snes,&it);
1161: if (!(it % 2)) {
1162: [compute Jacobian here]
1163: }
1164: .ve
1165: can be used in your ComputeJacobian() function to cause the Jacobian to be
1166: recomputed every second SNES iteration.
1168: After the SNES solve is complete this will return the number of nonlinear iterations used.
1170: Level: intermediate
1172: .keywords: SNES, nonlinear, get, iteration, number,
1174: .seealso: SNESGetLinearSolveIterations()
1175: @*/
1176: PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt *iter)
1177: {
1181: *iter = snes->iter;
1182: return(0);
1183: }
1185: /*@
1186: SNESSetIterationNumber - Sets the current iteration number.
1188: Not Collective
1190: Input Parameter:
1191: . snes - SNES context
1192: . iter - iteration number
1194: Level: developer
1196: .keywords: SNES, nonlinear, set, iteration, number,
1198: .seealso: SNESGetLinearSolveIterations()
1199: @*/
1200: PetscErrorCode SNESSetIterationNumber(SNES snes,PetscInt iter)
1201: {
1206: PetscObjectSAWsTakeAccess((PetscObject)snes);
1207: snes->iter = iter;
1208: PetscObjectSAWsGrantAccess((PetscObject)snes);
1209: return(0);
1210: }
1212: /*@
1213: SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1214: attempted by the nonlinear solver.
1216: Not Collective
1218: Input Parameter:
1219: . snes - SNES context
1221: Output Parameter:
1222: . nfails - number of unsuccessful steps attempted
1224: Notes:
1225: This counter is reset to zero for each successive call to SNESSolve().
1227: Level: intermediate
1229: .keywords: SNES, nonlinear, get, number, unsuccessful, steps
1231: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1232: SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1233: @*/
1234: PetscErrorCode SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1235: {
1239: *nfails = snes->numFailures;
1240: return(0);
1241: }
1243: /*@
1244: SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1245: attempted by the nonlinear solver before it gives up.
1247: Not Collective
1249: Input Parameters:
1250: + snes - SNES context
1251: - maxFails - maximum of unsuccessful steps
1253: Level: intermediate
1255: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
1257: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1258: SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1259: @*/
1260: PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1261: {
1264: snes->maxFailures = maxFails;
1265: return(0);
1266: }
1268: /*@
1269: SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1270: attempted by the nonlinear solver before it gives up.
1272: Not Collective
1274: Input Parameter:
1275: . snes - SNES context
1277: Output Parameter:
1278: . maxFails - maximum of unsuccessful steps
1280: Level: intermediate
1282: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1284: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1285: SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1287: @*/
1288: PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1289: {
1293: *maxFails = snes->maxFailures;
1294: return(0);
1295: }
1297: /*@
1298: SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1299: done by SNES.
1301: Not Collective
1303: Input Parameter:
1304: . snes - SNES context
1306: Output Parameter:
1307: . nfuncs - number of evaluations
1309: Level: intermediate
1311: Notes: Reset every time SNESSolve is called unless SNESSetCountersReset() is used.
1313: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1315: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1316: @*/
1317: PetscErrorCode SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1318: {
1322: *nfuncs = snes->nfuncs;
1323: return(0);
1324: }
1326: /*@
1327: SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1328: linear solvers.
1330: Not Collective
1332: Input Parameter:
1333: . snes - SNES context
1335: Output Parameter:
1336: . nfails - number of failed solves
1338: Level: intermediate
1340: Options Database Keys:
1341: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1343: Notes:
1344: This counter is reset to zero for each successive call to SNESSolve().
1346: .keywords: SNES, nonlinear, get, number, unsuccessful, steps
1348: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1349: @*/
1350: PetscErrorCode SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1351: {
1355: *nfails = snes->numLinearSolveFailures;
1356: return(0);
1357: }
1359: /*@
1360: SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1361: allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
1363: Logically Collective on SNES
1365: Input Parameters:
1366: + snes - SNES context
1367: - maxFails - maximum allowed linear solve failures
1369: Level: intermediate
1371: Options Database Keys:
1372: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1374: Notes: By default this is 0; that is SNES returns on the first failed linear solve
1376: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
1378: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1379: @*/
1380: PetscErrorCode SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1381: {
1385: snes->maxLinearSolveFailures = maxFails;
1386: return(0);
1387: }
1389: /*@
1390: SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1391: are allowed before SNES terminates
1393: Not Collective
1395: Input Parameter:
1396: . snes - SNES context
1398: Output Parameter:
1399: . maxFails - maximum of unsuccessful solves allowed
1401: Level: intermediate
1403: Notes: By default this is 1; that is SNES returns on the first failed linear solve
1405: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1407: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1408: @*/
1409: PetscErrorCode SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1410: {
1414: *maxFails = snes->maxLinearSolveFailures;
1415: return(0);
1416: }
1418: /*@
1419: SNESGetLinearSolveIterations - Gets the total number of linear iterations
1420: used by the nonlinear solver.
1422: Not Collective
1424: Input Parameter:
1425: . snes - SNES context
1427: Output Parameter:
1428: . lits - number of linear iterations
1430: Notes:
1431: This counter is reset to zero for each successive call to SNESSolve() unless SNESSetCountersReset() is used.
1433: If the linear solver fails inside the SNESSolve() the iterations for that call to the linear solver are not included. If you wish to count them
1434: then call KSPGetIterationNumber() after the failed solve.
1436: Level: intermediate
1438: .keywords: SNES, nonlinear, get, number, linear, iterations
1440: .seealso: SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1441: @*/
1442: PetscErrorCode SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1443: {
1447: *lits = snes->linear_its;
1448: return(0);
1449: }
1451: /*@
1452: SNESSetCountersReset - Sets whether or not the counters for linear iterations and function evaluations
1453: are reset every time SNESSolve() is called.
1455: Logically Collective on SNES
1457: Input Parameter:
1458: + snes - SNES context
1459: - reset - whether to reset the counters or not
1461: Notes:
1462: This defaults to PETSC_TRUE
1464: Level: developer
1466: .keywords: SNES, nonlinear, set, reset, number, linear, iterations
1468: .seealso: SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC()
1469: @*/
1470: PetscErrorCode SNESSetCountersReset(SNES snes,PetscBool reset)
1471: {
1475: snes->counters_reset = reset;
1476: return(0);
1477: }
1480: /*@
1481: SNESSetKSP - Sets a KSP context for the SNES object to use
1483: Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
1485: Input Parameters:
1486: + snes - the SNES context
1487: - ksp - the KSP context
1489: Notes:
1490: The SNES object already has its KSP object, you can obtain with SNESGetKSP()
1491: so this routine is rarely needed.
1493: The KSP object that is already in the SNES object has its reference count
1494: decreased by one.
1496: Level: developer
1498: .keywords: SNES, nonlinear, get, KSP, context
1500: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1501: @*/
1502: PetscErrorCode SNESSetKSP(SNES snes,KSP ksp)
1503: {
1510: PetscObjectReference((PetscObject)ksp);
1511: if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1512: snes->ksp = ksp;
1513: return(0);
1514: }
1516: /* -----------------------------------------------------------*/
1517: /*@
1518: SNESCreate - Creates a nonlinear solver context.
1520: Collective on MPI_Comm
1522: Input Parameters:
1523: . comm - MPI communicator
1525: Output Parameter:
1526: . outsnes - the new SNES context
1528: Options Database Keys:
1529: + -snes_mf - Activates default matrix-free Jacobian-vector products,
1530: and no preconditioning matrix
1531: . -snes_mf_operator - Activates default matrix-free Jacobian-vector
1532: products, and a user-provided preconditioning matrix
1533: as set by SNESSetJacobian()
1534: - -snes_fd - Uses (slow!) finite differences to compute Jacobian
1536: Level: beginner
1538: Developer Notes: SNES always creates a KSP object even though many SNES methods do not use it. This is
1539: unfortunate and should be fixed at some point. The flag snes->usesksp indicates if the
1540: particular method does use KSP and regulates if the information about the KSP is printed
1541: in SNESView(). TSSetFromOptions() does call SNESSetFromOptions() which can lead to users being confused
1542: by help messages about meaningless SNES options.
1544: SNES always creates the snes->kspconvctx even though it is used by only one type. This should
1545: be fixed.
1547: .keywords: SNES, nonlinear, create, context
1549: .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()
1551: @*/
1552: PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes)
1553: {
1555: SNES snes;
1556: SNESKSPEW *kctx;
1560: *outsnes = NULL;
1561: SNESInitializePackage();
1563: PetscHeaderCreate(snes,SNES_CLASSID,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);
1565: snes->ops->converged = SNESConvergedDefault;
1566: snes->usesksp = PETSC_TRUE;
1567: snes->tolerancesset = PETSC_FALSE;
1568: snes->max_its = 50;
1569: snes->max_funcs = 10000;
1570: snes->norm = 0.0;
1571: snes->normschedule = SNES_NORM_ALWAYS;
1572: snes->functype = SNES_FUNCTION_DEFAULT;
1573: #if defined(PETSC_USE_REAL_SINGLE)
1574: snes->rtol = 1.e-5;
1575: #else
1576: snes->rtol = 1.e-8;
1577: #endif
1578: snes->ttol = 0.0;
1579: #if defined(PETSC_USE_REAL_SINGLE)
1580: snes->abstol = 1.e-25;
1581: #else
1582: snes->abstol = 1.e-50;
1583: #endif
1584: #if defined(PETSC_USE_REAL_SINGLE)
1585: snes->stol = 1.e-5;
1586: #else
1587: snes->stol = 1.e-8;
1588: #endif
1589: #if defined(PETSC_USE_REAL_SINGLE)
1590: snes->deltatol = 1.e-6;
1591: #else
1592: snes->deltatol = 1.e-12;
1593: #endif
1594: snes->divtol = 1.e4;
1595: snes->rnorm0 = 0;
1596: snes->nfuncs = 0;
1597: snes->numFailures = 0;
1598: snes->maxFailures = 1;
1599: snes->linear_its = 0;
1600: snes->lagjacobian = 1;
1601: snes->jac_iter = 0;
1602: snes->lagjac_persist = PETSC_FALSE;
1603: snes->lagpreconditioner = 1;
1604: snes->pre_iter = 0;
1605: snes->lagpre_persist = PETSC_FALSE;
1606: snes->numbermonitors = 0;
1607: snes->data = 0;
1608: snes->setupcalled = PETSC_FALSE;
1609: snes->ksp_ewconv = PETSC_FALSE;
1610: snes->nwork = 0;
1611: snes->work = 0;
1612: snes->nvwork = 0;
1613: snes->vwork = 0;
1614: snes->conv_hist_len = 0;
1615: snes->conv_hist_max = 0;
1616: snes->conv_hist = NULL;
1617: snes->conv_hist_its = NULL;
1618: snes->conv_hist_reset = PETSC_TRUE;
1619: snes->counters_reset = PETSC_TRUE;
1620: snes->vec_func_init_set = PETSC_FALSE;
1621: snes->reason = SNES_CONVERGED_ITERATING;
1622: snes->npcside = PC_RIGHT;
1624: snes->mf = PETSC_FALSE;
1625: snes->mf_operator = PETSC_FALSE;
1626: snes->mf_version = 1;
1628: snes->numLinearSolveFailures = 0;
1629: snes->maxLinearSolveFailures = 1;
1631: snes->vizerotolerance = 1.e-8;
1633: /* Set this to true if the implementation of SNESSolve_XXX does compute the residual at the final solution. */
1634: snes->alwayscomputesfinalresidual = PETSC_FALSE;
1636: /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
1637: PetscNewLog(snes,&kctx);
1639: snes->kspconvctx = (void*)kctx;
1640: kctx->version = 2;
1641: kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1642: this was too large for some test cases */
1643: kctx->rtol_last = 0.0;
1644: kctx->rtol_max = .9;
1645: kctx->gamma = 1.0;
1646: kctx->alpha = .5*(1.0 + PetscSqrtReal(5.0));
1647: kctx->alpha2 = kctx->alpha;
1648: kctx->threshold = .1;
1649: kctx->lresid_last = 0.0;
1650: kctx->norm_last = 0.0;
1652: *outsnes = snes;
1653: return(0);
1654: }
1656: /*MC
1657: SNESFunction - Functional form used to convey the nonlinear function to be solved by SNES
1659: Synopsis:
1660: #include "petscsnes.h"
1661: PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);
1663: Input Parameters:
1664: + snes - the SNES context
1665: . x - state at which to evaluate residual
1666: - ctx - optional user-defined function context, passed in with SNESSetFunction()
1668: Output Parameter:
1669: . f - vector to put residual (function value)
1671: Level: intermediate
1673: .seealso: SNESSetFunction(), SNESGetFunction()
1674: M*/
1676: /*@C
1677: SNESSetFunction - Sets the function evaluation routine and function
1678: vector for use by the SNES routines in solving systems of nonlinear
1679: equations.
1681: Logically Collective on SNES
1683: Input Parameters:
1684: + snes - the SNES context
1685: . r - vector to store function value
1686: . f - function evaluation routine; see SNESFunction for calling sequence details
1687: - ctx - [optional] user-defined context for private data for the
1688: function evaluation routine (may be NULL)
1690: Notes:
1691: The Newton-like methods typically solve linear systems of the form
1692: $ f'(x) x = -f(x),
1693: where f'(x) denotes the Jacobian matrix and f(x) is the function.
1695: Level: beginner
1697: .keywords: SNES, nonlinear, set, function
1699: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1700: @*/
1701: PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1702: {
1704: DM dm;
1708: if (r) {
1711: PetscObjectReference((PetscObject)r);
1712: VecDestroy(&snes->vec_func);
1714: snes->vec_func = r;
1715: }
1716: SNESGetDM(snes,&dm);
1717: DMSNESSetFunction(dm,f,ctx);
1718: return(0);
1719: }
1722: /*@C
1723: SNESSetInitialFunction - Sets the function vector to be used as the
1724: function norm at the initialization of the method. In some
1725: instances, the user has precomputed the function before calling
1726: SNESSolve. This function allows one to avoid a redundant call
1727: to SNESComputeFunction in that case.
1729: Logically Collective on SNES
1731: Input Parameters:
1732: + snes - the SNES context
1733: - f - vector to store function value
1735: Notes:
1736: This should not be modified during the solution procedure.
1738: This is used extensively in the SNESFAS hierarchy and in nonlinear preconditioning.
1740: Level: developer
1742: .keywords: SNES, nonlinear, set, function
1744: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1745: @*/
1746: PetscErrorCode SNESSetInitialFunction(SNES snes, Vec f)
1747: {
1749: Vec vec_func;
1755: if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1756: snes->vec_func_init_set = PETSC_FALSE;
1757: return(0);
1758: }
1759: SNESGetFunction(snes,&vec_func,NULL,NULL);
1760: VecCopy(f, vec_func);
1762: snes->vec_func_init_set = PETSC_TRUE;
1763: return(0);
1764: }
1766: /*@
1767: SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1768: of the SNES method.
1770: Logically Collective on SNES
1772: Input Parameters:
1773: + snes - the SNES context
1774: - normschedule - the frequency of norm computation
1776: Options Database Key:
1777: . -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly>
1779: Notes:
1780: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1781: of the nonlinear function and the taking of its norm at every iteration to
1782: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1783: (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
1784: may either be monitored for convergence or not. As these are often used as nonlinear
1785: preconditioners, monitoring the norm of their error is not a useful enterprise within
1786: their solution.
1788: Level: developer
1790: .keywords: SNES, nonlinear, set, function, norm, type
1792: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1793: @*/
1794: PetscErrorCode SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1795: {
1798: snes->normschedule = normschedule;
1799: return(0);
1800: }
1803: /*@
1804: SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1805: of the SNES method.
1807: Logically Collective on SNES
1809: Input Parameters:
1810: + snes - the SNES context
1811: - normschedule - the type of the norm used
1813: Level: advanced
1815: .keywords: SNES, nonlinear, set, function, norm, type
1817: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1818: @*/
1819: PetscErrorCode SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1820: {
1823: *normschedule = snes->normschedule;
1824: return(0);
1825: }
1828: /*@
1829: SNESSetFunctionNorm - Sets the last computed residual norm.
1831: Logically Collective on SNES
1833: Input Parameters:
1834: + snes - the SNES context
1836: - normschedule - the frequency of norm computation
1838: Level: developer
1840: .keywords: SNES, nonlinear, set, function, norm, type
1841: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1842: @*/
1843: PetscErrorCode SNESSetFunctionNorm(SNES snes, PetscReal norm)
1844: {
1847: snes->norm = norm;
1848: return(0);
1849: }
1851: /*@
1852: SNESGetFunctionNorm - Gets the last computed norm of the residual
1854: Not Collective
1856: Input Parameter:
1857: . snes - the SNES context
1859: Output Parameter:
1860: . norm - the last computed residual norm
1862: Level: developer
1864: .keywords: SNES, nonlinear, set, function, norm, type
1865: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1866: @*/
1867: PetscErrorCode SNESGetFunctionNorm(SNES snes, PetscReal *norm)
1868: {
1872: *norm = snes->norm;
1873: return(0);
1874: }
1876: /*@C
1877: SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
1878: of the SNES method.
1880: Logically Collective on SNES
1882: Input Parameters:
1883: + snes - the SNES context
1884: - normschedule - the frequency of norm computation
1886: Notes:
1887: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1888: of the nonlinear function and the taking of its norm at every iteration to
1889: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1890: (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
1891: may either be monitored for convergence or not. As these are often used as nonlinear
1892: preconditioners, monitoring the norm of their error is not a useful enterprise within
1893: their solution.
1895: Level: developer
1897: .keywords: SNES, nonlinear, set, function, norm, type
1899: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1900: @*/
1901: PetscErrorCode SNESSetFunctionType(SNES snes, SNESFunctionType type)
1902: {
1905: snes->functype = type;
1906: return(0);
1907: }
1910: /*@C
1911: SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
1912: of the SNES method.
1914: Logically Collective on SNES
1916: Input Parameters:
1917: + snes - the SNES context
1918: - normschedule - the type of the norm used
1920: Level: advanced
1922: .keywords: SNES, nonlinear, set, function, norm, type
1924: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1925: @*/
1926: PetscErrorCode SNESGetFunctionType(SNES snes, SNESFunctionType *type)
1927: {
1930: *type = snes->functype;
1931: return(0);
1932: }
1934: /*MC
1935: SNESNGSFunction - function used to convey a Gauss-Seidel sweep on the nonlinear function
1937: Synopsis:
1938: #include <petscsnes.h>
1939: $ SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);
1941: + X - solution vector
1942: . B - RHS vector
1943: - ctx - optional user-defined Gauss-Seidel context
1945: Level: intermediate
1947: .seealso: SNESSetNGS(), SNESGetNGS()
1948: M*/
1950: /*@C
1951: SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
1952: use with composed nonlinear solvers.
1954: Input Parameters:
1955: + snes - the SNES context
1956: . f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
1957: - ctx - [optional] user-defined context for private data for the
1958: smoother evaluation routine (may be NULL)
1960: Notes:
1961: The NGS routines are used by the composed nonlinear solver to generate
1962: a problem appropriate update to the solution, particularly FAS.
1964: Level: intermediate
1966: .keywords: SNES, nonlinear, set, Gauss-Seidel
1968: .seealso: SNESGetFunction(), SNESComputeNGS()
1969: @*/
1970: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1971: {
1973: DM dm;
1977: SNESGetDM(snes,&dm);
1978: DMSNESSetNGS(dm,f,ctx);
1979: return(0);
1980: }
1982: PETSC_EXTERN PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
1983: {
1985: DM dm;
1986: DMSNES sdm;
1989: SNESGetDM(snes,&dm);
1990: DMGetDMSNES(dm,&sdm);
1991: /* A(x)*x - b(x) */
1992: if (sdm->ops->computepfunction) {
1993: (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
1994: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");
1996: if (sdm->ops->computepjacobian) {
1997: (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
1998: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard matrix.");
1999: VecScale(f,-1.0);
2000: MatMultAdd(snes->jacobian,x,f,f);
2001: return(0);
2002: }
2004: PETSC_EXTERN PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
2005: {
2007: /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
2008: return(0);
2009: }
2011: /*@C
2012: SNESSetPicard - Use SNES to solve the semilinear-system A(x) x = b(x) via a Picard type iteration (Picard linearization)
2014: Logically Collective on SNES
2016: Input Parameters:
2017: + snes - the SNES context
2018: . r - vector to store function value
2019: . b - function evaluation routine
2020: . Amat - matrix with which A(x) x - b(x) is to be computed
2021: . Pmat - matrix from which preconditioner is computed (usually the same as Amat)
2022: . J - function to compute matrix value, see SNESJacobianFunction for details on its calling sequence
2023: - ctx - [optional] user-defined context for private data for the
2024: function evaluation routine (may be NULL)
2026: Notes:
2027: 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
2028: 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.
2030: One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both
2032: $ 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}
2033: $ Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration.
2035: Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner.
2037: We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then
2038: the direct Picard iteration A(x^n) x^{n+1} = b(x^n)
2040: 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
2041: 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
2042: different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).
2044: Level: intermediate
2046: .keywords: SNES, nonlinear, set, function
2048: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
2049: @*/
2050: 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)
2051: {
2053: DM dm;
2057: SNESGetDM(snes, &dm);
2058: DMSNESSetPicard(dm,b,J,ctx);
2059: SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
2060: SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
2061: return(0);
2062: }
2064: /*@C
2065: SNESGetPicard - Returns the context for the Picard iteration
2067: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
2069: Input Parameter:
2070: . snes - the SNES context
2072: Output Parameter:
2073: + r - the function (or NULL)
2074: . f - the function (or NULL); see SNESFunction for calling sequence details
2075: . Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
2076: . Pmat - the matrix from which the preconditioner will be constructed (or NULL)
2077: . J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
2078: - ctx - the function context (or NULL)
2080: Level: advanced
2082: .keywords: SNES, nonlinear, get, function
2084: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
2085: @*/
2086: 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)
2087: {
2089: DM dm;
2093: SNESGetFunction(snes,r,NULL,NULL);
2094: SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
2095: SNESGetDM(snes,&dm);
2096: DMSNESGetPicard(dm,f,J,ctx);
2097: return(0);
2098: }
2100: /*@C
2101: SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem
2103: Logically Collective on SNES
2105: Input Parameters:
2106: + snes - the SNES context
2107: . func - function evaluation routine
2108: - ctx - [optional] user-defined context for private data for the
2109: function evaluation routine (may be NULL)
2111: Calling sequence of func:
2112: $ func (SNES snes,Vec x,void *ctx);
2114: . f - function vector
2115: - ctx - optional user-defined function context
2117: Level: intermediate
2119: .keywords: SNES, nonlinear, set, function
2121: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
2122: @*/
2123: PetscErrorCode SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
2124: {
2127: if (func) snes->ops->computeinitialguess = func;
2128: if (ctx) snes->initialguessP = ctx;
2129: return(0);
2130: }
2132: /* --------------------------------------------------------------- */
2133: /*@C
2134: SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
2135: it assumes a zero right hand side.
2137: Logically Collective on SNES
2139: Input Parameter:
2140: . snes - the SNES context
2142: Output Parameter:
2143: . rhs - the right hand side vector or NULL if the right hand side vector is null
2145: Level: intermediate
2147: .keywords: SNES, nonlinear, get, function, right hand side
2149: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
2150: @*/
2151: PetscErrorCode SNESGetRhs(SNES snes,Vec *rhs)
2152: {
2156: *rhs = snes->vec_rhs;
2157: return(0);
2158: }
2160: /*@
2161: SNESComputeFunction - Calls the function that has been set with SNESSetFunction().
2163: Collective on SNES
2165: Input Parameters:
2166: + snes - the SNES context
2167: - x - input vector
2169: Output Parameter:
2170: . y - function vector, as set by SNESSetFunction()
2172: Notes:
2173: SNESComputeFunction() is typically used within nonlinear solvers
2174: implementations, so most users would not generally call this routine
2175: themselves.
2177: Level: developer
2179: .keywords: SNES, nonlinear, compute, function
2181: .seealso: SNESSetFunction(), SNESGetFunction()
2182: @*/
2183: PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y)
2184: {
2186: DM dm;
2187: DMSNES sdm;
2195: VecValidValues(x,2,PETSC_TRUE);
2197: SNESGetDM(snes,&dm);
2198: DMGetDMSNES(dm,&sdm);
2199: if (sdm->ops->computefunction) {
2200: if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2201: PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2202: }
2203: VecLockPush(x);
2204: PetscStackPush("SNES user function");
2205: (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2206: PetscStackPop;
2207: VecLockPop(x);
2208: if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2209: PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2210: }
2211: } else if (snes->vec_rhs) {
2212: MatMult(snes->jacobian, x, y);
2213: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2214: if (snes->vec_rhs) {
2215: VecAXPY(y,-1.0,snes->vec_rhs);
2216: }
2217: snes->nfuncs++;
2218: /*
2219: domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2220: propagate the value to all processes
2221: */
2222: if (snes->domainerror) {
2223: VecSetInf(y);
2224: }
2225: return(0);
2226: }
2228: /*@
2229: SNESComputeNGS - Calls the Gauss-Seidel function that has been set with SNESSetNGS().
2231: Collective on SNES
2233: Input Parameters:
2234: + snes - the SNES context
2235: . x - input vector
2236: - b - rhs vector
2238: Output Parameter:
2239: . x - new solution vector
2241: Notes:
2242: SNESComputeNGS() is typically used within composed nonlinear solver
2243: implementations, so most users would not generally call this routine
2244: themselves.
2246: Level: developer
2248: .keywords: SNES, nonlinear, compute, function
2250: .seealso: SNESSetNGS(), SNESComputeFunction()
2251: @*/
2252: PetscErrorCode SNESComputeNGS(SNES snes,Vec b,Vec x)
2253: {
2255: DM dm;
2256: DMSNES sdm;
2264: if (b) {VecValidValues(b,2,PETSC_TRUE);}
2265: PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2266: SNESGetDM(snes,&dm);
2267: DMGetDMSNES(dm,&sdm);
2268: if (sdm->ops->computegs) {
2269: if (b) {VecLockPush(b);}
2270: PetscStackPush("SNES user NGS");
2271: (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2272: PetscStackPop;
2273: if (b) {VecLockPop(b);}
2274: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2275: PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2276: return(0);
2277: }
2279: PetscErrorCode SNESTestJacobian(SNES snes)
2280: {
2281: Mat A,B,C,D,jacobian;
2282: Vec x = snes->vec_sol,f = snes->vec_func;
2283: PetscErrorCode ierr;
2284: PetscReal nrm,gnorm;
2285: PetscReal threshold = 1.e-5;
2286: PetscInt m,n,M,N;
2287: void *functx;
2288: PetscBool complete_print = PETSC_FALSE,threshold_print = PETSC_FALSE,test = PETSC_FALSE,flg;
2289: PetscViewer viewer,mviewer;
2290: MPI_Comm comm;
2291: PetscInt tabs;
2292: static PetscBool directionsprinted = PETSC_FALSE;
2293: PetscViewerFormat format;
2296: PetscObjectOptionsBegin((PetscObject)snes);
2297: PetscOptionsName("-snes_test_jacobian","Compare hand-coded and finite difference Jacobians","None",&test);
2298: PetscOptionsReal("-snes_test_jacobian", "Threshold for element difference between hand-coded and finite difference being meaningful", "None", threshold, &threshold,NULL);
2299: PetscOptionsViewer("-snes_test_jacobian_view","View difference between hand-coded and finite difference Jacobians element entries","None",&mviewer,&format,&complete_print);
2300: if (!complete_print) {
2301: PetscOptionsViewer("-snes_test_jacobian_display","Display difference between hand-coded and finite difference Jacobians","None",&mviewer,&format,&complete_print);
2302: }
2303: /* for compatibility with PETSc 3.9 and older. */
2304: PetscOptionsReal("-snes_test_jacobian_display_threshold", "Display difference between hand-coded and finite difference Jacobians which exceed input threshold", "None", threshold, &threshold, &threshold_print);
2305: PetscOptionsEnd();
2306: if (!test) return(0);
2308: PetscObjectGetComm((PetscObject)snes,&comm);
2309: PetscViewerASCIIGetStdout(comm,&viewer);
2310: PetscViewerASCIIGetTab(viewer, &tabs);
2311: PetscViewerASCIISetTab(viewer, ((PetscObject)snes)->tablevel);
2312: PetscViewerASCIIPrintf(viewer," ---------- Testing Jacobian -------------\n");
2313: if (!complete_print && !directionsprinted) {
2314: PetscViewerASCIIPrintf(viewer," Run with -snes_test_jacobian_view and optionally -snes_test_jacobian <threshold> to show difference\n");
2315: PetscViewerASCIIPrintf(viewer," of hand-coded and finite difference Jacobian entries greater than <threshold>.\n");
2316: }
2317: if (!directionsprinted) {
2318: PetscViewerASCIIPrintf(viewer," Testing hand-coded Jacobian, if (for double precision runs) ||J - Jfd||_F/||J||_F is\n");
2319: PetscViewerASCIIPrintf(viewer," O(1.e-8), the hand-coded Jacobian is probably correct.\n");
2320: directionsprinted = PETSC_TRUE;
2321: }
2322: if (complete_print) {
2323: PetscViewerPushFormat(mviewer,format);
2324: }
2326: /* evaluate the function at this point because SNESComputeJacobianDefault() assumes that the function has been evaluated and put into snes->vec_func */
2327: SNESComputeFunction(snes,x,f);
2329: PetscObjectTypeCompare((PetscObject)snes->jacobian,MATMFFD,&flg);
2330: if (!flg) jacobian = snes->jacobian;
2331: else jacobian = snes->jacobian_pre;
2333: while (jacobian) {
2334: PetscObjectTypeCompareAny((PetscObject)jacobian,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2335: if (flg) {
2336: A = jacobian;
2337: PetscObjectReference((PetscObject)A);
2338: } else {
2339: MatComputeExplicitOperator(jacobian,&A);
2340: }
2342: MatCreate(PetscObjectComm((PetscObject)A),&B);
2343: MatGetSize(A,&M,&N);
2344: MatGetLocalSize(A,&m,&n);
2345: MatSetSizes(B,m,n,M,N);
2346: MatSetType(B,((PetscObject)A)->type_name);
2347: MatSetUp(B);
2348: MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
2350: SNESGetFunction(snes,NULL,NULL,&functx);
2351: SNESComputeJacobianDefault(snes,x,B,B,functx);
2353: MatDuplicate(B,MAT_COPY_VALUES,&D);
2354: MatAYPX(D,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2355: MatNorm(D,NORM_FROBENIUS,&nrm);
2356: MatNorm(A,NORM_FROBENIUS,&gnorm);
2357: MatDestroy(&D);
2358: if (!gnorm) gnorm = 1; /* just in case */
2359: PetscViewerASCIIPrintf(viewer," ||J - Jfd||_F/||J||_F = %g, ||J - Jfd||_F = %g\n",(double)(nrm/gnorm),(double)nrm);
2361: if (complete_print) {
2362: PetscViewerASCIIPrintf(viewer," Hand-coded Jacobian ----------\n");
2363: MatView(jacobian,mviewer);
2364: PetscViewerASCIIPrintf(viewer," Finite difference Jacobian ----------\n");
2365: MatView(B,mviewer);
2366: }
2368: if (threshold_print) {
2369: PetscInt Istart, Iend, *ccols, bncols, cncols, j, row;
2370: PetscScalar *cvals;
2371: const PetscInt *bcols;
2372: const PetscScalar *bvals;
2374: MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);
2375: MatCreate(PetscObjectComm((PetscObject)A),&C);
2376: MatSetSizes(C,m,n,M,N);
2377: MatSetType(C,((PetscObject)A)->type_name);
2378: MatSetUp(C);
2379: MatSetOption(C,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
2380: MatGetOwnershipRange(B,&Istart,&Iend);
2382: for (row = Istart; row < Iend; row++) {
2383: MatGetRow(B,row,&bncols,&bcols,&bvals);
2384: PetscMalloc2(bncols,&ccols,bncols,&cvals);
2385: for (j = 0, cncols = 0; j < bncols; j++) {
2386: if (PetscAbsScalar(bvals[j]) > threshold) {
2387: ccols[cncols] = bcols[j];
2388: cvals[cncols] = bvals[j];
2389: cncols += 1;
2390: }
2391: }
2392: if (cncols) {
2393: MatSetValues(C,1,&row,cncols,ccols,cvals,INSERT_VALUES);
2394: }
2395: MatRestoreRow(B,row,&bncols,&bcols,&bvals);
2396: PetscFree2(ccols,cvals);
2397: }
2398: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2399: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2400: PetscViewerASCIIPrintf(viewer," Hand-coded minus finite-difference Jacobian with tolerance %g ----------\n",(double)threshold);
2401: MatView(C,complete_print ? mviewer : viewer);
2402: MatDestroy(&C);
2403: }
2404: MatDestroy(&A);
2405: MatDestroy(&B);
2407: if (jacobian != snes->jacobian_pre) {
2408: jacobian = snes->jacobian_pre;
2409: PetscViewerASCIIPrintf(viewer," ---------- Testing Jacobian for preconditioner -------------\n");
2410: }
2411: else jacobian = NULL;
2412: }
2413: if (complete_print) {
2414: PetscViewerPopFormat(mviewer);
2415: }
2416: PetscViewerASCIISetTab(viewer,tabs);
2417: return(0);
2418: }
2420: /*@
2421: SNESComputeJacobian - Computes the Jacobian matrix that has been set with SNESSetJacobian().
2423: Collective on SNES and Mat
2425: Input Parameters:
2426: + snes - the SNES context
2427: - x - input vector
2429: Output Parameters:
2430: + A - Jacobian matrix
2431: - B - optional preconditioning matrix
2433: Options Database Keys:
2434: + -snes_lag_preconditioner <lag>
2435: . -snes_lag_jacobian <lag>
2436: . -snes_test_jacobian - compare the user provided Jacobian with one compute via finite differences to check for errors
2437: . -snes_test_jacobian_display - display the user provided Jacobian, the finite difference Jacobian and the difference between them to help users detect the location of errors in the user provided Jacobian
2438: . -snes_test_jacobian_display_threshold <numerical value> - display entries in the difference between the user provided Jacobian and finite difference Jacobian that are greater than a certain value to help users detect errors
2439: . -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2440: . -snes_compare_explicit_draw - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2441: . -snes_compare_explicit_contour - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2442: . -snes_compare_operator - Make the comparison options above use the operator instead of the preconditioning matrix
2443: . -snes_compare_coloring - Compute the finite difference Jacobian using coloring and display norms of difference
2444: . -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
2445: . -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2446: . -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2447: . -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2448: . -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
2449: - -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences
2452: Notes:
2453: Most users should not need to explicitly call this routine, as it
2454: is used internally within the nonlinear solvers.
2456: Developer Notes: This has duplicative ways of checking the accuracy of the user provided Jacobian (see the options above). This is for historical reasons, the routine SNESTestJacobian() use to used
2457: for with the SNESType of test that has been removed.
2459: Level: developer
2461: .keywords: SNES, compute, Jacobian, matrix
2463: .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2464: @*/
2465: PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2466: {
2468: PetscBool flag;
2469: DM dm;
2470: DMSNES sdm;
2471: KSP ksp;
2477: VecValidValues(X,2,PETSC_TRUE);
2478: SNESGetDM(snes,&dm);
2479: DMGetDMSNES(dm,&sdm);
2481: if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_USER,"Must call SNESSetJacobian(), DMSNESSetJacobian(), DMDASNESSetJacobianLocal(), etc");
2483: /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
2485: if (snes->lagjacobian == -2) {
2486: snes->lagjacobian = -1;
2488: PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2489: } else if (snes->lagjacobian == -1) {
2490: PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2491: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2492: if (flag) {
2493: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2494: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2495: }
2496: return(0);
2497: } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2498: PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2499: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2500: if (flag) {
2501: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2502: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2503: }
2504: return(0);
2505: }
2506: if (snes->npc && snes->npcside== PC_LEFT) {
2507: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2508: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2509: return(0);
2510: }
2512: PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2513: VecLockPush(X);
2514: PetscStackPush("SNES user Jacobian function");
2515: (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2516: PetscStackPop;
2517: VecLockPop(X);
2518: PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);
2520: /* the next line ensures that snes->ksp exists */
2521: SNESGetKSP(snes,&ksp);
2522: if (snes->lagpreconditioner == -2) {
2523: PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2524: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2525: snes->lagpreconditioner = -1;
2526: } else if (snes->lagpreconditioner == -1) {
2527: PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2528: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2529: } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2530: PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2531: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2532: } else {
2533: PetscInfo(snes,"Rebuilding preconditioner\n");
2534: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2535: }
2537: SNESTestJacobian(snes);
2538: /* make sure user returned a correct Jacobian and preconditioner */
2541: {
2542: PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
2543: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_explicit",NULL,NULL,&flag);
2544: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",NULL,NULL,&flag_draw);
2545: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",NULL,NULL,&flag_contour);
2546: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_operator",NULL,NULL,&flag_operator);
2547: if (flag || flag_draw || flag_contour) {
2548: Mat Bexp_mine = NULL,Bexp,FDexp;
2549: PetscViewer vdraw,vstdout;
2550: PetscBool flg;
2551: if (flag_operator) {
2552: MatComputeExplicitOperator(A,&Bexp_mine);
2553: Bexp = Bexp_mine;
2554: } else {
2555: /* See if the preconditioning matrix can be viewed and added directly */
2556: PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2557: if (flg) Bexp = B;
2558: else {
2559: /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
2560: MatComputeExplicitOperator(B,&Bexp_mine);
2561: Bexp = Bexp_mine;
2562: }
2563: }
2564: MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);
2565: SNESComputeJacobianDefault(snes,X,FDexp,FDexp,NULL);
2566: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2567: if (flag_draw || flag_contour) {
2568: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2569: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2570: } else vdraw = NULL;
2571: PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator ? "Jacobian" : "preconditioning Jacobian");
2572: if (flag) {MatView(Bexp,vstdout);}
2573: if (vdraw) {MatView(Bexp,vdraw);}
2574: PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");
2575: if (flag) {MatView(FDexp,vstdout);}
2576: if (vdraw) {MatView(FDexp,vdraw);}
2577: MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);
2578: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");
2579: if (flag) {MatView(FDexp,vstdout);}
2580: if (vdraw) { /* Always use contour for the difference */
2581: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2582: MatView(FDexp,vdraw);
2583: PetscViewerPopFormat(vdraw);
2584: }
2585: if (flag_contour) {PetscViewerPopFormat(vdraw);}
2586: PetscViewerDestroy(&vdraw);
2587: MatDestroy(&Bexp_mine);
2588: MatDestroy(&FDexp);
2589: }
2590: }
2591: {
2592: PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
2593: PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
2594: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring",NULL,NULL,&flag);
2595: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring_display",NULL,NULL,&flag_display);
2596: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",NULL,NULL,&flag_draw);
2597: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",NULL,NULL,&flag_contour);
2598: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",NULL,NULL,&flag_threshold);
2599: if (flag_threshold) {
2600: PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,NULL);
2601: PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,NULL);
2602: }
2603: if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
2604: Mat Bfd;
2605: PetscViewer vdraw,vstdout;
2606: MatColoring coloring;
2607: ISColoring iscoloring;
2608: MatFDColoring matfdcoloring;
2609: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2610: void *funcctx;
2611: PetscReal norm1,norm2,normmax;
2613: MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2614: MatColoringCreate(Bfd,&coloring);
2615: MatColoringSetType(coloring,MATCOLORINGSL);
2616: MatColoringSetFromOptions(coloring);
2617: MatColoringApply(coloring,&iscoloring);
2618: MatColoringDestroy(&coloring);
2619: MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2620: MatFDColoringSetFromOptions(matfdcoloring);
2621: MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2622: ISColoringDestroy(&iscoloring);
2624: /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2625: SNESGetFunction(snes,NULL,&func,&funcctx);
2626: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2627: PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2628: PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2629: MatFDColoringSetFromOptions(matfdcoloring);
2630: MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2631: MatFDColoringDestroy(&matfdcoloring);
2633: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2634: if (flag_draw || flag_contour) {
2635: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2636: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2637: } else vdraw = NULL;
2638: PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2639: if (flag_display) {MatView(B,vstdout);}
2640: if (vdraw) {MatView(B,vdraw);}
2641: PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2642: if (flag_display) {MatView(Bfd,vstdout);}
2643: if (vdraw) {MatView(Bfd,vdraw);}
2644: MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2645: MatNorm(Bfd,NORM_1,&norm1);
2646: MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2647: MatNorm(Bfd,NORM_MAX,&normmax);
2648: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2649: if (flag_display) {MatView(Bfd,vstdout);}
2650: if (vdraw) { /* Always use contour for the difference */
2651: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2652: MatView(Bfd,vdraw);
2653: PetscViewerPopFormat(vdraw);
2654: }
2655: if (flag_contour) {PetscViewerPopFormat(vdraw);}
2657: if (flag_threshold) {
2658: PetscInt bs,rstart,rend,i;
2659: MatGetBlockSize(B,&bs);
2660: MatGetOwnershipRange(B,&rstart,&rend);
2661: for (i=rstart; i<rend; i++) {
2662: const PetscScalar *ba,*ca;
2663: const PetscInt *bj,*cj;
2664: PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
2665: PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0;
2666: MatGetRow(B,i,&bn,&bj,&ba);
2667: MatGetRow(Bfd,i,&cn,&cj,&ca);
2668: if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
2669: for (j=0; j<bn; j++) {
2670: PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2671: if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
2672: maxentrycol = bj[j];
2673: maxentry = PetscRealPart(ba[j]);
2674: }
2675: if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
2676: maxdiffcol = bj[j];
2677: maxdiff = PetscRealPart(ca[j]);
2678: }
2679: if (rdiff > maxrdiff) {
2680: maxrdiffcol = bj[j];
2681: maxrdiff = rdiff;
2682: }
2683: }
2684: if (maxrdiff > 1) {
2685: 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);
2686: for (j=0; j<bn; j++) {
2687: PetscReal rdiff;
2688: rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2689: if (rdiff > 1) {
2690: PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
2691: }
2692: }
2693: PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
2694: }
2695: MatRestoreRow(B,i,&bn,&bj,&ba);
2696: MatRestoreRow(Bfd,i,&cn,&cj,&ca);
2697: }
2698: }
2699: PetscViewerDestroy(&vdraw);
2700: MatDestroy(&Bfd);
2701: }
2702: }
2703: return(0);
2704: }
2706: /*MC
2707: SNESJacobianFunction - Function used to convey the nonlinear Jacobian of the function to be solved by SNES
2709: Synopsis:
2710: #include "petscsnes.h"
2711: PetscErrorCode SNESJacobianFunction(SNES snes,Vec x,Mat Amat,Mat Pmat,void *ctx);
2713: + x - input vector
2714: . Amat - the matrix that defines the (approximate) Jacobian
2715: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2716: - ctx - [optional] user-defined Jacobian context
2718: Level: intermediate
2720: .seealso: SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2721: M*/
2723: /*@C
2724: SNESSetJacobian - Sets the function to compute Jacobian as well as the
2725: location to store the matrix.
2727: Logically Collective on SNES and Mat
2729: Input Parameters:
2730: + snes - the SNES context
2731: . Amat - the matrix that defines the (approximate) Jacobian
2732: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2733: . J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
2734: - ctx - [optional] user-defined context for private data for the
2735: Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)
2737: Notes:
2738: If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2739: each matrix.
2741: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
2742: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
2744: If using SNESComputeJacobianDefaultColor() to assemble a Jacobian, the ctx argument
2745: must be a MatFDColoring.
2747: Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian. One common
2748: example is to use the "Picard linearization" which only differentiates through the highest order parts of each term.
2750: Level: beginner
2752: .keywords: SNES, nonlinear, set, Jacobian, matrix
2754: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J,
2755: SNESSetPicard(), SNESJacobianFunction
2756: @*/
2757: PetscErrorCode SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2758: {
2760: DM dm;
2768: SNESGetDM(snes,&dm);
2769: DMSNESSetJacobian(dm,J,ctx);
2770: if (Amat) {
2771: PetscObjectReference((PetscObject)Amat);
2772: MatDestroy(&snes->jacobian);
2774: snes->jacobian = Amat;
2775: }
2776: if (Pmat) {
2777: PetscObjectReference((PetscObject)Pmat);
2778: MatDestroy(&snes->jacobian_pre);
2780: snes->jacobian_pre = Pmat;
2781: }
2782: return(0);
2783: }
2785: /*@C
2786: SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2787: provided context for evaluating the Jacobian.
2789: Not Collective, but Mat object will be parallel if SNES object is
2791: Input Parameter:
2792: . snes - the nonlinear solver context
2794: Output Parameters:
2795: + Amat - location to stash (approximate) Jacobian matrix (or NULL)
2796: . Pmat - location to stash matrix used to compute the preconditioner (or NULL)
2797: . J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
2798: - ctx - location to stash Jacobian ctx (or NULL)
2800: Level: advanced
2802: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
2803: @*/
2804: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2805: {
2807: DM dm;
2808: DMSNES sdm;
2812: if (Amat) *Amat = snes->jacobian;
2813: if (Pmat) *Pmat = snes->jacobian_pre;
2814: SNESGetDM(snes,&dm);
2815: DMGetDMSNES(dm,&sdm);
2816: if (J) *J = sdm->ops->computejacobian;
2817: if (ctx) *ctx = sdm->jacobianctx;
2818: return(0);
2819: }
2821: /*@
2822: SNESSetUp - Sets up the internal data structures for the later use
2823: of a nonlinear solver.
2825: Collective on SNES
2827: Input Parameters:
2828: . snes - the SNES context
2830: Notes:
2831: For basic use of the SNES solvers the user need not explicitly call
2832: SNESSetUp(), since these actions will automatically occur during
2833: the call to SNESSolve(). However, if one wishes to control this
2834: phase separately, SNESSetUp() should be called after SNESCreate()
2835: and optional routines of the form SNESSetXXX(), but before SNESSolve().
2837: Level: advanced
2839: .keywords: SNES, nonlinear, setup
2841: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
2842: @*/
2843: PetscErrorCode SNESSetUp(SNES snes)
2844: {
2846: DM dm;
2847: DMSNES sdm;
2848: SNESLineSearch linesearch, pclinesearch;
2849: void *lsprectx,*lspostctx;
2850: PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
2851: PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
2852: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2853: Vec f,fpc;
2854: void *funcctx;
2855: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2856: void *jacctx,*appctx;
2857: Mat j,jpre;
2861: if (snes->setupcalled) return(0);
2863: if (!((PetscObject)snes)->type_name) {
2864: SNESSetType(snes,SNESNEWTONLS);
2865: }
2867: SNESGetFunction(snes,&snes->vec_func,NULL,NULL);
2869: SNESGetDM(snes,&dm);
2870: DMGetDMSNES(dm,&sdm);
2871: if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
2872: if (!sdm->ops->computejacobian) {
2873: DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
2874: }
2875: if (!snes->vec_func) {
2876: DMCreateGlobalVector(dm,&snes->vec_func);
2877: }
2879: if (!snes->ksp) {
2880: SNESGetKSP(snes, &snes->ksp);
2881: }
2883: if (!snes->linesearch) {
2884: SNESGetLineSearch(snes, &snes->linesearch);
2885: }
2886: SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);
2888: if (snes->npc && (snes->npcside== PC_LEFT)) {
2889: snes->mf = PETSC_TRUE;
2890: snes->mf_operator = PETSC_FALSE;
2891: }
2893: if (snes->npc) {
2894: /* copy the DM over */
2895: SNESGetDM(snes,&dm);
2896: SNESSetDM(snes->npc,dm);
2898: SNESGetFunction(snes,&f,&func,&funcctx);
2899: VecDuplicate(f,&fpc);
2900: SNESSetFunction(snes->npc,fpc,func,funcctx);
2901: SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
2902: SNESSetJacobian(snes->npc,j,jpre,jac,jacctx);
2903: SNESGetApplicationContext(snes,&appctx);
2904: SNESSetApplicationContext(snes->npc,appctx);
2905: VecDestroy(&fpc);
2907: /* copy the function pointers over */
2908: PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)snes->npc);
2910: /* default to 1 iteration */
2911: SNESSetTolerances(snes->npc,0.0,0.0,0.0,1,snes->npc->max_funcs);
2912: if (snes->npcside==PC_RIGHT) {
2913: SNESSetNormSchedule(snes->npc,SNES_NORM_FINAL_ONLY);
2914: } else {
2915: SNESSetNormSchedule(snes->npc,SNES_NORM_NONE);
2916: }
2917: SNESSetFromOptions(snes->npc);
2919: /* copy the line search context over */
2920: SNESGetLineSearch(snes,&linesearch);
2921: SNESGetLineSearch(snes->npc,&pclinesearch);
2922: SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
2923: SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
2924: SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
2925: SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
2926: PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
2927: }
2928: if (snes->mf) {
2929: SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
2930: }
2931: if (snes->ops->usercompute && !snes->user) {
2932: (*snes->ops->usercompute)(snes,(void**)&snes->user);
2933: }
2935: snes->jac_iter = 0;
2936: snes->pre_iter = 0;
2938: if (snes->ops->setup) {
2939: (*snes->ops->setup)(snes);
2940: }
2942: if (snes->npc && (snes->npcside== PC_LEFT)) {
2943: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2944: SNESGetLineSearch(snes,&linesearch);
2945: SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
2946: }
2947: }
2949: snes->setupcalled = PETSC_TRUE;
2950: return(0);
2951: }
2953: /*@
2954: SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats
2956: Collective on SNES
2958: Input Parameter:
2959: . snes - iterative context obtained from SNESCreate()
2961: Level: intermediate
2963: Notes: Also calls the application context destroy routine set with SNESSetComputeApplicationContext()
2965: .keywords: SNES, destroy
2967: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
2968: @*/
2969: PetscErrorCode SNESReset(SNES snes)
2970: {
2975: if (snes->ops->userdestroy && snes->user) {
2976: (*snes->ops->userdestroy)((void**)&snes->user);
2977: snes->user = NULL;
2978: }
2979: if (snes->npc) {
2980: SNESReset(snes->npc);
2981: }
2983: if (snes->ops->reset) {
2984: (*snes->ops->reset)(snes);
2985: }
2986: if (snes->ksp) {
2987: KSPReset(snes->ksp);
2988: }
2990: if (snes->linesearch) {
2991: SNESLineSearchReset(snes->linesearch);
2992: }
2994: VecDestroy(&snes->vec_rhs);
2995: VecDestroy(&snes->vec_sol);
2996: VecDestroy(&snes->vec_sol_update);
2997: VecDestroy(&snes->vec_func);
2998: MatDestroy(&snes->jacobian);
2999: MatDestroy(&snes->jacobian_pre);
3000: VecDestroyVecs(snes->nwork,&snes->work);
3001: VecDestroyVecs(snes->nvwork,&snes->vwork);
3003: snes->alwayscomputesfinalresidual = PETSC_FALSE;
3005: snes->nwork = snes->nvwork = 0;
3006: snes->setupcalled = PETSC_FALSE;
3007: return(0);
3008: }
3010: /*@
3011: SNESDestroy - Destroys the nonlinear solver context that was created
3012: with SNESCreate().
3014: Collective on SNES
3016: Input Parameter:
3017: . snes - the SNES context
3019: Level: beginner
3021: .keywords: SNES, nonlinear, destroy
3023: .seealso: SNESCreate(), SNESSolve()
3024: @*/
3025: PetscErrorCode SNESDestroy(SNES *snes)
3026: {
3030: if (!*snes) return(0);
3032: if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; return(0);}
3034: SNESReset((*snes));
3035: SNESDestroy(&(*snes)->npc);
3037: /* if memory was published with SAWs then destroy it */
3038: PetscObjectSAWsViewOff((PetscObject)*snes);
3039: if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}
3041: if ((*snes)->dm) {DMCoarsenHookRemove((*snes)->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,*snes);}
3042: DMDestroy(&(*snes)->dm);
3043: KSPDestroy(&(*snes)->ksp);
3044: SNESLineSearchDestroy(&(*snes)->linesearch);
3046: PetscFree((*snes)->kspconvctx);
3047: if ((*snes)->ops->convergeddestroy) {
3048: (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
3049: }
3050: if ((*snes)->conv_malloc) {
3051: PetscFree((*snes)->conv_hist);
3052: PetscFree((*snes)->conv_hist_its);
3053: }
3054: SNESMonitorCancel((*snes));
3055: PetscHeaderDestroy(snes);
3056: return(0);
3057: }
3059: /* ----------- Routines to set solver parameters ---------- */
3061: /*@
3062: SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
3064: Logically Collective on SNES
3066: Input Parameters:
3067: + snes - the SNES context
3068: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3069: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
3071: Options Database Keys:
3072: . -snes_lag_preconditioner <lag>
3074: Notes:
3075: The default is 1
3076: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3077: If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
3079: Level: intermediate
3081: .keywords: SNES, nonlinear, set, convergence, tolerances
3083: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
3085: @*/
3086: PetscErrorCode SNESSetLagPreconditioner(SNES snes,PetscInt lag)
3087: {
3090: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3091: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3093: snes->lagpreconditioner = lag;
3094: return(0);
3095: }
3097: /*@
3098: SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does
3100: Logically Collective on SNES
3102: Input Parameters:
3103: + snes - the SNES context
3104: - steps - the number of refinements to do, defaults to 0
3106: Options Database Keys:
3107: . -snes_grid_sequence <steps>
3109: Level: intermediate
3111: Notes:
3112: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
3114: .keywords: SNES, nonlinear, set, convergence, tolerances
3116: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetGridSequence()
3118: @*/
3119: PetscErrorCode SNESSetGridSequence(SNES snes,PetscInt steps)
3120: {
3124: snes->gridsequence = steps;
3125: return(0);
3126: }
3128: /*@
3129: SNESGetGridSequence - gets the number of steps of grid sequencing that SNES does
3131: Logically Collective on SNES
3133: Input Parameter:
3134: . snes - the SNES context
3136: Output Parameter:
3137: . steps - the number of refinements to do, defaults to 0
3139: Options Database Keys:
3140: . -snes_grid_sequence <steps>
3142: Level: intermediate
3144: Notes:
3145: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
3147: .keywords: SNES, nonlinear, set, convergence, tolerances
3149: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESSetGridSequence()
3151: @*/
3152: PetscErrorCode SNESGetGridSequence(SNES snes,PetscInt *steps)
3153: {
3156: *steps = snes->gridsequence;
3157: return(0);
3158: }
3160: /*@
3161: SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
3163: Not Collective
3165: Input Parameter:
3166: . snes - the SNES context
3168: Output Parameter:
3169: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3170: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
3172: Options Database Keys:
3173: . -snes_lag_preconditioner <lag>
3175: Notes:
3176: The default is 1
3177: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3179: Level: intermediate
3181: .keywords: SNES, nonlinear, set, convergence, tolerances
3183: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
3185: @*/
3186: PetscErrorCode SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
3187: {
3190: *lag = snes->lagpreconditioner;
3191: return(0);
3192: }
3194: /*@
3195: SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
3196: often the preconditioner is rebuilt.
3198: Logically Collective on SNES
3200: Input Parameters:
3201: + snes - the SNES context
3202: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3203: the Jacobian is built etc. -2 means rebuild at next chance but then never again
3205: Options Database Keys:
3206: . -snes_lag_jacobian <lag>
3208: Notes:
3209: The default is 1
3210: The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3211: 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
3212: at the next Newton step but never again (unless it is reset to another value)
3214: Level: intermediate
3216: .keywords: SNES, nonlinear, set, convergence, tolerances
3218: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
3220: @*/
3221: PetscErrorCode SNESSetLagJacobian(SNES snes,PetscInt lag)
3222: {
3225: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3226: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3228: snes->lagjacobian = lag;
3229: return(0);
3230: }
3232: /*@
3233: SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
3235: Not Collective
3237: Input Parameter:
3238: . snes - the SNES context
3240: Output Parameter:
3241: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3242: the Jacobian is built etc.
3244: Options Database Keys:
3245: . -snes_lag_jacobian <lag>
3247: Notes:
3248: The default is 1
3249: The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3251: Level: intermediate
3253: .keywords: SNES, nonlinear, set, convergence, tolerances
3255: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
3257: @*/
3258: PetscErrorCode SNESGetLagJacobian(SNES snes,PetscInt *lag)
3259: {
3262: *lag = snes->lagjacobian;
3263: return(0);
3264: }
3266: /*@
3267: SNESSetLagJacobianPersists - Set whether or not the Jacobian lagging persists through multiple solves
3269: Logically collective on SNES
3271: Input Parameter:
3272: + snes - the SNES context
3273: - flg - jacobian lagging persists if true
3275: Options Database Keys:
3276: . -snes_lag_jacobian_persists <flg>
3278: Notes: This is useful both for nonlinear preconditioning, where it's appropriate to have the Jacobian be stale by
3279: several solves, and for implicit time-stepping, where Jacobian lagging in the inner nonlinear solve over several
3280: timesteps may present huge efficiency gains.
3282: Level: developer
3284: .keywords: SNES, nonlinear, lag
3286: .seealso: SNESSetLagPreconditionerPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()
3288: @*/
3289: PetscErrorCode SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
3290: {
3294: snes->lagjac_persist = flg;
3295: return(0);
3296: }
3298: /*@
3299: SNESSetLagPreconditionerPersists - Set whether or not the preconditioner lagging persists through multiple solves
3301: Logically Collective on SNES
3303: Input Parameter:
3304: + snes - the SNES context
3305: - flg - preconditioner lagging persists if true
3307: Options Database Keys:
3308: . -snes_lag_jacobian_persists <flg>
3310: Notes: This is useful both for nonlinear preconditioning, where it's appropriate to have the preconditioner be stale
3311: by several solves, and for implicit time-stepping, where preconditioner lagging in the inner nonlinear solve over
3312: several timesteps may present huge efficiency gains.
3314: Level: developer
3316: .keywords: SNES, nonlinear, lag
3318: .seealso: SNESSetLagJacobianPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()
3320: @*/
3321: PetscErrorCode SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3322: {
3326: snes->lagpre_persist = flg;
3327: return(0);
3328: }
3330: /*@
3331: SNESSetForceIteration - force SNESSolve() to take at least one iteration regardless of the initial residual norm
3333: Logically Collective on SNES
3335: Input Parameters:
3336: + snes - the SNES context
3337: - force - PETSC_TRUE require at least one iteration
3339: Options Database Keys:
3340: . -snes_force_iteration <force> - Sets forcing an iteration
3342: Notes:
3343: This is used sometimes with TS to prevent TS from detecting a false steady state solution
3345: Level: intermediate
3347: .keywords: SNES, nonlinear, set, convergence, tolerances
3349: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3350: @*/
3351: PetscErrorCode SNESSetForceIteration(SNES snes,PetscBool force)
3352: {
3355: snes->forceiteration = force;
3356: return(0);
3357: }
3359: /*@
3360: SNESGetForceIteration - Whether or not to force SNESSolve() take at least one iteration regardless of the initial residual norm
3362: Logically Collective on SNES
3364: Input Parameters:
3365: . snes - the SNES context
3367: Output Parameter:
3368: . force - PETSC_TRUE requires at least one iteration.
3370: .keywords: SNES, nonlinear, set, convergence, tolerances
3372: .seealso: SNESSetForceIteration(), SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3373: @*/
3374: PetscErrorCode SNESGetForceIteration(SNES snes,PetscBool *force)
3375: {
3378: *force = snes->forceiteration;
3379: return(0);
3380: }
3382: /*@
3383: SNESSetTolerances - Sets various parameters used in convergence tests.
3385: Logically Collective on SNES
3387: Input Parameters:
3388: + snes - the SNES context
3389: . abstol - absolute convergence tolerance
3390: . rtol - relative convergence tolerance
3391: . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x ||
3392: . maxit - maximum number of iterations
3393: - maxf - maximum number of function evaluations
3395: Options Database Keys:
3396: + -snes_atol <abstol> - Sets abstol
3397: . -snes_rtol <rtol> - Sets rtol
3398: . -snes_stol <stol> - Sets stol
3399: . -snes_max_it <maxit> - Sets maxit
3400: - -snes_max_funcs <maxf> - Sets maxf
3402: Notes:
3403: The default maximum number of iterations is 50.
3404: The default maximum number of function evaluations is 1000.
3406: Level: intermediate
3408: .keywords: SNES, nonlinear, set, convergence, tolerances
3410: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance(), SNESSetForceIteration()
3411: @*/
3412: PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3413: {
3422: if (abstol != PETSC_DEFAULT) {
3423: if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3424: snes->abstol = abstol;
3425: }
3426: if (rtol != PETSC_DEFAULT) {
3427: 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);
3428: snes->rtol = rtol;
3429: }
3430: if (stol != PETSC_DEFAULT) {
3431: if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3432: snes->stol = stol;
3433: }
3434: if (maxit != PETSC_DEFAULT) {
3435: if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3436: snes->max_its = maxit;
3437: }
3438: if (maxf != PETSC_DEFAULT) {
3439: if (maxf < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
3440: snes->max_funcs = maxf;
3441: }
3442: snes->tolerancesset = PETSC_TRUE;
3443: return(0);
3444: }
3446: /*@
3447: SNESSetDivergenceTolerance - Sets the divergence tolerance used for the SNES divergence test.
3449: Logically Collective on SNES
3451: Input Parameters:
3452: + snes - the SNES context
3453: - divtol - the divergence tolerance. Use -1 to deactivate the test.
3455: Options Database Keys:
3456: + -snes_divergence_tolerance <divtol> - Sets divtol
3458: Notes:
3459: The default divergence tolerance is 1e4.
3461: Level: intermediate
3463: .keywords: SNES, nonlinear, set, divergence, tolerance
3465: .seealso: SNESSetTolerances(), SNESGetDivergenceTolerance
3466: @*/
3467: PetscErrorCode SNESSetDivergenceTolerance(SNES snes,PetscReal divtol)
3468: {
3473: if (divtol != PETSC_DEFAULT) {
3474: snes->divtol = divtol;
3475: }
3476: else {
3477: snes->divtol = 1.0e4;
3478: }
3479: return(0);
3480: }
3482: /*@
3483: SNESGetTolerances - Gets various parameters used in convergence tests.
3485: Not Collective
3487: Input Parameters:
3488: + snes - the SNES context
3489: . atol - absolute convergence tolerance
3490: . rtol - relative convergence tolerance
3491: . stol - convergence tolerance in terms of the norm
3492: of the change in the solution between steps
3493: . maxit - maximum number of iterations
3494: - maxf - maximum number of function evaluations
3496: Notes:
3497: The user can specify NULL for any parameter that is not needed.
3499: Level: intermediate
3501: .keywords: SNES, nonlinear, get, convergence, tolerances
3503: .seealso: SNESSetTolerances()
3504: @*/
3505: PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3506: {
3509: if (atol) *atol = snes->abstol;
3510: if (rtol) *rtol = snes->rtol;
3511: if (stol) *stol = snes->stol;
3512: if (maxit) *maxit = snes->max_its;
3513: if (maxf) *maxf = snes->max_funcs;
3514: return(0);
3515: }
3517: /*@
3518: SNESGetDivergenceTolerance - Gets divergence tolerance used in divergence test.
3520: Not Collective
3522: Input Parameters:
3523: + snes - the SNES context
3524: - divtol - divergence tolerance
3526: Level: intermediate
3528: .keywords: SNES, nonlinear, get, divergence, tolerance
3530: .seealso: SNESSetDivergenceTolerance()
3531: @*/
3532: PetscErrorCode SNESGetDivergenceTolerance(SNES snes,PetscReal *divtol)
3533: {
3536: if (divtol) *divtol = snes->divtol;
3537: return(0);
3538: }
3540: /*@
3541: SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
3543: Logically Collective on SNES
3545: Input Parameters:
3546: + snes - the SNES context
3547: - tol - tolerance
3549: Options Database Key:
3550: . -snes_trtol <tol> - Sets tol
3552: Level: intermediate
3554: .keywords: SNES, nonlinear, set, trust region, tolerance
3556: .seealso: SNESSetTolerances()
3557: @*/
3558: PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3559: {
3563: snes->deltatol = tol;
3564: return(0);
3565: }
3567: /*
3568: Duplicate the lg monitors for SNES from KSP; for some reason with
3569: dynamic libraries things don't work under Sun4 if we just use
3570: macros instead of functions
3571: */
3572: PetscErrorCode SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,void *ctx)
3573: {
3578: KSPMonitorLGResidualNorm((KSP)snes,it,norm,ctx);
3579: return(0);
3580: }
3582: PetscErrorCode SNESMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)
3583: {
3587: KSPMonitorLGResidualNormCreate(comm,host,label,x,y,m,n,lgctx);
3588: return(0);
3589: }
3591: PETSC_INTERN PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3593: PetscErrorCode SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3594: {
3595: PetscDrawLG lg;
3596: PetscErrorCode ierr;
3597: PetscReal x,y,per;
3598: PetscViewer v = (PetscViewer)monctx;
3599: static PetscReal prev; /* should be in the context */
3600: PetscDraw draw;
3604: PetscViewerDrawGetDrawLG(v,0,&lg);
3605: if (!n) {PetscDrawLGReset(lg);}
3606: PetscDrawLGGetDraw(lg,&draw);
3607: PetscDrawSetTitle(draw,"Residual norm");
3608: x = (PetscReal)n;
3609: if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3610: else y = -15.0;
3611: PetscDrawLGAddPoint(lg,&x,&y);
3612: if (n < 20 || !(n % 5) || snes->reason) {
3613: PetscDrawLGDraw(lg);
3614: PetscDrawLGSave(lg);
3615: }
3617: PetscViewerDrawGetDrawLG(v,1,&lg);
3618: if (!n) {PetscDrawLGReset(lg);}
3619: PetscDrawLGGetDraw(lg,&draw);
3620: PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3621: SNESMonitorRange_Private(snes,n,&per);
3622: x = (PetscReal)n;
3623: y = 100.0*per;
3624: PetscDrawLGAddPoint(lg,&x,&y);
3625: if (n < 20 || !(n % 5) || snes->reason) {
3626: PetscDrawLGDraw(lg);
3627: PetscDrawLGSave(lg);
3628: }
3630: PetscViewerDrawGetDrawLG(v,2,&lg);
3631: if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3632: PetscDrawLGGetDraw(lg,&draw);
3633: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3634: x = (PetscReal)n;
3635: y = (prev - rnorm)/prev;
3636: PetscDrawLGAddPoint(lg,&x,&y);
3637: if (n < 20 || !(n % 5) || snes->reason) {
3638: PetscDrawLGDraw(lg);
3639: PetscDrawLGSave(lg);
3640: }
3642: PetscViewerDrawGetDrawLG(v,3,&lg);
3643: if (!n) {PetscDrawLGReset(lg);}
3644: PetscDrawLGGetDraw(lg,&draw);
3645: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3646: x = (PetscReal)n;
3647: y = (prev - rnorm)/(prev*per);
3648: if (n > 2) { /*skip initial crazy value */
3649: PetscDrawLGAddPoint(lg,&x,&y);
3650: }
3651: if (n < 20 || !(n % 5) || snes->reason) {
3652: PetscDrawLGDraw(lg);
3653: PetscDrawLGSave(lg);
3654: }
3655: prev = rnorm;
3656: return(0);
3657: }
3659: /*@
3660: SNESMonitor - runs the user provided monitor routines, if they exist
3662: Collective on SNES
3664: Input Parameters:
3665: + snes - nonlinear solver context obtained from SNESCreate()
3666: . iter - iteration number
3667: - rnorm - relative norm of the residual
3669: Notes:
3670: This routine is called by the SNES implementations.
3671: It does not typically need to be called by the user.
3673: Level: developer
3675: .seealso: SNESMonitorSet()
3676: @*/
3677: PetscErrorCode SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3678: {
3680: PetscInt i,n = snes->numbermonitors;
3683: VecLockPush(snes->vec_sol);
3684: for (i=0; i<n; i++) {
3685: (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3686: }
3687: VecLockPop(snes->vec_sol);
3688: return(0);
3689: }
3691: /* ------------ Routines to set performance monitoring options ----------- */
3693: /*MC
3694: SNESMonitorFunction - functional form passed to SNESMonitorSet() to monitor convergence of nonlinear solver
3696: Synopsis:
3697: #include <petscsnes.h>
3698: $ PetscErrorCode SNESMonitorFunction(SNES snes,PetscInt its, PetscReal norm,void *mctx)
3700: + snes - the SNES context
3701: . its - iteration number
3702: . norm - 2-norm function value (may be estimated)
3703: - mctx - [optional] monitoring context
3705: Level: advanced
3707: .seealso: SNESMonitorSet(), SNESMonitorGet()
3708: M*/
3710: /*@C
3711: SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3712: iteration of the nonlinear solver to display the iteration's
3713: progress.
3715: Logically Collective on SNES
3717: Input Parameters:
3718: + snes - the SNES context
3719: . f - the monitor function, see SNESMonitorFunction for the calling sequence
3720: . mctx - [optional] user-defined context for private data for the
3721: monitor routine (use NULL if no context is desired)
3722: - monitordestroy - [optional] routine that frees monitor context
3723: (may be NULL)
3725: Options Database Keys:
3726: + -snes_monitor - sets SNESMonitorDefault()
3727: . -snes_monitor_lg_residualnorm - sets line graph monitor,
3728: uses SNESMonitorLGCreate()
3729: - -snes_monitor_cancel - cancels all monitors that have
3730: been hardwired into a code by
3731: calls to SNESMonitorSet(), but
3732: does not cancel those set via
3733: the options database.
3735: Notes:
3736: Several different monitoring routines may be set by calling
3737: SNESMonitorSet() multiple times; all will be called in the
3738: order in which they were set.
3740: Fortran notes: Only a single monitor function can be set for each SNES object
3742: Level: intermediate
3744: .keywords: SNES, nonlinear, set, monitor
3746: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3747: @*/
3748: PetscErrorCode SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3749: {
3750: PetscInt i;
3752: PetscBool identical;
3756: for (i=0; i<snes->numbermonitors;i++) {
3757: PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))snes->monitor[i],snes->monitorcontext[i],snes->monitordestroy[i],&identical);
3758: if (identical) return(0);
3759: }
3760: if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3761: snes->monitor[snes->numbermonitors] = f;
3762: snes->monitordestroy[snes->numbermonitors] = monitordestroy;
3763: snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3764: return(0);
3765: }
3767: /*@
3768: SNESMonitorCancel - Clears all the monitor functions for a SNES object.
3770: Logically Collective on SNES
3772: Input Parameters:
3773: . snes - the SNES context
3775: Options Database Key:
3776: . -snes_monitor_cancel - cancels all monitors that have been hardwired
3777: into a code by calls to SNESMonitorSet(), but does not cancel those
3778: set via the options database
3780: Notes:
3781: There is no way to clear one specific monitor from a SNES object.
3783: Level: intermediate
3785: .keywords: SNES, nonlinear, set, monitor
3787: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3788: @*/
3789: PetscErrorCode SNESMonitorCancel(SNES snes)
3790: {
3792: PetscInt i;
3796: for (i=0; i<snes->numbermonitors; i++) {
3797: if (snes->monitordestroy[i]) {
3798: (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3799: }
3800: }
3801: snes->numbermonitors = 0;
3802: return(0);
3803: }
3805: /*MC
3806: SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver
3808: Synopsis:
3809: #include <petscsnes.h>
3810: $ PetscErrorCode SNESConvergenceTest(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
3812: + snes - the SNES context
3813: . it - current iteration (0 is the first and is before any Newton step)
3814: . cctx - [optional] convergence context
3815: . reason - reason for convergence/divergence
3816: . xnorm - 2-norm of current iterate
3817: . gnorm - 2-norm of current step
3818: - f - 2-norm of function
3820: Level: intermediate
3822: .seealso: SNESSetConvergenceTest(), SNESGetConvergenceTest()
3823: M*/
3825: /*@C
3826: SNESSetConvergenceTest - Sets the function that is to be used
3827: to test for convergence of the nonlinear iterative solution.
3829: Logically Collective on SNES
3831: Input Parameters:
3832: + snes - the SNES context
3833: . SNESConvergenceTestFunction - routine to test for convergence
3834: . cctx - [optional] context for private data for the convergence routine (may be NULL)
3835: - destroy - [optional] destructor for the context (may be NULL; NULL_FUNCTION in Fortran)
3837: Level: advanced
3839: .keywords: SNES, nonlinear, set, convergence, test
3841: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3842: @*/
3843: PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3844: {
3849: if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
3850: if (snes->ops->convergeddestroy) {
3851: (*snes->ops->convergeddestroy)(snes->cnvP);
3852: }
3853: snes->ops->converged = SNESConvergenceTestFunction;
3854: snes->ops->convergeddestroy = destroy;
3855: snes->cnvP = cctx;
3856: return(0);
3857: }
3859: /*@
3860: SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
3862: Not Collective
3864: Input Parameter:
3865: . snes - the SNES context
3867: Output Parameter:
3868: . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
3869: manual pages for the individual convergence tests for complete lists
3871: Options Database:
3872: . -snes_converged_reason - prints the reason to standard out
3874: Level: intermediate
3876: Notes: Should only be called after the call the SNESSolve() is complete, if it is called earlier it returns the value SNES__CONVERGED_ITERATING.
3878: .keywords: SNES, nonlinear, set, convergence, test
3880: .seealso: SNESSetConvergenceTest(), SNESSetConvergedReason(), SNESConvergedReason
3881: @*/
3882: PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
3883: {
3887: *reason = snes->reason;
3888: return(0);
3889: }
3891: /*@
3892: SNESSetConvergedReason - Sets the reason the SNES iteration was stopped.
3894: Not Collective
3896: Input Parameters:
3897: + snes - the SNES context
3898: - reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
3899: manual pages for the individual convergence tests for complete lists
3901: Level: intermediate
3903: .keywords: SNES, nonlinear, set, convergence, test
3904: .seealso: SNESGetConvergedReason(), SNESSetConvergenceTest(), SNESConvergedReason
3905: @*/
3906: PetscErrorCode SNESSetConvergedReason(SNES snes,SNESConvergedReason reason)
3907: {
3910: snes->reason = reason;
3911: return(0);
3912: }
3914: /*@
3915: SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
3917: Logically Collective on SNES
3919: Input Parameters:
3920: + snes - iterative context obtained from SNESCreate()
3921: . a - array to hold history, this array will contain the function norms computed at each step
3922: . its - integer array holds the number of linear iterations for each solve.
3923: . na - size of a and its
3924: - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
3925: else it continues storing new values for new nonlinear solves after the old ones
3927: Notes:
3928: If 'a' and 'its' are NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
3929: default array of length 10000 is allocated.
3931: This routine is useful, e.g., when running a code for purposes
3932: of accurate performance monitoring, when no I/O should be done
3933: during the section of code that is being timed.
3935: Level: intermediate
3937: .keywords: SNES, set, convergence, history
3939: .seealso: SNESGetConvergenceHistory()
3941: @*/
3942: PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
3943: {
3950: if (!a) {
3951: if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
3952: PetscCalloc1(na,&a);
3953: PetscCalloc1(na,&its);
3955: snes->conv_malloc = PETSC_TRUE;
3956: }
3957: snes->conv_hist = a;
3958: snes->conv_hist_its = its;
3959: snes->conv_hist_max = na;
3960: snes->conv_hist_len = 0;
3961: snes->conv_hist_reset = reset;
3962: return(0);
3963: }
3965: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3966: #include <engine.h> /* MATLAB include file */
3967: #include <mex.h> /* MATLAB include file */
3969: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
3970: {
3971: mxArray *mat;
3972: PetscInt i;
3973: PetscReal *ar;
3976: mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
3977: ar = (PetscReal*) mxGetData(mat);
3978: for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
3979: PetscFunctionReturn(mat);
3980: }
3981: #endif
3983: /*@C
3984: SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
3986: Not Collective
3988: Input Parameter:
3989: . snes - iterative context obtained from SNESCreate()
3991: Output Parameters:
3992: . a - array to hold history
3993: . its - integer array holds the number of linear iterations (or
3994: negative if not converged) for each solve.
3995: - na - size of a and its
3997: Notes:
3998: The calling sequence for this routine in Fortran is
3999: $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
4001: This routine is useful, e.g., when running a code for purposes
4002: of accurate performance monitoring, when no I/O should be done
4003: during the section of code that is being timed.
4005: Level: intermediate
4007: .keywords: SNES, get, convergence, history
4009: .seealso: SNESSetConvergencHistory()
4011: @*/
4012: PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
4013: {
4016: if (a) *a = snes->conv_hist;
4017: if (its) *its = snes->conv_hist_its;
4018: if (na) *na = snes->conv_hist_len;
4019: return(0);
4020: }
4022: /*@C
4023: SNESSetUpdate - Sets the general-purpose update function called
4024: at the beginning of every iteration of the nonlinear solve. Specifically
4025: it is called just before the Jacobian is "evaluated".
4027: Logically Collective on SNES
4029: Input Parameters:
4030: . snes - The nonlinear solver context
4031: . func - The function
4033: Calling sequence of func:
4034: . func (SNES snes, PetscInt step);
4036: . step - The current step of the iteration
4038: Level: advanced
4040: 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()
4041: This is not used by most users.
4043: .keywords: SNES, update
4045: .seealso SNESSetJacobian(), SNESSolve()
4046: @*/
4047: PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
4048: {
4051: snes->ops->update = func;
4052: return(0);
4053: }
4055: /*
4056: SNESScaleStep_Private - Scales a step so that its length is less than the
4057: positive parameter delta.
4059: Input Parameters:
4060: + snes - the SNES context
4061: . y - approximate solution of linear system
4062: . fnorm - 2-norm of current function
4063: - delta - trust region size
4065: Output Parameters:
4066: + gpnorm - predicted function norm at the new point, assuming local
4067: linearization. The value is zero if the step lies within the trust
4068: region, and exceeds zero otherwise.
4069: - ynorm - 2-norm of the step
4071: Note:
4072: For non-trust region methods such as SNESNEWTONLS, the parameter delta
4073: is set to be the maximum allowable step size.
4075: .keywords: SNES, nonlinear, scale, step
4076: */
4077: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
4078: {
4079: PetscReal nrm;
4080: PetscScalar cnorm;
4088: VecNorm(y,NORM_2,&nrm);
4089: if (nrm > *delta) {
4090: nrm = *delta/nrm;
4091: *gpnorm = (1.0 - nrm)*(*fnorm);
4092: cnorm = nrm;
4093: VecScale(y,cnorm);
4094: *ynorm = *delta;
4095: } else {
4096: *gpnorm = 0.0;
4097: *ynorm = nrm;
4098: }
4099: return(0);
4100: }
4102: /*@
4103: SNESReasonView - Displays the reason a SNES solve converged or diverged to a viewer
4105: Collective on SNES
4107: Parameter:
4108: + snes - iterative context obtained from SNESCreate()
4109: - viewer - the viewer to display the reason
4112: Options Database Keys:
4113: . -snes_converged_reason - print reason for converged or diverged, also prints number of iterations
4115: Level: beginner
4117: .keywords: SNES, solve, linear system
4119: .seealso: SNESCreate(), SNESSetUp(), SNESDestroy(), SNESSetTolerances(), SNESConvergedDefault()
4121: @*/
4122: PetscErrorCode SNESReasonView(SNES snes,PetscViewer viewer)
4123: {
4124: PetscViewerFormat format;
4125: PetscBool isAscii;
4126: PetscErrorCode ierr;
4129: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
4130: if (isAscii) {
4131: PetscViewerGetFormat(viewer, &format);
4132: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
4133: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
4134: DM dm;
4135: Vec u;
4136: PetscDS prob;
4137: PetscInt Nf, f;
4138: PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
4139: PetscReal error;
4141: SNESGetDM(snes, &dm);
4142: SNESGetSolution(snes, &u);
4143: DMGetDS(dm, &prob);
4144: PetscDSGetNumFields(prob, &Nf);
4145: PetscMalloc1(Nf, &exactFuncs);
4146: for (f = 0; f < Nf; ++f) {PetscDSGetExactSolution(prob, f, &exactFuncs[f]);}
4147: DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);
4148: PetscFree(exactFuncs);
4149: if (error < 1.0e-11) {PetscViewerASCIIPrintf(viewer, "L_2 Error: < 1.0e-11\n");}
4150: else {PetscViewerASCIIPrintf(viewer, "L_2 Error: %g\n", error);}
4151: }
4152: if (snes->reason > 0) {
4153: if (((PetscObject) snes)->prefix) {
4154: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4155: } else {
4156: PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4157: }
4158: } else {
4159: if (((PetscObject) snes)->prefix) {
4160: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
4161: } else {
4162: PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
4163: }
4164: }
4165: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
4166: }
4167: return(0);
4168: }
4170: /*@C
4171: SNESReasonViewFromOptions - Processes command line options to determine if/how a SNESReason is to be viewed.
4173: Collective on SNES
4175: Input Parameters:
4176: . snes - the SNES object
4178: Level: intermediate
4180: @*/
4181: PetscErrorCode SNESReasonViewFromOptions(SNES snes)
4182: {
4183: PetscErrorCode ierr;
4184: PetscViewer viewer;
4185: PetscBool flg;
4186: static PetscBool incall = PETSC_FALSE;
4187: PetscViewerFormat format;
4190: if (incall) return(0);
4191: incall = PETSC_TRUE;
4192: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
4193: if (flg) {
4194: PetscViewerPushFormat(viewer,format);
4195: SNESReasonView(snes,viewer);
4196: PetscViewerPopFormat(viewer);
4197: PetscViewerDestroy(&viewer);
4198: }
4199: incall = PETSC_FALSE;
4200: return(0);
4201: }
4203: /*@C
4204: SNESSolve - Solves a nonlinear system F(x) = b.
4205: Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
4207: Collective on SNES
4209: Input Parameters:
4210: + snes - the SNES context
4211: . b - the constant part of the equation F(x) = b, or NULL to use zero.
4212: - x - the solution vector.
4214: Notes:
4215: The user should initialize the vector,x, with the initial guess
4216: for the nonlinear solve prior to calling SNESSolve. In particular,
4217: to employ an initial guess of zero, the user should explicitly set
4218: this vector to zero by calling VecSet().
4220: Level: beginner
4222: .keywords: SNES, nonlinear, solve
4224: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
4225: @*/
4226: PetscErrorCode SNESSolve(SNES snes,Vec b,Vec x)
4227: {
4228: PetscErrorCode ierr;
4229: PetscBool flg;
4230: PetscInt grid;
4231: Vec xcreated = NULL;
4232: DM dm;
4241: /* High level operations using the nonlinear solver */
4242: {
4243: PetscViewer viewer;
4244: PetscViewerFormat format;
4245: PetscInt num;
4246: PetscBool flg;
4247: static PetscBool incall = PETSC_FALSE;
4249: if (!incall) {
4250: /* Estimate the convergence rate of the discretization */
4251: PetscOptionsGetViewer(PetscObjectComm((PetscObject) snes), ((PetscObject) snes)->prefix, "-snes_convergence_estimate", &viewer, &format, &flg);
4252: if (flg) {
4253: PetscConvEst conv;
4254: PetscReal alpha; /* Convergence rate of the solution error in the L_2 norm */
4256: incall = PETSC_TRUE;
4257: PetscConvEstCreate(PetscObjectComm((PetscObject) snes), &conv);
4258: PetscConvEstSetSolver(conv, snes);
4259: PetscConvEstSetFromOptions(conv);
4260: PetscConvEstSetUp(conv);
4261: PetscConvEstGetConvRate(conv, &alpha);
4262: PetscViewerPushFormat(viewer, format);
4263: PetscConvEstRateView(conv, alpha, viewer);
4264: PetscViewerPopFormat(viewer);
4265: PetscViewerDestroy(&viewer);
4266: PetscConvEstDestroy(&conv);
4267: incall = PETSC_FALSE;
4268: }
4269: /* Adaptively refine the initial grid */
4270: num = 1;
4271: PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_initial", &num, &flg);
4272: if (flg) {
4273: DMAdaptor adaptor;
4275: incall = PETSC_TRUE;
4276: DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);
4277: DMAdaptorSetSolver(adaptor, snes);
4278: DMAdaptorSetSequenceLength(adaptor, num);
4279: DMAdaptorSetFromOptions(adaptor);
4280: DMAdaptorSetUp(adaptor);
4281: DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_INITIAL, &dm, &x);
4282: DMAdaptorDestroy(&adaptor);
4283: incall = PETSC_FALSE;
4284: }
4285: /* Use grid sequencing to adapt */
4286: num = 0;
4287: PetscOptionsGetInt(NULL, ((PetscObject) snes)->prefix, "-snes_adapt_sequence", &num, NULL);
4288: if (num) {
4289: DMAdaptor adaptor;
4291: incall = PETSC_TRUE;
4292: DMAdaptorCreate(PETSC_COMM_WORLD, &adaptor);
4293: DMAdaptorSetSolver(adaptor, snes);
4294: DMAdaptorSetSequenceLength(adaptor, num);
4295: DMAdaptorSetFromOptions(adaptor);
4296: DMAdaptorSetUp(adaptor);
4297: DMAdaptorAdapt(adaptor, x, DM_ADAPTATION_SEQUENTIAL, &dm, &x);
4298: DMAdaptorDestroy(&adaptor);
4299: incall = PETSC_FALSE;
4300: }
4301: }
4302: }
4303: if (!x) {
4304: SNESGetDM(snes,&dm);
4305: DMCreateGlobalVector(dm,&xcreated);
4306: x = xcreated;
4307: }
4308: SNESViewFromOptions(snes,NULL,"-snes_view_pre");
4310: for (grid=0; grid<snes->gridsequence; grid++) {PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));}
4311: for (grid=0; grid<snes->gridsequence+1; grid++) {
4313: /* set solution vector */
4314: if (!grid) {PetscObjectReference((PetscObject)x);}
4315: VecDestroy(&snes->vec_sol);
4316: snes->vec_sol = x;
4317: SNESGetDM(snes,&dm);
4319: /* set affine vector if provided */
4320: if (b) { PetscObjectReference((PetscObject)b); }
4321: VecDestroy(&snes->vec_rhs);
4322: snes->vec_rhs = b;
4324: if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
4325: if (snes->vec_rhs == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
4326: if (!snes->vec_sol_update /* && snes->vec_sol */) {
4327: VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
4328: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
4329: }
4330: DMShellSetGlobalVector(dm,snes->vec_sol);
4331: SNESSetUp(snes);
4333: if (!grid) {
4334: if (snes->ops->computeinitialguess) {
4335: (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
4336: }
4337: }
4339: if (snes->conv_hist_reset) snes->conv_hist_len = 0;
4340: if (snes->counters_reset) {snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;}
4342: PetscLogEventBegin(SNES_Solve,snes,0,0,0);
4343: (*snes->ops->solve)(snes);
4344: PetscLogEventEnd(SNES_Solve,snes,0,0,0);
4345: if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
4346: snes->domainerror = PETSC_FALSE; /* clear the flag if it has been set */
4348: if (snes->lagjac_persist) snes->jac_iter += snes->iter;
4349: if (snes->lagpre_persist) snes->pre_iter += snes->iter;
4351: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_test_local_min",NULL,NULL,&flg);
4352: if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
4353: SNESReasonViewFromOptions(snes);
4355: if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
4356: if (snes->reason < 0) break;
4357: if (grid < snes->gridsequence) {
4358: DM fine;
4359: Vec xnew;
4360: Mat interp;
4362: DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
4363: if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
4364: DMCreateInterpolation(snes->dm,fine,&interp,NULL);
4365: DMCreateGlobalVector(fine,&xnew);
4366: MatInterpolate(interp,x,xnew);
4367: DMInterpolate(snes->dm,interp,fine);
4368: MatDestroy(&interp);
4369: x = xnew;
4371: SNESReset(snes);
4372: SNESSetDM(snes,fine);
4373: DMDestroy(&fine);
4374: PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
4375: }
4376: }
4377: SNESViewFromOptions(snes,NULL,"-snes_view");
4378: VecViewFromOptions(snes->vec_sol,(PetscObject)snes,"-snes_view_solution");
4380: VecDestroy(&xcreated);
4381: PetscObjectSAWsBlock((PetscObject)snes);
4382: return(0);
4383: }
4385: /* --------- Internal routines for SNES Package --------- */
4387: /*@C
4388: SNESSetType - Sets the method for the nonlinear solver.
4390: Collective on SNES
4392: Input Parameters:
4393: + snes - the SNES context
4394: - type - a known method
4396: Options Database Key:
4397: . -snes_type <type> - Sets the method; use -help for a list
4398: of available methods (for instance, newtonls or newtontr)
4400: Notes:
4401: See "petsc/include/petscsnes.h" for available methods (for instance)
4402: + SNESNEWTONLS - Newton's method with line search
4403: (systems of nonlinear equations)
4404: . SNESNEWTONTR - Newton's method with trust region
4405: (systems of nonlinear equations)
4407: Normally, it is best to use the SNESSetFromOptions() command and then
4408: set the SNES solver type from the options database rather than by using
4409: this routine. Using the options database provides the user with
4410: maximum flexibility in evaluating the many nonlinear solvers.
4411: The SNESSetType() routine is provided for those situations where it
4412: is necessary to set the nonlinear solver independently of the command
4413: line or options database. This might be the case, for example, when
4414: the choice of solver changes during the execution of the program,
4415: and the user's application is taking responsibility for choosing the
4416: appropriate method.
4418: Developer Notes: SNESRegister() adds a constructor for a new SNESType to SNESList, SNESSetType() locates
4419: the constructor in that list and calls it to create the spexific object.
4421: Level: intermediate
4423: .keywords: SNES, set, type
4425: .seealso: SNESType, SNESCreate(), SNESDestroy(), SNESGetType(), SNESSetFromOptions()
4427: @*/
4428: PetscErrorCode SNESSetType(SNES snes,SNESType type)
4429: {
4430: PetscErrorCode ierr,(*r)(SNES);
4431: PetscBool match;
4437: PetscObjectTypeCompare((PetscObject)snes,type,&match);
4438: if (match) return(0);
4440: PetscFunctionListFind(SNESList,type,&r);
4441: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
4442: /* Destroy the previous private SNES context */
4443: if (snes->ops->destroy) {
4444: (*(snes)->ops->destroy)(snes);
4445: snes->ops->destroy = NULL;
4446: }
4447: /* Reinitialize function pointers in SNESOps structure */
4448: snes->ops->setup = 0;
4449: snes->ops->solve = 0;
4450: snes->ops->view = 0;
4451: snes->ops->setfromoptions = 0;
4452: snes->ops->destroy = 0;
4453: /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
4454: snes->setupcalled = PETSC_FALSE;
4456: PetscObjectChangeTypeName((PetscObject)snes,type);
4457: (*r)(snes);
4458: return(0);
4459: }
4461: /*@C
4462: SNESGetType - Gets the SNES method type and name (as a string).
4464: Not Collective
4466: Input Parameter:
4467: . snes - nonlinear solver context
4469: Output Parameter:
4470: . type - SNES method (a character string)
4472: Level: intermediate
4474: .keywords: SNES, nonlinear, get, type, name
4475: @*/
4476: PetscErrorCode SNESGetType(SNES snes,SNESType *type)
4477: {
4481: *type = ((PetscObject)snes)->type_name;
4482: return(0);
4483: }
4485: /*@
4486: SNESSetSolution - Sets the solution vector for use by the SNES routines.
4488: Logically Collective on SNES and Vec
4490: Input Parameters:
4491: + snes - the SNES context obtained from SNESCreate()
4492: - u - the solution vector
4494: Level: beginner
4496: .keywords: SNES, set, solution
4497: @*/
4498: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4499: {
4500: DM dm;
4506: PetscObjectReference((PetscObject) u);
4507: VecDestroy(&snes->vec_sol);
4509: snes->vec_sol = u;
4511: SNESGetDM(snes, &dm);
4512: DMShellSetGlobalVector(dm, u);
4513: return(0);
4514: }
4516: /*@
4517: SNESGetSolution - Returns the vector where the approximate solution is
4518: stored. This is the fine grid solution when using SNESSetGridSequence().
4520: Not Collective, but Vec is parallel if SNES is parallel
4522: Input Parameter:
4523: . snes - the SNES context
4525: Output Parameter:
4526: . x - the solution
4528: Level: intermediate
4530: .keywords: SNES, nonlinear, get, solution
4532: .seealso: SNESGetSolutionUpdate(), SNESGetFunction()
4533: @*/
4534: PetscErrorCode SNESGetSolution(SNES snes,Vec *x)
4535: {
4539: *x = snes->vec_sol;
4540: return(0);
4541: }
4543: /*@
4544: SNESGetSolutionUpdate - Returns the vector where the solution update is
4545: stored.
4547: Not Collective, but Vec is parallel if SNES is parallel
4549: Input Parameter:
4550: . snes - the SNES context
4552: Output Parameter:
4553: . x - the solution update
4555: Level: advanced
4557: .keywords: SNES, nonlinear, get, solution, update
4559: .seealso: SNESGetSolution(), SNESGetFunction()
4560: @*/
4561: PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x)
4562: {
4566: *x = snes->vec_sol_update;
4567: return(0);
4568: }
4570: /*@C
4571: SNESGetFunction - Returns the vector where the function is stored.
4573: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
4575: Input Parameter:
4576: . snes - the SNES context
4578: Output Parameter:
4579: + r - the vector that is used to store residuals (or NULL if you don't want it)
4580: . f - the function (or NULL if you don't want it); see SNESFunction for calling sequence details
4581: - ctx - the function context (or NULL if you don't want it)
4583: Level: advanced
4585: Notes: The vector r DOES NOT, in general contain the current value of the SNES nonlinear function
4587: .keywords: SNES, nonlinear, get, function
4589: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4590: @*/
4591: PetscErrorCode SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4592: {
4594: DM dm;
4598: if (r) {
4599: if (!snes->vec_func) {
4600: if (snes->vec_rhs) {
4601: VecDuplicate(snes->vec_rhs,&snes->vec_func);
4602: } else if (snes->vec_sol) {
4603: VecDuplicate(snes->vec_sol,&snes->vec_func);
4604: } else if (snes->dm) {
4605: DMCreateGlobalVector(snes->dm,&snes->vec_func);
4606: }
4607: }
4608: *r = snes->vec_func;
4609: }
4610: SNESGetDM(snes,&dm);
4611: DMSNESGetFunction(dm,f,ctx);
4612: return(0);
4613: }
4615: /*@C
4616: SNESGetNGS - Returns the NGS function and context.
4618: Input Parameter:
4619: . snes - the SNES context
4621: Output Parameter:
4622: + f - the function (or NULL) see SNESNGSFunction for details
4623: - ctx - the function context (or NULL)
4625: Level: advanced
4627: .keywords: SNES, nonlinear, get, function
4629: .seealso: SNESSetNGS(), SNESGetFunction()
4630: @*/
4632: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4633: {
4635: DM dm;
4639: SNESGetDM(snes,&dm);
4640: DMSNESGetNGS(dm,f,ctx);
4641: return(0);
4642: }
4644: /*@C
4645: SNESSetOptionsPrefix - Sets the prefix used for searching for all
4646: SNES options in the database.
4648: Logically Collective on SNES
4650: Input Parameter:
4651: + snes - the SNES context
4652: - prefix - the prefix to prepend to all option names
4654: Notes:
4655: A hyphen (-) must NOT be given at the beginning of the prefix name.
4656: The first character of all runtime options is AUTOMATICALLY the hyphen.
4658: Level: advanced
4660: .keywords: SNES, set, options, prefix, database
4662: .seealso: SNESSetFromOptions()
4663: @*/
4664: PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[])
4665: {
4670: PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4671: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4672: if (snes->linesearch) {
4673: SNESGetLineSearch(snes,&snes->linesearch);
4674: PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4675: }
4676: KSPSetOptionsPrefix(snes->ksp,prefix);
4677: return(0);
4678: }
4680: /*@C
4681: SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4682: SNES options in the database.
4684: Logically Collective on SNES
4686: Input Parameters:
4687: + snes - the SNES context
4688: - prefix - the prefix to prepend to all option names
4690: Notes:
4691: A hyphen (-) must NOT be given at the beginning of the prefix name.
4692: The first character of all runtime options is AUTOMATICALLY the hyphen.
4694: Level: advanced
4696: .keywords: SNES, append, options, prefix, database
4698: .seealso: SNESGetOptionsPrefix()
4699: @*/
4700: PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4701: {
4706: PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4707: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4708: if (snes->linesearch) {
4709: SNESGetLineSearch(snes,&snes->linesearch);
4710: PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4711: }
4712: KSPAppendOptionsPrefix(snes->ksp,prefix);
4713: return(0);
4714: }
4716: /*@C
4717: SNESGetOptionsPrefix - Sets the prefix used for searching for all
4718: SNES options in the database.
4720: Not Collective
4722: Input Parameter:
4723: . snes - the SNES context
4725: Output Parameter:
4726: . prefix - pointer to the prefix string used
4728: Notes: On the fortran side, the user should pass in a string 'prefix' of
4729: sufficient length to hold the prefix.
4731: Level: advanced
4733: .keywords: SNES, get, options, prefix, database
4735: .seealso: SNESAppendOptionsPrefix()
4736: @*/
4737: PetscErrorCode SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4738: {
4743: PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4744: return(0);
4745: }
4748: /*@C
4749: SNESRegister - Adds a method to the nonlinear solver package.
4751: Not collective
4753: Input Parameters:
4754: + name_solver - name of a new user-defined solver
4755: - routine_create - routine to create method context
4757: Notes:
4758: SNESRegister() may be called multiple times to add several user-defined solvers.
4760: Sample usage:
4761: .vb
4762: SNESRegister("my_solver",MySolverCreate);
4763: .ve
4765: Then, your solver can be chosen with the procedural interface via
4766: $ SNESSetType(snes,"my_solver")
4767: or at runtime via the option
4768: $ -snes_type my_solver
4770: Level: advanced
4772: Note: If your function is not being put into a shared library then use SNESRegister() instead
4774: .keywords: SNES, nonlinear, register
4776: .seealso: SNESRegisterAll(), SNESRegisterDestroy()
4778: Level: advanced
4779: @*/
4780: PetscErrorCode SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4781: {
4785: PetscFunctionListAdd(&SNESList,sname,function);
4786: return(0);
4787: }
4789: PetscErrorCode SNESTestLocalMin(SNES snes)
4790: {
4792: PetscInt N,i,j;
4793: Vec u,uh,fh;
4794: PetscScalar value;
4795: PetscReal norm;
4798: SNESGetSolution(snes,&u);
4799: VecDuplicate(u,&uh);
4800: VecDuplicate(u,&fh);
4802: /* currently only works for sequential */
4803: PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4804: VecGetSize(u,&N);
4805: for (i=0; i<N; i++) {
4806: VecCopy(u,uh);
4807: PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4808: for (j=-10; j<11; j++) {
4809: value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4810: VecSetValue(uh,i,value,ADD_VALUES);
4811: SNESComputeFunction(snes,uh,fh);
4812: VecNorm(fh,NORM_2,&norm);
4813: PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);
4814: value = -value;
4815: VecSetValue(uh,i,value,ADD_VALUES);
4816: }
4817: }
4818: VecDestroy(&uh);
4819: VecDestroy(&fh);
4820: return(0);
4821: }
4823: /*@
4824: SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4825: computing relative tolerance for linear solvers within an inexact
4826: Newton method.
4828: Logically Collective on SNES
4830: Input Parameters:
4831: + snes - SNES context
4832: - flag - PETSC_TRUE or PETSC_FALSE
4834: Options Database:
4835: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4836: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
4837: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4838: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4839: . -snes_ksp_ew_gamma <gamma> - Sets gamma
4840: . -snes_ksp_ew_alpha <alpha> - Sets alpha
4841: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4842: - -snes_ksp_ew_threshold <threshold> - Sets threshold
4844: Notes:
4845: Currently, the default is to use a constant relative tolerance for
4846: the inner linear solvers. Alternatively, one can use the
4847: Eisenstat-Walker method, where the relative convergence tolerance
4848: is reset at each Newton iteration according progress of the nonlinear
4849: solver.
4851: Level: advanced
4853: Reference:
4854: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4855: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
4857: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
4859: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4860: @*/
4861: PetscErrorCode SNESKSPSetUseEW(SNES snes,PetscBool flag)
4862: {
4866: snes->ksp_ewconv = flag;
4867: return(0);
4868: }
4870: /*@
4871: SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
4872: for computing relative tolerance for linear solvers within an
4873: inexact Newton method.
4875: Not Collective
4877: Input Parameter:
4878: . snes - SNES context
4880: Output Parameter:
4881: . flag - PETSC_TRUE or PETSC_FALSE
4883: Notes:
4884: Currently, the default is to use a constant relative tolerance for
4885: the inner linear solvers. Alternatively, one can use the
4886: Eisenstat-Walker method, where the relative convergence tolerance
4887: is reset at each Newton iteration according progress of the nonlinear
4888: solver.
4890: Level: advanced
4892: Reference:
4893: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4894: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
4896: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
4898: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4899: @*/
4900: PetscErrorCode SNESKSPGetUseEW(SNES snes, PetscBool *flag)
4901: {
4905: *flag = snes->ksp_ewconv;
4906: return(0);
4907: }
4909: /*@
4910: SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
4911: convergence criteria for the linear solvers within an inexact
4912: Newton method.
4914: Logically Collective on SNES
4916: Input Parameters:
4917: + snes - SNES context
4918: . version - version 1, 2 (default is 2) or 3
4919: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4920: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4921: . gamma - multiplicative factor for version 2 rtol computation
4922: (0 <= gamma2 <= 1)
4923: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
4924: . alpha2 - power for safeguard
4925: - threshold - threshold for imposing safeguard (0 < threshold < 1)
4927: Note:
4928: Version 3 was contributed by Luis Chacon, June 2006.
4930: Use PETSC_DEFAULT to retain the default for any of the parameters.
4932: Level: advanced
4934: Reference:
4935: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4936: inexact Newton method", Utah State University Math. Stat. Dept. Res.
4937: Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
4939: .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
4941: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
4942: @*/
4943: PetscErrorCode SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
4944: {
4945: SNESKSPEW *kctx;
4949: kctx = (SNESKSPEW*)snes->kspconvctx;
4950: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4959: if (version != PETSC_DEFAULT) kctx->version = version;
4960: if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0;
4961: if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max;
4962: if (gamma != PETSC_DEFAULT) kctx->gamma = gamma;
4963: if (alpha != PETSC_DEFAULT) kctx->alpha = alpha;
4964: if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2;
4965: if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
4967: 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);
4968: 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);
4969: 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);
4970: 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);
4971: 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);
4972: 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);
4973: return(0);
4974: }
4976: /*@
4977: SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
4978: convergence criteria for the linear solvers within an inexact
4979: Newton method.
4981: Not Collective
4983: Input Parameters:
4984: snes - SNES context
4986: Output Parameters:
4987: + version - version 1, 2 (default is 2) or 3
4988: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4989: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4990: . gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
4991: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
4992: . alpha2 - power for safeguard
4993: - threshold - threshold for imposing safeguard (0 < threshold < 1)
4995: Level: advanced
4997: .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
4999: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
5000: @*/
5001: PetscErrorCode SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
5002: {
5003: SNESKSPEW *kctx;
5007: kctx = (SNESKSPEW*)snes->kspconvctx;
5008: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
5009: if (version) *version = kctx->version;
5010: if (rtol_0) *rtol_0 = kctx->rtol_0;
5011: if (rtol_max) *rtol_max = kctx->rtol_max;
5012: if (gamma) *gamma = kctx->gamma;
5013: if (alpha) *alpha = kctx->alpha;
5014: if (alpha2) *alpha2 = kctx->alpha2;
5015: if (threshold) *threshold = kctx->threshold;
5016: return(0);
5017: }
5019: PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5020: {
5022: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
5023: PetscReal rtol = PETSC_DEFAULT,stol;
5026: if (!snes->ksp_ewconv) return(0);
5027: if (!snes->iter) {
5028: rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
5029: VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
5030: }
5031: else {
5032: if (kctx->version == 1) {
5033: rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
5034: if (rtol < 0.0) rtol = -rtol;
5035: stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
5036: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5037: } else if (kctx->version == 2) {
5038: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5039: stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
5040: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
5041: } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
5042: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
5043: /* safeguard: avoid sharp decrease of rtol */
5044: stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
5045: stol = PetscMax(rtol,stol);
5046: rtol = PetscMin(kctx->rtol_0,stol);
5047: /* safeguard: avoid oversolving */
5048: stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
5049: stol = PetscMax(rtol,stol);
5050: rtol = PetscMin(kctx->rtol_0,stol);
5051: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
5052: }
5053: /* safeguard: avoid rtol greater than one */
5054: rtol = PetscMin(rtol,kctx->rtol_max);
5055: KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
5056: PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
5057: return(0);
5058: }
5060: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
5061: {
5063: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
5064: PCSide pcside;
5065: Vec lres;
5068: if (!snes->ksp_ewconv) return(0);
5069: KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
5070: kctx->norm_last = snes->norm;
5071: if (kctx->version == 1) {
5072: PC pc;
5073: PetscBool isNone;
5075: KSPGetPC(ksp, &pc);
5076: PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
5077: KSPGetPCSide(ksp,&pcside);
5078: if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
5079: /* KSP residual is true linear residual */
5080: KSPGetResidualNorm(ksp,&kctx->lresid_last);
5081: } else {
5082: /* KSP residual is preconditioned residual */
5083: /* compute true linear residual norm */
5084: VecDuplicate(b,&lres);
5085: MatMult(snes->jacobian,x,lres);
5086: VecAYPX(lres,-1.0,b);
5087: VecNorm(lres,NORM_2,&kctx->lresid_last);
5088: VecDestroy(&lres);
5089: }
5090: }
5091: return(0);
5092: }
5094: /*@
5095: SNESGetKSP - Returns the KSP context for a SNES solver.
5097: Not Collective, but if SNES object is parallel, then KSP object is parallel
5099: Input Parameter:
5100: . snes - the SNES context
5102: Output Parameter:
5103: . ksp - the KSP context
5105: Notes:
5106: The user can then directly manipulate the KSP context to set various
5107: options, etc. Likewise, the user can then extract and manipulate the
5108: PC contexts as well.
5110: Level: beginner
5112: .keywords: SNES, nonlinear, get, KSP, context
5114: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
5115: @*/
5116: PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp)
5117: {
5124: if (!snes->ksp) {
5125: PetscBool monitor = PETSC_FALSE;
5127: KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
5128: PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
5129: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);
5131: KSPSetPreSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPreSolve_SNESEW,snes);
5132: KSPSetPostSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPostSolve_SNESEW,snes);
5134: PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
5135: if (monitor) {
5136: KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
5137: }
5138: monitor = PETSC_FALSE;
5139: PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
5140: if (monitor) {
5141: PetscObject *objs;
5142: KSPMonitorSNESLGResidualNormCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
5143: objs[0] = (PetscObject) snes;
5144: KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
5145: }
5146: }
5147: *ksp = snes->ksp;
5148: return(0);
5149: }
5152: #include <petsc/private/dmimpl.h>
5153: /*@
5154: SNESSetDM - Sets the DM that may be used by some nonlinear solvers or their underlying preconditioners
5156: Logically Collective on SNES
5158: Input Parameters:
5159: + snes - the nonlinear solver context
5160: - dm - the dm, cannot be NULL
5162: Level: intermediate
5164: .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
5165: @*/
5166: PetscErrorCode SNESSetDM(SNES snes,DM dm)
5167: {
5169: KSP ksp;
5170: DMSNES sdm;
5175: PetscObjectReference((PetscObject)dm);
5176: if (snes->dm) { /* Move the DMSNES context over to the new DM unless the new DM already has one */
5177: if (snes->dm->dmsnes && !dm->dmsnes) {
5178: DMCopyDMSNES(snes->dm,dm);
5179: DMGetDMSNES(snes->dm,&sdm);
5180: if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
5181: }
5182: DMCoarsenHookRemove(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
5183: DMDestroy(&snes->dm);
5184: }
5185: snes->dm = dm;
5186: snes->dmAuto = PETSC_FALSE;
5188: SNESGetKSP(snes,&ksp);
5189: KSPSetDM(ksp,dm);
5190: KSPSetDMActive(ksp,PETSC_FALSE);
5191: if (snes->npc) {
5192: SNESSetDM(snes->npc, snes->dm);
5193: SNESSetNPCSide(snes,snes->npcside);
5194: }
5195: return(0);
5196: }
5198: /*@
5199: SNESGetDM - Gets the DM that may be used by some preconditioners
5201: Not Collective but DM obtained is parallel on SNES
5203: Input Parameter:
5204: . snes - the preconditioner context
5206: Output Parameter:
5207: . dm - the dm
5209: Level: intermediate
5211: .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
5212: @*/
5213: PetscErrorCode SNESGetDM(SNES snes,DM *dm)
5214: {
5219: if (!snes->dm) {
5220: DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
5221: snes->dmAuto = PETSC_TRUE;
5222: }
5223: *dm = snes->dm;
5224: return(0);
5225: }
5227: /*@
5228: SNESSetNPC - Sets the nonlinear preconditioner to be used.
5230: Collective on SNES
5232: Input Parameters:
5233: + snes - iterative context obtained from SNESCreate()
5234: - pc - the preconditioner object
5236: Notes:
5237: Use SNESGetNPC() to retrieve the preconditioner context (for example,
5238: to configure it using the API).
5240: Level: developer
5242: .keywords: SNES, set, precondition
5243: .seealso: SNESGetNPC(), SNESHasNPC()
5244: @*/
5245: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
5246: {
5253: PetscObjectReference((PetscObject) pc);
5254: SNESDestroy(&snes->npc);
5255: snes->npc = pc;
5256: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->npc);
5257: return(0);
5258: }
5260: /*@
5261: SNESGetNPC - Creates a nonlinear preconditioning solver (SNES) to be used to precondition the nonlinear solver.
5263: Not Collective
5265: Input Parameter:
5266: . snes - iterative context obtained from SNESCreate()
5268: Output Parameter:
5269: . pc - preconditioner context
5271: Notes: If a SNES was previously set with SNESSetNPC() then that SNES is returned.
5273: Level: developer
5275: .keywords: SNES, get, preconditioner
5276: .seealso: SNESSetNPC(), SNESHasNPC()
5277: @*/
5278: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
5279: {
5281: const char *optionsprefix;
5286: if (!snes->npc) {
5287: SNESCreate(PetscObjectComm((PetscObject)snes),&snes->npc);
5288: PetscObjectIncrementTabLevel((PetscObject)snes->npc,(PetscObject)snes,1);
5289: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->npc);
5290: SNESGetOptionsPrefix(snes,&optionsprefix);
5291: SNESSetOptionsPrefix(snes->npc,optionsprefix);
5292: SNESAppendOptionsPrefix(snes->npc,"npc_");
5293: SNESSetCountersReset(snes->npc,PETSC_FALSE);
5294: }
5295: *pc = snes->npc;
5296: return(0);
5297: }
5299: /*@
5300: SNESHasNPC - Returns whether a nonlinear preconditioner exists
5302: Not Collective
5304: Input Parameter:
5305: . snes - iterative context obtained from SNESCreate()
5307: Output Parameter:
5308: . has_npc - whether the SNES has an NPC or not
5310: Level: developer
5312: .keywords: SNES, has, preconditioner
5313: .seealso: SNESSetNPC(), SNESGetNPC()
5314: @*/
5315: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
5316: {
5319: *has_npc = (PetscBool) (snes->npc ? PETSC_TRUE : PETSC_FALSE);
5320: return(0);
5321: }
5323: /*@
5324: SNESSetNPCSide - Sets the preconditioning side.
5326: Logically Collective on SNES
5328: Input Parameter:
5329: . snes - iterative context obtained from SNESCreate()
5331: Output Parameter:
5332: . side - the preconditioning side, where side is one of
5333: .vb
5334: PC_LEFT - left preconditioning
5335: PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5336: .ve
5338: Options Database Keys:
5339: . -snes_pc_side <right,left>
5341: Notes: SNESNRICHARDSON and SNESNCG only support left preconditioning.
5343: Level: intermediate
5345: .keywords: SNES, set, right, left, side, preconditioner, flag
5347: .seealso: SNESGetNPCSide(), KSPSetPCSide()
5348: @*/
5349: PetscErrorCode SNESSetNPCSide(SNES snes,PCSide side)
5350: {
5354: snes->npcside= side;
5355: return(0);
5356: }
5358: /*@
5359: SNESGetNPCSide - Gets the preconditioning side.
5361: Not Collective
5363: Input Parameter:
5364: . snes - iterative context obtained from SNESCreate()
5366: Output Parameter:
5367: . side - the preconditioning side, where side is one of
5368: .vb
5369: PC_LEFT - left preconditioning
5370: PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5371: .ve
5373: Level: intermediate
5375: .keywords: SNES, get, right, left, side, preconditioner, flag
5377: .seealso: SNESSetNPCSide(), KSPGetPCSide()
5378: @*/
5379: PetscErrorCode SNESGetNPCSide(SNES snes,PCSide *side)
5380: {
5384: *side = snes->npcside;
5385: return(0);
5386: }
5388: /*@
5389: SNESSetLineSearch - Sets the linesearch on the SNES instance.
5391: Collective on SNES
5393: Input Parameters:
5394: + snes - iterative context obtained from SNESCreate()
5395: - linesearch - the linesearch object
5397: Notes:
5398: Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
5399: to configure it using the API).
5401: Level: developer
5403: .keywords: SNES, set, linesearch
5404: .seealso: SNESGetLineSearch()
5405: @*/
5406: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
5407: {
5414: PetscObjectReference((PetscObject) linesearch);
5415: SNESLineSearchDestroy(&snes->linesearch);
5417: snes->linesearch = linesearch;
5419: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5420: return(0);
5421: }
5423: /*@
5424: SNESGetLineSearch - Returns a pointer to the line search context set with SNESSetLineSearch()
5425: or creates a default line search instance associated with the SNES and returns it.
5427: Not Collective
5429: Input Parameter:
5430: . snes - iterative context obtained from SNESCreate()
5432: Output Parameter:
5433: . linesearch - linesearch context
5435: Level: beginner
5437: .keywords: SNES, get, linesearch
5438: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5439: @*/
5440: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5441: {
5443: const char *optionsprefix;
5448: if (!snes->linesearch) {
5449: SNESGetOptionsPrefix(snes, &optionsprefix);
5450: SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5451: SNESLineSearchSetSNES(snes->linesearch, snes);
5452: SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5453: PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5454: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5455: }
5456: *linesearch = snes->linesearch;
5457: return(0);
5458: }
5460: #if defined(PETSC_HAVE_MATLAB_ENGINE)
5461: #include <mex.h>
5463: typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
5465: /*
5466: SNESComputeFunction_Matlab - Calls the function that has been set with SNESSetFunctionMatlab().
5468: Collective on SNES
5470: Input Parameters:
5471: + snes - the SNES context
5472: - x - input vector
5474: Output Parameter:
5475: . y - function vector, as set by SNESSetFunction()
5477: Notes:
5478: SNESComputeFunction() is typically used within nonlinear solvers
5479: implementations, so most users would not generally call this routine
5480: themselves.
5482: Level: developer
5484: .keywords: SNES, nonlinear, compute, function
5486: .seealso: SNESSetFunction(), SNESGetFunction()
5487: */
5488: PetscErrorCode SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
5489: {
5490: PetscErrorCode ierr;
5491: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5492: int nlhs = 1,nrhs = 5;
5493: mxArray *plhs[1],*prhs[5];
5494: long long int lx = 0,ly = 0,ls = 0;
5503: /* call Matlab function in ctx with arguments x and y */
5505: PetscMemcpy(&ls,&snes,sizeof(snes));
5506: PetscMemcpy(&lx,&x,sizeof(x));
5507: PetscMemcpy(&ly,&y,sizeof(x));
5508: prhs[0] = mxCreateDoubleScalar((double)ls);
5509: prhs[1] = mxCreateDoubleScalar((double)lx);
5510: prhs[2] = mxCreateDoubleScalar((double)ly);
5511: prhs[3] = mxCreateString(sctx->funcname);
5512: prhs[4] = sctx->ctx;
5513: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");
5514: mxGetScalar(plhs[0]);
5515: mxDestroyArray(prhs[0]);
5516: mxDestroyArray(prhs[1]);
5517: mxDestroyArray(prhs[2]);
5518: mxDestroyArray(prhs[3]);
5519: mxDestroyArray(plhs[0]);
5520: return(0);
5521: }
5523: /*
5524: SNESSetFunctionMatlab - Sets the function evaluation routine and function
5525: vector for use by the SNES routines in solving systems of nonlinear
5526: equations from MATLAB. Here the function is a string containing the name of a MATLAB function
5528: Logically Collective on SNES
5530: Input Parameters:
5531: + snes - the SNES context
5532: . r - vector to store function value
5533: - f - function evaluation routine
5535: Notes:
5536: The Newton-like methods typically solve linear systems of the form
5537: $ f'(x) x = -f(x),
5538: where f'(x) denotes the Jacobian matrix and f(x) is the function.
5540: Level: beginner
5542: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5544: .keywords: SNES, nonlinear, set, function
5546: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5547: */
5548: PetscErrorCode SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx)
5549: {
5550: PetscErrorCode ierr;
5551: SNESMatlabContext *sctx;
5554: /* currently sctx is memory bleed */
5555: PetscNew(&sctx);
5556: PetscStrallocpy(f,&sctx->funcname);
5557: /*
5558: This should work, but it doesn't
5559: sctx->ctx = ctx;
5560: mexMakeArrayPersistent(sctx->ctx);
5561: */
5562: sctx->ctx = mxDuplicateArray(ctx);
5563: SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);
5564: return(0);
5565: }
5567: /*
5568: SNESComputeJacobian_Matlab - Calls the function that has been set with SNESSetJacobianMatlab().
5570: Collective on SNES
5572: Input Parameters:
5573: + snes - the SNES context
5574: . x - input vector
5575: . A, B - the matrices
5576: - ctx - user context
5578: Level: developer
5580: .keywords: SNES, nonlinear, compute, function
5582: .seealso: SNESSetFunction(), SNESGetFunction()
5583: @*/
5584: PetscErrorCode SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx)
5585: {
5586: PetscErrorCode ierr;
5587: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5588: int nlhs = 2,nrhs = 6;
5589: mxArray *plhs[2],*prhs[6];
5590: long long int lx = 0,lA = 0,ls = 0, lB = 0;
5596: /* call Matlab function in ctx with arguments x and y */
5598: PetscMemcpy(&ls,&snes,sizeof(snes));
5599: PetscMemcpy(&lx,&x,sizeof(x));
5600: PetscMemcpy(&lA,A,sizeof(x));
5601: PetscMemcpy(&lB,B,sizeof(x));
5602: prhs[0] = mxCreateDoubleScalar((double)ls);
5603: prhs[1] = mxCreateDoubleScalar((double)lx);
5604: prhs[2] = mxCreateDoubleScalar((double)lA);
5605: prhs[3] = mxCreateDoubleScalar((double)lB);
5606: prhs[4] = mxCreateString(sctx->funcname);
5607: prhs[5] = sctx->ctx;
5608: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");
5609: mxGetScalar(plhs[0]);
5610: mxDestroyArray(prhs[0]);
5611: mxDestroyArray(prhs[1]);
5612: mxDestroyArray(prhs[2]);
5613: mxDestroyArray(prhs[3]);
5614: mxDestroyArray(prhs[4]);
5615: mxDestroyArray(plhs[0]);
5616: mxDestroyArray(plhs[1]);
5617: return(0);
5618: }
5620: /*
5621: SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
5622: vector for use by the SNES routines in solving systems of nonlinear
5623: equations from MATLAB. Here the function is a string containing the name of a MATLAB function
5625: Logically Collective on SNES
5627: Input Parameters:
5628: + snes - the SNES context
5629: . A,B - Jacobian matrices
5630: . J - function evaluation routine
5631: - ctx - user context
5633: Level: developer
5635: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5637: .keywords: SNES, nonlinear, set, function
5639: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J
5640: */
5641: PetscErrorCode SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx)
5642: {
5643: PetscErrorCode ierr;
5644: SNESMatlabContext *sctx;
5647: /* currently sctx is memory bleed */
5648: PetscNew(&sctx);
5649: PetscStrallocpy(J,&sctx->funcname);
5650: /*
5651: This should work, but it doesn't
5652: sctx->ctx = ctx;
5653: mexMakeArrayPersistent(sctx->ctx);
5654: */
5655: sctx->ctx = mxDuplicateArray(ctx);
5656: SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);
5657: return(0);
5658: }
5660: /*
5661: SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().
5663: Collective on SNES
5665: .seealso: SNESSetFunction(), SNESGetFunction()
5666: @*/
5667: PetscErrorCode SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
5668: {
5669: PetscErrorCode ierr;
5670: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5671: int nlhs = 1,nrhs = 6;
5672: mxArray *plhs[1],*prhs[6];
5673: long long int lx = 0,ls = 0;
5674: Vec x = snes->vec_sol;
5679: PetscMemcpy(&ls,&snes,sizeof(snes));
5680: PetscMemcpy(&lx,&x,sizeof(x));
5681: prhs[0] = mxCreateDoubleScalar((double)ls);
5682: prhs[1] = mxCreateDoubleScalar((double)it);
5683: prhs[2] = mxCreateDoubleScalar((double)fnorm);
5684: prhs[3] = mxCreateDoubleScalar((double)lx);
5685: prhs[4] = mxCreateString(sctx->funcname);
5686: prhs[5] = sctx->ctx;
5687: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");
5688: mxGetScalar(plhs[0]);
5689: mxDestroyArray(prhs[0]);
5690: mxDestroyArray(prhs[1]);
5691: mxDestroyArray(prhs[2]);
5692: mxDestroyArray(prhs[3]);
5693: mxDestroyArray(prhs[4]);
5694: mxDestroyArray(plhs[0]);
5695: return(0);
5696: }
5698: /*
5699: SNESMonitorSetMatlab - Sets the monitor function from MATLAB
5701: Level: developer
5703: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5705: .keywords: SNES, nonlinear, set, function
5707: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5708: */
5709: PetscErrorCode SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx)
5710: {
5711: PetscErrorCode ierr;
5712: SNESMatlabContext *sctx;
5715: /* currently sctx is memory bleed */
5716: PetscNew(&sctx);
5717: PetscStrallocpy(f,&sctx->funcname);
5718: /*
5719: This should work, but it doesn't
5720: sctx->ctx = ctx;
5721: mexMakeArrayPersistent(sctx->ctx);
5722: */
5723: sctx->ctx = mxDuplicateArray(ctx);
5724: SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);
5725: return(0);
5726: }
5728: #endif