Actual source code: snesnoise.c
petsc-3.9.2 2018-05-20
2: #include <petsc/private/snesimpl.h>
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;
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);
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 */
28: PetscBool flg;
29: char noise_file[PETSC_MAX_PATH_LEN];
32: PetscNewLog(snes,&neP);
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;
41: /* Create work vectors */
42: VecDuplicateVecs(x,3,&neP->workv);
43: w = neP->workv[0];
45: /* Set components of vector w to random numbers */
46: PetscRandomCreate(PetscObjectComm((PetscObject)snes),&rctx);
47: PetscRandomSetFromOptions(rctx);
48: VecSetRandom(w,rctx);
49: PetscRandomDestroy(&rctx);
51: /* Open output file */
52: PetscOptionsGetString(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_noise_file",noise_file,PETSC_MAX_PATH_LEN,&flg);
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: PetscInfo(snes,"Creating Jorge's differencing parameter context\n");
58: *outneP = neP;
59: return(0);
60: }
62: PetscErrorCode SNESDiffParameterDestroy_More(void *nePv)
63: {
64: DIFFPAR_MORE *neP = (DIFFPAR_MORE*)nePv;
66: int err;
69: /* Destroy work vectors and close output file */
70: VecDestroyVecs(3,&neP->workv);
71: err = fclose(neP->fp);
72: if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file");
73: PetscFree(neP);
74: return(0);
75: }
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;
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;
94: PetscObjectGetComm((PetscObject)snes,&comm);
95: /* Call to SNESSetUp() just to set data structures in SNES context */
96: if (!snes->setupcalled) {SNESSetUp(snes);}
98: w = neP->workv[0];
99: xp = neP->workv[1];
100: fvec = neP->workv[2];
101: fp = neP->fp;
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;
111: /* We have 5 tries to attempt to compute a good hopt value */
112: SNESGetIterationNumber(snes,&i);
113: PetscFPrintf(comm,fp,"\n ------- SNES iteration %D ---------\n",i);
114: for (iter=0; iter<5; iter++) {
115: neP->h_first_try = h;
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: VecWAXPY(xp,alpha,p,x);
122: SNESComputeFunction(snes,xp,fvec);
123: neP->function_count++;
124: VecDot(fvec,w,&fval[k]);
125: }
126: f = fval[(nf+1)/2 - 1];
128: /* Construct the difference table */
129: for (i=0; i<nf; i++) tab[i][0] = fval[i];
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: }
137: /* Print the difference table */
138: PetscFPrintf(comm,fp,"Difference Table: iter = %D\n",iter);
139: for (i=0; i<nf; i++) {
140: for (j=0; j<nf-i; j++) {
141: PetscFPrintf(comm,fp," %10.2e ",tab[i][j]);
142: }
143: PetscFPrintf(comm,fp,"\n");
144: }
146: /* Call the noise estimator */
147: SNESNoise_dnest_(&nf,fval,&h,fnoise,&fder2,hopt,&info,eps);
149: /* Output statements */
150: rerrf = *fnoise/PetscAbsScalar(f);
151: if (info == 1) {PetscFPrintf(comm,fp,"%s\n","Noise detected");}
152: if (info == 2) {PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too small");}
153: if (info == 3) {PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too large");}
154: if (info == 4) {PetscFPrintf(comm,fp,"%s\n","Noise detected, but unreliable hopt");}
155: 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]);
156: 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);
158: /* Save fnoise and fder2. */
159: if (*fnoise) fnoise_s = *fnoise;
160: if (fder2) fder2_s = fder2;
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 {
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:
183: if (*fnoise < neP->fnoise_min) {
184: PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %g, fnoise_min = %g\n",(double)*fnoise,(double)neP->fnoise_min);
185: *fnoise = neP->fnoise_min;
186: neP->fnoise_resets++;
187: }
188: if (*hopt < neP->hopt_min) {
189: PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %g, hopt_min = %g\n",(double)*hopt,(double)neP->hopt_min);
190: *hopt = neP->hopt_min;
191: neP->hopt_resets++;
192: }
194: PetscFPrintf(comm,fp,"Errors in derivative:\n");
195: PetscFPrintf(comm,fp,"f = %g, fnoise = %g, fder2 = %g, hopt = %g\n",(double)f,(double)*fnoise,(double)fder2,(double)*hopt);
197: /* For now, compute h **each** MV Mult!! */
198: /*
199: PetscOptionsHasName(NULL,"-matrix_free_jorge_each_mvp",&flg);
200: if (!flg) {
201: Mat mat;
202: SNESGetJacobian(snes,&mat,NULL,NULL);
203: SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);
204: }
205: */
206: fcount = neP->function_count - fcount;
207: 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);
209: PetscOptionsGetBool(NULL,NULL,"-noise_test",&noise_test,NULL);
210: if (noise_test) {
211: JacMatMultCompare(snes,x,p,*hopt);
212: }
213: return(0);
214: }
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;
227: PetscInt i;
228: PetscBool printv = PETSC_FALSE;
229: char filename[32];
230: MPI_Comm comm;
233: PetscObjectGetComm((PetscObject)snes,&comm);
234: /* Compute function and analytic Jacobian at x */
235: SNESGetJacobian(snes,&Jmf,&J,NULL,NULL);
236: SNESComputeJacobian(snes,x,Jmf,J);
237: SNESGetFunction(snes,&f,NULL,NULL);
238: SNESComputeFunction(snes,x,f);
240: /* Duplicate work vectors */
241: VecDuplicate(x,&yy2);
242: VecDuplicate(x,&yy1);
244: /* Compute true matrix-vector product */
245: MatMult(J,p,yy1);
246: VecNorm(yy1,NORM_2,&yy1n);
248: /* View product vector if desired */
249: PetscOptionsGetBool(NULL,NULL,"-print_vecs",&printv,NULL);
250: if (printv) {
251: PetscViewerASCIIOpen(comm,"y1.out",&view2);
252: PetscViewerPushFormat(view2,PETSC_VIEWER_ASCII_COMMON);
253: VecView(yy1,view2);
254: PetscViewerPopFormat(view2);
255: PetscViewerDestroy(&view2);
256: }
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: SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);
265: /* Compute matrix-vector product via differencing approximation */
266: MatMult(Jmf,p,yy2);
267: VecNorm(yy2,NORM_2,&yy2n);
269: /* View product vector if desired */
270: if (printv) {
271: sprintf(filename,"y2.%d.out",(int)i);
272: PetscViewerASCIIOpen(comm,filename,&view2);
273: PetscViewerPushFormat(view2,PETSC_VIEWER_ASCII_COMMON);
274: VecView(yy2,view2);
275: PetscViewerPopFormat(view2);
276: PetscViewerDestroy(&view2);
277: }
279: /* Compute relative error */
280: VecAXPY(yy2,alpha,yy1);
281: VecNorm(yy2,NORM_2,&enorm);
282: enorm = enorm/yy1n;
283: PetscFPrintf(comm,stdout,"h = %g: relative error = %g\n",(double)h,(double)enorm);
284: h *= 10.0;
285: }
286: return(0);
287: }
289: static PetscInt lin_its_total = 0;
291: PetscErrorCode SNESNoiseMonitor(SNES snes,PetscInt its,double fnorm,void *dummy)
292: {
294: PetscInt lin_its;
297: SNESGetLinearSolveIterations(snes,&lin_its);
298: lin_its_total += lin_its;
299: 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);
301: SNESUnSetMatrixFreeParameter(snes);
302: return(0);
303: }