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