138774f0aSPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/ 238774f0aSPeter Brune #include <petscblaslapack.h> 338774f0aSPeter Brune 438774f0aSPeter Brune PetscErrorCode SNESNGMRESUpdateSubspace_Private(SNES snes,PetscInt ivec,PetscInt l,Vec F,PetscReal fnorm,Vec X) 538774f0aSPeter Brune { 638774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 738774f0aSPeter Brune Vec *Fdot = ngmres->Fdot; 838774f0aSPeter Brune Vec *Xdot = ngmres->Xdot; 938774f0aSPeter Brune 1038774f0aSPeter Brune PetscFunctionBegin; 1163a3b9bcSJacob Faibussowitsch PetscCheck(ivec <= l,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot update vector %" PetscInt_FMT " with space size %" PetscInt_FMT "!",ivec,l); 129566063dSJacob Faibussowitsch PetscCall(VecCopy(F,Fdot[ivec])); 139566063dSJacob Faibussowitsch PetscCall(VecCopy(X,Xdot[ivec])); 141aa26658SKarl Rupp 1538774f0aSPeter Brune ngmres->fnorms[ivec] = fnorm; 1638774f0aSPeter Brune PetscFunctionReturn(0); 1738774f0aSPeter Brune } 1838774f0aSPeter Brune 19b3c6a99cSPeter Brune PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES snes,PetscInt ivec,PetscInt l,Vec XM,Vec FM,PetscReal fMnorm,Vec X,Vec XA,Vec FA) 2038774f0aSPeter Brune { 2138774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 2238774f0aSPeter Brune PetscInt i,j; 2338774f0aSPeter Brune Vec *Fdot = ngmres->Fdot; 2438774f0aSPeter Brune Vec *Xdot = ngmres->Xdot; 2538774f0aSPeter Brune PetscScalar *beta = ngmres->beta; 2638774f0aSPeter Brune PetscScalar *xi = ngmres->xi; 2738774f0aSPeter Brune PetscScalar alph_total = 0.; 2838774f0aSPeter Brune PetscReal nu; 2938774f0aSPeter Brune Vec Y = snes->work[2]; 3038774f0aSPeter Brune PetscBool changed_y,changed_w; 3138774f0aSPeter Brune 3238774f0aSPeter Brune PetscFunctionBegin; 3338774f0aSPeter Brune nu = fMnorm*fMnorm; 3438774f0aSPeter Brune 3538774f0aSPeter Brune /* construct the right hand side and xi factors */ 36b3c6a99cSPeter Brune if (l > 0) { 379566063dSJacob Faibussowitsch PetscCall(VecMDotBegin(FM,l,Fdot,xi)); 389566063dSJacob Faibussowitsch PetscCall(VecMDotBegin(Fdot[ivec],l,Fdot,beta)); 399566063dSJacob Faibussowitsch PetscCall(VecMDotEnd(FM,l,Fdot,xi)); 409566063dSJacob Faibussowitsch PetscCall(VecMDotEnd(Fdot[ivec],l,Fdot,beta)); 41b3c6a99cSPeter Brune for (i = 0; i < l; i++) { 42b3c6a99cSPeter Brune Q(i,ivec) = beta[i]; 43b3c6a99cSPeter Brune Q(ivec,i) = beta[i]; 44b3c6a99cSPeter Brune } 45b3c6a99cSPeter Brune } else { 46b3c6a99cSPeter Brune Q(0,0) = ngmres->fnorms[ivec]*ngmres->fnorms[ivec]; 47b3c6a99cSPeter Brune } 48b3c6a99cSPeter Brune 491aa26658SKarl Rupp for (i = 0; i < l; i++) beta[i] = nu - xi[i]; 501aa26658SKarl Rupp 5138774f0aSPeter Brune /* construct h */ 5238774f0aSPeter Brune for (j = 0; j < l; j++) { 5338774f0aSPeter Brune for (i = 0; i < l; i++) { 5438774f0aSPeter Brune H(i,j) = Q(i,j)-xi[i]-xi[j]+nu; 5538774f0aSPeter Brune } 5638774f0aSPeter Brune } 5738774f0aSPeter Brune if (l == 1) { 5838774f0aSPeter Brune /* simply set alpha[0] = beta[0] / H[0, 0] */ 591aa26658SKarl Rupp if (H(0,0) != 0.) beta[0] = beta[0]/H(0,0); 601aa26658SKarl Rupp else beta[0] = 0.; 6138774f0aSPeter Brune } else { 629566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(l,&ngmres->m)); 639566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(l,&ngmres->n)); 6438774f0aSPeter Brune ngmres->info = 0; 6538774f0aSPeter Brune ngmres->rcond = -1.; 669566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 6738774f0aSPeter Brune #if defined(PETSC_USE_COMPLEX) 68*792fecdfSBarry Smith PetscCallBLAS("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)); 6938774f0aSPeter Brune #else 70*792fecdfSBarry Smith PetscCallBLAS("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)); 7138774f0aSPeter Brune #endif 729566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 7308401ef6SPierre Jolivet PetscCheck(ngmres->info >= 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS"); 7408401ef6SPierre Jolivet PetscCheck(ngmres->info <= 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge"); 7538774f0aSPeter Brune } 7638774f0aSPeter Brune for (i=0; i<l; i++) { 770b121fc5SBarry Smith PetscCheck(!PetscIsInfOrNanScalar(beta[i]),PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output"); 7838774f0aSPeter Brune } 7938774f0aSPeter Brune alph_total = 0.; 801aa26658SKarl Rupp for (i = 0; i < l; i++) alph_total += beta[i]; 811aa26658SKarl Rupp 829566063dSJacob Faibussowitsch PetscCall(VecCopy(XM,XA)); 839566063dSJacob Faibussowitsch PetscCall(VecScale(XA,1.-alph_total)); 849566063dSJacob Faibussowitsch PetscCall(VecMAXPY(XA,l,beta,Xdot)); 8538774f0aSPeter Brune /* check the validity of the step */ 869566063dSJacob Faibussowitsch PetscCall(VecCopy(XA,Y)); 879566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y,-1.0,X)); 889566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPostCheck(snes->linesearch,X,Y,XA,&changed_y,&changed_w)); 8946159c86SPeter Brune if (!ngmres->approxfunc) { 90efd4aadfSBarry Smith if (snes->npc && snes->npcside== PC_LEFT) { 919566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes,XA,NULL,FA)); 9246159c86SPeter Brune } else { 939566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,XA,FA)); 9446159c86SPeter Brune } 9535f895b4SBarry Smith } else { 969566063dSJacob Faibussowitsch PetscCall(VecCopy(FM,FA)); 979566063dSJacob Faibussowitsch PetscCall(VecScale(FA,1.-alph_total)); 989566063dSJacob Faibussowitsch PetscCall(VecMAXPY(FA,l,beta,Fdot)); 99077c4231SPeter Brune } 10038774f0aSPeter Brune PetscFunctionReturn(0); 10138774f0aSPeter Brune } 10238774f0aSPeter Brune 103b3c6a99cSPeter 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) 10438774f0aSPeter Brune { 10538774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 106b3c6a99cSPeter Brune PetscReal dcurnorm,dmin = -1.0; 10738774f0aSPeter Brune Vec *Xdot = ngmres->Xdot; 10838774f0aSPeter Brune PetscInt i; 10938774f0aSPeter Brune 11038774f0aSPeter Brune PetscFunctionBegin; 1111baa6e33SBarry Smith if (xMnorm) PetscCall(VecNormBegin(XM,NORM_2,xMnorm)); 1121baa6e33SBarry Smith if (fMnorm) PetscCall(VecNormBegin(FM,NORM_2,fMnorm)); 113b3c6a99cSPeter Brune if (yMnorm) { 1149566063dSJacob Faibussowitsch PetscCall(VecCopy(X,D)); 1159566063dSJacob Faibussowitsch PetscCall(VecAXPY(D,-1.0,XM)); 1169566063dSJacob Faibussowitsch PetscCall(VecNormBegin(D,NORM_2,yMnorm)); 117b3c6a99cSPeter Brune } 1181baa6e33SBarry Smith if (xAnorm) PetscCall(VecNormBegin(XA,NORM_2,xAnorm)); 1191baa6e33SBarry Smith if (fAnorm) PetscCall(VecNormBegin(FA,NORM_2,fAnorm)); 120b3c6a99cSPeter Brune if (yAnorm) { 1219566063dSJacob Faibussowitsch PetscCall(VecCopy(X,D)); 1229566063dSJacob Faibussowitsch PetscCall(VecAXPY(D,-1.0,XA)); 1239566063dSJacob Faibussowitsch PetscCall(VecNormBegin(D,NORM_2,yAnorm)); 124b3c6a99cSPeter Brune } 12538774f0aSPeter Brune if (dnorm) { 1269566063dSJacob Faibussowitsch PetscCall(VecCopy(XA,D)); 1279566063dSJacob Faibussowitsch PetscCall(VecAXPY(D,-1.0,XM)); 1289566063dSJacob Faibussowitsch PetscCall(VecNormBegin(D,NORM_2,dnorm)); 12938774f0aSPeter Brune } 13038774f0aSPeter Brune if (dminnorm) { 13138774f0aSPeter Brune for (i=0; i<l; i++) { 1329566063dSJacob Faibussowitsch PetscCall(VecCopy(Xdot[i],D)); 1339566063dSJacob Faibussowitsch PetscCall(VecAXPY(D,-1.0,XA)); 1349566063dSJacob Faibussowitsch PetscCall(VecNormBegin(D,NORM_2,&ngmres->xnorms[i])); 13538774f0aSPeter Brune } 13638774f0aSPeter Brune } 1379566063dSJacob Faibussowitsch if (xMnorm) PetscCall(VecNormEnd(XM,NORM_2,xMnorm)); 1389566063dSJacob Faibussowitsch if (fMnorm) PetscCall(VecNormEnd(FM,NORM_2,fMnorm)); 1399566063dSJacob Faibussowitsch if (yMnorm) PetscCall(VecNormEnd(D,NORM_2,yMnorm)); 1409566063dSJacob Faibussowitsch if (xAnorm) PetscCall(VecNormEnd(XA,NORM_2,xAnorm)); 1419566063dSJacob Faibussowitsch if (fAnorm) PetscCall(VecNormEnd(FA,NORM_2,fAnorm)); 1429566063dSJacob Faibussowitsch if (yAnorm) PetscCall(VecNormEnd(D,NORM_2,yAnorm)); 1439566063dSJacob Faibussowitsch if (dnorm) PetscCall(VecNormEnd(D,NORM_2,dnorm)); 14438774f0aSPeter Brune if (dminnorm) { 14538774f0aSPeter Brune for (i=0; i<l; i++) { 1469566063dSJacob Faibussowitsch PetscCall(VecNormEnd(D,NORM_2,&ngmres->xnorms[i])); 14738774f0aSPeter Brune dcurnorm = ngmres->xnorms[i]; 148b3c6a99cSPeter Brune if ((dcurnorm < dmin) || (dmin < 0.0)) dmin = dcurnorm; 14938774f0aSPeter Brune } 150b3c6a99cSPeter Brune *dminnorm = dmin; 15138774f0aSPeter Brune } 15238774f0aSPeter Brune PetscFunctionReturn(0); 15338774f0aSPeter Brune } 15438774f0aSPeter Brune 155b3c6a99cSPeter 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) 15638774f0aSPeter Brune { 15738774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data; 158422a814eSBarry Smith SNESLineSearchReason lssucceed; 159422a814eSBarry Smith PetscBool selectA; 16038774f0aSPeter Brune 16138774f0aSPeter Brune PetscFunctionBegin; 16238774f0aSPeter Brune if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 16338774f0aSPeter Brune /* X = X + \lambda(XA - X) */ 16438774f0aSPeter Brune if (ngmres->monitor) { 16563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"||F_A||_2 = %e, ||F_M||_2 = %e\n",(double)fAnorm,(double)fMnorm)); 16638774f0aSPeter Brune } 1679566063dSJacob Faibussowitsch PetscCall(VecCopy(FM,F)); 1689566063dSJacob Faibussowitsch PetscCall(VecCopy(XM,X)); 1699566063dSJacob Faibussowitsch PetscCall(VecCopy(XA,Y)); 1709566063dSJacob Faibussowitsch PetscCall(VecAYPX(Y,-1.0,X)); 17138774f0aSPeter Brune *fnorm = fMnorm; 1729566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(ngmres->additive_linesearch,X,F,fnorm,Y)); 1739566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(ngmres->additive_linesearch,&lssucceed)); 1749566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(ngmres->additive_linesearch,xnorm,fnorm,ynorm)); 175422a814eSBarry Smith if (lssucceed) { 17638774f0aSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 17738774f0aSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 17838774f0aSPeter Brune PetscFunctionReturn(0); 17938774f0aSPeter Brune } 18038774f0aSPeter Brune } 18138774f0aSPeter Brune if (ngmres->monitor) { 18263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"Additive solution: ||F||_2 = %e\n",(double)*fnorm)); 18338774f0aSPeter Brune } 18438774f0aSPeter Brune } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 18538774f0aSPeter Brune selectA = PETSC_TRUE; 18638774f0aSPeter Brune /* Conditions for choosing the accelerated answer */ 18738774f0aSPeter Brune /* Criterion A -- the norm of the function isn't increased above the minimum by too much */ 1881aa26658SKarl Rupp if (fAnorm >= ngmres->gammaA*fminnorm) selectA = PETSC_FALSE; 1891aa26658SKarl Rupp 19038774f0aSPeter Brune /* Criterion B -- the choice of x^A isn't too close to some other choice */ 19138774f0aSPeter Brune if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(*fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) { 1921aa26658SKarl Rupp } else selectA=PETSC_FALSE; 1931aa26658SKarl Rupp 19438774f0aSPeter Brune if (selectA) { 19538774f0aSPeter Brune if (ngmres->monitor) { 19663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n",(double)fAnorm,(double)fMnorm)); 19738774f0aSPeter Brune } 19838774f0aSPeter Brune /* copy it over */ 199b3c6a99cSPeter Brune *xnorm = xAnorm; 20038774f0aSPeter Brune *fnorm = fAnorm; 201b3c6a99cSPeter Brune *ynorm = yAnorm; 2029566063dSJacob Faibussowitsch PetscCall(VecCopy(FA,F)); 2039566063dSJacob Faibussowitsch PetscCall(VecCopy(XA,X)); 20438774f0aSPeter Brune } else { 20538774f0aSPeter Brune if (ngmres->monitor) { 20663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n",(double)fAnorm,(double)fMnorm)); 20738774f0aSPeter Brune } 208b3c6a99cSPeter Brune *xnorm = xMnorm; 20938774f0aSPeter Brune *fnorm = fMnorm; 210b3c6a99cSPeter Brune *ynorm = yMnorm; 2119566063dSJacob Faibussowitsch PetscCall(VecCopy(XM,Y)); 2129566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y,-1.0,X)); 2139566063dSJacob Faibussowitsch PetscCall(VecCopy(FM,F)); 2149566063dSJacob Faibussowitsch PetscCall(VecCopy(XM,X)); 21538774f0aSPeter Brune } 21638774f0aSPeter Brune } else { /* none */ 217b3c6a99cSPeter Brune *xnorm = xAnorm; 21838774f0aSPeter Brune *fnorm = fAnorm; 219b3c6a99cSPeter Brune *ynorm = yAnorm; 2209566063dSJacob Faibussowitsch PetscCall(VecCopy(FA,F)); 2219566063dSJacob Faibussowitsch PetscCall(VecCopy(XA,X)); 22238774f0aSPeter Brune } 22338774f0aSPeter Brune PetscFunctionReturn(0); 22438774f0aSPeter Brune } 22538774f0aSPeter Brune 22623b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes,PetscInt l,PetscReal fMnorm, PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,PetscBool *selectRestart) 22738774f0aSPeter Brune { 22838774f0aSPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data; 22938774f0aSPeter Brune 23038774f0aSPeter Brune PetscFunctionBegin; 23138774f0aSPeter Brune *selectRestart = PETSC_FALSE; 23238774f0aSPeter Brune /* difference stagnation restart */ 23321687c63SPeter Brune if ((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm)) && l > 0) { 23438774f0aSPeter Brune if (ngmres->monitor) { 23563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"difference restart: %e > %e\n",(double)(ngmres->epsilonB*dnorm),(double)dminnorm)); 23638774f0aSPeter Brune } 23738774f0aSPeter Brune *selectRestart = PETSC_TRUE; 23838774f0aSPeter Brune } 23938774f0aSPeter Brune /* residual stagnation restart */ 24038774f0aSPeter Brune if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) { 24138774f0aSPeter Brune if (ngmres->monitor) { 24263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"residual restart: %e > %e\n",(double)PetscSqrtReal(fAnorm),(double)(ngmres->gammaC*PetscSqrtReal(fminnorm)))); 24338774f0aSPeter Brune } 24438774f0aSPeter Brune *selectRestart = PETSC_TRUE; 24538774f0aSPeter Brune } 24623b3e82cSAsbjørn Nilsen Riseth 24723b3e82cSAsbjørn Nilsen Riseth /* F_M stagnation restart */ 24823b3e82cSAsbjørn Nilsen Riseth if (ngmres->restart_fm_rise && fMnorm > snes->norm) { 24923b3e82cSAsbjørn Nilsen Riseth if (ngmres->monitor) { 25063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"F_M rise restart: %e > %e\n",(double)fMnorm,(double)snes->norm)); 25123b3e82cSAsbjørn Nilsen Riseth } 25223b3e82cSAsbjørn Nilsen Riseth *selectRestart = PETSC_TRUE; 25323b3e82cSAsbjørn Nilsen Riseth } 25423b3e82cSAsbjørn Nilsen Riseth 25538774f0aSPeter Brune PetscFunctionReturn(0); 25638774f0aSPeter Brune } 257