Actual source code: ngmres.c

  2: #include <../src/ksp/ksp/impls/ngmres/ngmresimpl.h>       /*I "petscksp.h" I*/

  4: static PetscErrorCode    BuildNGmresSoln(PetscScalar*,Vec,KSP,PetscInt,PetscInt);
  5: /*
  6:      KSPSetUp_NGMRES - Sets up the workspace needed by the NGMRES method. 

  8:       This is called once, usually automatically by KSPSolve() or KSPSetUp()
  9: */
 12: PetscErrorCode KSPSetUp_NGMRES(KSP ksp)
 13: {

 15:   PetscInt       hh;
 16:   PetscInt       msize;
 17:   KSP_NGMRES     *ngmres = (KSP_NGMRES*)ksp->data;
 19:   PetscBool      diagonalscale;

 22:   /* 
 23:        This implementation of NGMRES only handles left preconditioning
 24:      so generate an error otherwise.
 25:   */
 26:   if (ksp->pc_side == PC_RIGHT) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"No right preconditioning for KSPNGMRES");
 27:   else if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"No symmetric preconditioning for KSPNGMRES");
 28:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 29:   if (diagonalscale) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
 30:   ngmres->beta  = 1.0;
 31:   msize         = ngmres->msize;  /* restart size */
 32:   hh            = msize * msize;
 33: 
 34:   //  PetscMalloc1(hh,PetscScalar,&ngmres->hh_origin,bb,PetscScalar,&ngmres->bb_origin,rs,PetscScalar,&ngmres->rs_origin,cc,PetscScalar,&ngmres->cc_origin,cc,PetscScalar,&ngmres->ss_origin);
 35:   PetscMalloc2(hh,PetscScalar,&ngmres->hh_origin,msize,PetscScalar,&ngmres->nrs);

 37:   PetscMemzero(ngmres->hh_origin,hh*sizeof(PetscScalar));
 38:   PetscMemzero(ngmres->nrs,msize*sizeof(PetscScalar));

 40:   PetscLogObjectMemory(ksp,(hh+msize)*sizeof(PetscScalar));

 42:   KSPGetVecs(ksp,ngmres->msize,&ngmres->v,ngmres->msize*2,&ngmres->w);
 43:   VecDuplicateVecs(ngmres->w[0], ngmres->msize, &ngmres->q);
 44: 
 45:   KSPDefaultGetWork(ksp,3);
 46:   return(0);
 47: }

 49: /*
 50:        KSPSolve_NGMRES - This routine actually applies the method


 53:    Input Parameter:
 54: .     ksp - the Krylov space object that was set to use conjugate gradient, by, for 
 55:             example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPNGMRES);
 56: */
 59: PetscErrorCode  KSPSolve_NGMRES(KSP ksp)
 60: {
 62:   PetscInt       i,ivec,j,k,l,flag,it;
 63:   KSP_NGMRES    *ngmres = (KSP_NGMRES*)ksp->data;
 64:   Mat            Amat;
 65:   //  Vec            X,F,R, B,Fold, Xold,temp,*dX = ngmres->w,*dF = ngmres->w+ngmres->msize;
 66:   Vec            X,F,R, B,Fold, Xold,temp,*dX = ngmres->w,*dF = ngmres->w+ngmres->msize;
 67:   PetscScalar    *nrs=ngmres->nrs;
 68:   PetscReal      gnorm;
 69: 

 72:   X             = ksp->vec_sol;
 73:   B             = ksp->vec_rhs;
 74:   F             = ksp->work[0];
 75:   Fold          = ksp->work[1];
 76:   R             = ksp->work[2];
 77:   VecDuplicate(X,&Fold);
 78:   VecDuplicate(X,&Xold);
 79:   VecDuplicate(X,&temp);


 82:   PCGetOperators(ksp->pc,&Amat,PETSC_NULL,PETSC_NULL);

 84:   if (!ksp->guess_zero) {
 85:     KSP_MatMult(ksp,Amat,X,R);            /*     r <- b - Ax     */
 86:     VecAYPX(R,-1.0,B);
 87:   } else {
 88:     VecCopy(B,R);                         /*     r <- b (x is 0) */
 89:   }
 90:   PCApply(ksp->pc,R,F);                /*     F = B r */
 91:   if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
 92:     VecNorm(R,NORM_2,&gnorm);
 93:   } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
 94:     VecNorm(F,NORM_2,&gnorm);
 95:   } else SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"NormType not supported");
 96:   KSPLogResidualHistory(ksp,gnorm);
 97:   KSPMonitor(ksp,0,gnorm);
 98:   (*ksp->converged)(ksp,0,gnorm,&ksp->reason,ksp->cnvP);


