xref: /petsc/src/snes/interface/noise/snesnoise.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
269371c9d4SSatish Balay PetscErrorCode SNESDiffParameterCreate_More(SNES snes, Vec x, void **outneP) {
277736150eSBarry Smith   DIFFPAR_MORE *neP;
287736150eSBarry Smith   Vec           w;
297736150eSBarry Smith   PetscRandom   rctx; /* random number generator context */
30ace3abfcSBarry Smith   PetscBool     flg;
31e2d1d2b7SBarry Smith   char          noise_file[PETSC_MAX_PATH_LEN];
327736150eSBarry Smith 
337736150eSBarry Smith   PetscFunctionBegin;
349566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &neP));
35f5af7f23SKarl Rupp 
367736150eSBarry Smith   neP->function_count = 0;
377736150eSBarry Smith   neP->fnoise_min     = 1.0e-20;
387736150eSBarry Smith   neP->hopt_min       = 1.0e-8;
397736150eSBarry Smith   neP->h_first_try    = 1.0e-3;
407736150eSBarry Smith   neP->fnoise_resets  = 0;
417736150eSBarry Smith   neP->hopt_resets    = 0;
427736150eSBarry Smith 
437736150eSBarry Smith   /* Create work vectors */
449566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(x, 3, &neP->workv));
457736150eSBarry Smith   w = neP->workv[0];
467736150eSBarry Smith 
477736150eSBarry Smith   /* Set components of vector w to random numbers */
489566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)snes), &rctx));
499566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rctx));
509566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(w, rctx));
519566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rctx));
527736150eSBarry Smith 
537736150eSBarry Smith   /* Open output file */
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_mf_noise_file", noise_file, sizeof(noise_file), &flg));
557736150eSBarry Smith   if (flg) neP->fp = fopen(noise_file, "w");
567736150eSBarry Smith   else neP->fp = fopen("noise.out", "w");
5728b400f6SJacob Faibussowitsch   PetscCheck(neP->fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file");
589566063dSJacob Faibussowitsch   PetscCall(PetscInfo(snes, "Creating Jorge's differencing parameter context\n"));
597736150eSBarry Smith 
607736150eSBarry Smith   *outneP = neP;
617736150eSBarry Smith   PetscFunctionReturn(0);
627736150eSBarry Smith }
637736150eSBarry Smith 
649371c9d4SSatish Balay PetscErrorCode SNESDiffParameterDestroy_More(void *nePv) {
657736150eSBarry Smith   DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv;
66ed9cf6e9SBarry Smith   int           err;
677736150eSBarry Smith 
687736150eSBarry Smith   PetscFunctionBegin;
697736150eSBarry Smith   /* Destroy work vectors and close output file */
709566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(3, &neP->workv));
71ed9cf6e9SBarry Smith   err = fclose(neP->fp);
7228b400f6SJacob Faibussowitsch   PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_SYS, "fclose() failed on file");
739566063dSJacob Faibussowitsch   PetscCall(PetscFree(neP));
747736150eSBarry Smith   PetscFunctionReturn(0);
757736150eSBarry Smith }
767736150eSBarry Smith 
779371c9d4SSatish Balay PetscErrorCode SNESDiffParameterCompute_More(SNES snes, void *nePv, Vec x, Vec p, double *fnoise, double *hopt) {
787736150eSBarry Smith   DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv;
797736150eSBarry Smith   Vec           w, xp, fvec; /* work vectors to use in computing h */
8090c26489SBarry Smith   double        zero = 0.0, hl, hu, h, fnoise_s, fder2_s;
8190c26489SBarry Smith   PetscScalar   alpha;
823964eb88SJed Brown   PetscScalar   fval[7], tab[7][7], eps[7], f = -1;
833964eb88SJed Brown   double        rerrf = -1., fder2;
8477431f27SBarry Smith   PetscInt      iter, k, i, j, info;
8577431f27SBarry Smith   PetscInt      nf = 7; /* number of function evaluations */
8677431f27SBarry Smith   PetscInt      fcount;
87ce94432eSBarry Smith   MPI_Comm      comm;
887736150eSBarry Smith   FILE         *fp;
89ace3abfcSBarry Smith   PetscBool     noise_test = PETSC_FALSE;
907736150eSBarry Smith 
917736150eSBarry Smith   PetscFunctionBegin;
929566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
937736150eSBarry Smith   /* Call to SNESSetUp() just to set data structures in SNES context */
949566063dSJacob Faibussowitsch   if (!snes->setupcalled) PetscCall(SNESSetUp(snes));
957736150eSBarry Smith 
967736150eSBarry Smith   w    = neP->workv[0];
977736150eSBarry Smith   xp   = neP->workv[1];
987736150eSBarry Smith   fvec = neP->workv[2];
997736150eSBarry Smith   fp   = neP->fp;
1007736150eSBarry Smith 
1017736150eSBarry Smith   /* Initialize parameters */
1027736150eSBarry Smith   hl       = zero;
1037736150eSBarry Smith   hu       = zero;
1047736150eSBarry Smith   h        = neP->h_first_try;
1057736150eSBarry Smith   fnoise_s = zero;
1067736150eSBarry Smith   fder2_s  = zero;
1077736150eSBarry Smith   fcount   = neP->function_count;
1087736150eSBarry Smith 
1097736150eSBarry Smith   /* We have 5 tries to attempt to compute a good hopt value */
1109566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &i));
11163a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "\n ------- SNES iteration %" PetscInt_FMT " ---------\n", i));
1127736150eSBarry Smith   for (iter = 0; iter < 5; iter++) {
1137736150eSBarry Smith     neP->h_first_try = h;
1147736150eSBarry Smith 
1157736150eSBarry Smith     /* Compute the nf function values needed to estimate the noise from
1167736150eSBarry Smith        the difference table */
1177736150eSBarry Smith     for (k = 0; k < nf; k++) {
1187736150eSBarry Smith       alpha = h * (k + 1 - (nf + 1) / 2);
1199566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(xp, alpha, p, x));
1209566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, xp, fvec));
1217736150eSBarry Smith       neP->function_count++;
1229566063dSJacob Faibussowitsch       PetscCall(VecDot(fvec, w, &fval[k]));
1237736150eSBarry Smith     }
1247736150eSBarry Smith     f = fval[(nf + 1) / 2 - 1];
1257736150eSBarry Smith 
1267736150eSBarry Smith     /* Construct the difference table */
127f5af7f23SKarl Rupp     for (i = 0; i < nf; i++) tab[i][0] = fval[i];
128f5af7f23SKarl Rupp 
129e108cb99SStefano Zampini     for (j = 0; j < nf - 1; j++) {
1309371c9d4SSatish Balay       for (i = 0; i < nf - j - 1; i++) { tab[i][j + 1] = tab[i + 1][j] - tab[i][j]; }
1317736150eSBarry Smith     }
1327736150eSBarry Smith 
1337736150eSBarry Smith     /* Print the difference table */
13463a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "Difference Table: iter = %" PetscInt_FMT "\n", iter));
1357736150eSBarry Smith     for (i = 0; i < nf; i++) {
136*48a46eb9SPierre Jolivet       for (j = 0; j < nf - i; j++) PetscCall(PetscFPrintf(comm, fp, " %10.2e ", tab[i][j]));
1379566063dSJacob Faibussowitsch       PetscCall(PetscFPrintf(comm, fp, "\n"));
1387736150eSBarry Smith     }
1397736150eSBarry Smith 
1407736150eSBarry Smith     /* Call the noise estimator */
1419566063dSJacob Faibussowitsch     PetscCall(SNESNoise_dnest_(&nf, fval, &h, fnoise, &fder2, hopt, &info, eps));
1427736150eSBarry Smith 
1437736150eSBarry Smith     /* Output statements */
1447736150eSBarry Smith     rerrf = *fnoise / PetscAbsScalar(f);
1459566063dSJacob Faibussowitsch     if (info == 1) PetscCall(PetscFPrintf(comm, fp, "%s\n", "Noise detected"));
146*48a46eb9SPierre Jolivet     if (info == 2) PetscCall(PetscFPrintf(comm, fp, "%s\n", "Noise not detected; h is too small"));
147*48a46eb9SPierre Jolivet     if (info == 3) PetscCall(PetscFPrintf(comm, fp, "%s\n", "Noise not detected; h is too large"));
1489566063dSJacob Faibussowitsch     if (info == 4) PetscCall(PetscFPrintf(comm, fp, "%s\n", "Noise detected, but unreliable hopt"));
1499566063dSJacob Faibussowitsch     PetscCall(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]));
1509566063dSJacob Faibussowitsch     PetscCall(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));
1517736150eSBarry Smith 
1527736150eSBarry Smith     /* Save fnoise and fder2. */
1537736150eSBarry Smith     if (*fnoise) fnoise_s = *fnoise;
1547736150eSBarry Smith     if (fder2) fder2_s = fder2;
1557736150eSBarry Smith 
1567736150eSBarry Smith     /* Check for noise detection. */
1577736150eSBarry Smith     if (fnoise_s && fder2_s) {
1587736150eSBarry Smith       *fnoise = fnoise_s;
1597736150eSBarry Smith       fder2   = fder2_s;
1607736150eSBarry Smith       *hopt   = 1.68 * sqrt(*fnoise / PetscAbsScalar(fder2));
1617736150eSBarry Smith       goto theend;
1627736150eSBarry Smith     } else {
1637736150eSBarry Smith       /* Update hl and hu, and determine new h */
1647736150eSBarry Smith       if (info == 2 || info == 4) {
1657736150eSBarry Smith         hl = h;
1667736150eSBarry Smith         if (hu == zero) h = 100 * h;
1677736150eSBarry Smith         else h = PetscMin(100 * h, 0.1 * hu);
1687736150eSBarry Smith       } else if (info == 3) {
1697736150eSBarry Smith         hu = h;
1707736150eSBarry Smith         h  = PetscMax(1.0e-3, sqrt(hl / hu)) * hu;
1717736150eSBarry Smith       }
1727736150eSBarry Smith     }
1737736150eSBarry Smith   }
1747736150eSBarry Smith theend:
1757736150eSBarry Smith 
1767736150eSBarry Smith   if (*fnoise < neP->fnoise_min) {
1779566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "Resetting fnoise: fnoise1 = %g, fnoise_min = %g\n", (double)*fnoise, (double)neP->fnoise_min));
1787736150eSBarry Smith     *fnoise = neP->fnoise_min;
1797736150eSBarry Smith     neP->fnoise_resets++;
1807736150eSBarry Smith   }
1817736150eSBarry Smith   if (*hopt < neP->hopt_min) {
1829566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "Resetting hopt: hopt1 = %g, hopt_min = %g\n", (double)*hopt, (double)neP->hopt_min));
1837736150eSBarry Smith     *hopt = neP->hopt_min;
1847736150eSBarry Smith     neP->hopt_resets++;
1857736150eSBarry Smith   }
1867736150eSBarry Smith 
1879566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "Errors in derivative:\n"));
1889566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "f = %g, fnoise = %g, fder2 = %g, hopt = %g\n", (double)f, (double)*fnoise, (double)fder2, (double)*hopt));
1897736150eSBarry Smith 
1907736150eSBarry Smith   /* For now, compute h **each** MV Mult!! */
1917736150eSBarry Smith   /*
1929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,"-matrix_free_jorge_each_mvp",&flg));
1937736150eSBarry Smith   if (!flg) {
1947736150eSBarry Smith     Mat mat;
1959566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes,&mat,NULL,NULL));
1969566063dSJacob Faibussowitsch     PetscCall(SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt));
1977736150eSBarry Smith   }
1987736150eSBarry Smith   */
1997736150eSBarry Smith   fcount = neP->function_count - fcount;
20063a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(snes, "fct_now = %" PetscInt_FMT ", fct_cum = %" PetscInt_FMT ", rerrf=%g, sqrt(noise)=%g, h_more=%g\n", fcount, neP->function_count, (double)rerrf, (double)PetscSqrtReal(*fnoise), (double)*hopt));
2017736150eSBarry Smith 
2029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-noise_test", &noise_test, NULL));
2031baa6e33SBarry Smith   if (noise_test) PetscCall(JacMatMultCompare(snes, x, p, *hopt));
2047736150eSBarry Smith   PetscFunctionReturn(0);
2057736150eSBarry Smith }
2067736150eSBarry Smith 
2079371c9d4SSatish Balay PetscErrorCode JacMatMultCompare(SNES snes, Vec x, Vec p, double hopt) {
2087736150eSBarry Smith   Vec         yy1, yy2; /* work vectors */
20952ed857cSBarry Smith   PetscViewer view2;    /* viewer */
2107736150eSBarry Smith   Mat         J;        /* analytic Jacobian (set as preconditioner matrix) */
2117736150eSBarry Smith   Mat         Jmf;      /* matrix-free Jacobian (set as true system matrix) */
2127736150eSBarry Smith   double      h;        /* differencing parameter */
2137736150eSBarry Smith   Vec         f;
21490c26489SBarry Smith   PetscScalar alpha;
21590c26489SBarry Smith   PetscReal   yy1n, yy2n, enorm;
216a7cc72afSBarry Smith   PetscInt    i;
217ace3abfcSBarry Smith   PetscBool   printv = PETSC_FALSE;
2187736150eSBarry Smith   char        filename[32];
219ce94432eSBarry Smith   MPI_Comm    comm;
2207736150eSBarry Smith 
2217736150eSBarry Smith   PetscFunctionBegin;
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
2237736150eSBarry Smith   /* Compute function and analytic Jacobian at x */
2249566063dSJacob Faibussowitsch   PetscCall(SNESGetJacobian(snes, &Jmf, &J, NULL, NULL));
2259566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, x, Jmf, J));
2269566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &f, NULL, NULL));
2279566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, x, f));
2287736150eSBarry Smith 
2297736150eSBarry Smith   /* Duplicate work vectors */
2309566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &yy2));
2319566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &yy1));
2327736150eSBarry Smith 
2337736150eSBarry Smith   /* Compute true matrix-vector product */
2349566063dSJacob Faibussowitsch   PetscCall(MatMult(J, p, yy1));
2359566063dSJacob Faibussowitsch   PetscCall(VecNorm(yy1, NORM_2, &yy1n));
2367736150eSBarry Smith 
2377736150eSBarry Smith   /* View product vector if desired */
2389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-print_vecs", &printv, NULL));
2397736150eSBarry Smith   if (printv) {
2409566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm, "y1.out", &view2));
2419566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(view2, PETSC_VIEWER_ASCII_COMMON));
2429566063dSJacob Faibussowitsch     PetscCall(VecView(yy1, view2));
2439566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(view2));
2449566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&view2));
2457736150eSBarry Smith   }
2467736150eSBarry Smith 
2477736150eSBarry Smith   /* Test Jacobian-vector product computation */
2487736150eSBarry Smith   alpha = -1.0;
2497736150eSBarry Smith   h     = 0.01 * hopt;
2507736150eSBarry Smith   for (i = 0; i < 5; i++) {
2517736150eSBarry Smith     /* Set differencing parameter for matrix-free multiplication */
2529566063dSJacob Faibussowitsch     PetscCall(SNESDefaultMatrixFreeSetParameters2(Jmf, PETSC_DEFAULT, PETSC_DEFAULT, h));
2537736150eSBarry Smith 
2547736150eSBarry Smith     /* Compute matrix-vector product via differencing approximation */
2559566063dSJacob Faibussowitsch     PetscCall(MatMult(Jmf, p, yy2));
2569566063dSJacob Faibussowitsch     PetscCall(VecNorm(yy2, NORM_2, &yy2n));
2577736150eSBarry Smith 
2587736150eSBarry Smith     /* View product vector if desired */
2597736150eSBarry Smith     if (printv) {
26077431f27SBarry Smith       sprintf(filename, "y2.%d.out", (int)i);
2619566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIOpen(comm, filename, &view2));
2629566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(view2, PETSC_VIEWER_ASCII_COMMON));
2639566063dSJacob Faibussowitsch       PetscCall(VecView(yy2, view2));
2649566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(view2));
2659566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&view2));
2667736150eSBarry Smith     }
2677736150eSBarry Smith 
2687736150eSBarry Smith     /* Compute relative error */
2699566063dSJacob Faibussowitsch     PetscCall(VecAXPY(yy2, alpha, yy1));
2709566063dSJacob Faibussowitsch     PetscCall(VecNorm(yy2, NORM_2, &enorm));
2717736150eSBarry Smith     enorm = enorm / yy1n;
2729566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, stdout, "h = %g: relative error = %g\n", (double)h, (double)enorm));
2737736150eSBarry Smith     h *= 10.0;
2747736150eSBarry Smith   }
2757736150eSBarry Smith   PetscFunctionReturn(0);
2767736150eSBarry Smith }
2777736150eSBarry Smith 
278a7cc72afSBarry Smith static PetscInt lin_its_total = 0;
2797736150eSBarry Smith 
2809371c9d4SSatish Balay PetscErrorCode SNESNoiseMonitor(SNES snes, PetscInt its, double fnorm, void *dummy) {
28177431f27SBarry Smith   PetscInt lin_its;
2827736150eSBarry Smith 
2837736150eSBarry Smith   PetscFunctionBegin;
2849566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(snes, &lin_its));
2857736150eSBarry Smith   lin_its_total += lin_its;
28663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), "iter = %" PetscInt_FMT ", SNES Function norm = %g, lin_its = %" PetscInt_FMT ", total_lin_its = %" PetscInt_FMT "\n", its, (double)fnorm, lin_its, lin_its_total));
2877736150eSBarry Smith 
2889566063dSJacob Faibussowitsch   PetscCall(SNESUnSetMatrixFreeParameter(snes));
2897736150eSBarry Smith   PetscFunctionReturn(0);
2907736150eSBarry Smith }
291