138774f0aSPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 238774f0aSPeter Brune #include <petscblaslapack.h> 338774f0aSPeter Brune 438774f0aSPeter Brune #undef __FUNCT__ 538774f0aSPeter Brune #define __FUNCT__ "SNESNGMRESUpdateSubspace_Private" 638774f0aSPeter Brune PetscErrorCode SNESNGMRESUpdateSubspace_Private(SNES snes,PetscInt ivec,PetscInt l,Vec F,PetscReal fnorm,Vec X) 738774f0aSPeter Brune { 838774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 938774f0aSPeter Brune Vec *Fdot = ngmres->Fdot; 1038774f0aSPeter Brune Vec *Xdot = ngmres->Xdot; 1138774f0aSPeter Brune PetscErrorCode ierr; 1238774f0aSPeter Brune 1338774f0aSPeter Brune PetscFunctionBegin; 14ce94432eSBarry Smith if (ivec > l) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot update vector %d with space size %d!",ivec,l); 1538774f0aSPeter Brune ierr = VecCopy(F,Fdot[ivec]);CHKERRQ(ierr); 1638774f0aSPeter Brune ierr = VecCopy(X,Xdot[ivec]);CHKERRQ(ierr); 171aa26658SKarl Rupp 1838774f0aSPeter Brune ngmres->fnorms[ivec] = fnorm; 1938774f0aSPeter Brune PetscFunctionReturn(0); 2038774f0aSPeter Brune } 2138774f0aSPeter Brune 2238774f0aSPeter Brune #undef __FUNCT__ 2338774f0aSPeter Brune #define __FUNCT__ "SNESNGMRESFormCombinedSolution_Private" 24*b3c6a99cSPeter Brune PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES snes,PetscInt ivec,PetscInt l,Vec XM,Vec FM,PetscReal fMnorm,Vec X,Vec XA,Vec FA) 2538774f0aSPeter Brune { 2638774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 2738774f0aSPeter Brune PetscInt i,j; 2838774f0aSPeter Brune Vec *Fdot = ngmres->Fdot; 2938774f0aSPeter Brune Vec *Xdot = ngmres->Xdot; 3038774f0aSPeter Brune PetscScalar *beta = ngmres->beta; 3138774f0aSPeter Brune PetscScalar *xi = ngmres->xi; 3238774f0aSPeter Brune PetscScalar alph_total = 0.; 3338774f0aSPeter Brune PetscErrorCode ierr; 3438774f0aSPeter Brune PetscReal nu; 3538774f0aSPeter Brune Vec Y = snes->work[2]; 3638774f0aSPeter Brune PetscBool changed_y,changed_w; 3738774f0aSPeter Brune 3838774f0aSPeter Brune PetscFunctionBegin; 3938774f0aSPeter Brune nu = fMnorm*fMnorm; 4038774f0aSPeter Brune 4138774f0aSPeter Brune /* construct the right hand side and xi factors */ 42*b3c6a99cSPeter Brune if (l > 0) { 43*b3c6a99cSPeter Brune ierr = VecMDotBegin(FM,l,Fdot,xi);CHKERRQ(ierr); 44*b3c6a99cSPeter Brune ierr = VecMDotBegin(Fdot[ivec],l,Fdot,beta);CHKERRQ(ierr); 45*b3c6a99cSPeter Brune ierr = VecMDotEnd(FM,l,Fdot,xi);CHKERRQ(ierr); 46*b3c6a99cSPeter Brune ierr = VecMDotEnd(Fdot[ivec],l,Fdot,beta);CHKERRQ(ierr); 47*b3c6a99cSPeter Brune for (i = 0; i < l; i++) { 48*b3c6a99cSPeter Brune Q(i,ivec) = beta[i]; 49*b3c6a99cSPeter Brune Q(ivec,i) = beta[i]; 50*b3c6a99cSPeter Brune } 51*b3c6a99cSPeter Brune } else { 52*b3c6a99cSPeter Brune Q(0,0) = ngmres->fnorms[ivec]*ngmres->fnorms[ivec]; 53*b3c6a99cSPeter Brune } 54*b3c6a99cSPeter Brune 551aa26658SKarl Rupp for (i = 0; i < l; i++) beta[i] = nu - xi[i]; 561aa26658SKarl Rupp 5738774f0aSPeter Brune /* construct h */ 5838774f0aSPeter Brune for (j = 0; j < l; j++) { 5938774f0aSPeter Brune for (i = 0; i < l; i++) { 6038774f0aSPeter Brune H(i,j) = Q(i,j)-xi[i]-xi[j]+nu; 6138774f0aSPeter Brune } 6238774f0aSPeter Brune } 6338774f0aSPeter Brune if (l == 1) { 6438774f0aSPeter Brune /* simply set alpha[0] = beta[0] / H[0, 0] */ 651aa26658SKarl Rupp if (H(0,0) != 0.) beta[0] = beta[0]/H(0,0); 661aa26658SKarl Rupp else beta[0] = 0.; 6738774f0aSPeter Brune } else { 6838774f0aSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS) 69ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"NGMRES with LS requires the LAPACK GELSS routine."); 7038774f0aSPeter Brune #else 7138774f0aSPeter Brune ierr = PetscBLASIntCast(l,&ngmres->m);CHKERRQ(ierr); 7238774f0aSPeter Brune ierr = PetscBLASIntCast(l,&ngmres->n);CHKERRQ(ierr); 7338774f0aSPeter Brune ngmres->info = 0; 7438774f0aSPeter Brune ngmres->rcond = -1.; 7538774f0aSPeter Brune ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 7638774f0aSPeter Brune #if defined(PETSC_USE_COMPLEX) 778b83055fSJed Brown PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&ngmres->m,&ngmres->n,&ngmres->nrhs,ngmres->h,&ngmres->lda,ngmres->beta,&ngmres->ldb,ngmres->s,&ngmres->rcond,&ngmres->rank,ngmres->work,&ngmres->lwork,ngmres->rwork,&ngmres->info)); 7838774f0aSPeter Brune #else 798b83055fSJed Brown PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&ngmres->m,&ngmres->n,&ngmres->nrhs,ngmres->h,&ngmres->lda,ngmres->beta,&ngmres->ldb,ngmres->s,&ngmres->rcond,&ngmres->rank,ngmres->work,&ngmres->lwork,&ngmres->info)); 8038774f0aSPeter Brune #endif 8138774f0aSPeter Brune ierr = PetscFPTrapPop();CHKERRQ(ierr); 82ce94432eSBarry Smith if (ngmres->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS"); 83ce94432eSBarry Smith if (ngmres->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge"); 8438774f0aSPeter Brune #endif 8538774f0aSPeter Brune } 8638774f0aSPeter Brune for (i=0; i<l; i++) { 87ce94432eSBarry Smith if (PetscIsInfOrNanScalar(beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output"); 8838774f0aSPeter Brune } 8938774f0aSPeter Brune alph_total = 0.; 901aa26658SKarl Rupp for (i = 0; i < l; i++) alph_total += beta[i]; 911aa26658SKarl Rupp 9238774f0aSPeter Brune ierr = VecCopy(XM,XA);CHKERRQ(ierr); 9338774f0aSPeter Brune ierr = VecScale(XA,1.-alph_total);CHKERRQ(ierr); 9438774f0aSPeter Brune ierr = VecMAXPY(XA,l,beta,Xdot);CHKERRQ(ierr); 9538774f0aSPeter Brune /* check the validity of the step */ 9638774f0aSPeter Brune ierr = VecCopy(XA,Y);CHKERRQ(ierr); 9738774f0aSPeter Brune ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); 9838774f0aSPeter Brune ierr = SNESLineSearchPostCheck(snes->linesearch,X,Y,XA,&changed_y,&changed_w);CHKERRQ(ierr); 9946159c86SPeter Brune if (!ngmres->approxfunc) { 10046159c86SPeter Brune if (snes->pc && snes->pcside == PC_LEFT) { 10146159c86SPeter Brune ierr = SNESApplyPC(snes,XA,NULL,NULL,FA);CHKERRQ(ierr); 10246159c86SPeter Brune } else { 10346159c86SPeter Brune ierr =SNESComputeFunction(snes,XA,FA);CHKERRQ(ierr); 10446159c86SPeter Brune } 10546159c86SPeter Brune } 106077c4231SPeter Brune else { 107077c4231SPeter Brune ierr = VecCopy(FM,FA);CHKERRQ(ierr); 108077c4231SPeter Brune ierr = VecScale(FA,1.-alph_total);CHKERRQ(ierr); 109077c4231SPeter Brune ierr = VecMAXPY(FA,l,beta,Fdot);CHKERRQ(ierr); 110077c4231SPeter Brune } 11138774f0aSPeter Brune PetscFunctionReturn(0); 11238774f0aSPeter Brune } 11338774f0aSPeter Brune 11438774f0aSPeter Brune #undef __FUNCT__ 115*b3c6a99cSPeter Brune #define __FUNCT__ "SNESNGMRESNorms_Private" 116*b3c6a99cSPeter Brune PetscErrorCode SNESNGMRESNorms_Private(SNES snes,PetscInt l,Vec X,Vec F,Vec XM,Vec FM,Vec XA,Vec FA,Vec D,PetscReal *dnorm,PetscReal *dminnorm,PetscReal *xMnorm,PetscReal *fMnorm,PetscReal *yMnorm, PetscReal *xAnorm,PetscReal *fAnorm,PetscReal *yAnorm) 11738774f0aSPeter Brune { 11838774f0aSPeter Brune PetscErrorCode ierr; 11938774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 120*b3c6a99cSPeter Brune PetscReal dcurnorm,dmin = -1.0; 12138774f0aSPeter Brune Vec *Xdot = ngmres->Xdot; 12238774f0aSPeter Brune PetscInt i; 12338774f0aSPeter Brune 12438774f0aSPeter Brune PetscFunctionBegin; 125*b3c6a99cSPeter Brune if (xMnorm) { 126*b3c6a99cSPeter Brune ierr = VecNormBegin(XM,NORM_2,xMnorm);CHKERRQ(ierr); 127*b3c6a99cSPeter Brune } 128*b3c6a99cSPeter Brune if (fMnorm) { 129*b3c6a99cSPeter Brune ierr = VecNormBegin(FM,NORM_2,fMnorm);CHKERRQ(ierr); 130*b3c6a99cSPeter Brune } 131*b3c6a99cSPeter Brune if (yMnorm) { 132*b3c6a99cSPeter Brune ierr = VecCopy(X,D);CHKERRQ(ierr); 133*b3c6a99cSPeter Brune ierr = VecAXPY(D,-1.0,XM);CHKERRQ(ierr); 134*b3c6a99cSPeter Brune ierr = VecNormBegin(D,NORM_2,yMnorm);CHKERRQ(ierr); 135*b3c6a99cSPeter Brune } 136*b3c6a99cSPeter Brune if (xAnorm) { 137*b3c6a99cSPeter Brune ierr = VecNormBegin(XA,NORM_2,xAnorm);CHKERRQ(ierr); 138*b3c6a99cSPeter Brune } 13938774f0aSPeter Brune if (fAnorm) { 14038774f0aSPeter Brune ierr = VecNormBegin(FA,NORM_2,fAnorm);CHKERRQ(ierr); 14138774f0aSPeter Brune } 142*b3c6a99cSPeter Brune if (yAnorm) { 143*b3c6a99cSPeter Brune ierr = VecCopy(X,D);CHKERRQ(ierr); 144*b3c6a99cSPeter Brune ierr = VecAXPY(D,-1.0,XA);CHKERRQ(ierr); 145*b3c6a99cSPeter Brune ierr = VecNormBegin(D,NORM_2,yAnorm);CHKERRQ(ierr); 146*b3c6a99cSPeter Brune } 14738774f0aSPeter Brune if (dnorm) { 14838774f0aSPeter Brune ierr = VecCopy(XA,D);CHKERRQ(ierr); 14938774f0aSPeter Brune ierr = VecAXPY(D,-1.0,XM);CHKERRQ(ierr); 15038774f0aSPeter Brune ierr = VecNormBegin(D,NORM_2,dnorm);CHKERRQ(ierr); 15138774f0aSPeter Brune } 15238774f0aSPeter Brune if (dminnorm) { 15338774f0aSPeter Brune for (i=0; i<l; i++) { 154*b3c6a99cSPeter Brune ierr = VecCopy(Xdot[i],D);CHKERRQ(ierr); 155*b3c6a99cSPeter Brune ierr = VecAXPY(D,-1.0,XA);CHKERRQ(ierr); 156*b3c6a99cSPeter Brune ierr = VecNormBegin(D,NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr); 15738774f0aSPeter Brune } 15838774f0aSPeter Brune } 159*b3c6a99cSPeter Brune if (xMnorm) {ierr = VecNormEnd(XM,NORM_2,xMnorm);CHKERRQ(ierr);} 160*b3c6a99cSPeter Brune if (fMnorm) {ierr = VecNormEnd(FM,NORM_2,fMnorm);CHKERRQ(ierr);} 161*b3c6a99cSPeter Brune if (yMnorm) {ierr = VecNormEnd(D,NORM_2,yMnorm);CHKERRQ(ierr);} 162*b3c6a99cSPeter Brune if (xAnorm) {ierr = VecNormEnd(XA,NORM_2,xAnorm);CHKERRQ(ierr);} 16338774f0aSPeter Brune if (fAnorm) {ierr = VecNormEnd(FA,NORM_2,fAnorm);CHKERRQ(ierr);} 164*b3c6a99cSPeter Brune if (yAnorm) {ierr = VecNormEnd(D,NORM_2,yAnorm);CHKERRQ(ierr);} 16538774f0aSPeter Brune if (dnorm) {ierr = VecNormEnd(D,NORM_2,dnorm);CHKERRQ(ierr);} 16638774f0aSPeter Brune if (dminnorm) { 16738774f0aSPeter Brune for (i=0; i<l; i++) { 168*b3c6a99cSPeter Brune ierr = VecNormEnd(D,NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr); 16938774f0aSPeter Brune dcurnorm = ngmres->xnorms[i]; 170*b3c6a99cSPeter Brune if ((dcurnorm < dmin) || (dmin < 0.0)) dmin = dcurnorm; 17138774f0aSPeter Brune } 172*b3c6a99cSPeter Brune *dminnorm = dmin; 17338774f0aSPeter Brune } 17438774f0aSPeter Brune PetscFunctionReturn(0); 17538774f0aSPeter Brune } 17638774f0aSPeter Brune 17738774f0aSPeter Brune #undef __FUNCT__ 17838774f0aSPeter Brune #define __FUNCT__ "SNESNGMRESSelect_Private" 179*b3c6a99cSPeter Brune PetscErrorCode SNESNGMRESSelect_Private(SNES snes,PetscInt k_restart,Vec XM,Vec FM,PetscReal xMnorm,PetscReal fMnorm,PetscReal yMnorm,Vec XA,Vec FA,PetscReal xAnorm,PetscReal fAnorm,PetscReal yAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,Vec X,Vec F,Vec Y,PetscReal *xnorm,PetscReal *fnorm,PetscReal *ynorm) 18038774f0aSPeter Brune { 18138774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 18238774f0aSPeter Brune PetscErrorCode ierr; 18338774f0aSPeter Brune PetscBool lssucceed,selectA; 18438774f0aSPeter Brune 18538774f0aSPeter Brune PetscFunctionBegin; 18638774f0aSPeter Brune if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 18738774f0aSPeter Brune /* X = X + \lambda(XA - X) */ 18838774f0aSPeter Brune if (ngmres->monitor) { 18938774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr); 19038774f0aSPeter Brune } 19138774f0aSPeter Brune ierr = VecCopy(FM,F);CHKERRQ(ierr); 19238774f0aSPeter Brune ierr = VecCopy(XM,X);CHKERRQ(ierr); 19338774f0aSPeter Brune ierr = VecCopy(XA,Y);CHKERRQ(ierr); 19438774f0aSPeter Brune ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 19538774f0aSPeter Brune *fnorm = fMnorm; 19638774f0aSPeter Brune ierr = SNESLineSearchApply(ngmres->additive_linesearch,X,F,fnorm,Y);CHKERRQ(ierr); 19738774f0aSPeter Brune ierr = SNESLineSearchGetSuccess(ngmres->additive_linesearch,&lssucceed);CHKERRQ(ierr); 19838774f0aSPeter Brune if (!lssucceed) { 19938774f0aSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 20038774f0aSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 20138774f0aSPeter Brune PetscFunctionReturn(0); 20238774f0aSPeter Brune } 20338774f0aSPeter Brune } 204*b3c6a99cSPeter Brune ierr = SNESLineSearchGetNorms(ngmres->additive_linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr); 20538774f0aSPeter Brune if (ngmres->monitor) { 20638774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Additive solution: ||F||_2 = %e\n",*fnorm);CHKERRQ(ierr); 20738774f0aSPeter Brune } 20838774f0aSPeter Brune } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 20938774f0aSPeter Brune selectA = PETSC_TRUE; 21038774f0aSPeter Brune /* Conditions for choosing the accelerated answer */ 21138774f0aSPeter Brune /* Criterion A -- the norm of the function isn't increased above the minimum by too much */ 2121aa26658SKarl Rupp if (fAnorm >= ngmres->gammaA*fminnorm) selectA = PETSC_FALSE; 2131aa26658SKarl Rupp 21438774f0aSPeter Brune /* Criterion B -- the choice of x^A isn't too close to some other choice */ 21538774f0aSPeter Brune if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(*fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) { 2161aa26658SKarl Rupp } else selectA=PETSC_FALSE; 2171aa26658SKarl Rupp 21838774f0aSPeter Brune if (selectA) { 21938774f0aSPeter Brune if (ngmres->monitor) { 22038774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr); 22138774f0aSPeter Brune } 22238774f0aSPeter Brune /* copy it over */ 223*b3c6a99cSPeter Brune *xnorm = xAnorm; 22438774f0aSPeter Brune *fnorm = fAnorm; 225*b3c6a99cSPeter Brune *ynorm = yAnorm; 22638774f0aSPeter Brune ierr = VecCopy(FA,F);CHKERRQ(ierr); 22738774f0aSPeter Brune ierr = VecCopy(XA,X);CHKERRQ(ierr); 22838774f0aSPeter Brune } else { 22938774f0aSPeter Brune if (ngmres->monitor) { 23038774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr); 23138774f0aSPeter Brune } 232*b3c6a99cSPeter Brune *xnorm = xMnorm; 23338774f0aSPeter Brune *fnorm = fMnorm; 234*b3c6a99cSPeter Brune *ynorm = yMnorm; 23538774f0aSPeter Brune ierr = VecCopy(XM,Y);CHKERRQ(ierr); 23638774f0aSPeter Brune ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); 23738774f0aSPeter Brune ierr = VecCopy(FM,F);CHKERRQ(ierr); 23838774f0aSPeter Brune ierr = VecCopy(XM,X);CHKERRQ(ierr); 23938774f0aSPeter Brune } 24038774f0aSPeter Brune } else { /* none */ 241*b3c6a99cSPeter Brune *xnorm = xAnorm; 24238774f0aSPeter Brune *fnorm = fAnorm; 243*b3c6a99cSPeter Brune *ynorm = yAnorm; 24438774f0aSPeter Brune ierr = VecCopy(FA,F);CHKERRQ(ierr); 24538774f0aSPeter Brune ierr = VecCopy(XA,X);CHKERRQ(ierr); 24638774f0aSPeter Brune } 24738774f0aSPeter Brune PetscFunctionReturn(0); 24838774f0aSPeter Brune } 24938774f0aSPeter Brune 25038774f0aSPeter Brune #undef __FUNCT__ 25138774f0aSPeter Brune #define __FUNCT__ "SNESNGMRESSelectRestart_Private" 25221687c63SPeter Brune PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes,PetscInt l,PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,PetscBool *selectRestart) 25338774f0aSPeter Brune { 25438774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 25538774f0aSPeter Brune PetscErrorCode ierr; 25638774f0aSPeter Brune 25738774f0aSPeter Brune PetscFunctionBegin; 25838774f0aSPeter Brune *selectRestart = PETSC_FALSE; 25938774f0aSPeter Brune /* difference stagnation restart */ 26021687c63SPeter Brune if ((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm)) && l > 0) { 26138774f0aSPeter Brune if (ngmres->monitor) { 26238774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"difference restart: %e > %e\n",ngmres->epsilonB*dnorm,dminnorm);CHKERRQ(ierr); 26338774f0aSPeter Brune } 26438774f0aSPeter Brune *selectRestart = PETSC_TRUE; 26538774f0aSPeter Brune } 26638774f0aSPeter Brune /* residual stagnation restart */ 26738774f0aSPeter Brune if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 26838774f0aSPeter Brune if (ngmres->monitor) { 26938774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"residual restart: %e > %e\n",PetscSqrtReal(fAnorm),ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr); 27038774f0aSPeter Brune } 27138774f0aSPeter Brune *selectRestart = PETSC_TRUE; 27238774f0aSPeter Brune } 27338774f0aSPeter Brune PetscFunctionReturn(0); 27438774f0aSPeter Brune } 275