xref: /petsc/src/snes/interface/noise/snesnoise.c (revision d2fd7bfc6f0fd2e1d083decbb7cc7d77e16824f0)
1 
2 #include <petsc/private/snesimpl.h>
3 
4 PETSC_INTERN PetscErrorCode SNESDiffParameterCreate_More(SNES,Vec,void**);
5 PETSC_INTERN PetscErrorCode SNESDiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*);
6 PETSC_INTERN PetscErrorCode SNESDiffParameterDestroy_More(void*);
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   PetscInt 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   PetscInt fnoise_resets;    /* number of times we've reset the noise estimate */
17   PetscInt hopt_resets;      /* number of times we've reset the hopt estimate */
18 } DIFFPAR_MORE;
19 
20 
21 PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat,double,double,double);
22 PETSC_INTERN PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes);
23 PETSC_INTERN PetscErrorCode SNESNoise_dnest_(PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*);
24 
25 static PetscErrorCode JacMatMultCompare(SNES,Vec,Vec,double);
26 
27 PetscErrorCode SNESDiffParameterCreate_More(SNES snes,Vec x,void **outneP)
28 {
29   DIFFPAR_MORE   *neP;
30   Vec            w;
31   PetscRandom    rctx;  /* random number generator context */
32   PetscErrorCode ierr;
33   PetscBool      flg;
34   char           noise_file[PETSC_MAX_PATH_LEN];
35 
36   PetscFunctionBegin;
37   ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr);
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(PetscObjectComm((PetscObject)snes),&rctx);CHKERRQ(ierr);
52   ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr);
53   ierr = VecSetRandom(w,rctx);CHKERRQ(ierr);
54   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
55 
56   /* Open output file */
57   ierr = PetscOptionsGetString(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_noise_file",noise_file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
58   if (flg) neP->fp = fopen(noise_file,"w");
59   else     neP->fp = fopen("noise.out","w");
60   if (!neP->fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file");
61   ierr = PetscInfo(snes,"Creating Jorge's differencing parameter context\n");CHKERRQ(ierr);
62 
63   *outneP = neP;
64   PetscFunctionReturn(0);
65 }
66 
67 PetscErrorCode SNESDiffParameterDestroy_More(void *nePv)
68 {
69   DIFFPAR_MORE   *neP = (DIFFPAR_MORE*)nePv;
70   PetscErrorCode ierr;
71   int            err;
72 
73   PetscFunctionBegin;
74   /* Destroy work vectors and close output file */
75   ierr = VecDestroyVecs(3,&neP->workv);CHKERRQ(ierr);
76   err  = fclose(neP->fp);
77   if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
78   ierr = PetscFree(neP);CHKERRQ(ierr);
79   PetscFunctionReturn(0);
80 }
81 
82 PetscErrorCode SNESDiffParameterCompute_More(SNES snes,void *nePv,Vec x,Vec p,double *fnoise,double *hopt)
83 {
84   DIFFPAR_MORE   *neP = (DIFFPAR_MORE*)nePv;
85   Vec            w, xp, fvec;    /* work vectors to use in computing h */
86   double         zero = 0.0, hl, hu, h, fnoise_s, fder2_s;
87   PetscScalar    alpha;
88   PetscScalar    fval[7], tab[7][7], eps[7], f = -1;
89   double         rerrf = -1., fder2;
90   PetscErrorCode ierr;
91   PetscInt       iter, k, i, j,  info;
92   PetscInt       nf = 7;         /* number of function evaluations */
93   PetscInt       fcount;
94   MPI_Comm       comm;
95   FILE           *fp;
96   PetscBool      noise_test = PETSC_FALSE;
97 
98   PetscFunctionBegin;
99   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
100   /* Call to SNESSetUp() just to set data structures in SNES context */
101   if (!snes->setupcalled) {ierr = SNESSetUp(snes);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(xp,alpha,p,x);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++) tab[i][0] = fval[i];
135 
136     for (j=0; j<nf-1; j++) {
137       for (i=0; i<nf-j-1; i++) {
138         tab[i][j+1] = tab[i+1][j] - tab[i][j];
139       }
140     }
141 
142     /* Print the difference table */
143     ierr = PetscFPrintf(comm,fp,"Difference Table: iter = %D\n",iter);CHKERRQ(ierr);
144     for (i=0; i<nf; i++) {
145       for (j=0; j<nf-i; j++) {
146         ierr = PetscFPrintf(comm,fp," %10.2e ",tab[i][j]);CHKERRQ(ierr);
147       }
148       ierr = PetscFPrintf(comm,fp,"\n");CHKERRQ(ierr);
149     }
150 
151     /* Call the noise estimator */
152     ierr = SNESNoise_dnest_(&nf,fval,&h,fnoise,&fder2,hopt,&info,eps);CHKERRQ(ierr);
153 
154     /* Output statements */
155     rerrf = *fnoise/PetscAbsScalar(f);
156     if (info == 1) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected");CHKERRQ(ierr);}
157     if (info == 2) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too small");CHKERRQ(ierr);}
158     if (info == 3) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too large");CHKERRQ(ierr);}
159     if (info == 4) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected, but unreliable hopt");CHKERRQ(ierr);}
160     ierr = PetscFPrintf(comm,fp,"Approximate epsfcn %g  %g  %g  %g  %g  %g\n",(double)eps[0],(double)eps[1],(double)eps[2],(double)eps[3],(double)eps[4],(double)eps[5]);CHKERRQ(ierr);
161     ierr = PetscFPrintf(comm,fp,"h = %g, fnoise = %g, fder2 = %g, rerrf = %g, hopt = %g\n\n",(double)h, (double)*fnoise, (double)fder2, (double)rerrf, (double)*hopt);CHKERRQ(ierr);
162 
163     /* Save fnoise and fder2. */
164     if (*fnoise) fnoise_s = *fnoise;
165     if (fder2) fder2_s = fder2;
166 
167     /* Check for noise detection. */
168     if (fnoise_s && fder2_s) {
169       *fnoise = fnoise_s;
170       fder2   = fder2_s;
171       *hopt   = 1.68*sqrt(*fnoise/PetscAbsScalar(fder2));
172       goto theend;
173     } else {
174 
175       /* Update hl and hu, and determine new h */
176       if (info == 2 || info == 4) {
177         hl = h;
178         if (hu == zero) h = 100*h;
179         else            h = PetscMin(100*h,0.1*hu);
180       } else if (info == 3) {
181         hu = h;
182         h  = PetscMax(1.0e-3,sqrt(hl/hu))*hu;
183       }
184     }
185   }
186 theend:
187 
188   if (*fnoise < neP->fnoise_min) {
189     ierr    = PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %g, fnoise_min = %g\n",(double)*fnoise,(double)neP->fnoise_min);CHKERRQ(ierr);
190     *fnoise = neP->fnoise_min;
191     neP->fnoise_resets++;
192   }
193   if (*hopt < neP->hopt_min) {
194     ierr  = PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %g, hopt_min = %g\n",(double)*hopt,(double)neP->hopt_min);CHKERRQ(ierr);
195     *hopt = neP->hopt_min;
196     neP->hopt_resets++;
197   }
198 
199   ierr = PetscFPrintf(comm,fp,"Errors in derivative:\n");CHKERRQ(ierr);
200   ierr = PetscFPrintf(comm,fp,"f = %g, fnoise = %g, fder2 = %g, hopt = %g\n",(double)f,(double)*fnoise,(double)fder2,(double)*hopt);CHKERRQ(ierr);
201 
202   /* For now, compute h **each** MV Mult!! */
203   /*
204   ierr = PetscOptionsHasName(NULL,"-matrix_free_jorge_each_mvp",&flg);CHKERRQ(ierr);
205   if (!flg) {
206     Mat mat;
207     ierr = SNESGetJacobian(snes,&mat,NULL,NULL);CHKERRQ(ierr);
208     ierr = SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);CHKERRQ(ierr);
209   }
210   */
211   fcount = neP->function_count - fcount;
212   ierr   = PetscInfo5(snes,"fct_now = %D, fct_cum = %D, rerrf=%g, sqrt(noise)=%g, h_more=%g\n",fcount,neP->function_count,(double)rerrf,(double)PetscSqrtReal(*fnoise),(double)*hopt);CHKERRQ(ierr);
213 
214   ierr = PetscOptionsGetBool(NULL,NULL,"-noise_test",&noise_test,NULL);CHKERRQ(ierr);
215   if (noise_test) {
216     ierr = JacMatMultCompare(snes,x,p,*hopt);CHKERRQ(ierr);
217   }
218   PetscFunctionReturn(0);
219 }
220 
221 PetscErrorCode JacMatMultCompare(SNES snes,Vec x,Vec p,double hopt)
222 {
223   Vec            yy1, yy2; /* work vectors */
224   PetscViewer    view2;    /* viewer */
225   Mat            J;        /* analytic Jacobian (set as preconditioner matrix) */
226   Mat            Jmf;      /* matrix-free Jacobian (set as true system matrix) */
227   double         h;        /* differencing parameter */
228   Vec            f;
229   PetscScalar    alpha;
230   PetscReal      yy1n,yy2n,enorm;
231   PetscErrorCode ierr;
232   PetscInt       i;
233   PetscBool      printv = PETSC_FALSE;
234   char           filename[32];
235   MPI_Comm       comm;
236 
237   PetscFunctionBegin;
238   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
239   /* Compute function and analytic Jacobian at x */
240   ierr = SNESGetJacobian(snes,&Jmf,&J,NULL,NULL);CHKERRQ(ierr);
241   ierr = SNESComputeJacobian(snes,x,Jmf,J);CHKERRQ(ierr);
242   ierr = SNESGetFunction(snes,&f,NULL,NULL);CHKERRQ(ierr);
243   ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
244 
245   /* Duplicate work vectors */
246   ierr = VecDuplicate(x,&yy2);CHKERRQ(ierr);
247   ierr = VecDuplicate(x,&yy1);CHKERRQ(ierr);
248 
249   /* Compute true matrix-vector product */
250   ierr = MatMult(J,p,yy1);CHKERRQ(ierr);
251   ierr = VecNorm(yy1,NORM_2,&yy1n);CHKERRQ(ierr);
252 
253   /* View product vector if desired */
254   ierr = PetscOptionsGetBool(NULL,NULL,"-print_vecs",&printv,NULL);CHKERRQ(ierr);
255   if (printv) {
256     ierr = PetscViewerASCIIOpen(comm,"y1.out",&view2);CHKERRQ(ierr);
257     ierr = PetscViewerPushFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
258     ierr = VecView(yy1,view2);CHKERRQ(ierr);
259     ierr = PetscViewerPopFormat(view2);CHKERRQ(ierr);
260     ierr = PetscViewerDestroy(&view2);CHKERRQ(ierr);
261   }
262 
263   /* Test Jacobian-vector product computation */
264   alpha = -1.0;
265   h     = 0.01 * hopt;
266   for (i=0; i<5; i++) {
267     /* Set differencing parameter for matrix-free multiplication */
268     ierr = SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);CHKERRQ(ierr);
269 
270     /* Compute matrix-vector product via differencing approximation */
271     ierr = MatMult(Jmf,p,yy2);CHKERRQ(ierr);
272     ierr = VecNorm(yy2,NORM_2,&yy2n);CHKERRQ(ierr);
273 
274     /* View product vector if desired */
275     if (printv) {
276       sprintf(filename,"y2.%d.out",(int)i);
277       ierr = PetscViewerASCIIOpen(comm,filename,&view2);CHKERRQ(ierr);
278       ierr = PetscViewerPushFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
279       ierr = VecView(yy2,view2);CHKERRQ(ierr);
280       ierr = PetscViewerPopFormat(view2);CHKERRQ(ierr);
281       ierr = PetscViewerDestroy(&view2);CHKERRQ(ierr);
282     }
283 
284     /* Compute relative error */
285     ierr  = VecAXPY(yy2,alpha,yy1);CHKERRQ(ierr);
286     ierr  = VecNorm(yy2,NORM_2,&enorm);CHKERRQ(ierr);
287     enorm = enorm/yy1n;
288     ierr  = PetscFPrintf(comm,stdout,"h = %g: relative error = %g\n",(double)h,(double)enorm);CHKERRQ(ierr);
289     h    *= 10.0;
290   }
291   PetscFunctionReturn(0);
292 }
293 
294 static PetscInt lin_its_total = 0;
295 
296 PetscErrorCode SNESNoiseMonitor(SNES snes,PetscInt its,double fnorm,void *dummy)
297 {
298   PetscErrorCode ierr;
299   PetscInt       lin_its;
300 
301   PetscFunctionBegin;
302   ierr           = SNESGetLinearSolveIterations(snes,&lin_its);CHKERRQ(ierr);
303   lin_its_total += lin_its;
304   ierr           = PetscPrintf(PetscObjectComm((PetscObject)snes), "iter = %D, SNES Function norm = %g, lin_its = %D, total_lin_its = %D\n",its,(double)fnorm,lin_its,lin_its_total);CHKERRQ(ierr);
305 
306   ierr = SNESUnSetMatrixFreeParameter(snes);CHKERRQ(ierr);
307   PetscFunctionReturn(0);
308 }
309