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 PetscScalar *xi = ngmres->xi; 1238774f0aSPeter Brune PetscInt i; 1338774f0aSPeter Brune PetscReal nu; 1438774f0aSPeter Brune PetscErrorCode ierr; 1538774f0aSPeter Brune 1638774f0aSPeter Brune PetscFunctionBegin; 17ce94432eSBarry Smith if (ivec > l) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot update vector %d with space size %d!",ivec,l); 1838774f0aSPeter Brune ierr = VecCopy(F,Fdot[ivec]);CHKERRQ(ierr); 1938774f0aSPeter Brune ierr = VecCopy(X,Xdot[ivec]);CHKERRQ(ierr); 201aa26658SKarl Rupp 2138774f0aSPeter Brune ngmres->fnorms[ivec] = fnorm; 2238774f0aSPeter Brune if (l > 0) { 2338774f0aSPeter Brune ierr = VecMDot(F,l,Fdot,xi);CHKERRQ(ierr); 2438774f0aSPeter Brune for (i = 0; i < l; i++) { 2538774f0aSPeter Brune Q(i,ivec) = xi[i]; 2638774f0aSPeter Brune Q(ivec,i) = xi[i]; 2738774f0aSPeter Brune } 2838774f0aSPeter Brune } else { 2938774f0aSPeter Brune nu = fnorm*fnorm; 3038774f0aSPeter Brune Q(0,0) = nu; 3138774f0aSPeter Brune } 3238774f0aSPeter Brune PetscFunctionReturn(0); 3338774f0aSPeter Brune } 3438774f0aSPeter Brune 3538774f0aSPeter Brune #undef __FUNCT__ 3638774f0aSPeter Brune #define __FUNCT__ "SNESNGMRESFormCombinedSolution_Private" 3738774f0aSPeter Brune PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES snes,PetscInt l,Vec XM,Vec FM,PetscReal fMnorm,Vec X,Vec XA,Vec FA) 3838774f0aSPeter Brune { 3938774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 4038774f0aSPeter Brune PetscInt i,j; 4138774f0aSPeter Brune Vec *Fdot = ngmres->Fdot; 4238774f0aSPeter Brune Vec *Xdot = ngmres->Xdot; 4338774f0aSPeter Brune PetscScalar *beta = ngmres->beta; 4438774f0aSPeter Brune PetscScalar *xi = ngmres->xi; 4538774f0aSPeter Brune PetscScalar alph_total = 0.; 4638774f0aSPeter Brune PetscErrorCode ierr; 4738774f0aSPeter Brune PetscReal nu; 4838774f0aSPeter Brune Vec Y = snes->work[2]; 4938774f0aSPeter Brune PetscBool changed_y,changed_w; 5038774f0aSPeter Brune 5138774f0aSPeter Brune PetscFunctionBegin; 5238774f0aSPeter Brune nu = fMnorm*fMnorm; 5338774f0aSPeter Brune 5438774f0aSPeter Brune /* construct the right hand side and xi factors */ 5538774f0aSPeter Brune ierr = VecMDot(FM,l,Fdot,xi);CHKERRQ(ierr); 561aa26658SKarl Rupp for (i = 0; i < l; i++) beta[i] = nu - xi[i]; 571aa26658SKarl Rupp 5838774f0aSPeter Brune /* construct h */ 5938774f0aSPeter Brune for (j = 0; j < l; j++) { 6038774f0aSPeter Brune for (i = 0; i < l; i++) { 6138774f0aSPeter Brune H(i,j) = Q(i,j)-xi[i]-xi[j]+nu; 6238774f0aSPeter Brune } 6338774f0aSPeter Brune } 6438774f0aSPeter Brune if (l == 1) { 6538774f0aSPeter Brune /* simply set alpha[0] = beta[0] / H[0, 0] */ 661aa26658SKarl Rupp if (H(0,0) != 0.) beta[0] = beta[0]/H(0,0); 671aa26658SKarl Rupp else beta[0] = 0.; 6838774f0aSPeter Brune } else { 6938774f0aSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS) 70ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"NGMRES with LS requires the LAPACK GELSS routine."); 7138774f0aSPeter Brune #else 7238774f0aSPeter Brune ierr = PetscBLASIntCast(l,&ngmres->m);CHKERRQ(ierr); 7338774f0aSPeter Brune ierr = PetscBLASIntCast(l,&ngmres->n);CHKERRQ(ierr); 7438774f0aSPeter Brune ngmres->info = 0; 7538774f0aSPeter Brune ngmres->rcond = -1.; 7638774f0aSPeter Brune ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 7738774f0aSPeter Brune #if defined(PETSC_USE_COMPLEX) 78*8b83055fSJed 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)); 7938774f0aSPeter Brune #else 80*8b83055fSJed 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)); 8138774f0aSPeter Brune #endif 8238774f0aSPeter Brune ierr = PetscFPTrapPop();CHKERRQ(ierr); 83ce94432eSBarry Smith if (ngmres->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS"); 84ce94432eSBarry Smith if (ngmres->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge"); 8538774f0aSPeter Brune #endif 8638774f0aSPeter Brune } 8738774f0aSPeter Brune for (i=0; i<l; i++) { 88ce94432eSBarry Smith if (PetscIsInfOrNanScalar(beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output"); 8938774f0aSPeter Brune } 9038774f0aSPeter Brune alph_total = 0.; 911aa26658SKarl Rupp for (i = 0; i < l; i++) alph_total += beta[i]; 921aa26658SKarl Rupp 9338774f0aSPeter Brune ierr = VecCopy(XM,XA);CHKERRQ(ierr); 9438774f0aSPeter Brune ierr = VecScale(XA,1.-alph_total);CHKERRQ(ierr); 9538774f0aSPeter Brune ierr = VecMAXPY(XA,l,beta,Xdot);CHKERRQ(ierr); 9638774f0aSPeter Brune /* check the validity of the step */ 9738774f0aSPeter Brune ierr = VecCopy(XA,Y);CHKERRQ(ierr); 9838774f0aSPeter Brune ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); 9938774f0aSPeter Brune ierr = SNESLineSearchPostCheck(snes->linesearch,X,Y,XA,&changed_y,&changed_w);CHKERRQ(ierr); 100077c4231SPeter Brune if (!ngmres->approxfunc) {ierr = SNESComputeFunction(snes,XA,FA);CHKERRQ(ierr);} 101077c4231SPeter Brune else { 102077c4231SPeter Brune ierr = VecCopy(FM,FA);CHKERRQ(ierr); 103077c4231SPeter Brune ierr = VecScale(FA,1.-alph_total);CHKERRQ(ierr); 104077c4231SPeter Brune ierr = VecMAXPY(FA,l,beta,Fdot);CHKERRQ(ierr); 105077c4231SPeter Brune } 10638774f0aSPeter Brune PetscFunctionReturn(0); 10738774f0aSPeter Brune } 10838774f0aSPeter Brune 10938774f0aSPeter Brune #undef __FUNCT__ 11038774f0aSPeter Brune #define __FUNCT__ "SNESNGMRESCalculateDifferences_Private" 11138774f0aSPeter Brune PetscErrorCode SNESNGMRESCalculateDifferences_Private(SNES snes,PetscInt l,Vec X,Vec F,Vec XM,Vec FM,Vec XA,Vec FA,Vec D,PetscReal *dnorm,PetscReal *dminnorm,PetscReal *fAnorm) 11238774f0aSPeter Brune { 11338774f0aSPeter Brune PetscErrorCode ierr; 11438774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 11538774f0aSPeter Brune PetscReal dcurnorm; 11638774f0aSPeter Brune Vec *Xdot = ngmres->Xdot; 11738774f0aSPeter Brune PetscInt i; 11838774f0aSPeter Brune 11938774f0aSPeter Brune PetscFunctionBegin; 12038774f0aSPeter Brune if (ngmres->singlereduction) { 12138774f0aSPeter Brune *dminnorm = -1.0; 12238774f0aSPeter Brune if (fAnorm) { 12338774f0aSPeter Brune ierr = VecNormBegin(FA,NORM_2,fAnorm);CHKERRQ(ierr); 12438774f0aSPeter Brune } 12538774f0aSPeter Brune if (dnorm) { 12638774f0aSPeter Brune ierr = VecCopy(XA,D);CHKERRQ(ierr); 12738774f0aSPeter Brune ierr = VecAXPY(D,-1.0,XM);CHKERRQ(ierr); 12838774f0aSPeter Brune ierr = VecNormBegin(D,NORM_2,dnorm);CHKERRQ(ierr); 12938774f0aSPeter Brune } 13038774f0aSPeter Brune if (dminnorm) { 13138774f0aSPeter Brune for (i=0; i<l; i++) { 13238774f0aSPeter Brune ierr=VecAXPY(Xdot[i],-1.0,XA);CHKERRQ(ierr); 13338774f0aSPeter Brune } 13438774f0aSPeter Brune for (i=0; i<l; i++) { 13538774f0aSPeter Brune ierr = VecNormBegin(Xdot[i],NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr); 13638774f0aSPeter Brune } 13738774f0aSPeter Brune } 13838774f0aSPeter Brune if (fAnorm) {ierr = VecNormEnd(FA,NORM_2,fAnorm);CHKERRQ(ierr);} 13938774f0aSPeter Brune if (dnorm) {ierr = VecNormEnd(D,NORM_2,dnorm);CHKERRQ(ierr);} 14038774f0aSPeter Brune if (dminnorm) { 14138774f0aSPeter Brune for (i=0; i<l; i++) { 14238774f0aSPeter Brune ierr = VecNormEnd(Xdot[i],NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr); 14338774f0aSPeter Brune } 14438774f0aSPeter Brune for (i=0; i<l; i++) { 14538774f0aSPeter Brune dcurnorm = ngmres->xnorms[i]; 14638774f0aSPeter Brune if ((dcurnorm < *dminnorm) || (*dminnorm < 0.0)) *dminnorm = dcurnorm; 14738774f0aSPeter Brune ierr=VecAXPY(Xdot[i],1.0,XA);CHKERRQ(ierr); 14838774f0aSPeter Brune } 14938774f0aSPeter Brune } 15038774f0aSPeter Brune } else { 15138774f0aSPeter Brune if (dnorm) { 15238774f0aSPeter Brune ierr=VecCopy(XA,D);CHKERRQ(ierr); 15338774f0aSPeter Brune ierr=VecAXPY(D,-1.0,XM);CHKERRQ(ierr); 15438774f0aSPeter Brune ierr=VecNormBegin(D,NORM_2,dnorm);CHKERRQ(ierr); 15538774f0aSPeter Brune } 15638774f0aSPeter Brune if (fAnorm) { 15738774f0aSPeter Brune ierr=VecNormBegin(FA,NORM_2,fAnorm);CHKERRQ(ierr); 15838774f0aSPeter Brune } 15938774f0aSPeter Brune if (dnorm) { 16038774f0aSPeter Brune ierr=VecNormEnd(D,NORM_2,dnorm);CHKERRQ(ierr); 16138774f0aSPeter Brune } 16238774f0aSPeter Brune if (fAnorm) { 16338774f0aSPeter Brune ierr=VecNormEnd(FA,NORM_2,fAnorm);CHKERRQ(ierr); 16438774f0aSPeter Brune } 16538774f0aSPeter Brune if (dminnorm) { 16638774f0aSPeter Brune *dminnorm = -1.0; 16738774f0aSPeter Brune for (i=0; i<l; i++) { 16838774f0aSPeter Brune ierr=VecCopy(XA,D);CHKERRQ(ierr); 16938774f0aSPeter Brune ierr=VecAXPY(D,-1.0,Xdot[i]);CHKERRQ(ierr); 17038774f0aSPeter Brune ierr=VecNorm(D,NORM_2,&dcurnorm);CHKERRQ(ierr); 17138774f0aSPeter Brune if ((dcurnorm < *dminnorm) || (*dminnorm < 0.0)) *dminnorm = dcurnorm; 17238774f0aSPeter Brune } 17338774f0aSPeter Brune } 17438774f0aSPeter Brune } 17538774f0aSPeter Brune PetscFunctionReturn(0); 17638774f0aSPeter Brune } 17738774f0aSPeter Brune 17838774f0aSPeter Brune #undef __FUNCT__ 17938774f0aSPeter Brune #define __FUNCT__ "SNESNGMRESSelect_Private" 18038774f0aSPeter Brune PetscErrorCode SNESNGMRESSelect_Private(SNES snes,PetscInt k_restart,Vec XM,Vec FM,PetscReal fMnorm,Vec XA,Vec FA,PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,Vec X,Vec F,Vec Y,PetscReal *fnorm) 18138774f0aSPeter Brune { 18238774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 18338774f0aSPeter Brune PetscErrorCode ierr; 18438774f0aSPeter Brune PetscBool lssucceed,selectA; 18538774f0aSPeter Brune 18638774f0aSPeter Brune PetscFunctionBegin; 18738774f0aSPeter Brune if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 18838774f0aSPeter Brune /* X = X + \lambda(XA - X) */ 18938774f0aSPeter Brune if (ngmres->monitor) { 19038774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr); 19138774f0aSPeter Brune } 19238774f0aSPeter Brune ierr = VecCopy(FM,F);CHKERRQ(ierr); 19338774f0aSPeter Brune ierr = VecCopy(XM,X);CHKERRQ(ierr); 19438774f0aSPeter Brune ierr = VecCopy(XA,Y);CHKERRQ(ierr); 19538774f0aSPeter Brune ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr); 19638774f0aSPeter Brune *fnorm = fMnorm; 19738774f0aSPeter Brune ierr = SNESLineSearchApply(ngmres->additive_linesearch,X,F,fnorm,Y);CHKERRQ(ierr); 19838774f0aSPeter Brune ierr = SNESLineSearchGetSuccess(ngmres->additive_linesearch,&lssucceed);CHKERRQ(ierr); 19938774f0aSPeter Brune if (!lssucceed) { 20038774f0aSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 20138774f0aSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 20238774f0aSPeter Brune PetscFunctionReturn(0); 20338774f0aSPeter Brune } 20438774f0aSPeter Brune } 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 */ 22338774f0aSPeter Brune *fnorm = fAnorm; 22438774f0aSPeter Brune ierr = VecCopy(FA,F);CHKERRQ(ierr); 22538774f0aSPeter Brune ierr = VecCopy(XA,X);CHKERRQ(ierr); 22638774f0aSPeter Brune } else { 22738774f0aSPeter Brune if (ngmres->monitor) { 22838774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr); 22938774f0aSPeter Brune } 23038774f0aSPeter Brune *fnorm = fMnorm; 23138774f0aSPeter Brune ierr = VecCopy(XM,Y);CHKERRQ(ierr); 23238774f0aSPeter Brune ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); 23338774f0aSPeter Brune ierr = VecCopy(FM,F);CHKERRQ(ierr); 23438774f0aSPeter Brune ierr = VecCopy(XM,X);CHKERRQ(ierr); 23538774f0aSPeter Brune } 23638774f0aSPeter Brune } else { /* none */ 23738774f0aSPeter Brune *fnorm = fAnorm; 23838774f0aSPeter Brune ierr = VecCopy(FA,F);CHKERRQ(ierr); 23938774f0aSPeter Brune ierr = VecCopy(XA,X);CHKERRQ(ierr); 24038774f0aSPeter Brune } 24138774f0aSPeter Brune PetscFunctionReturn(0); 24238774f0aSPeter Brune } 24338774f0aSPeter Brune 24438774f0aSPeter Brune #undef __FUNCT__ 24538774f0aSPeter Brune #define __FUNCT__ "SNESNGMRESSelectRestart_Private" 24621687c63SPeter Brune PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes,PetscInt l,PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,PetscBool *selectRestart) 24738774f0aSPeter Brune { 24838774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 24938774f0aSPeter Brune PetscErrorCode ierr; 25038774f0aSPeter Brune 25138774f0aSPeter Brune PetscFunctionBegin; 25238774f0aSPeter Brune *selectRestart = PETSC_FALSE; 25338774f0aSPeter Brune /* difference stagnation restart */ 25421687c63SPeter Brune if ((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm)) && l > 0) { 25538774f0aSPeter Brune if (ngmres->monitor) { 25638774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"difference restart: %e > %e\n",ngmres->epsilonB*dnorm,dminnorm);CHKERRQ(ierr); 25738774f0aSPeter Brune } 25838774f0aSPeter Brune *selectRestart = PETSC_TRUE; 25938774f0aSPeter Brune } 26038774f0aSPeter Brune /* residual stagnation restart */ 26138774f0aSPeter Brune if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 26238774f0aSPeter Brune if (ngmres->monitor) { 26338774f0aSPeter Brune ierr = PetscViewerASCIIPrintf(ngmres->monitor,"residual restart: %e > %e\n",PetscSqrtReal(fAnorm),ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr); 26438774f0aSPeter Brune } 26538774f0aSPeter Brune *selectRestart = PETSC_TRUE; 26638774f0aSPeter Brune } 26738774f0aSPeter Brune PetscFunctionReturn(0); 26838774f0aSPeter Brune } 269