Actual source code: trlanczos.c

  1: /*                       

  3:    SLEPc singular value solver: "trlanczos"

  5:    Method: Golub-Kahan-Lanczos bidiagonalization with thick-restart

  7:    Last update: Jun 2007

  9:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 10:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 11:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

 13:    This file is part of SLEPc.
 14:       
 15:    SLEPc is free software: you can redistribute it and/or modify it under  the
 16:    terms of version 3 of the GNU Lesser General Public License as published by
 17:    the Free Software Foundation.

 19:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 20:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 21:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 22:    more details.

 24:    You  should have received a copy of the GNU Lesser General  Public  License
 25:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 26:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 27: */

 29: #include <private/svdimpl.h>                /*I "slepcsvd.h" I*/
 30: #include <private/ipimpl.h>
 31: #include <slepcblaslapack.h>

 33: typedef struct {
 34:   PetscBool oneside;
 35: } SVD_TRLANCZOS;

 39: PetscErrorCode SVDSetUp_TRLanczos(SVD svd)
 40: {
 42:   PetscInt       N;

 45:   SVDMatGetSize(svd,PETSC_NULL,&N);
 46:   if (svd->ncv) { /* ncv set */
 47:     if (svd->ncv<svd->nsv) SETERRQ(((PetscObject)svd)->comm,1,"The value of ncv must be at least nsv");
 48:   }
 49:   else if (svd->mpd) { /* mpd set */
 50:     svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
 51:   }
 52:   else { /* neither set: defaults depend on nsv being small or large */
 53:     if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
 54:     else { svd->mpd = 500; svd->ncv = PetscMin(N,svd->nsv+svd->mpd); }
 55:   }
 56:   if (!svd->mpd) svd->mpd = svd->ncv;
 57:   if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(((PetscObject)svd)->comm,1,"The value of ncv must not be larger than nev+mpd");
 58:   if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
 59:   if (svd->ncv!=svd->n) {
 60:     VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
 61:   }
 62:   return(0);
 63: }

 67: static PetscErrorCode SVDOneSideTRLanczosMGS(SVD svd,PetscReal *alpha,PetscReal *beta,PetscScalar* bb,Vec *V,Vec v,Vec* U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
 68: {
 70:   PetscReal      a,b;
 71:   PetscInt       i,k=nconv+l;

 74:   SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
 75:   if (l>0) {
 76:     SlepcVecMAXPBY(U[k],1.0,-1.0,l,bb,U+nconv);
 77:   }
 78:   IPNorm(svd->ip,U[k],&a);
 79:   VecScale(U[k],1.0/a);
 80:   alpha[0] = a;
 81:   for (i=k+1;i<n;i++) {
 82:     SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
 83:     IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,&b,PETSC_NULL);
 84:     VecScale(V[i],1.0/b);
 85:     beta[i-k-1] = b;
 86: 
 87:     SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
 88:     VecAXPY(U[i],-b,U[i-1]);
 89:     IPNorm(svd->ip,U[i],&a);
 90:     VecScale(U[i],1.0/a);
 91:     alpha[i-k] = a;
 92:   }
 93:   SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
 94:   IPOrthogonalize(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,&b,PETSC_NULL);
 95:   beta[n-k-1] = b;
 96:   return(0);
 97: }

