xref: /petsc/src/snes/impls/ngmres/ngmresfunc.c (revision 647c55fd8ae18e0846d428867ca031dc04cbd941)
138774f0aSPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
238774f0aSPeter Brune #include <petscblaslapack.h>
338774f0aSPeter Brune 
44936d809SStefano Zampini PetscErrorCode SNESNGMRESGetAdditiveLineSearch_Private(SNES snes, SNESLineSearch *linesearch)
54936d809SStefano Zampini {
64936d809SStefano Zampini   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
74936d809SStefano Zampini 
84936d809SStefano Zampini   PetscFunctionBegin;
94936d809SStefano Zampini   if (!ngmres->additive_linesearch) {
104936d809SStefano Zampini     const char *optionsprefix;
114936d809SStefano Zampini     PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
124936d809SStefano Zampini     PetscCall(SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &ngmres->additive_linesearch));
134936d809SStefano Zampini     PetscCall(SNESLineSearchSetSNES(ngmres->additive_linesearch, snes));
144936d809SStefano Zampini     PetscCall(SNESLineSearchSetType(ngmres->additive_linesearch, SNESLINESEARCHL2));
154936d809SStefano Zampini     PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, "snes_ngmres_additive_"));
164936d809SStefano Zampini     PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, optionsprefix));
174936d809SStefano Zampini     PetscCall(PetscObjectIncrementTabLevel((PetscObject)ngmres->additive_linesearch, (PetscObject)snes, 1));
184936d809SStefano Zampini   }
194936d809SStefano Zampini   *linesearch = ngmres->additive_linesearch;
204936d809SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
214936d809SStefano Zampini }
224936d809SStefano Zampini 
23d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESUpdateSubspace_Private(SNES snes, PetscInt ivec, PetscInt l, Vec F, PetscReal fnorm, Vec X)
24d71ae5a4SJacob Faibussowitsch {
2538774f0aSPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
2638774f0aSPeter Brune   Vec         *Fdot   = ngmres->Fdot;
2738774f0aSPeter Brune   Vec         *Xdot   = ngmres->Xdot;
2838774f0aSPeter Brune 
2938774f0aSPeter Brune   PetscFunctionBegin;
3063a3b9bcSJacob Faibussowitsch   PetscCheck(ivec <= l, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Cannot update vector %" PetscInt_FMT " with space size %" PetscInt_FMT "!", ivec, l);
319566063dSJacob Faibussowitsch   PetscCall(VecCopy(F, Fdot[ivec]));
329566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, Xdot[ivec]));
331aa26658SKarl Rupp 
3438774f0aSPeter Brune   ngmres->fnorms[ivec] = fnorm;
353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3638774f0aSPeter Brune }
3738774f0aSPeter Brune 
38d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES snes, PetscInt ivec, PetscInt l, Vec XM, Vec FM, PetscReal fMnorm, Vec X, Vec XA, Vec FA)
39d71ae5a4SJacob Faibussowitsch {
4038774f0aSPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
4138774f0aSPeter Brune   PetscInt     i, j;
4238774f0aSPeter Brune   Vec         *Fdot       = ngmres->Fdot;
4338774f0aSPeter Brune   Vec         *Xdot       = ngmres->Xdot;
4438774f0aSPeter Brune   PetscScalar *beta       = ngmres->beta;
4538774f0aSPeter Brune   PetscScalar *xi         = ngmres->xi;
4638774f0aSPeter Brune   PetscScalar  alph_total = 0.;
4738774f0aSPeter Brune   PetscReal    nu;
485ef867c2SRené Chenard   Vec          Y = snes->vec_sol_update;
4938774f0aSPeter Brune   PetscBool    changed_y, changed_w;
5038774f0aSPeter Brune 
5138774f0aSPeter Brune   PetscFunctionBegin;
5238774f0aSPeter Brune   nu = fMnorm * fMnorm;
5338774f0aSPeter Brune 
54dd8e379bSPierre Jolivet   /* construct the right-hand side and xi factors */
55b3c6a99cSPeter Brune   if (l > 0) {
569566063dSJacob Faibussowitsch     PetscCall(VecMDotBegin(FM, l, Fdot, xi));
579566063dSJacob Faibussowitsch     PetscCall(VecMDotBegin(Fdot[ivec], l, Fdot, beta));
589566063dSJacob Faibussowitsch     PetscCall(VecMDotEnd(FM, l, Fdot, xi));
599566063dSJacob Faibussowitsch     PetscCall(VecMDotEnd(Fdot[ivec], l, Fdot, beta));
60b3c6a99cSPeter Brune     for (i = 0; i < l; i++) {
61b3c6a99cSPeter Brune       Q(i, ivec) = beta[i];
62b3c6a99cSPeter Brune       Q(ivec, i) = beta[i];
63b3c6a99cSPeter Brune     }
64b3c6a99cSPeter Brune   } else {
65b3c6a99cSPeter Brune     Q(0, 0) = ngmres->fnorms[ivec] * ngmres->fnorms[ivec];
66b3c6a99cSPeter Brune   }
67b3c6a99cSPeter Brune 
681aa26658SKarl Rupp   for (i = 0; i < l; i++) beta[i] = nu - xi[i];
691aa26658SKarl Rupp 
7038774f0aSPeter Brune   /* construct h */
7138774f0aSPeter Brune   for (j = 0; j < l; j++) {
72ad540459SPierre Jolivet     for (i = 0; i < l; i++) H(i, j) = Q(i, j) - xi[i] - xi[j] + nu;
7338774f0aSPeter Brune   }
7438774f0aSPeter Brune   if (l == 1) {
7538774f0aSPeter Brune     /* simply set alpha[0] = beta[0] / H[0, 0] */
761aa26658SKarl Rupp     if (H(0, 0) != 0.) beta[0] = beta[0] / H(0, 0);
771aa26658SKarl Rupp     else beta[0] = 0.;
7838774f0aSPeter Brune   } else {
799566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(l, &ngmres->m));
809566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(l, &ngmres->n));
8138774f0aSPeter Brune     ngmres->info  = 0;
8238774f0aSPeter Brune     ngmres->rcond = -1.;
839566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
8438774f0aSPeter Brune #if defined(PETSC_USE_COMPLEX)
85792fecdfSBarry 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));
8638774f0aSPeter Brune #else
87792fecdfSBarry 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));
8838774f0aSPeter Brune #endif
899566063dSJacob Faibussowitsch     PetscCall(PetscFPTrapPop());
9008401ef6SPierre Jolivet     PetscCheck(ngmres->info >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "Bad argument to GELSS");
9108401ef6SPierre Jolivet     PetscCheck(ngmres->info <= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD failed to converge");
9238774f0aSPeter Brune   }
93ad540459SPierre Jolivet   for (i = 0; i < l; i++) PetscCheck(!PetscIsInfOrNanScalar(beta[i]), PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD generated inconsistent output");
9438774f0aSPeter Brune   alph_total = 0.;
951aa26658SKarl Rupp   for (i = 0; i < l; i++) alph_total += beta[i];
961aa26658SKarl Rupp 
97*647c55fdSStefano Zampini   PetscCall(VecAXPBY(XA, 1.0 - alph_total, 0.0, XM));
989566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(XA, l, beta, Xdot));
9938774f0aSPeter Brune   /* check the validity of the step */
100*647c55fdSStefano Zampini   PetscCall(VecWAXPY(Y, -1.0, X, XA));
1019566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPostCheck(snes->linesearch, X, Y, XA, &changed_y, &changed_w));
10246159c86SPeter Brune   if (!ngmres->approxfunc) {
103efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_LEFT) {
1049566063dSJacob Faibussowitsch       PetscCall(SNESApplyNPC(snes, XA, NULL, FA));
10546159c86SPeter Brune     } else {
1069566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, XA, FA));
10746159c86SPeter Brune     }
10835f895b4SBarry Smith   } else {
109*647c55fdSStefano Zampini     PetscCall(VecAXPBY(FA, 1.0 - alph_total, 0.0, FM));
1109566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(FA, l, beta, Fdot));
111077c4231SPeter Brune   }
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11338774f0aSPeter Brune }
11438774f0aSPeter Brune 
115d71ae5a4SJacob 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)
116d71ae5a4SJacob Faibussowitsch {
11738774f0aSPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
118b3c6a99cSPeter Brune   PetscReal    dcurnorm, dmin = -1.0;
11938774f0aSPeter Brune   Vec         *Xdot = ngmres->Xdot;
12038774f0aSPeter Brune   PetscInt     i;
12138774f0aSPeter Brune 
12238774f0aSPeter Brune   PetscFunctionBegin;
1231baa6e33SBarry Smith   if (xMnorm) PetscCall(VecNormBegin(XM, NORM_2, xMnorm));
1241baa6e33SBarry Smith   if (fMnorm) PetscCall(VecNormBegin(FM, NORM_2, fMnorm));
125b3c6a99cSPeter Brune   if (yMnorm) {
126*647c55fdSStefano Zampini     PetscCall(VecWAXPY(D, -1.0, XM, X));
1279566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(D, NORM_2, yMnorm));
128b3c6a99cSPeter Brune   }
1291baa6e33SBarry Smith   if (xAnorm) PetscCall(VecNormBegin(XA, NORM_2, xAnorm));
1301baa6e33SBarry Smith   if (fAnorm) PetscCall(VecNormBegin(FA, NORM_2, fAnorm));
131b3c6a99cSPeter Brune   if (yAnorm) {
132*647c55fdSStefano Zampini     PetscCall(VecWAXPY(D, -1.0, XA, X));
1339566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(D, NORM_2, yAnorm));
134b3c6a99cSPeter Brune   }
13538774f0aSPeter Brune   if (dnorm) {
136*647c55fdSStefano Zampini     PetscCall(VecWAXPY(D, -1.0, XM, XA));
1379566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(D, NORM_2, dnorm));
13838774f0aSPeter Brune   }
13938774f0aSPeter Brune   if (dminnorm) {
14038774f0aSPeter Brune     for (i = 0; i < l; i++) {
141*647c55fdSStefano Zampini       PetscCall(VecWAXPY(D, -1.0, XA, Xdot[i]));
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   }
1603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16138774f0aSPeter Brune }
16238774f0aSPeter Brune 
1634936d809SStefano Zampini PetscErrorCode SNESNGMRESSelect_Private(SNES snes, PetscInt k_restart, Vec XM, Vec FM, PetscReal xMnorm, PetscReal fMnorm, PetscReal yMnorm, PetscReal objM, Vec XA, Vec FA, PetscReal xAnorm, PetscReal fAnorm, PetscReal yAnorm, PetscReal objA, PetscReal dnorm, PetscReal objmin, PetscReal dminnorm, Vec X, Vec F, Vec Y, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm)
164d71ae5a4SJacob Faibussowitsch {
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) */
1724936d809SStefano Zampini     if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "obj(X_A) = %e, ||F_A||_2 = %e, obj(X_M) = %e, ||F_M||_2 = %e\n", (double)objA, (double)fAnorm, (double)objM, (double)fMnorm));
1734936d809SStefano Zampini     /* Test if is XA - XM is a descent direction: we want < F(XM), XA - XM > not positive
1744936d809SStefano Zampini        If positive, GMRES will be restarted see https://epubs.siam.org/doi/pdf/10.1137/110835530 */
1759566063dSJacob Faibussowitsch     PetscCall(VecCopy(FM, F));
1769566063dSJacob Faibussowitsch     PetscCall(VecCopy(XM, X));
1774936d809SStefano Zampini     PetscCall(VecWAXPY(Y, -1.0, XA, X));                        /* minus sign since linesearch expects to find Xnew = X - lambda * Y */
1784936d809SStefano Zampini     PetscCall(VecDotRealPart(FM, Y, &ngmres->descent_ls_test)); /* this is actually < F(XM), XM - XA > */
17938774f0aSPeter Brune     *fnorm = fMnorm;
1804936d809SStefano Zampini     if (ngmres->descent_ls_test < 0) { /* XA - XM is not a descent direction, select XM */
1814936d809SStefano Zampini       *xnorm = xMnorm;
1824936d809SStefano Zampini       *fnorm = fMnorm;
1834936d809SStefano Zampini       *ynorm = yMnorm;
1844936d809SStefano Zampini       PetscCall(VecWAXPY(Y, -1.0, X, XM));
1854936d809SStefano Zampini       PetscCall(VecCopy(FM, F));
1864936d809SStefano Zampini       PetscCall(VecCopy(XM, X));
1874936d809SStefano Zampini     } else {
1884936d809SStefano Zampini       PetscCall(SNESNGMRESGetAdditiveLineSearch_Private(snes, &ngmres->additive_linesearch));
1899566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchApply(ngmres->additive_linesearch, X, F, fnorm, Y));
1909566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetReason(ngmres->additive_linesearch, &lssucceed));
1919566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetNorms(ngmres->additive_linesearch, xnorm, fnorm, ynorm));
192422a814eSBarry Smith       if (lssucceed) {
19338774f0aSPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
19438774f0aSPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
1953ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
19638774f0aSPeter Brune         }
19738774f0aSPeter Brune       }
1984936d809SStefano Zampini     }
1994936d809SStefano Zampini     if (ngmres->monitor) {
2004936d809SStefano Zampini       PetscReal objT = *fnorm;
2014936d809SStefano Zampini       PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
2021aa26658SKarl Rupp 
2034936d809SStefano Zampini       PetscCall(SNESGetObjective(snes, &objective, NULL));
2044936d809SStefano Zampini       if (objective) PetscCall(SNESComputeObjective(snes, X, &objT));
2054936d809SStefano Zampini       PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Additive solution: objective = %e\n", (double)objT));
2064936d809SStefano Zampini     }
2074936d809SStefano Zampini   } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
2084936d809SStefano Zampini     /* Conditions for choosing the accelerated answer:
2094936d809SStefano Zampini           Criterion A -- the objective function isn't increased above the minimum by too much
2104936d809SStefano Zampini           Criterion B -- the choice of x^A isn't too close to some other choice
2114936d809SStefano Zampini     */
2124936d809SStefano Zampini     selectA = (PetscBool)(/* A */ (objA < ngmres->gammaA * objmin) && /* B */ (ngmres->epsilonB * dnorm < dminnorm || objA < ngmres->deltaB * objmin));
2131aa26658SKarl Rupp 
21438774f0aSPeter Brune     if (selectA) {
2154936d809SStefano Zampini       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "picked X_A, obj(X_A) = %e, ||F_A||_2 = %e, obj(X_M) = %e, ||F_M||_2 = %e\n", (double)objA, (double)fAnorm, (double)objM, (double)fMnorm));
21638774f0aSPeter Brune       /* copy it over */
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     } else {
2234936d809SStefano Zampini       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "picked X_M, obj(X_A) = %e, ||F_A||_2 = %e, obj(X_M) = %e, ||F_M||_2 = %e\n", (double)objA, (double)fAnorm, (double)objM, (double)fMnorm));
224b3c6a99cSPeter Brune       *xnorm = xMnorm;
22538774f0aSPeter Brune       *fnorm = fMnorm;
226b3c6a99cSPeter Brune       *ynorm = yMnorm;
2274936d809SStefano Zampini       PetscCall(VecWAXPY(Y, -1.0, X, XM));
2289566063dSJacob Faibussowitsch       PetscCall(VecCopy(FM, F));
2299566063dSJacob Faibussowitsch       PetscCall(VecCopy(XM, X));
23038774f0aSPeter Brune     }
23138774f0aSPeter Brune   } else { /* none */
232b3c6a99cSPeter Brune     *xnorm = xAnorm;
23338774f0aSPeter Brune     *fnorm = fAnorm;
234b3c6a99cSPeter Brune     *ynorm = yAnorm;
2359566063dSJacob Faibussowitsch     PetscCall(VecCopy(FA, F));
2369566063dSJacob Faibussowitsch     PetscCall(VecCopy(XA, X));
23738774f0aSPeter Brune   }
2383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23938774f0aSPeter Brune }
24038774f0aSPeter Brune 
2414936d809SStefano Zampini PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes, PetscInt l, PetscReal obj, PetscReal objM, PetscReal objA, PetscReal dnorm, PetscReal objmin, PetscReal dminnorm, PetscBool *selectRestart)
242d71ae5a4SJacob Faibussowitsch {
24338774f0aSPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
24438774f0aSPeter Brune 
24538774f0aSPeter Brune   PetscFunctionBegin;
24638774f0aSPeter Brune   *selectRestart = PETSC_FALSE;
2474936d809SStefano Zampini   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
2484936d809SStefano Zampini     if (ngmres->descent_ls_test < 0) { /* XA - XM is not a descent direction */
2494936d809SStefano Zampini       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "ascent restart: %e > 0\n", (double)-ngmres->descent_ls_test));
2504936d809SStefano Zampini       *selectRestart = PETSC_TRUE;
2514936d809SStefano Zampini     }
2524936d809SStefano Zampini   } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
25338774f0aSPeter Brune     /* difference stagnation restart */
2544936d809SStefano Zampini     if (ngmres->epsilonB * dnorm > dminnorm && objA > ngmres->deltaB * objmin && l > 0) {
25548a46eb9SPierre Jolivet       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "difference restart: %e > %e\n", (double)(ngmres->epsilonB * dnorm), (double)dminnorm));
25638774f0aSPeter Brune       *selectRestart = PETSC_TRUE;
25738774f0aSPeter Brune     }
25838774f0aSPeter Brune     /* residual stagnation restart */
2594936d809SStefano Zampini     if (objA > ngmres->gammaC * objmin) {
2604936d809SStefano Zampini       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "residual restart: %e > %e\n", (double)objA, (double)(ngmres->gammaC * objmin)));
26138774f0aSPeter Brune       *selectRestart = PETSC_TRUE;
26238774f0aSPeter Brune     }
26323b3e82cSAsbjørn Nilsen Riseth 
26423b3e82cSAsbjørn Nilsen Riseth     /* F_M stagnation restart */
2654936d809SStefano Zampini     if (ngmres->restart_fm_rise && objM > obj) {
2664936d809SStefano Zampini       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "F_M rise restart: %e > %e\n", (double)objM, (double)obj));
26723b3e82cSAsbjørn Nilsen Riseth       *selectRestart = PETSC_TRUE;
26823b3e82cSAsbjørn Nilsen Riseth     }
2694936d809SStefano Zampini   }
2703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27138774f0aSPeter Brune }
272