xref: /petsc/src/snes/impls/ngmres/ngmresfunc.c (revision a7fac7c2b47dbcd8ece0cc2dfdfe0e63be1bb7b5)
1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2 #include <petscblaslapack.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESNGMRESUpdateSubspace_Private"
6 PetscErrorCode SNESNGMRESUpdateSubspace_Private(SNES snes,PetscInt ivec,PetscInt l,Vec F,PetscReal fnorm,Vec X)
7 {
8   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
9   Vec            *Fdot   = ngmres->Fdot;
10   Vec            *Xdot   = ngmres->Xdot;
11   PetscErrorCode ierr;
12 
13   PetscFunctionBegin;
14   if (ivec > l) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot update vector %d with space size %d!",ivec,l);
15   ierr = VecCopy(F,Fdot[ivec]);CHKERRQ(ierr);
16   ierr = VecCopy(X,Xdot[ivec]);CHKERRQ(ierr);
17 
18   ngmres->fnorms[ivec] = fnorm;
19   PetscFunctionReturn(0);
20 }
21 
22 #undef __FUNCT__
23 #define __FUNCT__ "SNESNGMRESFormCombinedSolution_Private"
24 PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES snes,PetscInt ivec,PetscInt l,Vec XM,Vec FM,PetscReal fMnorm,Vec X,Vec XA,Vec FA)
25 {
26   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
27   PetscInt       i,j;
28   Vec            *Fdot      = ngmres->Fdot;
29   Vec            *Xdot      = ngmres->Xdot;
30   PetscScalar    *beta      = ngmres->beta;
31   PetscScalar    *xi        = ngmres->xi;
32   PetscScalar    alph_total = 0.;
33   PetscErrorCode ierr;
34   PetscReal      nu;
35   Vec            Y = snes->work[2];
36   PetscBool      changed_y,changed_w;
37 
38   PetscFunctionBegin;
39   nu = fMnorm*fMnorm;
40 
41   /* construct the right hand side and xi factors */
42   if (l > 0) {
43     ierr = VecMDotBegin(FM,l,Fdot,xi);CHKERRQ(ierr);
44     ierr = VecMDotBegin(Fdot[ivec],l,Fdot,beta);CHKERRQ(ierr);
45     ierr = VecMDotEnd(FM,l,Fdot,xi);CHKERRQ(ierr);
46     ierr = VecMDotEnd(Fdot[ivec],l,Fdot,beta);CHKERRQ(ierr);
47     for (i = 0; i < l; i++) {
48       Q(i,ivec) = beta[i];
49       Q(ivec,i) = beta[i];
50     }
51   } else {
52     Q(0,0) = ngmres->fnorms[ivec]*ngmres->fnorms[ivec];
53   }
54 
55   for (i = 0; i < l; i++) beta[i] = nu - xi[i];
56 
57   /* construct h */
58   for (j = 0; j < l; j++) {
59     for (i = 0; i < l; i++) {
60       H(i,j) = Q(i,j)-xi[i]-xi[j]+nu;
61     }
62   }
63   if (l == 1) {
64     /* simply set alpha[0] = beta[0] / H[0, 0] */
65     if (H(0,0) != 0.) beta[0] = beta[0]/H(0,0);
66     else beta[0] = 0.;
67   } else {
68 #if defined(PETSC_MISSING_LAPACK_GELSS)
69     SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"NGMRES with LS requires the LAPACK GELSS routine.");
70 #else
71     ierr          = PetscBLASIntCast(l,&ngmres->m);CHKERRQ(ierr);
72     ierr          = PetscBLASIntCast(l,&ngmres->n);CHKERRQ(ierr);
73     ngmres->info  = 0;
74     ngmres->rcond = -1.;
75     ierr          = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
76 #if defined(PETSC_USE_COMPLEX)
77     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));
78 #else
79     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));
80 #endif
81     ierr = PetscFPTrapPop();CHKERRQ(ierr);
82     if (ngmres->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
83     if (ngmres->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
84 #endif
85   }
86   for (i=0; i<l; i++) {
87     if (PetscIsInfOrNanScalar(beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
88   }
89   alph_total = 0.;
90   for (i = 0; i < l; i++) alph_total += beta[i];
91 
92   ierr = VecCopy(XM,XA);CHKERRQ(ierr);
93   ierr = VecScale(XA,1.-alph_total);CHKERRQ(ierr);
94   ierr = VecMAXPY(XA,l,beta,Xdot);CHKERRQ(ierr);
95   /* check the validity of the step */
96   ierr = VecCopy(XA,Y);CHKERRQ(ierr);
97   ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr);
98   ierr = SNESLineSearchPostCheck(snes->linesearch,X,Y,XA,&changed_y,&changed_w);CHKERRQ(ierr);
99   if (!ngmres->approxfunc) {
100     if (snes->pc && snes->pcside == PC_LEFT) {
101       ierr = SNESApplyNPC(snes,XA,NULL,FA);CHKERRQ(ierr);
102     } else {
103       ierr =SNESComputeFunction(snes,XA,FA);CHKERRQ(ierr);
104     }
105   }
106   else {
107     ierr = VecCopy(FM,FA);CHKERRQ(ierr);
108     ierr = VecScale(FA,1.-alph_total);CHKERRQ(ierr);
109     ierr = VecMAXPY(FA,l,beta,Fdot);CHKERRQ(ierr);
110   }
111   PetscFunctionReturn(0);
112 }
113 
114 #undef __FUNCT__
115 #define __FUNCT__ "SNESNGMRESNorms_Private"
116 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)
117 {
118   PetscErrorCode ierr;
119   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
120   PetscReal      dcurnorm,dmin = -1.0;
121   Vec            *Xdot = ngmres->Xdot;
122   PetscInt       i;
123 
124   PetscFunctionBegin;
125   if (xMnorm) {
126     ierr = VecNormBegin(XM,NORM_2,xMnorm);CHKERRQ(ierr);
127   }
128   if (fMnorm) {
129     ierr = VecNormBegin(FM,NORM_2,fMnorm);CHKERRQ(ierr);
130   }
131   if (yMnorm) {
132     ierr = VecCopy(X,D);CHKERRQ(ierr);
133     ierr = VecAXPY(D,-1.0,XM);CHKERRQ(ierr);
134     ierr = VecNormBegin(D,NORM_2,yMnorm);CHKERRQ(ierr);
135   }
136   if (xAnorm) {
137     ierr = VecNormBegin(XA,NORM_2,xAnorm);CHKERRQ(ierr);
138   }
139   if (fAnorm) {
140     ierr = VecNormBegin(FA,NORM_2,fAnorm);CHKERRQ(ierr);
141   }
142   if (yAnorm) {
143     ierr = VecCopy(X,D);CHKERRQ(ierr);
144     ierr = VecAXPY(D,-1.0,XA);CHKERRQ(ierr);
145     ierr = VecNormBegin(D,NORM_2,yAnorm);CHKERRQ(ierr);
146   }
147   if (dnorm) {
148     ierr = VecCopy(XA,D);CHKERRQ(ierr);
149     ierr = VecAXPY(D,-1.0,XM);CHKERRQ(ierr);
150     ierr = VecNormBegin(D,NORM_2,dnorm);CHKERRQ(ierr);
151   }
152   if (dminnorm) {
153     for (i=0; i<l; i++) {
154       ierr = VecCopy(Xdot[i],D);CHKERRQ(ierr);
155       ierr = VecAXPY(D,-1.0,XA);CHKERRQ(ierr);
156       ierr = VecNormBegin(D,NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr);
157     }
158   }
159   if (xMnorm) {ierr = VecNormEnd(XM,NORM_2,xMnorm);CHKERRQ(ierr);}
160   if (fMnorm) {ierr = VecNormEnd(FM,NORM_2,fMnorm);CHKERRQ(ierr);}
161   if (yMnorm) {ierr = VecNormEnd(D,NORM_2,yMnorm);CHKERRQ(ierr);}
162   if (xAnorm) {ierr = VecNormEnd(XA,NORM_2,xAnorm);CHKERRQ(ierr);}
163   if (fAnorm) {ierr = VecNormEnd(FA,NORM_2,fAnorm);CHKERRQ(ierr);}
164   if (yAnorm) {ierr = VecNormEnd(D,NORM_2,yAnorm);CHKERRQ(ierr);}
165   if (dnorm) {ierr = VecNormEnd(D,NORM_2,dnorm);CHKERRQ(ierr);}
166   if (dminnorm) {
167     for (i=0; i<l; i++) {
168       ierr = VecNormEnd(D,NORM_2,&ngmres->xnorms[i]);CHKERRQ(ierr);
169       dcurnorm = ngmres->xnorms[i];
170       if ((dcurnorm < dmin) || (dmin < 0.0)) dmin = dcurnorm;
171     }
172     *dminnorm = dmin;
173   }
174   PetscFunctionReturn(0);
175 }
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "SNESNGMRESSelect_Private"
179 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)
180 {
181   SNES_NGMRES          *ngmres = (SNES_NGMRES*) snes->data;
182   PetscErrorCode       ierr;
183   SNESLineSearchReason lssucceed;
184   PetscBool            selectA;
185 
186   PetscFunctionBegin;
187   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
188     /* X = X + \lambda(XA - X) */
189     if (ngmres->monitor) {
190       ierr = PetscViewerASCIIPrintf(ngmres->monitor,"||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr);
191     }
192     ierr   = VecCopy(FM,F);CHKERRQ(ierr);
193     ierr   = VecCopy(XM,X);CHKERRQ(ierr);
194     ierr   = VecCopy(XA,Y);CHKERRQ(ierr);
195     ierr   = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
196     *fnorm = fMnorm;
197     ierr   = SNESLineSearchApply(ngmres->additive_linesearch,X,F,fnorm,Y);CHKERRQ(ierr);
198     ierr   = SNESLineSearchGetReason(ngmres->additive_linesearch,&lssucceed);CHKERRQ(ierr);
199     ierr   = SNESLineSearchGetNorms(ngmres->additive_linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
200     if (lssucceed) {
201       if (++snes->numFailures >= snes->maxFailures) {
202         snes->reason = SNES_DIVERGED_LINE_SEARCH;
203         PetscFunctionReturn(0);
204       }
205     }
206     if (ngmres->monitor) {
207       ierr = PetscViewerASCIIPrintf(ngmres->monitor,"Additive solution: ||F||_2 = %e\n",*fnorm);CHKERRQ(ierr);
208     }
209   } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
210     selectA = PETSC_TRUE;
211     /* Conditions for choosing the accelerated answer */
212     /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
213     if (fAnorm >= ngmres->gammaA*fminnorm) selectA = PETSC_FALSE;
214 
215     /* Criterion B -- the choice of x^A isn't too close to some other choice */
216     if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(*fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) {
217     } else selectA=PETSC_FALSE;
218 
219     if (selectA) {
220       if (ngmres->monitor) {
221         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr);
222       }
223       /* copy it over */
224       *xnorm = xAnorm;
225       *fnorm = fAnorm;
226       *ynorm = yAnorm;
227       ierr   = VecCopy(FA,F);CHKERRQ(ierr);
228       ierr   = VecCopy(XA,X);CHKERRQ(ierr);
229     } else {
230       if (ngmres->monitor) {
231         ierr = PetscViewerASCIIPrintf(ngmres->monitor,"picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);CHKERRQ(ierr);
232       }
233       *xnorm = xMnorm;
234       *fnorm = fMnorm;
235       *ynorm = yMnorm;
236       ierr   = VecCopy(XM,Y);CHKERRQ(ierr);
237       ierr   = VecAXPY(Y,-1.0,X);CHKERRQ(ierr);
238       ierr   = VecCopy(FM,F);CHKERRQ(ierr);
239       ierr   = VecCopy(XM,X);CHKERRQ(ierr);
240     }
241   } else { /* none */
242     *xnorm = xAnorm;
243     *fnorm = fAnorm;
244     *ynorm = yAnorm;
245     ierr   = VecCopy(FA,F);CHKERRQ(ierr);
246     ierr   = VecCopy(XA,X);CHKERRQ(ierr);
247   }
248   PetscFunctionReturn(0);
249 }
250 
251 #undef __FUNCT__
252 #define __FUNCT__ "SNESNGMRESSelectRestart_Private"
253 PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes,PetscInt l,PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,PetscBool *selectRestart)
254 {
255   SNES_NGMRES    *ngmres = (SNES_NGMRES*)snes->data;
256   PetscErrorCode ierr;
257 
258   PetscFunctionBegin;
259   *selectRestart = PETSC_FALSE;
260   /* difference stagnation restart */
261   if ((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm)) && l > 0) {
262     if (ngmres->monitor) {
263       ierr = PetscViewerASCIIPrintf(ngmres->monitor,"difference restart: %e > %e\n",ngmres->epsilonB*dnorm,dminnorm);CHKERRQ(ierr);
264     }
265     *selectRestart = PETSC_TRUE;
266   }
267   /* residual stagnation restart */
268   if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) {
269     if (ngmres->monitor) {
270       ierr = PetscViewerASCIIPrintf(ngmres->monitor,"residual restart: %e > %e\n",PetscSqrtReal(fAnorm),ngmres->gammaC*PetscSqrtReal(fminnorm));CHKERRQ(ierr);
271     }
272     *selectRestart = PETSC_TRUE;
273   }
274   PetscFunctionReturn(0);
275 }
276