101: static PetscErrorCode SVDOneSideTRLanczosCGS(SVD svd,PetscReal *alpha,PetscReal *beta,PetscScalar* bb,Vec *V,Vec v,Vec* U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
102: {
104:   PetscReal      a,b,sum,onorm;
105:   PetscScalar    dot;
106:   PetscInt       i,j,k=nconv+l;

109:   SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
110:   if (l>0) {
111:     SlepcVecMAXPBY(U[k],1.0,-1.0,l,bb,U+nconv);
112:   }
113:   for (i=k+1;i<n;i++) {
114:     SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
115:     IPNormBegin(svd->ip,U[i-1],&a);
116:     if (svd->ip->orthog_ref == IP_ORTHOG_REFINE_IFNEEDED) {
117:       IPInnerProductBegin(svd->ip,V[i],V[i],&dot);
118:     }
119:     IPMInnerProductBegin(svd->ip,V[i],i,V,work);
120:     IPNormEnd(svd->ip,U[i-1],&a);
121:     if (svd->ip->orthog_ref == IP_ORTHOG_REFINE_IFNEEDED) {
122:       IPInnerProductEnd(svd->ip,V[i],V[i],&dot);
123:     }
124:     IPMInnerProductEnd(svd->ip,V[i],i,V,work);
125: 
126:     VecScale(U[i-1],1.0/a);
127:     for (j=0;j<i;j++) work[j] = work[j] / a;
128:     SlepcVecMAXPBY(V[i],1.0/a,-1.0,i,work,V);

130:     switch (svd->ip->orthog_ref) {
131:     case IP_ORTHOG_REFINE_NEVER:
132:       IPNorm(svd->ip,V[i],&b);
133:       break;
134:     case IP_ORTHOG_REFINE_ALWAYS:
135:       IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b);
136:       break;
137:     case IP_ORTHOG_REFINE_IFNEEDED:
138:       onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
139:       sum = 0.0;
140:       for (j=0;j<i;j++) {
141:         sum += PetscRealPart(work[j] * PetscConj(work[j]));
142:       }
143:       b = PetscRealPart(dot)/(a*a) - sum;
144:       if (b>0.0) b = PetscSqrtReal(b);
145:       else {
146:         IPNorm(svd->ip,V[i],&b);
147:       }
148:       if (b < svd->ip->orthog_eta * onorm) {
149:         IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b);
150:       }
151:       break;
152:     }
153: 
154:     VecScale(V[i],1.0/b);
155: 
156:     SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
157:     VecAXPY(U[i],-b,U[i-1]);

159:     alpha[i-k-1] = a;
160:     beta[i-k-1] = b;
161:   }
162:   SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
163:   IPNormBegin(svd->ip,U[n-1],&a);
164:   if (svd->ip->orthog_ref == IP_ORTHOG_REFINE_IFNEEDED) {
165:     IPInnerProductBegin(svd->ip,v,v,&dot);
166:   }
167:   IPMInnerProductBegin(svd->ip,v,n,V,work);
168:   IPNormEnd(svd->ip,U[n-1],&a);
169:   if (svd->ip->orthog_ref == IP_ORTHOG_REFINE_IFNEEDED) {
170:     IPInnerProductEnd(svd->ip,v,v,&dot);
171:   }
172:   IPMInnerProductEnd(svd->ip,v,n,V,work);
173: 
174:   VecScale(U[n-1],1.0/a);
175:   for (j=0;j<n;j++) work[j] = work[j] / a;
176:   SlepcVecMAXPBY(v,1.0/a,-1.0,n,work,V);

178:   switch (svd->ip->orthog_ref) {
179:   case IP_ORTHOG_REFINE_NEVER:
180:     IPNorm(svd->ip,v,&b);
181:     break;
182:   case IP_ORTHOG_REFINE_ALWAYS:
183:     IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,PETSC_NULL,&b);
184:     break;
185:   case IP_ORTHOG_REFINE_IFNEEDED:
186:     onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
187:     sum = 0.0;
188:     for (j=0;j<i;j++) {
189:       sum += PetscRealPart(work[j] * PetscConj(work[j]));
190:     }
191:     b = PetscRealPart(dot)/(a*a) - sum;
192:     if (b>0.0) b = PetscSqrtReal(b);
193:     else {
194:       IPNorm(svd->ip,v,&b);
195:     }
196:     if (b < svd->ip->orthog_eta * onorm) {
197:       IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,PETSC_NULL,&b);
198:     }
199:     break;
200:   }
201: 
202:   alpha[n-k-1] = a;
203:   beta[n-k-1] = b;
204:   return(0);
205: }

