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