101:  /* for k=0 */

103:   k=0;
104:   ierr= VecCopy(X,Xold);  /* Xbar_0= X_0 */
105:   ierr= VecCopy(F,Fold);  /* Fbar_0 = f_0= B(b-Ax_0) */
106:   ierr= VecWAXPY(X, ngmres->beta, Fold, Xold);        /*X_1 = X_bar + beta * Fbar */

108:   /* to calculate f_1 */
109:   KSP_MatMult(ksp,Amat,X,R);            /*     r <- b - Ax     */
110:   VecAYPX(R,-1.0,B);
111:   PCApply(ksp->pc,R,F);                /*     F= f_1 = B(b-Ax) */

113:   /* calculate dX and dF for k=0 */
114:   ierr= VecWAXPY(dX[k],-1.0, Xold, X);  /* dX= X_1 - X_0 */
115:   ierr= VecWAXPY(dF[k],-1.0, Fold, F);  /* dF= f_1 - f_0 */
116: 
117: 
118:   ierr= VecCopy(X,Xold);  /* Xbar_0= X_0 */
119:   ierr= VecCopy(F,Fold);  /* Fbar_0 = f_0= B(b-Ax_0) */

121:   flag=0;

123:   for (k=1; k<ksp->max_it; k += 1) {  /* begin the iteration */
124:     l=ngmres->msize;
125:     if(k<l) {
126:       l=k;
127:       ivec=l;
128:     }else{
129:       ivec=l-1;
130:     }
131:     it=l-1;
132:     BuildNGmresSoln(nrs,Fold,ksp,it,flag);
133: 
134:     /* to obtain the solution at k+1 step */
135:     ierr= VecCopy(Xold,X);  /* X=Xold+Fold-(dX + dF) *nrd */
136:      ierr= VecAXPY(X,1.0,Fold);  /* X= Xold+Fold */
137:     for(i=0;i<l;i++){      /*X= Xold+Fold- (dX+dF*beta) *innerb */
138:       ierr= VecAXPY(X,-nrs[i], dX[i]);
139:       ierr= VecAXPY(X,-nrs[i]*ngmres->beta, dF[i]);
140:     }


143:     /* to test with GMRES */
144:     //ierr= VecCopy(Xold,temp);  /* X=Xold-(dX ) *nrd */
145:     /*for(i=0;i<l;i++){      
146:       ierr= VecAXPY(temp,-nrs[i], dX[i]);
147:     }
148:     KSP_MatMult(ksp,Amat,temp,R);            
149:     VecAYPX(R,-1.0,B); 
150:     PCApply(ksp->pc,R,F);                
151:      VecNorm(R,NORM_2,&gnorm); 
152:      printf("gmres residual norm=%g\n",gnorm);
153:     */

155:     /* to calculate f_k+1 */
156:     KSP_MatMult(ksp,Amat,X,R);            /*     r <- b - Ax     */
157:     VecAYPX(R,-1.0,B);
158:     PCApply(ksp->pc,R,F);                /*     F= f_1 = B(b-Ax) */
159: 

161:     /* check the convegence */
162: 
163:     if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
164:         VecNorm(R,NORM_2,&gnorm);
165:       } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
166:         VecNorm(F,NORM_2,&gnorm);
167:       } else SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"NormType not supported");
168:       KSPLogResidualHistory(ksp,gnorm);

170:       KSPMonitor(ksp,k,gnorm);
171:       ksp->its=k;
172:       (*ksp->converged)(ksp,k,gnorm,&ksp->reason,ksp->cnvP);
173:       if (ksp->reason) return(0);


176:     /* calculate dX and dF for k=0 */
177:       if( k>ivec) {/* we need to replace the old vectors */
178:         flag=1;
179:         for(i=0;i<it;i++){
180:           ierr= VecCopy(dX[i+1],dX[i]);  /* X=Xold+Fold-(dX + dF) *nrd */
181:           ierr= VecCopy(dF[i+1],dF[i]);  /* X=Xold+Fold-(dX + dF) *nrd */
182:           for(j=0;j<l;j++)
183:             *HH(j,i)=*HH(j,i+1);
184:         }
185:       }
186:       ierr= VecWAXPY(dX[ivec],-1.0, Xold, X);  /* dX= X_1 - X_0 */
187:       ierr= VecWAXPY(dF[ivec],-1.0, Fold, F);  /* dF= f_1 - f_0 */
188:       ierr= VecCopy(X,Xold);
189:       ierr= VecCopy(F,Fold);
190: 
191:   }
192:   ksp->reason = KSP_DIVERGED_ITS;
193:   return(0);
194: }
197: PetscErrorCode KSPReset_NGMRES(KSP ksp)
198: {
199:   KSP_NGMRES    *ngmres = (KSP_NGMRES*)ksp->data;

203:   VecDestroyVecs(ngmres->msize,&ngmres->v);
204:   VecDestroyVecs(ngmres->msize,&ngmres->w);
205:   VecDestroyVecs(ngmres->msize,&ngmres->q);
206:   return(0);
207: }