209: PetscErrorCode SVDSolve_TRLanczos(SVD svd)
210: {
212:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS *)svd->data;
213:   PetscReal      *alpha,*beta,norm;
214:   PetscScalar    *b,*Q,*PT,*swork;
215:   PetscInt       i,j,k,l,m,n,nv;
216:   Vec            v;
217:   PetscBool      conv;
218:   IPOrthogType   orthog;
219: 
221:   /* allocate working space */
222:   PetscMalloc(sizeof(PetscReal)*svd->n,&alpha);
223:   PetscMalloc(sizeof(PetscReal)*svd->n,&beta);
224:   PetscMalloc(sizeof(PetscScalar)*svd->n,&b);
225:   PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&Q);
226:   PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&PT);
227: #if !defined(PETSC_USE_COMPLEX)
228:   if (svd->which == SVD_SMALLEST) {
229: #endif
230:     PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&swork);
231: #if !defined(PETSC_USE_COMPLEX)
232:   } else {
233:     PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
234:   }
235: #endif
236:   VecDuplicate(svd->V[0],&v);
237:   IPGetOrthogonalization(svd->ip,&orthog,PETSC_NULL,PETSC_NULL);
238: 
239:   /* normalize start vector */
240:   if (svd->nini==0) {
241:     SlepcVecSetRandom(svd->V[0],svd->rand);
242:   }
243:   VecNormalize(svd->V[0],&norm);
244: 
245:   l = 0;
246:   while (svd->reason == SVD_CONVERGED_ITERATING) {
247:     svd->its++;

249:     /* inner loop */
250:     nv = PetscMin(svd->nconv+svd->mpd,svd->n);
251:     if (lanczos->oneside) {
252:       if (orthog == IP_ORTHOG_MGS) {
253:         SVDOneSideTRLanczosMGS(svd,alpha,beta,b+svd->nconv,svd->V,v,svd->U,svd->nconv,l,nv,swork);
254:       } else {
255:         SVDOneSideTRLanczosCGS(svd,alpha,beta,b+svd->nconv,svd->V,v,svd->U,svd->nconv,l,nv,swork);
256:       }
257:     } else {
258:       SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv+l,nv,swork);
259:     }
260:     VecScale(v,1.0/beta[nv-svd->nconv-l-1]);
261: 
262:     /* compute SVD of general matrix */
263:     n = nv - svd->nconv;
264:     /* first l columns */
265:     for (j=0;j<l;j++) {
266:       for (i=0;i<j;i++) Q[j*n+i] = 0.0;
267:       Q[j*n+j] = svd->sigma[svd->nconv+j];
268:       for (i=j+1;i<n;i++) Q[j*n+i] = 0.0;
269:     }
270:     /* l+1 column */
271:     for (i=0;i<l;i++) Q[l*n+i] = b[i+svd->nconv];
272:     Q[l*n+l] = alpha[0];
273:     for (i=l+1;i<n;i++) Q[l*n+i] = 0.0;
274:     /* rest of matrix */
275:     for (j=l+1;j<n;j++) {
276:       for (i=0;i<j-1;i++) Q[j*n+i] = 0.0;
277:       Q[j*n+j-1] = beta[j-l-1];
278:       Q[j*n+j] = alpha[j-l];
279:       for (i=j+1;i<n;i++) Q[j*n+i] = 0.0;
280:     }
281:     SVDDense(n,n,Q,alpha,PETSC_NULL,PT);

283:     /* compute error estimates */
284:     k = 0;
285:     conv = PETSC_TRUE;
286:     for (i=svd->nconv;i<nv;i++) {
287:       if (svd->which == SVD_SMALLEST) j = n-i+svd->nconv-1;
288:       else j = i-svd->nconv;
289:       svd->sigma[i] = alpha[j];
290:       b[i] = Q[j*n+n-1]*beta[n-l-1];
291:       svd->errest[i] = PetscAbsScalar(b[i]);
292:       if (alpha[j] > svd->tol) svd->errest[i] /= alpha[j];
293:       if (conv) {
294:         if (svd->errest[i] < svd->tol) k++;
295:         else conv = PETSC_FALSE;
296:       }
297:     }
298: 
299:     /* check convergence and update l */
300:     if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
301:     if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
302:     if (svd->reason != SVD_CONVERGED_ITERATING) l = 0;
303:     else l = PetscMax((nv - svd->nconv - k) / 2,0);
304: 
305:     /* compute converged singular vectors and restart vectors*/
306: #if !defined(PETSC_USE_COMPLEX)
307:     if (svd->which == SVD_SMALLEST) {
308: #endif
309:     for (i=0;i<k+l;i++) {
310:       if (svd->which == SVD_SMALLEST) j = n-i-1;
311:       else j = i;
312:       for (m=0;m<n;m++) swork[j*n+m] = PT[m*n+j];
313:     }
314:     SlepcUpdateVectors(n,svd->V+svd->nconv,0,k+l,swork,n,PETSC_FALSE);
315:     for (i=0;i<k+l;i++) {
316:       if (svd->which == SVD_SMALLEST) j = n-i-1;
317:       else j = i;
318:       for (m=0;m<n;m++) swork[j*n+m] = Q[j*n+m];
319:     }
320:     SlepcUpdateVectors(n,svd->U+svd->nconv,0,k+l,swork,n,PETSC_FALSE);
321: #if !defined(PETSC_USE_COMPLEX)
322:     } else {
323:       SlepcUpdateVectors(n,svd->V+svd->nconv,0,k+l,PT,n,PETSC_TRUE);
324:       SlepcUpdateVectors(n,svd->U+svd->nconv,0,k+l,Q,n,PETSC_FALSE);
325:     }
326: #endif
327: 
328:     /* copy the last vector to be the next initial vector */
329:     if (svd->reason == SVD_CONVERGED_ITERATING) {
330:       VecCopy(v,svd->V[svd->nconv+k+l]);
331:     }
332: 
333:     svd->nconv += k;
334:     SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
335:   }
336: 
337:   /* orthonormalize U columns in one side method */
338:   if (lanczos->oneside) {
339:     for (i=0;i<svd->nconv;i++) {
340:       IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,svd->U,svd->U[i],PETSC_NULL,&norm,PETSC_NULL);
341:       VecScale(svd->U[i],1.0/norm);
342:     }
343:   }
344: 
345:   /* free working space */
346:   VecDestroy(&v);

