xref: /petsc/src/snes/impls/ngmres/ngmresfunc.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode ierr;
1038774f0aSPeter Brune 
1138774f0aSPeter Brune   PetscFunctionBegin;
12*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(ivec > l,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot update vector %D with space size %D!",ivec,l);
1338774f0aSPeter Brune   ierr = VecCopy(F,Fdot[ivec]);CHKERRQ(ierr);
1438774f0aSPeter Brune   ierr = VecCopy(X,Xdot[ivec]);CHKERRQ(ierr);
151aa26658SKarl Rupp 
1638774f0aSPeter Brune   ngmres->fnorms[ivec] = fnorm;
1738774f0aSPeter Brune   PetscFunctionReturn(0);
1838774f0aSPeter Brune }
1938774f0aSPeter Brune 
20b3c6a99cSPeter Brune PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES snes,PetscInt ivec,PetscInt l,Vec XM,Vec FM,PetscReal fMnorm,Vec X,Vec XA,Vec FA)
2138774f0aSPeter Brune {
2238774f0aSPeter Brune   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
2338774f0aSPeter Brune   PetscInt       i,j;
2438774f0aSPeter Brune   Vec            *Fdot      = ngmres->Fdot;
2538774f0aSPeter Brune   Vec            *Xdot      = ngmres->Xdot;
2638774f0aSPeter Brune   PetscScalar    *beta      = ngmres->beta;
2738774f0aSPeter Brune   PetscScalar    *xi        = ngmres->xi;
2838774f0aSPeter Brune   PetscScalar    alph_total = 0.;
2938774f0aSPeter Brune   PetscErrorCode ierr;
3038774f0aSPeter Brune   PetscReal      nu;
3138774f0aSPeter Brune   Vec            Y = snes->work[2];
3238774f0aSPeter Brune   PetscBool      changed_y,changed_w;
3338774f0aSPeter Brune 
3438774f0aSPeter Brune   PetscFunctionBegin;
3538774f0aSPeter Brune   nu = fMnorm*fMnorm;
3638774f0aSPeter Brune 
3738774f0aSPeter Brune   /* construct the right hand side and xi factors */
38b3c6a99cSPeter Brune   if (l > 0) {
39b3c6a99cSPeter Brune     ierr = VecMDotBegin(FM,l,Fdot,xi);CHKERRQ(ierr);
40b3c6a99cSPeter Brune     ierr = VecMDotBegin(Fdot[ivec],l,Fdot,beta);CHKERRQ(ierr);
41b3c6a99cSPeter Brune     ierr = VecMDotEnd(FM,l,Fdot,xi);CHKERRQ(ierr);
42b3c6a99cSPeter Brune     ierr = VecMDotEnd(Fdot[ivec],l,Fdot,beta);CHKERRQ(ierr);
43b3c6a99cSPeter Brune     for (i = 0; i < l; i++) {
44b3c6a99cSPeter Brune       Q(i,ivec) = beta[i];
45b3c6a99cSPeter Brune       Q(ivec,i) = beta[i];
46b3c6a99cSPeter Brune     }
47b3c6a99cSPeter Brune   } else {
48b3c6a99cSPeter Brune     Q(0,0) = ngmres->fnorms[ivec]*ngmres->fnorms[ivec];
49b3c6a99cSPeter Brune   }
50b3c6a99cSPeter Brune 
511aa26658SKarl Rupp   for (i = 0; i < l; i++) beta[i] = nu - xi[i];
521aa26658SKarl Rupp 
5338774f0aSPeter Brune   /* construct h */
5438774f0aSPeter Brune   for (j = 0; j < l; j++) {
5538774f0aSPeter Brune     for (i = 0; i < l; i++) {
5638774f0aSPeter Brune       H(i,j) = Q(i,j)-xi[i]-xi[j]+nu;
5738774f0aSPeter Brune     }
5838774f0aSPeter Brune   }
5938774f0aSPeter Brune   if (l == 1) {
6038774f0aSPeter Brune     /* simply set alpha[0] = beta[0] / H[0, 0] */
611aa26658SKarl Rupp     if (H(0,0) != 0.) beta[0] = beta[0]/H(0,0);
621aa26658SKarl Rupp     else beta[0] = 0.;
6338774f0aSPeter Brune   } else {
6438774f0aSPeter Brune     ierr          = PetscBLASIntCast(l,&ngmres->m);CHKERRQ(ierr);
6538774f0aSPeter Brune     ierr          = PetscBLASIntCast(l,&ngmres->n);CHKERRQ(ierr);
6638774f0aSPeter Brune     ngmres->info  = 0;
6738774f0aSPeter Brune     ngmres->rcond = -1.;
6838774f0aSPeter Brune     ierr          = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
6938774f0aSPeter Brune #if defined(PETSC_USE_COMPLEX)
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->rwork,&ngmres->info));
7138774f0aSPeter Brune #else
728b83055fSJed 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));
7338774f0aSPeter Brune #endif
7438774f0aSPeter Brune     ierr = PetscFPTrapPop();CHKERRQ(ierr);
75*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(ngmres->info < 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
76*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(ngmres->info > 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
7738774f0aSPeter Brune   }
7838774f0aSPeter Brune   for (i=0; i<l; i++) {
79*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PetscIsInfOrNanScalar(beta[i]),PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
8038774f0aSPeter Brune   }
8138774f0aSPeter Brune   alph_total = 0.;
821aa26658SKarl Rupp   for (i = 0; i < l; i++) alph_total += beta[i];
831aa26658SKarl Rupp 
8438774f0aSPeter Brune   ierr = VecCopy(XM,XA);CHKERRQ(ierr);
8538774f0aSPeter Brune   ierr = VecScale(XA,1.-alph_total);CHKERRQ(ierr);
8638774f0aSPeter Brune   ierr = VecMAXPY(XA,l,beta,Xdot);CHKERRQ(ierr);
8738774f0aSPeter Brune   /* check the validity of the step */
8838774f0aSPeter Brune   ierr = VecCopy(XA,Y);CHKERRQ(ierr);
8938774f0aSPeter Brune   ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr);
9038774f0aSPeter Brune   ierr = SNESLineSearchPostCheck(snes->linesearch,X,Y,XA,&changed_y,&changed_w);CHKERRQ(ierr);
9146159c86SPeter Brune   if (!ngmres->approxfunc) {
92efd4aadfSBarry Smith     if (snes->npc && snes->npcside== PC_LEFT) {
93be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,XA,NULL,FA);CHKERRQ(ierr);
9446159c86SPeter Brune     } else {
9546159c86SPeter Brune       ierr = SNESComputeFunction(snes,XA,FA);CHKERRQ(ierr);
9646159c86SPeter Brune     }
9735f895b4SBarry Smith   } else {
98077c4231SPeter Brune     ierr = VecCopy(FM,FA);CHKERRQ(ierr);
99077c4231SPeter Brune     ierr = VecScale(FA,1.-alph_total);CHKERRQ(ierr);
100077c4231SPeter Brune     ierr = VecMAXPY(FA,l,beta,Fdot);CHKERRQ(ierr);
101077c4231SPeter Brune   }
10238774f0aSPeter Brune   PetscFunctionReturn(0);
10338774f0aSPeter Brune }
10438774f0aSPeter Brune 
105b3c6a99cSPeter 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)
10638774f0aSPeter Brune {
10738774f0aSPeter Brune   PetscErrorCode ierr;
10838774f0aSPeter Brune   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
109b3c6a99cSPeter Brune   PetscReal      dcurnorm,dmin = -1.0;
11038774f0aSPeter Brune   Vec            *Xdot = ngmres->Xdot;
11138774f0aSPeter Brune   PetscInt       i;
11238774f0aSPeter Brune 
11338774f0aSPeter Brune   PetscFunctionBegin;
114b3c6a99cSPeter Brune   if (xMnorm) {
115b3c6a99cSPeter Brune     ierr = VecNormBegin(XM,NORM_2,xMnorm);CHKERRQ(ierr);
116b3c6a99cSPeter Brune   }
117b3c6a99cSPeter Brune   if (fMnorm) {
118b3c6a99cSPeter Brune     ierr = VecNormBegin(FM,NORM_2,fMnorm);CHKERRQ(ierr);
119b3c6a99cSPeter Brune   }
120b3c6a99cSPeter Brune   if (yMnorm) {
121b3c6a99cSPeter Brune     ierr = VecCopy(X,D);CHKERRQ(ierr);
122b3c6a99cSPeter Brune     ierr = VecAXPY(D,-1.0,XM);CHKERRQ(ierr);
123b3c6a99cSPeter Brune     ierr = VecNormBegin(D,NORM_2,yMnorm);CHKERRQ(ierr);
124b3c6a99cSPeter Brune   }
125b3c6a99cSPeter Brune   if (xAnorm) {
126b3c6a99cSPeter Brune     ierr = VecNormBegin(XA,NORM_2,xAnorm);CHKERRQ(ierr);
127b3c6a99cSPeter Brune   }
12838774f0aSPeter Brune   if (fAnorm) {
12938774f0aSPeter Brune     ierr = VecNormBegin(FA,NORM_2,fAnorm);CHKERRQ(ierr);
13038774f0aSPeter Brune   }
131b3c6a99cSPeter Brune   if (yAnorm) {
132b3c6a99cSPeter Brune     ierr = VecCopy(X,D);CHKERRQ(ierr);
133b3c6a99cSPeter Brune     ierr = VecAXPY(D,-1.0,XA);CHKERRQ(ierr);
134b3c6a99cSPeter Brune     ierr = VecNormBegin(D,NORM_2,yAnorm);CHKERRQ(ierr);
135b3c6a99cSPeter Brune   }
13638774f0aSPeter Brune   if (dnorm) {
13738774f0aSPeter Brune     ierr = VecCopy(XA,D);CHKERRQ(ierr);
13838774f0aSPeter Brune     ierr = VecAXPY(D,-1.0,XM);CHKERRQ(ierr);
13938774f0aSPeter Brune     ierr = VecNormBegin(D,NORM_2,dnorm);CHKERRQ(ierr);
14038774f0aSPeter Brune   }
14138774f0aSPeter Brune   if (dminnorm) {
14238774f0aSPeter Brune     for (i=0; i<l; i++) {
143b3c6a99cSPeter Brune       ierr = VecCopy(Xdot[i],D);CHKERRQ(ierr);
144b3c6a99cSPeter Brune       ierr = VecAXPY(D,-1.0,XA);CHKERRQ(ierr);
145b3c6a99cSPeter Brune       ierr = VecNormBegin(D,NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr);
14638774f0aSPeter Brune     }
14738774f0aSPeter Brune   }
148b3c6a99cSPeter Brune   if (xMnorm) {ierr = VecNormEnd(XM,NORM_2,xMnorm);CHKERRQ(ierr);}
149b3c6a99cSPeter Brune   if (fMnorm) {ierr = VecNormEnd(FM,NORM_2,fMnorm);CHKERRQ(ierr);}
150b3c6a99cSPeter Brune   if (yMnorm) {ierr = VecNormEnd(D,NORM_2,yMnorm);CHKERRQ(ierr);}
151b3c6a99cSPeter Brune   if (xAnorm) {ierr = VecNormEnd(XA,NORM_2,xAnorm);CHKERRQ(ierr);}
15238774f0aSPeter Brune   if (fAnorm) {ierr = VecNormEnd(FA,NORM_2,fAnorm);CHKERRQ(ierr);}
153b3c6a99cSPeter Brune   if (yAnorm) {ierr = VecNormEnd(D,NORM_2,yAnorm);CHKERRQ(ierr);}
15438774f0aSPeter Brune   if (dnorm) {ierr = VecNormEnd(D,NORM_2,dnorm);CHKERRQ(ierr);}
15538774f0aSPeter Brune   if (dminnorm) {
15638774f0aSPeter Brune     for (i=0; i<l; i++) {
157b3c6a99cSPeter Brune       ierr = VecNormEnd(D,NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr);
15838774f0aSPeter Brune       dcurnorm = ngmres->xnorms[i];
159b3c6a99cSPeter Brune       if ((dcurnorm < dmin) || (dmin < 0.0)) dmin = dcurnorm;
16038774f0aSPeter Brune     }
161b3c6a99cSPeter Brune     *dminnorm = dmin;
16238774f0aSPeter Brune   }
16338774f0aSPeter Brune   PetscFunctionReturn(0);
16438774f0aSPeter Brune }
16538774f0aSPeter Brune 
166b3c6a99cSPeter 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)
16738774f0aSPeter Brune {
16838774f0aSPeter Brune   SNES_NGMRES          *ngmres = (SNES_NGMRES*) snes->data;
16938774f0aSPeter Brune   PetscErrorCode       ierr;
170422a814eSBarry Smith   SNESLineSearchReason lssucceed;
171422a814eSBarry Smith   PetscBool            selectA;
17238774f0aSPeter Brune 
17338774f0aSPeter Brune   PetscFunctionBegin;
17438774f0aSPeter Brune   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
17538774f0aSPeter Brune     /* X = X + \lambda(XA - X) */
17638774f0aSPeter Brune     if (ngmres->monitor) {
17738774f0aSPeter Brune       ierr = PetscViewerASCIIPrintf(ngmres->monitor,"||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr);
17838774f0aSPeter Brune     }
17938774f0aSPeter Brune     ierr   = VecCopy(FM,F);CHKERRQ(ierr);
18038774f0aSPeter Brune     ierr   = VecCopy(XM,X);CHKERRQ(ierr);
18138774f0aSPeter Brune     ierr   = VecCopy(XA,Y);CHKERRQ(ierr);
18238774f0aSPeter Brune     ierr   = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
18338774f0aSPeter Brune     *fnorm = fMnorm;
18438774f0aSPeter Brune     ierr   = SNESLineSearchApply(ngmres->additive_linesearch,X,F,fnorm,Y);CHKERRQ(ierr);
185422a814eSBarry Smith     ierr   = SNESLineSearchGetReason(ngmres->additive_linesearch,&lssucceed);CHKERRQ(ierr);
186422a814eSBarry Smith     ierr   = SNESLineSearchGetNorms(ngmres->additive_linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
187422a814eSBarry Smith     if (lssucceed) {
18838774f0aSPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
18938774f0aSPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
19038774f0aSPeter Brune         PetscFunctionReturn(0);
19138774f0aSPeter Brune       }
19238774f0aSPeter Brune     }
19338774f0aSPeter Brune     if (ngmres->monitor) {
19438774f0aSPeter Brune       ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Additive solution: ||F||_2 = %e\n",*fnorm);CHKERRQ(ierr);
19538774f0aSPeter Brune     }
19638774f0aSPeter Brune   } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
19738774f0aSPeter Brune     selectA = PETSC_TRUE;
19838774f0aSPeter Brune     /* Conditions for choosing the accelerated answer */
19938774f0aSPeter Brune     /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
2001aa26658SKarl Rupp     if (fAnorm >= ngmres->gammaA*fminnorm) selectA = PETSC_FALSE;
2011aa26658SKarl Rupp 
20238774f0aSPeter Brune     /* Criterion B -- the choice of x^A isn't too close to some other choice */
20338774f0aSPeter Brune     if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(*fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) {
2041aa26658SKarl Rupp     } else selectA=PETSC_FALSE;
2051aa26658SKarl Rupp 
20638774f0aSPeter Brune     if (selectA) {
20738774f0aSPeter Brune       if (ngmres->monitor) {
20838774f0aSPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr);
20938774f0aSPeter Brune       }
21038774f0aSPeter Brune       /* copy it over */
211b3c6a99cSPeter Brune       *xnorm = xAnorm;
21238774f0aSPeter Brune       *fnorm = fAnorm;
213b3c6a99cSPeter Brune       *ynorm = yAnorm;
21438774f0aSPeter Brune       ierr   = VecCopy(FA,F);CHKERRQ(ierr);
21538774f0aSPeter Brune       ierr   = VecCopy(XA,X);CHKERRQ(ierr);
21638774f0aSPeter Brune     } else {
21738774f0aSPeter Brune       if (ngmres->monitor) {
21838774f0aSPeter Brune         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr);
21938774f0aSPeter Brune       }
220b3c6a99cSPeter Brune       *xnorm = xMnorm;
22138774f0aSPeter Brune       *fnorm = fMnorm;
222b3c6a99cSPeter Brune       *ynorm = yMnorm;
22338774f0aSPeter Brune       ierr   = VecCopy(XM,Y);CHKERRQ(ierr);
22438774f0aSPeter Brune       ierr   = VecAXPY(Y,-1.0,X);CHKERRQ(ierr);
22538774f0aSPeter Brune       ierr   = VecCopy(FM,F);CHKERRQ(ierr);
22638774f0aSPeter Brune       ierr   = VecCopy(XM,X);CHKERRQ(ierr);
22738774f0aSPeter Brune     }
22838774f0aSPeter Brune   } else { /* none */
229b3c6a99cSPeter Brune     *xnorm = xAnorm;
23038774f0aSPeter Brune     *fnorm = fAnorm;
231b3c6a99cSPeter Brune     *ynorm = yAnorm;
23238774f0aSPeter Brune     ierr   = VecCopy(FA,F);CHKERRQ(ierr);
23338774f0aSPeter Brune     ierr   = VecCopy(XA,X);CHKERRQ(ierr);
23438774f0aSPeter Brune   }
23538774f0aSPeter Brune   PetscFunctionReturn(0);
23638774f0aSPeter Brune }
23738774f0aSPeter Brune 
23823b3e82cSAsbjørn Nilsen Riseth PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes,PetscInt l,PetscReal fMnorm, PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,PetscBool *selectRestart)
23938774f0aSPeter Brune {
24038774f0aSPeter Brune   SNES_NGMRES    *ngmres = (SNES_NGMRES*)snes->data;
24138774f0aSPeter Brune   PetscErrorCode ierr;
24238774f0aSPeter Brune 
24338774f0aSPeter Brune   PetscFunctionBegin;
24438774f0aSPeter Brune   *selectRestart = PETSC_FALSE;
24538774f0aSPeter Brune   /* difference stagnation restart */
24621687c63SPeter Brune   if ((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm)) && l > 0) {
24738774f0aSPeter Brune     if (ngmres->monitor) {
24838774f0aSPeter Brune       ierr = PetscViewerASCIIPrintf(ngmres->monitor,"difference restart: %e > %e\n",ngmres->epsilonB*dnorm,dminnorm);CHKERRQ(ierr);
24938774f0aSPeter Brune     }
25038774f0aSPeter Brune     *selectRestart = PETSC_TRUE;
25138774f0aSPeter Brune   }
25238774f0aSPeter Brune   /* residual stagnation restart */
25338774f0aSPeter Brune   if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) {
25438774f0aSPeter Brune     if (ngmres->monitor) {
25538774f0aSPeter Brune       ierr = PetscViewerASCIIPrintf(ngmres->monitor,"residual restart: %e > %e\n",PetscSqrtReal(fAnorm),ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr);
25638774f0aSPeter Brune     }
25738774f0aSPeter Brune     *selectRestart = PETSC_TRUE;
25838774f0aSPeter Brune   }
25923b3e82cSAsbjørn Nilsen Riseth 
26023b3e82cSAsbjørn Nilsen Riseth   /* F_M stagnation restart */
26123b3e82cSAsbjørn Nilsen Riseth   if (ngmres->restart_fm_rise && fMnorm > snes->norm) {
26223b3e82cSAsbjørn Nilsen Riseth     if (ngmres->monitor) {
26323b3e82cSAsbjørn Nilsen Riseth       ierr = PetscViewerASCIIPrintf(ngmres->monitor,"F_M rise restart: %e > %e\n",fMnorm,snes->norm);CHKERRQ(ierr);
26423b3e82cSAsbjørn Nilsen Riseth     }
26523b3e82cSAsbjørn Nilsen Riseth     *selectRestart = PETSC_TRUE;
26623b3e82cSAsbjørn Nilsen Riseth   }
26723b3e82cSAsbjørn Nilsen Riseth 
26838774f0aSPeter Brune   PetscFunctionReturn(0);
26938774f0aSPeter Brune }
270