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