Actual source code: ngmresfunc.c

petsc-3.4.2 2013-07-02
  1: #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
  2: #include <petscblaslapack.h>

  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:   PetscScalar    *xi     = ngmres->xi;
 12:   PetscInt       i;
 13:   PetscReal      nu;

 17:   if (ivec > l) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot update vector %d with space size %d!",ivec,l);
 18:   VecCopy(F,Fdot[ivec]);
 19:   VecCopy(X,Xdot[ivec]);

 21:   ngmres->fnorms[ivec] = fnorm;
 22:   if (l > 0) {
 23:     VecMDot(F,l,Fdot,xi);
 24:     for (i = 0; i < l; i++) {
 25:       Q(i,ivec) = xi[i];
 26:       Q(ivec,i) = xi[i];
 27:     }
 28:   } else {
 29:     nu     = fnorm*fnorm;
 30:     Q(0,0) = nu;
 31:   }
 32:   return(0);
 33: }

 37: PetscErrorCode SNESNGMRESFormCombinedSolution_Private(SNES snes,PetscInt l,Vec XM,Vec FM,PetscReal fMnorm,Vec X,Vec XA,Vec FA)
 38: {
 39:   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
 40:   PetscInt       i,j;
 41:   Vec            *Fdot      = ngmres->Fdot;
 42:   Vec            *Xdot      = ngmres->Xdot;
 43:   PetscScalar    *beta      = ngmres->beta;
 44:   PetscScalar    *xi        = ngmres->xi;
 45:   PetscScalar    alph_total = 0.;
 47:   PetscReal      nu;
 48:   Vec            Y = snes->work[2];
 49:   PetscBool      changed_y,changed_w;

 52:   nu = fMnorm*fMnorm;

 54:   /* construct the right hand side and xi factors */
 55:   VecMDot(FM,l,Fdot,xi);
 56:   for (i = 0; i < l; i++) beta[i] = nu - xi[i];

 58:   /* construct h */
 59:   for (j = 0; j < l; j++) {
 60:     for (i = 0; i < l; i++) {
 61:       H(i,j) = Q(i,j)-xi[i]-xi[j]+nu;
 62:     }
 63:   }
 64:   if (l == 1) {
 65:     /* simply set alpha[0] = beta[0] / H[0, 0] */
 66:     if (H(0,0) != 0.) beta[0] = beta[0]/H(0,0);
 67:     else beta[0] = 0.;
 68:   } else {
 69: #if defined(PETSC_MISSING_LAPACK_GELSS)
 70:     SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"NGMRES with LS requires the LAPACK GELSS routine.");
 71: #else
 72:     PetscBLASIntCast(l,&ngmres->m);
 73:     PetscBLASIntCast(l,&ngmres->n);
 74:     ngmres->info  = 0;
 75:     ngmres->rcond = -1.;
 76:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
 77: #if defined(PETSC_USE_COMPLEX)
 78:     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));
 79: #else
 80:     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));
 81: #endif
 82:     PetscFPTrapPop();
 83:     if (ngmres->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
 84:     if (ngmres->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
 85: #endif
 86:   }
 87:   for (i=0; i<l; i++) {
 88:     if (PetscIsInfOrNanScalar(beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
 89:   }
 90:   alph_total = 0.;
 91:   for (i = 0; i < l; i++) alph_total += beta[i];

 93:   VecCopy(XM,XA);
 94:   VecScale(XA,1.-alph_total);
 95:   VecMAXPY(XA,l,beta,Xdot);
 96:   /* check the validity of the step */
 97:   VecCopy(XA,Y);
 98:   VecAXPY(Y,-1.0,X);
 99:   SNESLineSearchPostCheck(snes->linesearch,X,Y,XA,&changed_y,&changed_w);
100:   if (!ngmres->approxfunc) {SNESComputeFunction(snes,XA,FA);}
101:   else {
102:     VecCopy(FM,FA);
103:     VecScale(FA,1.-alph_total);
104:     VecMAXPY(FA,l,beta,Fdot);
105:   }
106:   return(0);
107: }

111: PetscErrorCode SNESNGMRESCalculateDifferences_Private(SNES snes,PetscInt l,Vec X,Vec F,Vec XM,Vec FM,Vec XA,Vec FA,Vec D,PetscReal *dnorm,PetscReal *dminnorm,PetscReal *fAnorm)
112: {
114:   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
115:   PetscReal      dcurnorm;
116:   Vec            *Xdot = ngmres->Xdot;
117:   PetscInt       i;

120:   if (ngmres->singlereduction) {
121:     *dminnorm = -1.0;
122:     if (fAnorm) {
123:       VecNormBegin(FA,NORM_2,fAnorm);
124:     }
125:     if (dnorm) {
126:       VecCopy(XA,D);
127:       VecAXPY(D,-1.0,XM);
128:       VecNormBegin(D,NORM_2,dnorm);
129:     }
130:     if (dminnorm) {
131:       for (i=0; i<l; i++) {
132:         ierr=VecAXPY(Xdot[i],-1.0,XA);
133:       }
134:       for (i=0; i<l; i++) {
135:         VecNormBegin(Xdot[i],NORM_2,&ngmres->xnorms[i]);
136:       }
137:     }
138:     if (fAnorm) {VecNormEnd(FA,NORM_2,fAnorm);}
139:     if (dnorm) {VecNormEnd(D,NORM_2,dnorm);}
140:     if (dminnorm) {
141:       for (i=0; i<l; i++) {
142:         VecNormEnd(Xdot[i],NORM_2,&ngmres->xnorms[i]);
143:       }
144:       for (i=0; i<l; i++) {
145:         dcurnorm = ngmres->xnorms[i];
146:         if ((dcurnorm < *dminnorm) || (*dminnorm < 0.0)) *dminnorm = dcurnorm;
147:         ierr=VecAXPY(Xdot[i],1.0,XA);
148:       }
149:     }
150:   } else {
151:     if (dnorm) {
152:       ierr=VecCopy(XA,D);
153:       ierr=VecAXPY(D,-1.0,XM);
154:       ierr=VecNormBegin(D,NORM_2,dnorm);
155:     }
156:     if (fAnorm) {
157:       ierr=VecNormBegin(FA,NORM_2,fAnorm);
158:     }
159:     if (dnorm) {
160:       ierr=VecNormEnd(D,NORM_2,dnorm);
161:     }
162:     if (fAnorm) {
163:       ierr=VecNormEnd(FA,NORM_2,fAnorm);
164:     }
165:     if (dminnorm) {
166:       *dminnorm = -1.0;
167:       for (i=0; i<l; i++) {
168:         ierr=VecCopy(XA,D);
169:         ierr=VecAXPY(D,-1.0,Xdot[i]);
170:         ierr=VecNorm(D,NORM_2,&dcurnorm);
171:         if ((dcurnorm < *dminnorm) || (*dminnorm < 0.0)) *dminnorm = dcurnorm;
172:       }
173:     }
174:   }
175:   return(0);
176: }

180: PetscErrorCode SNESNGMRESSelect_Private(SNES snes,PetscInt k_restart,Vec XM,Vec FM,PetscReal fMnorm,Vec XA,Vec FA,PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,Vec X,Vec F,Vec Y,PetscReal *fnorm)
181: {
182:   SNES_NGMRES    *ngmres = (SNES_NGMRES*) snes->data;
184:   PetscBool      lssucceed,selectA;

187:   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
188:     /* X = X + \lambda(XA - X) */
189:     if (ngmres->monitor) {
190:       PetscViewerASCIIPrintf(ngmres->monitor,"||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);
191:     }
192:     VecCopy(FM,F);
193:     VecCopy(XM,X);
194:     VecCopy(XA,Y);
195:     VecAYPX(Y,-1.0,X);
196:     *fnorm = fMnorm;
197:     SNESLineSearchApply(ngmres->additive_linesearch,X,F,fnorm,Y);
198:     SNESLineSearchGetSuccess(ngmres->additive_linesearch,&lssucceed);
199:     if (!lssucceed) {
200:       if (++snes->numFailures >= snes->maxFailures) {
201:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
202:         return(0);
203:       }
204:     }
205:     if (ngmres->monitor) {
206:       PetscViewerASCIIPrintf(ngmres->monitor,"Additive solution: ||F||_2 = %e\n",*fnorm);
207:     }
208:   } else if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
209:     selectA = PETSC_TRUE;
210:     /* Conditions for choosing the accelerated answer */
211:     /* Criterion A -- the norm of the function isn't increased above the minimum by too much */
212:     if (fAnorm >= ngmres->gammaA*fminnorm) selectA = PETSC_FALSE;

214:     /* Criterion B -- the choice of x^A isn't too close to some other choice */
215:     if (ngmres->epsilonB*dnorm<dminnorm || PetscSqrtReal(*fnorm)<ngmres->deltaB*PetscSqrtReal(fminnorm)) {
216:     } else selectA=PETSC_FALSE;

218:     if (selectA) {
219:       if (ngmres->monitor) {
220:         PetscViewerASCIIPrintf(ngmres->monitor,"picked X_A, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);
221:       }
222:       /* copy it over */
223:       *fnorm = fAnorm;
224:       VecCopy(FA,F);
225:       VecCopy(XA,X);
226:     } else {
227:       if (ngmres->monitor) {
228:         PetscViewerASCIIPrintf(ngmres->monitor,"picked X_M, ||F_A||_2 = %e, ||F_M||_2 = %e\n",fAnorm,fMnorm);
229:       }
230:       *fnorm = fMnorm;
231:       VecCopy(XM,Y);
232:       VecAXPY(Y,-1.0,X);
233:       VecCopy(FM,F);
234:       VecCopy(XM,X);
235:     }
236:   } else { /* none */
237:     *fnorm = fAnorm;
238:     VecCopy(FA,F);
239:     VecCopy(XA,X);
240:   }
241:   return(0);
242: }

246: PetscErrorCode SNESNGMRESSelectRestart_Private(SNES snes,PetscInt l,PetscReal fAnorm,PetscReal dnorm,PetscReal fminnorm,PetscReal dminnorm,PetscBool *selectRestart)
247: {
248:   SNES_NGMRES    *ngmres = (SNES_NGMRES*)snes->data;

252:   *selectRestart = PETSC_FALSE;
253:   /* difference stagnation restart */
254:   if ((ngmres->epsilonB*dnorm > dminnorm) && (PetscSqrtReal(fAnorm) > ngmres->deltaB*PetscSqrtReal(fminnorm)) && l > 0) {
255:     if (ngmres->monitor) {
256:       PetscViewerASCIIPrintf(ngmres->monitor,"difference restart: %e > %e\n",ngmres->epsilonB*dnorm,dminnorm);
257:     }
258:     *selectRestart = PETSC_TRUE;
259:   }
260:   /* residual stagnation restart */
261:   if (PetscSqrtReal(fAnorm) > ngmres->gammaC*PetscSqrtReal(fminnorm)) {
262:     if (ngmres->monitor) {
263:       PetscViewerASCIIPrintf(ngmres->monitor,"residual restart: %e > %e\n",PetscSqrtReal(fAnorm),ngmres->gammaC*PetscSqrtReal(fminnorm));
264:     }
265:     *selectRestart = PETSC_TRUE;
266:   }
267:   return(0);
268: }