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