209: /*
210:        KSPDestroy_NGMRES - Frees all memory space used by the Krylov method

212: */
215: PetscErrorCode KSPDestroy_NGMRES(KSP ksp)
216: {

220:   KSPReset_NGMRES(ksp);
221:   KSPDefaultDestroy(ksp);
222:   return(0);
223: }

225: /*
226:      KSPView_NGMRES - Prints information about the current Krylov method being used

228:       Currently this only prints information to a file (or stdout) about the 
229:       symmetry of the problem. If your Krylov method has special options or 
230:       flags that information should be printed here.

232: */
235: PetscErrorCode KSPView_NGMRES(KSP ksp,PetscViewer viewer)
236: {
237:   KSP_NGMRES    *ngmres = (KSP_NGMRES *)ksp->data;
239:   PetscBool      iascii;

242:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
243:   if (iascii) {
244:     PetscViewerASCIIPrintf(viewer,"  Size of space %d\n",ngmres->msize);
245:   } else {
246:     SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for KSP ngmres",((PetscObject)viewer)->type_name);
247:   }
248:   return(0);
249: }

251: /*
252:     KSPSetFromOptions_NGMRES - Checks the options database for options related to the  method
253: */
256: PetscErrorCode KSPSetFromOptions_NGMRES(KSP ksp)
257: {
259:   KSP_NGMRES    *ngmres = (KSP_NGMRES *)ksp->data;

262:   PetscOptionsHead("KSP NGMRES options");
263:     PetscOptionsInt("-ksp_ngmres_restart","Number of directions","None",ngmres->msize,&ngmres->msize,PETSC_NULL);
264:   PetscOptionsTail();
265:   return(0);
266: }

268: /*
269:     KSPCreate_NGMRES - Creates the data structure for the Krylov method NGMRES and sets the 
270:        function pointers for all the routines it needs to call (KSPSolve_NGMRES() etc)

273: */
274: /*MC
275:      KSPNGMRES - The nonlinear generalized minimum residual (NGMRES) iterative method

277:    Level: beginner

279:    Notes: Supports only left preconditioning

281: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP

283: M*/
287: PetscErrorCode  KSPCreate_NGMRES(KSP ksp)
288: {
290:   KSP_NGMRES     *ngmres;

293:   PetscNewLog(ksp,KSP_NGMRES,&ngmres);
294:   ksp->data = (void*)ngmres;
295:   ngmres->msize = 30;
296:   ngmres->csize = 0;

298:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
299:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);

301:   /*
302:        Sets the functions that are associated with this data structure 
303:        (in C++ this is the same as defining virtual functions)
304:   */
305:   ksp->ops->setup                = KSPSetUp_NGMRES;
306:   ksp->ops->solve                = KSPSolve_NGMRES;
307:   ksp->ops->reset                = KSPReset_NGMRES;
308:   ksp->ops->destroy              = KSPDestroy_NGMRES;
309:   ksp->ops->view                 = KSPView_NGMRES;
310:   ksp->ops->setfromoptions       = KSPSetFromOptions_NGMRES;
311:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
312:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
313:   return(0);
314: }




