17736150eSBarry Smith #ifdef PETSC_RCS_HEADER 2*051ac159SSatish Balay static char vcid[] = "$Id: snesnoise.c,v 1.1 2001/01/19 21:11:04 bsmith Exp balay $"; 37736150eSBarry Smith #endif 47736150eSBarry Smith 5*051ac159SSatish Balay 67736150eSBarry Smith #include "src/snes/snesimpl.h" 77736150eSBarry Smith 87736150eSBarry Smith /* Data used by Jorge's diff parameter computation method */ 97736150eSBarry Smith typedef struct { 107736150eSBarry Smith Vec *workv; /* work vectors */ 117736150eSBarry Smith FILE *fp; /* output file */ 127736150eSBarry Smith int function_count; /* count of function evaluations for diff param estimation */ 137736150eSBarry Smith double fnoise_min; /* minimim allowable noise */ 147736150eSBarry Smith double hopt_min; /* minimum allowable hopt */ 157736150eSBarry Smith double h_first_try; /* first try for h used in diff parameter estimate */ 167736150eSBarry Smith int fnoise_resets; /* number of times we've reset the noise estimate */ 177736150eSBarry Smith int hopt_resets; /* number of times we've reset the hopt estimate */ 187736150eSBarry Smith } DIFFPAR_MORE; 197736150eSBarry Smith 207736150eSBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS) 217736150eSBarry Smith #define dnest_ DNEST 227736150eSBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 237736150eSBarry Smith #define dnest_ dnest 247736150eSBarry Smith #endif 257736150eSBarry Smith 267736150eSBarry Smith EXTERN_C_BEGIN 277736150eSBarry Smith extern void dnest_(int*,Scalar*,Scalar*,Scalar*,Scalar*,Scalar*,int*,Scalar*); 287736150eSBarry Smith EXTERN_C_END 297736150eSBarry Smith extern int JacMatMultCompare(SNES,Vec,Vec,double); 307736150eSBarry Smith extern int SNESDefaultMatrixFreeSetParameters2(Mat,double,double,double); 317736150eSBarry Smith extern int SNESUnSetMatrixFreeParameter(SNES snes); 327736150eSBarry Smith 337736150eSBarry Smith #undef __FUNC__ 347736150eSBarry Smith #define __FUNC__ "DiffParameterCreate_More" 357736150eSBarry Smith int DiffParameterCreate_More(SNES snes,Vec x,void **outneP) 367736150eSBarry Smith { 377736150eSBarry Smith DIFFPAR_MORE *neP; 387736150eSBarry Smith Vec w; 397736150eSBarry Smith PetscRandom rctx; /* random number generator context */ 407736150eSBarry Smith int ierr, flg; 417736150eSBarry Smith char noise_file[128]; 427736150eSBarry Smith 437736150eSBarry Smith PetscFunctionBegin; 447736150eSBarry Smith 457736150eSBarry Smith neP = PetscNew(DIFFPAR_MORE);CHKPTRQ(neP); 467736150eSBarry Smith ierr = PetscMemzero(neP,sizeof(DIFFPAR_MORE));CHKERRQ(ierr); 477736150eSBarry Smith PLogObjectMemory(snes,sizeof(DIFFPAR_MORE)); 487736150eSBarry Smith 497736150eSBarry Smith neP->function_count = 0; 507736150eSBarry Smith neP->fnoise_min = 1.0e-20; 517736150eSBarry Smith neP->hopt_min = 1.0e-8; 527736150eSBarry Smith neP->h_first_try = 1.0e-3; 537736150eSBarry Smith neP->fnoise_resets = 0; 547736150eSBarry Smith neP->hopt_resets = 0; 557736150eSBarry Smith 567736150eSBarry Smith /* Create work vectors */ 577736150eSBarry Smith ierr = VecDuplicateVecs(x,3,&neP->workv);CHKERRQ(ierr); 587736150eSBarry Smith w = neP->workv[0]; 597736150eSBarry Smith 607736150eSBarry Smith /* Set components of vector w to random numbers */ 617736150eSBarry Smith ierr = PetscRandomCreate(snes->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr); 627736150eSBarry Smith ierr = VecSetRandom(rctx,w);CHKERRQ(ierr); 637736150eSBarry Smith ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr); 647736150eSBarry Smith 657736150eSBarry Smith /* Open output file */ 667736150eSBarry Smith ierr = OptionsGetString(snes->prefix,"-snes_mf_noise_file",noise_file,128,&flg);CHKERRQ(ierr); 677736150eSBarry Smith if (flg) neP->fp = fopen(noise_file,"w"); 687736150eSBarry Smith else neP->fp = fopen("noise.out","w"); 697736150eSBarry Smith if (!neP->fp) SETERRQ(PETSC_ERR_FILE_OPEN,0,"Cannot open file"); 707736150eSBarry Smith PLogInfo(snes,"DiffParameterCreate_More: Creating Jorge's differencing parameter context\n"); 717736150eSBarry Smith 727736150eSBarry Smith *outneP = neP; 737736150eSBarry Smith PetscFunctionReturn(0); 747736150eSBarry Smith } 757736150eSBarry Smith 767736150eSBarry Smith #undef __FUNC__ 777736150eSBarry Smith #define __FUNC__ "DiffParameterDestroy_More" 787736150eSBarry Smith int DiffParameterDestroy_More(void *nePv) 797736150eSBarry Smith { 807736150eSBarry Smith DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv; 817736150eSBarry Smith int ierr; 827736150eSBarry Smith 837736150eSBarry Smith PetscFunctionBegin; 847736150eSBarry Smith /* Destroy work vectors and close output file */ 857736150eSBarry Smith ierr = VecDestroyVecs(neP->workv,3);CHKERRQ(ierr); 867736150eSBarry Smith fclose(neP->fp); 877736150eSBarry Smith ierr = PetscFree(neP);CHKERRQ(ierr); 887736150eSBarry Smith PetscFunctionReturn(0); 897736150eSBarry Smith } 907736150eSBarry Smith 917736150eSBarry Smith #undef __FUNC__ 927736150eSBarry Smith #define __FUNC__ "DiffParameterCompute_More" 937736150eSBarry Smith int DiffParameterCompute_More(SNES snes,void *nePv,Vec x,Vec p,double *fnoise,double *hopt) 947736150eSBarry Smith { 957736150eSBarry Smith DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv; 967736150eSBarry Smith Vec w, xp, fvec; /* work vectors to use in computing h */ 977736150eSBarry Smith double zero = 0.0, hl, hu, h, fnoise_s, fder2_s, alpha; 987736150eSBarry Smith double fval[7], tab[7][7], eps[7]; 997736150eSBarry Smith double f, rerrf, fder2; 1007736150eSBarry Smith int iter, k, i, j, ierr, info; 1017736150eSBarry Smith int nf = 7; /* number of function evaluations */ 1027736150eSBarry Smith int noise_test, fcount; 1037736150eSBarry Smith MPI_Comm comm = snes->comm; 1047736150eSBarry Smith FILE *fp; 1057736150eSBarry Smith 1067736150eSBarry Smith PetscFunctionBegin; 1077736150eSBarry Smith /* Call to SNESSetUp() just to set data structures in SNES context */ 1087736150eSBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1097736150eSBarry Smith 1107736150eSBarry Smith w = neP->workv[0]; 1117736150eSBarry Smith xp = neP->workv[1]; 1127736150eSBarry Smith fvec = neP->workv[2]; 1137736150eSBarry Smith fp = neP->fp; 1147736150eSBarry Smith 1157736150eSBarry Smith /* Initialize parameters */ 1167736150eSBarry Smith hl = zero; 1177736150eSBarry Smith hu = zero; 1187736150eSBarry Smith h = neP->h_first_try; 1197736150eSBarry Smith fnoise_s = zero; 1207736150eSBarry Smith fder2_s = zero; 1217736150eSBarry Smith fcount = neP->function_count; 1227736150eSBarry Smith 1237736150eSBarry Smith /* We have 5 tries to attempt to compute a good hopt value */ 1247736150eSBarry Smith ierr = SNESGetIterationNumber(snes,&i);CHKERRQ(ierr); 1257736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"\n ------- SNES iteration %d ---------\n",i);CHKERRQ(ierr); 1267736150eSBarry Smith for (iter=0; iter<5; iter++) { 1277736150eSBarry Smith neP->h_first_try = h; 1287736150eSBarry Smith 1297736150eSBarry Smith /* Compute the nf function values needed to estimate the noise from 1307736150eSBarry Smith the difference table */ 1317736150eSBarry Smith for (k=0; k<nf; k++) { 1327736150eSBarry Smith alpha = h * ( k+1 - (nf+1)/2 ); 1337736150eSBarry Smith ierr = VecWAXPY(&alpha,p,x,xp);CHKERRQ(ierr); 1347736150eSBarry Smith ierr = SNESComputeFunction(snes,xp,fvec);CHKERRQ(ierr); 1357736150eSBarry Smith neP->function_count++; 1367736150eSBarry Smith ierr = VecDot(fvec,w,&fval[k]);CHKERRQ(ierr); 1377736150eSBarry Smith } 1387736150eSBarry Smith f = fval[(nf+1)/2 - 1]; 1397736150eSBarry Smith 1407736150eSBarry Smith /* Construct the difference table */ 1417736150eSBarry Smith for (i=0; i<nf; i++) { 1427736150eSBarry Smith tab[i][0] = fval[i]; 1437736150eSBarry Smith } 1447736150eSBarry Smith for (j=0; j<6; j++) { 1457736150eSBarry Smith for (i=0; i<nf-j; i++) { 1467736150eSBarry Smith tab[i][j+1] = tab[i+1][j] - tab[i][j]; 1477736150eSBarry Smith } 1487736150eSBarry Smith } 1497736150eSBarry Smith 1507736150eSBarry Smith /* Print the difference table */ 1517736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"Difference Table: iter = %d\n",iter);CHKERRQ(ierr); 1527736150eSBarry Smith for (i=0; i<nf; i++) { 1537736150eSBarry Smith for (j=0; j<nf-i; j++) { 1547736150eSBarry Smith ierr = PetscFPrintf(comm,fp," %10.2e ",tab[i][j]);CHKERRQ(ierr); 1557736150eSBarry Smith } 1567736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"\n");CHKERRQ(ierr); 1577736150eSBarry Smith } 1587736150eSBarry Smith 1597736150eSBarry Smith /* Call the noise estimator */ 1607736150eSBarry Smith dnest_(&nf,fval,&h,fnoise,&fder2,hopt,&info,eps); 1617736150eSBarry Smith 1627736150eSBarry Smith /* Output statements */ 1637736150eSBarry Smith rerrf = *fnoise/PetscAbsScalar(f); 1647736150eSBarry Smith if (info == 1) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected");CHKERRQ(ierr);} 1657736150eSBarry Smith if (info == 2) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too small");CHKERRQ(ierr);} 1667736150eSBarry Smith if (info == 3) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too large");CHKERRQ(ierr);} 1677736150eSBarry Smith if (info == 4) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected, but unreliable hopt");CHKERRQ(ierr);} 1687736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"Approximate epsfcn %g %g %g %g %g %g\n", 1697736150eSBarry Smith eps[0],eps[1],eps[2],eps[3],eps[4],eps[5]);CHKERRQ(ierr); 1707736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"h = %g, fnoise = %g, fder2 = %g, rerrf = %g, hopt = %g\n\n", 1717736150eSBarry Smith h, *fnoise, fder2, rerrf, *hopt);CHKERRQ(ierr); 1727736150eSBarry Smith 1737736150eSBarry Smith /* Save fnoise and fder2. */ 1747736150eSBarry Smith if (*fnoise) fnoise_s = *fnoise; 1757736150eSBarry Smith if (fder2) fder2_s = fder2; 1767736150eSBarry Smith 1777736150eSBarry Smith /* Check for noise detection. */ 1787736150eSBarry Smith if (fnoise_s && fder2_s) { 1797736150eSBarry Smith *fnoise = fnoise_s; 1807736150eSBarry Smith fder2 = fder2_s; 1817736150eSBarry Smith *hopt = 1.68*sqrt(*fnoise/PetscAbsScalar(fder2)); 1827736150eSBarry Smith goto theend; 1837736150eSBarry Smith } else { 1847736150eSBarry Smith 1857736150eSBarry Smith /* Update hl and hu, and determine new h */ 1867736150eSBarry Smith if (info == 2 || info == 4) { 1877736150eSBarry Smith hl = h; 1887736150eSBarry Smith if (hu == zero) h = 100*h; 1897736150eSBarry Smith else h = PetscMin(100*h,0.1*hu); 1907736150eSBarry Smith } else if (info == 3) { 1917736150eSBarry Smith hu = h; 1927736150eSBarry Smith h = PetscMax(1.0e-3,sqrt(hl/hu))*hu; 1937736150eSBarry Smith } 1947736150eSBarry Smith } 1957736150eSBarry Smith } 1967736150eSBarry Smith theend: 1977736150eSBarry Smith 1987736150eSBarry Smith if (*fnoise < neP->fnoise_min) { 1997736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %g, fnoise_min = %g\n",*fnoise,neP->fnoise_min);CHKERRQ(ierr); 2007736150eSBarry Smith *fnoise = neP->fnoise_min; 2017736150eSBarry Smith neP->fnoise_resets++; 2027736150eSBarry Smith } 2037736150eSBarry Smith if (*hopt < neP->hopt_min) { 2047736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %g, hopt_min = %g\n",*hopt,neP->hopt_min);CHKERRQ(ierr); 2057736150eSBarry Smith *hopt = neP->hopt_min; 2067736150eSBarry Smith neP->hopt_resets++; 2077736150eSBarry Smith } 2087736150eSBarry Smith 2097736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"Errors in derivative:\n");CHKERRQ(ierr); 2107736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"f = %g, fnoise = %g, fder2 = %g, hopt = %g\n",f,*fnoise,fder2,*hopt);CHKERRQ(ierr); 2117736150eSBarry Smith 2127736150eSBarry Smith /* For now, compute h **each** MV Mult!! */ 2137736150eSBarry Smith /* 2147736150eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-matrix_free_jorge_each_mvp",&flg);CHKERRQ(ierr); 2157736150eSBarry Smith if (!flg) { 2167736150eSBarry Smith Mat mat; 2177736150eSBarry Smith ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2187736150eSBarry Smith ierr = SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);CHKERRQ(ierr); 2197736150eSBarry Smith } 2207736150eSBarry Smith */ 2217736150eSBarry Smith fcount = neP->function_count - fcount; 2227736150eSBarry Smith PLogInfo(snes,"DiffParameterCompute_More: fct_now = %d, fct_cum = %d, rerrf=%g, sqrt(noise)=%g, h_more=%g\n", 2237736150eSBarry Smith fcount,neP->function_count,rerrf,sqrt(*fnoise),*hopt); 2247736150eSBarry Smith 2257736150eSBarry Smith 2267736150eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-noise_test",&noise_test);CHKERRQ(ierr); 2277736150eSBarry Smith if (noise_test) { 2287736150eSBarry Smith ierr = JacMatMultCompare(snes,x,p,*hopt);CHKERRQ(ierr); 2297736150eSBarry Smith } 2307736150eSBarry Smith PetscFunctionReturn(0); 2317736150eSBarry Smith } 2327736150eSBarry Smith 2337736150eSBarry Smith #undef __FUNC__ 2347736150eSBarry Smith #define __FUNC__ "JacMatMultCompare" 2357736150eSBarry Smith int JacMatMultCompare(SNES snes,Vec x,Vec p,double hopt) 2367736150eSBarry Smith { 2377736150eSBarry Smith Vec yy1, yy2; /* work vectors */ 2387736150eSBarry Smith Viewer view2; /* viewer */ 2397736150eSBarry Smith Mat J; /* analytic Jacobian (set as preconditioner matrix) */ 2407736150eSBarry Smith Mat Jmf; /* matrix-free Jacobian (set as true system matrix) */ 2417736150eSBarry Smith double h; /* differencing parameter */ 2427736150eSBarry Smith Vec f; 2437736150eSBarry Smith MatStructure sparsity = DIFFERENT_NONZERO_PATTERN; 2447736150eSBarry Smith Scalar alpha, yy1n, yy2n, enorm; 2457736150eSBarry Smith int i, ierr, printv; 2467736150eSBarry Smith char filename[32]; 2477736150eSBarry Smith MPI_Comm comm = snes->comm; 2487736150eSBarry Smith 2497736150eSBarry Smith PetscFunctionBegin; 2507736150eSBarry Smith 2517736150eSBarry Smith /* Compute function and analytic Jacobian at x */ 2527736150eSBarry Smith ierr = SNESGetJacobian(snes,&Jmf,&J,PETSC_NULL);CHKERRQ(ierr); 2537736150eSBarry Smith ierr = SNESComputeJacobian(snes,x,&Jmf,&J,&sparsity);CHKERRQ(ierr); 2547736150eSBarry Smith ierr = SNESGetFunction(snes,&f,PETSC_NULL);CHKERRQ(ierr); 2557736150eSBarry Smith ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); 2567736150eSBarry Smith 2577736150eSBarry Smith /* Duplicate work vectors */ 2587736150eSBarry Smith ierr = VecDuplicate(x,&yy2);CHKERRQ(ierr); 2597736150eSBarry Smith ierr = VecDuplicate(x,&yy1);CHKERRQ(ierr); 2607736150eSBarry Smith 2617736150eSBarry Smith /* Compute true matrix-vector product */ 2627736150eSBarry Smith ierr = MatMult(J,p,yy1);CHKERRQ(ierr); 2637736150eSBarry Smith ierr = VecNorm(yy1,NORM_2,&yy1n);CHKERRQ(ierr); 2647736150eSBarry Smith 2657736150eSBarry Smith /* View product vector if desired */ 2667736150eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-print_vecs",&printv);CHKERRQ(ierr); 2677736150eSBarry Smith if (printv) { 2687736150eSBarry Smith ierr = ViewerASCIIOpen(comm,"y1.out",&view2);CHKERRQ(ierr); 2697736150eSBarry Smith ierr = ViewerSetFormat(view2,VIEWER_FORMAT_ASCII_COMMON,PETSC_NULL);CHKERRQ(ierr); 2707736150eSBarry Smith ierr = VecView(yy1,view2);CHKERRQ(ierr); 2717736150eSBarry Smith ierr = ViewerDestroy(view2);CHKERRQ(ierr); 2727736150eSBarry Smith } 2737736150eSBarry Smith 2747736150eSBarry Smith /* Test Jacobian-vector product computation */ 2757736150eSBarry Smith alpha = -1.0; 2767736150eSBarry Smith h = 0.01 * hopt; 2777736150eSBarry Smith for (i=0; i<5; i++) { 2787736150eSBarry Smith /* Set differencing parameter for matrix-free multiplication */ 2797736150eSBarry Smith ierr = SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);CHKERRQ(ierr); 2807736150eSBarry Smith 2817736150eSBarry Smith /* Compute matrix-vector product via differencing approximation */ 2827736150eSBarry Smith ierr = MatMult(Jmf,p,yy2);CHKERRQ(ierr); 2837736150eSBarry Smith ierr = VecNorm(yy2,NORM_2,&yy2n);CHKERRQ(ierr); 2847736150eSBarry Smith 2857736150eSBarry Smith /* View product vector if desired */ 2867736150eSBarry Smith if (printv) { 2877736150eSBarry Smith sprintf(filename,"y2.%d.out",i); 2887736150eSBarry Smith ierr = ViewerASCIIOpen(comm,filename,&view2);CHKERRQ(ierr); 2897736150eSBarry Smith ierr = ViewerSetFormat(view2,VIEWER_FORMAT_ASCII_COMMON,PETSC_NULL);CHKERRQ(ierr); 2907736150eSBarry Smith ierr = VecView(yy2,view2);CHKERRQ(ierr); 2917736150eSBarry Smith ierr = ViewerDestroy(view2);CHKERRQ(ierr); 2927736150eSBarry Smith } 2937736150eSBarry Smith 2947736150eSBarry Smith /* Compute relative error */ 2957736150eSBarry Smith ierr = VecAXPY(&alpha,yy1,yy2);CHKERRQ(ierr); 2967736150eSBarry Smith ierr = VecNorm(yy2,NORM_2,&enorm);CHKERRQ(ierr); 2977736150eSBarry Smith enorm = enorm/yy1n; 2987736150eSBarry Smith ierr = PetscFPrintf(comm,stdout,"h = %g: relative error = %g\n",h,enorm);CHKERRQ(ierr); 2997736150eSBarry Smith h *= 10.0; 3007736150eSBarry Smith } 3017736150eSBarry Smith PetscFunctionReturn(0); 3027736150eSBarry Smith } 3037736150eSBarry Smith 3047736150eSBarry Smith static int lin_its_total = 0; 3057736150eSBarry Smith 3067736150eSBarry Smith int MyMonitor(SNES snes,int its,double fnorm,void *dummy) 3077736150eSBarry Smith { 3087736150eSBarry Smith int ierr, lin_its; 3097736150eSBarry Smith 3107736150eSBarry Smith PetscFunctionBegin; 3117736150eSBarry Smith ierr = SNESGetNumberLinearIterations(snes,&lin_its);CHKERRQ(ierr); 3127736150eSBarry Smith lin_its_total += lin_its; 3137736150eSBarry Smith 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); 3147736150eSBarry Smith 3157736150eSBarry Smith ierr = SNESUnSetMatrixFreeParameter(snes);CHKERRQ(ierr); 3167736150eSBarry Smith PetscFunctionReturn(0); 3177736150eSBarry Smith } 318