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-2010, Universidad 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
 30:  #include private/ipimpl.h
 31:  #include slepcblaslapack.h

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

 39: PetscErrorCode SVDSetUp_TRLANCZOS(SVD svd)
 40: {
 41:   PetscErrorCode  ierr;
 42:   PetscInt        i,N,nloc;
 43:   PetscScalar     *pU;

 46:   SVDMatGetSize(svd,PETSC_NULL,&N);
 47:   if (svd->ncv) { /* ncv set */
 48:     if (svd->ncv<svd->nsv) SETERRQ(1,"The value of ncv must be at least nsv");
 49:   }
 50:   else if (svd->mpd) { /* mpd set */
 51:     svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
 52:   }
 53:   else { /* neither set: defaults depend on nsv being small or large */
 54:     if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
 55:     else { svd->mpd = 500; svd->ncv = PetscMin(N,svd->nsv+svd->mpd); }
 56:   }
 57:   if (!svd->mpd) svd->mpd = svd->ncv;
 58:   if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(1,"The value of ncv must not be larger than nev+mpd");
 59:   if (!svd->max_it)
 60:     svd->max_it = PetscMax(N/svd->ncv,100);
 61:   if (svd->ncv!=svd->n) {
 62:     if (svd->U) {
 63:       VecGetArray(svd->U[0],&pU);
 64:       for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]);  }
 65:       PetscFree(pU);
 66:       PetscFree(svd->U);
 67:     }
 68:     PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
 69:     SVDMatGetLocalSize(svd,&nloc,PETSC_NULL);
 70:     PetscMalloc(svd->ncv*nloc*sizeof(PetscScalar),&pU);
 71:     for (i=0;i<svd->ncv;i++) {
 72:       VecCreateMPIWithArray(((PetscObject)svd)->comm,nloc,PETSC_DECIDE,pU+i*nloc,&svd->U[i]);
 73:     }
 74:   }
 75:   return(0);
 76: }

 80: 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)
 81: {
 83:   PetscReal      a,b;
 84:   PetscInt       i,k=nconv+l;

 87:   SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
 88:   if (l>0) {
 89:     SlepcVecMAXPBY(U[k],1.0,-1.0,l,bb,U+nconv);
 90:   }
 91:   IPNorm(svd->ip,U[k],&a);
 92:   VecScale(U[k],1.0/a);
 93:   alpha[0] = a;
 94:   for (i=k+1;i<n;i++) {
 95:     SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
 96:     IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,&b,PETSC_NULL);
 97:     VecScale(V[i],1.0/b);
 98:     beta[i-k-1] = b;
 99: 
100:     SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
101:     VecAXPY(U[i],-b,U[i-1]);
102:     IPNorm(svd->ip,U[i],&a);
103:     VecScale(U[i],1.0/a);
104:     alpha[i-k] = a;
105:   }
106:   SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
107:   IPOrthogonalize(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,&b,PETSC_NULL);
108:   beta[n-k-1] = b;
109:   return(0);
110: }