348:   PetscFree(alpha);
349:   PetscFree(beta);
350:   PetscFree(b);
351:   PetscFree(Q);
352:   PetscFree(PT);
353:   PetscFree(swork);
354:   return(0);
355: }

359: PetscErrorCode SVDSetFromOptions_TRLanczos(SVD svd)
360: {
362:   PetscBool      set,val;
363:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS *)svd->data;

366:   PetscOptionsHead("SVD TRLanczos Options");
367:   PetscOptionsBool("-svd_trlanczos_oneside","Lanczos one-side reorthogonalization","SVDTRLanczosSetOneSide",lanczos->oneside,&val,&set);
368:   if (set) {
369:     SVDTRLanczosSetOneSide(svd,val);
370:   }
371:   PetscOptionsTail();
372:   return(0);
373: }

375: EXTERN_C_BEGIN
378: PetscErrorCode SVDTRLanczosSetOneSide_TRLanczos(SVD svd,PetscBool oneside)
379: {
380:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;

383:   lanczos->oneside = oneside;
384:   return(0);
385: }
386: EXTERN_C_END

390: /*@
391:    SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method 
392:    to be used is one-sided or two-sided.

394:    Logically Collective on SVD

396:    Input Parameters:
397: +  svd     - singular value solver
398: -  oneside - boolean flag indicating if the method is one-sided or not

400:    Options Database Key:
401: .  -svd_trlanczos_oneside <boolean> - Indicates the boolean flag

403:    Note:
404:    By default, a two-sided variant is selected, which is sometimes slightly
405:    more robust. However, the one-sided variant is faster because it avoids 
406:    the orthogonalization associated to left singular vectors. 

408:    Level: advanced

410: .seealso: SVDLanczosSetOneSide()
411: @*/
412: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscBool oneside)
413: {

419:   PetscTryMethod(svd,"SVDTRLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
420:   return(0);
421: }

425: /*@
426:    SVDTRLanczosGetOneSide - Gets if the variant of the Lanczos method 
427:    to be used is one-sided or two-sided.

429:    Not Collective

431:    Input Parameters:
432: .  svd     - singular value solver

434:    Output Parameters:
435: .  oneside - boolean flag indicating if the method is one-sided or not

437:    Level: advanced

439: .seealso: SVDTRLanczosSetOneSide()
440: @*/
441: PetscErrorCode SVDTRLanczosGetOneSide(SVD svd,PetscBool *oneside)
442: {

448:   PetscTryMethod(svd,"SVDTRLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
449:   return(0);
450: }

452: EXTERN_C_BEGIN
455: PetscErrorCode SVDTRLanczosGetOneSide_TRLanczos(SVD svd,PetscBool *oneside)
456: {
457:   SVD_TRLANCZOS    *lanczos = (SVD_TRLANCZOS *)svd->data;

460:   *oneside = lanczos->oneside;
461:   return(0);
462: }
463: EXTERN_C_END

467: PetscErrorCode SVDReset_TRLanczos(SVD svd)
468: {

472:   VecDestroyVecs(svd->n,&svd->U);
473:   return(0);
474: }

478: PetscErrorCode SVDDestroy_TRLanczos(SVD svd)
479: {

483:   PetscFree(svd->data);
484:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosSetOneSide_C","",PETSC_NULL);
485:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosGetOneSide_C","",PETSC_NULL);
486:   return(0);
487: }

491: PetscErrorCode SVDView_TRLanczos(SVD svd,PetscViewer viewer)
492: {
494:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS *)svd->data;

497:   PetscViewerASCIIPrintf(viewer,"  TRLanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
498:   return(0);
499: }

501: EXTERN_C_BEGIN
504: PetscErrorCode SVDCreate_TRLanczos(SVD svd)
505: {

509:   PetscNewLog(svd,SVD_TRLANCZOS,&svd->data);
510:   svd->ops->setup          = SVDSetUp_TRLanczos;
511:   svd->ops->solve          = SVDSolve_TRLanczos;
512:   svd->ops->destroy        = SVDDestroy_TRLanczos;
513:   svd->ops->reset          = SVDReset_TRLanczos;
514:   svd->ops->setfromoptions = SVDSetFromOptions_TRLanczos;
515:   svd->ops->view           = SVDView_TRLanczos;
516:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosSetOneSide_C","SVDTRLanczosSetOneSide_TRLanczos",SVDTRLanczosSetOneSide_TRLanczos);
517:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosGetOneSide_C","SVDTRLanczosGetOneSide_TRLanczos",SVDTRLanczosGetOneSide_TRLanczos);
518:   if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
519:   return(0);
520: }
521: EXTERN_C_END