1051ac159SSatish Balay 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> 37736150eSBarry Smith 425acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCreate_More(SNES,Vec,void**); 525acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterCompute_More(SNES,void*,Vec,Vec,PetscReal*,PetscReal*); 625acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDiffParameterDestroy_More(void*); 725acbd8eSLisandro Dalcin 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 */ 125bd1e576SStefano Zampini PetscInt 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 */ 16a7cc72afSBarry Smith PetscInt fnoise_resets; /* number of times we've reset the noise estimate */ 17a7cc72afSBarry Smith PetscInt hopt_resets; /* number of times we've reset the hopt estimate */ 187736150eSBarry Smith } DIFFPAR_MORE; 197736150eSBarry Smith 2025acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat,double,double,double); 2125acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes); 2225acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESNoise_dnest_(PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*); 2325acbd8eSLisandro Dalcin 24ca8e9878SJed Brown static PetscErrorCode JacMatMultCompare(SNES,Vec,Vec,double); 257736150eSBarry Smith 268b5db460SBarry Smith PetscErrorCode SNESDiffParameterCreate_More(SNES snes,Vec x,void **outneP) 277736150eSBarry Smith { 287736150eSBarry Smith DIFFPAR_MORE *neP; 297736150eSBarry Smith Vec w; 307736150eSBarry Smith PetscRandom rctx; /* random number generator context */ 31ace3abfcSBarry Smith PetscBool flg; 32e2d1d2b7SBarry Smith char noise_file[PETSC_MAX_PATH_LEN]; 337736150eSBarry Smith 347736150eSBarry Smith PetscFunctionBegin; 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(snes,&neP)); 36f5af7f23SKarl Rupp 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 */ 455f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicateVecs(x,3,&neP->workv)); 467736150eSBarry Smith w = neP->workv[0]; 477736150eSBarry Smith 487736150eSBarry Smith /* Set components of vector w to random numbers */ 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject)snes),&rctx)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rctx)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(w,rctx)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rctx)); 537736150eSBarry Smith 547736150eSBarry Smith /* Open output file */ 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_noise_file",noise_file,sizeof(noise_file),&flg)); 567736150eSBarry Smith if (flg) neP->fp = fopen(noise_file,"w"); 577736150eSBarry Smith else neP->fp = fopen("noise.out","w"); 58*28b400f6SJacob Faibussowitsch PetscCheck(neP->fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file"); 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(snes,"Creating Jorge's differencing parameter context\n")); 607736150eSBarry Smith 617736150eSBarry Smith *outneP = neP; 627736150eSBarry Smith PetscFunctionReturn(0); 637736150eSBarry Smith } 647736150eSBarry Smith 658b5db460SBarry Smith PetscErrorCode SNESDiffParameterDestroy_More(void *nePv) 667736150eSBarry Smith { 677736150eSBarry Smith DIFFPAR_MORE *neP = (DIFFPAR_MORE*)nePv; 68ed9cf6e9SBarry Smith int err; 697736150eSBarry Smith 707736150eSBarry Smith PetscFunctionBegin; 717736150eSBarry Smith /* Destroy work vectors and close output file */ 725f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroyVecs(3,&neP->workv)); 73ed9cf6e9SBarry Smith err = fclose(neP->fp); 74*28b400f6SJacob Faibussowitsch PetscCheck(!err,PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(neP)); 767736150eSBarry Smith PetscFunctionReturn(0); 777736150eSBarry Smith } 787736150eSBarry Smith 798b5db460SBarry Smith PetscErrorCode SNESDiffParameterCompute_More(SNES snes,void *nePv,Vec x,Vec p,double *fnoise,double *hopt) 807736150eSBarry Smith { 817736150eSBarry Smith DIFFPAR_MORE *neP = (DIFFPAR_MORE*)nePv; 827736150eSBarry Smith Vec w, xp, fvec; /* work vectors to use in computing h */ 8390c26489SBarry Smith double zero = 0.0, hl, hu, h, fnoise_s, fder2_s; 8490c26489SBarry Smith PetscScalar alpha; 853964eb88SJed Brown PetscScalar fval[7], tab[7][7], eps[7], f = -1; 863964eb88SJed Brown double rerrf = -1., fder2; 8777431f27SBarry Smith PetscInt iter, k, i, j, info; 8877431f27SBarry Smith PetscInt nf = 7; /* number of function evaluations */ 8977431f27SBarry Smith PetscInt fcount; 90ce94432eSBarry Smith MPI_Comm comm; 917736150eSBarry Smith FILE *fp; 92ace3abfcSBarry Smith PetscBool noise_test = PETSC_FALSE; 937736150eSBarry Smith 947736150eSBarry Smith PetscFunctionBegin; 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)snes,&comm)); 967736150eSBarry Smith /* Call to SNESSetUp() just to set data structures in SNES context */ 975f80ce2aSJacob Faibussowitsch if (!snes->setupcalled) CHKERRQ(SNESSetUp(snes)); 987736150eSBarry Smith 997736150eSBarry Smith w = neP->workv[0]; 1007736150eSBarry Smith xp = neP->workv[1]; 1017736150eSBarry Smith fvec = neP->workv[2]; 1027736150eSBarry Smith fp = neP->fp; 1037736150eSBarry Smith 1047736150eSBarry Smith /* Initialize parameters */ 1057736150eSBarry Smith hl = zero; 1067736150eSBarry Smith hu = zero; 1077736150eSBarry Smith h = neP->h_first_try; 1087736150eSBarry Smith fnoise_s = zero; 1097736150eSBarry Smith fder2_s = zero; 1107736150eSBarry Smith fcount = neP->function_count; 1117736150eSBarry Smith 1127736150eSBarry Smith /* We have 5 tries to attempt to compute a good hopt value */ 1135f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&i)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp,"\n ------- SNES iteration %D ---------\n",i)); 1157736150eSBarry Smith for (iter=0; iter<5; iter++) { 1167736150eSBarry Smith neP->h_first_try = h; 1177736150eSBarry Smith 1187736150eSBarry Smith /* Compute the nf function values needed to estimate the noise from 1197736150eSBarry Smith the difference table */ 1207736150eSBarry Smith for (k=0; k<nf; k++) { 1217736150eSBarry Smith alpha = h * (k+1 - (nf+1)/2); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(xp,alpha,p,x)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeFunction(snes,xp,fvec)); 1247736150eSBarry Smith neP->function_count++; 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecDot(fvec,w,&fval[k])); 1267736150eSBarry Smith } 1277736150eSBarry Smith f = fval[(nf+1)/2 - 1]; 1287736150eSBarry Smith 1297736150eSBarry Smith /* Construct the difference table */ 130f5af7f23SKarl Rupp for (i=0; i<nf; i++) tab[i][0] = fval[i]; 131f5af7f23SKarl Rupp 132e108cb99SStefano Zampini for (j=0; j<nf-1; j++) { 133e108cb99SStefano Zampini for (i=0; i<nf-j-1; i++) { 1347736150eSBarry Smith tab[i][j+1] = tab[i+1][j] - tab[i][j]; 1357736150eSBarry Smith } 1367736150eSBarry Smith } 1377736150eSBarry Smith 1387736150eSBarry Smith /* Print the difference table */ 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp,"Difference Table: iter = %D\n",iter)); 1407736150eSBarry Smith for (i=0; i<nf; i++) { 1417736150eSBarry Smith for (j=0; j<nf-i; j++) { 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp," %10.2e ",tab[i][j])); 1437736150eSBarry Smith } 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp,"\n")); 1457736150eSBarry Smith } 1467736150eSBarry Smith 1477736150eSBarry Smith /* Call the noise estimator */ 1485f80ce2aSJacob Faibussowitsch CHKERRQ(SNESNoise_dnest_(&nf,fval,&h,fnoise,&fder2,hopt,&info,eps)); 1497736150eSBarry Smith 1507736150eSBarry Smith /* Output statements */ 1517736150eSBarry Smith rerrf = *fnoise/PetscAbsScalar(f); 1525f80ce2aSJacob Faibussowitsch if (info == 1) CHKERRQ(PetscFPrintf(comm,fp,"%s\n","Noise detected")); 1535f80ce2aSJacob Faibussowitsch if (info == 2) {CHKERRQ(PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too small"));} 1545f80ce2aSJacob Faibussowitsch if (info == 3) {CHKERRQ(PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too large"));} 1555f80ce2aSJacob Faibussowitsch if (info == 4) CHKERRQ(PetscFPrintf(comm,fp,"%s\n","Noise detected, but unreliable hopt")); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp,"Approximate epsfcn %g %g %g %g %g %g\n",(double)eps[0],(double)eps[1],(double)eps[2],(double)eps[3],(double)eps[4],(double)eps[5])); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp,"h = %g, fnoise = %g, fder2 = %g, rerrf = %g, hopt = %g\n\n",(double)h, (double)*fnoise, (double)fder2, (double)rerrf, (double)*hopt)); 1587736150eSBarry Smith 1597736150eSBarry Smith /* Save fnoise and fder2. */ 1607736150eSBarry Smith if (*fnoise) fnoise_s = *fnoise; 1617736150eSBarry Smith if (fder2) fder2_s = fder2; 1627736150eSBarry Smith 1637736150eSBarry Smith /* Check for noise detection. */ 1647736150eSBarry Smith if (fnoise_s && fder2_s) { 1657736150eSBarry Smith *fnoise = fnoise_s; 1667736150eSBarry Smith fder2 = fder2_s; 1677736150eSBarry Smith *hopt = 1.68*sqrt(*fnoise/PetscAbsScalar(fder2)); 1687736150eSBarry Smith goto theend; 1697736150eSBarry Smith } else { 1707736150eSBarry Smith 1717736150eSBarry Smith /* Update hl and hu, and determine new h */ 1727736150eSBarry Smith if (info == 2 || info == 4) { 1737736150eSBarry Smith hl = h; 1747736150eSBarry Smith if (hu == zero) h = 100*h; 1757736150eSBarry Smith else h = PetscMin(100*h,0.1*hu); 1767736150eSBarry Smith } else if (info == 3) { 1777736150eSBarry Smith hu = h; 1787736150eSBarry Smith h = PetscMax(1.0e-3,sqrt(hl/hu))*hu; 1797736150eSBarry Smith } 1807736150eSBarry Smith } 1817736150eSBarry Smith } 1827736150eSBarry Smith theend: 1837736150eSBarry Smith 1847736150eSBarry Smith if (*fnoise < neP->fnoise_min) { 1855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %g, fnoise_min = %g\n",(double)*fnoise,(double)neP->fnoise_min)); 1867736150eSBarry Smith *fnoise = neP->fnoise_min; 1877736150eSBarry Smith neP->fnoise_resets++; 1887736150eSBarry Smith } 1897736150eSBarry Smith if (*hopt < neP->hopt_min) { 1905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %g, hopt_min = %g\n",(double)*hopt,(double)neP->hopt_min)); 1917736150eSBarry Smith *hopt = neP->hopt_min; 1927736150eSBarry Smith neP->hopt_resets++; 1937736150eSBarry Smith } 1947736150eSBarry Smith 1955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp,"Errors in derivative:\n")); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,fp,"f = %g, fnoise = %g, fder2 = %g, hopt = %g\n",(double)f,(double)*fnoise,(double)fder2,(double)*hopt)); 1977736150eSBarry Smith 1987736150eSBarry Smith /* For now, compute h **each** MV Mult!! */ 1997736150eSBarry Smith /* 2005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,"-matrix_free_jorge_each_mvp",&flg)); 2017736150eSBarry Smith if (!flg) { 2027736150eSBarry Smith Mat mat; 2035f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetJacobian(snes,&mat,NULL,NULL)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt)); 2057736150eSBarry Smith } 2067736150eSBarry Smith */ 2077736150eSBarry Smith fcount = neP->function_count - fcount; 2085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(snes,"fct_now = %D, fct_cum = %D, rerrf=%g, sqrt(noise)=%g, h_more=%g\n",fcount,neP->function_count,(double)rerrf,(double)PetscSqrtReal(*fnoise),(double)*hopt)); 2097736150eSBarry Smith 2105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-noise_test",&noise_test,NULL)); 2117736150eSBarry Smith if (noise_test) { 2125f80ce2aSJacob Faibussowitsch CHKERRQ(JacMatMultCompare(snes,x,p,*hopt)); 2137736150eSBarry Smith } 2147736150eSBarry Smith PetscFunctionReturn(0); 2157736150eSBarry Smith } 2167736150eSBarry Smith 217dfbe8321SBarry Smith PetscErrorCode JacMatMultCompare(SNES snes,Vec x,Vec p,double hopt) 2187736150eSBarry Smith { 2197736150eSBarry Smith Vec yy1, yy2; /* work vectors */ 22052ed857cSBarry Smith PetscViewer view2; /* viewer */ 2217736150eSBarry Smith Mat J; /* analytic Jacobian (set as preconditioner matrix) */ 2227736150eSBarry Smith Mat Jmf; /* matrix-free Jacobian (set as true system matrix) */ 2237736150eSBarry Smith double h; /* differencing parameter */ 2247736150eSBarry Smith Vec f; 22590c26489SBarry Smith PetscScalar alpha; 22690c26489SBarry Smith PetscReal yy1n,yy2n,enorm; 227a7cc72afSBarry Smith PetscInt i; 228ace3abfcSBarry Smith PetscBool printv = PETSC_FALSE; 2297736150eSBarry Smith char filename[32]; 230ce94432eSBarry Smith MPI_Comm comm; 2317736150eSBarry Smith 2327736150eSBarry Smith PetscFunctionBegin; 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)snes,&comm)); 2347736150eSBarry Smith /* Compute function and analytic Jacobian at x */ 2355f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetJacobian(snes,&Jmf,&J,NULL,NULL)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeJacobian(snes,x,Jmf,J)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes,&f,NULL,NULL)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(SNESComputeFunction(snes,x,f)); 2397736150eSBarry Smith 2407736150eSBarry Smith /* Duplicate work vectors */ 2415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&yy2)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&yy1)); 2437736150eSBarry Smith 2447736150eSBarry Smith /* Compute true matrix-vector product */ 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(J,p,yy1)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(yy1,NORM_2,&yy1n)); 2477736150eSBarry Smith 2487736150eSBarry Smith /* View product vector if desired */ 2495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-print_vecs",&printv,NULL)); 2507736150eSBarry Smith if (printv) { 2515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(comm,"y1.out",&view2)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(view2,PETSC_VIEWER_ASCII_COMMON)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(yy1,view2)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(view2)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view2)); 2567736150eSBarry Smith } 2577736150eSBarry Smith 2587736150eSBarry Smith /* Test Jacobian-vector product computation */ 2597736150eSBarry Smith alpha = -1.0; 2607736150eSBarry Smith h = 0.01 * hopt; 2617736150eSBarry Smith for (i=0; i<5; i++) { 2627736150eSBarry Smith /* Set differencing parameter for matrix-free multiplication */ 2635f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h)); 2647736150eSBarry Smith 2657736150eSBarry Smith /* Compute matrix-vector product via differencing approximation */ 2665f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Jmf,p,yy2)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(yy2,NORM_2,&yy2n)); 2687736150eSBarry Smith 2697736150eSBarry Smith /* View product vector if desired */ 2707736150eSBarry Smith if (printv) { 27177431f27SBarry Smith sprintf(filename,"y2.%d.out",(int)i); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIOpen(comm,filename,&view2)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(view2,PETSC_VIEWER_ASCII_COMMON)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(yy2,view2)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(view2)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view2)); 2777736150eSBarry Smith } 2787736150eSBarry Smith 2797736150eSBarry Smith /* Compute relative error */ 2805f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(yy2,alpha,yy1)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(yy2,NORM_2,&enorm)); 2827736150eSBarry Smith enorm = enorm/yy1n; 2835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,stdout,"h = %g: relative error = %g\n",(double)h,(double)enorm)); 2847736150eSBarry Smith h *= 10.0; 2857736150eSBarry Smith } 2867736150eSBarry Smith PetscFunctionReturn(0); 2877736150eSBarry Smith } 2887736150eSBarry Smith 289a7cc72afSBarry Smith static PetscInt lin_its_total = 0; 2907736150eSBarry Smith 291b12ffc70SBarry Smith PetscErrorCode SNESNoiseMonitor(SNES snes,PetscInt its,double fnorm,void *dummy) 2927736150eSBarry Smith { 29377431f27SBarry Smith PetscInt lin_its; 2947736150eSBarry Smith 2957736150eSBarry Smith PetscFunctionBegin; 2965f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetLinearSolveIterations(snes,&lin_its)); 2977736150eSBarry Smith lin_its_total += lin_its; 2985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)snes), "iter = %D, SNES Function norm = %g, lin_its = %D, total_lin_its = %D\n",its,(double)fnorm,lin_its,lin_its_total)); 2997736150eSBarry Smith 3005f80ce2aSJacob Faibussowitsch CHKERRQ(SNESUnSetMatrixFreeParameter(snes)); 3017736150eSBarry Smith PetscFunctionReturn(0); 3027736150eSBarry Smith } 303