114: 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)
115: {
117:   PetscReal      a,b,sum,onorm;
118:   PetscScalar    dot;
119:   PetscInt       i,j,k=nconv+l;

122:   SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
123:   if (l>0) {
124:     SlepcVecMAXPBY(U[k],1.0,-1.0,l,bb,U+nconv);
125:   }
126:   for (i=k+1;i<n;i++) {
127:     SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
128:     IPNormBegin(svd->ip,U[i-1],&a);
129:     if (svd->ip->orthog_ref == IP_ORTH_REFINE_IFNEEDED) {
130:       IPInnerProductBegin(svd->ip,V[i],V[i],&dot);
131:     }
132:     IPMInnerProductBegin(svd->ip,V[i],i,V,work);
133:     IPNormEnd(svd->ip,U[i-1],&a);
134:     if (svd->ip->orthog_ref == IP_ORTH_REFINE_IFNEEDED) {
135:       IPInnerProductEnd(svd->ip,V[i],V[i],&dot);
136:     }
137:     IPMInnerProductEnd(svd->ip,V[i],i,V,work);
138: 
139:     VecScale(U[i-1],1.0/a);
140:     for (j=0;j<i;j++) work[j] = work[j] / a;
141:     SlepcVecMAXPBY(V[i],1.0/a,-1.0,i,work,V);

143:     switch (svd->ip->orthog_ref) {
144:     case IP_ORTH_REFINE_NEVER:
145:       IPNorm(svd->ip,V[i],&b);
146:       break;
147:     case IP_ORTH_REFINE_ALWAYS:
148:       IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b);
149:       break;
150:     case IP_ORTH_REFINE_IFNEEDED:
151:       onorm = sqrt(PetscRealPart(dot)) / a;
152:       sum = 0.0;
153:       for (j=0;j<i;j++) {
154:         sum += PetscRealPart(work[j] * PetscConj(work[j]));
155:       }
156:       b = PetscRealPart(dot)/(a*a) - sum;
157:       if (b>0.0) b = sqrt(b);
158:       else {
159:         IPNorm(svd->ip,V[i],&b);
160:       }
161:       if (b < svd->ip->orthog_eta * onorm) {
162:         IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b);
163:       }
164:       break;
165:     }
166: 
167:     VecScale(V[i],1.0/b);
168: 
169:     SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
170:     VecAXPY(U[i],-b,U[i-1]);

172:     alpha[i-k-1] = a;
173:     beta[i-k-1] = b;
174:   }
175:   SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
176:   IPNormBegin(svd->ip,U[n-1],&a);
177:   if (svd->ip->orthog_ref == IP_ORTH_REFINE_IFNEEDED) {
178:     IPInnerProductBegin(svd->ip,v,v,&dot);
179:   }
180:   IPMInnerProductBegin(svd->ip,v,n,V,work);
181:   IPNormEnd(svd->ip,U[n-1],&a);
182:   if (svd->ip->orthog_ref == IP_ORTH_REFINE_IFNEEDED) {
183:     IPInnerProductEnd(svd->ip,v,v,&dot);
184:   }
185:   IPMInnerProductEnd(svd->ip,v,n,V,work);
186: 
187:   VecScale(U[n-1],1.0/a);
188:   for (j=0;j<n;j++) work[j] = work[j] / a;
189:   SlepcVecMAXPBY(v,1.0/a,-1.0,n,work,V);

191:   switch (svd->ip->orthog_ref) {
192:   case IP_ORTH_REFINE_NEVER:
193:     IPNorm(svd->ip,v,&b);
194:     break;
195:   case IP_ORTH_REFINE_ALWAYS:
196:     IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,PETSC_NULL,&b);
197:     break;
198:   case IP_ORTH_REFINE_IFNEEDED:
199:     onorm = sqrt(PetscRealPart(dot)) / a;
200:     sum = 0.0;
201:     for (j=0;j<i;j++) {
202:       sum += PetscRealPart(work[j] * PetscConj(work[j]));
203:     }
204:     b = PetscRealPart(dot)/(a*a) - sum;
205:     if (b>0.0) b = sqrt(b);
206:     else {
207:       IPNorm(svd->ip,v,&b);
208:     }
209:     if (b < svd->ip->orthog_eta * onorm) {
210:       IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,PETSC_NULL,&b);
211:     }
212:     break;
213:   }
214: 
215:   alpha[n-k-1] = a;
216:   beta[n-k-1] = b;
217:   return(0);
218: }

