Actual source code: snes.c
petsc-3.8.4 2018-03-24
2: #include <petsc/private/snesimpl.h>
3: #include <petscdmshell.h>
4: #include <petscdraw.h>
5: #include <petscds.h>
6: #include <petscconvest.h>
8: PetscBool SNESRegisterAllCalled = PETSC_FALSE;
9: PetscFunctionList SNESList = NULL;
11: /* Logging support */
12: PetscClassId SNES_CLASSID, DMSNES_CLASSID;
13: PetscLogEvent SNES_Solve, SNES_FunctionEval, SNES_JacobianEval, SNES_NGSEval, SNES_NGSFuncEval, SNES_NPCSolve, SNES_ObjectiveEval;
15: /*@
16: SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged.
18: Logically Collective on SNES
20: Input Parameters:
21: + snes - iterative context obtained from SNESCreate()
22: - flg - PETSC_TRUE indicates you want the error generated
24: Options database keys:
25: . -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)
27: Level: intermediate
29: Notes:
30: Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve()
31: to determine if it has converged.
33: .keywords: SNES, set, initial guess, nonzero
35: .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
36: @*/
37: PetscErrorCode SNESSetErrorIfNotConverged(SNES snes,PetscBool flg)
38: {
42: snes->errorifnotconverged = flg;
43: return(0);
44: }
46: /*@
47: SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge?
49: Not Collective
51: Input Parameter:
52: . snes - iterative context obtained from SNESCreate()
54: Output Parameter:
55: . flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE
57: Level: intermediate
59: .keywords: SNES, set, initial guess, nonzero
61: .seealso: SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
62: @*/
63: PetscErrorCode SNESGetErrorIfNotConverged(SNES snes,PetscBool *flag)
64: {
68: *flag = snes->errorifnotconverged;
69: return(0);
70: }
72: /*@
73: SNESSetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?
75: Logically Collective on SNES
77: Input Parameters:
78: + snes - the shell SNES
79: - flg - is the residual computed?
81: Level: advanced
83: .seealso: SNESGetAlwaysComputesFinalResidual()
84: @*/
85: PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES snes, PetscBool flg)
86: {
89: snes->alwayscomputesfinalresidual = flg;
90: return(0);
91: }
93: /*@
94: SNESGetAlwaysComputesFinalResidual - does the SNES always compute the residual at the final solution?
96: Logically Collective on SNES
98: Input Parameter:
99: . snes - the shell SNES
101: Output Parameter:
102: . flg - is the residual computed?
104: Level: advanced
106: .seealso: SNESSetAlwaysComputesFinalResidual()
107: @*/
108: PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES snes, PetscBool *flg)
109: {
112: *flg = snes->alwayscomputesfinalresidual;
113: return(0);
114: }
116: /*@
117: SNESSetFunctionDomainError - tells SNES that the input vector to your SNESFunction is not
118: in the functions domain. For example, negative pressure.
120: Logically Collective on SNES
122: Input Parameters:
123: . snes - the SNES context
125: Level: advanced
127: .keywords: SNES, view
129: .seealso: SNESCreate(), SNESSetFunction(), SNESFunction
130: @*/
131: PetscErrorCode SNESSetFunctionDomainError(SNES snes)
132: {
135: if (snes->errorifnotconverged) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"User code indicates input vector is not in the function domain");
136: snes->domainerror = PETSC_TRUE;
137: return(0);
138: }
140: /*@
141: SNESGetFunctionDomainError - Gets the status of the domain error after a call to SNESComputeFunction;
143: Logically Collective on SNES
145: Input Parameters:
146: . snes - the SNES context
148: Output Parameters:
149: . domainerror - Set to PETSC_TRUE if there's a domain error; PETSC_FALSE otherwise.
151: Level: advanced
153: .keywords: SNES, view
155: .seealso: SNESSetFunctionDomainError(), SNESComputeFunction()
156: @*/
157: PetscErrorCode SNESGetFunctionDomainError(SNES snes, PetscBool *domainerror)
158: {
162: *domainerror = snes->domainerror;
163: return(0);
164: }
166: /*@C
167: SNESLoad - Loads a SNES that has been stored in binary with SNESView().
169: Collective on PetscViewer
171: Input Parameters:
172: + newdm - the newly loaded SNES, this needs to have been created with SNESCreate() or
173: some related function before a call to SNESLoad().
174: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
176: Level: intermediate
178: Notes:
179: The type is determined by the data in the file, any type set into the SNES before this call is ignored.
181: Notes for advanced users:
182: Most users should not need to know the details of the binary storage
183: format, since SNESLoad() and TSView() completely hide these details.
184: But for anyone who's interested, the standard binary matrix storage
185: format is
186: .vb
187: has not yet been determined
188: .ve
190: .seealso: PetscViewerBinaryOpen(), SNESView(), MatLoad(), VecLoad()
191: @*/
192: PetscErrorCode SNESLoad(SNES snes, PetscViewer viewer)
193: {
195: PetscBool isbinary;
196: PetscInt classid;
197: char type[256];
198: KSP ksp;
199: DM dm;
200: DMSNES dmsnes;
205: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
206: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
208: PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);
209: if (classid != SNES_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONG,"Not SNES next in file");
210: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
211: SNESSetType(snes, type);
212: if (snes->ops->load) {
213: (*snes->ops->load)(snes,viewer);
214: }
215: SNESGetDM(snes,&dm);
216: DMGetDMSNES(dm,&dmsnes);
217: DMSNESLoad(dmsnes,viewer);
218: SNESGetKSP(snes,&ksp);
219: KSPLoad(ksp,viewer);
220: return(0);
221: }
223: #include <petscdraw.h>
224: #if defined(PETSC_HAVE_SAWS)
225: #include <petscviewersaws.h>
226: #endif
227: /*@C
228: SNESView - Prints the SNES data structure.
230: Collective on SNES
232: Input Parameters:
233: + SNES - the SNES context
234: - viewer - visualization context
236: Options Database Key:
237: . -snes_view - Calls SNESView() at end of SNESSolve()
239: Notes:
240: The available visualization contexts include
241: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
242: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
243: output where only the first processor opens
244: the file. All other processors send their
245: data to the first processor to print.
247: The user can open an alternative visualization context with
248: PetscViewerASCIIOpen() - output to a specified file.
250: Level: beginner
252: .keywords: SNES, view
254: .seealso: PetscViewerASCIIOpen()
255: @*/
256: PetscErrorCode SNESView(SNES snes,PetscViewer viewer)
257: {
258: SNESKSPEW *kctx;
260: KSP ksp;
261: SNESLineSearch linesearch;
262: PetscBool iascii,isstring,isbinary,isdraw;
263: DMSNES dmsnes;
264: #if defined(PETSC_HAVE_SAWS)
265: PetscBool issaws;
266: #endif
270: if (!viewer) {
271: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);
272: }
276: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
277: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
278: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
279: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
280: #if defined(PETSC_HAVE_SAWS)
281: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSAWS,&issaws);
282: #endif
283: if (iascii) {
284: SNESNormSchedule normschedule;
286: PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer);
287: if (!snes->setupcalled) {
288: PetscViewerASCIIPrintf(viewer," SNES has not been set up so information may be incomplete\n");
289: }
290: if (snes->ops->view) {
291: PetscViewerASCIIPushTab(viewer);
292: (*snes->ops->view)(snes,viewer);
293: PetscViewerASCIIPopTab(viewer);
294: }
295: PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);
296: PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n",(double)snes->rtol,(double)snes->abstol,(double)snes->stol);
297: if (snes->usesksp) {
298: PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);
299: }
300: PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);
301: SNESGetNormSchedule(snes, &normschedule);
302: if (normschedule > 0) {PetscViewerASCIIPrintf(viewer," norm schedule %s\n",SNESNormSchedules[normschedule]);}
303: if (snes->gridsequence) {
304: PetscViewerASCIIPrintf(viewer," total number of grid sequence refinements=%D\n",snes->gridsequence);
305: }
306: if (snes->ksp_ewconv) {
307: kctx = (SNESKSPEW*)snes->kspconvctx;
308: if (kctx) {
309: PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);
310: PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",(double)kctx->rtol_0,(double)kctx->rtol_max,(double)kctx->threshold);
311: PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",(double)kctx->gamma,(double)kctx->alpha,(double)kctx->alpha2);
312: }
313: }
314: if (snes->lagpreconditioner == -1) {
315: PetscViewerASCIIPrintf(viewer," Preconditioned is never rebuilt\n");
316: } else if (snes->lagpreconditioner > 1) {
317: PetscViewerASCIIPrintf(viewer," Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);
318: }
319: if (snes->lagjacobian == -1) {
320: PetscViewerASCIIPrintf(viewer," Jacobian is never rebuilt\n");
321: } else if (snes->lagjacobian > 1) {
322: PetscViewerASCIIPrintf(viewer," Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);
323: }
324: } else if (isstring) {
325: const char *type;
326: SNESGetType(snes,&type);
327: PetscViewerStringSPrintf(viewer," %-3.3s",type);
328: } else if (isbinary) {
329: PetscInt classid = SNES_FILE_CLASSID;
330: MPI_Comm comm;
331: PetscMPIInt rank;
332: char type[256];
334: PetscObjectGetComm((PetscObject)snes,&comm);
335: MPI_Comm_rank(comm,&rank);
336: if (!rank) {
337: PetscViewerBinaryWrite(viewer,&classid,1,PETSC_INT,PETSC_FALSE);
338: PetscStrncpy(type,((PetscObject)snes)->type_name,sizeof(type));
339: PetscViewerBinaryWrite(viewer,type,sizeof(type),PETSC_CHAR,PETSC_FALSE);
340: }
341: if (snes->ops->view) {
342: (*snes->ops->view)(snes,viewer);
343: }
344: } else if (isdraw) {
345: PetscDraw draw;
346: char str[36];
347: PetscReal x,y,bottom,h;
349: PetscViewerDrawGetDraw(viewer,0,&draw);
350: PetscDrawGetCurrentPoint(draw,&x,&y);
351: PetscStrcpy(str,"SNES: ");
352: PetscStrcat(str,((PetscObject)snes)->type_name);
353: PetscDrawStringBoxed(draw,x,y,PETSC_DRAW_BLUE,PETSC_DRAW_BLACK,str,NULL,&h);
354: bottom = y - h;
355: PetscDrawPushCurrentPoint(draw,x,bottom);
356: if (snes->ops->view) {
357: (*snes->ops->view)(snes,viewer);
358: }
359: #if defined(PETSC_HAVE_SAWS)
360: } else if (issaws) {
361: PetscMPIInt rank;
362: const char *name;
364: PetscObjectGetName((PetscObject)snes,&name);
365: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
366: if (!((PetscObject)snes)->amsmem && !rank) {
367: char dir[1024];
369: PetscObjectViewSAWs((PetscObject)snes,viewer);
370: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/its",name);
371: PetscStackCallSAWs(SAWs_Register,(dir,&snes->iter,1,SAWs_READ,SAWs_INT));
372: if (!snes->conv_hist) {
373: SNESSetConvergenceHistory(snes,NULL,NULL,PETSC_DECIDE,PETSC_TRUE);
374: }
375: PetscSNPrintf(dir,1024,"/PETSc/Objects/%s/conv_hist",name);
376: PetscStackCallSAWs(SAWs_Register,(dir,snes->conv_hist,10,SAWs_READ,SAWs_DOUBLE));
377: }
378: #endif
379: }
380: if (snes->linesearch) {
381: PetscViewerASCIIPushTab(viewer);
382: SNESGetLineSearch(snes, &linesearch);
383: SNESLineSearchView(linesearch, viewer);
384: PetscViewerASCIIPopTab(viewer);
385: }
386: if (snes->npc && snes->usesnpc) {
387: PetscViewerASCIIPushTab(viewer);
388: SNESView(snes->npc, viewer);
389: PetscViewerASCIIPopTab(viewer);
390: }
391: PetscViewerASCIIPushTab(viewer);
392: DMGetDMSNES(snes->dm,&dmsnes);
393: DMSNESView(dmsnes, viewer);
394: PetscViewerASCIIPopTab(viewer);
395: if (snes->usesksp) {
396: SNESGetKSP(snes,&ksp);
397: PetscViewerASCIIPushTab(viewer);
398: KSPView(ksp,viewer);
399: PetscViewerASCIIPopTab(viewer);
400: }
401: if (isdraw) {
402: PetscDraw draw;
403: PetscViewerDrawGetDraw(viewer,0,&draw);
404: PetscDrawPopCurrentPoint(draw);
405: }
406: return(0);
407: }
409: /*
410: We retain a list of functions that also take SNES command
411: line options. These are called at the end SNESSetFromOptions()
412: */
413: #define MAXSETFROMOPTIONS 5
414: static PetscInt numberofsetfromoptions;
415: static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
417: /*@C
418: SNESAddOptionsChecker - Adds an additional function to check for SNES options.
420: Not Collective
422: Input Parameter:
423: . snescheck - function that checks for options
425: Level: developer
427: .seealso: SNESSetFromOptions()
428: @*/
429: PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
430: {
432: if (numberofsetfromoptions >= MAXSETFROMOPTIONS) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
433: othersetfromoptions[numberofsetfromoptions++] = snescheck;
434: return(0);
435: }
437: extern PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
439: static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool hasOperator, PetscInt version)
440: {
441: Mat J;
442: KSP ksp;
443: PC pc;
444: PetscBool match;
446: MatNullSpace nullsp;
451: if (!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
452: Mat A = snes->jacobian, B = snes->jacobian_pre;
453: MatCreateVecs(A ? A : B, NULL,&snes->vec_func);
454: }
456: if (version == 1) {
457: MatCreateSNESMF(snes,&J);
458: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
459: MatSetFromOptions(J);
460: } else if (version == 2) {
461: if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
462: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
463: SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);
464: #else
465: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)");
466: #endif
467: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2");
469: /* attach any user provided null space that was on Amat to the newly created matrix free matrix */
470: if (snes->jacobian) {
471: MatGetNullSpace(snes->jacobian,&nullsp);
472: if (nullsp) {
473: MatSetNullSpace(J,nullsp);
474: }
475: }
477: PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);
478: if (hasOperator) {
480: /* This version replaces the user provided Jacobian matrix with a
481: matrix-free version but still employs the user-provided preconditioner matrix. */
482: SNESSetJacobian(snes,J,0,0,0);
483: } else {
484: /* This version replaces both the user-provided Jacobian and the user-
485: provided preconditioner Jacobian with the default matrix free version. */
486: if ((snes->npcside== PC_LEFT) && snes->npc) {
487: if (!snes->jacobian){SNESSetJacobian(snes,J,0,0,0);}
488: } else {
489: SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,0);
490: }
491: /* Force no preconditioner */
492: SNESGetKSP(snes,&ksp);
493: KSPGetPC(ksp,&pc);
494: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&match);
495: if (!match) {
496: PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");
497: PCSetType(pc,PCNONE);
498: }
499: }
500: MatDestroy(&J);
501: return(0);
502: }
504: static PetscErrorCode DMRestrictHook_SNESVecSol(DM dmfine,Mat Restrict,Vec Rscale,Mat Inject,DM dmcoarse,void *ctx)
505: {
506: SNES snes = (SNES)ctx;
508: Vec Xfine,Xfine_named = NULL,Xcoarse;
511: if (PetscLogPrintInfo) {
512: PetscInt finelevel,coarselevel,fineclevel,coarseclevel;
513: DMGetRefineLevel(dmfine,&finelevel);
514: DMGetCoarsenLevel(dmfine,&fineclevel);
515: DMGetRefineLevel(dmcoarse,&coarselevel);
516: DMGetCoarsenLevel(dmcoarse,&coarseclevel);
517: PetscInfo4(dmfine,"Restricting SNES solution vector from level %D-%D to level %D-%D\n",finelevel,fineclevel,coarselevel,coarseclevel);
518: }
519: if (dmfine == snes->dm) Xfine = snes->vec_sol;
520: else {
521: DMGetNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);
522: Xfine = Xfine_named;
523: }
524: DMGetNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
525: if (Inject) {
526: MatRestrict(Inject,Xfine,Xcoarse);
527: } else {
528: MatRestrict(Restrict,Xfine,Xcoarse);
529: VecPointwiseMult(Xcoarse,Xcoarse,Rscale);
530: }
531: DMRestoreNamedGlobalVector(dmcoarse,"SNESVecSol",&Xcoarse);
532: if (Xfine_named) {DMRestoreNamedGlobalVector(dmfine,"SNESVecSol",&Xfine_named);}
533: return(0);
534: }
536: static PetscErrorCode DMCoarsenHook_SNESVecSol(DM dm,DM dmc,void *ctx)
537: {
541: DMCoarsenHookAdd(dmc,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,ctx);
542: return(0);
543: }
545: /* This may be called to rediscretize the operator on levels of linear multigrid. The DM shuffle is so the user can
546: * safely call SNESGetDM() in their residual evaluation routine. */
547: static PetscErrorCode KSPComputeOperators_SNES(KSP ksp,Mat A,Mat B,void *ctx)
548: {
549: SNES snes = (SNES)ctx;
551: Mat Asave = A,Bsave = B;
552: Vec X,Xnamed = NULL;
553: DM dmsave;
554: void *ctxsave;
555: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
558: dmsave = snes->dm;
559: KSPGetDM(ksp,&snes->dm);
560: if (dmsave == snes->dm) X = snes->vec_sol; /* We are on the finest level */
561: else { /* We are on a coarser level, this vec was initialized using a DM restrict hook */
562: DMGetNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
563: X = Xnamed;
564: SNESGetJacobian(snes,NULL,NULL,&jac,&ctxsave);
565: /* If the DM's don't match up, the MatFDColoring context needed for the jacobian won't match up either -- fixit. */
566: if (jac == SNESComputeJacobianDefaultColor) {
567: SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefaultColor,0);
568: }
569: }
570: /* put the previous context back */
572: SNESComputeJacobian(snes,X,A,B);
573: if (snes->dm != dmsave && jac == SNESComputeJacobianDefaultColor) {
574: SNESSetJacobian(snes,NULL,NULL,jac,ctxsave);
575: }
577: if (A != Asave || B != Bsave) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"No support for changing matrices at this time");
578: if (Xnamed) {
579: DMRestoreNamedGlobalVector(snes->dm,"SNESVecSol",&Xnamed);
580: }
581: snes->dm = dmsave;
582: return(0);
583: }
585: /*@
586: SNESSetUpMatrices - ensures that matrices are available for SNES, to be called by SNESSetUp_XXX()
588: Collective
590: Input Arguments:
591: . snes - snes to configure
593: Level: developer
595: .seealso: SNESSetUp()
596: @*/
597: PetscErrorCode SNESSetUpMatrices(SNES snes)
598: {
600: DM dm;
601: DMSNES sdm;
604: SNESGetDM(snes,&dm);
605: DMGetDMSNES(dm,&sdm);
606: if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"DMSNES not properly configured");
607: else if (!snes->jacobian && snes->mf) {
608: Mat J;
609: void *functx;
610: MatCreateSNESMF(snes,&J);
611: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
612: MatSetFromOptions(J);
613: SNESGetFunction(snes,NULL,NULL,&functx);
614: SNESSetJacobian(snes,J,J,0,0);
615: MatDestroy(&J);
616: } else if (snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
617: Mat J,B;
618: MatCreateSNESMF(snes,&J);
619: MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);
620: MatSetFromOptions(J);
621: DMCreateMatrix(snes->dm,&B);
622: /* sdm->computejacobian was already set to reach here */
623: SNESSetJacobian(snes,J,B,NULL,NULL);
624: MatDestroy(&J);
625: MatDestroy(&B);
626: } else if (!snes->jacobian_pre) {
627: PetscDS prob;
628: Mat J, B;
629: PetscBool hasPrec = PETSC_FALSE;
631: J = snes->jacobian;
632: DMGetDS(dm, &prob);
633: if (prob) {PetscDSHasJacobianPreconditioner(prob, &hasPrec);}
634: if (J) {PetscObjectReference((PetscObject) J);}
635: else if (hasPrec) {DMCreateMatrix(snes->dm, &J);}
636: DMCreateMatrix(snes->dm, &B);
637: SNESSetJacobian(snes, J ? J : B, B, NULL, NULL);
638: MatDestroy(&J);
639: MatDestroy(&B);
640: }
641: {
642: KSP ksp;
643: SNESGetKSP(snes,&ksp);
644: KSPSetComputeOperators(ksp,KSPComputeOperators_SNES,snes);
645: DMCoarsenHookAdd(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
646: }
647: return(0);
648: }
650: /*@C
651: SNESMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
653: Collective on SNES
655: Input Parameters:
656: + snes - SNES object you wish to monitor
657: . name - the monitor type one is seeking
658: . help - message indicating what monitoring is done
659: . manual - manual page for the monitor
660: . monitor - the monitor function
661: - 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
663: Level: developer
665: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
666: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
667: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
668: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
669: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
670: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
671: PetscOptionsFList(), PetscOptionsEList()
672: @*/
673: PetscErrorCode SNESMonitorSetFromOptions(SNES snes,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(SNES,PetscViewerAndFormat*))
674: {
675: PetscErrorCode ierr;
676: PetscViewer viewer;
677: PetscViewerFormat format;
678: PetscBool flg;
681: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,name,&viewer,&format,&flg);
682: if (flg) {
683: PetscViewerAndFormat *vf;
684: PetscViewerAndFormatCreate(viewer,format,&vf);
685: PetscObjectDereference((PetscObject)viewer);
686: if (monitorsetup) {
687: (*monitorsetup)(snes,vf);
688: }
689: SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
690: }
691: return(0);
692: }
694: /*@
695: SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
697: Collective on SNES
699: Input Parameter:
700: . snes - the SNES context
702: Options Database Keys:
703: + -snes_type <type> - newtonls, newtontr, ngmres, ncg, nrichardson, qn, vi, fas, SNESType for complete list
704: . -snes_stol - convergence tolerance in terms of the norm
705: of the change in the solution between steps
706: . -snes_atol <abstol> - absolute tolerance of residual norm
707: . -snes_rtol <rtol> - relative decrease in tolerance norm from initial
708: . -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
709: . -snes_force_iteration <force> - force SNESSolve() to take at least one iteration
710: . -snes_max_it <max_it> - maximum number of iterations
711: . -snes_max_funcs <max_funcs> - maximum number of function evaluations
712: . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
713: . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
714: . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
715: . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
716: . -snes_trtol <trtol> - trust region tolerance
717: . -snes_no_convergence_test - skip convergence test in nonlinear
718: solver; hence iterations will continue until max_it
719: or some other criterion is reached. Saves expense
720: of convergence test
721: . -snes_monitor [ascii][:filename][:viewer format] - prints residual norm at each iteration. if no filename given prints to stdout
722: . -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration
723: . -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration
724: . -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration
725: . -snes_monitor_lg_residualnorm - plots residual norm at each iteration
726: . -snes_monitor_lg_range - plots residual norm at each iteration
727: . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
728: . -snes_fd_color - use finite differences with coloring to compute Jacobian
729: . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
730: - -snes_converged_reason - print the reason for convergence/divergence after each solve
732: Options Database for Eisenstat-Walker method:
733: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
734: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
735: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
736: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
737: . -snes_ksp_ew_gamma <gamma> - Sets gamma
738: . -snes_ksp_ew_alpha <alpha> - Sets alpha
739: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
740: - -snes_ksp_ew_threshold <threshold> - Sets threshold
742: Notes:
743: To see all options, run your program with the -help option or consult
744: Users-Manual: ch_snes
746: Level: beginner
748: .keywords: SNES, nonlinear, set, options, database
750: .seealso: SNESSetOptionsPrefix()
751: @*/
752: PetscErrorCode SNESSetFromOptions(SNES snes)
753: {
754: PetscBool flg,pcset,persist,set;
755: PetscInt i,indx,lag,grids;
756: const char *deft = SNESNEWTONLS;
757: const char *convtests[] = {"default","skip"};
758: SNESKSPEW *kctx = NULL;
759: char type[256], monfilename[PETSC_MAX_PATH_LEN];
761: PCSide pcside;
762: const char *optionsprefix;
766: SNESRegisterAll();
767: PetscObjectOptionsBegin((PetscObject)snes);
768: if (((PetscObject)snes)->type_name) deft = ((PetscObject)snes)->type_name;
769: PetscOptionsFList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);
770: if (flg) {
771: SNESSetType(snes,type);
772: } else if (!((PetscObject)snes)->type_name) {
773: SNESSetType(snes,deft);
774: }
775: PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->stol,&snes->stol,NULL);
776: PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,NULL);
778: PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,NULL);
779: PetscOptionsReal("-snes_divergence_tolerance","Stop if residual norm increases by this factor","SNESSetDivergenceTolerance",snes->divtol,&snes->divtol,NULL);
780: PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,NULL);
781: PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,NULL);
782: PetscOptionsInt("-snes_max_fail","Maximum nonlinear step failures","SNESSetMaxNonlinearStepFailures",snes->maxFailures,&snes->maxFailures,NULL);
783: PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,NULL);
784: PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,NULL);
785: PetscOptionsBool("-snes_force_iteration","Force SNESSolve() to take at least one iteration","SNESForceIteration",snes->forceiteration,&snes->forceiteration,NULL);
787: PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);
788: if (flg) {
789: SNESSetLagPreconditioner(snes,lag);
790: }
791: PetscOptionsBool("-snes_lag_preconditioner_persists","Preconditioner lagging through multiple solves","SNESSetLagPreconditionerPersists",snes->lagjac_persist,&persist,&flg);
792: if (flg) {
793: SNESSetLagPreconditionerPersists(snes,persist);
794: }
795: PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);
796: if (flg) {
797: SNESSetLagJacobian(snes,lag);
798: }
799: PetscOptionsBool("-snes_lag_jacobian_persists","Jacobian lagging through multiple solves","SNESSetLagJacobianPersists",snes->lagjac_persist,&persist,&flg);
800: if (flg) {
801: SNESSetLagJacobianPersists(snes,persist);
802: }
804: PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);
805: if (flg) {
806: SNESSetGridSequence(snes,grids);
807: }
809: PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);
810: if (flg) {
811: switch (indx) {
812: case 0: SNESSetConvergenceTest(snes,SNESConvergedDefault,NULL,NULL); break;
813: case 1: SNESSetConvergenceTest(snes,SNESConvergedSkip,NULL,NULL); break;
814: }
815: }
817: PetscOptionsEList("-snes_norm_schedule","SNES Norm schedule","SNESSetNormSchedule",SNESNormSchedules,5,"function",&indx,&flg);
818: if (flg) { SNESSetNormSchedule(snes,(SNESNormSchedule)indx); }
820: PetscOptionsEList("-snes_function_type","SNES Norm schedule","SNESSetFunctionType",SNESFunctionTypes,2,"unpreconditioned",&indx,&flg);
821: if (flg) { SNESSetFunctionType(snes,(SNESFunctionType)indx); }
823: kctx = (SNESKSPEW*)snes->kspconvctx;
825: PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,NULL);
827: PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,NULL);
828: PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,NULL);
829: PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,NULL);
830: PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,NULL);
831: PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,NULL);
832: PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,NULL);
833: PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,NULL);
835: flg = PETSC_FALSE;
836: PetscOptionsBool("-snes_check_jacobian","Check each Jacobian with a differenced one","SNESUpdateCheckJacobian",flg,&flg,&set);
837: if (set && flg) {
838: SNESSetUpdate(snes,SNESUpdateCheckJacobian);
839: }
841: flg = PETSC_FALSE;
842: PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,&set);
843: if (set && flg) {SNESMonitorCancel(snes);}
845: SNESMonitorSetFromOptions(snes,"-snes_monitor","Monitor norm of function","SNESMonitorDefault",SNESMonitorDefault,NULL);
846: SNESMonitorSetFromOptions(snes,"-snes_monitor_short","Monitor norm of function with fewer digits","SNESMonitorDefaultShort",SNESMonitorDefaultShort,NULL);
847: SNESMonitorSetFromOptions(snes,"-snes_monitor_range","Monitor range of elements of function","SNESMonitorRange",SNESMonitorRange,NULL);
849: SNESMonitorSetFromOptions(snes,"-snes_monitor_ratio","Monitor ratios of the norm of function for consecutive steps","SNESMonitorRatio",SNESMonitorRatio,SNESMonitorRatioSetUp);
850: SNESMonitorSetFromOptions(snes,"-snes_monitor_field","Monitor norm of function (split into fields)","SNESMonitorDefaultField",SNESMonitorDefaultField,NULL);
851: SNESMonitorSetFromOptions(snes,"-snes_monitor_solution","View solution at each iteration","SNESMonitorSolution",SNESMonitorSolution,NULL);
852: SNESMonitorSetFromOptions(snes,"-snes_monitor_solution_update","View correction at each iteration","SNESMonitorSolutionUpdate",SNESMonitorSolutionUpdate,NULL);
853: SNESMonitorSetFromOptions(snes,"-snes_monitor_residual","View residual at each iteration","SNESMonitorResidual",SNESMonitorResidual,NULL);
854: SNESMonitorSetFromOptions(snes,"-snes_monitor_jacupdate_spectrum","Print the change in the spectrum of the Jacobian","SNESMonitorJacUpdateSpectrum",SNESMonitorJacUpdateSpectrum,NULL);
855: SNESMonitorSetFromOptions(snes,"-snes_monitor_fields","Monitor norm of function per field","SNESMonitorSet",SNESMonitorFields,NULL);
857: PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);
858: if (flg) {PetscPythonMonitorSet((PetscObject)snes,monfilename);}
861: flg = PETSC_FALSE;
862: PetscOptionsBool("-snes_monitor_lg_residualnorm","Plot function norm at each iteration","SNESMonitorLGResidualNorm",flg,&flg,NULL);
863: if (flg) {
864: PetscDrawLG ctx;
866: SNESMonitorLGCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
867: SNESMonitorSet(snes,SNESMonitorLGResidualNorm,ctx,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
868: }
869: flg = PETSC_FALSE;
870: PetscOptionsBool("-snes_monitor_lg_range","Plot function range at each iteration","SNESMonitorLGRange",flg,&flg,NULL);
871: if (flg) {
872: PetscViewer ctx;
874: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,400,300,&ctx);
875: SNESMonitorSet(snes,SNESMonitorLGRange,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);
876: }
880: flg = PETSC_FALSE;
881: PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESComputeJacobianDefault",flg,&flg,NULL);
882: if (flg) {
883: void *functx;
884: DM dm;
885: DMSNES sdm;
886: SNESGetDM(snes,&dm);
887: DMGetDMSNES(dm,&sdm);
888: sdm->jacobianctx = NULL;
889: SNESGetFunction(snes,NULL,NULL,&functx);
890: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefault,functx);
891: PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");
892: }
894: flg = PETSC_FALSE;
895: PetscOptionsBool("-snes_fd_function","Use finite differences (slow) to compute function from user objective","SNESObjectiveComputeFunctionDefaultFD",flg,&flg,NULL);
896: if (flg) {
897: SNESSetFunction(snes,NULL,SNESObjectiveComputeFunctionDefaultFD,NULL);
898: }
900: flg = PETSC_FALSE;
901: PetscOptionsBool("-snes_fd_color","Use finite differences with coloring to compute Jacobian","SNESComputeJacobianDefaultColor",flg,&flg,NULL);
902: if (flg) {
903: DM dm;
904: DMSNES sdm;
905: SNESGetDM(snes,&dm);
906: DMGetDMSNES(dm,&sdm);
907: sdm->jacobianctx = NULL;
908: SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESComputeJacobianDefaultColor,0);
909: PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
910: }
912: flg = PETSC_FALSE;
913: PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf_operator,&flg);
914: if (flg && snes->mf_operator) {
915: snes->mf_operator = PETSC_TRUE;
916: snes->mf = PETSC_TRUE;
917: }
918: flg = PETSC_FALSE;
919: PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","SNESSetUseMatrixFree",PETSC_FALSE,&snes->mf,&flg);
920: if (!flg && snes->mf_operator) snes->mf = PETSC_TRUE;
921: PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",snes->mf_version,&snes->mf_version,0);
923: flg = PETSC_FALSE;
924: SNESGetNPCSide(snes,&pcside);
925: PetscOptionsEnum("-snes_npc_side","SNES nonlinear preconditioner side","SNESSetNPCSide",PCSides,(PetscEnum)pcside,(PetscEnum*)&pcside,&flg);
926: if (flg) {SNESSetNPCSide(snes,pcside);}
928: #if defined(PETSC_HAVE_SAWS)
929: /*
930: Publish convergence information using SAWs
931: */
932: flg = PETSC_FALSE;
933: PetscOptionsBool("-snes_monitor_saws","Publish SNES progress using SAWs","SNESMonitorSet",flg,&flg,NULL);
934: if (flg) {
935: void *ctx;
936: SNESMonitorSAWsCreate(snes,&ctx);
937: SNESMonitorSet(snes,SNESMonitorSAWs,ctx,SNESMonitorSAWsDestroy);
938: }
939: #endif
940: #if defined(PETSC_HAVE_SAWS)
941: {
942: PetscBool set;
943: flg = PETSC_FALSE;
944: PetscOptionsBool("-snes_saws_block","Block for SAWs at end of SNESSolve","PetscObjectSAWsBlock",((PetscObject)snes)->amspublishblock,&flg,&set);
945: if (set) {
946: PetscObjectSAWsSetBlock((PetscObject)snes,flg);
947: }
948: }
949: #endif
951: for (i = 0; i < numberofsetfromoptions; i++) {
952: (*othersetfromoptions[i])(snes);
953: }
955: if (snes->ops->setfromoptions) {
956: (*snes->ops->setfromoptions)(PetscOptionsObject,snes);
957: }
959: /* process any options handlers added with PetscObjectAddOptionsHandler() */
960: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)snes);
961: PetscOptionsEnd();
963: if (!snes->linesearch) {
964: SNESGetLineSearch(snes, &snes->linesearch);
965: }
966: SNESLineSearchSetFromOptions(snes->linesearch);
968: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
969: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
970: KSPSetFromOptions(snes->ksp);
972: /* if someone has set the SNES NPC type, create it. */
973: SNESGetOptionsPrefix(snes, &optionsprefix);
974: PetscOptionsHasName(((PetscObject)snes)->options,optionsprefix, "-npc_snes_type", &pcset);
975: if (pcset && (!snes->npc)) {
976: SNESGetNPC(snes, &snes->npc);
977: }
978: return(0);
979: }
981: /*@C
982: SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
983: the nonlinear solvers.
985: Logically Collective on SNES
987: Input Parameters:
988: + snes - the SNES context
989: . compute - function to compute the context
990: - destroy - function to destroy the context
992: Level: intermediate
994: Notes:
995: This function is currently not available from Fortran.
997: .keywords: SNES, nonlinear, set, application, context
999: .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
1000: @*/
1001: PetscErrorCode SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
1002: {
1005: snes->ops->usercompute = compute;
1006: snes->ops->userdestroy = destroy;
1007: return(0);
1008: }
1010: /*@
1011: SNESSetApplicationContext - Sets the optional user-defined context for
1012: the nonlinear solvers.
1014: Logically Collective on SNES
1016: Input Parameters:
1017: + snes - the SNES context
1018: - usrP - optional user context
1020: Level: intermediate
1022: Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
1023: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1025: .keywords: SNES, nonlinear, set, application, context
1027: .seealso: SNESGetApplicationContext()
1028: @*/
1029: PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP)
1030: {
1032: KSP ksp;
1036: SNESGetKSP(snes,&ksp);
1037: KSPSetApplicationContext(ksp,usrP);
1038: snes->user = usrP;
1039: return(0);
1040: }
1042: /*@
1043: SNESGetApplicationContext - Gets the user-defined context for the
1044: nonlinear solvers.
1046: Not Collective
1048: Input Parameter:
1049: . snes - SNES context
1051: Output Parameter:
1052: . usrP - user context
1054: Fortran Notes: To use this from Fortran you must write a Fortran interface definition for this
1055: function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.
1057: Level: intermediate
1059: .keywords: SNES, nonlinear, get, application, context
1061: .seealso: SNESSetApplicationContext()
1062: @*/
1063: PetscErrorCode SNESGetApplicationContext(SNES snes,void *usrP)
1064: {
1067: *(void**)usrP = snes->user;
1068: return(0);
1069: }
1071: /*@
1072: SNESSetUseMatrixFree - indicates that SNES should use matrix free finite difference matrix vector products internally to apply
1073: the Jacobian.
1075: Collective on SNES
1077: Input Parameters:
1078: + snes - SNES context
1079: . mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1080: - mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used
1082: Options Database:
1083: + -snes_mf - use matrix free for both the mat and pmat operator
1084: - -snes_mf_operator - use matrix free only for the mat operator
1086: Level: intermediate
1088: .keywords: SNES, nonlinear, get, iteration, number,
1090: .seealso: SNESGetUseMatrixFree(), MatCreateSNESMF()
1091: @*/
1092: PetscErrorCode SNESSetUseMatrixFree(SNES snes,PetscBool mf_operator,PetscBool mf)
1093: {
1098: if (mf && !mf_operator) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"If using mf must also use mf_operator");
1099: snes->mf = mf;
1100: snes->mf_operator = mf_operator;
1101: return(0);
1102: }
1104: /*@
1105: SNESGetUseMatrixFree - indicates if the SNES uses matrix free finite difference matrix vector products to apply
1106: the Jacobian.
1108: Collective on SNES
1110: Input Parameter:
1111: . snes - SNES context
1113: Output Parameters:
1114: + mf - use matrix-free for both the Amat and Pmat used by SNESSetJacobian(), both the Amat and Pmat set in SNESSetJacobian() will be ignored
1115: - mf_operator - use matrix-free only for the Amat used by SNESSetJacobian(), this means the user provided Pmat will continue to be used
1117: Options Database:
1118: + -snes_mf - use matrix free for both the mat and pmat operator
1119: - -snes_mf_operator - use matrix free only for the mat operator
1121: Level: intermediate
1123: .keywords: SNES, nonlinear, get, iteration, number,
1125: .seealso: SNESSetUseMatrixFree(), MatCreateSNESMF()
1126: @*/
1127: PetscErrorCode SNESGetUseMatrixFree(SNES snes,PetscBool *mf_operator,PetscBool *mf)
1128: {
1131: if (mf) *mf = snes->mf;
1132: if (mf_operator) *mf_operator = snes->mf_operator;
1133: return(0);
1134: }
1136: /*@
1137: SNESGetIterationNumber - Gets the number of nonlinear iterations completed
1138: at this time.
1140: Not Collective
1142: Input Parameter:
1143: . snes - SNES context
1145: Output Parameter:
1146: . iter - iteration number
1148: Notes:
1149: For example, during the computation of iteration 2 this would return 1.
1151: This is useful for using lagged Jacobians (where one does not recompute the
1152: Jacobian at each SNES iteration). For example, the code
1153: .vb
1154: SNESGetIterationNumber(snes,&it);
1155: if (!(it % 2)) {
1156: [compute Jacobian here]
1157: }
1158: .ve
1159: can be used in your ComputeJacobian() function to cause the Jacobian to be
1160: recomputed every second SNES iteration.
1162: Level: intermediate
1164: .keywords: SNES, nonlinear, get, iteration, number,
1166: .seealso: SNESGetLinearSolveIterations()
1167: @*/
1168: PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt *iter)
1169: {
1173: *iter = snes->iter;
1174: return(0);
1175: }
1177: /*@
1178: SNESSetIterationNumber - Sets the current iteration number.
1180: Not Collective
1182: Input Parameter:
1183: . snes - SNES context
1184: . iter - iteration number
1186: Level: developer
1188: .keywords: SNES, nonlinear, set, iteration, number,
1190: .seealso: SNESGetLinearSolveIterations()
1191: @*/
1192: PetscErrorCode SNESSetIterationNumber(SNES snes,PetscInt iter)
1193: {
1198: PetscObjectSAWsTakeAccess((PetscObject)snes);
1199: snes->iter = iter;
1200: PetscObjectSAWsGrantAccess((PetscObject)snes);
1201: return(0);
1202: }
1204: /*@
1205: SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
1206: attempted by the nonlinear solver.
1208: Not Collective
1210: Input Parameter:
1211: . snes - SNES context
1213: Output Parameter:
1214: . nfails - number of unsuccessful steps attempted
1216: Notes:
1217: This counter is reset to zero for each successive call to SNESSolve().
1219: Level: intermediate
1221: .keywords: SNES, nonlinear, get, number, unsuccessful, steps
1223: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1224: SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
1225: @*/
1226: PetscErrorCode SNESGetNonlinearStepFailures(SNES snes,PetscInt *nfails)
1227: {
1231: *nfails = snes->numFailures;
1232: return(0);
1233: }
1235: /*@
1236: SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
1237: attempted by the nonlinear solver before it gives up.
1239: Not Collective
1241: Input Parameters:
1242: + snes - SNES context
1243: - maxFails - maximum of unsuccessful steps
1245: Level: intermediate
1247: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
1249: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1250: SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1251: @*/
1252: PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
1253: {
1256: snes->maxFailures = maxFails;
1257: return(0);
1258: }
1260: /*@
1261: SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
1262: attempted by the nonlinear solver before it gives up.
1264: Not Collective
1266: Input Parameter:
1267: . snes - SNES context
1269: Output Parameter:
1270: . maxFails - maximum of unsuccessful steps
1272: Level: intermediate
1274: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1276: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
1277: SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
1279: @*/
1280: PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
1281: {
1285: *maxFails = snes->maxFailures;
1286: return(0);
1287: }
1289: /*@
1290: SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
1291: done by SNES.
1293: Not Collective
1295: Input Parameter:
1296: . snes - SNES context
1298: Output Parameter:
1299: . nfuncs - number of evaluations
1301: Level: intermediate
1303: Notes: Reset every time SNESSolve is called unless SNESSetCountersReset() is used.
1305: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1307: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), SNESSetCountersReset()
1308: @*/
1309: PetscErrorCode SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
1310: {
1314: *nfuncs = snes->nfuncs;
1315: return(0);
1316: }
1318: /*@
1319: SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
1320: linear solvers.
1322: Not Collective
1324: Input Parameter:
1325: . snes - SNES context
1327: Output Parameter:
1328: . nfails - number of failed solves
1330: Level: intermediate
1332: Options Database Keys:
1333: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1335: Notes:
1336: This counter is reset to zero for each successive call to SNESSolve().
1338: .keywords: SNES, nonlinear, get, number, unsuccessful, steps
1340: .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
1341: @*/
1342: PetscErrorCode SNESGetLinearSolveFailures(SNES snes,PetscInt *nfails)
1343: {
1347: *nfails = snes->numLinearSolveFailures;
1348: return(0);
1349: }
1351: /*@
1352: SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
1353: allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
1355: Logically Collective on SNES
1357: Input Parameters:
1358: + snes - SNES context
1359: - maxFails - maximum allowed linear solve failures
1361: Level: intermediate
1363: Options Database Keys:
1364: . -snes_max_linear_solve_fail <num> - The number of failures before the solve is terminated
1366: Notes: By default this is 0; that is SNES returns on the first failed linear solve
1368: .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
1370: .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
1371: @*/
1372: PetscErrorCode SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
1373: {
1377: snes->maxLinearSolveFailures = maxFails;
1378: return(0);
1379: }
1381: /*@
1382: SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
1383: are allowed before SNES terminates
1385: Not Collective
1387: Input Parameter:
1388: . snes - SNES context
1390: Output Parameter:
1391: . maxFails - maximum of unsuccessful solves allowed
1393: Level: intermediate
1395: Notes: By default this is 1; that is SNES returns on the first failed linear solve
1397: .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
1399: .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
1400: @*/
1401: PetscErrorCode SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
1402: {
1406: *maxFails = snes->maxLinearSolveFailures;
1407: return(0);
1408: }
1410: /*@
1411: SNESGetLinearSolveIterations - Gets the total number of linear iterations
1412: used by the nonlinear solver.
1414: Not Collective
1416: Input Parameter:
1417: . snes - SNES context
1419: Output Parameter:
1420: . lits - number of linear iterations
1422: Notes:
1423: This counter is reset to zero for each successive call to SNESSolve() unless SNESSetCountersReset() is used.
1425: 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
1426: then call KSPGetIterationNumber() after the failed solve.
1428: Level: intermediate
1430: .keywords: SNES, nonlinear, get, number, linear, iterations
1432: .seealso: SNESGetIterationNumber(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESSetCountersReset()
1433: @*/
1434: PetscErrorCode SNESGetLinearSolveIterations(SNES snes,PetscInt *lits)
1435: {
1439: *lits = snes->linear_its;
1440: return(0);
1441: }
1443: /*@
1444: SNESSetCountersReset - Sets whether or not the counters for linear iterations and function evaluations
1445: are reset every time SNESSolve() is called.
1447: Logically Collective on SNES
1449: Input Parameter:
1450: + snes - SNES context
1451: - reset - whether to reset the counters or not
1453: Notes:
1454: This defaults to PETSC_TRUE
1456: Level: developer
1458: .keywords: SNES, nonlinear, set, reset, number, linear, iterations
1460: .seealso: SNESGetNumberFunctionEvals(), SNESGetLinearSolveIterations(), SNESGetNPC()
1461: @*/
1462: PetscErrorCode SNESSetCountersReset(SNES snes,PetscBool reset)
1463: {
1467: snes->counters_reset = reset;
1468: return(0);
1469: }
1472: /*@
1473: SNESSetKSP - Sets a KSP context for the SNES object to use
1475: Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
1477: Input Parameters:
1478: + snes - the SNES context
1479: - ksp - the KSP context
1481: Notes:
1482: The SNES object already has its KSP object, you can obtain with SNESGetKSP()
1483: so this routine is rarely needed.
1485: The KSP object that is already in the SNES object has its reference count
1486: decreased by one.
1488: Level: developer
1490: .keywords: SNES, nonlinear, get, KSP, context
1492: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1493: @*/
1494: PetscErrorCode SNESSetKSP(SNES snes,KSP ksp)
1495: {
1502: PetscObjectReference((PetscObject)ksp);
1503: if (snes->ksp) {PetscObjectDereference((PetscObject)snes->ksp);}
1504: snes->ksp = ksp;
1505: return(0);
1506: }
1508: /* -----------------------------------------------------------*/
1509: /*@
1510: SNESCreate - Creates a nonlinear solver context.
1512: Collective on MPI_Comm
1514: Input Parameters:
1515: . comm - MPI communicator
1517: Output Parameter:
1518: . outsnes - the new SNES context
1520: Options Database Keys:
1521: + -snes_mf - Activates default matrix-free Jacobian-vector products,
1522: and no preconditioning matrix
1523: . -snes_mf_operator - Activates default matrix-free Jacobian-vector
1524: products, and a user-provided preconditioning matrix
1525: as set by SNESSetJacobian()
1526: - -snes_fd - Uses (slow!) finite differences to compute Jacobian
1528: Level: beginner
1530: Developer Notes: SNES always creates a KSP object even though many SNES methods do not use it. This is
1531: unfortunate and should be fixed at some point. The flag snes->usesksp indicates if the
1532: particular method does use KSP and regulates if the information about the KSP is printed
1533: in SNESView(). TSSetFromOptions() does call SNESSetFromOptions() which can lead to users being confused
1534: by help messages about meaningless SNES options.
1536: SNES always creates the snes->kspconvctx even though it is used by only one type. This should
1537: be fixed.
1539: .keywords: SNES, nonlinear, create, context
1541: .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()
1543: @*/
1544: PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes)
1545: {
1547: SNES snes;
1548: SNESKSPEW *kctx;
1552: *outsnes = NULL;
1553: SNESInitializePackage();
1555: PetscHeaderCreate(snes,SNES_CLASSID,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);
1557: snes->ops->converged = SNESConvergedDefault;
1558: snes->usesksp = PETSC_TRUE;
1559: snes->tolerancesset = PETSC_FALSE;
1560: snes->max_its = 50;
1561: snes->max_funcs = 10000;
1562: snes->norm = 0.0;
1563: snes->normschedule = SNES_NORM_ALWAYS;
1564: snes->functype = SNES_FUNCTION_DEFAULT;
1565: #if defined(PETSC_USE_REAL_SINGLE)
1566: snes->rtol = 1.e-5;
1567: #else
1568: snes->rtol = 1.e-8;
1569: #endif
1570: snes->ttol = 0.0;
1571: #if defined(PETSC_USE_REAL_SINGLE)
1572: snes->abstol = 1.e-25;
1573: #else
1574: snes->abstol = 1.e-50;
1575: #endif
1576: #if defined(PETSC_USE_REAL_SINGLE)
1577: snes->stol = 1.e-5;
1578: #else
1579: snes->stol = 1.e-8;
1580: #endif
1581: #if defined(PETSC_USE_REAL_SINGLE)
1582: snes->deltatol = 1.e-6;
1583: #else
1584: snes->deltatol = 1.e-12;
1585: #endif
1586: snes->divtol = 1.e4;
1587: snes->rnorm0 = 0;
1588: snes->nfuncs = 0;
1589: snes->numFailures = 0;
1590: snes->maxFailures = 1;
1591: snes->linear_its = 0;
1592: snes->lagjacobian = 1;
1593: snes->jac_iter = 0;
1594: snes->lagjac_persist = PETSC_FALSE;
1595: snes->lagpreconditioner = 1;
1596: snes->pre_iter = 0;
1597: snes->lagpre_persist = PETSC_FALSE;
1598: snes->numbermonitors = 0;
1599: snes->data = 0;
1600: snes->setupcalled = PETSC_FALSE;
1601: snes->ksp_ewconv = PETSC_FALSE;
1602: snes->nwork = 0;
1603: snes->work = 0;
1604: snes->nvwork = 0;
1605: snes->vwork = 0;
1606: snes->conv_hist_len = 0;
1607: snes->conv_hist_max = 0;
1608: snes->conv_hist = NULL;
1609: snes->conv_hist_its = NULL;
1610: snes->conv_hist_reset = PETSC_TRUE;
1611: snes->counters_reset = PETSC_TRUE;
1612: snes->vec_func_init_set = PETSC_FALSE;
1613: snes->reason = SNES_CONVERGED_ITERATING;
1614: snes->npcside = PC_RIGHT;
1616: snes->mf = PETSC_FALSE;
1617: snes->mf_operator = PETSC_FALSE;
1618: snes->mf_version = 1;
1620: snes->numLinearSolveFailures = 0;
1621: snes->maxLinearSolveFailures = 1;
1623: snes->vizerotolerance = 1.e-8;
1625: /* Set this to true if the implementation of SNESSolve_XXX does compute the residual at the final solution. */
1626: snes->alwayscomputesfinalresidual = PETSC_FALSE;
1628: /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
1629: PetscNewLog(snes,&kctx);
1631: snes->kspconvctx = (void*)kctx;
1632: kctx->version = 2;
1633: kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1634: this was too large for some test cases */
1635: kctx->rtol_last = 0.0;
1636: kctx->rtol_max = .9;
1637: kctx->gamma = 1.0;
1638: kctx->alpha = .5*(1.0 + PetscSqrtReal(5.0));
1639: kctx->alpha2 = kctx->alpha;
1640: kctx->threshold = .1;
1641: kctx->lresid_last = 0.0;
1642: kctx->norm_last = 0.0;
1644: *outsnes = snes;
1645: return(0);
1646: }
1648: /*MC
1649: SNESFunction - Functional form used to convey the nonlinear function to be solved by SNES
1651: Synopsis:
1652: #include "petscsnes.h"
1653: PetscErrorCode SNESFunction(SNES snes,Vec x,Vec f,void *ctx);
1655: Input Parameters:
1656: + snes - the SNES context
1657: . x - state at which to evaluate residual
1658: - ctx - optional user-defined function context, passed in with SNESSetFunction()
1660: Output Parameter:
1661: . f - vector to put residual (function value)
1663: Level: intermediate
1665: .seealso: SNESSetFunction(), SNESGetFunction()
1666: M*/
1668: /*@C
1669: SNESSetFunction - Sets the function evaluation routine and function
1670: vector for use by the SNES routines in solving systems of nonlinear
1671: equations.
1673: Logically Collective on SNES
1675: Input Parameters:
1676: + snes - the SNES context
1677: . r - vector to store function value
1678: . f - function evaluation routine; see SNESFunction for calling sequence details
1679: - ctx - [optional] user-defined context for private data for the
1680: function evaluation routine (may be NULL)
1682: Notes:
1683: The Newton-like methods typically solve linear systems of the form
1684: $ f'(x) x = -f(x),
1685: where f'(x) denotes the Jacobian matrix and f(x) is the function.
1687: Level: beginner
1689: .keywords: SNES, nonlinear, set, function
1691: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard(), SNESFunction
1692: @*/
1693: PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1694: {
1696: DM dm;
1700: if (r) {
1703: PetscObjectReference((PetscObject)r);
1704: VecDestroy(&snes->vec_func);
1706: snes->vec_func = r;
1707: }
1708: SNESGetDM(snes,&dm);
1709: DMSNESSetFunction(dm,f,ctx);
1710: return(0);
1711: }
1714: /*@C
1715: SNESSetInitialFunction - Sets the function vector to be used as the
1716: function norm at the initialization of the method. In some
1717: instances, the user has precomputed the function before calling
1718: SNESSolve. This function allows one to avoid a redundant call
1719: to SNESComputeFunction in that case.
1721: Logically Collective on SNES
1723: Input Parameters:
1724: + snes - the SNES context
1725: - f - vector to store function value
1727: Notes:
1728: This should not be modified during the solution procedure.
1730: This is used extensively in the SNESFAS hierarchy and in nonlinear preconditioning.
1732: Level: developer
1734: .keywords: SNES, nonlinear, set, function
1736: .seealso: SNESSetFunction(), SNESComputeFunction(), SNESSetInitialFunctionNorm()
1737: @*/
1738: PetscErrorCode SNESSetInitialFunction(SNES snes, Vec f)
1739: {
1741: Vec vec_func;
1747: if (snes->npcside== PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1748: snes->vec_func_init_set = PETSC_FALSE;
1749: return(0);
1750: }
1751: SNESGetFunction(snes,&vec_func,NULL,NULL);
1752: VecCopy(f, vec_func);
1754: snes->vec_func_init_set = PETSC_TRUE;
1755: return(0);
1756: }
1758: /*@
1759: SNESSetNormSchedule - Sets the SNESNormSchedule used in covergence and monitoring
1760: of the SNES method.
1762: Logically Collective on SNES
1764: Input Parameters:
1765: + snes - the SNES context
1766: - normschedule - the frequency of norm computation
1768: Options Database Key:
1769: . -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly>
1771: Notes:
1772: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1773: of the nonlinear function and the taking of its norm at every iteration to
1774: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1775: (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
1776: may either be monitored for convergence or not. As these are often used as nonlinear
1777: preconditioners, monitoring the norm of their error is not a useful enterprise within
1778: their solution.
1780: Level: developer
1782: .keywords: SNES, nonlinear, set, function, norm, type
1784: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1785: @*/
1786: PetscErrorCode SNESSetNormSchedule(SNES snes, SNESNormSchedule normschedule)
1787: {
1790: snes->normschedule = normschedule;
1791: return(0);
1792: }
1795: /*@
1796: SNESGetNormSchedule - Gets the SNESNormSchedule used in covergence and monitoring
1797: of the SNES method.
1799: Logically Collective on SNES
1801: Input Parameters:
1802: + snes - the SNES context
1803: - normschedule - the type of the norm used
1805: Level: advanced
1807: .keywords: SNES, nonlinear, set, function, norm, type
1809: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1810: @*/
1811: PetscErrorCode SNESGetNormSchedule(SNES snes, SNESNormSchedule *normschedule)
1812: {
1815: *normschedule = snes->normschedule;
1816: return(0);
1817: }
1820: /*@
1821: SNESSetFunctionNorm - Sets the last computed residual norm.
1823: Logically Collective on SNES
1825: Input Parameters:
1826: + snes - the SNES context
1828: - normschedule - the frequency of norm computation
1830: Level: developer
1832: .keywords: SNES, nonlinear, set, function, norm, type
1833: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1834: @*/
1835: PetscErrorCode SNESSetFunctionNorm(SNES snes, PetscReal norm)
1836: {
1839: snes->norm = norm;
1840: return(0);
1841: }
1843: /*@
1844: SNESGetFunctionNorm - Gets the last computed norm of the residual
1846: Not Collective
1848: Input Parameter:
1849: . snes - the SNES context
1851: Output Parameter:
1852: . norm - the last computed residual norm
1854: Level: developer
1856: .keywords: SNES, nonlinear, set, function, norm, type
1857: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1858: @*/
1859: PetscErrorCode SNESGetFunctionNorm(SNES snes, PetscReal *norm)
1860: {
1864: *norm = snes->norm;
1865: return(0);
1866: }
1868: /*@C
1869: SNESSetFunctionType - Sets the SNESNormSchedule used in covergence and monitoring
1870: of the SNES method.
1872: Logically Collective on SNES
1874: Input Parameters:
1875: + snes - the SNES context
1876: - normschedule - the frequency of norm computation
1878: Notes:
1879: Only certain SNES methods support certain SNESNormSchedules. Most require evaluation
1880: of the nonlinear function and the taking of its norm at every iteration to
1881: even ensure convergence at all. However, methods such as custom Gauss-Seidel methods
1882: (SNESNGS) and the like do not require the norm of the function to be computed, and therfore
1883: may either be monitored for convergence or not. As these are often used as nonlinear
1884: preconditioners, monitoring the norm of their error is not a useful enterprise within
1885: their solution.
1887: Level: developer
1889: .keywords: SNES, nonlinear, set, function, norm, type
1891: .seealso: SNESGetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1892: @*/
1893: PetscErrorCode SNESSetFunctionType(SNES snes, SNESFunctionType type)
1894: {
1897: snes->functype = type;
1898: return(0);
1899: }
1902: /*@C
1903: SNESGetFunctionType - Gets the SNESNormSchedule used in covergence and monitoring
1904: of the SNES method.
1906: Logically Collective on SNES
1908: Input Parameters:
1909: + snes - the SNES context
1910: - normschedule - the type of the norm used
1912: Level: advanced
1914: .keywords: SNES, nonlinear, set, function, norm, type
1916: .seealso: SNESSetNormSchedule(), SNESComputeFunction(), VecNorm(), SNESSetFunction(), SNESSetInitialFunction(), SNESNormSchedule
1917: @*/
1918: PetscErrorCode SNESGetFunctionType(SNES snes, SNESFunctionType *type)
1919: {
1922: *type = snes->functype;
1923: return(0);
1924: }
1926: /*MC
1927: SNESNGSFunction - function used to convey a Gauss-Seidel sweep on the nonlinear function
1929: Synopsis:
1930: #include <petscsnes.h>
1931: $ SNESNGSFunction(SNES snes,Vec x,Vec b,void *ctx);
1933: + X - solution vector
1934: . B - RHS vector
1935: - ctx - optional user-defined Gauss-Seidel context
1937: Level: intermediate
1939: .seealso: SNESSetNGS(), SNESGetNGS()
1940: M*/
1942: /*@C
1943: SNESSetNGS - Sets the user nonlinear Gauss-Seidel routine for
1944: use with composed nonlinear solvers.
1946: Input Parameters:
1947: + snes - the SNES context
1948: . f - function evaluation routine to apply Gauss-Seidel see SNESNGSFunction
1949: - ctx - [optional] user-defined context for private data for the
1950: smoother evaluation routine (may be NULL)
1952: Notes:
1953: The NGS routines are used by the composed nonlinear solver to generate
1954: a problem appropriate update to the solution, particularly FAS.
1956: Level: intermediate
1958: .keywords: SNES, nonlinear, set, Gauss-Seidel
1960: .seealso: SNESGetFunction(), SNESComputeNGS()
1961: @*/
1962: PetscErrorCode SNESSetNGS(SNES snes,PetscErrorCode (*f)(SNES,Vec,Vec,void*),void *ctx)
1963: {
1965: DM dm;
1969: SNESGetDM(snes,&dm);
1970: DMSNESSetNGS(dm,f,ctx);
1971: return(0);
1972: }
1974: PETSC_EXTERN PetscErrorCode SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
1975: {
1977: DM dm;
1978: DMSNES sdm;
1981: SNESGetDM(snes,&dm);
1982: DMGetDMSNES(dm,&sdm);
1983: /* A(x)*x - b(x) */
1984: if (sdm->ops->computepfunction) {
1985: (*sdm->ops->computepfunction)(snes,x,f,sdm->pctx);
1986: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard function.");
1988: if (sdm->ops->computepjacobian) {
1989: (*sdm->ops->computepjacobian)(snes,x,snes->jacobian,snes->jacobian_pre,sdm->pctx);
1990: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetPicard() to provide Picard matrix.");
1991: VecScale(f,-1.0);
1992: MatMultAdd(snes->jacobian,x,f,f);
1993: return(0);
1994: }
1996: PETSC_EXTERN PetscErrorCode SNESPicardComputeJacobian(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
1997: {
1999: /* the jacobian matrix should be pre-filled in SNESPicardComputeFunction */
2000: return(0);
2001: }
2003: /*@C
2004: SNESSetPicard - Use SNES to solve the semilinear-system A(x) x = b(x) via a Picard type iteration (Picard linearization)
2006: Logically Collective on SNES
2008: Input Parameters:
2009: + snes - the SNES context
2010: . r - vector to store function value
2011: . b - function evaluation routine
2012: . Amat - matrix with which A(x) x - b(x) is to be computed
2013: . Pmat - matrix from which preconditioner is computed (usually the same as Amat)
2014: . J - function to compute matrix value, see SNESJacobianFunction for details on its calling sequence
2015: - ctx - [optional] user-defined context for private data for the
2016: function evaluation routine (may be NULL)
2018: Notes:
2019: 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
2020: 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.
2022: One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both
2024: $ 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}
2025: $ Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration.
2027: Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner.
2029: We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then
2030: the direct Picard iteration A(x^n) x^{n+1} = b(x^n)
2032: 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
2033: 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
2034: different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).
2036: Level: intermediate
2038: .keywords: SNES, nonlinear, set, function
2040: .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard(), SNESJacobianFunction
2041: @*/
2042: 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)
2043: {
2045: DM dm;
2049: SNESGetDM(snes, &dm);
2050: DMSNESSetPicard(dm,b,J,ctx);
2051: SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);
2052: SNESSetJacobian(snes,Amat,Pmat,SNESPicardComputeJacobian,ctx);
2053: return(0);
2054: }
2056: /*@C
2057: SNESGetPicard - Returns the context for the Picard iteration
2059: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
2061: Input Parameter:
2062: . snes - the SNES context
2064: Output Parameter:
2065: + r - the function (or NULL)
2066: . f - the function (or NULL); see SNESFunction for calling sequence details
2067: . Amat - the matrix used to defined the operation A(x) x - b(x) (or NULL)
2068: . Pmat - the matrix from which the preconditioner will be constructed (or NULL)
2069: . J - the function for matrix evaluation (or NULL); see SNESJacobianFunction for calling sequence details
2070: - ctx - the function context (or NULL)
2072: Level: advanced
2074: .keywords: SNES, nonlinear, get, function
2076: .seealso: SNESSetPicard(), SNESGetFunction(), SNESGetJacobian(), SNESGetDM(), SNESFunction, SNESJacobianFunction
2077: @*/
2078: 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)
2079: {
2081: DM dm;
2085: SNESGetFunction(snes,r,NULL,NULL);
2086: SNESGetJacobian(snes,Amat,Pmat,NULL,NULL);
2087: SNESGetDM(snes,&dm);
2088: DMSNESGetPicard(dm,f,J,ctx);
2089: return(0);
2090: }
2092: /*@C
2093: SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem
2095: Logically Collective on SNES
2097: Input Parameters:
2098: + snes - the SNES context
2099: . func - function evaluation routine
2100: - ctx - [optional] user-defined context for private data for the
2101: function evaluation routine (may be NULL)
2103: Calling sequence of func:
2104: $ func (SNES snes,Vec x,void *ctx);
2106: . f - function vector
2107: - ctx - optional user-defined function context
2109: Level: intermediate
2111: .keywords: SNES, nonlinear, set, function
2113: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
2114: @*/
2115: PetscErrorCode SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
2116: {
2119: if (func) snes->ops->computeinitialguess = func;
2120: if (ctx) snes->initialguessP = ctx;
2121: return(0);
2122: }
2124: /* --------------------------------------------------------------- */
2125: /*@C
2126: SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
2127: it assumes a zero right hand side.
2129: Logically Collective on SNES
2131: Input Parameter:
2132: . snes - the SNES context
2134: Output Parameter:
2135: . rhs - the right hand side vector or NULL if the right hand side vector is null
2137: Level: intermediate
2139: .keywords: SNES, nonlinear, get, function, right hand side
2141: .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
2142: @*/
2143: PetscErrorCode SNESGetRhs(SNES snes,Vec *rhs)
2144: {
2148: *rhs = snes->vec_rhs;
2149: return(0);
2150: }
2152: /*@
2153: SNESComputeFunction - Calls the function that has been set with SNESSetFunction().
2155: Collective on SNES
2157: Input Parameters:
2158: + snes - the SNES context
2159: - x - input vector
2161: Output Parameter:
2162: . y - function vector, as set by SNESSetFunction()
2164: Notes:
2165: SNESComputeFunction() is typically used within nonlinear solvers
2166: implementations, so most users would not generally call this routine
2167: themselves.
2169: Level: developer
2171: .keywords: SNES, nonlinear, compute, function
2173: .seealso: SNESSetFunction(), SNESGetFunction()
2174: @*/
2175: PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y)
2176: {
2178: DM dm;
2179: DMSNES sdm;
2187: VecValidValues(x,2,PETSC_TRUE);
2189: SNESGetDM(snes,&dm);
2190: DMGetDMSNES(dm,&sdm);
2191: if (sdm->ops->computefunction) {
2192: if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2193: PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);
2194: }
2195: VecLockPush(x);
2196: PetscStackPush("SNES user function");
2197: (*sdm->ops->computefunction)(snes,x,y,sdm->functionctx);
2198: PetscStackPop;
2199: VecLockPop(x);
2200: if (sdm->ops->computefunction != SNESObjectiveComputeFunctionDefaultFD) {
2201: PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);
2202: }
2203: } else if (snes->vec_rhs) {
2204: MatMult(snes->jacobian, x, y);
2205: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
2206: if (snes->vec_rhs) {
2207: VecAXPY(y,-1.0,snes->vec_rhs);
2208: }
2209: snes->nfuncs++;
2210: /*
2211: domainerror might not be set on all processes; so we tag vector locally with Inf and the next inner product or norm will
2212: propagate the value to all processes
2213: */
2214: if (snes->domainerror) {
2215: VecSetInf(y);
2216: }
2217: return(0);
2218: }
2220: /*@
2221: SNESComputeNGS - Calls the Gauss-Seidel function that has been set with SNESSetNGS().
2223: Collective on SNES
2225: Input Parameters:
2226: + snes - the SNES context
2227: . x - input vector
2228: - b - rhs vector
2230: Output Parameter:
2231: . x - new solution vector
2233: Notes:
2234: SNESComputeNGS() is typically used within composed nonlinear solver
2235: implementations, so most users would not generally call this routine
2236: themselves.
2238: Level: developer
2240: .keywords: SNES, nonlinear, compute, function
2242: .seealso: SNESSetNGS(), SNESComputeFunction()
2243: @*/
2244: PetscErrorCode SNESComputeNGS(SNES snes,Vec b,Vec x)
2245: {
2247: DM dm;
2248: DMSNES sdm;
2256: if (b) {VecValidValues(b,2,PETSC_TRUE);}
2257: PetscLogEventBegin(SNES_NGSEval,snes,x,b,0);
2258: SNESGetDM(snes,&dm);
2259: DMGetDMSNES(dm,&sdm);
2260: if (sdm->ops->computegs) {
2261: if (b) {VecLockPush(b);}
2262: PetscStackPush("SNES user NGS");
2263: (*sdm->ops->computegs)(snes,x,b,sdm->gsctx);
2264: PetscStackPop;
2265: if (b) {VecLockPop(b);}
2266: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetNGS() before SNESComputeNGS(), likely called from SNESSolve().");
2267: PetscLogEventEnd(SNES_NGSEval,snes,x,b,0);
2268: return(0);
2269: }
2271: /*@
2272: SNESComputeJacobian - Computes the Jacobian matrix that has been set with SNESSetJacobian().
2274: Collective on SNES and Mat
2276: Input Parameters:
2277: + snes - the SNES context
2278: - x - input vector
2280: Output Parameters:
2281: + A - Jacobian matrix
2282: - B - optional preconditioning matrix
2284: Options Database Keys:
2285: + -snes_lag_preconditioner <lag>
2286: . -snes_lag_jacobian <lag>
2287: . -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
2288: . -snes_compare_explicit_draw - Compare the computed Jacobian to the finite difference Jacobian and draw the result
2289: . -snes_compare_explicit_contour - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
2290: . -snes_compare_operator - Make the comparison options above use the operator instead of the preconditioning matrix
2291: . -snes_compare_coloring - Compute the finite difference Jacobian using coloring and display norms of difference
2292: . -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
2293: . -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
2294: . -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2295: . -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
2296: . -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
2297: - -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences
2300: Notes:
2301: Most users should not need to explicitly call this routine, as it
2302: is used internally within the nonlinear solvers.
2304: Level: developer
2306: .keywords: SNES, compute, Jacobian, matrix
2308: .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
2309: @*/
2310: PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat A,Mat B)
2311: {
2313: PetscBool flag;
2314: DM dm;
2315: DMSNES sdm;
2316: KSP ksp;
2322: VecValidValues(X,2,PETSC_TRUE);
2323: SNESGetDM(snes,&dm);
2324: DMGetDMSNES(dm,&sdm);
2326: if (!sdm->ops->computejacobian) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_USER,"Must call SNESSetJacobian(), DMSNESSetJacobian(), DMDASNESSetJacobianLocal(), etc");
2328: /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
2330: if (snes->lagjacobian == -2) {
2331: snes->lagjacobian = -1;
2333: PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");
2334: } else if (snes->lagjacobian == -1) {
2335: PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");
2336: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2337: if (flag) {
2338: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2339: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2340: }
2341: return(0);
2342: } else if (snes->lagjacobian > 1 && (snes->iter + snes->jac_iter) % snes->lagjacobian) {
2343: PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);
2344: PetscObjectTypeCompare((PetscObject)A,MATMFFD,&flag);
2345: if (flag) {
2346: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2347: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2348: }
2349: return(0);
2350: }
2351: if (snes->npc && snes->npcside== PC_LEFT) {
2352: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2353: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2354: return(0);
2355: }
2357: PetscLogEventBegin(SNES_JacobianEval,snes,X,A,B);
2358: VecLockPush(X);
2359: PetscStackPush("SNES user Jacobian function");
2360: (*sdm->ops->computejacobian)(snes,X,A,B,sdm->jacobianctx);
2361: PetscStackPop;
2362: VecLockPop(X);
2363: PetscLogEventEnd(SNES_JacobianEval,snes,X,A,B);
2365: /* the next line ensures that snes->ksp exists */
2366: SNESGetKSP(snes,&ksp);
2367: if (snes->lagpreconditioner == -2) {
2368: PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");
2369: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2370: snes->lagpreconditioner = -1;
2371: } else if (snes->lagpreconditioner == -1) {
2372: PetscInfo(snes,"Reusing preconditioner because lag is -1\n");
2373: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2374: } else if (snes->lagpreconditioner > 1 && (snes->iter + snes->pre_iter) % snes->lagpreconditioner) {
2375: PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);
2376: KSPSetReusePreconditioner(snes->ksp,PETSC_TRUE);
2377: } else {
2378: PetscInfo(snes,"Rebuilding preconditioner\n");
2379: KSPSetReusePreconditioner(snes->ksp,PETSC_FALSE);
2380: }
2382: /* make sure user returned a correct Jacobian and preconditioner */
2385: {
2386: PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
2387: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_explicit",NULL,NULL,&flag);
2388: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",NULL,NULL,&flag_draw);
2389: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",NULL,NULL,&flag_contour);
2390: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_operator",NULL,NULL,&flag_operator);
2391: if (flag || flag_draw || flag_contour) {
2392: Mat Bexp_mine = NULL,Bexp,FDexp;
2393: PetscViewer vdraw,vstdout;
2394: PetscBool flg;
2395: if (flag_operator) {
2396: MatComputeExplicitOperator(A,&Bexp_mine);
2397: Bexp = Bexp_mine;
2398: } else {
2399: /* See if the preconditioning matrix can be viewed and added directly */
2400: PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");
2401: if (flg) Bexp = B;
2402: else {
2403: /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
2404: MatComputeExplicitOperator(B,&Bexp_mine);
2405: Bexp = Bexp_mine;
2406: }
2407: }
2408: MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);
2409: SNESComputeJacobianDefault(snes,X,FDexp,FDexp,NULL);
2410: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2411: if (flag_draw || flag_contour) {
2412: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2413: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2414: } else vdraw = NULL;
2415: PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator ? "Jacobian" : "preconditioning Jacobian");
2416: if (flag) {MatView(Bexp,vstdout);}
2417: if (vdraw) {MatView(Bexp,vdraw);}
2418: PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");
2419: if (flag) {MatView(FDexp,vstdout);}
2420: if (vdraw) {MatView(FDexp,vdraw);}
2421: MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);
2422: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");
2423: if (flag) {MatView(FDexp,vstdout);}
2424: if (vdraw) { /* Always use contour for the difference */
2425: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2426: MatView(FDexp,vdraw);
2427: PetscViewerPopFormat(vdraw);
2428: }
2429: if (flag_contour) {PetscViewerPopFormat(vdraw);}
2430: PetscViewerDestroy(&vdraw);
2431: MatDestroy(&Bexp_mine);
2432: MatDestroy(&FDexp);
2433: }
2434: }
2435: {
2436: PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
2437: PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
2438: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring",NULL,NULL,&flag);
2439: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring_display",NULL,NULL,&flag_display);
2440: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",NULL,NULL,&flag_draw);
2441: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",NULL,NULL,&flag_contour);
2442: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",NULL,NULL,&flag_threshold);
2443: if (flag_threshold) {
2444: PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,NULL);
2445: PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,NULL);
2446: }
2447: if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
2448: Mat Bfd;
2449: PetscViewer vdraw,vstdout;
2450: MatColoring coloring;
2451: ISColoring iscoloring;
2452: MatFDColoring matfdcoloring;
2453: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2454: void *funcctx;
2455: PetscReal norm1,norm2,normmax;
2457: MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&Bfd);
2458: MatColoringCreate(Bfd,&coloring);
2459: MatColoringSetType(coloring,MATCOLORINGSL);
2460: MatColoringSetFromOptions(coloring);
2461: MatColoringApply(coloring,&iscoloring);
2462: MatColoringDestroy(&coloring);
2463: MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);
2464: MatFDColoringSetFromOptions(matfdcoloring);
2465: MatFDColoringSetUp(Bfd,iscoloring,matfdcoloring);
2466: ISColoringDestroy(&iscoloring);
2468: /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
2469: SNESGetFunction(snes,NULL,&func,&funcctx);
2470: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))func,funcctx);
2471: PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);
2472: PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");
2473: MatFDColoringSetFromOptions(matfdcoloring);
2474: MatFDColoringApply(Bfd,matfdcoloring,X,snes);
2475: MatFDColoringDestroy(&matfdcoloring);
2477: PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&vstdout);
2478: if (flag_draw || flag_contour) {
2479: PetscViewerDrawOpen(PetscObjectComm((PetscObject)snes),0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);
2480: if (flag_contour) {PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);}
2481: } else vdraw = NULL;
2482: PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");
2483: if (flag_display) {MatView(B,vstdout);}
2484: if (vdraw) {MatView(B,vdraw);}
2485: PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");
2486: if (flag_display) {MatView(Bfd,vstdout);}
2487: if (vdraw) {MatView(Bfd,vdraw);}
2488: MatAYPX(Bfd,-1.0,B,SAME_NONZERO_PATTERN);
2489: MatNorm(Bfd,NORM_1,&norm1);
2490: MatNorm(Bfd,NORM_FROBENIUS,&norm2);
2491: MatNorm(Bfd,NORM_MAX,&normmax);
2492: PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%g normFrob=%g normmax=%g\n",(double)norm1,(double)norm2,(double)normmax);
2493: if (flag_display) {MatView(Bfd,vstdout);}
2494: if (vdraw) { /* Always use contour for the difference */
2495: PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);
2496: MatView(Bfd,vdraw);
2497: PetscViewerPopFormat(vdraw);
2498: }
2499: if (flag_contour) {PetscViewerPopFormat(vdraw);}
2501: if (flag_threshold) {
2502: PetscInt bs,rstart,rend,i;
2503: MatGetBlockSize(B,&bs);
2504: MatGetOwnershipRange(B,&rstart,&rend);
2505: for (i=rstart; i<rend; i++) {
2506: const PetscScalar *ba,*ca;
2507: const PetscInt *bj,*cj;
2508: PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
2509: PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0;
2510: MatGetRow(B,i,&bn,&bj,&ba);
2511: MatGetRow(Bfd,i,&cn,&cj,&ca);
2512: if (bn != cn) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
2513: for (j=0; j<bn; j++) {
2514: PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2515: if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
2516: maxentrycol = bj[j];
2517: maxentry = PetscRealPart(ba[j]);
2518: }
2519: if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
2520: maxdiffcol = bj[j];
2521: maxdiff = PetscRealPart(ca[j]);
2522: }
2523: if (rdiff > maxrdiff) {
2524: maxrdiffcol = bj[j];
2525: maxrdiff = rdiff;
2526: }
2527: }
2528: if (maxrdiff > 1) {
2529: 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);
2530: for (j=0; j<bn; j++) {
2531: PetscReal rdiff;
2532: rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
2533: if (rdiff > 1) {
2534: PetscViewerASCIIPrintf(vstdout," (%D,%g:%g)",bj[j],(double)PetscRealPart(ba[j]),(double)PetscRealPart(ca[j]));
2535: }
2536: }
2537: PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);
2538: }
2539: MatRestoreRow(B,i,&bn,&bj,&ba);
2540: MatRestoreRow(Bfd,i,&cn,&cj,&ca);
2541: }
2542: }
2543: PetscViewerDestroy(&vdraw);
2544: MatDestroy(&Bfd);
2545: }
2546: }
2547: return(0);
2548: }
2550: /*MC
2551: SNESJacobianFunction - Function used to convey the nonlinear Jacobian of the function to be solved by SNES
2553: Synopsis:
2554: #include "petscsnes.h"
2555: PetscErrorCode SNESJacobianFunction(SNES snes,Vec x,Mat Amat,Mat Pmat,void *ctx);
2557: + x - input vector
2558: . Amat - the matrix that defines the (approximate) Jacobian
2559: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2560: - ctx - [optional] user-defined Jacobian context
2562: Level: intermediate
2564: .seealso: SNESSetFunction(), SNESGetFunction(), SNESSetJacobian(), SNESGetJacobian()
2565: M*/
2567: /*@C
2568: SNESSetJacobian - Sets the function to compute Jacobian as well as the
2569: location to store the matrix.
2571: Logically Collective on SNES and Mat
2573: Input Parameters:
2574: + snes - the SNES context
2575: . Amat - the matrix that defines the (approximate) Jacobian
2576: . Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat.
2577: . J - Jacobian evaluation routine (if NULL then SNES retains any previously set value), see SNESJacobianFunction for details
2578: - ctx - [optional] user-defined context for private data for the
2579: Jacobian evaluation routine (may be NULL) (if NULL then SNES retains any previously set value)
2581: Notes:
2582: If the Amat matrix and Pmat matrix are different you must call MatAssemblyBegin/End() on
2583: each matrix.
2585: If you know the operator Amat has a null space you can use MatSetNullSpace() and MatSetTransposeNullSpace() to supply the null
2586: space to Amat and the KSP solvers will automatically use that null space as needed during the solution process.
2588: If using SNESComputeJacobianDefaultColor() to assemble a Jacobian, the ctx argument
2589: must be a MatFDColoring.
2591: Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian. One common
2592: example is to use the "Picard linearization" which only differentiates through the highest order parts of each term.
2594: Level: beginner
2596: .keywords: SNES, nonlinear, set, Jacobian, matrix
2598: .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESComputeJacobianDefaultColor(), MatStructure, J,
2599: SNESSetPicard(), SNESJacobianFunction
2600: @*/
2601: PetscErrorCode SNESSetJacobian(SNES snes,Mat Amat,Mat Pmat,PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*),void *ctx)
2602: {
2604: DM dm;
2612: SNESGetDM(snes,&dm);
2613: DMSNESSetJacobian(dm,J,ctx);
2614: if (Amat) {
2615: PetscObjectReference((PetscObject)Amat);
2616: MatDestroy(&snes->jacobian);
2618: snes->jacobian = Amat;
2619: }
2620: if (Pmat) {
2621: PetscObjectReference((PetscObject)Pmat);
2622: MatDestroy(&snes->jacobian_pre);
2624: snes->jacobian_pre = Pmat;
2625: }
2626: return(0);
2627: }
2629: /*@C
2630: SNESGetJacobian - Returns the Jacobian matrix and optionally the user
2631: provided context for evaluating the Jacobian.
2633: Not Collective, but Mat object will be parallel if SNES object is
2635: Input Parameter:
2636: . snes - the nonlinear solver context
2638: Output Parameters:
2639: + Amat - location to stash (approximate) Jacobian matrix (or NULL)
2640: . Pmat - location to stash matrix used to compute the preconditioner (or NULL)
2641: . J - location to put Jacobian function (or NULL), see SNESJacobianFunction for details on its calling sequence
2642: - ctx - location to stash Jacobian ctx (or NULL)
2644: Level: advanced
2646: .seealso: SNESSetJacobian(), SNESComputeJacobian(), SNESJacobianFunction, SNESGetFunction()
2647: @*/
2648: PetscErrorCode SNESGetJacobian(SNES snes,Mat *Amat,Mat *Pmat,PetscErrorCode (**J)(SNES,Vec,Mat,Mat,void*),void **ctx)
2649: {
2651: DM dm;
2652: DMSNES sdm;
2656: if (Amat) *Amat = snes->jacobian;
2657: if (Pmat) *Pmat = snes->jacobian_pre;
2658: SNESGetDM(snes,&dm);
2659: DMGetDMSNES(dm,&sdm);
2660: if (J) *J = sdm->ops->computejacobian;
2661: if (ctx) *ctx = sdm->jacobianctx;
2662: return(0);
2663: }
2665: /*@
2666: SNESSetUp - Sets up the internal data structures for the later use
2667: of a nonlinear solver.
2669: Collective on SNES
2671: Input Parameters:
2672: . snes - the SNES context
2674: Notes:
2675: For basic use of the SNES solvers the user need not explicitly call
2676: SNESSetUp(), since these actions will automatically occur during
2677: the call to SNESSolve(). However, if one wishes to control this
2678: phase separately, SNESSetUp() should be called after SNESCreate()
2679: and optional routines of the form SNESSetXXX(), but before SNESSolve().
2681: Level: advanced
2683: .keywords: SNES, nonlinear, setup
2685: .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
2686: @*/
2687: PetscErrorCode SNESSetUp(SNES snes)
2688: {
2690: DM dm;
2691: DMSNES sdm;
2692: SNESLineSearch linesearch, pclinesearch;
2693: void *lsprectx,*lspostctx;
2694: PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
2695: PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
2696: PetscErrorCode (*func)(SNES,Vec,Vec,void*);
2697: Vec f,fpc;
2698: void *funcctx;
2699: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*);
2700: void *jacctx,*appctx;
2701: Mat j,jpre;
2705: if (snes->setupcalled) return(0);
2707: if (!((PetscObject)snes)->type_name) {
2708: SNESSetType(snes,SNESNEWTONLS);
2709: }
2711: SNESGetFunction(snes,&snes->vec_func,NULL,NULL);
2713: SNESGetDM(snes,&dm);
2714: DMGetDMSNES(dm,&sdm);
2715: if (!sdm->ops->computefunction) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Function never provided to SNES object");
2716: if (!sdm->ops->computejacobian) {
2717: DMSNESSetJacobian(dm,SNESComputeJacobianDefaultColor,NULL);
2718: }
2719: if (!snes->vec_func) {
2720: DMCreateGlobalVector(dm,&snes->vec_func);
2721: }
2723: if (!snes->ksp) {
2724: SNESGetKSP(snes, &snes->ksp);
2725: }
2727: if (!snes->linesearch) {
2728: SNESGetLineSearch(snes, &snes->linesearch);
2729: }
2730: SNESLineSearchSetFunction(snes->linesearch,SNESComputeFunction);
2732: if (snes->npc && (snes->npcside== PC_LEFT)) {
2733: snes->mf = PETSC_TRUE;
2734: snes->mf_operator = PETSC_FALSE;
2735: }
2737: if (snes->npc) {
2738: /* copy the DM over */
2739: SNESGetDM(snes,&dm);
2740: SNESSetDM(snes->npc,dm);
2742: SNESGetFunction(snes,&f,&func,&funcctx);
2743: VecDuplicate(f,&fpc);
2744: SNESSetFunction(snes->npc,fpc,func,funcctx);
2745: SNESGetJacobian(snes,&j,&jpre,&jac,&jacctx);
2746: SNESSetJacobian(snes->npc,j,jpre,jac,jacctx);
2747: SNESGetApplicationContext(snes,&appctx);
2748: SNESSetApplicationContext(snes->npc,appctx);
2749: VecDestroy(&fpc);
2751: /* copy the function pointers over */
2752: PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)snes->npc);
2754: /* default to 1 iteration */
2755: SNESSetTolerances(snes->npc,0.0,0.0,0.0,1,snes->npc->max_funcs);
2756: if (snes->npcside==PC_RIGHT) {
2757: SNESSetNormSchedule(snes->npc,SNES_NORM_FINAL_ONLY);
2758: } else {
2759: SNESSetNormSchedule(snes->npc,SNES_NORM_NONE);
2760: }
2761: SNESSetFromOptions(snes->npc);
2763: /* copy the line search context over */
2764: SNESGetLineSearch(snes,&linesearch);
2765: SNESGetLineSearch(snes->npc,&pclinesearch);
2766: SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);
2767: SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);
2768: SNESLineSearchSetPreCheck(pclinesearch,precheck,lsprectx);
2769: SNESLineSearchSetPostCheck(pclinesearch,postcheck,lspostctx);
2770: PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)pclinesearch);
2771: }
2772: if (snes->mf) {
2773: SNESSetUpMatrixFree_Private(snes, snes->mf_operator, snes->mf_version);
2774: }
2775: if (snes->ops->usercompute && !snes->user) {
2776: (*snes->ops->usercompute)(snes,(void**)&snes->user);
2777: }
2779: snes->jac_iter = 0;
2780: snes->pre_iter = 0;
2782: if (snes->ops->setup) {
2783: (*snes->ops->setup)(snes);
2784: }
2786: if (snes->npc && (snes->npcside== PC_LEFT)) {
2787: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
2788: SNESGetLineSearch(snes,&linesearch);
2789: SNESLineSearchSetFunction(linesearch,SNESComputeFunctionDefaultNPC);
2790: }
2791: }
2793: snes->setupcalled = PETSC_TRUE;
2794: return(0);
2795: }
2797: /*@
2798: SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats
2800: Collective on SNES
2802: Input Parameter:
2803: . snes - iterative context obtained from SNESCreate()
2805: Level: intermediate
2807: Notes: Also calls the application context destroy routine set with SNESSetComputeApplicationContext()
2809: .keywords: SNES, destroy
2811: .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
2812: @*/
2813: PetscErrorCode SNESReset(SNES snes)
2814: {
2819: if (snes->ops->userdestroy && snes->user) {
2820: (*snes->ops->userdestroy)((void**)&snes->user);
2821: snes->user = NULL;
2822: }
2823: if (snes->npc) {
2824: SNESReset(snes->npc);
2825: }
2827: if (snes->ops->reset) {
2828: (*snes->ops->reset)(snes);
2829: }
2830: if (snes->ksp) {
2831: KSPReset(snes->ksp);
2832: }
2834: if (snes->linesearch) {
2835: SNESLineSearchReset(snes->linesearch);
2836: }
2838: VecDestroy(&snes->vec_rhs);
2839: VecDestroy(&snes->vec_sol);
2840: VecDestroy(&snes->vec_sol_update);
2841: VecDestroy(&snes->vec_func);
2842: MatDestroy(&snes->jacobian);
2843: MatDestroy(&snes->jacobian_pre);
2844: VecDestroyVecs(snes->nwork,&snes->work);
2845: VecDestroyVecs(snes->nvwork,&snes->vwork);
2847: snes->alwayscomputesfinalresidual = PETSC_FALSE;
2849: snes->nwork = snes->nvwork = 0;
2850: snes->setupcalled = PETSC_FALSE;
2851: return(0);
2852: }
2854: /*@
2855: SNESDestroy - Destroys the nonlinear solver context that was created
2856: with SNESCreate().
2858: Collective on SNES
2860: Input Parameter:
2861: . snes - the SNES context
2863: Level: beginner
2865: .keywords: SNES, nonlinear, destroy
2867: .seealso: SNESCreate(), SNESSolve()
2868: @*/
2869: PetscErrorCode SNESDestroy(SNES *snes)
2870: {
2874: if (!*snes) return(0);
2876: if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; return(0);}
2878: SNESReset((*snes));
2879: SNESDestroy(&(*snes)->npc);
2881: /* if memory was published with SAWs then destroy it */
2882: PetscObjectSAWsViewOff((PetscObject)*snes);
2883: if ((*snes)->ops->destroy) {(*((*snes))->ops->destroy)((*snes));}
2885: if ((*snes)->dm) {DMCoarsenHookRemove((*snes)->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,*snes);}
2886: DMDestroy(&(*snes)->dm);
2887: KSPDestroy(&(*snes)->ksp);
2888: SNESLineSearchDestroy(&(*snes)->linesearch);
2890: PetscFree((*snes)->kspconvctx);
2891: if ((*snes)->ops->convergeddestroy) {
2892: (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);
2893: }
2894: if ((*snes)->conv_malloc) {
2895: PetscFree((*snes)->conv_hist);
2896: PetscFree((*snes)->conv_hist_its);
2897: }
2898: SNESMonitorCancel((*snes));
2899: PetscHeaderDestroy(snes);
2900: return(0);
2901: }
2903: /* ----------- Routines to set solver parameters ---------- */
2905: /*@
2906: SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
2908: Logically Collective on SNES
2910: Input Parameters:
2911: + snes - the SNES context
2912: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2913: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
2915: Options Database Keys:
2916: . -snes_lag_preconditioner <lag>
2918: Notes:
2919: The default is 1
2920: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2921: If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
2923: Level: intermediate
2925: .keywords: SNES, nonlinear, set, convergence, tolerances
2927: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
2929: @*/
2930: PetscErrorCode SNESSetLagPreconditioner(SNES snes,PetscInt lag)
2931: {
2934: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2935: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2937: snes->lagpreconditioner = lag;
2938: return(0);
2939: }
2941: /*@
2942: SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does
2944: Logically Collective on SNES
2946: Input Parameters:
2947: + snes - the SNES context
2948: - steps - the number of refinements to do, defaults to 0
2950: Options Database Keys:
2951: . -snes_grid_sequence <steps>
2953: Level: intermediate
2955: Notes:
2956: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
2958: .keywords: SNES, nonlinear, set, convergence, tolerances
2960: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetGridSequence()
2962: @*/
2963: PetscErrorCode SNESSetGridSequence(SNES snes,PetscInt steps)
2964: {
2968: snes->gridsequence = steps;
2969: return(0);
2970: }
2972: /*@
2973: SNESGetGridSequence - gets the number of steps of grid sequencing that SNES does
2975: Logically Collective on SNES
2977: Input Parameter:
2978: . snes - the SNES context
2980: Output Parameter:
2981: . steps - the number of refinements to do, defaults to 0
2983: Options Database Keys:
2984: . -snes_grid_sequence <steps>
2986: Level: intermediate
2988: Notes:
2989: Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
2991: .keywords: SNES, nonlinear, set, convergence, tolerances
2993: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESSetGridSequence()
2995: @*/
2996: PetscErrorCode SNESGetGridSequence(SNES snes,PetscInt *steps)
2997: {
3000: *steps = snes->gridsequence;
3001: return(0);
3002: }
3004: /*@
3005: SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
3007: Not Collective
3009: Input Parameter:
3010: . snes - the SNES context
3012: Output Parameter:
3013: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3014: the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
3016: Options Database Keys:
3017: . -snes_lag_preconditioner <lag>
3019: Notes:
3020: The default is 1
3021: The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3023: Level: intermediate
3025: .keywords: SNES, nonlinear, set, convergence, tolerances
3027: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
3029: @*/
3030: PetscErrorCode SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
3031: {
3034: *lag = snes->lagpreconditioner;
3035: return(0);
3036: }
3038: /*@
3039: SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
3040: often the preconditioner is rebuilt.
3042: Logically Collective on SNES
3044: Input Parameters:
3045: + snes - the SNES context
3046: - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3047: the Jacobian is built etc. -2 means rebuild at next chance but then never again
3049: Options Database Keys:
3050: . -snes_lag_jacobian <lag>
3052: Notes:
3053: The default is 1
3054: The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3055: 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
3056: at the next Newton step but never again (unless it is reset to another value)
3058: Level: intermediate
3060: .keywords: SNES, nonlinear, set, convergence, tolerances
3062: .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
3064: @*/
3065: PetscErrorCode SNESSetLagJacobian(SNES snes,PetscInt lag)
3066: {
3069: if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
3070: if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
3072: snes->lagjacobian = lag;
3073: return(0);
3074: }
3076: /*@
3077: SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
3079: Not Collective
3081: Input Parameter:
3082: . snes - the SNES context
3084: Output Parameter:
3085: . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
3086: the Jacobian is built etc.
3088: Options Database Keys:
3089: . -snes_lag_jacobian <lag>
3091: Notes:
3092: The default is 1
3093: The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
3095: Level: intermediate
3097: .keywords: SNES, nonlinear, set, convergence, tolerances
3099: .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
3101: @*/
3102: PetscErrorCode SNESGetLagJacobian(SNES snes,PetscInt *lag)
3103: {
3106: *lag = snes->lagjacobian;
3107: return(0);
3108: }
3110: /*@
3111: SNESSetLagJacobianPersists - Set whether or not the Jacobian lagging persists through multiple solves
3113: Logically collective on SNES
3115: Input Parameter:
3116: + snes - the SNES context
3117: - flg - jacobian lagging persists if true
3119: Options Database Keys:
3120: . -snes_lag_jacobian_persists <flg>
3122: Notes: This is useful both for nonlinear preconditioning, where it's appropriate to have the Jacobian be stale by
3123: several solves, and for implicit time-stepping, where Jacobian lagging in the inner nonlinear solve over several
3124: timesteps may present huge efficiency gains.
3126: Level: developer
3128: .keywords: SNES, nonlinear, lag
3130: .seealso: SNESSetLagPreconditionerPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()
3132: @*/
3133: PetscErrorCode SNESSetLagJacobianPersists(SNES snes,PetscBool flg)
3134: {
3138: snes->lagjac_persist = flg;
3139: return(0);
3140: }
3142: /*@
3143: SNESSetLagPreconditionerPersists - Set whether or not the preconditioner lagging persists through multiple solves
3145: Logically Collective on SNES
3147: Input Parameter:
3148: + snes - the SNES context
3149: - flg - preconditioner lagging persists if true
3151: Options Database Keys:
3152: . -snes_lag_jacobian_persists <flg>
3154: Notes: This is useful both for nonlinear preconditioning, where it's appropriate to have the preconditioner be stale
3155: by several solves, and for implicit time-stepping, where preconditioner lagging in the inner nonlinear solve over
3156: several timesteps may present huge efficiency gains.
3158: Level: developer
3160: .keywords: SNES, nonlinear, lag
3162: .seealso: SNESSetLagJacobianPersists(), SNESSetLagJacobian(), SNESGetLagJacobian(), SNESGetNPC()
3164: @*/
3165: PetscErrorCode SNESSetLagPreconditionerPersists(SNES snes,PetscBool flg)
3166: {
3170: snes->lagpre_persist = flg;
3171: return(0);
3172: }
3174: /*@
3175: SNESSetForceIteration - force SNESSolve() to take at least one iteration regardless of the initial residual norm
3177: Logically Collective on SNES
3179: Input Parameters:
3180: + snes - the SNES context
3181: - force - PETSC_TRUE require at least one iteration
3183: Options Database Keys:
3184: . -snes_force_iteration <force> - Sets forcing an iteration
3186: Notes:
3187: This is used sometimes with TS to prevent TS from detecting a false steady state solution
3189: Level: intermediate
3191: .keywords: SNES, nonlinear, set, convergence, tolerances
3193: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3194: @*/
3195: PetscErrorCode SNESSetForceIteration(SNES snes,PetscBool force)
3196: {
3199: snes->forceiteration = force;
3200: return(0);
3201: }
3203: /*@
3204: SNESGetForceIteration - Whether or not to force SNESSolve() take at least one iteration regardless of the initial residual norm
3206: Logically Collective on SNES
3208: Input Parameters:
3209: . snes - the SNES context
3211: Output Parameter:
3212: . force - PETSC_TRUE requires at least one iteration.
3214: .keywords: SNES, nonlinear, set, convergence, tolerances
3216: .seealso: SNESSetForceIteration(), SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance()
3217: @*/
3218: PetscErrorCode SNESGetForceIteration(SNES snes,PetscBool *force)
3219: {
3222: *force = snes->forceiteration;
3223: return(0);
3224: }
3226: /*@
3227: SNESSetTolerances - Sets various parameters used in convergence tests.
3229: Logically Collective on SNES
3231: Input Parameters:
3232: + snes - the SNES context
3233: . abstol - absolute convergence tolerance
3234: . rtol - relative convergence tolerance
3235: . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x ||
3236: . maxit - maximum number of iterations
3237: - maxf - maximum number of function evaluations
3239: Options Database Keys:
3240: + -snes_atol <abstol> - Sets abstol
3241: . -snes_rtol <rtol> - Sets rtol
3242: . -snes_stol <stol> - Sets stol
3243: . -snes_max_it <maxit> - Sets maxit
3244: - -snes_max_funcs <maxf> - Sets maxf
3246: Notes:
3247: The default maximum number of iterations is 50.
3248: The default maximum number of function evaluations is 1000.
3250: Level: intermediate
3252: .keywords: SNES, nonlinear, set, convergence, tolerances
3254: .seealso: SNESSetTrustRegionTolerance(), SNESSetDivergenceTolerance(), SNESSetForceIteration()
3255: @*/
3256: PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
3257: {
3266: if (abstol != PETSC_DEFAULT) {
3267: if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
3268: snes->abstol = abstol;
3269: }
3270: if (rtol != PETSC_DEFAULT) {
3271: 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);
3272: snes->rtol = rtol;
3273: }
3274: if (stol != PETSC_DEFAULT) {
3275: if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
3276: snes->stol = stol;
3277: }
3278: if (maxit != PETSC_DEFAULT) {
3279: if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
3280: snes->max_its = maxit;
3281: }
3282: if (maxf != PETSC_DEFAULT) {
3283: if (maxf < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
3284: snes->max_funcs = maxf;
3285: }
3286: snes->tolerancesset = PETSC_TRUE;
3287: return(0);
3288: }
3290: /*@
3291: SNESSetDivergenceTolerance - Sets the divergence tolerance used for the SNES divergence test.
3293: Logically Collective on SNES
3295: Input Parameters:
3296: + snes - the SNES context
3297: - divtol - the divergence tolerance. Use -1 to deactivate the test.
3299: Options Database Keys:
3300: + -snes_divergence_tolerance <divtol> - Sets divtol
3302: Notes:
3303: The default divergence tolerance is 1e4.
3305: Level: intermediate
3307: .keywords: SNES, nonlinear, set, divergence, tolerance
3309: .seealso: SNESSetTolerances(), SNESGetDivergenceTolerance
3310: @*/
3311: PetscErrorCode SNESSetDivergenceTolerance(SNES snes,PetscReal divtol)
3312: {
3317: if (divtol != PETSC_DEFAULT) {
3318: snes->divtol = divtol;
3319: }
3320: else {
3321: snes->divtol = 1.0e4;
3322: }
3323: return(0);
3324: }
3326: /*@
3327: SNESGetTolerances - Gets various parameters used in convergence tests.
3329: Not Collective
3331: Input Parameters:
3332: + snes - the SNES context
3333: . atol - absolute convergence tolerance
3334: . rtol - relative convergence tolerance
3335: . stol - convergence tolerance in terms of the norm
3336: of the change in the solution between steps
3337: . maxit - maximum number of iterations
3338: - maxf - maximum number of function evaluations
3340: Notes:
3341: The user can specify NULL for any parameter that is not needed.
3343: Level: intermediate
3345: .keywords: SNES, nonlinear, get, convergence, tolerances
3347: .seealso: SNESSetTolerances()
3348: @*/
3349: PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
3350: {
3353: if (atol) *atol = snes->abstol;
3354: if (rtol) *rtol = snes->rtol;
3355: if (stol) *stol = snes->stol;
3356: if (maxit) *maxit = snes->max_its;
3357: if (maxf) *maxf = snes->max_funcs;
3358: return(0);
3359: }
3361: /*@
3362: SNESGetDivergenceTolerance - Gets divergence tolerance used in divergence test.
3364: Not Collective
3366: Input Parameters:
3367: + snes - the SNES context
3368: - divtol - divergence tolerance
3370: Level: intermediate
3372: .keywords: SNES, nonlinear, get, divergence, tolerance
3374: .seealso: SNESSetDivergenceTolerance()
3375: @*/
3376: PetscErrorCode SNESGetDivergenceTolerance(SNES snes,PetscReal *divtol)
3377: {
3380: if (divtol) *divtol = snes->divtol;
3381: return(0);
3382: }
3384: /*@
3385: SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
3387: Logically Collective on SNES
3389: Input Parameters:
3390: + snes - the SNES context
3391: - tol - tolerance
3393: Options Database Key:
3394: . -snes_trtol <tol> - Sets tol
3396: Level: intermediate
3398: .keywords: SNES, nonlinear, set, trust region, tolerance
3400: .seealso: SNESSetTolerances()
3401: @*/
3402: PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
3403: {
3407: snes->deltatol = tol;
3408: return(0);
3409: }
3411: /*
3412: Duplicate the lg monitors for SNES from KSP; for some reason with
3413: dynamic libraries things don't work under Sun4 if we just use
3414: macros instead of functions
3415: */
3416: PetscErrorCode SNESMonitorLGResidualNorm(SNES snes,PetscInt it,PetscReal norm,void *ctx)
3417: {
3422: KSPMonitorLGResidualNorm((KSP)snes,it,norm,ctx);
3423: return(0);
3424: }
3426: PetscErrorCode SNESMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)
3427: {
3431: KSPMonitorLGResidualNormCreate(comm,host,label,x,y,m,n,lgctx);
3432: return(0);
3433: }
3435: PETSC_INTERN PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
3437: PetscErrorCode SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
3438: {
3439: PetscDrawLG lg;
3440: PetscErrorCode ierr;
3441: PetscReal x,y,per;
3442: PetscViewer v = (PetscViewer)monctx;
3443: static PetscReal prev; /* should be in the context */
3444: PetscDraw draw;
3448: PetscViewerDrawGetDrawLG(v,0,&lg);
3449: if (!n) {PetscDrawLGReset(lg);}
3450: PetscDrawLGGetDraw(lg,&draw);
3451: PetscDrawSetTitle(draw,"Residual norm");
3452: x = (PetscReal)n;
3453: if (rnorm > 0.0) y = PetscLog10Real(rnorm);
3454: else y = -15.0;
3455: PetscDrawLGAddPoint(lg,&x,&y);
3456: if (n < 20 || !(n % 5) || snes->reason) {
3457: PetscDrawLGDraw(lg);
3458: PetscDrawLGSave(lg);
3459: }
3461: PetscViewerDrawGetDrawLG(v,1,&lg);
3462: if (!n) {PetscDrawLGReset(lg);}
3463: PetscDrawLGGetDraw(lg,&draw);
3464: PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
3465: SNESMonitorRange_Private(snes,n,&per);
3466: x = (PetscReal)n;
3467: y = 100.0*per;
3468: PetscDrawLGAddPoint(lg,&x,&y);
3469: if (n < 20 || !(n % 5) || snes->reason) {
3470: PetscDrawLGDraw(lg);
3471: PetscDrawLGSave(lg);
3472: }
3474: PetscViewerDrawGetDrawLG(v,2,&lg);
3475: if (!n) {prev = rnorm;PetscDrawLGReset(lg);}
3476: PetscDrawLGGetDraw(lg,&draw);
3477: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");
3478: x = (PetscReal)n;
3479: y = (prev - rnorm)/prev;
3480: PetscDrawLGAddPoint(lg,&x,&y);
3481: if (n < 20 || !(n % 5) || snes->reason) {
3482: PetscDrawLGDraw(lg);
3483: PetscDrawLGSave(lg);
3484: }
3486: PetscViewerDrawGetDrawLG(v,3,&lg);
3487: if (!n) {PetscDrawLGReset(lg);}
3488: PetscDrawLGGetDraw(lg,&draw);
3489: PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
3490: x = (PetscReal)n;
3491: y = (prev - rnorm)/(prev*per);
3492: if (n > 2) { /*skip initial crazy value */
3493: PetscDrawLGAddPoint(lg,&x,&y);
3494: }
3495: if (n < 20 || !(n % 5) || snes->reason) {
3496: PetscDrawLGDraw(lg);
3497: PetscDrawLGSave(lg);
3498: }
3499: prev = rnorm;
3500: return(0);
3501: }
3503: /*@
3504: SNESMonitor - runs the user provided monitor routines, if they exist
3506: Collective on SNES
3508: Input Parameters:
3509: + snes - nonlinear solver context obtained from SNESCreate()
3510: . iter - iteration number
3511: - rnorm - relative norm of the residual
3513: Notes:
3514: This routine is called by the SNES implementations.
3515: It does not typically need to be called by the user.
3517: Level: developer
3519: .seealso: SNESMonitorSet()
3520: @*/
3521: PetscErrorCode SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
3522: {
3524: PetscInt i,n = snes->numbermonitors;
3527: VecLockPush(snes->vec_sol);
3528: for (i=0; i<n; i++) {
3529: (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);
3530: }
3531: VecLockPop(snes->vec_sol);
3532: return(0);
3533: }
3535: /* ------------ Routines to set performance monitoring options ----------- */
3537: /*MC
3538: SNESMonitorFunction - functional form passed to SNESMonitorSet() to monitor convergence of nonlinear solver
3540: Synopsis:
3541: #include <petscsnes.h>
3542: $ PetscErrorCode SNESMonitorFunction(SNES snes,PetscInt its, PetscReal norm,void *mctx)
3544: + snes - the SNES context
3545: . its - iteration number
3546: . norm - 2-norm function value (may be estimated)
3547: - mctx - [optional] monitoring context
3549: Level: advanced
3551: .seealso: SNESMonitorSet(), SNESMonitorGet()
3552: M*/
3554: /*@C
3555: SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
3556: iteration of the nonlinear solver to display the iteration's
3557: progress.
3559: Logically Collective on SNES
3561: Input Parameters:
3562: + snes - the SNES context
3563: . f - the monitor function, see SNESMonitorFunction for the calling sequence
3564: . mctx - [optional] user-defined context for private data for the
3565: monitor routine (use NULL if no context is desired)
3566: - monitordestroy - [optional] routine that frees monitor context
3567: (may be NULL)
3569: Options Database Keys:
3570: + -snes_monitor - sets SNESMonitorDefault()
3571: . -snes_monitor_lg_residualnorm - sets line graph monitor,
3572: uses SNESMonitorLGCreate()
3573: - -snes_monitor_cancel - cancels all monitors that have
3574: been hardwired into a code by
3575: calls to SNESMonitorSet(), but
3576: does not cancel those set via
3577: the options database.
3579: Notes:
3580: Several different monitoring routines may be set by calling
3581: SNESMonitorSet() multiple times; all will be called in the
3582: order in which they were set.
3584: Fortran notes: Only a single monitor function can be set for each SNES object
3586: Level: intermediate
3588: .keywords: SNES, nonlinear, set, monitor
3590: .seealso: SNESMonitorDefault(), SNESMonitorCancel(), SNESMonitorFunction
3591: @*/
3592: PetscErrorCode SNESMonitorSet(SNES snes,PetscErrorCode (*f)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
3593: {
3594: PetscInt i;
3596: PetscBool identical;
3600: for (i=0; i<snes->numbermonitors;i++) {
3601: PetscMonitorCompare((PetscErrorCode (*)(void))f,mctx,monitordestroy,(PetscErrorCode (*)(void))snes->monitor[i],snes->monitorcontext[i],snes->monitordestroy[i],&identical);
3602: if (identical) return(0);
3603: }
3604: if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
3605: snes->monitor[snes->numbermonitors] = f;
3606: snes->monitordestroy[snes->numbermonitors] = monitordestroy;
3607: snes->monitorcontext[snes->numbermonitors++] = (void*)mctx;
3608: return(0);
3609: }
3611: /*@
3612: SNESMonitorCancel - Clears all the monitor functions for a SNES object.
3614: Logically Collective on SNES
3616: Input Parameters:
3617: . snes - the SNES context
3619: Options Database Key:
3620: . -snes_monitor_cancel - cancels all monitors that have been hardwired
3621: into a code by calls to SNESMonitorSet(), but does not cancel those
3622: set via the options database
3624: Notes:
3625: There is no way to clear one specific monitor from a SNES object.
3627: Level: intermediate
3629: .keywords: SNES, nonlinear, set, monitor
3631: .seealso: SNESMonitorDefault(), SNESMonitorSet()
3632: @*/
3633: PetscErrorCode SNESMonitorCancel(SNES snes)
3634: {
3636: PetscInt i;
3640: for (i=0; i<snes->numbermonitors; i++) {
3641: if (snes->monitordestroy[i]) {
3642: (*snes->monitordestroy[i])(&snes->monitorcontext[i]);
3643: }
3644: }
3645: snes->numbermonitors = 0;
3646: return(0);
3647: }
3649: /*MC
3650: SNESConvergenceTestFunction - functional form used for testing of convergence of nonlinear solver
3652: Synopsis:
3653: #include <petscsnes.h>
3654: $ PetscErrorCode SNESConvergenceTest(SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
3656: + snes - the SNES context
3657: . it - current iteration (0 is the first and is before any Newton step)
3658: . cctx - [optional] convergence context
3659: . reason - reason for convergence/divergence
3660: . xnorm - 2-norm of current iterate
3661: . gnorm - 2-norm of current step
3662: - f - 2-norm of function
3664: Level: intermediate
3666: .seealso: SNESSetConvergenceTest(), SNESGetConvergenceTest()
3667: M*/
3669: /*@C
3670: SNESSetConvergenceTest - Sets the function that is to be used
3671: to test for convergence of the nonlinear iterative solution.
3673: Logically Collective on SNES
3675: Input Parameters:
3676: + snes - the SNES context
3677: . SNESConvergenceTestFunction - routine to test for convergence
3678: . cctx - [optional] context for private data for the convergence routine (may be NULL)
3679: - destroy - [optional] destructor for the context (may be NULL; NULL_FUNCTION in Fortran)
3681: Level: advanced
3683: .keywords: SNES, nonlinear, set, convergence, test
3685: .seealso: SNESConvergedDefault(), SNESConvergedSkip(), SNESConvergenceTestFunction
3686: @*/
3687: PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*SNESConvergenceTestFunction)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
3688: {
3693: if (!SNESConvergenceTestFunction) SNESConvergenceTestFunction = SNESConvergedSkip;
3694: if (snes->ops->convergeddestroy) {
3695: (*snes->ops->convergeddestroy)(snes->cnvP);
3696: }
3697: snes->ops->converged = SNESConvergenceTestFunction;
3698: snes->ops->convergeddestroy = destroy;
3699: snes->cnvP = cctx;
3700: return(0);
3701: }
3703: /*@
3704: SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
3706: Not Collective
3708: Input Parameter:
3709: . snes - the SNES context
3711: Output Parameter:
3712: . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
3713: manual pages for the individual convergence tests for complete lists
3715: Options Database:
3716: . -snes_converged_reason - prints the reason to standard out
3718: Level: intermediate
3720: Notes: Should only be called after the call the SNESSolve() is complete, if it is called earlier it returns the value SNES__CONVERGED_ITERATING.
3722: .keywords: SNES, nonlinear, set, convergence, test
3724: .seealso: SNESSetConvergenceTest(), SNESSetConvergedReason(), SNESConvergedReason
3725: @*/
3726: PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
3727: {
3731: *reason = snes->reason;
3732: return(0);
3733: }
3735: /*@
3736: SNESSetConvergedReason - Sets the reason the SNES iteration was stopped.
3738: Not Collective
3740: Input Parameters:
3741: + snes - the SNES context
3742: - reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
3743: manual pages for the individual convergence tests for complete lists
3745: Level: intermediate
3747: .keywords: SNES, nonlinear, set, convergence, test
3748: .seealso: SNESGetConvergedReason(), SNESSetConvergenceTest(), SNESConvergedReason
3749: @*/
3750: PetscErrorCode SNESSetConvergedReason(SNES snes,SNESConvergedReason reason)
3751: {
3754: snes->reason = reason;
3755: return(0);
3756: }
3758: /*@
3759: SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
3761: Logically Collective on SNES
3763: Input Parameters:
3764: + snes - iterative context obtained from SNESCreate()
3765: . a - array to hold history, this array will contain the function norms computed at each step
3766: . its - integer array holds the number of linear iterations for each solve.
3767: . na - size of a and its
3768: - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
3769: else it continues storing new values for new nonlinear solves after the old ones
3771: Notes:
3772: If 'a' and 'its' are NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
3773: default array of length 10000 is allocated.
3775: This routine is useful, e.g., when running a code for purposes
3776: of accurate performance monitoring, when no I/O should be done
3777: during the section of code that is being timed.
3779: Level: intermediate
3781: .keywords: SNES, set, convergence, history
3783: .seealso: SNESGetConvergenceHistory()
3785: @*/
3786: PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool reset)
3787: {
3794: if (!a) {
3795: if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
3796: PetscCalloc1(na,&a);
3797: PetscCalloc1(na,&its);
3799: snes->conv_malloc = PETSC_TRUE;
3800: }
3801: snes->conv_hist = a;
3802: snes->conv_hist_its = its;
3803: snes->conv_hist_max = na;
3804: snes->conv_hist_len = 0;
3805: snes->conv_hist_reset = reset;
3806: return(0);
3807: }
3809: #if defined(PETSC_HAVE_MATLAB_ENGINE)
3810: #include <engine.h> /* MATLAB include file */
3811: #include <mex.h> /* MATLAB include file */
3813: PETSC_EXTERN mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
3814: {
3815: mxArray *mat;
3816: PetscInt i;
3817: PetscReal *ar;
3820: mat = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
3821: ar = (PetscReal*) mxGetData(mat);
3822: for (i=0; i<snes->conv_hist_len; i++) ar[i] = snes->conv_hist[i];
3823: PetscFunctionReturn(mat);
3824: }
3825: #endif
3827: /*@C
3828: SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
3830: Not Collective
3832: Input Parameter:
3833: . snes - iterative context obtained from SNESCreate()
3835: Output Parameters:
3836: . a - array to hold history
3837: . its - integer array holds the number of linear iterations (or
3838: negative if not converged) for each solve.
3839: - na - size of a and its
3841: Notes:
3842: The calling sequence for this routine in Fortran is
3843: $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
3845: This routine is useful, e.g., when running a code for purposes
3846: of accurate performance monitoring, when no I/O should be done
3847: during the section of code that is being timed.
3849: Level: intermediate
3851: .keywords: SNES, get, convergence, history
3853: .seealso: SNESSetConvergencHistory()
3855: @*/
3856: PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
3857: {
3860: if (a) *a = snes->conv_hist;
3861: if (its) *its = snes->conv_hist_its;
3862: if (na) *na = snes->conv_hist_len;
3863: return(0);
3864: }
3866: /*@C
3867: SNESSetUpdate - Sets the general-purpose update function called
3868: at the beginning of every iteration of the nonlinear solve. Specifically
3869: it is called just before the Jacobian is "evaluated".
3871: Logically Collective on SNES
3873: Input Parameters:
3874: . snes - The nonlinear solver context
3875: . func - The function
3877: Calling sequence of func:
3878: . func (SNES snes, PetscInt step);
3880: . step - The current step of the iteration
3882: Level: advanced
3884: 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()
3885: This is not used by most users.
3887: .keywords: SNES, update
3889: .seealso SNESSetJacobian(), SNESSolve()
3890: @*/
3891: PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
3892: {
3895: snes->ops->update = func;
3896: return(0);
3897: }
3899: /*
3900: SNESScaleStep_Private - Scales a step so that its length is less than the
3901: positive parameter delta.
3903: Input Parameters:
3904: + snes - the SNES context
3905: . y - approximate solution of linear system
3906: . fnorm - 2-norm of current function
3907: - delta - trust region size
3909: Output Parameters:
3910: + gpnorm - predicted function norm at the new point, assuming local
3911: linearization. The value is zero if the step lies within the trust
3912: region, and exceeds zero otherwise.
3913: - ynorm - 2-norm of the step
3915: Note:
3916: For non-trust region methods such as SNESNEWTONLS, the parameter delta
3917: is set to be the maximum allowable step size.
3919: .keywords: SNES, nonlinear, scale, step
3920: */
3921: PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
3922: {
3923: PetscReal nrm;
3924: PetscScalar cnorm;
3932: VecNorm(y,NORM_2,&nrm);
3933: if (nrm > *delta) {
3934: nrm = *delta/nrm;
3935: *gpnorm = (1.0 - nrm)*(*fnorm);
3936: cnorm = nrm;
3937: VecScale(y,cnorm);
3938: *ynorm = *delta;
3939: } else {
3940: *gpnorm = 0.0;
3941: *ynorm = nrm;
3942: }
3943: return(0);
3944: }
3946: /*@
3947: SNESReasonView - Displays the reason a SNES solve converged or diverged to a viewer
3949: Collective on SNES
3951: Parameter:
3952: + snes - iterative context obtained from SNESCreate()
3953: - viewer - the viewer to display the reason
3956: Options Database Keys:
3957: . -snes_converged_reason - print reason for converged or diverged, also prints number of iterations
3959: Level: beginner
3961: .keywords: SNES, solve, linear system
3963: .seealso: SNESCreate(), SNESSetUp(), SNESDestroy(), SNESSetTolerances(), SNESConvergedDefault()
3965: @*/
3966: PetscErrorCode SNESReasonView(SNES snes,PetscViewer viewer)
3967: {
3969: PetscBool isAscii;
3972: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
3973: if (isAscii) {
3974: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
3975: if (snes->reason > 0) {
3976: if (((PetscObject) snes)->prefix) {
3977: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve converged due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3978: } else {
3979: PetscViewerASCIIPrintf(viewer,"Nonlinear solve converged due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3980: }
3981: } else {
3982: if (((PetscObject) snes)->prefix) {
3983: PetscViewerASCIIPrintf(viewer,"Nonlinear %s solve did not converge due to %s iterations %D\n",((PetscObject) snes)->prefix,SNESConvergedReasons[snes->reason],snes->iter);
3984: } else {
3985: PetscViewerASCIIPrintf(viewer,"Nonlinear solve did not converge due to %s iterations %D\n",SNESConvergedReasons[snes->reason],snes->iter);
3986: }
3987: }
3988: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
3989: }
3990: return(0);
3991: }
3993: /*@C
3994: SNESReasonViewFromOptions - Processes command line options to determine if/how a SNESReason is to be viewed.
3996: Collective on SNES
3998: Input Parameters:
3999: . snes - the SNES object
4001: Level: intermediate
4003: @*/
4004: PetscErrorCode SNESReasonViewFromOptions(SNES snes)
4005: {
4006: PetscErrorCode ierr;
4007: PetscViewer viewer;
4008: PetscBool flg;
4009: static PetscBool incall = PETSC_FALSE;
4010: PetscViewerFormat format;
4013: if (incall) return(0);
4014: incall = PETSC_TRUE;
4015: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_converged_reason",&viewer,&format,&flg);
4016: if (flg) {
4017: PetscViewerPushFormat(viewer,format);
4018: SNESReasonView(snes,viewer);
4019: PetscViewerPopFormat(viewer);
4020: PetscViewerDestroy(&viewer);
4021: }
4022: incall = PETSC_FALSE;
4023: return(0);
4024: }
4026: /*@C
4027: SNESSolve - Solves a nonlinear system F(x) = b.
4028: Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
4030: Collective on SNES
4032: Input Parameters:
4033: + snes - the SNES context
4034: . b - the constant part of the equation F(x) = b, or NULL to use zero.
4035: - x - the solution vector.
4037: Notes:
4038: The user should initialize the vector,x, with the initial guess
4039: for the nonlinear solve prior to calling SNESSolve. In particular,
4040: to employ an initial guess of zero, the user should explicitly set
4041: this vector to zero by calling VecSet().
4043: Level: beginner
4045: .keywords: SNES, nonlinear, solve
4047: .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
4048: @*/
4049: PetscErrorCode SNESSolve(SNES snes,Vec b,Vec x)
4050: {
4051: PetscErrorCode ierr;
4052: PetscBool flg;
4053: PetscInt grid;
4054: Vec xcreated = NULL;
4055: DM dm;
4064: {
4065: PetscViewer viewer;
4066: PetscViewerFormat format;
4067: PetscBool flg;
4068: static PetscBool incall = PETSC_FALSE;
4070: if (!incall) {
4071: PetscOptionsGetViewer(PetscObjectComm((PetscObject) snes), ((PetscObject) snes)->prefix, "-snes_convergence_estimate", &viewer, &format, &flg);
4072: if (flg) {
4073: PetscConvEst conv;
4074: PetscReal alpha; /* Convergence rate of the solution error in the L_2 norm */
4076: incall = PETSC_TRUE;
4077: PetscConvEstCreate(PetscObjectComm((PetscObject) snes), &conv);
4078: PetscConvEstSetSolver(conv, snes);
4079: PetscConvEstSetFromOptions(conv);
4080: PetscConvEstSetUp(conv);
4081: PetscConvEstGetConvRate(conv, &alpha);
4082: PetscViewerPushFormat(viewer, format);
4083: PetscConvEstRateView(conv, alpha, viewer);
4084: PetscViewerPopFormat(viewer);
4085: PetscViewerDestroy(&viewer);
4086: PetscConvEstDestroy(&conv);
4087: incall = PETSC_FALSE;
4088: }
4089: }
4090: }
4091: if (!x) {
4092: SNESGetDM(snes,&dm);
4093: DMCreateGlobalVector(dm,&xcreated);
4094: x = xcreated;
4095: }
4096: SNESViewFromOptions(snes,NULL,"-snes_view_pre");
4098: for (grid=0; grid<snes->gridsequence; grid++) {PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));}
4099: for (grid=0; grid<snes->gridsequence+1; grid++) {
4101: /* set solution vector */
4102: if (!grid) {PetscObjectReference((PetscObject)x);}
4103: VecDestroy(&snes->vec_sol);
4104: snes->vec_sol = x;
4105: SNESGetDM(snes,&dm);
4107: /* set affine vector if provided */
4108: if (b) { PetscObjectReference((PetscObject)b); }
4109: VecDestroy(&snes->vec_rhs);
4110: snes->vec_rhs = b;
4112: if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
4113: if (snes->vec_rhs == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
4114: if (!snes->vec_sol_update /* && snes->vec_sol */) {
4115: VecDuplicate(snes->vec_sol,&snes->vec_sol_update);
4116: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->vec_sol_update);
4117: }
4118: DMShellSetGlobalVector(dm,snes->vec_sol);
4119: SNESSetUp(snes);
4121: if (!grid) {
4122: if (snes->ops->computeinitialguess) {
4123: (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);
4124: }
4125: }
4127: if (snes->conv_hist_reset) snes->conv_hist_len = 0;
4128: if (snes->counters_reset) {snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;}
4130: PetscLogEventBegin(SNES_Solve,snes,0,0,0);
4131: (*snes->ops->solve)(snes);
4132: PetscLogEventEnd(SNES_Solve,snes,0,0,0);
4133: if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
4134: snes->domainerror = PETSC_FALSE; /* clear the flag if it has been set */
4136: if (snes->lagjac_persist) snes->jac_iter += snes->iter;
4137: if (snes->lagpre_persist) snes->pre_iter += snes->iter;
4139: PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"-snes_test_local_min",NULL,NULL,&flg);
4140: if (flg && !PetscPreLoadingOn) { SNESTestLocalMin(snes); }
4141: SNESReasonViewFromOptions(snes);
4143: if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
4144: if (snes->reason < 0) break;
4145: if (grid < snes->gridsequence) {
4146: DM fine;
4147: Vec xnew;
4148: Mat interp;
4150: DMRefine(snes->dm,PetscObjectComm((PetscObject)snes),&fine);
4151: if (!fine) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"DMRefine() did not perform any refinement, cannot continue grid sequencing");
4152: DMCreateInterpolation(snes->dm,fine,&interp,NULL);
4153: DMCreateGlobalVector(fine,&xnew);
4154: MatInterpolate(interp,x,xnew);
4155: DMInterpolate(snes->dm,interp,fine);
4156: MatDestroy(&interp);
4157: x = xnew;
4159: SNESReset(snes);
4160: SNESSetDM(snes,fine);
4161: DMDestroy(&fine);
4162: PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)));
4163: }
4164: }
4165: SNESViewFromOptions(snes,NULL,"-snes_view");
4166: VecViewFromOptions(snes->vec_sol,(PetscObject)snes,"-snes_view_solution");
4168: VecDestroy(&xcreated);
4169: PetscObjectSAWsBlock((PetscObject)snes);
4170: return(0);
4171: }
4173: /* --------- Internal routines for SNES Package --------- */
4175: /*@C
4176: SNESSetType - Sets the method for the nonlinear solver.
4178: Collective on SNES
4180: Input Parameters:
4181: + snes - the SNES context
4182: - type - a known method
4184: Options Database Key:
4185: . -snes_type <type> - Sets the method; use -help for a list
4186: of available methods (for instance, newtonls or newtontr)
4188: Notes:
4189: See "petsc/include/petscsnes.h" for available methods (for instance)
4190: + SNESNEWTONLS - Newton's method with line search
4191: (systems of nonlinear equations)
4192: . SNESNEWTONTR - Newton's method with trust region
4193: (systems of nonlinear equations)
4195: Normally, it is best to use the SNESSetFromOptions() command and then
4196: set the SNES solver type from the options database rather than by using
4197: this routine. Using the options database provides the user with
4198: maximum flexibility in evaluating the many nonlinear solvers.
4199: The SNESSetType() routine is provided for those situations where it
4200: is necessary to set the nonlinear solver independently of the command
4201: line or options database. This might be the case, for example, when
4202: the choice of solver changes during the execution of the program,
4203: and the user's application is taking responsibility for choosing the
4204: appropriate method.
4206: Developer Notes: SNESRegister() adds a constructor for a new SNESType to SNESList, SNESSetType() locates
4207: the constructor in that list and calls it to create the spexific object.
4209: Level: intermediate
4211: .keywords: SNES, set, type
4213: .seealso: SNESType, SNESCreate(), SNESDestroy(), SNESGetType(), SNESSetFromOptions()
4215: @*/
4216: PetscErrorCode SNESSetType(SNES snes,SNESType type)
4217: {
4218: PetscErrorCode ierr,(*r)(SNES);
4219: PetscBool match;
4225: PetscObjectTypeCompare((PetscObject)snes,type,&match);
4226: if (match) return(0);
4228: PetscFunctionListFind(SNESList,type,&r);
4229: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
4230: /* Destroy the previous private SNES context */
4231: if (snes->ops->destroy) {
4232: (*(snes)->ops->destroy)(snes);
4233: snes->ops->destroy = NULL;
4234: }
4235: /* Reinitialize function pointers in SNESOps structure */
4236: snes->ops->setup = 0;
4237: snes->ops->solve = 0;
4238: snes->ops->view = 0;
4239: snes->ops->setfromoptions = 0;
4240: snes->ops->destroy = 0;
4241: /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
4242: snes->setupcalled = PETSC_FALSE;
4244: PetscObjectChangeTypeName((PetscObject)snes,type);
4245: (*r)(snes);
4246: return(0);
4247: }
4249: /*@C
4250: SNESGetType - Gets the SNES method type and name (as a string).
4252: Not Collective
4254: Input Parameter:
4255: . snes - nonlinear solver context
4257: Output Parameter:
4258: . type - SNES method (a character string)
4260: Level: intermediate
4262: .keywords: SNES, nonlinear, get, type, name
4263: @*/
4264: PetscErrorCode SNESGetType(SNES snes,SNESType *type)
4265: {
4269: *type = ((PetscObject)snes)->type_name;
4270: return(0);
4271: }
4273: /*@
4274: SNESSetSolution - Sets the solution vector for use by the SNES routines.
4276: Logically Collective on SNES and Vec
4278: Input Parameters:
4279: + snes - the SNES context obtained from SNESCreate()
4280: - u - the solution vector
4282: Level: beginner
4284: .keywords: SNES, set, solution
4285: @*/
4286: PetscErrorCode SNESSetSolution(SNES snes, Vec u)
4287: {
4288: DM dm;
4294: PetscObjectReference((PetscObject) u);
4295: VecDestroy(&snes->vec_sol);
4297: snes->vec_sol = u;
4299: SNESGetDM(snes, &dm);
4300: DMShellSetGlobalVector(dm, u);
4301: return(0);
4302: }
4304: /*@
4305: SNESGetSolution - Returns the vector where the approximate solution is
4306: stored. This is the fine grid solution when using SNESSetGridSequence().
4308: Not Collective, but Vec is parallel if SNES is parallel
4310: Input Parameter:
4311: . snes - the SNES context
4313: Output Parameter:
4314: . x - the solution
4316: Level: intermediate
4318: .keywords: SNES, nonlinear, get, solution
4320: .seealso: SNESGetSolutionUpdate(), SNESGetFunction()
4321: @*/
4322: PetscErrorCode SNESGetSolution(SNES snes,Vec *x)
4323: {
4327: *x = snes->vec_sol;
4328: return(0);
4329: }
4331: /*@
4332: SNESGetSolutionUpdate - Returns the vector where the solution update is
4333: stored.
4335: Not Collective, but Vec is parallel if SNES is parallel
4337: Input Parameter:
4338: . snes - the SNES context
4340: Output Parameter:
4341: . x - the solution update
4343: Level: advanced
4345: .keywords: SNES, nonlinear, get, solution, update
4347: .seealso: SNESGetSolution(), SNESGetFunction()
4348: @*/
4349: PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x)
4350: {
4354: *x = snes->vec_sol_update;
4355: return(0);
4356: }
4358: /*@C
4359: SNESGetFunction - Returns the vector where the function is stored.
4361: Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
4363: Input Parameter:
4364: . snes - the SNES context
4366: Output Parameter:
4367: + r - the vector that is used to store residuals (or NULL if you don't want it)
4368: . f - the function (or NULL if you don't want it); see SNESFunction for calling sequence details
4369: - ctx - the function context (or NULL if you don't want it)
4371: Level: advanced
4373: .keywords: SNES, nonlinear, get, function
4375: .seealso: SNESSetFunction(), SNESGetSolution(), SNESFunction
4376: @*/
4377: PetscErrorCode SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**f)(SNES,Vec,Vec,void*),void **ctx)
4378: {
4380: DM dm;
4384: if (r) {
4385: if (!snes->vec_func) {
4386: if (snes->vec_rhs) {
4387: VecDuplicate(snes->vec_rhs,&snes->vec_func);
4388: } else if (snes->vec_sol) {
4389: VecDuplicate(snes->vec_sol,&snes->vec_func);
4390: } else if (snes->dm) {
4391: DMCreateGlobalVector(snes->dm,&snes->vec_func);
4392: }
4393: }
4394: *r = snes->vec_func;
4395: }
4396: SNESGetDM(snes,&dm);
4397: DMSNESGetFunction(dm,f,ctx);
4398: return(0);
4399: }
4401: /*@C
4402: SNESGetNGS - Returns the NGS function and context.
4404: Input Parameter:
4405: . snes - the SNES context
4407: Output Parameter:
4408: + f - the function (or NULL) see SNESNGSFunction for details
4409: - ctx - the function context (or NULL)
4411: Level: advanced
4413: .keywords: SNES, nonlinear, get, function
4415: .seealso: SNESSetNGS(), SNESGetFunction()
4416: @*/
4418: PetscErrorCode SNESGetNGS (SNES snes, PetscErrorCode (**f)(SNES, Vec, Vec, void*), void ** ctx)
4419: {
4421: DM dm;
4425: SNESGetDM(snes,&dm);
4426: DMSNESGetNGS(dm,f,ctx);
4427: return(0);
4428: }
4430: /*@C
4431: SNESSetOptionsPrefix - Sets the prefix used for searching for all
4432: SNES options in the database.
4434: Logically Collective on SNES
4436: Input Parameter:
4437: + snes - the SNES context
4438: - prefix - the prefix to prepend to all option names
4440: Notes:
4441: A hyphen (-) must NOT be given at the beginning of the prefix name.
4442: The first character of all runtime options is AUTOMATICALLY the hyphen.
4444: Level: advanced
4446: .keywords: SNES, set, options, prefix, database
4448: .seealso: SNESSetFromOptions()
4449: @*/
4450: PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[])
4451: {
4456: PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);
4457: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4458: if (snes->linesearch) {
4459: SNESGetLineSearch(snes,&snes->linesearch);
4460: PetscObjectSetOptionsPrefix((PetscObject)snes->linesearch,prefix);
4461: }
4462: KSPSetOptionsPrefix(snes->ksp,prefix);
4463: return(0);
4464: }
4466: /*@C
4467: SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
4468: SNES options in the database.
4470: Logically Collective on SNES
4472: Input Parameters:
4473: + snes - the SNES context
4474: - prefix - the prefix to prepend to all option names
4476: Notes:
4477: A hyphen (-) must NOT be given at the beginning of the prefix name.
4478: The first character of all runtime options is AUTOMATICALLY the hyphen.
4480: Level: advanced
4482: .keywords: SNES, append, options, prefix, database
4484: .seealso: SNESGetOptionsPrefix()
4485: @*/
4486: PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[])
4487: {
4492: PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);
4493: if (!snes->ksp) {SNESGetKSP(snes,&snes->ksp);}
4494: if (snes->linesearch) {
4495: SNESGetLineSearch(snes,&snes->linesearch);
4496: PetscObjectAppendOptionsPrefix((PetscObject)snes->linesearch,prefix);
4497: }
4498: KSPAppendOptionsPrefix(snes->ksp,prefix);
4499: return(0);
4500: }
4502: /*@C
4503: SNESGetOptionsPrefix - Sets the prefix used for searching for all
4504: SNES options in the database.
4506: Not Collective
4508: Input Parameter:
4509: . snes - the SNES context
4511: Output Parameter:
4512: . prefix - pointer to the prefix string used
4514: Notes: On the fortran side, the user should pass in a string 'prefix' of
4515: sufficient length to hold the prefix.
4517: Level: advanced
4519: .keywords: SNES, get, options, prefix, database
4521: .seealso: SNESAppendOptionsPrefix()
4522: @*/
4523: PetscErrorCode SNESGetOptionsPrefix(SNES snes,const char *prefix[])
4524: {
4529: PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);
4530: return(0);
4531: }
4534: /*@C
4535: SNESRegister - Adds a method to the nonlinear solver package.
4537: Not collective
4539: Input Parameters:
4540: + name_solver - name of a new user-defined solver
4541: - routine_create - routine to create method context
4543: Notes:
4544: SNESRegister() may be called multiple times to add several user-defined solvers.
4546: Sample usage:
4547: .vb
4548: SNESRegister("my_solver",MySolverCreate);
4549: .ve
4551: Then, your solver can be chosen with the procedural interface via
4552: $ SNESSetType(snes,"my_solver")
4553: or at runtime via the option
4554: $ -snes_type my_solver
4556: Level: advanced
4558: Note: If your function is not being put into a shared library then use SNESRegister() instead
4560: .keywords: SNES, nonlinear, register
4562: .seealso: SNESRegisterAll(), SNESRegisterDestroy()
4564: Level: advanced
4565: @*/
4566: PetscErrorCode SNESRegister(const char sname[],PetscErrorCode (*function)(SNES))
4567: {
4571: PetscFunctionListAdd(&SNESList,sname,function);
4572: return(0);
4573: }
4575: PetscErrorCode SNESTestLocalMin(SNES snes)
4576: {
4578: PetscInt N,i,j;
4579: Vec u,uh,fh;
4580: PetscScalar value;
4581: PetscReal norm;
4584: SNESGetSolution(snes,&u);
4585: VecDuplicate(u,&uh);
4586: VecDuplicate(u,&fh);
4588: /* currently only works for sequential */
4589: PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
4590: VecGetSize(u,&N);
4591: for (i=0; i<N; i++) {
4592: VecCopy(u,uh);
4593: PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);
4594: for (j=-10; j<11; j++) {
4595: value = PetscSign(j)*PetscExpReal(PetscAbs(j)-10.0);
4596: VecSetValue(uh,i,value,ADD_VALUES);
4597: SNESComputeFunction(snes,uh,fh);
4598: VecNorm(fh,NORM_2,&norm);
4599: PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);
4600: value = -value;
4601: VecSetValue(uh,i,value,ADD_VALUES);
4602: }
4603: }
4604: VecDestroy(&uh);
4605: VecDestroy(&fh);
4606: return(0);
4607: }
4609: /*@
4610: SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
4611: computing relative tolerance for linear solvers within an inexact
4612: Newton method.
4614: Logically Collective on SNES
4616: Input Parameters:
4617: + snes - SNES context
4618: - flag - PETSC_TRUE or PETSC_FALSE
4620: Options Database:
4621: + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
4622: . -snes_ksp_ew_version ver - version of Eisenstat-Walker method
4623: . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
4624: . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
4625: . -snes_ksp_ew_gamma <gamma> - Sets gamma
4626: . -snes_ksp_ew_alpha <alpha> - Sets alpha
4627: . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
4628: - -snes_ksp_ew_threshold <threshold> - Sets threshold
4630: Notes:
4631: Currently, the default is to use a constant relative tolerance for
4632: the inner linear solvers. Alternatively, one can use the
4633: Eisenstat-Walker method, where the relative convergence tolerance
4634: is reset at each Newton iteration according progress of the nonlinear
4635: solver.
4637: Level: advanced
4639: Reference:
4640: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4641: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
4643: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
4645: .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4646: @*/
4647: PetscErrorCode SNESKSPSetUseEW(SNES snes,PetscBool flag)
4648: {
4652: snes->ksp_ewconv = flag;
4653: return(0);
4654: }
4656: /*@
4657: SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
4658: for computing relative tolerance for linear solvers within an
4659: inexact Newton method.
4661: Not Collective
4663: Input Parameter:
4664: . snes - SNES context
4666: Output Parameter:
4667: . flag - PETSC_TRUE or PETSC_FALSE
4669: Notes:
4670: Currently, the default is to use a constant relative tolerance for
4671: the inner linear solvers. Alternatively, one can use the
4672: Eisenstat-Walker method, where the relative convergence tolerance
4673: is reset at each Newton iteration according progress of the nonlinear
4674: solver.
4676: Level: advanced
4678: Reference:
4679: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4680: inexact Newton method", SISC 17 (1), pp.16-32, 1996.
4682: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
4684: .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
4685: @*/
4686: PetscErrorCode SNESKSPGetUseEW(SNES snes, PetscBool *flag)
4687: {
4691: *flag = snes->ksp_ewconv;
4692: return(0);
4693: }
4695: /*@
4696: SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
4697: convergence criteria for the linear solvers within an inexact
4698: Newton method.
4700: Logically Collective on SNES
4702: Input Parameters:
4703: + snes - SNES context
4704: . version - version 1, 2 (default is 2) or 3
4705: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4706: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4707: . gamma - multiplicative factor for version 2 rtol computation
4708: (0 <= gamma2 <= 1)
4709: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
4710: . alpha2 - power for safeguard
4711: - threshold - threshold for imposing safeguard (0 < threshold < 1)
4713: Note:
4714: Version 3 was contributed by Luis Chacon, June 2006.
4716: Use PETSC_DEFAULT to retain the default for any of the parameters.
4718: Level: advanced
4720: Reference:
4721: S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
4722: inexact Newton method", Utah State University Math. Stat. Dept. Res.
4723: Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
4725: .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
4727: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
4728: @*/
4729: PetscErrorCode SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
4730: {
4731: SNESKSPEW *kctx;
4735: kctx = (SNESKSPEW*)snes->kspconvctx;
4736: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4745: if (version != PETSC_DEFAULT) kctx->version = version;
4746: if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0;
4747: if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max;
4748: if (gamma != PETSC_DEFAULT) kctx->gamma = gamma;
4749: if (alpha != PETSC_DEFAULT) kctx->alpha = alpha;
4750: if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2;
4751: if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
4753: 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);
4754: 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);
4755: 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);
4756: 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);
4757: 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);
4758: 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);
4759: return(0);
4760: }
4762: /*@
4763: SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
4764: convergence criteria for the linear solvers within an inexact
4765: Newton method.
4767: Not Collective
4769: Input Parameters:
4770: snes - SNES context
4772: Output Parameters:
4773: + version - version 1, 2 (default is 2) or 3
4774: . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
4775: . rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
4776: . gamma - multiplicative factor for version 2 rtol computation (0 <= gamma2 <= 1)
4777: . alpha - power for version 2 rtol computation (1 < alpha <= 2)
4778: . alpha2 - power for safeguard
4779: - threshold - threshold for imposing safeguard (0 < threshold < 1)
4781: Level: advanced
4783: .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
4785: .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
4786: @*/
4787: PetscErrorCode SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
4788: {
4789: SNESKSPEW *kctx;
4793: kctx = (SNESKSPEW*)snes->kspconvctx;
4794: if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
4795: if (version) *version = kctx->version;
4796: if (rtol_0) *rtol_0 = kctx->rtol_0;
4797: if (rtol_max) *rtol_max = kctx->rtol_max;
4798: if (gamma) *gamma = kctx->gamma;
4799: if (alpha) *alpha = kctx->alpha;
4800: if (alpha2) *alpha2 = kctx->alpha2;
4801: if (threshold) *threshold = kctx->threshold;
4802: return(0);
4803: }
4805: PetscErrorCode KSPPreSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4806: {
4808: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
4809: PetscReal rtol = PETSC_DEFAULT,stol;
4812: if (!snes->ksp_ewconv) return(0);
4813: if (!snes->iter) {
4814: rtol = kctx->rtol_0; /* first time in, so use the original user rtol */
4815: VecNorm(snes->vec_func,NORM_2,&kctx->norm_first);
4816: }
4817: else {
4818: if (kctx->version == 1) {
4819: rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
4820: if (rtol < 0.0) rtol = -rtol;
4821: stol = PetscPowReal(kctx->rtol_last,kctx->alpha2);
4822: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4823: } else if (kctx->version == 2) {
4824: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4825: stol = kctx->gamma * PetscPowReal(kctx->rtol_last,kctx->alpha);
4826: if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
4827: } else if (kctx->version == 3) { /* contributed by Luis Chacon, June 2006. */
4828: rtol = kctx->gamma * PetscPowReal(snes->norm/kctx->norm_last,kctx->alpha);
4829: /* safeguard: avoid sharp decrease of rtol */
4830: stol = kctx->gamma*PetscPowReal(kctx->rtol_last,kctx->alpha);
4831: stol = PetscMax(rtol,stol);
4832: rtol = PetscMin(kctx->rtol_0,stol);
4833: /* safeguard: avoid oversolving */
4834: stol = kctx->gamma*(kctx->norm_first*snes->rtol)/snes->norm;
4835: stol = PetscMax(rtol,stol);
4836: rtol = PetscMin(kctx->rtol_0,stol);
4837: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
4838: }
4839: /* safeguard: avoid rtol greater than one */
4840: rtol = PetscMin(rtol,kctx->rtol_max);
4841: KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
4842: PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%g\n",snes->iter,kctx->version,(double)rtol);
4843: return(0);
4844: }
4846: PetscErrorCode KSPPostSolve_SNESEW(KSP ksp, Vec b, Vec x, SNES snes)
4847: {
4849: SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx;
4850: PCSide pcside;
4851: Vec lres;
4854: if (!snes->ksp_ewconv) return(0);
4855: KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);
4856: kctx->norm_last = snes->norm;
4857: if (kctx->version == 1) {
4858: PC pc;
4859: PetscBool isNone;
4861: KSPGetPC(ksp, &pc);
4862: PetscObjectTypeCompare((PetscObject) pc, PCNONE, &isNone);
4863: KSPGetPCSide(ksp,&pcside);
4864: if (pcside == PC_RIGHT || isNone) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
4865: /* KSP residual is true linear residual */
4866: KSPGetResidualNorm(ksp,&kctx->lresid_last);
4867: } else {
4868: /* KSP residual is preconditioned residual */
4869: /* compute true linear residual norm */
4870: VecDuplicate(b,&lres);
4871: MatMult(snes->jacobian,x,lres);
4872: VecAYPX(lres,-1.0,b);
4873: VecNorm(lres,NORM_2,&kctx->lresid_last);
4874: VecDestroy(&lres);
4875: }
4876: }
4877: return(0);
4878: }
4880: /*@
4881: SNESGetKSP - Returns the KSP context for a SNES solver.
4883: Not Collective, but if SNES object is parallel, then KSP object is parallel
4885: Input Parameter:
4886: . snes - the SNES context
4888: Output Parameter:
4889: . ksp - the KSP context
4891: Notes:
4892: The user can then directly manipulate the KSP context to set various
4893: options, etc. Likewise, the user can then extract and manipulate the
4894: PC contexts as well.
4896: Level: beginner
4898: .keywords: SNES, nonlinear, get, KSP, context
4900: .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
4901: @*/
4902: PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp)
4903: {
4910: if (!snes->ksp) {
4911: PetscBool monitor = PETSC_FALSE;
4913: KSPCreate(PetscObjectComm((PetscObject)snes),&snes->ksp);
4914: PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);
4915: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->ksp);
4917: KSPSetPreSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPreSolve_SNESEW,snes);
4918: KSPSetPostSolve(snes->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPostSolve_SNESEW,snes);
4920: PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes",&monitor,NULL);
4921: if (monitor) {
4922: KSPMonitorSet(snes->ksp,KSPMonitorSNES,snes,NULL);
4923: }
4924: monitor = PETSC_FALSE;
4925: PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-ksp_monitor_snes_lg",&monitor,NULL);
4926: if (monitor) {
4927: PetscObject *objs;
4928: KSPMonitorSNESLGResidualNormCreate(PetscObjectComm((PetscObject)snes),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,600,600,&objs);
4929: objs[0] = (PetscObject) snes;
4930: KSPMonitorSet(snes->ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSNESLGResidualNorm,objs,(PetscErrorCode (*)(void**))KSPMonitorSNESLGResidualNormDestroy);
4931: }
4932: }
4933: *ksp = snes->ksp;
4934: return(0);
4935: }
4938: #include <petsc/private/dmimpl.h>
4939: /*@
4940: SNESSetDM - Sets the DM that may be used by some nonlinear solvers or their underlying preconditioners
4942: Logically Collective on SNES
4944: Input Parameters:
4945: + snes - the nonlinear solver context
4946: - dm - the dm, cannot be NULL
4948: Level: intermediate
4950: .seealso: SNESGetDM(), SNESHasDM(), KSPSetDM(), KSPGetDM()
4951: @*/
4952: PetscErrorCode SNESSetDM(SNES snes,DM dm)
4953: {
4955: KSP ksp;
4956: DMSNES sdm;
4961: PetscObjectReference((PetscObject)dm);
4962: if (snes->dm) { /* Move the DMSNES context over to the new DM unless the new DM already has one */
4963: if (snes->dm->dmsnes && !dm->dmsnes) {
4964: DMCopyDMSNES(snes->dm,dm);
4965: DMGetDMSNES(snes->dm,&sdm);
4966: if (sdm->originaldm == snes->dm) sdm->originaldm = dm; /* Grant write privileges to the replacement DM */
4967: }
4968: DMCoarsenHookRemove(snes->dm,DMCoarsenHook_SNESVecSol,DMRestrictHook_SNESVecSol,snes);
4969: DMDestroy(&snes->dm);
4970: }
4971: snes->dm = dm;
4972: snes->dmAuto = PETSC_FALSE;
4974: SNESGetKSP(snes,&ksp);
4975: KSPSetDM(ksp,dm);
4976: KSPSetDMActive(ksp,PETSC_FALSE);
4977: if (snes->npc) {
4978: SNESSetDM(snes->npc, snes->dm);
4979: SNESSetNPCSide(snes,snes->npcside);
4980: }
4981: return(0);
4982: }
4984: /*@
4985: SNESGetDM - Gets the DM that may be used by some preconditioners
4987: Not Collective but DM obtained is parallel on SNES
4989: Input Parameter:
4990: . snes - the preconditioner context
4992: Output Parameter:
4993: . dm - the dm
4995: Level: intermediate
4997: .seealso: SNESSetDM(), SNESHasDM(), KSPSetDM(), KSPGetDM()
4998: @*/
4999: PetscErrorCode SNESGetDM(SNES snes,DM *dm)
5000: {
5005: if (!snes->dm) {
5006: DMShellCreate(PetscObjectComm((PetscObject)snes),&snes->dm);
5007: snes->dmAuto = PETSC_TRUE;
5008: }
5009: *dm = snes->dm;
5010: return(0);
5011: }
5014: /*@
5015: SNESHasDM - Whether snes has dm
5017: Not collective but all processes must return the same value
5019: Input Parameter:
5020: . snes - the nonlinear solver object
5022: Output Parameter:
5023: . hasdm - a flag indicates whether there is dm in snes
5025: Level: intermediate
5027: .seealso: SNESGetDM(), SNESSetDM(), KSPSetDM(), KSPGetDM()
5028: @*/
5029: PetscErrorCode SNESHasDM(SNES snes,PetscBool *hasdm)
5030: {
5034: if (snes->dm) *hasdm = PETSC_TRUE;
5035: else *hasdm = PETSC_FALSE;
5036: return(0);
5037: }
5039: /*@
5040: SNESSetNPC - Sets the nonlinear preconditioner to be used.
5042: Collective on SNES
5044: Input Parameters:
5045: + snes - iterative context obtained from SNESCreate()
5046: - pc - the preconditioner object
5048: Notes:
5049: Use SNESGetNPC() to retrieve the preconditioner context (for example,
5050: to configure it using the API).
5052: Level: developer
5054: .keywords: SNES, set, precondition
5055: .seealso: SNESGetNPC(), SNESHasNPC()
5056: @*/
5057: PetscErrorCode SNESSetNPC(SNES snes, SNES pc)
5058: {
5065: PetscObjectReference((PetscObject) pc);
5066: SNESDestroy(&snes->npc);
5067: snes->npc = pc;
5068: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->npc);
5069: return(0);
5070: }
5072: /*@
5073: SNESGetNPC - Creates a nonlinear preconditioning solver (SNES) to be used to precondition the nonlinear solver.
5075: Not Collective
5077: Input Parameter:
5078: . snes - iterative context obtained from SNESCreate()
5080: Output Parameter:
5081: . pc - preconditioner context
5083: Notes: If a SNES was previously set with SNESSetNPC() then that SNES is returned.
5085: Level: developer
5087: .keywords: SNES, get, preconditioner
5088: .seealso: SNESSetNPC(), SNESHasNPC()
5089: @*/
5090: PetscErrorCode SNESGetNPC(SNES snes, SNES *pc)
5091: {
5093: const char *optionsprefix;
5098: if (!snes->npc) {
5099: SNESCreate(PetscObjectComm((PetscObject)snes),&snes->npc);
5100: PetscObjectIncrementTabLevel((PetscObject)snes->npc,(PetscObject)snes,1);
5101: PetscLogObjectParent((PetscObject)snes,(PetscObject)snes->npc);
5102: SNESGetOptionsPrefix(snes,&optionsprefix);
5103: SNESSetOptionsPrefix(snes->npc,optionsprefix);
5104: SNESAppendOptionsPrefix(snes->npc,"npc_");
5105: SNESSetCountersReset(snes->npc,PETSC_FALSE);
5106: }
5107: *pc = snes->npc;
5108: return(0);
5109: }
5111: /*@
5112: SNESHasNPC - Returns whether a nonlinear preconditioner exists
5114: Not Collective
5116: Input Parameter:
5117: . snes - iterative context obtained from SNESCreate()
5119: Output Parameter:
5120: . has_npc - whether the SNES has an NPC or not
5122: Level: developer
5124: .keywords: SNES, has, preconditioner
5125: .seealso: SNESSetNPC(), SNESGetNPC()
5126: @*/
5127: PetscErrorCode SNESHasNPC(SNES snes, PetscBool *has_npc)
5128: {
5131: *has_npc = (PetscBool) (snes->npc ? PETSC_TRUE : PETSC_FALSE);
5132: return(0);
5133: }
5135: /*@
5136: SNESSetNPCSide - Sets the preconditioning side.
5138: Logically Collective on SNES
5140: Input Parameter:
5141: . snes - iterative context obtained from SNESCreate()
5143: Output Parameter:
5144: . side - the preconditioning side, where side is one of
5145: .vb
5146: PC_LEFT - left preconditioning
5147: PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5148: .ve
5150: Options Database Keys:
5151: . -snes_pc_side <right,left>
5153: Notes: SNESNRICHARDSON and SNESNCG only support left preconditioning.
5155: Level: intermediate
5157: .keywords: SNES, set, right, left, side, preconditioner, flag
5159: .seealso: SNESGetNPCSide(), KSPSetPCSide()
5160: @*/
5161: PetscErrorCode SNESSetNPCSide(SNES snes,PCSide side)
5162: {
5166: snes->npcside= side;
5167: return(0);
5168: }
5170: /*@
5171: SNESGetNPCSide - Gets the preconditioning side.
5173: Not Collective
5175: Input Parameter:
5176: . snes - iterative context obtained from SNESCreate()
5178: Output Parameter:
5179: . side - the preconditioning side, where side is one of
5180: .vb
5181: PC_LEFT - left preconditioning
5182: PC_RIGHT - right preconditioning (default for most nonlinear solvers)
5183: .ve
5185: Level: intermediate
5187: .keywords: SNES, get, right, left, side, preconditioner, flag
5189: .seealso: SNESSetNPCSide(), KSPGetPCSide()
5190: @*/
5191: PetscErrorCode SNESGetNPCSide(SNES snes,PCSide *side)
5192: {
5196: *side = snes->npcside;
5197: return(0);
5198: }
5200: /*@
5201: SNESSetLineSearch - Sets the linesearch on the SNES instance.
5203: Collective on SNES
5205: Input Parameters:
5206: + snes - iterative context obtained from SNESCreate()
5207: - linesearch - the linesearch object
5209: Notes:
5210: Use SNESGetLineSearch() to retrieve the preconditioner context (for example,
5211: to configure it using the API).
5213: Level: developer
5215: .keywords: SNES, set, linesearch
5216: .seealso: SNESGetLineSearch()
5217: @*/
5218: PetscErrorCode SNESSetLineSearch(SNES snes, SNESLineSearch linesearch)
5219: {
5226: PetscObjectReference((PetscObject) linesearch);
5227: SNESLineSearchDestroy(&snes->linesearch);
5229: snes->linesearch = linesearch;
5231: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5232: return(0);
5233: }
5235: /*@
5236: SNESGetLineSearch - Returns a pointer to the line search context set with SNESSetLineSearch()
5237: or creates a default line search instance associated with the SNES and returns it.
5239: Not Collective
5241: Input Parameter:
5242: . snes - iterative context obtained from SNESCreate()
5244: Output Parameter:
5245: . linesearch - linesearch context
5247: Level: beginner
5249: .keywords: SNES, get, linesearch
5250: .seealso: SNESSetLineSearch(), SNESLineSearchCreate()
5251: @*/
5252: PetscErrorCode SNESGetLineSearch(SNES snes, SNESLineSearch *linesearch)
5253: {
5255: const char *optionsprefix;
5260: if (!snes->linesearch) {
5261: SNESGetOptionsPrefix(snes, &optionsprefix);
5262: SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &snes->linesearch);
5263: SNESLineSearchSetSNES(snes->linesearch, snes);
5264: SNESLineSearchAppendOptionsPrefix(snes->linesearch, optionsprefix);
5265: PetscObjectIncrementTabLevel((PetscObject) snes->linesearch, (PetscObject) snes, 1);
5266: PetscLogObjectParent((PetscObject)snes, (PetscObject)snes->linesearch);
5267: }
5268: *linesearch = snes->linesearch;
5269: return(0);
5270: }
5272: #if defined(PETSC_HAVE_MATLAB_ENGINE)
5273: #include <mex.h>
5275: typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
5277: /*
5278: SNESComputeFunction_Matlab - Calls the function that has been set with SNESSetFunctionMatlab().
5280: Collective on SNES
5282: Input Parameters:
5283: + snes - the SNES context
5284: - x - input vector
5286: Output Parameter:
5287: . y - function vector, as set by SNESSetFunction()
5289: Notes:
5290: SNESComputeFunction() is typically used within nonlinear solvers
5291: implementations, so most users would not generally call this routine
5292: themselves.
5294: Level: developer
5296: .keywords: SNES, nonlinear, compute, function
5298: .seealso: SNESSetFunction(), SNESGetFunction()
5299: */
5300: PetscErrorCode SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
5301: {
5302: PetscErrorCode ierr;
5303: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5304: int nlhs = 1,nrhs = 5;
5305: mxArray *plhs[1],*prhs[5];
5306: long long int lx = 0,ly = 0,ls = 0;
5315: /* call Matlab function in ctx with arguments x and y */
5317: PetscMemcpy(&ls,&snes,sizeof(snes));
5318: PetscMemcpy(&lx,&x,sizeof(x));
5319: PetscMemcpy(&ly,&y,sizeof(x));
5320: prhs[0] = mxCreateDoubleScalar((double)ls);
5321: prhs[1] = mxCreateDoubleScalar((double)lx);
5322: prhs[2] = mxCreateDoubleScalar((double)ly);
5323: prhs[3] = mxCreateString(sctx->funcname);
5324: prhs[4] = sctx->ctx;
5325: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");
5326: mxGetScalar(plhs[0]);
5327: mxDestroyArray(prhs[0]);
5328: mxDestroyArray(prhs[1]);
5329: mxDestroyArray(prhs[2]);
5330: mxDestroyArray(prhs[3]);
5331: mxDestroyArray(plhs[0]);
5332: return(0);
5333: }
5335: /*
5336: SNESSetFunctionMatlab - Sets the function evaluation routine and function
5337: vector for use by the SNES routines in solving systems of nonlinear
5338: equations from MATLAB. Here the function is a string containing the name of a MATLAB function
5340: Logically Collective on SNES
5342: Input Parameters:
5343: + snes - the SNES context
5344: . r - vector to store function value
5345: - f - function evaluation routine
5347: Notes:
5348: The Newton-like methods typically solve linear systems of the form
5349: $ f'(x) x = -f(x),
5350: where f'(x) denotes the Jacobian matrix and f(x) is the function.
5352: Level: beginner
5354: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5356: .keywords: SNES, nonlinear, set, function
5358: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5359: */
5360: PetscErrorCode SNESSetFunctionMatlab(SNES snes,Vec r,const char *f,mxArray *ctx)
5361: {
5362: PetscErrorCode ierr;
5363: SNESMatlabContext *sctx;
5366: /* currently sctx is memory bleed */
5367: PetscNew(&sctx);
5368: PetscStrallocpy(f,&sctx->funcname);
5369: /*
5370: This should work, but it doesn't
5371: sctx->ctx = ctx;
5372: mexMakeArrayPersistent(sctx->ctx);
5373: */
5374: sctx->ctx = mxDuplicateArray(ctx);
5375: SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);
5376: return(0);
5377: }
5379: /*
5380: SNESComputeJacobian_Matlab - Calls the function that has been set with SNESSetJacobianMatlab().
5382: Collective on SNES
5384: Input Parameters:
5385: + snes - the SNES context
5386: . x - input vector
5387: . A, B - the matrices
5388: - ctx - user context
5390: Level: developer
5392: .keywords: SNES, nonlinear, compute, function
5394: .seealso: SNESSetFunction(), SNESGetFunction()
5395: @*/
5396: PetscErrorCode SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat A,Mat B,void *ctx)
5397: {
5398: PetscErrorCode ierr;
5399: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5400: int nlhs = 2,nrhs = 6;
5401: mxArray *plhs[2],*prhs[6];
5402: long long int lx = 0,lA = 0,ls = 0, lB = 0;
5408: /* call Matlab function in ctx with arguments x and y */
5410: PetscMemcpy(&ls,&snes,sizeof(snes));
5411: PetscMemcpy(&lx,&x,sizeof(x));
5412: PetscMemcpy(&lA,A,sizeof(x));
5413: PetscMemcpy(&lB,B,sizeof(x));
5414: prhs[0] = mxCreateDoubleScalar((double)ls);
5415: prhs[1] = mxCreateDoubleScalar((double)lx);
5416: prhs[2] = mxCreateDoubleScalar((double)lA);
5417: prhs[3] = mxCreateDoubleScalar((double)lB);
5418: prhs[4] = mxCreateString(sctx->funcname);
5419: prhs[5] = sctx->ctx;
5420: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");
5421: mxGetScalar(plhs[0]);
5422: mxDestroyArray(prhs[0]);
5423: mxDestroyArray(prhs[1]);
5424: mxDestroyArray(prhs[2]);
5425: mxDestroyArray(prhs[3]);
5426: mxDestroyArray(prhs[4]);
5427: mxDestroyArray(plhs[0]);
5428: mxDestroyArray(plhs[1]);
5429: return(0);
5430: }
5432: /*
5433: SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
5434: vector for use by the SNES routines in solving systems of nonlinear
5435: equations from MATLAB. Here the function is a string containing the name of a MATLAB function
5437: Logically Collective on SNES
5439: Input Parameters:
5440: + snes - the SNES context
5441: . A,B - Jacobian matrices
5442: . J - function evaluation routine
5443: - ctx - user context
5445: Level: developer
5447: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5449: .keywords: SNES, nonlinear, set, function
5451: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction(), J
5452: */
5453: PetscErrorCode SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *J,mxArray *ctx)
5454: {
5455: PetscErrorCode ierr;
5456: SNESMatlabContext *sctx;
5459: /* currently sctx is memory bleed */
5460: PetscNew(&sctx);
5461: PetscStrallocpy(J,&sctx->funcname);
5462: /*
5463: This should work, but it doesn't
5464: sctx->ctx = ctx;
5465: mexMakeArrayPersistent(sctx->ctx);
5466: */
5467: sctx->ctx = mxDuplicateArray(ctx);
5468: SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);
5469: return(0);
5470: }
5472: /*
5473: SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().
5475: Collective on SNES
5477: .seealso: SNESSetFunction(), SNESGetFunction()
5478: @*/
5479: PetscErrorCode SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
5480: {
5481: PetscErrorCode ierr;
5482: SNESMatlabContext *sctx = (SNESMatlabContext*)ctx;
5483: int nlhs = 1,nrhs = 6;
5484: mxArray *plhs[1],*prhs[6];
5485: long long int lx = 0,ls = 0;
5486: Vec x = snes->vec_sol;
5491: PetscMemcpy(&ls,&snes,sizeof(snes));
5492: PetscMemcpy(&lx,&x,sizeof(x));
5493: prhs[0] = mxCreateDoubleScalar((double)ls);
5494: prhs[1] = mxCreateDoubleScalar((double)it);
5495: prhs[2] = mxCreateDoubleScalar((double)fnorm);
5496: prhs[3] = mxCreateDoubleScalar((double)lx);
5497: prhs[4] = mxCreateString(sctx->funcname);
5498: prhs[5] = sctx->ctx;
5499: mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");
5500: mxGetScalar(plhs[0]);
5501: mxDestroyArray(prhs[0]);
5502: mxDestroyArray(prhs[1]);
5503: mxDestroyArray(prhs[2]);
5504: mxDestroyArray(prhs[3]);
5505: mxDestroyArray(prhs[4]);
5506: mxDestroyArray(plhs[0]);
5507: return(0);
5508: }
5510: /*
5511: SNESMonitorSetMatlab - Sets the monitor function from MATLAB
5513: Level: developer
5515: Developer Note: This bleeds the allocated memory SNESMatlabContext *sctx;
5517: .keywords: SNES, nonlinear, set, function
5519: .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
5520: */
5521: PetscErrorCode SNESMonitorSetMatlab(SNES snes,const char *f,mxArray *ctx)
5522: {
5523: PetscErrorCode ierr;
5524: SNESMatlabContext *sctx;
5527: /* currently sctx is memory bleed */
5528: PetscNew(&sctx);
5529: PetscStrallocpy(f,&sctx->funcname);
5530: /*
5531: This should work, but it doesn't
5532: sctx->ctx = ctx;
5533: mexMakeArrayPersistent(sctx->ctx);
5534: */
5535: sctx->ctx = mxDuplicateArray(ctx);
5536: SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,NULL);
5537: return(0);
5538: }
5540: #endif