xref: /petsc/src/snes/interface/noise/snesnoise.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
1051ac159SSatish Balay 
2b45d2f2cSJed Brown #include <petsc-private/snesimpl.h>
37736150eSBarry Smith 
47736150eSBarry Smith /* Data used by Jorge's diff parameter computation method */
57736150eSBarry Smith typedef struct {
67736150eSBarry Smith   Vec      *workv;           /* work vectors */
77736150eSBarry Smith   FILE     *fp;              /* output file */
87736150eSBarry Smith   int      function_count;   /* count of function evaluations for diff param estimation */
97736150eSBarry Smith   double   fnoise_min;       /* minimim allowable noise */
107736150eSBarry Smith   double   hopt_min;         /* minimum allowable hopt */
117736150eSBarry Smith   double   h_first_try;      /* first try for h used in diff parameter estimate */
12a7cc72afSBarry Smith   PetscInt fnoise_resets;    /* number of times we've reset the noise estimate */
13a7cc72afSBarry Smith   PetscInt hopt_resets;      /* number of times we've reset the hopt estimate */
147736150eSBarry Smith } DIFFPAR_MORE;
157736150eSBarry Smith 
167736150eSBarry Smith 
17ca8e9878SJed Brown extern PetscErrorCode SNESNoise_dnest_(PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*);
18ca8e9878SJed Brown static PetscErrorCode JacMatMultCompare(SNES,Vec,Vec,double);
1909573ac7SBarry Smith extern PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat,double,double,double);
2009573ac7SBarry Smith extern PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes);
217736150eSBarry Smith 
224a2ae208SSatish Balay #undef __FUNCT__
238b5db460SBarry Smith #define __FUNCT__ "SNESDiffParameterCreate_More"
248b5db460SBarry Smith PetscErrorCode SNESDiffParameterCreate_More(SNES snes,Vec x,void **outneP)
257736150eSBarry Smith {
267736150eSBarry Smith   DIFFPAR_MORE   *neP;
277736150eSBarry Smith   Vec            w;
287736150eSBarry Smith   PetscRandom    rctx;  /* random number generator context */
29dfbe8321SBarry Smith   PetscErrorCode ierr;
30ace3abfcSBarry Smith   PetscBool      flg;
31e2d1d2b7SBarry Smith   char           noise_file[PETSC_MAX_PATH_LEN];
327736150eSBarry Smith 
337736150eSBarry Smith   PetscFunctionBegin;
3438f2d2fdSLisandro Dalcin   ierr = PetscNewLog(snes,DIFFPAR_MORE,&neP);CHKERRQ(ierr);
35f5af7f23SKarl Rupp 
367736150eSBarry Smith   neP->function_count = 0;
377736150eSBarry Smith   neP->fnoise_min     = 1.0e-20;
387736150eSBarry Smith   neP->hopt_min       = 1.0e-8;
397736150eSBarry Smith   neP->h_first_try    = 1.0e-3;
407736150eSBarry Smith   neP->fnoise_resets  = 0;
417736150eSBarry Smith   neP->hopt_resets    = 0;
427736150eSBarry Smith 
437736150eSBarry Smith   /* Create work vectors */
447736150eSBarry Smith   ierr = VecDuplicateVecs(x,3,&neP->workv);CHKERRQ(ierr);
457736150eSBarry Smith   w    = neP->workv[0];
467736150eSBarry Smith 
477736150eSBarry Smith   /* Set components of vector w to random numbers */
48*ce94432eSBarry Smith   ierr = PetscRandomCreate(PetscObjectComm((PetscObject)snes),&rctx);CHKERRQ(ierr);
49c77d6671SHong Zhang   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
50abb0e124SMatthew Knepley   ierr = VecSetRandom(w,rctx);CHKERRQ(ierr);
516bf464f9SBarry Smith   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
527736150eSBarry Smith 
537736150eSBarry Smith   /* Open output file */
54bcbf2dc5SJed Brown   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_mf_noise_file",noise_file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
557736150eSBarry Smith   if (flg) neP->fp = fopen(noise_file,"w");
567736150eSBarry Smith   else     neP->fp = fopen("noise.out","w");
57e32f2f54SBarry Smith   if (!neP->fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file");
58ae15b995SBarry Smith   ierr = PetscInfo(snes,"Creating Jorge's differencing parameter context\n");CHKERRQ(ierr);
597736150eSBarry Smith 
607736150eSBarry Smith   *outneP = neP;
617736150eSBarry Smith   PetscFunctionReturn(0);
627736150eSBarry Smith }
637736150eSBarry Smith 
644a2ae208SSatish Balay #undef __FUNCT__
658b5db460SBarry Smith #define __FUNCT__ "SNESDiffParameterDestroy_More"
668b5db460SBarry Smith PetscErrorCode SNESDiffParameterDestroy_More(void *nePv)
677736150eSBarry Smith {
687736150eSBarry Smith   DIFFPAR_MORE   *neP = (DIFFPAR_MORE*)nePv;
69dfbe8321SBarry Smith   PetscErrorCode ierr;
70ed9cf6e9SBarry Smith   int            err;
717736150eSBarry Smith 
727736150eSBarry Smith   PetscFunctionBegin;
737736150eSBarry Smith   /* Destroy work vectors and close output file */
74574deadeSBarry Smith   ierr = VecDestroyVecs(3,&neP->workv);CHKERRQ(ierr);
75ed9cf6e9SBarry Smith   err  = fclose(neP->fp);
76e32f2f54SBarry Smith   if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
777736150eSBarry Smith   ierr = PetscFree(neP);CHKERRQ(ierr);
787736150eSBarry Smith   PetscFunctionReturn(0);
797736150eSBarry Smith }
807736150eSBarry Smith 
814a2ae208SSatish Balay #undef __FUNCT__
828b5db460SBarry Smith #define __FUNCT__ "SNESDiffParameterCompute_More"
838b5db460SBarry Smith PetscErrorCode SNESDiffParameterCompute_More(SNES snes,void *nePv,Vec x,Vec p,double *fnoise,double *hopt)
847736150eSBarry Smith {
857736150eSBarry Smith   DIFFPAR_MORE   *neP = (DIFFPAR_MORE*)nePv;
867736150eSBarry Smith   Vec            w, xp, fvec;    /* work vectors to use in computing h */
8790c26489SBarry Smith   double         zero = 0.0, hl, hu, h, fnoise_s, fder2_s;
8890c26489SBarry Smith   PetscScalar    alpha;
8990c26489SBarry Smith   PetscScalar    fval[7], tab[7][7], eps[7], f;
9090c26489SBarry Smith   double         rerrf, fder2;
916849ba73SBarry Smith   PetscErrorCode ierr;
9277431f27SBarry Smith   PetscInt       iter, k, i, j,  info;
9377431f27SBarry Smith   PetscInt       nf = 7;         /* number of function evaluations */
9477431f27SBarry Smith   PetscInt       fcount;
95*ce94432eSBarry Smith   MPI_Comm       comm;
967736150eSBarry Smith   FILE           *fp;
97ace3abfcSBarry Smith   PetscBool      noise_test = PETSC_FALSE;
987736150eSBarry Smith 
997736150eSBarry Smith   PetscFunctionBegin;
100*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1017736150eSBarry Smith   /* Call to SNESSetUp() just to set data structures in SNES context */
10270e92668SMatthew Knepley   if (!snes->setupcalled) {ierr = SNESSetUp(snes);CHKERRQ(ierr);}
1037736150eSBarry Smith 
1047736150eSBarry Smith   w    = neP->workv[0];
1057736150eSBarry Smith   xp   = neP->workv[1];
1067736150eSBarry Smith   fvec = neP->workv[2];
1077736150eSBarry Smith   fp   = neP->fp;
1087736150eSBarry Smith 
1097736150eSBarry Smith   /* Initialize parameters */
1107736150eSBarry Smith   hl       = zero;
1117736150eSBarry Smith   hu       = zero;
1127736150eSBarry Smith   h        = neP->h_first_try;
1137736150eSBarry Smith   fnoise_s = zero;
1147736150eSBarry Smith   fder2_s  = zero;
1157736150eSBarry Smith   fcount   = neP->function_count;
1167736150eSBarry Smith 
1177736150eSBarry Smith   /* We have 5 tries to attempt to compute a good hopt value */
1187736150eSBarry Smith   ierr = SNESGetIterationNumber(snes,&i);CHKERRQ(ierr);
11977431f27SBarry Smith   ierr = PetscFPrintf(comm,fp,"\n ------- SNES iteration %D ---------\n",i);CHKERRQ(ierr);
1207736150eSBarry Smith   for (iter=0; iter<5; iter++) {
1217736150eSBarry Smith     neP->h_first_try = h;
1227736150eSBarry Smith 
1237736150eSBarry Smith     /* Compute the nf function values needed to estimate the noise from
1247736150eSBarry Smith        the difference table */
1257736150eSBarry Smith     for (k=0; k<nf; k++) {
1267736150eSBarry Smith       alpha = h * (k+1 - (nf+1)/2);
1272dcb1b2aSMatthew Knepley       ierr  = VecWAXPY(xp,alpha,p,x);CHKERRQ(ierr);
1287736150eSBarry Smith       ierr  = SNESComputeFunction(snes,xp,fvec);CHKERRQ(ierr);
1297736150eSBarry Smith       neP->function_count++;
1307736150eSBarry Smith       ierr = VecDot(fvec,w,&fval[k]);CHKERRQ(ierr);
1317736150eSBarry Smith     }
1327736150eSBarry Smith     f = fval[(nf+1)/2 - 1];
1337736150eSBarry Smith 
1347736150eSBarry Smith     /* Construct the difference table */
135f5af7f23SKarl Rupp     for (i=0; i<nf; i++) tab[i][0] = fval[i];
136f5af7f23SKarl Rupp 
1377736150eSBarry Smith     for (j=0; j<6; j++) {
1387736150eSBarry Smith       for (i=0; i<nf-j; i++) {
1397736150eSBarry Smith         tab[i][j+1] = tab[i+1][j] - tab[i][j];
1407736150eSBarry Smith       }
1417736150eSBarry Smith     }
1427736150eSBarry Smith 
1437736150eSBarry Smith     /* Print the difference table */
14477431f27SBarry Smith     ierr = PetscFPrintf(comm,fp,"Difference Table: iter = %D\n",iter);CHKERRQ(ierr);
1457736150eSBarry Smith     for (i=0; i<nf; i++) {
1467736150eSBarry Smith       for (j=0; j<nf-i; j++) {
1477736150eSBarry Smith         ierr = PetscFPrintf(comm,fp," %10.2e ",tab[i][j]);CHKERRQ(ierr);
1487736150eSBarry Smith       }
1497736150eSBarry Smith       ierr = PetscFPrintf(comm,fp,"\n");CHKERRQ(ierr);
1507736150eSBarry Smith     }
1517736150eSBarry Smith 
1527736150eSBarry Smith     /* Call the noise estimator */
153ca8e9878SJed Brown     ierr = SNESNoise_dnest_(&nf,fval,&h,fnoise,&fder2,hopt,&info,eps);CHKERRQ(ierr);
1547736150eSBarry Smith 
1557736150eSBarry Smith     /* Output statements */
1567736150eSBarry Smith     rerrf = *fnoise/PetscAbsScalar(f);
1577736150eSBarry Smith     if (info == 1) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected");CHKERRQ(ierr);}
1587736150eSBarry Smith     if (info == 2) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too small");CHKERRQ(ierr);}
1597736150eSBarry Smith     if (info == 3) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too large");CHKERRQ(ierr);}
1607736150eSBarry Smith     if (info == 4) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected, but unreliable hopt");CHKERRQ(ierr);}
161f5af7f23SKarl Rupp     ierr = PetscFPrintf(comm,fp,"Approximate epsfcn %G  %G  %G  %G  %G  %G\n",eps[0],eps[1],eps[2],eps[3],eps[4],eps[5]);CHKERRQ(ierr);
162f5af7f23SKarl Rupp     ierr = PetscFPrintf(comm,fp,"h = %G, fnoise = %G, fder2 = %G, rerrf = %G, hopt = %G\n\n",h, *fnoise, fder2, rerrf, *hopt);CHKERRQ(ierr);
1637736150eSBarry Smith 
1647736150eSBarry Smith     /* Save fnoise and fder2. */
1657736150eSBarry Smith     if (*fnoise) fnoise_s = *fnoise;
1667736150eSBarry Smith     if (fder2) fder2_s = fder2;
1677736150eSBarry Smith 
1687736150eSBarry Smith     /* Check for noise detection. */
1697736150eSBarry Smith     if (fnoise_s && fder2_s) {
1707736150eSBarry Smith       *fnoise = fnoise_s;
1717736150eSBarry Smith       fder2   = fder2_s;
1727736150eSBarry Smith       *hopt   = 1.68*sqrt(*fnoise/PetscAbsScalar(fder2));
1737736150eSBarry Smith       goto theend;
1747736150eSBarry Smith     } else {
1757736150eSBarry Smith 
1767736150eSBarry Smith       /* Update hl and hu, and determine new h */
1777736150eSBarry Smith       if (info == 2 || info == 4) {
1787736150eSBarry Smith         hl = h;
1797736150eSBarry Smith         if (hu == zero) h = 100*h;
1807736150eSBarry Smith         else            h = PetscMin(100*h,0.1*hu);
1817736150eSBarry Smith       } else if (info == 3) {
1827736150eSBarry Smith         hu = h;
1837736150eSBarry Smith         h  = PetscMax(1.0e-3,sqrt(hl/hu))*hu;
1847736150eSBarry Smith       }
1857736150eSBarry Smith     }
1867736150eSBarry Smith   }
1877736150eSBarry Smith theend:
1887736150eSBarry Smith 
1897736150eSBarry Smith   if (*fnoise < neP->fnoise_min) {
190a83599f4SBarry Smith     ierr    = PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %G, fnoise_min = %G\n",*fnoise,neP->fnoise_min);CHKERRQ(ierr);
1917736150eSBarry Smith     *fnoise = neP->fnoise_min;
1927736150eSBarry Smith     neP->fnoise_resets++;
1937736150eSBarry Smith   }
1947736150eSBarry Smith   if (*hopt < neP->hopt_min) {
195a83599f4SBarry Smith     ierr  = PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %G, hopt_min = %G\n",*hopt,neP->hopt_min);CHKERRQ(ierr);
1967736150eSBarry Smith     *hopt = neP->hopt_min;
1977736150eSBarry Smith     neP->hopt_resets++;
1987736150eSBarry Smith   }
1997736150eSBarry Smith 
2007736150eSBarry Smith   ierr = PetscFPrintf(comm,fp,"Errors in derivative:\n");CHKERRQ(ierr);
201a83599f4SBarry Smith   ierr = PetscFPrintf(comm,fp,"f = %G, fnoise = %G, fder2 = %G, hopt = %G\n",f,*fnoise,fder2,*hopt);CHKERRQ(ierr);
2027736150eSBarry Smith 
2037736150eSBarry Smith   /* For now, compute h **each** MV Mult!! */
2047736150eSBarry Smith   /*
2050298fd71SBarry Smith   ierr = PetscOptionsHasName(NULL,"-matrix_free_jorge_each_mvp",&flg);CHKERRQ(ierr);
2067736150eSBarry Smith   if (!flg) {
2077736150eSBarry Smith     Mat mat;
2080298fd71SBarry Smith     ierr = SNESGetJacobian(snes,&mat,NULL,NULL);CHKERRQ(ierr);
2097736150eSBarry Smith     ierr = SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);CHKERRQ(ierr);
2107736150eSBarry Smith   }
2117736150eSBarry Smith   */
2127736150eSBarry Smith   fcount = neP->function_count - fcount;
213ae15b995SBarry Smith   ierr   = PetscInfo5(snes,"fct_now = %D, fct_cum = %D, rerrf=%G, sqrt(noise)=%G, h_more=%G\n",fcount,neP->function_count,rerrf,sqrt(*fnoise),*hopt);CHKERRQ(ierr);
2147736150eSBarry Smith 
2150298fd71SBarry Smith   ierr = PetscOptionsGetBool(NULL,"-noise_test",&noise_test,NULL);CHKERRQ(ierr);
2167736150eSBarry Smith   if (noise_test) {
2177736150eSBarry Smith     ierr = JacMatMultCompare(snes,x,p,*hopt);CHKERRQ(ierr);
2187736150eSBarry Smith   }
2197736150eSBarry Smith   PetscFunctionReturn(0);
2207736150eSBarry Smith }
2217736150eSBarry Smith 
2224a2ae208SSatish Balay #undef __FUNCT__
2234a2ae208SSatish Balay #define __FUNCT__ "JacMatMultCompare"
224dfbe8321SBarry Smith PetscErrorCode JacMatMultCompare(SNES snes,Vec x,Vec p,double hopt)
2257736150eSBarry Smith {
2267736150eSBarry Smith   Vec            yy1, yy2; /* work vectors */
22752ed857cSBarry Smith   PetscViewer    view2;    /* viewer */
2287736150eSBarry Smith   Mat            J;        /* analytic Jacobian (set as preconditioner matrix) */
2297736150eSBarry Smith   Mat            Jmf;      /* matrix-free Jacobian (set as true system matrix) */
2307736150eSBarry Smith   double         h;        /* differencing parameter */
2317736150eSBarry Smith   Vec            f;
2327736150eSBarry Smith   MatStructure   sparsity = DIFFERENT_NONZERO_PATTERN;
23390c26489SBarry Smith   PetscScalar    alpha;
23490c26489SBarry Smith   PetscReal      yy1n,yy2n,enorm;
2356849ba73SBarry Smith   PetscErrorCode ierr;
236a7cc72afSBarry Smith   PetscInt       i;
237ace3abfcSBarry Smith   PetscBool      printv = PETSC_FALSE;
2387736150eSBarry Smith   char           filename[32];
239*ce94432eSBarry Smith   MPI_Comm       comm;
2407736150eSBarry Smith 
2417736150eSBarry Smith   PetscFunctionBegin;
242*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
2437736150eSBarry Smith   /* Compute function and analytic Jacobian at x */
2440298fd71SBarry Smith   ierr = SNESGetJacobian(snes,&Jmf,&J,NULL,NULL);CHKERRQ(ierr);
2457736150eSBarry Smith   ierr = SNESComputeJacobian(snes,x,&Jmf,&J,&sparsity);CHKERRQ(ierr);
2460298fd71SBarry Smith   ierr = SNESGetFunction(snes,&f,NULL,NULL);CHKERRQ(ierr);
2477736150eSBarry Smith   ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
2487736150eSBarry Smith 
2497736150eSBarry Smith   /* Duplicate work vectors */
2507736150eSBarry Smith   ierr = VecDuplicate(x,&yy2);CHKERRQ(ierr);
2517736150eSBarry Smith   ierr = VecDuplicate(x,&yy1);CHKERRQ(ierr);
2527736150eSBarry Smith 
2537736150eSBarry Smith   /* Compute true matrix-vector product */
2547736150eSBarry Smith   ierr = MatMult(J,p,yy1);CHKERRQ(ierr);
2557736150eSBarry Smith   ierr = VecNorm(yy1,NORM_2,&yy1n);CHKERRQ(ierr);
2567736150eSBarry Smith 
2577736150eSBarry Smith   /* View product vector if desired */
2580298fd71SBarry Smith   ierr = PetscOptionsGetBool(NULL,"-print_vecs",&printv,NULL);CHKERRQ(ierr);
2597736150eSBarry Smith   if (printv) {
26052ed857cSBarry Smith     ierr = PetscViewerASCIIOpen(comm,"y1.out",&view2);CHKERRQ(ierr);
26152ed857cSBarry Smith     ierr = PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
2627736150eSBarry Smith     ierr = VecView(yy1,view2);CHKERRQ(ierr);
2636bf464f9SBarry Smith     ierr = PetscViewerDestroy(&view2);CHKERRQ(ierr);
2647736150eSBarry Smith   }
2657736150eSBarry Smith 
2667736150eSBarry Smith   /* Test Jacobian-vector product computation */
2677736150eSBarry Smith   alpha = -1.0;
2687736150eSBarry Smith   h     = 0.01 * hopt;
2697736150eSBarry Smith   for (i=0; i<5; i++) {
2707736150eSBarry Smith     /* Set differencing parameter for matrix-free multiplication */
2717736150eSBarry Smith     ierr = SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);CHKERRQ(ierr);
2727736150eSBarry Smith 
2737736150eSBarry Smith     /* Compute matrix-vector product via differencing approximation */
2747736150eSBarry Smith     ierr = MatMult(Jmf,p,yy2);CHKERRQ(ierr);
2757736150eSBarry Smith     ierr = VecNorm(yy2,NORM_2,&yy2n);CHKERRQ(ierr);
2767736150eSBarry Smith 
2777736150eSBarry Smith     /* View product vector if desired */
2787736150eSBarry Smith     if (printv) {
27977431f27SBarry Smith       sprintf(filename,"y2.%d.out",(int)i);
28052ed857cSBarry Smith       ierr = PetscViewerASCIIOpen(comm,filename,&view2);CHKERRQ(ierr);
28152ed857cSBarry Smith       ierr = PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
2827736150eSBarry Smith       ierr = VecView(yy2,view2);CHKERRQ(ierr);
2836bf464f9SBarry Smith       ierr = PetscViewerDestroy(&view2);CHKERRQ(ierr);
2847736150eSBarry Smith     }
2857736150eSBarry Smith 
2867736150eSBarry Smith     /* Compute relative error */
2872dcb1b2aSMatthew Knepley     ierr  = VecAXPY(yy2,alpha,yy1);CHKERRQ(ierr);
2887736150eSBarry Smith     ierr  = VecNorm(yy2,NORM_2,&enorm);CHKERRQ(ierr);
2897736150eSBarry Smith     enorm = enorm/yy1n;
290a83599f4SBarry Smith     ierr  = PetscFPrintf(comm,stdout,"h = %G: relative error = %G\n",h,enorm);CHKERRQ(ierr);
2917736150eSBarry Smith     h    *= 10.0;
2927736150eSBarry Smith   }
2937736150eSBarry Smith   PetscFunctionReturn(0);
2947736150eSBarry Smith }
2957736150eSBarry Smith 
296a7cc72afSBarry Smith static PetscInt lin_its_total = 0;
2977736150eSBarry Smith 
298b12ffc70SBarry Smith PetscErrorCode SNESNoiseMonitor(SNES snes,PetscInt its,double fnorm,void *dummy)
2997736150eSBarry Smith {
3006849ba73SBarry Smith   PetscErrorCode ierr;
30177431f27SBarry Smith   PetscInt       lin_its;
3027736150eSBarry Smith 
3037736150eSBarry Smith   PetscFunctionBegin;
304b850b91aSLisandro Dalcin   ierr           = SNESGetLinearSolveIterations(snes,&lin_its);CHKERRQ(ierr);
3057736150eSBarry Smith   lin_its_total += lin_its;
306*ce94432eSBarry Smith   ierr           = PetscPrintf(PetscObjectComm((PetscObject)snes), "iter = %D, SNES Function norm = %G, lin_its = %D, total_lin_its = %D\n",its,fnorm,lin_its,lin_its_total);CHKERRQ(ierr);
3077736150eSBarry Smith 
3087736150eSBarry Smith   ierr = SNESUnSetMatrixFreeParameter(snes);CHKERRQ(ierr);
3097736150eSBarry Smith   PetscFunctionReturn(0);
3107736150eSBarry Smith }
311