xref: /petsc/src/snes/interface/noise/snesnoise.c (revision 7736150ef7de677169cd1cf3bf9c31f522b43d24)
1*7736150eSBarry Smith #ifdef PETSC_RCS_HEADER
2*7736150eSBarry Smith static char vcid[] = "$Id: snesnoise.c,v 1.13 1999/09/27 21:31:40 bsmith Exp $";
3*7736150eSBarry Smith #endif
4*7736150eSBarry Smith 
5*7736150eSBarry Smith #include "src/snes/snesimpl.h"
6*7736150eSBarry Smith 
7*7736150eSBarry Smith /* Data used by Jorge's diff parameter computation method */
8*7736150eSBarry Smith typedef struct {
9*7736150eSBarry Smith   Vec    *workv;          /* work vectors */
10*7736150eSBarry Smith   FILE   *fp;             /* output file */
11*7736150eSBarry Smith   int    function_count;  /* count of function evaluations for diff param estimation */
12*7736150eSBarry Smith   double fnoise_min;      /* minimim allowable noise */
13*7736150eSBarry Smith   double hopt_min;        /* minimum allowable hopt */
14*7736150eSBarry Smith   double h_first_try;     /* first try for h used in diff parameter estimate */
15*7736150eSBarry Smith   int    fnoise_resets;   /* number of times we've reset the noise estimate */
16*7736150eSBarry Smith   int    hopt_resets;     /* number of times we've reset the hopt estimate */
17*7736150eSBarry Smith } DIFFPAR_MORE;
18*7736150eSBarry Smith 
19*7736150eSBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
20*7736150eSBarry Smith #define dnest_ DNEST
21*7736150eSBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
22*7736150eSBarry Smith #define dnest_ dnest
23*7736150eSBarry Smith #endif
24*7736150eSBarry Smith 
25*7736150eSBarry Smith EXTERN_C_BEGIN
26*7736150eSBarry Smith extern void dnest_(int*,Scalar*,Scalar*,Scalar*,Scalar*,Scalar*,int*,Scalar*);
27*7736150eSBarry Smith EXTERN_C_END
28*7736150eSBarry Smith extern int JacMatMultCompare(SNES,Vec,Vec,double);
29*7736150eSBarry Smith extern int SNESDefaultMatrixFreeSetParameters2(Mat,double,double,double);
30*7736150eSBarry Smith extern int SNESUnSetMatrixFreeParameter(SNES snes);
31*7736150eSBarry Smith 
32*7736150eSBarry Smith #undef __FUNC__
33*7736150eSBarry Smith #define __FUNC__ "DiffParameterCreate_More"
34*7736150eSBarry Smith int DiffParameterCreate_More(SNES snes,Vec x,void **outneP)
35*7736150eSBarry Smith {
36*7736150eSBarry Smith   DIFFPAR_MORE *neP;
37*7736150eSBarry Smith   Vec          w;
38*7736150eSBarry Smith   PetscRandom  rctx;  /* random number generator context */
39*7736150eSBarry Smith   int          ierr, flg;
40*7736150eSBarry Smith   char         noise_file[128];
41*7736150eSBarry Smith 
42*7736150eSBarry Smith   PetscFunctionBegin;
43*7736150eSBarry Smith 
44*7736150eSBarry Smith   neP  = PetscNew(DIFFPAR_MORE);CHKPTRQ(neP);
45*7736150eSBarry Smith   ierr = PetscMemzero(neP,sizeof(DIFFPAR_MORE));CHKERRQ(ierr);
46*7736150eSBarry Smith   PLogObjectMemory(snes,sizeof(DIFFPAR_MORE));
47*7736150eSBarry Smith 
48*7736150eSBarry Smith   neP->function_count = 0;
49*7736150eSBarry Smith   neP->fnoise_min     = 1.0e-20;
50*7736150eSBarry Smith   neP->hopt_min       = 1.0e-8;
51*7736150eSBarry Smith   neP->h_first_try    = 1.0e-3;
52*7736150eSBarry Smith   neP->fnoise_resets  = 0;
53*7736150eSBarry Smith   neP->hopt_resets    = 0;
54*7736150eSBarry Smith 
55*7736150eSBarry Smith   /* Create work vectors */
56*7736150eSBarry Smith   ierr = VecDuplicateVecs(x,3,&neP->workv);CHKERRQ(ierr);
57*7736150eSBarry Smith   w = neP->workv[0];
58*7736150eSBarry Smith 
59*7736150eSBarry Smith   /* Set components of vector w to random numbers */
60*7736150eSBarry Smith   ierr = PetscRandomCreate(snes->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
61*7736150eSBarry Smith   ierr = VecSetRandom(rctx,w);CHKERRQ(ierr);
62*7736150eSBarry Smith   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
63*7736150eSBarry Smith 
64*7736150eSBarry Smith   /* Open output file */
65*7736150eSBarry Smith   ierr = OptionsGetString(snes->prefix,"-snes_mf_noise_file",noise_file,128,&flg);CHKERRQ(ierr);
66*7736150eSBarry Smith   if (flg) neP->fp = fopen(noise_file,"w");
67*7736150eSBarry Smith   else     neP->fp = fopen("noise.out","w");
68*7736150eSBarry Smith   if (!neP->fp) SETERRQ(PETSC_ERR_FILE_OPEN,0,"Cannot open file");
69*7736150eSBarry Smith   PLogInfo(snes,"DiffParameterCreate_More: Creating Jorge's differencing parameter context\n");
70*7736150eSBarry Smith 
71*7736150eSBarry Smith   *outneP = neP;
72*7736150eSBarry Smith   PetscFunctionReturn(0);
73*7736150eSBarry Smith }
74*7736150eSBarry Smith 
75*7736150eSBarry Smith #undef __FUNC__
76*7736150eSBarry Smith #define __FUNC__ "DiffParameterDestroy_More"
77*7736150eSBarry Smith int DiffParameterDestroy_More(void *nePv)
78*7736150eSBarry Smith {
79*7736150eSBarry Smith   DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv;
80*7736150eSBarry Smith   int          ierr;
81*7736150eSBarry Smith 
82*7736150eSBarry Smith   PetscFunctionBegin;
83*7736150eSBarry Smith   /* Destroy work vectors and close output file */
84*7736150eSBarry Smith   ierr = VecDestroyVecs(neP->workv,3);CHKERRQ(ierr);
85*7736150eSBarry Smith   fclose(neP->fp);
86*7736150eSBarry Smith   ierr = PetscFree(neP);CHKERRQ(ierr);
87*7736150eSBarry Smith   PetscFunctionReturn(0);
88*7736150eSBarry Smith }
89*7736150eSBarry Smith 
90*7736150eSBarry Smith #undef __FUNC__
91*7736150eSBarry Smith #define __FUNC__ "DiffParameterCompute_More"
92*7736150eSBarry Smith int DiffParameterCompute_More(SNES snes,void *nePv,Vec x,Vec p,double *fnoise,double *hopt)
93*7736150eSBarry Smith {
94*7736150eSBarry Smith   DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv;
95*7736150eSBarry Smith   Vec         w, xp, fvec;    /* work vectors to use in computing h */
96*7736150eSBarry Smith   double      zero = 0.0, hl, hu, h, fnoise_s, fder2_s, alpha;
97*7736150eSBarry Smith   double      fval[7], tab[7][7], eps[7];
98*7736150eSBarry Smith   double      f, rerrf, fder2;
99*7736150eSBarry Smith   int         iter, k, i, j, ierr, info;
100*7736150eSBarry Smith   int         nf = 7;         /* number of function evaluations */
101*7736150eSBarry Smith   int         noise_test, fcount;
102*7736150eSBarry Smith   MPI_Comm    comm = snes->comm;
103*7736150eSBarry Smith   FILE        *fp;
104*7736150eSBarry Smith 
105*7736150eSBarry Smith   PetscFunctionBegin;
106*7736150eSBarry Smith   /* Call to SNESSetUp() just to set data structures in SNES context */
107*7736150eSBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
108*7736150eSBarry Smith 
109*7736150eSBarry Smith   w    = neP->workv[0];
110*7736150eSBarry Smith   xp   = neP->workv[1];
111*7736150eSBarry Smith   fvec = neP->workv[2];
112*7736150eSBarry Smith   fp   = neP->fp;
113*7736150eSBarry Smith 
114*7736150eSBarry Smith   /* Initialize parameters */
115*7736150eSBarry Smith   hl       = zero;
116*7736150eSBarry Smith   hu       = zero;
117*7736150eSBarry Smith   h        = neP->h_first_try;
118*7736150eSBarry Smith   fnoise_s = zero;
119*7736150eSBarry Smith   fder2_s  = zero;
120*7736150eSBarry Smith   fcount   = neP->function_count;
121*7736150eSBarry Smith 
122*7736150eSBarry Smith   /* We have 5 tries to attempt to compute a good hopt value */
123*7736150eSBarry Smith   ierr = SNESGetIterationNumber(snes,&i);CHKERRQ(ierr);
124*7736150eSBarry Smith   ierr = PetscFPrintf(comm,fp,"\n ------- SNES iteration %d ---------\n",i);CHKERRQ(ierr);
125*7736150eSBarry Smith   for (iter=0; iter<5; iter++) {
126*7736150eSBarry Smith     neP->h_first_try = h;
127*7736150eSBarry Smith 
128*7736150eSBarry Smith     /* Compute the nf function values needed to estimate the noise from
129*7736150eSBarry Smith        the difference table */
130*7736150eSBarry Smith     for (k=0; k<nf; k++) {
131*7736150eSBarry Smith       alpha = h * ( k+1 - (nf+1)/2 );
132*7736150eSBarry Smith       ierr = VecWAXPY(&alpha,p,x,xp);CHKERRQ(ierr);
133*7736150eSBarry Smith       ierr = SNESComputeFunction(snes,xp,fvec);CHKERRQ(ierr);
134*7736150eSBarry Smith       neP->function_count++;
135*7736150eSBarry Smith       ierr = VecDot(fvec,w,&fval[k]);CHKERRQ(ierr);
136*7736150eSBarry Smith     }
137*7736150eSBarry Smith     f = fval[(nf+1)/2 - 1];
138*7736150eSBarry Smith 
139*7736150eSBarry Smith     /* Construct the difference table */
140*7736150eSBarry Smith     for (i=0; i<nf; i++) {
141*7736150eSBarry Smith       tab[i][0] = fval[i];
142*7736150eSBarry Smith     }
143*7736150eSBarry Smith     for (j=0; j<6; j++) {
144*7736150eSBarry Smith       for (i=0; i<nf-j; i++) {
145*7736150eSBarry Smith         tab[i][j+1] = tab[i+1][j] - tab[i][j];
146*7736150eSBarry Smith       }
147*7736150eSBarry Smith     }
148*7736150eSBarry Smith 
149*7736150eSBarry Smith     /* Print the difference table */
150*7736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"Difference Table: iter = %d\n",iter);CHKERRQ(ierr);
151*7736150eSBarry Smith     for (i=0; i<nf; i++) {
152*7736150eSBarry Smith       for (j=0; j<nf-i; j++) {
153*7736150eSBarry Smith         ierr = PetscFPrintf(comm,fp," %10.2e ",tab[i][j]);CHKERRQ(ierr);
154*7736150eSBarry Smith       }
155*7736150eSBarry Smith       ierr = PetscFPrintf(comm,fp,"\n");CHKERRQ(ierr);
156*7736150eSBarry Smith     }
157*7736150eSBarry Smith 
158*7736150eSBarry Smith     /* Call the noise estimator */
159*7736150eSBarry Smith     dnest_(&nf,fval,&h,fnoise,&fder2,hopt,&info,eps);
160*7736150eSBarry Smith 
161*7736150eSBarry Smith     /* Output statements */
162*7736150eSBarry Smith     rerrf = *fnoise/PetscAbsScalar(f);
163*7736150eSBarry Smith     if (info == 1) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected");CHKERRQ(ierr);}
164*7736150eSBarry Smith     if (info == 2) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too small");CHKERRQ(ierr);}
165*7736150eSBarry Smith     if (info == 3) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too large");CHKERRQ(ierr);}
166*7736150eSBarry Smith     if (info == 4) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected, but unreliable hopt");CHKERRQ(ierr);}
167*7736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"Approximate epsfcn %g  %g  %g  %g  %g  %g\n",
168*7736150eSBarry Smith         eps[0],eps[1],eps[2],eps[3],eps[4],eps[5]);CHKERRQ(ierr);
169*7736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"h = %g, fnoise = %g, fder2 = %g, rerrf = %g, hopt = %g\n\n",
170*7736150eSBarry Smith             h, *fnoise, fder2, rerrf, *hopt);CHKERRQ(ierr);
171*7736150eSBarry Smith 
172*7736150eSBarry Smith     /* Save fnoise and fder2. */
173*7736150eSBarry Smith     if (*fnoise) fnoise_s = *fnoise;
174*7736150eSBarry Smith     if (fder2)  fder2_s = fder2;
175*7736150eSBarry Smith 
176*7736150eSBarry Smith     /* Check for noise detection. */
177*7736150eSBarry Smith     if (fnoise_s && fder2_s) {
178*7736150eSBarry Smith       *fnoise = fnoise_s;
179*7736150eSBarry Smith       fder2 = fder2_s;
180*7736150eSBarry Smith       *hopt = 1.68*sqrt(*fnoise/PetscAbsScalar(fder2));
181*7736150eSBarry Smith       goto theend;
182*7736150eSBarry Smith     } else {
183*7736150eSBarry Smith 
184*7736150eSBarry Smith       /* Update hl and hu, and determine new h */
185*7736150eSBarry Smith       if (info == 2 || info == 4) {
186*7736150eSBarry Smith         hl = h;
187*7736150eSBarry Smith         if (hu == zero) h = 100*h;
188*7736150eSBarry Smith         else            h = PetscMin(100*h,0.1*hu);
189*7736150eSBarry Smith        } else if (info == 3) {
190*7736150eSBarry Smith         hu = h;
191*7736150eSBarry Smith         h = PetscMax(1.0e-3,sqrt(hl/hu))*hu;
192*7736150eSBarry Smith       }
193*7736150eSBarry Smith     }
194*7736150eSBarry Smith   }
195*7736150eSBarry Smith   theend:
196*7736150eSBarry Smith 
197*7736150eSBarry Smith   if (*fnoise < neP->fnoise_min) {
198*7736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %g, fnoise_min = %g\n",*fnoise,neP->fnoise_min);CHKERRQ(ierr);
199*7736150eSBarry Smith     *fnoise = neP->fnoise_min;
200*7736150eSBarry Smith     neP->fnoise_resets++;
201*7736150eSBarry Smith   }
202*7736150eSBarry Smith   if (*hopt < neP->hopt_min) {
203*7736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %g, hopt_min = %g\n",*hopt,neP->hopt_min);CHKERRQ(ierr);
204*7736150eSBarry Smith     *hopt = neP->hopt_min;
205*7736150eSBarry Smith     neP->hopt_resets++;
206*7736150eSBarry Smith   }
207*7736150eSBarry Smith 
208*7736150eSBarry Smith   ierr = PetscFPrintf(comm,fp,"Errors in derivative:\n");CHKERRQ(ierr);
209*7736150eSBarry Smith   ierr = PetscFPrintf(comm,fp,"f = %g, fnoise = %g, fder2 = %g, hopt = %g\n",f,*fnoise,fder2,*hopt);CHKERRQ(ierr);
210*7736150eSBarry Smith 
211*7736150eSBarry Smith   /* For now, compute h **each** MV Mult!! */
212*7736150eSBarry Smith   /*
213*7736150eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-matrix_free_jorge_each_mvp",&flg);CHKERRQ(ierr);
214*7736150eSBarry Smith   if (!flg) {
215*7736150eSBarry Smith     Mat mat;
216*7736150eSBarry Smith     ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
217*7736150eSBarry Smith     ierr = SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);CHKERRQ(ierr);
218*7736150eSBarry Smith   }
219*7736150eSBarry Smith   */
220*7736150eSBarry Smith   fcount = neP->function_count - fcount;
221*7736150eSBarry Smith   PLogInfo(snes,"DiffParameterCompute_More: fct_now = %d, fct_cum = %d, rerrf=%g, sqrt(noise)=%g, h_more=%g\n",
222*7736150eSBarry Smith            fcount,neP->function_count,rerrf,sqrt(*fnoise),*hopt);
223*7736150eSBarry Smith 
224*7736150eSBarry Smith 
225*7736150eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-noise_test",&noise_test);CHKERRQ(ierr);
226*7736150eSBarry Smith   if (noise_test) {
227*7736150eSBarry Smith     ierr = JacMatMultCompare(snes,x,p,*hopt);CHKERRQ(ierr);
228*7736150eSBarry Smith   }
229*7736150eSBarry Smith   PetscFunctionReturn(0);
230*7736150eSBarry Smith }
231*7736150eSBarry Smith 
232*7736150eSBarry Smith #undef __FUNC__
233*7736150eSBarry Smith #define __FUNC__ "JacMatMultCompare"
234*7736150eSBarry Smith int JacMatMultCompare(SNES snes,Vec x,Vec p,double hopt)
235*7736150eSBarry Smith {
236*7736150eSBarry Smith   Vec          yy1, yy2; /* work vectors */
237*7736150eSBarry Smith   Viewer       view2;    /* viewer */
238*7736150eSBarry Smith   Mat          J;        /* analytic Jacobian (set as preconditioner matrix) */
239*7736150eSBarry Smith   Mat          Jmf;      /* matrix-free Jacobian (set as true system matrix) */
240*7736150eSBarry Smith   double       h;        /* differencing parameter */
241*7736150eSBarry Smith   Vec          f;
242*7736150eSBarry Smith   MatStructure sparsity = DIFFERENT_NONZERO_PATTERN;
243*7736150eSBarry Smith   Scalar       alpha, yy1n, yy2n, enorm;
244*7736150eSBarry Smith   int          i, ierr, printv;
245*7736150eSBarry Smith   char         filename[32];
246*7736150eSBarry Smith   MPI_Comm     comm = snes->comm;
247*7736150eSBarry Smith 
248*7736150eSBarry Smith   PetscFunctionBegin;
249*7736150eSBarry Smith 
250*7736150eSBarry Smith   /* Compute function and analytic Jacobian at x */
251*7736150eSBarry Smith   ierr = SNESGetJacobian(snes,&Jmf,&J,PETSC_NULL);CHKERRQ(ierr);
252*7736150eSBarry Smith   ierr = SNESComputeJacobian(snes,x,&Jmf,&J,&sparsity);CHKERRQ(ierr);
253*7736150eSBarry Smith   ierr = SNESGetFunction(snes,&f,PETSC_NULL);CHKERRQ(ierr);
254*7736150eSBarry Smith   ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
255*7736150eSBarry Smith 
256*7736150eSBarry Smith   /* Duplicate work vectors */
257*7736150eSBarry Smith   ierr = VecDuplicate(x,&yy2);CHKERRQ(ierr);
258*7736150eSBarry Smith   ierr = VecDuplicate(x,&yy1);CHKERRQ(ierr);
259*7736150eSBarry Smith 
260*7736150eSBarry Smith   /* Compute true matrix-vector product */
261*7736150eSBarry Smith   ierr = MatMult(J,p,yy1);CHKERRQ(ierr);
262*7736150eSBarry Smith   ierr = VecNorm(yy1,NORM_2,&yy1n);CHKERRQ(ierr);
263*7736150eSBarry Smith 
264*7736150eSBarry Smith   /* View product vector if desired */
265*7736150eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-print_vecs",&printv);CHKERRQ(ierr);
266*7736150eSBarry Smith   if (printv) {
267*7736150eSBarry Smith     ierr = ViewerASCIIOpen(comm,"y1.out",&view2);CHKERRQ(ierr);
268*7736150eSBarry Smith     ierr = ViewerSetFormat(view2,VIEWER_FORMAT_ASCII_COMMON,PETSC_NULL);CHKERRQ(ierr);
269*7736150eSBarry Smith     ierr = VecView(yy1,view2);CHKERRQ(ierr);
270*7736150eSBarry Smith     ierr = ViewerDestroy(view2);CHKERRQ(ierr);
271*7736150eSBarry Smith   }
272*7736150eSBarry Smith 
273*7736150eSBarry Smith   /* Test Jacobian-vector product computation */
274*7736150eSBarry Smith   alpha = -1.0;
275*7736150eSBarry Smith   h = 0.01 * hopt;
276*7736150eSBarry Smith   for (i=0; i<5; i++) {
277*7736150eSBarry Smith     /* Set differencing parameter for matrix-free multiplication */
278*7736150eSBarry Smith     ierr = SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);CHKERRQ(ierr);
279*7736150eSBarry Smith 
280*7736150eSBarry Smith     /* Compute matrix-vector product via differencing approximation */
281*7736150eSBarry Smith     ierr = MatMult(Jmf,p,yy2);CHKERRQ(ierr);
282*7736150eSBarry Smith     ierr = VecNorm(yy2,NORM_2,&yy2n);CHKERRQ(ierr);
283*7736150eSBarry Smith 
284*7736150eSBarry Smith     /* View product vector if desired */
285*7736150eSBarry Smith     if (printv) {
286*7736150eSBarry Smith       sprintf(filename,"y2.%d.out",i);
287*7736150eSBarry Smith       ierr = ViewerASCIIOpen(comm,filename,&view2);CHKERRQ(ierr);
288*7736150eSBarry Smith       ierr = ViewerSetFormat(view2,VIEWER_FORMAT_ASCII_COMMON,PETSC_NULL);CHKERRQ(ierr);
289*7736150eSBarry Smith       ierr = VecView(yy2,view2);CHKERRQ(ierr);
290*7736150eSBarry Smith       ierr = ViewerDestroy(view2);CHKERRQ(ierr);
291*7736150eSBarry Smith     }
292*7736150eSBarry Smith 
293*7736150eSBarry Smith     /* Compute relative error */
294*7736150eSBarry Smith     ierr = VecAXPY(&alpha,yy1,yy2);CHKERRQ(ierr);
295*7736150eSBarry Smith     ierr = VecNorm(yy2,NORM_2,&enorm);CHKERRQ(ierr);
296*7736150eSBarry Smith     enorm = enorm/yy1n;
297*7736150eSBarry Smith     ierr = PetscFPrintf(comm,stdout,"h = %g: relative error = %g\n",h,enorm);CHKERRQ(ierr);
298*7736150eSBarry Smith     h *= 10.0;
299*7736150eSBarry Smith   }
300*7736150eSBarry Smith   PetscFunctionReturn(0);
301*7736150eSBarry Smith }
302*7736150eSBarry Smith 
303*7736150eSBarry Smith static int lin_its_total = 0;
304*7736150eSBarry Smith 
305*7736150eSBarry Smith int MyMonitor(SNES snes,int its,double fnorm,void *dummy)
306*7736150eSBarry Smith {
307*7736150eSBarry Smith   int ierr, lin_its;
308*7736150eSBarry Smith 
309*7736150eSBarry Smith   PetscFunctionBegin;
310*7736150eSBarry Smith   ierr = SNESGetNumberLinearIterations(snes,&lin_its);CHKERRQ(ierr);
311*7736150eSBarry Smith   lin_its_total += lin_its;
312*7736150eSBarry Smith   ierr = PetscPrintf(snes->comm, "iter = %d, SNES Function norm = %g, lin_its = %d, total_lin_its = %d\n",its,fnorm,lin_its,lin_its_total);CHKERRQ(ierr);
313*7736150eSBarry Smith 
314*7736150eSBarry Smith   ierr = SNESUnSetMatrixFreeParameter(snes);CHKERRQ(ierr);
315*7736150eSBarry Smith   PetscFunctionReturn(0);
316*7736150eSBarry Smith }
317