1 2 #include <petsc-private/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 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; 15 16 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); 21 22 #undef __FUNCT__ 23 #define __FUNCT__ "SNESDiffParameterCreate_More" 24 PetscErrorCode SNESDiffParameterCreate_More(SNES snes,Vec x,void **outneP) 25 { 26 DIFFPAR_MORE *neP; 27 Vec w; 28 PetscRandom rctx; /* random number generator context */ 29 PetscErrorCode ierr; 30 PetscBool flg; 31 char noise_file[PETSC_MAX_PATH_LEN]; 32 33 PetscFunctionBegin; 34 ierr = PetscNewLog(snes,DIFFPAR_MORE,&neP);CHKERRQ(ierr); 35 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; 42 43 /* Create work vectors */ 44 ierr = VecDuplicateVecs(x,3,&neP->workv);CHKERRQ(ierr); 45 w = neP->workv[0]; 46 47 /* Set components of vector w to random numbers */ 48 ierr = PetscRandomCreate(((PetscObject)snes)->comm,&rctx);CHKERRQ(ierr); 49 ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 50 ierr = VecSetRandom(w,rctx);CHKERRQ(ierr); 51 ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 52 53 /* Open output file */ 54 ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_mf_noise_file",noise_file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 55 if (flg) neP->fp = fopen(noise_file,"w"); 56 else neP->fp = fopen("noise.out","w"); 57 if (!neP->fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file"); 58 ierr = PetscInfo(snes,"Creating Jorge's differencing parameter context\n");CHKERRQ(ierr); 59 60 *outneP = neP; 61 PetscFunctionReturn(0); 62 } 63 64 #undef __FUNCT__ 65 #define __FUNCT__ "SNESDiffParameterDestroy_More" 66 PetscErrorCode SNESDiffParameterDestroy_More(void *nePv) 67 { 68 DIFFPAR_MORE *neP = (DIFFPAR_MORE*)nePv; 69 PetscErrorCode ierr; 70 int err; 71 72 PetscFunctionBegin; 73 /* Destroy work vectors and close output file */ 74 ierr = VecDestroyVecs(3,&neP->workv);CHKERRQ(ierr); 75 err = fclose(neP->fp); 76 if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 77 ierr = PetscFree(neP);CHKERRQ(ierr); 78 PetscFunctionReturn(0); 79 } 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "SNESDiffParameterCompute_More" 83 PetscErrorCode SNESDiffParameterCompute_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 PetscErrorCode ierr; 92 PetscInt iter, k, i, j, info; 93 PetscInt nf = 7; /* number of function evaluations */ 94 PetscInt fcount; 95 MPI_Comm comm = ((PetscObject)snes)->comm; 96 FILE *fp; 97 PetscBool noise_test = PETSC_FALSE; 98 99 PetscFunctionBegin; 100 /* Call to SNESSetUp() just to set data structures in SNES context */ 101 if (!snes->setupcalled) {ierr = SNESSetUp(snes);CHKERRQ(ierr);} 102 103 w = neP->workv[0]; 104 xp = neP->workv[1]; 105 fvec = neP->workv[2]; 106 fp = neP->fp; 107 108 /* Initialize parameters */ 109 hl = zero; 110 hu = zero; 111 h = neP->h_first_try; 112 fnoise_s = zero; 113 fder2_s = zero; 114 fcount = neP->function_count; 115 116 /* We have 5 tries to attempt to compute a good hopt value */ 117 ierr = SNESGetIterationNumber(snes,&i);CHKERRQ(ierr); 118 ierr = PetscFPrintf(comm,fp,"\n ------- SNES iteration %D ---------\n",i);CHKERRQ(ierr); 119 for (iter=0; iter<5; iter++) { 120 neP->h_first_try = h; 121 122 /* Compute the nf function values needed to estimate the noise from 123 the difference table */ 124 for (k=0; k<nf; k++) { 125 alpha = h * (k+1 - (nf+1)/2); 126 ierr = VecWAXPY(xp,alpha,p,x);CHKERRQ(ierr); 127 ierr = SNESComputeFunction(snes,xp,fvec);CHKERRQ(ierr); 128 neP->function_count++; 129 ierr = VecDot(fvec,w,&fval[k]);CHKERRQ(ierr); 130 } 131 f = fval[(nf+1)/2 - 1]; 132 133 /* Construct the difference table */ 134 for (i=0; i<nf; i++) 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 = SNESNoise_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",eps[0],eps[1],eps[2],eps[3],eps[4],eps[5]);CHKERRQ(ierr); 161 ierr = PetscFPrintf(comm,fp,"h = %G, fnoise = %G, fder2 = %G, rerrf = %G, hopt = %G\n\n",h, *fnoise, fder2, rerrf, *hopt);CHKERRQ(ierr); 162 163 /* Save fnoise and fder2. */ 164 if (*fnoise) fnoise_s = *fnoise; 165 if (fder2) fder2_s = fder2; 166 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 { 174 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: 187 188 if (*fnoise < neP->fnoise_min) { 189 ierr = PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %G, fnoise_min = %G\n",*fnoise,neP->fnoise_min);CHKERRQ(ierr); 190 *fnoise = neP->fnoise_min; 191 neP->fnoise_resets++; 192 } 193 if (*hopt < neP->hopt_min) { 194 ierr = PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %G, hopt_min = %G\n",*hopt,neP->hopt_min);CHKERRQ(ierr); 195 *hopt = neP->hopt_min; 196 neP->hopt_resets++; 197 } 198 199 ierr = PetscFPrintf(comm,fp,"Errors in derivative:\n");CHKERRQ(ierr); 200 ierr = PetscFPrintf(comm,fp,"f = %G, fnoise = %G, fder2 = %G, hopt = %G\n",f,*fnoise,fder2,*hopt);CHKERRQ(ierr); 201 202 /* For now, compute h **each** MV Mult!! */ 203 /* 204 ierr = PetscOptionsHasName(PETSC_NULL,"-matrix_free_jorge_each_mvp",&flg);CHKERRQ(ierr); 205 if (!flg) { 206 Mat mat; 207 ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 208 ierr = SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);CHKERRQ(ierr); 209 } 210 */ 211 fcount = neP->function_count - fcount; 212 ierr = PetscInfo5(snes,"fct_now = %D, fct_cum = %D, rerrf=%G, sqrt(noise)=%G, h_more=%G\n",fcount,neP->function_count,rerrf,sqrt(*fnoise),*hopt);CHKERRQ(ierr); 213 214 ierr = PetscOptionsGetBool(PETSC_NULL,"-noise_test",&noise_test,PETSC_NULL);CHKERRQ(ierr); 215 if (noise_test) { 216 ierr = JacMatMultCompare(snes,x,p,*hopt);CHKERRQ(ierr); 217 } 218 PetscFunctionReturn(0); 219 } 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "JacMatMultCompare" 223 PetscErrorCode 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 PetscErrorCode ierr; 235 PetscInt i; 236 PetscBool printv = PETSC_FALSE; 237 char filename[32]; 238 MPI_Comm comm = ((PetscObject)snes)->comm; 239 240 PetscFunctionBegin; 241 /* Compute function and analytic Jacobian at x */ 242 ierr = SNESGetJacobian(snes,&Jmf,&J,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 243 ierr = SNESComputeJacobian(snes,x,&Jmf,&J,&sparsity);CHKERRQ(ierr); 244 ierr = SNESGetFunction(snes,&f,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 245 ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); 246 247 /* Duplicate work vectors */ 248 ierr = VecDuplicate(x,&yy2);CHKERRQ(ierr); 249 ierr = VecDuplicate(x,&yy1);CHKERRQ(ierr); 250 251 /* Compute true matrix-vector product */ 252 ierr = MatMult(J,p,yy1);CHKERRQ(ierr); 253 ierr = VecNorm(yy1,NORM_2,&yy1n);CHKERRQ(ierr); 254 255 /* View product vector if desired */ 256 ierr = PetscOptionsGetBool(PETSC_NULL,"-print_vecs",&printv,PETSC_NULL);CHKERRQ(ierr); 257 if (printv) { 258 ierr = PetscViewerASCIIOpen(comm,"y1.out",&view2);CHKERRQ(ierr); 259 ierr = PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr); 260 ierr = VecView(yy1,view2);CHKERRQ(ierr); 261 ierr = PetscViewerDestroy(&view2);CHKERRQ(ierr); 262 } 263 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 ierr = SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);CHKERRQ(ierr); 270 271 /* Compute matrix-vector product via differencing approximation */ 272 ierr = MatMult(Jmf,p,yy2);CHKERRQ(ierr); 273 ierr = VecNorm(yy2,NORM_2,&yy2n);CHKERRQ(ierr); 274 275 /* View product vector if desired */ 276 if (printv) { 277 sprintf(filename,"y2.%d.out",(int)i); 278 ierr = PetscViewerASCIIOpen(comm,filename,&view2);CHKERRQ(ierr); 279 ierr = PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr); 280 ierr = VecView(yy2,view2);CHKERRQ(ierr); 281 ierr = PetscViewerDestroy(&view2);CHKERRQ(ierr); 282 } 283 284 /* Compute relative error */ 285 ierr = VecAXPY(yy2,alpha,yy1);CHKERRQ(ierr); 286 ierr = VecNorm(yy2,NORM_2,&enorm);CHKERRQ(ierr); 287 enorm = enorm/yy1n; 288 ierr = PetscFPrintf(comm,stdout,"h = %G: relative error = %G\n",h,enorm);CHKERRQ(ierr); 289 h *= 10.0; 290 } 291 PetscFunctionReturn(0); 292 } 293 294 static PetscInt lin_its_total = 0; 295 296 PetscErrorCode SNESNoiseMonitor(SNES snes,PetscInt its,double fnorm,void *dummy) 297 { 298 PetscErrorCode ierr; 299 PetscInt lin_its; 300 301 PetscFunctionBegin; 302 ierr = SNESGetLinearSolveIterations(snes,&lin_its);CHKERRQ(ierr); 303 lin_its_total += lin_its; 304 ierr = PetscPrintf(((PetscObject)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); 305 306 ierr = SNESUnSetMatrixFreeParameter(snes);CHKERRQ(ierr); 307 PetscFunctionReturn(0); 308 } 309