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