222: PetscErrorCode SVDSolve_TRLANCZOS(SVD svd)
223: {
225:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS *)svd->data;
226:   PetscReal      *alpha,*beta,norm;
227:   PetscScalar    *b,*Q,*PT,*swork;
228:   PetscInt       i,j,k,l,m,n,nv;
229:   Vec            v;
230:   PetscTruth     conv;
231:   IPOrthogonalizationType orthog;
232: 
234:   /* allocate working space */
235:   PetscMalloc(sizeof(PetscReal)*svd->n,&alpha);
236:   PetscMalloc(sizeof(PetscReal)*svd->n,&beta);
237:   PetscMalloc(sizeof(PetscScalar)*svd->n,&b);
238:   PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&Q);
239:   PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&PT);
240: #if !defined(PETSC_USE_COMPLEX)
241:   if (svd->which == SVD_SMALLEST) {
242: #endif
243:     PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&swork);
244: #if !defined(PETSC_USE_COMPLEX)
245:   } else {
246:     PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
247:   }
248: #endif
249:   VecDuplicate(svd->V[0],&v);
250:   IPGetOrthogonalization(svd->ip,&orthog,PETSC_NULL,PETSC_NULL);
251: 
252:   /* normalize start vector */
253:   if (svd->nini==0) {
254:     SlepcVecSetRandom(svd->V[0],svd->rand);
255:   }
256:   VecNormalize(svd->V[0],&norm);
257: 
258:   l = 0;
259:   while (svd->reason == SVD_CONVERGED_ITERATING) {
260:     svd->its++;

262:     /* inner loop */
263:     nv = PetscMin(svd->nconv+svd->mpd,svd->n);
264:     if (lanczos->oneside) {
265:       if (orthog == IP_ORTH_MGS) {
266:         SVDOneSideTRLanczosMGS(svd,alpha,beta,b+svd->nconv,svd->V,v,svd->U,svd->nconv,l,nv,swork);
267:       } else {
268:         SVDOneSideTRLanczosCGS(svd,alpha,beta,b+svd->nconv,svd->V,v,svd->U,svd->nconv,l,nv,swork);
269:       }
270:     } else {
271:       SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv+l,nv,swork);
272:     }
273:     VecScale(v,1.0/beta[nv-svd->nconv-l-1]);
274: 
275:     /* compute SVD of general matrix */
276:     n = nv - svd->nconv;
277:     /* first l columns */
278:     for (j=0;j<l;j++) {
279:       for (i=0;i<j;i++) Q[j*n+i] = 0.0;
280:       Q[j*n+j] = svd->sigma[svd->nconv+j];
281:       for (i=j+1;i<n;i++) Q[j*n+i] = 0.0;
282:     }
283:     /* l+1 column */
284:     for (i=0;i<l;i++) Q[l*n+i] = b[i+svd->nconv];
285:     Q[l*n+l] = alpha[0];
286:     for (i=l+1;i<n;i++) Q[l*n+i] = 0.0;
287:     /* rest of matrix */
288:     for (j=l+1;j<n;j++) {
289:       for (i=0;i<j-1;i++) Q[j*n+i] = 0.0;
290:       Q[j*n+j-1] = beta[j-l-1];
291:       Q[j*n+j] = alpha[j-l];
292:       for (i=j+1;i<n;i++) Q[j*n+i] = 0.0;
293:     }
294:     SVDDense(n,n,Q,alpha,PETSC_NULL,PT);

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

361:   PetscFree(alpha);
362:   PetscFree(beta);
363:   PetscFree(b);
364:   PetscFree(Q);
365:   PetscFree(PT);
366:   PetscFree(swork);
367:   return(0);
368: }

372: PetscErrorCode SVDSetFromOptions_TRLANCZOS(SVD svd)
373: {
375:   PetscTruth     set,val;
376:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS *)svd->data;

379:   PetscOptionsBegin(((PetscObject)svd)->comm,((PetscObject)svd)->prefix,"TRLANCZOS Singular Value Solver Options","SVD");
380:   PetscOptionsTruth("-svd_trlanczos_oneside","Lanczos one-side reorthogonalization","SVDTRLanczosSetOneSide",lanczos->oneside,&val,&set);
381:   if (set) {
382:     SVDTRLanczosSetOneSide(svd,val);
383:   }
384:   PetscOptionsEnd();
385:   return(0);
386: }

