xref: /petsc/src/snes/interface/noise/snesnoise.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1051ac159SSatish Balay 
27736150eSBarry Smith #include "src/snes/snesimpl.h"
37736150eSBarry Smith 
47736150eSBarry Smith /* Data used by Jorge's diff parameter computation method */
57736150eSBarry Smith typedef struct {
67736150eSBarry Smith   Vec    *workv;          /* work vectors */
77736150eSBarry Smith   FILE   *fp;             /* output file */
87736150eSBarry Smith   int    function_count;  /* count of function evaluations for diff param estimation */
97736150eSBarry Smith   double fnoise_min;      /* minimim allowable noise */
107736150eSBarry Smith   double hopt_min;        /* minimum allowable hopt */
117736150eSBarry Smith   double h_first_try;     /* first try for h used in diff parameter estimate */
127736150eSBarry Smith   int    fnoise_resets;   /* number of times we've reset the noise estimate */
137736150eSBarry Smith   int    hopt_resets;     /* number of times we've reset the hopt estimate */
147736150eSBarry Smith } DIFFPAR_MORE;
157736150eSBarry Smith 
167736150eSBarry Smith 
17f6428e95SBarry Smith extern int dnest_(int*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,int*,PetscScalar*);
18*6849ba73SBarry Smith EXTERN PetscErrorCode JacMatMultCompare(SNES,Vec,Vec,double);
19dfbe8321SBarry Smith EXTERN PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat,double,double,double);
20dfbe8321SBarry Smith EXTERN PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes);
217736150eSBarry Smith 
224a2ae208SSatish Balay #undef __FUNCT__
234a2ae208SSatish Balay #define __FUNCT__ "DiffParameterCreate_More"
24dfbe8321SBarry Smith PetscErrorCode DiffParameterCreate_More(SNES snes,Vec x,void **outneP)
257736150eSBarry Smith {
267736150eSBarry Smith   DIFFPAR_MORE *neP;
277736150eSBarry Smith   Vec          w;
287736150eSBarry Smith   PetscRandom  rctx;  /* random number generator context */
29dfbe8321SBarry Smith   PetscErrorCode ierr;
3052ed857cSBarry Smith   PetscTruth   flg;
31e2d1d2b7SBarry Smith   char         noise_file[PETSC_MAX_PATH_LEN];
327736150eSBarry Smith 
337736150eSBarry Smith   PetscFunctionBegin;
347736150eSBarry Smith 
3552ed857cSBarry Smith   ierr = PetscNew(DIFFPAR_MORE,&neP);CHKERRQ(ierr);
367736150eSBarry Smith   ierr = PetscMemzero(neP,sizeof(DIFFPAR_MORE));CHKERRQ(ierr);
3752ed857cSBarry Smith   PetscLogObjectMemory(snes,sizeof(DIFFPAR_MORE));
387736150eSBarry Smith 
397736150eSBarry Smith   neP->function_count = 0;
407736150eSBarry Smith   neP->fnoise_min     = 1.0e-20;
417736150eSBarry Smith   neP->hopt_min       = 1.0e-8;
427736150eSBarry Smith   neP->h_first_try    = 1.0e-3;
437736150eSBarry Smith   neP->fnoise_resets  = 0;
447736150eSBarry Smith   neP->hopt_resets    = 0;
457736150eSBarry Smith 
467736150eSBarry Smith   /* Create work vectors */
477736150eSBarry Smith   ierr = VecDuplicateVecs(x,3,&neP->workv);CHKERRQ(ierr);
487736150eSBarry Smith   w = neP->workv[0];
497736150eSBarry Smith 
507736150eSBarry Smith   /* Set components of vector w to random numbers */
517736150eSBarry Smith   ierr = PetscRandomCreate(snes->comm,RANDOM_DEFAULT,&rctx);CHKERRQ(ierr);
527736150eSBarry Smith   ierr = VecSetRandom(rctx,w);CHKERRQ(ierr);
537736150eSBarry Smith   ierr = PetscRandomDestroy(rctx);CHKERRQ(ierr);
547736150eSBarry Smith 
557736150eSBarry Smith   /* Open output file */
56e2d1d2b7SBarry Smith   ierr = PetscOptionsGetString(snes->prefix,"-snes_mf_noise_file",noise_file,PETSC_MAX_PATH_LEN-1,&flg);CHKERRQ(ierr);
577736150eSBarry Smith   if (flg) neP->fp = fopen(noise_file,"w");
587736150eSBarry Smith   else     neP->fp = fopen("noise.out","w");
5952ed857cSBarry Smith   if (!neP->fp) SETERRQ(PETSC_ERR_FILE_OPEN,"Cannot open file");
6052ed857cSBarry Smith   PetscLogInfo(snes,"DiffParameterCreate_More: Creating Jorge's differencing parameter context\n");
617736150eSBarry Smith 
627736150eSBarry Smith   *outneP = neP;
637736150eSBarry Smith   PetscFunctionReturn(0);
647736150eSBarry Smith }
657736150eSBarry Smith 
664a2ae208SSatish Balay #undef __FUNCT__
674a2ae208SSatish Balay #define __FUNCT__ "DiffParameterDestroy_More"
68dfbe8321SBarry Smith PetscErrorCode DiffParameterDestroy_More(void *nePv)
697736150eSBarry Smith {
707736150eSBarry Smith   DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv;
71dfbe8321SBarry Smith   PetscErrorCode ierr;
727736150eSBarry Smith 
737736150eSBarry Smith   PetscFunctionBegin;
747736150eSBarry Smith   /* Destroy work vectors and close output file */
757736150eSBarry Smith   ierr = VecDestroyVecs(neP->workv,3);CHKERRQ(ierr);
767736150eSBarry Smith   fclose(neP->fp);
777736150eSBarry Smith   ierr = PetscFree(neP);CHKERRQ(ierr);
787736150eSBarry Smith   PetscFunctionReturn(0);
797736150eSBarry Smith }
807736150eSBarry Smith 
814a2ae208SSatish Balay #undef __FUNCT__
824a2ae208SSatish Balay #define __FUNCT__ "DiffParameterCompute_More"
83dfbe8321SBarry Smith PetscErrorCode DiffParameterCompute_More(SNES snes,void *nePv,Vec x,Vec p,double *fnoise,double *hopt)
847736150eSBarry Smith {
857736150eSBarry Smith   DIFFPAR_MORE *neP = (DIFFPAR_MORE *)nePv;
867736150eSBarry Smith   Vec         w, xp, fvec;    /* work vectors to use in computing h */
8790c26489SBarry Smith   double      zero = 0.0, hl, hu, h, fnoise_s, fder2_s;
8890c26489SBarry Smith   PetscScalar alpha;
8990c26489SBarry Smith   PetscScalar fval[7], tab[7][7], eps[7], f;
9090c26489SBarry Smith   double      rerrf, fder2;
91*6849ba73SBarry Smith   PetscErrorCode ierr;
92*6849ba73SBarry Smith   int         iter, k, i, j,  info;
937736150eSBarry Smith   int         nf = 7;         /* number of function evaluations */
9452ed857cSBarry Smith   int         fcount;
957736150eSBarry Smith   MPI_Comm    comm = snes->comm;
967736150eSBarry Smith   FILE        *fp;
9752ed857cSBarry Smith   PetscTruth  noise_test;
987736150eSBarry Smith 
997736150eSBarry Smith   PetscFunctionBegin;
1007736150eSBarry Smith   /* Call to SNESSetUp() just to set data structures in SNES context */
1017736150eSBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
1027736150eSBarry Smith 
1037736150eSBarry Smith   w    = neP->workv[0];
1047736150eSBarry Smith   xp   = neP->workv[1];
1057736150eSBarry Smith   fvec = neP->workv[2];
1067736150eSBarry Smith   fp   = neP->fp;
1077736150eSBarry Smith 
1087736150eSBarry Smith   /* Initialize parameters */
1097736150eSBarry Smith   hl       = zero;
1107736150eSBarry Smith   hu       = zero;
1117736150eSBarry Smith   h        = neP->h_first_try;
1127736150eSBarry Smith   fnoise_s = zero;
1137736150eSBarry Smith   fder2_s  = zero;
1147736150eSBarry Smith   fcount   = neP->function_count;
1157736150eSBarry Smith 
1167736150eSBarry Smith   /* We have 5 tries to attempt to compute a good hopt value */
1177736150eSBarry Smith   ierr = SNESGetIterationNumber(snes,&i);CHKERRQ(ierr);
1187736150eSBarry Smith   ierr = PetscFPrintf(comm,fp,"\n ------- SNES iteration %d ---------\n",i);CHKERRQ(ierr);
1197736150eSBarry Smith   for (iter=0; iter<5; iter++) {
1207736150eSBarry Smith     neP->h_first_try = h;
1217736150eSBarry Smith 
1227736150eSBarry Smith     /* Compute the nf function values needed to estimate the noise from
1237736150eSBarry Smith        the difference table */
1247736150eSBarry Smith     for (k=0; k<nf; k++) {
1257736150eSBarry Smith       alpha = h * ( k+1 - (nf+1)/2 );
1267736150eSBarry Smith       ierr = VecWAXPY(&alpha,p,x,xp);CHKERRQ(ierr);
1277736150eSBarry Smith       ierr = SNESComputeFunction(snes,xp,fvec);CHKERRQ(ierr);
1287736150eSBarry Smith       neP->function_count++;
1297736150eSBarry Smith       ierr = VecDot(fvec,w,&fval[k]);CHKERRQ(ierr);
1307736150eSBarry Smith     }
1317736150eSBarry Smith     f = fval[(nf+1)/2 - 1];
1327736150eSBarry Smith 
1337736150eSBarry Smith     /* Construct the difference table */
1347736150eSBarry Smith     for (i=0; i<nf; i++) {
1357736150eSBarry Smith       tab[i][0] = fval[i];
1367736150eSBarry Smith     }
1377736150eSBarry Smith     for (j=0; j<6; j++) {
1387736150eSBarry Smith       for (i=0; i<nf-j; i++) {
1397736150eSBarry Smith         tab[i][j+1] = tab[i+1][j] - tab[i][j];
1407736150eSBarry Smith       }
1417736150eSBarry Smith     }
1427736150eSBarry Smith 
1437736150eSBarry Smith     /* Print the difference table */
1447736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"Difference Table: iter = %d\n",iter);CHKERRQ(ierr);
1457736150eSBarry Smith     for (i=0; i<nf; i++) {
1467736150eSBarry Smith       for (j=0; j<nf-i; j++) {
1477736150eSBarry Smith         ierr = PetscFPrintf(comm,fp," %10.2e ",tab[i][j]);CHKERRQ(ierr);
1487736150eSBarry Smith       }
1497736150eSBarry Smith       ierr = PetscFPrintf(comm,fp,"\n");CHKERRQ(ierr);
1507736150eSBarry Smith     }
1517736150eSBarry Smith 
1527736150eSBarry Smith     /* Call the noise estimator */
153f6428e95SBarry Smith     ierr = dnest_(&nf,fval,&h,fnoise,&fder2,hopt,&info,eps);CHKERRQ(ierr);
1547736150eSBarry Smith 
1557736150eSBarry Smith     /* Output statements */
1567736150eSBarry Smith     rerrf = *fnoise/PetscAbsScalar(f);
1577736150eSBarry Smith     if (info == 1) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected");CHKERRQ(ierr);}
1587736150eSBarry Smith     if (info == 2) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too small");CHKERRQ(ierr);}
1597736150eSBarry Smith     if (info == 3) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise not detected; h is too large");CHKERRQ(ierr);}
1607736150eSBarry Smith     if (info == 4) {ierr = PetscFPrintf(comm,fp,"%s\n","Noise detected, but unreliable hopt");CHKERRQ(ierr);}
1617736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"Approximate epsfcn %g  %g  %g  %g  %g  %g\n",
1627736150eSBarry Smith         eps[0],eps[1],eps[2],eps[3],eps[4],eps[5]);CHKERRQ(ierr);
1637736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"h = %g, fnoise = %g, fder2 = %g, rerrf = %g, hopt = %g\n\n",
1647736150eSBarry Smith             h, *fnoise, fder2, rerrf, *hopt);CHKERRQ(ierr);
1657736150eSBarry Smith 
1667736150eSBarry Smith     /* Save fnoise and fder2. */
1677736150eSBarry Smith     if (*fnoise) fnoise_s = *fnoise;
1687736150eSBarry Smith     if (fder2)  fder2_s = fder2;
1697736150eSBarry Smith 
1707736150eSBarry Smith     /* Check for noise detection. */
1717736150eSBarry Smith     if (fnoise_s && fder2_s) {
1727736150eSBarry Smith       *fnoise = fnoise_s;
1737736150eSBarry Smith       fder2 = fder2_s;
1747736150eSBarry Smith       *hopt = 1.68*sqrt(*fnoise/PetscAbsScalar(fder2));
1757736150eSBarry Smith       goto theend;
1767736150eSBarry Smith     } else {
1777736150eSBarry Smith 
1787736150eSBarry Smith       /* Update hl and hu, and determine new h */
1797736150eSBarry Smith       if (info == 2 || info == 4) {
1807736150eSBarry Smith         hl = h;
1817736150eSBarry Smith         if (hu == zero) h = 100*h;
1827736150eSBarry Smith         else            h = PetscMin(100*h,0.1*hu);
1837736150eSBarry Smith        } else if (info == 3) {
1847736150eSBarry Smith         hu = h;
1857736150eSBarry Smith         h = PetscMax(1.0e-3,sqrt(hl/hu))*hu;
1867736150eSBarry Smith       }
1877736150eSBarry Smith     }
1887736150eSBarry Smith   }
1897736150eSBarry Smith   theend:
1907736150eSBarry Smith 
1917736150eSBarry Smith   if (*fnoise < neP->fnoise_min) {
1927736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %g, fnoise_min = %g\n",*fnoise,neP->fnoise_min);CHKERRQ(ierr);
1937736150eSBarry Smith     *fnoise = neP->fnoise_min;
1947736150eSBarry Smith     neP->fnoise_resets++;
1957736150eSBarry Smith   }
1967736150eSBarry Smith   if (*hopt < neP->hopt_min) {
1977736150eSBarry Smith     ierr = PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %g, hopt_min = %g\n",*hopt,neP->hopt_min);CHKERRQ(ierr);
1987736150eSBarry Smith     *hopt = neP->hopt_min;
1997736150eSBarry Smith     neP->hopt_resets++;
2007736150eSBarry Smith   }
2017736150eSBarry Smith 
2027736150eSBarry Smith   ierr = PetscFPrintf(comm,fp,"Errors in derivative:\n");CHKERRQ(ierr);
2037736150eSBarry Smith   ierr = PetscFPrintf(comm,fp,"f = %g, fnoise = %g, fder2 = %g, hopt = %g\n",f,*fnoise,fder2,*hopt);CHKERRQ(ierr);
2047736150eSBarry Smith 
2057736150eSBarry Smith   /* For now, compute h **each** MV Mult!! */
2067736150eSBarry Smith   /*
20752ed857cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-matrix_free_jorge_each_mvp",&flg);CHKERRQ(ierr);
2087736150eSBarry Smith   if (!flg) {
2097736150eSBarry Smith     Mat mat;
2107736150eSBarry Smith     ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2117736150eSBarry Smith     ierr = SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);CHKERRQ(ierr);
2127736150eSBarry Smith   }
2137736150eSBarry Smith   */
2147736150eSBarry Smith   fcount = neP->function_count - fcount;
21552ed857cSBarry Smith   PetscLogInfo(snes,"DiffParameterCompute_More: fct_now = %d, fct_cum = %d, rerrf=%g, sqrt(noise)=%g, h_more=%g\n",
2167736150eSBarry Smith            fcount,neP->function_count,rerrf,sqrt(*fnoise),*hopt);
2177736150eSBarry Smith 
2187736150eSBarry Smith 
21952ed857cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-noise_test",&noise_test);CHKERRQ(ierr);
2207736150eSBarry Smith   if (noise_test) {
2217736150eSBarry Smith     ierr = JacMatMultCompare(snes,x,p,*hopt);CHKERRQ(ierr);
2227736150eSBarry Smith   }
2237736150eSBarry Smith   PetscFunctionReturn(0);
2247736150eSBarry Smith }
2257736150eSBarry Smith 
2264a2ae208SSatish Balay #undef __FUNCT__
2274a2ae208SSatish Balay #define __FUNCT__ "JacMatMultCompare"
228dfbe8321SBarry Smith PetscErrorCode JacMatMultCompare(SNES snes,Vec x,Vec p,double hopt)
2297736150eSBarry Smith {
2307736150eSBarry Smith   Vec          yy1, yy2; /* work vectors */
23152ed857cSBarry Smith   PetscViewer  view2;    /* viewer */
2327736150eSBarry Smith   Mat          J;        /* analytic Jacobian (set as preconditioner matrix) */
2337736150eSBarry Smith   Mat          Jmf;      /* matrix-free Jacobian (set as true system matrix) */
2347736150eSBarry Smith   double       h;        /* differencing parameter */
2357736150eSBarry Smith   Vec          f;
2367736150eSBarry Smith   MatStructure sparsity = DIFFERENT_NONZERO_PATTERN;
23790c26489SBarry Smith   PetscScalar  alpha;
23890c26489SBarry Smith   PetscReal    yy1n,yy2n,enorm;
239*6849ba73SBarry Smith   PetscErrorCode ierr;
240*6849ba73SBarry Smith   int          i;
24152ed857cSBarry Smith   PetscTruth   printv;
2427736150eSBarry Smith   char         filename[32];
2437736150eSBarry Smith   MPI_Comm     comm = snes->comm;
2447736150eSBarry Smith 
2457736150eSBarry Smith   PetscFunctionBegin;
2467736150eSBarry Smith 
2477736150eSBarry Smith   /* Compute function and analytic Jacobian at x */
24852ed857cSBarry Smith   ierr = SNESGetJacobian(snes,&Jmf,&J,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2497736150eSBarry Smith   ierr = SNESComputeJacobian(snes,x,&Jmf,&J,&sparsity);CHKERRQ(ierr);
25052ed857cSBarry Smith   ierr = SNESGetFunction(snes,&f,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2517736150eSBarry Smith   ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
2527736150eSBarry Smith 
2537736150eSBarry Smith   /* Duplicate work vectors */
2547736150eSBarry Smith   ierr = VecDuplicate(x,&yy2);CHKERRQ(ierr);
2557736150eSBarry Smith   ierr = VecDuplicate(x,&yy1);CHKERRQ(ierr);
2567736150eSBarry Smith 
2577736150eSBarry Smith   /* Compute true matrix-vector product */
2587736150eSBarry Smith   ierr = MatMult(J,p,yy1);CHKERRQ(ierr);
2597736150eSBarry Smith   ierr = VecNorm(yy1,NORM_2,&yy1n);CHKERRQ(ierr);
2607736150eSBarry Smith 
2617736150eSBarry Smith   /* View product vector if desired */
26252ed857cSBarry Smith   ierr = PetscOptionsHasName(PETSC_NULL,"-print_vecs",&printv);CHKERRQ(ierr);
2637736150eSBarry Smith   if (printv) {
26452ed857cSBarry Smith     ierr = PetscViewerASCIIOpen(comm,"y1.out",&view2);CHKERRQ(ierr);
26552ed857cSBarry Smith     ierr = PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
2667736150eSBarry Smith     ierr = VecView(yy1,view2);CHKERRQ(ierr);
26752ed857cSBarry Smith     ierr = PetscViewerDestroy(view2);CHKERRQ(ierr);
2687736150eSBarry Smith   }
2697736150eSBarry Smith 
2707736150eSBarry Smith   /* Test Jacobian-vector product computation */
2717736150eSBarry Smith   alpha = -1.0;
2727736150eSBarry Smith   h = 0.01 * hopt;
2737736150eSBarry Smith   for (i=0; i<5; i++) {
2747736150eSBarry Smith     /* Set differencing parameter for matrix-free multiplication */
2757736150eSBarry Smith     ierr = SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);CHKERRQ(ierr);
2767736150eSBarry Smith 
2777736150eSBarry Smith     /* Compute matrix-vector product via differencing approximation */
2787736150eSBarry Smith     ierr = MatMult(Jmf,p,yy2);CHKERRQ(ierr);
2797736150eSBarry Smith     ierr = VecNorm(yy2,NORM_2,&yy2n);CHKERRQ(ierr);
2807736150eSBarry Smith 
2817736150eSBarry Smith     /* View product vector if desired */
2827736150eSBarry Smith     if (printv) {
2837736150eSBarry Smith       sprintf(filename,"y2.%d.out",i);
28452ed857cSBarry Smith       ierr = PetscViewerASCIIOpen(comm,filename,&view2);CHKERRQ(ierr);
28552ed857cSBarry Smith       ierr = PetscViewerSetFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
2867736150eSBarry Smith       ierr = VecView(yy2,view2);CHKERRQ(ierr);
28752ed857cSBarry Smith       ierr = PetscViewerDestroy(view2);CHKERRQ(ierr);
2887736150eSBarry Smith     }
2897736150eSBarry Smith 
2907736150eSBarry Smith     /* Compute relative error */
2917736150eSBarry Smith     ierr = VecAXPY(&alpha,yy1,yy2);CHKERRQ(ierr);
2927736150eSBarry Smith     ierr = VecNorm(yy2,NORM_2,&enorm);CHKERRQ(ierr);
2937736150eSBarry Smith     enorm = enorm/yy1n;
2947736150eSBarry Smith     ierr = PetscFPrintf(comm,stdout,"h = %g: relative error = %g\n",h,enorm);CHKERRQ(ierr);
2957736150eSBarry Smith     h *= 10.0;
2967736150eSBarry Smith   }
2977736150eSBarry Smith   PetscFunctionReturn(0);
2987736150eSBarry Smith }
2997736150eSBarry Smith 
3007736150eSBarry Smith static int lin_its_total = 0;
3017736150eSBarry Smith 
302dfbe8321SBarry Smith PetscErrorCode MyMonitor(SNES snes,int its,double fnorm,void *dummy)
3037736150eSBarry Smith {
304*6849ba73SBarry Smith   PetscErrorCode ierr;
305*6849ba73SBarry Smith   int  lin_its;
3067736150eSBarry Smith 
3077736150eSBarry Smith   PetscFunctionBegin;
3087736150eSBarry Smith   ierr = SNESGetNumberLinearIterations(snes,&lin_its);CHKERRQ(ierr);
3097736150eSBarry Smith   lin_its_total += lin_its;
3107736150eSBarry 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);
3117736150eSBarry Smith 
3127736150eSBarry Smith   ierr = SNESUnSetMatrixFreeParameter(snes);CHKERRQ(ierr);
3137736150eSBarry Smith   PetscFunctionReturn(0);
3147736150eSBarry Smith }
315