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