163dd3a1aSKris Buschelman #define PETSCSNES_DLL 2051ac159SSatish Balay 3b9147fbbSdalcinl #include "include/private/snesimpl.h" 47736150eSBarry Smith 57736150eSBarry Smith /* Data used by Jorge's diff parameter computation method */ 67736150eSBarry Smith typedef struct { 77736150eSBarry Smith Vec *workv; /* work vectors */ 87736150eSBarry Smith FILE *fp; /* output file */ 97736150eSBarry Smith int function_count; /* count of function evaluations for diff param estimation */ 107736150eSBarry Smith double fnoise_min; /* minimim allowable noise */ 117736150eSBarry Smith double hopt_min; /* minimum allowable hopt */ 127736150eSBarry Smith double h_first_try; /* first try for h used in diff parameter estimate */ 13a7cc72afSBarry Smith PetscInt fnoise_resets; /* number of times we've reset the noise estimate */ 14a7cc72afSBarry Smith PetscInt hopt_resets; /* number of times we've reset the hopt estimate */ 157736150eSBarry Smith } DIFFPAR_MORE; 167736150eSBarry Smith 177736150eSBarry Smith 18a7cc72afSBarry Smith extern PetscErrorCode dnest_(PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*); 196849ba73SBarry Smith EXTERN PetscErrorCode JacMatMultCompare(SNES,Vec,Vec,double); 20dfbe8321SBarry Smith EXTERN PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat,double,double,double); 21dfbe8321SBarry Smith EXTERN PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes); 227736150eSBarry Smith 234a2ae208SSatish Balay #undef __FUNCT__ 244a2ae208SSatish Balay #define __FUNCT__ "DiffParameterCreate_More" 25dfbe8321SBarry Smith PetscErrorCode DiffParameterCreate_More(SNES snes,Vec x,void **outneP) 267736150eSBarry Smith { 277736150eSBarry Smith DIFFPAR_MORE *neP; 287736150eSBarry Smith Vec w; 297736150eSBarry Smith PetscRandom rctx; /* random number generator context */ 30dfbe8321SBarry Smith PetscErrorCode ierr; 3152ed857cSBarry Smith PetscTruth flg; 32e2d1d2b7SBarry Smith char noise_file[PETSC_MAX_PATH_LEN]; 337736150eSBarry Smith 347736150eSBarry Smith PetscFunctionBegin; 357736150eSBarry Smith 36*38f2d2fdSLisandro Dalcin ierr = PetscNewLog(snes,DIFFPAR_MORE,&neP);CHKERRQ(ierr); 377736150eSBarry Smith neP->function_count = 0; 387736150eSBarry Smith neP->fnoise_min = 1.0e-20; 397736150eSBarry Smith neP->hopt_min = 1.0e-8; 407736150eSBarry Smith neP->h_first_try = 1.0e-3; 417736150eSBarry Smith neP->fnoise_resets = 0; 427736150eSBarry Smith neP->hopt_resets = 0; 437736150eSBarry Smith 447736150eSBarry Smith /* Create work vectors */ 457736150eSBarry Smith ierr = VecDuplicateVecs(x,3,&neP->workv);CHKERRQ(ierr); 467736150eSBarry Smith w = neP->workv[0]; 477736150eSBarry Smith 487736150eSBarry Smith /* Set components of vector w to random numbers */ 49c77d6671SHong Zhang ierr = PetscRandomCreate(snes->comm,&rctx);CHKERRQ(ierr); 50c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 51abb0e124SMatthew Knepley ierr = VecSetRandom(w,rctx);CHKERRQ(ierr); 527736150eSBarry Smith ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr); 537736150eSBarry Smith 547736150eSBarry Smith /* Open output file */ 55e2d1d2b7SBarry Smith ierr = PetscOptionsGetString(snes->prefix,"-snes_mf_noise_file",noise_file,PETSC_MAX_PATH_LEN-1,&flg);CHKERRQ(ierr); 567736150eSBarry Smith if (flg) neP->fp = fopen(noise_file,"w"); 577736150eSBarry Smith else neP->fp = fopen("noise.out","w"); 5852ed857cSBarry Smith if (!neP->fp) SETERRQ(PETSC_ERR_FILE_OPEN,"Cannot open file"); 59ae15b995SBarry Smith ierr = PetscInfo(snes,"Creating Jorge's differencing parameter context\n");CHKERRQ(ierr); 607736150eSBarry Smith 617736150eSBarry Smith *outneP = neP; 627736150eSBarry Smith PetscFunctionReturn(0); 637736150eSBarry Smith } 647736150eSBarry Smith 654a2ae208SSatish Balay #undef __FUNCT__ 664a2ae208SSatish Balay #define __FUNCT__ "DiffParameterDestroy_More" 67dfbe8321SBarry Smith PetscErrorCode DiffParameterDestroy_More(void *nePv) 687736150eSBarry Smith { 697736150eSBarry Smith DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv; 70dfbe8321SBarry Smith PetscErrorCode ierr; 717736150eSBarry Smith 727736150eSBarry Smith PetscFunctionBegin; 737736150eSBarry Smith /* Destroy work vectors and close output file */ 747736150eSBarry Smith ierr = VecDestroyVecs(neP->workv,3);CHKERRQ(ierr); 757736150eSBarry Smith fclose(neP->fp); 767736150eSBarry Smith ierr = PetscFree(neP);CHKERRQ(ierr); 777736150eSBarry Smith PetscFunctionReturn(0); 787736150eSBarry Smith } 797736150eSBarry Smith 804a2ae208SSatish Balay #undef __FUNCT__ 814a2ae208SSatish Balay #define __FUNCT__ "DiffParameterCompute_More" 82dfbe8321SBarry Smith PetscErrorCode DiffParameterCompute_More(SNES snes,void *nePv,Vec x,Vec p,double *fnoise,double *hopt) 837736150eSBarry Smith { 847736150eSBarry Smith DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv; 857736150eSBarry Smith Vec w, xp, fvec; /* work vectors to use in computing h */ 8690c26489SBarry Smith double zero = 0.0, hl, hu, h, fnoise_s, fder2_s; 8790c26489SBarry Smith PetscScalar alpha; 8890c26489SBarry Smith PetscScalar fval[7], tab[7][7], eps[7], f; 8990c26489SBarry Smith double rerrf, fder2; 906849ba73SBarry Smith PetscErrorCode ierr; 9177431f27SBarry Smith PetscInt iter, k, i, j, info; 9277431f27SBarry Smith PetscInt nf = 7; /* number of function evaluations */ 9377431f27SBarry Smith PetscInt fcount; 947736150eSBarry Smith MPI_Comm comm = snes->comm; 957736150eSBarry Smith FILE *fp; 9652ed857cSBarry Smith PetscTruth noise_test; 977736150eSBarry Smith 987736150eSBarry Smith PetscFunctionBegin; 997736150eSBarry Smith /* Call to SNESSetUp() just to set data structures in SNES context */ 10070e92668SMatthew Knepley if (!snes->setupcalled) {ierr = SNESSetUp(snes);CHKERRQ(ierr);} 1017736150eSBarry Smith 1027736150eSBarry Smith w = neP->workv[0]; 1037736150eSBarry Smith xp = neP->workv[1]; 1047736150eSBarry Smith fvec = neP->workv[2]; 1057736150eSBarry Smith fp = neP->fp; 1067736150eSBarry Smith 1077736150eSBarry Smith /* Initialize parameters */ 1087736150eSBarry Smith hl = zero; 1097736150eSBarry Smith hu = zero; 1107736150eSBarry Smith h = neP->h_first_try; 1117736150eSBarry Smith fnoise_s = zero; 1127736150eSBarry Smith fder2_s = zero; 1137736150eSBarry Smith fcount = neP->function_count; 1147736150eSBarry Smith 1157736150eSBarry Smith /* We have 5 tries to attempt to compute a good hopt value */ 1167736150eSBarry Smith ierr = SNESGetIterationNumber(snes,&i);CHKERRQ(ierr); 11777431f27SBarry Smith ierr = PetscFPrintf(comm,fp,"\n ------- SNES iteration %D ---------\n",i);CHKERRQ(ierr); 1187736150eSBarry Smith for (iter=0; iter<5; iter++) { 1197736150eSBarry Smith neP->h_first_try = h; 1207736150eSBarry Smith 1217736150eSBarry Smith /* Compute the nf function values needed to estimate the noise from 1227736150eSBarry Smith the difference table */ 1237736150eSBarry Smith for (k=0; k<nf; k++) { 1247736150eSBarry Smith alpha = h * ( k+1 - (nf+1)/2 ); 1252dcb1b2aSMatthew Knepley ierr = VecWAXPY(xp,alpha,p,x);CHKERRQ(ierr); 1267736150eSBarry Smith ierr = SNESComputeFunction(snes,xp,fvec);CHKERRQ(ierr); 1277736150eSBarry Smith neP->function_count++; 1287736150eSBarry Smith ierr = VecDot(fvec,w,&fval[k]);CHKERRQ(ierr); 1297736150eSBarry Smith } 1307736150eSBarry Smith f = fval[(nf+1)/2 - 1]; 1317736150eSBarry Smith 1327736150eSBarry Smith /* Construct the difference table */ 1337736150eSBarry Smith for (i=0; i<nf; i++) { 1347736150eSBarry Smith tab[i][0] = fval[i]; 1357736150eSBarry Smith } 1367736150eSBarry Smith for (j=0; j<6; j++) { 1377736150eSBarry Smith for (i=0; i<nf-j; i++) { 1387736150eSBarry Smith tab[i][j+1] = tab[i+1][j] - tab[i][j]; 1397736150eSBarry Smith } 1407736150eSBarry Smith } 1417736150eSBarry Smith 1427736150eSBarry Smith /* Print the difference table */ 14377431f27SBarry Smith ierr = PetscFPrintf(comm,fp,"Difference Table: iter = %D\n",iter);CHKERRQ(ierr); 1447736150eSBarry Smith for (i=0; i<nf; i++) { 1457736150eSBarry Smith for (j=0; j<nf-i; j++) { 1467736150eSBarry Smith ierr = PetscFPrintf(comm,fp," %10.2e ",tab[i][j]);CHKERRQ(ierr); 1477736150eSBarry Smith } 1487736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"\n");CHKERRQ(ierr); 1497736150eSBarry Smith } 1507736150eSBarry Smith 1517736150eSBarry Smith /* Call the noise estimator */ 152f6428e95SBarry Smith ierr = dnest_(&nf,fval,&h,fnoise,&fder2,hopt,&info,eps);CHKERRQ(ierr); 1537736150eSBarry Smith 1547736150eSBarry Smith /* Output statements */ 1557736150eSBarry Smith rerrf = *fnoise/PetscAbsScalar(f); 1567736150eSBarry Smith if (info == 1) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected");CHKERRQ(ierr);} 1577736150eSBarry Smith if (info == 2) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too small");CHKERRQ(ierr);} 1587736150eSBarry Smith if (info == 3) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too large");CHKERRQ(ierr);} 1597736150eSBarry Smith if (info == 4) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected, but unreliable hopt");CHKERRQ(ierr);} 160a83599f4SBarry Smith ierr = PetscFPrintf(comm,fp,"Approximate epsfcn %G %G %G %G %G %G\n", 1617736150eSBarry Smith eps[0],eps[1],eps[2],eps[3],eps[4],eps[5]);CHKERRQ(ierr); 162a83599f4SBarry Smith ierr = PetscFPrintf(comm,fp,"h = %G, fnoise = %G, fder2 = %G, rerrf = %G, hopt = %G\n\n", 1637736150eSBarry Smith h, *fnoise, fder2, rerrf, *hopt);CHKERRQ(ierr); 1647736150eSBarry Smith 1657736150eSBarry Smith /* Save fnoise and fder2. */ 1667736150eSBarry Smith if (*fnoise) fnoise_s = *fnoise; 1677736150eSBarry Smith if (fder2) fder2_s = fder2; 1687736150eSBarry Smith 1697736150eSBarry Smith /* Check for noise detection. */ 1707736150eSBarry Smith if (fnoise_s && fder2_s) { 1717736150eSBarry Smith *fnoise = fnoise_s; 1727736150eSBarry Smith fder2 = fder2_s; 1737736150eSBarry Smith *hopt = 1.68*sqrt(*fnoise/PetscAbsScalar(fder2)); 1747736150eSBarry Smith goto theend; 1757736150eSBarry Smith } else { 1767736150eSBarry Smith 1777736150eSBarry Smith /* Update hl and hu, and determine new h */ 1787736150eSBarry Smith if (info == 2 || info == 4) { 1797736150eSBarry Smith hl = h; 1807736150eSBarry Smith if (hu == zero) h = 100*h; 1817736150eSBarry Smith else h = PetscMin(100*h,0.1*hu); 1827736150eSBarry Smith } else if (info == 3) { 1837736150eSBarry Smith hu = h; 1847736150eSBarry Smith h = PetscMax(1.0e-3,sqrt(hl/hu))*hu; 1857736150eSBarry Smith } 1867736150eSBarry Smith } 1877736150eSBarry Smith } 1887736150eSBarry Smith theend: 1897736150eSBarry Smith 1907736150eSBarry Smith if (*fnoise < neP->fnoise_min) { 191a83599f4SBarry Smith ierr = PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %G, fnoise_min = %G\n",*fnoise,neP->fnoise_min);CHKERRQ(ierr); 1927736150eSBarry Smith *fnoise = neP->fnoise_min; 1937736150eSBarry Smith neP->fnoise_resets++; 1947736150eSBarry Smith } 1957736150eSBarry Smith if (*hopt < neP->hopt_min) { 196a83599f4SBarry Smith ierr = PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %G, hopt_min = %G\n",*hopt,neP->hopt_min);CHKERRQ(ierr); 1977736150eSBarry Smith *hopt = neP->hopt_min; 1987736150eSBarry Smith neP->hopt_resets++; 1997736150eSBarry Smith } 2007736150eSBarry Smith 2017736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"Errors in derivative:\n");CHKERRQ(ierr); 202a83599f4SBarry Smith ierr = PetscFPrintf(comm,fp,"f = %G, fnoise = %G, fder2 = %G, hopt = %G\n",f,*fnoise,fder2,*hopt);CHKERRQ(ierr); 2037736150eSBarry Smith 2047736150eSBarry Smith /* For now, compute h **each** MV Mult!! */ 2057736150eSBarry Smith /* 20652ed857cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-matrix_free_jorge_each_mvp",&flg);CHKERRQ(ierr); 2077736150eSBarry Smith if (!flg) { 2087736150eSBarry Smith Mat mat; 2097736150eSBarry Smith ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2107736150eSBarry Smith ierr = SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);CHKERRQ(ierr); 2117736150eSBarry Smith } 2127736150eSBarry Smith */ 2137736150eSBarry Smith fcount = neP->function_count - fcount; 214ae15b995SBarry Smith 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); 2157736150eSBarry Smith 21652ed857cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-noise_test",&noise_test);CHKERRQ(ierr); 2177736150eSBarry Smith if (noise_test) { 2187736150eSBarry Smith ierr = JacMatMultCompare(snes,x,p,*hopt);CHKERRQ(ierr); 2197736150eSBarry Smith } 2207736150eSBarry Smith PetscFunctionReturn(0); 2217736150eSBarry Smith } 2227736150eSBarry Smith 2234a2ae208SSatish Balay #undef __FUNCT__ 2244a2ae208SSatish Balay #define __FUNCT__ "JacMatMultCompare" 225dfbe8321SBarry Smith PetscErrorCode JacMatMultCompare(SNES snes,Vec x,Vec p,double hopt) 2267736150eSBarry Smith { 2277736150eSBarry Smith Vec yy1, yy2; /* work vectors */ 22852ed857cSBarry Smith PetscViewer view2; /* viewer */ 2297736150eSBarry Smith Mat J; /* analytic Jacobian (set as preconditioner matrix) */ 2307736150eSBarry Smith Mat Jmf; /* matrix-free Jacobian (set as true system matrix) */ 2317736150eSBarry Smith double h; /* differencing parameter */ 2327736150eSBarry Smith Vec f; 2337736150eSBarry Smith MatStructure sparsity = DIFFERENT_NONZERO_PATTERN; 23490c26489SBarry Smith PetscScalar alpha; 23590c26489SBarry Smith PetscReal yy1n,yy2n,enorm; 2366849ba73SBarry Smith PetscErrorCode ierr; 237a7cc72afSBarry Smith PetscInt i; 23852ed857cSBarry Smith PetscTruth printv; 2397736150eSBarry Smith char filename[32]; 2407736150eSBarry Smith MPI_Comm comm = snes->comm; 2417736150eSBarry Smith 2427736150eSBarry Smith PetscFunctionBegin; 2437736150eSBarry Smith 2447736150eSBarry Smith /* Compute function and analytic Jacobian at x */ 24552ed857cSBarry Smith ierr = SNESGetJacobian(snes,&Jmf,&J,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2467736150eSBarry Smith ierr = SNESComputeJacobian(snes,x,&Jmf,&J,&sparsity);CHKERRQ(ierr); 24752ed857cSBarry Smith ierr = SNESGetFunction(snes,&f,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2487736150eSBarry Smith ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); 2497736150eSBarry Smith 2507736150eSBarry Smith /* Duplicate work vectors */ 2517736150eSBarry Smith ierr = VecDuplicate(x,&yy2);CHKERRQ(ierr); 2527736150eSBarry Smith ierr = VecDuplicate(x,&yy1);CHKERRQ(ierr); 2537736150eSBarry Smith 2547736150eSBarry Smith /* Compute true matrix-vector product */ 2557736150eSBarry Smith ierr = MatMult(J,p,yy1);CHKERRQ(ierr); 2567736150eSBarry Smith ierr = VecNorm(yy1,NORM_2,&yy1n);CHKERRQ(ierr); 2577736150eSBarry Smith 2587736150eSBarry Smith /* View product vector if desired */ 25952ed857cSBarry Smith ierr = PetscOptionsHasName(PETSC_NULL,"-print_vecs",&printv);CHKERRQ(ierr); 2607736150eSBarry Smith if (printv) { 26152ed857cSBarry Smith ierr = PetscViewerASCIIOpen(comm,"y1.out",&view2);CHKERRQ(ierr); 26252ed857cSBarry Smith ierr = PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr); 2637736150eSBarry Smith ierr = VecView(yy1,view2);CHKERRQ(ierr); 26452ed857cSBarry Smith ierr = PetscViewerDestroy(view2);CHKERRQ(ierr); 2657736150eSBarry Smith } 2667736150eSBarry Smith 2677736150eSBarry Smith /* Test Jacobian-vector product computation */ 2687736150eSBarry Smith alpha = -1.0; 2697736150eSBarry Smith h = 0.01 * hopt; 2707736150eSBarry Smith for (i=0; i<5; i++) { 2717736150eSBarry Smith /* Set differencing parameter for matrix-free multiplication */ 2727736150eSBarry Smith ierr = SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);CHKERRQ(ierr); 2737736150eSBarry Smith 2747736150eSBarry Smith /* Compute matrix-vector product via differencing approximation */ 2757736150eSBarry Smith ierr = MatMult(Jmf,p,yy2);CHKERRQ(ierr); 2767736150eSBarry Smith ierr = VecNorm(yy2,NORM_2,&yy2n);CHKERRQ(ierr); 2777736150eSBarry Smith 2787736150eSBarry Smith /* View product vector if desired */ 2797736150eSBarry Smith if (printv) { 28077431f27SBarry Smith sprintf(filename,"y2.%d.out",(int)i); 28152ed857cSBarry Smith ierr = PetscViewerASCIIOpen(comm,filename,&view2);CHKERRQ(ierr); 28252ed857cSBarry Smith ierr = PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr); 2837736150eSBarry Smith ierr = VecView(yy2,view2);CHKERRQ(ierr); 28452ed857cSBarry Smith ierr = PetscViewerDestroy(view2);CHKERRQ(ierr); 2857736150eSBarry Smith } 2867736150eSBarry Smith 2877736150eSBarry Smith /* Compute relative error */ 2882dcb1b2aSMatthew Knepley ierr = VecAXPY(yy2,alpha,yy1);CHKERRQ(ierr); 2897736150eSBarry Smith ierr = VecNorm(yy2,NORM_2,&enorm);CHKERRQ(ierr); 2907736150eSBarry Smith enorm = enorm/yy1n; 291a83599f4SBarry Smith ierr = PetscFPrintf(comm,stdout,"h = %G: relative error = %G\n",h,enorm);CHKERRQ(ierr); 2927736150eSBarry Smith h *= 10.0; 2937736150eSBarry Smith } 2947736150eSBarry Smith PetscFunctionReturn(0); 2957736150eSBarry Smith } 2967736150eSBarry Smith 297a7cc72afSBarry Smith static PetscInt lin_its_total = 0; 2987736150eSBarry Smith 29977431f27SBarry Smith PetscErrorCode MyMonitor(SNES snes,PetscInt its,double fnorm,void *dummy) 3007736150eSBarry Smith { 3016849ba73SBarry Smith PetscErrorCode ierr; 30277431f27SBarry Smith PetscInt lin_its; 3037736150eSBarry Smith 3047736150eSBarry Smith PetscFunctionBegin; 305b850b91aSLisandro Dalcin ierr = SNESGetLinearSolveIterations(snes,&lin_its);CHKERRQ(ierr); 3067736150eSBarry Smith lin_its_total += lin_its; 307a83599f4SBarry 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); 3087736150eSBarry Smith 3097736150eSBarry Smith ierr = SNESUnSetMatrixFreeParameter(snes);CHKERRQ(ierr); 3107736150eSBarry Smith PetscFunctionReturn(0); 3117736150eSBarry Smith } 312