xref: /petsc/src/snes/impls/ngmres/ngmresfunc.c (revision 5ef867c24bdfdb15299ef27085fe449976d0dbcc)
138774f0aSPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
238774f0aSPeter Brune #include <petscblaslapack.h>
338774f0aSPeter Brune 
4d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESUpdateSubspace_Private(SNES snes, PetscInt ivec, PetscInt l, Vec F, PetscReal fnorm, Vec X)
5d71ae5a4SJacob Faibussowitsch {
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;
163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1738774f0aSPeter Brune }
1838774f0aSPeter Brune 
19d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES snes, PetscInt ivec, PetscInt l, Vec XM, Vec FM, PetscReal fMnorm, Vec X, Vec XA, Vec FA)
20d71ae5a4SJacob Faibussowitsch {
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;
29*5ef867c2SRené Chenard   Vec          Y = snes->vec_sol_update;
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++) {
53ad540459SPierre Jolivet     for (i = 0; i < l; i++) H(i, j) = Q(i, j) - xi[i] - xi[j] + nu;
5438774f0aSPeter Brune   }
5538774f0aSPeter Brune   if (l == 1) {
5638774f0aSPeter Brune     /* simply set alpha[0] = beta[0] / H[0, 0] */
571aa26658SKarl Rupp     if (H(0, 0) != 0.) beta[0] = beta[0] / H(0, 0);
581aa26658SKarl Rupp     else beta[0] = 0.;
5938774f0aSPeter Brune   } else {
609566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(l, &ngmres->m));
619566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(l, &ngmres->n));
6238774f0aSPeter Brune     ngmres->info  = 0;
6338774f0aSPeter Brune     ngmres->rcond = -1.;
649566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
6538774f0aSPeter Brune #if defined(PETSC_USE_COMPLEX)
66792fecdfSBarry 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));
6738774f0aSPeter Brune #else
68792fecdfSBarry 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));
6938774f0aSPeter Brune #endif
709566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
7108401ef6SPierre Jolivet     PetscCheck(ngmres->info >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "Bad argument to GELSS");
7208401ef6SPierre Jolivet     PetscCheck(ngmres->info <= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD failed to converge");
7338774f0aSPeter Brune   }
74ad540459SPierre Jolivet   for (i = 0; i < l; i++) PetscCheck(!PetscIsInfOrNanScalar(beta[i]), PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD generated inconsistent output");
7538774f0aSPeter Brune   alph_total = 0.;
761aa26658SKarl Rupp   for (i = 0; i < l; i++) alph_total += beta[i];
771aa26658SKarl Rupp 
789566063dSJacob Faibussowitsch   PetscCall(VecCopy(XM, XA));
799566063dSJacob Faibussowitsch   PetscCall(VecScale(XA, 1. - alph_total));
809566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(XA, l, beta, Xdot));
8138774f0aSPeter Brune   /* check the validity of the step */
829566063dSJacob Faibussowitsch   PetscCall(VecCopy(XA, Y));
839566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Y, -1.0, X));
849566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPostCheck(snes->linesearch, X, Y, XA, &changed_y, &changed_w));
8546159c86SPeter Brune   if (!ngmres->approxfunc) {
86efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_LEFT) {
879566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, XA, NULL, FA));
8846159c86SPeter Brune     } else {
899566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, XA, FA));
9046159c86SPeter Brune     }
9135f895b4SBarry Smith   } else {
929566063dSJacob Faibussowitsch     PetscCall(VecCopy(FM, FA));
939566063dSJacob Faibussowitsch     PetscCall(VecScale(FA, 1. - alph_total));
949566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(FA, l, beta, Fdot));
95077c4231SPeter Brune   }
963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9738774f0aSPeter Brune }
9838774f0aSPeter Brune 
99d71ae5a4SJacob Faibussowitsch 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)
100d71ae5a4SJacob Faibussowitsch {
10138774f0aSPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
102b3c6a99cSPeter Brune   PetscReal    dcurnorm, dmin = -1.0;
10338774f0aSPeter Brune   Vec         *Xdot = ngmres->Xdot;
10438774f0aSPeter Brune   PetscInt     i;
10538774f0aSPeter Brune 
10638774f0aSPeter Brune   PetscFunctionBegin;
1071baa6e33SBarry Smith   if (xMnorm) PetscCall(VecNormBegin(XM, NORM_2, xMnorm));
1081baa6e33SBarry Smith   if (fMnorm) PetscCall(VecNormBegin(FM, NORM_2, fMnorm));
109b3c6a99cSPeter Brune   if (yMnorm) {
1109566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, D));
1119566063dSJacob Faibussowitsch     PetscCall(VecAXPY(D, -1.0, XM));
1129566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(D, NORM_2, yMnorm));
113b3c6a99cSPeter Brune   }
1141baa6e33SBarry Smith   if (xAnorm) PetscCall(VecNormBegin(XA, NORM_2, xAnorm));
1151baa6e33SBarry Smith   if (fAnorm) PetscCall(VecNormBegin(FA, NORM_2, fAnorm));
116b3c6a99cSPeter Brune   if (yAnorm) {
1179566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, D));
1189566063dSJacob Faibussowitsch     PetscCall(VecAXPY(D, -1.0, XA));
1199566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(D, NORM_2, yAnorm));
120b3c6a99cSPeter Brune   }
12138774f0aSPeter Brune   if (dnorm) {
1229566063dSJacob Faibussowitsch     PetscCall(VecCopy(XA, D));
1239566063dSJacob Faibussowitsch     PetscCall(VecAXPY(D, -1.0, XM));
1249566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(D, NORM_2, dnorm));
12538774f0aSPeter Brune   }
12638774f0aSPeter Brune   if (dminnorm) {
12738774f0aSPeter Brune     for (i = 0; i < l; i++) {
1289566063dSJacob Faibussowitsch       PetscCall(VecCopy(Xdot[i], D));
1299566063dSJacob Faibussowitsch       PetscCall(VecAXPY(D, -1.0, XA));
1309566063dSJacob Faibussowitsch       PetscCall(VecNormBegin(D, NORM_2, &ngmres->xnorms[i]));
13138774f0aSPeter Brune     }
13238774f0aSPeter Brune   }
1339566063dSJacob Faibussowitsch   if (xMnorm) PetscCall(VecNormEnd(XM, NORM_2, xMnorm));
1349566063dSJacob Faibussowitsch   if (fMnorm) PetscCall(VecNormEnd(FM, NORM_2, fMnorm));
1359566063dSJacob Faibussowitsch   if (yMnorm) PetscCall(VecNormEnd(D, NORM_2, yMnorm));
1369566063dSJacob Faibussowitsch   if (xAnorm) PetscCall(VecNormEnd(XA, NORM_2, xAnorm));
1379566063dSJacob Faibussowitsch   if (fAnorm) PetscCall(VecNormEnd(FA, NORM_2, fAnorm));
1389566063dSJacob Faibussowitsch   if (yAnorm) PetscCall(VecNormEnd(D, NORM_2, yAnorm));
1399566063dSJacob Faibussowitsch   if (dnorm) PetscCall(VecNormEnd(D, NORM_2, dnorm));
14038774f0aSPeter Brune   if (dminnorm) {
14138774f0aSPeter Brune     for (i = 0; i < l; i++) {
1429566063dSJacob Faibussowitsch       PetscCall(VecNormEnd(D, NORM_2, &ngmres->xnorms[i]));
14338774f0aSPeter Brune       dcurnorm = ngmres->xnorms[i];
144b3c6a99cSPeter Brune       if ((dcurnorm < dmin) || (dmin < 0.0)) dmin = dcurnorm;
14538774f0aSPeter Brune     }
146b3c6a99cSPeter Brune     *dminnorm = dmin;
14738774f0aSPeter Brune   }
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14938774f0aSPeter Brune }
15038774f0aSPeter Brune 
151d71ae5a4SJacob Faibussowitsch 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)
152d71ae5a4SJacob Faibussowitsch {
15338774f0aSPeter Brune   SNES_NGMRES         *ngmres = (SNES_NGMRES *)snes->data;
154422a814eSBarry Smith   SNESLineSearchReason lssucceed;
155422a814eSBarry Smith   PetscBool            selectA;
15638774f0aSPeter Brune 
15738774f0aSPeter Brune   PetscFunctionBegin;
15838774f0aSPeter Brune   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
15938774f0aSPeter Brune     /* X = X + \lambda(XA - X) */
16048a46eb9SPierre Jolivet     if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "||F_A||_2 = %e, ||F_M||_2 = %e\n", (double)fAnorm, (double)fMnorm));
1619566063dSJacob Faibussowitsch     PetscCall(VecCopy(FM, F));
1629566063dSJacob Faibussowitsch     PetscCall(VecCopy(XM, X));
1639566063dSJacob Faibussowitsch     PetscCall(VecCopy(XA, Y));
1649566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Y, -1.0, X));
16538774f0aSPeter Brune     *fnorm = fMnorm;
1669566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(ngmres->additive_linesearch, X, F, fnorm, Y));
1679566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(ngmres->additive_linesearch, &lssucceed));
1689566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(ngmres->additive_linesearch, xnorm, fnorm, ynorm));
169422a814eSBarry Smith     if (lssucceed) {
17038774f0aSPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
17138774f0aSPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
1723ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
17338774f0aSPeter Brune       }
17438774f0aSPeter Brune     }
17548a46eb9SPierre Jolivet     if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: ||F||_2 = %e\n", (double)*fnorm));
17638774f0aSPeter Brune   } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
17738774f0aSPeter Brune     selectA = PETSC_TRUE;
17838774f0aSPeter Brune     /* Conditions for choosing the accelerated answer */
17938774f0aSPeter Brune     /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
1801aa26658SKarl Rupp     if (fAnorm >= ngmres->gammaA * fminnorm) selectA = PETSC_FALSE;
1811aa26658SKarl Rupp 
18238774f0aSPeter Brune     /* Criterion B -- the choice of x^A isn't too close to some other choice */
18338774f0aSPeter Brune     if (ngmres->epsilonB * dnorm < dminnorm || PetscSqrtReal(*fnorm) < ngmres->deltaB * PetscSqrtReal(fminnorm)) {
1841aa26658SKarl Rupp     } else selectA = PETSC_FALSE;
1851aa26658SKarl Rupp 
18638774f0aSPeter Brune     if (selectA) {
18748a46eb9SPierre Jolivet       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n", (double)fAnorm, (double)fMnorm));
18838774f0aSPeter Brune       /* copy it over */
189b3c6a99cSPeter Brune       *xnorm = xAnorm;
19038774f0aSPeter Brune       *fnorm = fAnorm;
191b3c6a99cSPeter Brune       *ynorm = yAnorm;
1929566063dSJacob Faibussowitsch       PetscCall(VecCopy(FA, F));
1939566063dSJacob Faibussowitsch       PetscCall(VecCopy(XA, X));
19438774f0aSPeter Brune     } else {
19548a46eb9SPierre Jolivet       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n", (double)fAnorm, (double)fMnorm));
196b3c6a99cSPeter Brune       *xnorm = xMnorm;
19738774f0aSPeter Brune       *fnorm = fMnorm;
198b3c6a99cSPeter Brune       *ynorm = yMnorm;
1999566063dSJacob Faibussowitsch       PetscCall(VecCopy(XM, Y));
2009566063dSJacob Faibussowitsch       PetscCall(VecAXPY(Y, -1.0, X));
2019566063dSJacob Faibussowitsch       PetscCall(VecCopy(FM, F));
2029566063dSJacob Faibussowitsch       PetscCall(VecCopy(XM, X));
20338774f0aSPeter Brune     }
20438774f0aSPeter Brune   } else { /* none */
205b3c6a99cSPeter Brune     *xnorm = xAnorm;
20638774f0aSPeter Brune     *fnorm = fAnorm;
207b3c6a99cSPeter Brune     *ynorm = yAnorm;
2089566063dSJacob Faibussowitsch     PetscCall(VecCopy(FA, F));
2099566063dSJacob Faibussowitsch     PetscCall(VecCopy(XA, X));
21038774f0aSPeter Brune   }
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21238774f0aSPeter Brune }
21338774f0aSPeter Brune 
214d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes, PetscInt l, PetscReal fMnorm, PetscReal fAnorm, PetscReal dnorm, PetscReal fminnorm, PetscReal dminnorm, PetscBool *selectRestart)
215d71ae5a4SJacob Faibussowitsch {
21638774f0aSPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
21738774f0aSPeter Brune 
21838774f0aSPeter Brune   PetscFunctionBegin;
21938774f0aSPeter Brune   *selectRestart = PETSC_FALSE;
22038774f0aSPeter Brune   /* difference stagnation restart */
22121687c63SPeter Brune   if ((ngmres->epsilonB * dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB * PetscSqrtReal(fminnorm)) && l > 0) {
22248a46eb9SPierre Jolivet     if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", (double)(ngmres->epsilonB * dnorm), (double)dminnorm));
22338774f0aSPeter Brune     *selectRestart = PETSC_TRUE;
22438774f0aSPeter Brune   }
22538774f0aSPeter Brune   /* residual stagnation restart */
22638774f0aSPeter Brune   if (PetscSqrtReal(fAnorm) > ngmres->gammaC * PetscSqrtReal(fminnorm)) {
22748a46eb9SPierre Jolivet     if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", (double)PetscSqrtReal(fAnorm), (double)(ngmres->gammaC * PetscSqrtReal(fminnorm))));
22838774f0aSPeter Brune     *selectRestart = PETSC_TRUE;
22938774f0aSPeter Brune   }
23023b3e82cSAsbjørn Nilsen Riseth 
23123b3e82cSAsbjørn Nilsen Riseth   /* F_M stagnation restart */
23223b3e82cSAsbjørn Nilsen Riseth   if (ngmres->restart_fm_rise && fMnorm > snes->norm) {
23348a46eb9SPierre Jolivet     if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "F_M rise restart: %e > %e\n", (double)fMnorm, (double)snes->norm));
23423b3e82cSAsbjørn Nilsen Riseth     *selectRestart = PETSC_TRUE;
23523b3e82cSAsbjørn Nilsen Riseth   }
23623b3e82cSAsbjørn Nilsen Riseth 
2373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23838774f0aSPeter Brune }
239