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