xref: /petsc/src/snes/impls/ngmres/ngmresfunc.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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;
11*08401ef6SPierre Jolivet   PetscCheck(ivec <= l,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot update vector %D with space size %D!",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)
688b83055fSJed 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));
6938774f0aSPeter Brune #else
708b83055fSJed 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));
7138774f0aSPeter Brune #endif
729566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
73*08401ef6SPierre Jolivet     PetscCheck(ngmres->info >= 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
74*08401ef6SPierre 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++) {
772c71b3e2SJacob Faibussowitsch     PetscCheckFalse(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;
111b3c6a99cSPeter Brune   if (xMnorm) {
1129566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(XM,NORM_2,xMnorm));
113b3c6a99cSPeter Brune   }
114b3c6a99cSPeter Brune   if (fMnorm) {
1159566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(FM,NORM_2,fMnorm));
116b3c6a99cSPeter Brune   }
117b3c6a99cSPeter Brune   if (yMnorm) {
1189566063dSJacob Faibussowitsch     PetscCall(VecCopy(X,D));
1199566063dSJacob Faibussowitsch     PetscCall(VecAXPY(D,-1.0,XM));
1209566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(D,NORM_2,yMnorm));
121b3c6a99cSPeter Brune   }
122b3c6a99cSPeter Brune   if (xAnorm) {
1239566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(XA,NORM_2,xAnorm));
124b3c6a99cSPeter Brune   }
12538774f0aSPeter Brune   if (fAnorm) {
1269566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(FA,NORM_2,fAnorm));
12738774f0aSPeter Brune   }
128b3c6a99cSPeter Brune   if (yAnorm) {
1299566063dSJacob Faibussowitsch     PetscCall(VecCopy(X,D));
1309566063dSJacob Faibussowitsch     PetscCall(VecAXPY(D,-1.0,XA));
1319566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(D,NORM_2,yAnorm));
132b3c6a99cSPeter Brune   }
13338774f0aSPeter Brune   if (dnorm) {
1349566063dSJacob Faibussowitsch     PetscCall(VecCopy(XA,D));
1359566063dSJacob Faibussowitsch     PetscCall(VecAXPY(D,-1.0,XM));
1369566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(D,NORM_2,dnorm));
13738774f0aSPeter Brune   }
13838774f0aSPeter Brune   if (dminnorm) {
13938774f0aSPeter Brune     for (i=0; i<l; i++) {
1409566063dSJacob Faibussowitsch       PetscCall(VecCopy(Xdot[i],D));
1419566063dSJacob Faibussowitsch       PetscCall(VecAXPY(D,-1.0,XA));
1429566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(D,NORM_2,&ngmres->xnorms[i]));
14338774f0aSPeter Brune     }
14438774f0aSPeter Brune   }
1459566063dSJacob Faibussowitsch   if (xMnorm) PetscCall(VecNormEnd(XM,NORM_2,xMnorm));
1469566063dSJacob Faibussowitsch   if (fMnorm) PetscCall(VecNormEnd(FM,NORM_2,fMnorm));
1479566063dSJacob Faibussowitsch   if (yMnorm) PetscCall(VecNormEnd(D,NORM_2,yMnorm));
1489566063dSJacob Faibussowitsch   if (xAnorm) PetscCall(VecNormEnd(XA,NORM_2,xAnorm));
1499566063dSJacob Faibussowitsch   if (fAnorm) PetscCall(VecNormEnd(FA,NORM_2,fAnorm));
1509566063dSJacob Faibussowitsch   if (yAnorm) PetscCall(VecNormEnd(D,NORM_2,yAnorm));
1519566063dSJacob Faibussowitsch   if (dnorm) PetscCall(VecNormEnd(D,NORM_2,dnorm));
15238774f0aSPeter Brune   if (dminnorm) {
15338774f0aSPeter Brune     for (i=0; i<l; i++) {
1549566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(D,NORM_2,&ngmres->xnorms[i]));
15538774f0aSPeter Brune       dcurnorm = ngmres->xnorms[i];
156b3c6a99cSPeter Brune       if ((dcurnorm < dmin) || (dmin < 0.0)) dmin = dcurnorm;
15738774f0aSPeter Brune     }
158b3c6a99cSPeter Brune     *dminnorm = dmin;
15938774f0aSPeter Brune   }
16038774f0aSPeter Brune   PetscFunctionReturn(0);
16138774f0aSPeter Brune }
16238774f0aSPeter Brune 
163b3c6a99cSPeter 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)
16438774f0aSPeter Brune {
16538774f0aSPeter Brune   SNES_NGMRES          *ngmres = (SNES_NGMRES*) snes->data;
166422a814eSBarry Smith   SNESLineSearchReason lssucceed;
167422a814eSBarry Smith   PetscBool            selectA;
16838774f0aSPeter Brune 
16938774f0aSPeter Brune   PetscFunctionBegin;
17038774f0aSPeter Brune   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
17138774f0aSPeter Brune     /* X = X + \lambda(XA - X) */
17238774f0aSPeter Brune     if (ngmres->monitor) {
1739566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm));
17438774f0aSPeter Brune     }
1759566063dSJacob Faibussowitsch     PetscCall(VecCopy(FM,F));
1769566063dSJacob Faibussowitsch     PetscCall(VecCopy(XM,X));
1779566063dSJacob Faibussowitsch     PetscCall(VecCopy(XA,Y));
1789566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Y,-1.0,X));
17938774f0aSPeter Brune     *fnorm = fMnorm;
1809566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(ngmres->additive_linesearch,X,F,fnorm,Y));
1819566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(ngmres->additive_linesearch,&lssucceed));
1829566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(ngmres->additive_linesearch,xnorm,fnorm,ynorm));
183422a814eSBarry Smith     if (lssucceed) {
18438774f0aSPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
18538774f0aSPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
18638774f0aSPeter Brune         PetscFunctionReturn(0);
18738774f0aSPeter Brune       }
18838774f0aSPeter Brune     }
18938774f0aSPeter Brune     if (ngmres->monitor) {
1909566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"Additive solution: ||F||_2 = %e\n",*fnorm));
19138774f0aSPeter Brune     }
19238774f0aSPeter Brune   } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
19338774f0aSPeter Brune     selectA = PETSC_TRUE;
19438774f0aSPeter Brune     /* Conditions for choosing the accelerated answer */
19538774f0aSPeter Brune     /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
1961aa26658SKarl Rupp     if (fAnorm >= ngmres->gammaA*fminnorm) selectA = PETSC_FALSE;
1971aa26658SKarl Rupp 
19838774f0aSPeter Brune     /* Criterion B -- the choice of x^A isn't too close to some other choice */
19938774f0aSPeter Brune     if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(*fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) {
2001aa26658SKarl Rupp     } else selectA=PETSC_FALSE;
2011aa26658SKarl Rupp 
20238774f0aSPeter Brune     if (selectA) {
20338774f0aSPeter Brune       if (ngmres->monitor) {
2049566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm));
20538774f0aSPeter Brune       }
20638774f0aSPeter Brune       /* copy it over */
207b3c6a99cSPeter Brune       *xnorm = xAnorm;
20838774f0aSPeter Brune       *fnorm = fAnorm;
209b3c6a99cSPeter Brune       *ynorm = yAnorm;
2109566063dSJacob Faibussowitsch       PetscCall(VecCopy(FA,F));
2119566063dSJacob Faibussowitsch       PetscCall(VecCopy(XA,X));
21238774f0aSPeter Brune     } else {
21338774f0aSPeter Brune       if (ngmres->monitor) {
2149566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm));
21538774f0aSPeter Brune       }
216b3c6a99cSPeter Brune       *xnorm = xMnorm;
21738774f0aSPeter Brune       *fnorm = fMnorm;
218b3c6a99cSPeter Brune       *ynorm = yMnorm;
2199566063dSJacob Faibussowitsch       PetscCall(VecCopy(XM,Y));
2209566063dSJacob Faibussowitsch       PetscCall(VecAXPY(Y,-1.0,X));
2219566063dSJacob Faibussowitsch       PetscCall(VecCopy(FM,F));
2229566063dSJacob Faibussowitsch       PetscCall(VecCopy(XM,X));
22338774f0aSPeter Brune     }
22438774f0aSPeter Brune   } else { /* none */
225b3c6a99cSPeter Brune     *xnorm = xAnorm;
22638774f0aSPeter Brune     *fnorm = fAnorm;
227b3c6a99cSPeter Brune     *ynorm = yAnorm;
2289566063dSJacob Faibussowitsch     PetscCall(VecCopy(FA,F));
2299566063dSJacob Faibussowitsch     PetscCall(VecCopy(XA,X));
23038774f0aSPeter Brune   }
23138774f0aSPeter Brune   PetscFunctionReturn(0);
23238774f0aSPeter Brune }
23338774f0aSPeter Brune 
23423b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes,PetscInt l,PetscReal fMnorm, PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,PetscBool *selectRestart)
23538774f0aSPeter Brune {
23638774f0aSPeter Brune   SNES_NGMRES    *ngmres = (SNES_NGMRES*)snes->data;
23738774f0aSPeter Brune 
23838774f0aSPeter Brune   PetscFunctionBegin;
23938774f0aSPeter Brune   *selectRestart = PETSC_FALSE;
24038774f0aSPeter Brune   /* difference stagnation restart */
24121687c63SPeter Brune   if ((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm)) && l > 0) {
24238774f0aSPeter Brune     if (ngmres->monitor) {
2439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"difference restart: %e > %e\n",ngmres->epsilonB*dnorm,dminnorm));
24438774f0aSPeter Brune     }
24538774f0aSPeter Brune     *selectRestart = PETSC_TRUE;
24638774f0aSPeter Brune   }
24738774f0aSPeter Brune   /* residual stagnation restart */
24838774f0aSPeter Brune   if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) {
24938774f0aSPeter Brune     if (ngmres->monitor) {
2509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"residual restart: %e > %e\n",PetscSqrtReal(fAnorm),ngmres->gammaC*PetscSqrtReal(fminnorm)));
25138774f0aSPeter Brune     }
25238774f0aSPeter Brune     *selectRestart = PETSC_TRUE;
25338774f0aSPeter Brune   }
25423b3e82cSAsbjørn Nilsen Riseth 
25523b3e82cSAsbjørn Nilsen Riseth   /* F_M stagnation restart */
25623b3e82cSAsbjørn Nilsen Riseth   if (ngmres->restart_fm_rise && fMnorm > snes->norm) {
25723b3e82cSAsbjørn Nilsen Riseth     if (ngmres->monitor) {
2589566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(ngmres->monitor,"F_M rise restart: %e > %e\n",fMnorm,snes->norm));
25923b3e82cSAsbjørn Nilsen Riseth     }
26023b3e82cSAsbjørn Nilsen Riseth     *selectRestart = PETSC_TRUE;
26123b3e82cSAsbjørn Nilsen Riseth   }
26223b3e82cSAsbjørn Nilsen Riseth 
26338774f0aSPeter Brune   PetscFunctionReturn(0);
26438774f0aSPeter Brune }
265