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