391: PetscErrorCode SVDTRLanczosSetOneSide_TRLANCZOS(SVD svd,PetscTruth oneside)
392: {
393:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;

396:   lanczos->oneside = oneside;
397:   return(0);
398: }

403: /*@
404:    SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method 
405:    to be used is one-sided or two-sided.

407:    Collective on SVD

409:    Input Parameters:
410: +  svd     - singular value solver
411: -  oneside - boolean flag indicating if the method is one-sided or not

413:    Options Database Key:
414: .  -svd_trlanczos_oneside <boolean> - Indicates the boolean flag

416:    Note:
417:    By default, a two-sided variant is selected, which is sometimes slightly
418:    more robust. However, the one-sided variant is faster because it avoids 
419:    the orthogonalization associated to left singular vectors. 

421:    Level: advanced

423: .seealso: SVDLanczosSetOneSide()
424: @*/
425: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscTruth oneside)
426: {
427:   PetscErrorCode ierr, (*f)(SVD,PetscTruth);

431:   PetscObjectQueryFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",(void (**)())&f);
432:   if (f) {
433:     (*f)(svd,oneside);
434:   }
435:   return(0);
436: }

440: /*@
441:    SVDTRLanczosGetOneSide - Gets if the variant of the Lanczos method 
442:    to be used is one-sided or two-sided.

444:    Collective on SVD

446:    Input Parameters:
447: .  svd     - singular value solver

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

452:    Level: advanced

454: .seealso: SVDTRLanczosSetOneSide()
455: @*/
456: PetscErrorCode SVDTRLanczosGetOneSide(SVD svd,PetscTruth *oneside)
457: {
458:   PetscErrorCode ierr, (*f)(SVD,PetscTruth*);

462:   PetscObjectQueryFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",(void (**)())&f);
463:   if (f) {
464:     (*f)(svd,oneside);
465:   }
466:   return(0);
467: }

472: PetscErrorCode SVDTRLanczosGetOneSide_TRLANCZOS(SVD svd,PetscTruth *oneside)
473: {
474:   SVD_TRLANCZOS    *lanczos = (SVD_TRLANCZOS *)svd->data;

478:   *oneside = lanczos->oneside;
479:   return(0);
480: }


486: PetscErrorCode SVDDestroy_TRLANCZOS(SVD svd)
487: {

492:   SVDDestroy_Default(svd);
493:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosSetOneSide_C","",PETSC_NULL);
494:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosGetOneSide_C","",PETSC_NULL);
495:   return(0);
496: }

500: PetscErrorCode SVDView_TRLANCZOS(SVD svd,PetscViewer viewer)
501: {
503:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS *)svd->data;

506:   PetscViewerASCIIPrintf(viewer,"Lanczos reorthogonalization: %s\n",lanczos->oneside ? "one-side" : "two-side");
507:   return(0);
508: }

513: PetscErrorCode SVDCreate_TRLANCZOS(SVD svd)
514: {
516:   SVD_TRLANCZOS  *lanczos;

519:   PetscNew(SVD_TRLANCZOS,&lanczos);
520:   PetscLogObjectMemory(svd,sizeof(SVD_TRLANCZOS));
521:   svd->data                = (void *)lanczos;
522:   svd->ops->setup          = SVDSetUp_TRLANCZOS;
523:   svd->ops->solve          = SVDSolve_TRLANCZOS;
524:   svd->ops->destroy        = SVDDestroy_TRLANCZOS;
525:   svd->ops->setfromoptions = SVDSetFromOptions_TRLANCZOS;
526:   svd->ops->view           = SVDView_TRLANCZOS;
527:   lanczos->oneside         = PETSC_FALSE;
528:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosSetOneSide_C","SVDTRLanczosSetOneSide_TRLANCZOS",SVDTRLanczosSetOneSide_TRLANCZOS);
529:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosGetOneSide_C","SVDTRLanczosGetOneSide_TRLANCZOS",SVDTRLanczosGetOneSide_TRLANCZOS);
530:   return(0);
531: }