320: /*
321:     BuildNGmresSoln - create the solution of the least square problem to determine the coefficients

323:     Input parameters:
324:         nrs - the solution
325:         Fold  - the RHS vector
326:         flag - =1, we need to replace the old vector by new one, =0, we still add new vector 
327:      This is an internal routine that knows about the NGMRES internals.
328:  */
331: static PetscErrorCode BuildNGmresSoln(PetscScalar* nrs, Vec Fold, KSP ksp,PetscInt it, PetscInt flag)
332: {
333:   PetscScalar    tt,temps;
335:   PetscInt       i,ii,j,l;
336:   KSP_NGMRES      *ngmres = (KSP_NGMRES *)(ksp->data);
337:   //Vec *dF=ngmres->w+ngmres->msize, *Q=ngmres->w+ngmres->msize*2,temp;
338:   Vec *dF=ngmres->w+ngmres->msize,*Q=ngmres->q,temp;
339:   PetscReal      gam,areal;
340:   PetscScalar    a,b,c,s;
341: 
343:   VecDuplicate(Fold,&temp);
344:   l=it+1;
345: 
346:   /* Solve for solution vector that minimizes the residual */

348:   if(flag==1) { // we need to replace the old vector and need to modify the QR factors, use Givens rotation
349:       for(i=0;i<it;i++){
350:         /* calculate the Givens rotation */
351:         a=*HH(i,i);
352:         b=*HH(i+1,i);
353:         gam=1.0/PetscRealPart(PetscSqrtScalar(PetscConj(a)*a + PetscConj(b)*b));
354:         c= a*gam;
355:         s= b*gam;
356: 
357: #if defined(PETSC_USE_COMPLEX)
358:         /* update the Q factor */
359:         ierr= VecCopy(Q[i],temp);
360:         VecAXPBY(temp,s,PetscConj(c),Q[i+1]); /*temp= c*Q[i]+s*Q[i+1] */
361:         VecAXPBY(Q[i+1],-s,c,Q[i]); /* Q[i+1]= -s*Q[i] + c*Q[i+1] */
362:         ierr= VecCopy(temp,Q[i]);    /* Q[i]= c*Q[i] + s*Q[i+1] */
363:         /* update the R factor */
364:         for(j=0;j<l;j++){
365:           a= *HH(i,j);
366:           b=*HH(i+1,j);
367:           temps=PetscConj(c)* a+s* b;
368:           *HH(i+1,j)=-s*a+c*b;
369:           *HH(i,j)=temps;
370:         }
371: #else
372:         /* update the Q factor */
373:         ierr= VecCopy(Q[i],temp);
374:         VecAXPBY(temp,s,c,Q[i+1]); /*temp= c*Q[i]+s*Q[i+1] */
375:         VecAXPBY(Q[i+1],-s,c,Q[i]); /* Q[i+1]= -s*Q[i] + c*Q[i+1] */
376:         ierr= VecCopy(temp,Q[i]);    /* Q[i]= c*Q[i] + s*Q[i+1] */
377:         /* update the R factor */
378:         for(j=0;j<l;j++){
379:           a= *HH(i,j);
380:           b=*HH(i+1,j);
381:           temps=c* a+s* b;
382:           *HH(i+1,j)=-s*a+c*b;
383:           *HH(i,j)=temps;
384:         }
385: #endif 
386:       }
387:     }

389:     // add a new vector, use modified Gram-Schmidt
390:     ierr= VecCopy(dF[it],temp);
391:     for(i=0;i<it;i++){
392:       ierr=VecDot(temp,Q[i],HH(i,it)); /* h(i,l-1)= dF[l-1]'*Q[i] */
393:       VecAXPBY(temp,-*HH(i,it),1.0,Q[i]); /* temp= temp- h(i,l-1)*Q[i] */
394:     }
395:     ierr=VecCopy(temp,Q[it]);
396:     ierr=VecNormalize(Q[it],&areal);
397:     *HH(it,it) = areal;
398: 


401:     /* modify the RHS with Q'*Fold*/
402: 
403:   for(i=0;i<l;i++)
404:           ierr=VecDot(Fold,Q[i],ngmres->nrs+i); /* nrs= Fold'*Q[i] */
405: 
406:     /* start backsubstitution to solve the least square problem */

408:      if (*HH(it,it) != 0.0) {
409:       nrs[it] =  nrs[it]/ *HH(it,it);
410:     } else {
411:     ksp->reason = KSP_DIVERGED_BREAKDOWN;
412:     PetscInfo2(ksp,"Likely your matrix or preconditioner is singular. HH(it,it) is identically zero; it = %D nrs(it) = %G",it,10);
413:     return(0);
414:   }
415:   for (ii=1; ii<=it; ii++) {
416:     i   = it - ii;
417:     tt  = nrs[i];
418:     for (j=i+1; j<=it; j++) tt  = tt - *HH(i,j) * nrs[j];
419:     if (*HH(i,i) == 0.0) {
420:       ksp->reason = KSP_DIVERGED_BREAKDOWN;
421:       PetscInfo1(ksp,"Likely your matrix or preconditioner is singular. HH(k,k) is identically zero; i = %D",i);
422:       return(0);
423:     }
424:     nrs[i]   = tt / *HH(i,i);
425:   }

427: 
428:   return(0);
429: }