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 */ 12*5bd1e576SStefano 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 207736150eSBarry Smith 2125acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESDefaultMatrixFreeSetParameters2(Mat,double,double,double); 2225acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESUnSetMatrixFreeParameter(SNES snes); 2325acbd8eSLisandro Dalcin PETSC_INTERN PetscErrorCode SNESNoise_dnest_(PetscInt*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscScalar*,PetscInt*,PetscScalar*); 2425acbd8eSLisandro Dalcin 25ca8e9878SJed Brown static PetscErrorCode JacMatMultCompare(SNES,Vec,Vec,double); 267736150eSBarry Smith 278b5db460SBarry Smith PetscErrorCode SNESDiffParameterCreate_More(SNES snes,Vec x,void **outneP) 287736150eSBarry Smith { 297736150eSBarry Smith DIFFPAR_MORE *neP; 307736150eSBarry Smith Vec w; 317736150eSBarry Smith PetscRandom rctx; /* random number generator context */ 32dfbe8321SBarry Smith PetscErrorCode ierr; 33ace3abfcSBarry Smith PetscBool flg; 34e2d1d2b7SBarry Smith char noise_file[PETSC_MAX_PATH_LEN]; 357736150eSBarry Smith 367736150eSBarry Smith PetscFunctionBegin; 37b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 38f5af7f23SKarl Rupp 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 */ 51ce94432eSBarry Smith ierr = PetscRandomCreate(PetscObjectComm((PetscObject)snes),&rctx);CHKERRQ(ierr); 52c77d6671SHong Zhang ierr = PetscRandomSetFromOptions(rctx);CHKERRQ(ierr); 53abb0e124SMatthew Knepley ierr = VecSetRandom(w,rctx);CHKERRQ(ierr); 546bf464f9SBarry Smith ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr); 557736150eSBarry Smith 567736150eSBarry Smith /* Open output file */ 57c5929fdfSBarry Smith ierr = PetscOptionsGetString(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_mf_noise_file",noise_file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 587736150eSBarry Smith if (flg) neP->fp = fopen(noise_file,"w"); 597736150eSBarry Smith else neP->fp = fopen("noise.out","w"); 60e32f2f54SBarry Smith if (!neP->fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file"); 61ae15b995SBarry Smith ierr = PetscInfo(snes,"Creating Jorge's differencing parameter context\n");CHKERRQ(ierr); 627736150eSBarry Smith 637736150eSBarry Smith *outneP = neP; 647736150eSBarry Smith PetscFunctionReturn(0); 657736150eSBarry Smith } 667736150eSBarry Smith 678b5db460SBarry Smith PetscErrorCode SNESDiffParameterDestroy_More(void *nePv) 687736150eSBarry Smith { 697736150eSBarry Smith DIFFPAR_MORE *neP = (DIFFPAR_MORE*)nePv; 70dfbe8321SBarry Smith PetscErrorCode ierr; 71ed9cf6e9SBarry Smith int err; 727736150eSBarry Smith 737736150eSBarry Smith PetscFunctionBegin; 747736150eSBarry Smith /* Destroy work vectors and close output file */ 75574deadeSBarry Smith ierr = VecDestroyVecs(3,&neP->workv);CHKERRQ(ierr); 76ed9cf6e9SBarry Smith err = fclose(neP->fp); 77e32f2f54SBarry Smith if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fclose() failed on file"); 787736150eSBarry Smith ierr = PetscFree(neP);CHKERRQ(ierr); 797736150eSBarry Smith PetscFunctionReturn(0); 807736150eSBarry Smith } 817736150eSBarry Smith 828b5db460SBarry Smith PetscErrorCode SNESDiffParameterCompute_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; 883964eb88SJed Brown PetscScalar fval[7], tab[7][7], eps[7], f = -1; 893964eb88SJed Brown double rerrf = -1., 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; 94ce94432eSBarry Smith MPI_Comm comm; 957736150eSBarry Smith FILE *fp; 96ace3abfcSBarry Smith PetscBool noise_test = PETSC_FALSE; 977736150eSBarry Smith 987736150eSBarry Smith PetscFunctionBegin; 99ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 1007736150eSBarry Smith /* Call to SNESSetUp() just to set data structures in SNES context */ 10170e92668SMatthew Knepley if (!snes->setupcalled) {ierr = SNESSetUp(snes);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); 11877431f27SBarry 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); 1262dcb1b2aSMatthew Knepley ierr = VecWAXPY(xp,alpha,p,x);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 */ 134f5af7f23SKarl Rupp for (i=0; i<nf; i++) tab[i][0] = fval[i]; 135f5af7f23SKarl Rupp 136e108cb99SStefano Zampini for (j=0; j<nf-1; j++) { 137e108cb99SStefano Zampini for (i=0; i<nf-j-1; 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 */ 152ca8e9878SJed Brown ierr = SNESNoise_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);} 16057622a8eSBarry Smith ierr = 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]);CHKERRQ(ierr); 16157622a8eSBarry Smith ierr = 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);CHKERRQ(ierr); 1627736150eSBarry Smith 1637736150eSBarry Smith /* Save fnoise and fder2. */ 1647736150eSBarry Smith if (*fnoise) fnoise_s = *fnoise; 1657736150eSBarry Smith if (fder2) fder2_s = fder2; 1667736150eSBarry Smith 1677736150eSBarry Smith /* Check for noise detection. */ 1687736150eSBarry Smith if (fnoise_s && fder2_s) { 1697736150eSBarry Smith *fnoise = fnoise_s; 1707736150eSBarry Smith fder2 = fder2_s; 1717736150eSBarry Smith *hopt = 1.68*sqrt(*fnoise/PetscAbsScalar(fder2)); 1727736150eSBarry Smith goto theend; 1737736150eSBarry Smith } else { 1747736150eSBarry Smith 1757736150eSBarry Smith /* Update hl and hu, and determine new h */ 1767736150eSBarry Smith if (info == 2 || info == 4) { 1777736150eSBarry Smith hl = h; 1787736150eSBarry Smith if (hu == zero) h = 100*h; 1797736150eSBarry Smith else h = PetscMin(100*h,0.1*hu); 1807736150eSBarry Smith } else if (info == 3) { 1817736150eSBarry Smith hu = h; 1827736150eSBarry Smith h = PetscMax(1.0e-3,sqrt(hl/hu))*hu; 1837736150eSBarry Smith } 1847736150eSBarry Smith } 1857736150eSBarry Smith } 1867736150eSBarry Smith theend: 1877736150eSBarry Smith 1887736150eSBarry Smith if (*fnoise < neP->fnoise_min) { 18957622a8eSBarry Smith ierr = PetscFPrintf(comm,fp,"Resetting fnoise: fnoise1 = %g, fnoise_min = %g\n",(double)*fnoise,(double)neP->fnoise_min);CHKERRQ(ierr); 1907736150eSBarry Smith *fnoise = neP->fnoise_min; 1917736150eSBarry Smith neP->fnoise_resets++; 1927736150eSBarry Smith } 1937736150eSBarry Smith if (*hopt < neP->hopt_min) { 19457622a8eSBarry Smith ierr = PetscFPrintf(comm,fp,"Resetting hopt: hopt1 = %g, hopt_min = %g\n",(double)*hopt,(double)neP->hopt_min);CHKERRQ(ierr); 1957736150eSBarry Smith *hopt = neP->hopt_min; 1967736150eSBarry Smith neP->hopt_resets++; 1977736150eSBarry Smith } 1987736150eSBarry Smith 1997736150eSBarry Smith ierr = PetscFPrintf(comm,fp,"Errors in derivative:\n");CHKERRQ(ierr); 20057622a8eSBarry Smith ierr = PetscFPrintf(comm,fp,"f = %g, fnoise = %g, fder2 = %g, hopt = %g\n",(double)f,(double)*fnoise,(double)fder2,(double)*hopt);CHKERRQ(ierr); 2017736150eSBarry Smith 2027736150eSBarry Smith /* For now, compute h **each** MV Mult!! */ 2037736150eSBarry Smith /* 2040298fd71SBarry Smith ierr = PetscOptionsHasName(NULL,"-matrix_free_jorge_each_mvp",&flg);CHKERRQ(ierr); 2057736150eSBarry Smith if (!flg) { 2067736150eSBarry Smith Mat mat; 2070298fd71SBarry Smith ierr = SNESGetJacobian(snes,&mat,NULL,NULL);CHKERRQ(ierr); 2087736150eSBarry Smith ierr = SNESDefaultMatrixFreeSetParameters2(mat,PETSC_DEFAULT,PETSC_DEFAULT,*hopt);CHKERRQ(ierr); 2097736150eSBarry Smith } 2107736150eSBarry Smith */ 2117736150eSBarry Smith fcount = neP->function_count - fcount; 21257622a8eSBarry Smith ierr = PetscInfo5(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);CHKERRQ(ierr); 2137736150eSBarry Smith 214c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-noise_test",&noise_test,NULL);CHKERRQ(ierr); 2157736150eSBarry Smith if (noise_test) { 2167736150eSBarry Smith ierr = JacMatMultCompare(snes,x,p,*hopt);CHKERRQ(ierr); 2177736150eSBarry Smith } 2187736150eSBarry Smith PetscFunctionReturn(0); 2197736150eSBarry Smith } 2207736150eSBarry Smith 221dfbe8321SBarry Smith PetscErrorCode JacMatMultCompare(SNES snes,Vec x,Vec p,double hopt) 2227736150eSBarry Smith { 2237736150eSBarry Smith Vec yy1, yy2; /* work vectors */ 22452ed857cSBarry Smith PetscViewer view2; /* viewer */ 2257736150eSBarry Smith Mat J; /* analytic Jacobian (set as preconditioner matrix) */ 2267736150eSBarry Smith Mat Jmf; /* matrix-free Jacobian (set as true system matrix) */ 2277736150eSBarry Smith double h; /* differencing parameter */ 2287736150eSBarry Smith Vec f; 22990c26489SBarry Smith PetscScalar alpha; 23090c26489SBarry Smith PetscReal yy1n,yy2n,enorm; 2316849ba73SBarry Smith PetscErrorCode ierr; 232a7cc72afSBarry Smith PetscInt i; 233ace3abfcSBarry Smith PetscBool printv = PETSC_FALSE; 2347736150eSBarry Smith char filename[32]; 235ce94432eSBarry Smith MPI_Comm comm; 2367736150eSBarry Smith 2377736150eSBarry Smith PetscFunctionBegin; 238ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 2397736150eSBarry Smith /* Compute function and analytic Jacobian at x */ 2400298fd71SBarry Smith ierr = SNESGetJacobian(snes,&Jmf,&J,NULL,NULL);CHKERRQ(ierr); 241d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,x,Jmf,J);CHKERRQ(ierr); 2420298fd71SBarry Smith ierr = SNESGetFunction(snes,&f,NULL,NULL);CHKERRQ(ierr); 2437736150eSBarry Smith ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); 2447736150eSBarry Smith 2457736150eSBarry Smith /* Duplicate work vectors */ 2467736150eSBarry Smith ierr = VecDuplicate(x,&yy2);CHKERRQ(ierr); 2477736150eSBarry Smith ierr = VecDuplicate(x,&yy1);CHKERRQ(ierr); 2487736150eSBarry Smith 2497736150eSBarry Smith /* Compute true matrix-vector product */ 2507736150eSBarry Smith ierr = MatMult(J,p,yy1);CHKERRQ(ierr); 2517736150eSBarry Smith ierr = VecNorm(yy1,NORM_2,&yy1n);CHKERRQ(ierr); 2527736150eSBarry Smith 2537736150eSBarry Smith /* View product vector if desired */ 254c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-print_vecs",&printv,NULL);CHKERRQ(ierr); 2557736150eSBarry Smith if (printv) { 25652ed857cSBarry Smith ierr = PetscViewerASCIIOpen(comm,"y1.out",&view2);CHKERRQ(ierr); 2576a9046bcSBarry Smith ierr = PetscViewerPushFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr); 2587736150eSBarry Smith ierr = VecView(yy1,view2);CHKERRQ(ierr); 2596a9046bcSBarry Smith ierr = PetscViewerPopFormat(view2);CHKERRQ(ierr); 2606bf464f9SBarry Smith ierr = PetscViewerDestroy(&view2);CHKERRQ(ierr); 2617736150eSBarry Smith } 2627736150eSBarry Smith 2637736150eSBarry Smith /* Test Jacobian-vector product computation */ 2647736150eSBarry Smith alpha = -1.0; 2657736150eSBarry Smith h = 0.01 * hopt; 2667736150eSBarry Smith for (i=0; i<5; i++) { 2677736150eSBarry Smith /* Set differencing parameter for matrix-free multiplication */ 2687736150eSBarry Smith ierr = SNESDefaultMatrixFreeSetParameters2(Jmf,PETSC_DEFAULT,PETSC_DEFAULT,h);CHKERRQ(ierr); 2697736150eSBarry Smith 2707736150eSBarry Smith /* Compute matrix-vector product via differencing approximation */ 2717736150eSBarry Smith ierr = MatMult(Jmf,p,yy2);CHKERRQ(ierr); 2727736150eSBarry Smith ierr = VecNorm(yy2,NORM_2,&yy2n);CHKERRQ(ierr); 2737736150eSBarry Smith 2747736150eSBarry Smith /* View product vector if desired */ 2757736150eSBarry Smith if (printv) { 27677431f27SBarry Smith sprintf(filename,"y2.%d.out",(int)i); 27752ed857cSBarry Smith ierr = PetscViewerASCIIOpen(comm,filename,&view2);CHKERRQ(ierr); 2786a9046bcSBarry Smith ierr = PetscViewerPushFormat(view2,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr); 2797736150eSBarry Smith ierr = VecView(yy2,view2);CHKERRQ(ierr); 2806a9046bcSBarry Smith ierr = PetscViewerPopFormat(view2);CHKERRQ(ierr); 2816bf464f9SBarry Smith ierr = PetscViewerDestroy(&view2);CHKERRQ(ierr); 2827736150eSBarry Smith } 2837736150eSBarry Smith 2847736150eSBarry Smith /* Compute relative error */ 2852dcb1b2aSMatthew Knepley ierr = VecAXPY(yy2,alpha,yy1);CHKERRQ(ierr); 2867736150eSBarry Smith ierr = VecNorm(yy2,NORM_2,&enorm);CHKERRQ(ierr); 2877736150eSBarry Smith enorm = enorm/yy1n; 28857622a8eSBarry Smith ierr = PetscFPrintf(comm,stdout,"h = %g: relative error = %g\n",(double)h,(double)enorm);CHKERRQ(ierr); 2897736150eSBarry Smith h *= 10.0; 2907736150eSBarry Smith } 2917736150eSBarry Smith PetscFunctionReturn(0); 2927736150eSBarry Smith } 2937736150eSBarry Smith 294a7cc72afSBarry Smith static PetscInt lin_its_total = 0; 2957736150eSBarry Smith 296b12ffc70SBarry Smith PetscErrorCode SNESNoiseMonitor(SNES snes,PetscInt its,double fnorm,void *dummy) 2977736150eSBarry Smith { 2986849ba73SBarry Smith PetscErrorCode ierr; 29977431f27SBarry Smith PetscInt lin_its; 3007736150eSBarry Smith 3017736150eSBarry Smith PetscFunctionBegin; 302b850b91aSLisandro Dalcin ierr = SNESGetLinearSolveIterations(snes,&lin_its);CHKERRQ(ierr); 3037736150eSBarry Smith lin_its_total += lin_its; 30457622a8eSBarry Smith ierr = 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);CHKERRQ(ierr); 3057736150eSBarry Smith 3067736150eSBarry Smith ierr = SNESUnSetMatrixFreeParameter(snes);CHKERRQ(ierr); 3077736150eSBarry Smith PetscFunctionReturn(0); 3087736150eSBarry Smith } 309