14b11644fSPeter Brune 24b11644fSPeter Brune #include <private/snesimpl.h> 34b11644fSPeter Brune 44b11644fSPeter Brune typedef struct { 540012acdSPeter Brune PetscReal lambda; /* The default step length for the update */ 64b11644fSPeter Brune Vec * dX; /* The change in X */ 74b11644fSPeter Brune Vec * dF; /* The change in F */ 84b11644fSPeter Brune PetscInt m; /* the number of kept previous steps */ 9*bd052dfeSPeter Brune PetscScalar * alpha; 10*bd052dfeSPeter Brune PetscScalar * beta; 11*bd052dfeSPeter Brune PetscScalar * rho; 124b11644fSPeter Brune } QNContext; 134b11644fSPeter Brune 144b11644fSPeter Brune #undef __FUNCT__ 154b11644fSPeter Brune #define __FUNCT__ "LBGFSApplyJinv_Private" 164b11644fSPeter Brune PetscErrorCode LBGFSApplyJinv_Private(SNES snes, PetscInt it, Vec g, Vec z) { 174b11644fSPeter Brune 184b11644fSPeter Brune PetscErrorCode ierr; 194b11644fSPeter Brune 204b11644fSPeter Brune QNContext * qn = (QNContext *)snes->data; 214b11644fSPeter Brune 224b11644fSPeter Brune Vec * dX = qn->dX; 234b11644fSPeter Brune Vec * dF = qn->dF; 244b11644fSPeter Brune 25*bd052dfeSPeter Brune PetscScalar * alpha = qn->alpha; 26*bd052dfeSPeter Brune PetscScalar * beta = qn->beta; 27*bd052dfeSPeter Brune PetscScalar * rho = qn->rho; 28*bd052dfeSPeter Brune 294b11644fSPeter Brune PetscInt k, i; 304b11644fSPeter Brune PetscInt m = qn->m; 314b11644fSPeter Brune PetscScalar t; 324b11644fSPeter Brune PetscInt l = m; 334b11644fSPeter Brune 344b11644fSPeter Brune PetscFunctionBegin; 354b11644fSPeter Brune 36*bd052dfeSPeter Brune if (it < m) l = it; 374b11644fSPeter Brune 384b11644fSPeter Brune ierr = VecCopy(g, z);CHKERRQ(ierr); 394b11644fSPeter Brune 404b11644fSPeter Brune /* outward recursion starting at iteration k's update and working back */ 414b11644fSPeter Brune for (i = 0; i < l; i++) { 42*bd052dfeSPeter Brune /* k = (it - i - 1) % m; */ 43*bd052dfeSPeter Brune k = (it + i - l) % m; 444b11644fSPeter Brune ierr = VecDot(dX[k], z, &t);CHKERRQ(ierr); 45*bd052dfeSPeter Brune alpha[k] = t*rho[k]; 464b11644fSPeter Brune ierr = VecAXPY(z, -alpha[k], dF[k]);CHKERRQ(ierr); 474b11644fSPeter Brune } 484b11644fSPeter Brune 494b11644fSPeter Brune /* inner application of the initial inverse jacobian approximation */ 504b11644fSPeter Brune /* right now it's just the identity. Nothing needs to go here. */ 514b11644fSPeter Brune 524b11644fSPeter Brune /* inward recursion starting at the first update and working forward*/ 53*bd052dfeSPeter Brune for (i = 0; i < l; i++) { 54*bd052dfeSPeter Brune k = (it - i - 1) % m; 55*bd052dfeSPeter Brune /* k = (it + i - l) % m; */ 564b11644fSPeter Brune ierr = VecDot(dF[k], z, &t);CHKERRQ(ierr); 574b11644fSPeter Brune beta[k] = rho[k]*t; 584b11644fSPeter Brune ierr = VecAXPY(z, (alpha[k] - beta[k]), dX[k]); 594b11644fSPeter Brune } 60*bd052dfeSPeter Brune ierr = VecScale(z, 1.0);CHKERRQ(ierr); 614b11644fSPeter Brune 624b11644fSPeter Brune PetscFunctionReturn(0); 634b11644fSPeter Brune } 644b11644fSPeter Brune 65*bd052dfeSPeter Brune 66*bd052dfeSPeter Brune #undef __FUNCT__ 67*bd052dfeSPeter Brune #define __FUNCT__ "QNLineSearchQuadratic" 68*bd052dfeSPeter Brune PetscErrorCode QNLineSearchQuadratic(SNES snes,void *lsctx,Vec X,Vec F,Vec G,Vec Y,Vec W,PetscReal fnorm,PetscReal dummyXnorm,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag) 69*bd052dfeSPeter Brune { 70*bd052dfeSPeter Brune PetscInt i; 71*bd052dfeSPeter Brune PetscReal alphas[3] = {0.0, 0.5, 1.0}; 72*bd052dfeSPeter Brune PetscReal norms[3]; 73*bd052dfeSPeter Brune PetscReal alpha,a,b; 74*bd052dfeSPeter Brune PetscErrorCode ierr; 75*bd052dfeSPeter Brune PetscFunctionBegin; 76*bd052dfeSPeter Brune norms[0] = fnorm; 77*bd052dfeSPeter Brune /* Calculate trial solutions */ 78*bd052dfeSPeter Brune for(i=1; i < 3; ++i) { 79*bd052dfeSPeter Brune /* Calculate X^{n+1} = (1 - \alpha) X^n + \alpha Y */ 80*bd052dfeSPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 81*bd052dfeSPeter Brune ierr = VecAXPBY(W, alphas[i], 1 - alphas[i], Y);CHKERRQ(ierr); 82*bd052dfeSPeter Brune ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr); 83*bd052dfeSPeter Brune ierr = VecNorm(F, NORM_2, &norms[i]);CHKERRQ(ierr); 84*bd052dfeSPeter Brune } 85*bd052dfeSPeter Brune for(i = 0; i < 3; ++i) { 86*bd052dfeSPeter Brune norms[i] = PetscSqr(norms[i]); 87*bd052dfeSPeter Brune } 88*bd052dfeSPeter Brune /* Fit a quadratic: 89*bd052dfeSPeter Brune If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2} 90*bd052dfeSPeter Brune a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1) 91*bd052dfeSPeter Brune b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1) 92*bd052dfeSPeter Brune c = y_0 93*bd052dfeSPeter Brune x_min = -b/2a 94*bd052dfeSPeter Brune 95*bd052dfeSPeter Brune If we let x_{0,1,2} = 0, 0.5, 1.0 96*bd052dfeSPeter Brune a = 2 y_2 - 4 y_1 + 2 y_0 97*bd052dfeSPeter Brune b = -y_2 + 4 y_1 - 3 y_0 98*bd052dfeSPeter Brune c = y_0 99*bd052dfeSPeter Brune */ 100*bd052dfeSPeter Brune a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1])); 101*bd052dfeSPeter Brune b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]); 102*bd052dfeSPeter Brune /* Check for positive a (concave up) */ 103*bd052dfeSPeter Brune if (a >= 0.0) { 104*bd052dfeSPeter Brune alpha = -b/(2.0*a); 105*bd052dfeSPeter Brune alpha = PetscMin(alpha, alphas[2]); 106*bd052dfeSPeter Brune alpha = PetscMax(alpha, alphas[0]); 107*bd052dfeSPeter Brune } else { 108*bd052dfeSPeter Brune alpha = 1.0; 109*bd052dfeSPeter Brune } 110*bd052dfeSPeter Brune ierr = VecAXPBY(X, alpha, 1 - alpha, Y);CHKERRQ(ierr); 111*bd052dfeSPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 112*bd052dfeSPeter Brune if (alpha != 1.0) { 113*bd052dfeSPeter Brune ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr); 114*bd052dfeSPeter Brune } else { 115*bd052dfeSPeter Brune *gnorm = PetscSqrtReal(norms[2]); 116*bd052dfeSPeter Brune } 117*bd052dfeSPeter Brune *flag = PETSC_TRUE; 118*bd052dfeSPeter Brune PetscFunctionReturn(0); 119*bd052dfeSPeter Brune } 120*bd052dfeSPeter Brune 121*bd052dfeSPeter Brune 1224b11644fSPeter Brune #undef __FUNCT__ 1234b11644fSPeter Brune #define __FUNCT__ "SNESSolve_QN" 1244b11644fSPeter Brune static PetscErrorCode SNESSolve_QN(SNES snes) 1254b11644fSPeter Brune { 1264b11644fSPeter Brune 1274b11644fSPeter Brune PetscErrorCode ierr; 1284b11644fSPeter Brune QNContext * qn = (QNContext*) snes->data; 1294b11644fSPeter Brune 130*bd052dfeSPeter Brune Vec x, xold; 131*bd052dfeSPeter Brune Vec f, fold; 132*bd052dfeSPeter Brune Vec w, y; 1334b11644fSPeter Brune 134*bd052dfeSPeter Brune PetscInt i, k; 1354b11644fSPeter Brune 136*bd052dfeSPeter Brune PetscReal fnorm, xnorm; 1374b11644fSPeter Brune PetscInt m = qn->m; 138*bd052dfeSPeter Brune PetscBool ls_OK; 1394b11644fSPeter Brune 140*bd052dfeSPeter Brune PetscScalar rhosc; 141*bd052dfeSPeter Brune 142*bd052dfeSPeter Brune Vec * dX = qn->dX; 143*bd052dfeSPeter Brune Vec * dF = qn->dF; 144*bd052dfeSPeter Brune PetscScalar * rho = qn->rho; 1454b11644fSPeter Brune 1464b11644fSPeter Brune /* basically just a regular newton's method except for the application of the jacobian */ 1474b11644fSPeter Brune PetscFunctionBegin; 1484b11644fSPeter Brune 1494b11644fSPeter Brune x = snes->vec_sol; 150*bd052dfeSPeter Brune xold = snes->vec_sol_update; /* dX_k */ 151*bd052dfeSPeter Brune w = snes->work[1]; 1524b11644fSPeter Brune f = snes->vec_func; 153*bd052dfeSPeter Brune fold = snes->work[0]; 154*bd052dfeSPeter Brune y = snes->work[2]; 1554b11644fSPeter Brune 1564b11644fSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1574b11644fSPeter Brune 1584b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1594b11644fSPeter Brune snes->iter = 0; 1604b11644fSPeter Brune snes->norm = 0.; 1614b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1624b11644fSPeter Brune ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); 1634b11644fSPeter Brune if (snes->domainerror) { 1644b11644fSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1654b11644fSPeter Brune PetscFunctionReturn(0); 1664b11644fSPeter Brune } 1674b11644fSPeter Brune ierr = VecNorm(f, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1684b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1694b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1704b11644fSPeter Brune snes->norm = fnorm; 1714b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1724b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,0); 1734b11644fSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1744b11644fSPeter Brune 1754b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 1764b11644fSPeter Brune snes->ttol = fnorm*snes->rtol; 1774b11644fSPeter Brune /* test convergence */ 1784b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1794b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 180*bd052dfeSPeter Brune ierr = VecCopy(f, fold);CHKERRQ(ierr); 181*bd052dfeSPeter Brune ierr = VecCopy(x, xold);CHKERRQ(ierr); 1824b11644fSPeter Brune for(i = 0; i < snes->max_its; i++) { 183*bd052dfeSPeter Brune /* general purpose update */ 184*bd052dfeSPeter Brune if (snes->ops->update) { 185*bd052dfeSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 186*bd052dfeSPeter Brune } 1874b11644fSPeter Brune 188*bd052dfeSPeter Brune /* apply the current iteration of the approximate jacobian */ 189*bd052dfeSPeter Brune ierr = LBGFSApplyJinv_Private(snes, i, f, y);CHKERRQ(ierr); 190*bd052dfeSPeter Brune 191*bd052dfeSPeter Brune /* line search for lambda */ 192*bd052dfeSPeter Brune ierr = VecAYPX(y,-1.0,x);CHKERRQ(ierr); 193*bd052dfeSPeter Brune ierr = QNLineSearchQuadratic(snes, PETSC_NULL, x, f, 0, y, w, fnorm, xnorm, &xnorm, &fnorm, &ls_OK);CHKERRQ(ierr); 1944b11644fSPeter Brune ierr = SNESComputeFunction(snes, x, f);CHKERRQ(ierr); 1954b11644fSPeter Brune ierr = VecNorm(f, NORM_2, &fnorm);CHKERRQ(ierr); 1964b11644fSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1974b11644fSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1984b11644fSPeter Brune snes->norm = fnorm; 1994b11644fSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2004b11644fSPeter Brune SNESLogConvHistory(snes,fnorm,i+1); 2014b11644fSPeter Brune ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 2024b11644fSPeter Brune /* set parameter for default relative tolerance convergence test */ 2034b11644fSPeter Brune ierr = (*snes->ops->converged)(snes,i+1,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2044b11644fSPeter Brune if (snes->reason) PetscFunctionReturn(0); 2054b11644fSPeter Brune 206*bd052dfeSPeter Brune /* set the differences */ 207*bd052dfeSPeter Brune k = i % m; 208*bd052dfeSPeter Brune ierr = VecCopy(f, dF[k]);CHKERRQ(ierr); 209*bd052dfeSPeter Brune ierr = VecAXPY(dF[k], -1.0, fold);CHKERRQ(ierr); 210*bd052dfeSPeter Brune ierr = VecCopy(x, dX[k]);CHKERRQ(ierr); 211*bd052dfeSPeter Brune ierr = VecAXPY(dX[k], -1.0, xold);CHKERRQ(ierr); 212*bd052dfeSPeter Brune ierr = VecDot(dX[k], dF[k], &rhosc);CHKERRQ(ierr); 213*bd052dfeSPeter Brune rho[k] = 1. / rhosc; 214*bd052dfeSPeter Brune ierr = VecCopy(f, fold);CHKERRQ(ierr); 215*bd052dfeSPeter Brune ierr = VecCopy(x, xold);CHKERRQ(ierr); 2164b11644fSPeter Brune } 2174b11644fSPeter Brune if (i == snes->max_its) { 2184b11644fSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 2194b11644fSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 2204b11644fSPeter Brune } 2214b11644fSPeter Brune PetscFunctionReturn(0); 2224b11644fSPeter Brune } 2234b11644fSPeter Brune 2244b11644fSPeter Brune 2254b11644fSPeter Brune #undef __FUNCT__ 2264b11644fSPeter Brune #define __FUNCT__ "SNESSetUp_QN" 2274b11644fSPeter Brune static PetscErrorCode SNESSetUp_QN(SNES snes) 2284b11644fSPeter Brune { 2294b11644fSPeter Brune QNContext * qn = (QNContext *)snes->data; 2304b11644fSPeter Brune PetscErrorCode ierr; 2314b11644fSPeter Brune PetscFunctionBegin; 2324b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dX);CHKERRQ(ierr); 2334b11644fSPeter Brune ierr = VecDuplicateVecs(snes->vec_sol, qn->m, &qn->dF);CHKERRQ(ierr); 234*bd052dfeSPeter Brune ierr = PetscMalloc3(qn->m, PetscScalar, &qn->alpha, qn->m, PetscScalar, &qn->beta, qn->m, PetscScalar, &qn->rho);CHKERRQ(ierr); 235*bd052dfeSPeter Brune ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr); 2364b11644fSPeter Brune PetscFunctionReturn(0); 2374b11644fSPeter Brune } 2384b11644fSPeter Brune 2394b11644fSPeter Brune #undef __FUNCT__ 2404b11644fSPeter Brune #define __FUNCT__ "SNESReset_QN" 2414b11644fSPeter Brune static PetscErrorCode SNESReset_QN(SNES snes) 2424b11644fSPeter Brune { 2434b11644fSPeter Brune PetscErrorCode ierr; 2444b11644fSPeter Brune QNContext * qn; 2454b11644fSPeter Brune PetscFunctionBegin; 2464b11644fSPeter Brune if (snes->data) { 2474b11644fSPeter Brune qn = (QNContext *)snes->data; 2484b11644fSPeter Brune if (qn->dX) { 2494b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dX);CHKERRQ(ierr); 2504b11644fSPeter Brune } 2514b11644fSPeter Brune if (qn->dF) { 2524b11644fSPeter Brune ierr = VecDestroyVecs(qn->m, &qn->dF);CHKERRQ(ierr); 2534b11644fSPeter Brune } 254*bd052dfeSPeter Brune ierr = PetscFree3(qn->alpha, qn->beta, qn->rho);CHKERRQ(ierr); 2554b11644fSPeter Brune } 2564b11644fSPeter Brune if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 2574b11644fSPeter Brune PetscFunctionReturn(0); 2584b11644fSPeter Brune } 2594b11644fSPeter Brune 2604b11644fSPeter Brune #undef __FUNCT__ 2614b11644fSPeter Brune #define __FUNCT__ "SNESDestroy_QN" 2624b11644fSPeter Brune static PetscErrorCode SNESDestroy_QN(SNES snes) 2634b11644fSPeter Brune { 2644b11644fSPeter Brune PetscErrorCode ierr; 2654b11644fSPeter Brune PetscFunctionBegin; 2664b11644fSPeter Brune ierr = SNESReset_QN(snes);CHKERRQ(ierr); 2674b11644fSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 2684b11644fSPeter Brune PetscFunctionReturn(0); 2694b11644fSPeter Brune } 2704b11644fSPeter Brune 2714b11644fSPeter Brune #undef __FUNCT__ 2724b11644fSPeter Brune #define __FUNCT__ "SNESSetFromOptions_QN" 2734b11644fSPeter Brune static PetscErrorCode SNESSetFromOptions_QN(SNES snes) 2744b11644fSPeter Brune { 2754b11644fSPeter Brune 2764b11644fSPeter Brune PetscErrorCode ierr; 2774b11644fSPeter Brune QNContext * qn; 2784b11644fSPeter Brune 2794b11644fSPeter Brune PetscFunctionBegin; 2804b11644fSPeter Brune 2814b11644fSPeter Brune qn = (QNContext *)snes->data; 2824b11644fSPeter Brune 2834b11644fSPeter Brune ierr = PetscOptionsHead("SNES QN options");CHKERRQ(ierr); 2844b11644fSPeter Brune ierr = PetscOptionsReal("-snes_qn_lambda", "SOR mixing parameter", "SNES", qn->lambda, &qn->lambda, PETSC_NULL);CHKERRQ(ierr); 2854b11644fSPeter Brune ierr = PetscOptionsInt("-snes_qn_m", "Number of past states saved for L-Broyden methods", "SNES", qn->m, &qn->m, PETSC_NULL);CHKERRQ(ierr); 2864b11644fSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 2874b11644fSPeter Brune PetscFunctionReturn(0); 2884b11644fSPeter Brune } 2894b11644fSPeter Brune 2904b11644fSPeter Brune /* -------------------------------------------------------------------------- */ 2914b11644fSPeter Brune /*MC 2924b11644fSPeter Brune SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems. 2934b11644fSPeter Brune 29435f34c5cSPeter Brune Implements a limited-memory "good" Broyden update method. 2954b11644fSPeter Brune 2964b11644fSPeter Brune Level: beginner 2974b11644fSPeter Brune 2984b11644fSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR 2994b11644fSPeter Brune M*/ 3004b11644fSPeter Brune EXTERN_C_BEGIN 3014b11644fSPeter Brune #undef __FUNCT__ 3024b11644fSPeter Brune #define __FUNCT__ "SNESCreate_QN" 3034b11644fSPeter Brune PetscErrorCode SNESCreate_QN(SNES snes) 3044b11644fSPeter Brune { 3054b11644fSPeter Brune 3064b11644fSPeter Brune PetscErrorCode ierr; 3074b11644fSPeter Brune QNContext * qn; 3084b11644fSPeter Brune 3094b11644fSPeter Brune PetscFunctionBegin; 3104b11644fSPeter Brune snes->ops->setup = SNESSetUp_QN; 3114b11644fSPeter Brune snes->ops->solve = SNESSolve_QN; 3124b11644fSPeter Brune snes->ops->destroy = SNESDestroy_QN; 3134b11644fSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_QN; 3144b11644fSPeter Brune snes->ops->view = 0; 3154b11644fSPeter Brune snes->ops->reset = SNESReset_QN; 3164b11644fSPeter Brune 3174b11644fSPeter Brune ierr = PetscNewLog(snes, QNContext, &qn);CHKERRQ(ierr); 3184b11644fSPeter Brune snes->data = (void *) qn; 319*bd052dfeSPeter Brune qn->m = 100; 3204b11644fSPeter Brune qn->lambda = 1.; 3214b11644fSPeter Brune qn->dX = PETSC_NULL; 3224b11644fSPeter Brune qn->dF = PETSC_NULL; 3234b11644fSPeter Brune PetscFunctionReturn(0); 3244b11644fSPeter Brune } 3254b11644fSPeter Brune EXTERN_C_END 326