Actual source code: trlanczos.c

slepc-3.18.1 2022-11-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc singular value solver: "trlanczos"

 13:    Method: Thick-restart Lanczos

 15:    Algorithm:

 17:        Golub-Kahan-Lanczos bidiagonalization with thick-restart.

 19:    References:

 21:        [1] G.H. Golub and W. Kahan, "Calculating the singular values
 22:            and pseudo-inverse of a matrix", SIAM J. Numer. Anal. Ser.
 23:            B 2:205-224, 1965.

 25:        [2] V. Hernandez, J.E. Roman, and A. Tomas, "A robust and
 26:            efficient parallel SVD solver based on restarted Lanczos
 27:            bidiagonalization", Elec. Trans. Numer. Anal. 31:68-85,
 28:            2008.
 29: */

 31: #include <slepc/private/svdimpl.h>

 33: static PetscBool  cited = PETSC_FALSE,citedg = PETSC_FALSE;
 34: static const char citation[] =
 35:   "@Article{slepc-svd,\n"
 36:   "   author = \"V. Hern{\\'a}ndez and J. E. Rom{\\'a}n and A. Tom{\\'a}s\",\n"
 37:   "   title = \"A robust and efficient parallel {SVD} solver based on restarted {Lanczos} bidiagonalization\",\n"
 38:   "   journal = \"Electron. Trans. Numer. Anal.\",\n"
 39:   "   volume = \"31\",\n"
 40:   "   pages = \"68--85\",\n"
 41:   "   year = \"2008\"\n"
 42:   "}\n";
 43: static const char citationg[] =
 44:   "@Article{slepc-gsvd,\n"
 45:   "   author = \"F. Alvarruiz and C. Campos and J. E. Roman\",\n"
 46:   "   title = \"Thick-restarted {Lanczos} bidiagonalization methods for the {GSVD}\",\n"
 47:   "   note = \"arXiv:2206.03768\",\n"
 48:   "   year = \"2022\"\n"
 49:   "}\n";

 51: typedef struct {
 52:   /* user parameters */
 53:   PetscBool           oneside;   /* one-sided variant */
 54:   PetscReal           keep;      /* restart parameter */
 55:   PetscBool           lock;      /* locking/non-locking variant */
 56:   KSP                 ksp;       /* solver for least-squares problem in GSVD */
 57:   SVDTRLanczosGBidiag bidiag;    /* bidiagonalization variant for GSVD */
 58:   PetscReal           scalef;    /* scale factor for matrix B */
 59:   PetscReal           scaleth;   /* scale threshold for automatic scaling */
 60:   PetscBool           explicitmatrix;
 61:   /* auxiliary variables */
 62:   Mat                 Z;         /* aux matrix for GSVD, Z=[A;B] */
 63: } SVD_TRLANCZOS;

 65: /* Context for shell matrix [A; B] */
 66: typedef struct {
 67:   Mat       A,B,AT,BT;
 68:   Vec       y1,y2,y;
 69:   PetscInt  m;
 70:   PetscReal scalef;
 71: } MatZData;

 73: static PetscErrorCode MatZCreateContext(SVD svd,MatZData **zdata)
 74: {
 75:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;

 77:   PetscNew(zdata);
 78:   (*zdata)->A = svd->A;
 79:   (*zdata)->B = svd->B;
 80:   (*zdata)->AT = svd->AT;
 81:   (*zdata)->BT = svd->BT;
 82:   (*zdata)->scalef = lanczos->scalef;
 83:   MatCreateVecsEmpty(svd->A,NULL,&(*zdata)->y1);
 84:   MatCreateVecsEmpty(svd->B,NULL,&(*zdata)->y2);
 85:   VecGetLocalSize((*zdata)->y1,&(*zdata)->m);
 86:   BVCreateVec(svd->U,&(*zdata)->y);
 87:   return 0;
 88: }

 90: /* Update scale factor for B in Z=[A;B]
 91:    If matrices are swapped, the scale factor is inverted.*/
 92: static PetscErrorCode MatZUpdateScale(SVD svd)
 93: {
 94:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
 95:   MatZData       *zdata;
 96:   Mat            mats[2],normal;
 97:   MatType        Atype;
 98:   PetscBool      sametype;
 99:   PetscReal      scalef = svd->swapped? 1.0/lanczos->scalef : lanczos->scalef;

101:   if (lanczos->explicitmatrix) {
102:     /* Destroy the matrix Z and create it again */
103:     MatDestroy(&lanczos->Z);
104:     mats[0] = svd->A;
105:     if (scalef == 1.0) {
106:       mats[1] = svd->B;
107:     } else {
108:       MatDuplicate(svd->B,MAT_COPY_VALUES,&mats[1]);
109:       MatScale(mats[1],scalef);
110:     }
111:     MatCreateNest(PetscObjectComm((PetscObject)svd),2,NULL,1,NULL,mats,&lanczos->Z);
112:     MatGetType(svd->A,&Atype);
113:     PetscObjectTypeCompare((PetscObject)svd->B,Atype,&sametype);
114:     if (!sametype) Atype = MATAIJ;
115:     MatConvert(lanczos->Z,Atype,MAT_INPLACE_MATRIX,&lanczos->Z);
116:     if (scalef != 1.0) MatDestroy(&mats[1]);
117:   } else {
118:     MatShellGetContext(lanczos->Z,&zdata);
119:     zdata->scalef = scalef;
120:   }

122:   /* create normal equations matrix, to build the preconditioner in LSQR */
123:   MatCreateNormalHermitian(lanczos->Z,&normal);

125:   if (!lanczos->ksp) SVDTRLanczosGetKSP(svd,&lanczos->ksp);
126:   SVD_KSPSetOperators(lanczos->ksp,lanczos->Z,normal);
127:   KSPSetUp(lanczos->ksp);
128:   MatDestroy(&normal);
129:   return 0;
130: }

132: static PetscErrorCode MatDestroy_Z(Mat Z)
133: {
134:   MatZData       *zdata;

136:   MatShellGetContext(Z,&zdata);
137:   VecDestroy(&zdata->y1);
138:   VecDestroy(&zdata->y2);
139:   VecDestroy(&zdata->y);
140:   PetscFree(zdata);
141:   return 0;
142: }

144: static PetscErrorCode MatMult_Z(Mat Z,Vec x,Vec y)
145: {
146:   MatZData       *zdata;
147:   PetscScalar    *y_elems;

149:   MatShellGetContext(Z,&zdata);
150:   VecGetArray(y,&y_elems);
151:   VecPlaceArray(zdata->y1,y_elems);
152:   VecPlaceArray(zdata->y2,y_elems+zdata->m);

154:   MatMult(zdata->A,x,zdata->y1);
155:   MatMult(zdata->B,x,zdata->y2);
156:   VecScale(zdata->y2,zdata->scalef);

158:   VecResetArray(zdata->y1);
159:   VecResetArray(zdata->y2);
160:   VecRestoreArray(y,&y_elems);
161:   return 0;
162: }

164: static PetscErrorCode MatMultTranspose_Z(Mat Z,Vec y,Vec x)
165: {
166:   MatZData          *zdata;
167:   const PetscScalar *y_elems;

169:   MatShellGetContext(Z,&zdata);
170:   VecGetArrayRead(y,&y_elems);
171:   VecPlaceArray(zdata->y1,y_elems);
172:   VecPlaceArray(zdata->y2,y_elems+zdata->m);

174:   MatMult(zdata->BT,zdata->y2,x);
175:   VecScale(x,zdata->scalef);
176:   MatMultAdd(zdata->AT,zdata->y1,x,x);

178:   VecResetArray(zdata->y1);
179:   VecResetArray(zdata->y2);
180:   VecRestoreArrayRead(y,&y_elems);
181:   return 0;
182: }

184: static PetscErrorCode MatCreateVecs_Z(Mat Z,Vec *right,Vec *left)
185: {
186:   MatZData       *zdata;

188:   MatShellGetContext(Z,&zdata);
189:   if (right) MatCreateVecs(zdata->A,right,NULL);
190:   if (left) VecDuplicate(zdata->y,left);
191:   return 0;
192: }

194: #define SWAP(a,b,t) {t=a;a=b;b=t;}

196: PetscErrorCode SVDSetUp_TRLanczos(SVD svd)
197: {
198:   PetscInt       M,N,P,m,n,p;
199:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
200:   DSType         dstype;
201:   MatZData       *zdata;
202:   Mat            aux;

204:   SVDCheckDefinite(svd);
205:   MatGetSize(svd->A,&M,&N);
206:   SVDSetDimensions_Default(svd);
209:   if (svd->max_it==PETSC_DEFAULT) svd->max_it = PetscMax(N/svd->ncv,100);
210:   if (!lanczos->keep) lanczos->keep = 0.5;
211:   svd->leftbasis = PETSC_TRUE;
212:   SVDAllocateSolution(svd,1);
213:   dstype = DSSVD;
214:   if (svd->isgeneralized) {
215:     MatGetSize(svd->B,&P,NULL);
216:     if (lanczos->bidiag == SVD_TRLANCZOS_GBIDIAG_LOWER && ((svd->which==SVD_LARGEST && P<=N) || (svd->which==SVD_SMALLEST && M>N && P<=N))) {
217:       SWAP(svd->A,svd->B,aux);
218:       SWAP(svd->AT,svd->BT,aux);
219:       svd->swapped = PETSC_TRUE;
220:     } else svd->swapped = PETSC_FALSE;

222:     if (lanczos->bidiag==SVD_TRLANCZOS_GBIDIAG_UPPER || lanczos->bidiag==SVD_TRLANCZOS_GBIDIAG_LOWER) dstype = DSGSVD;
223:     SVDSetWorkVecs(svd,1,1);

225:     if (svd->conv==SVD_CONV_ABS) {  /* Residual norms are multiplied by matrix norms */
226:       if (!svd->nrma) MatNorm(svd->A,NORM_INFINITY,&svd->nrma);
227:       if (!svd->nrmb) MatNorm(svd->B,NORM_INFINITY,&svd->nrmb);
228:     }

230:     /* Create the matrix Z=[A;B] */
231:     MatGetLocalSize(svd->A,&m,&n);
232:     MatGetLocalSize(svd->B,&p,NULL);
233:     if (!lanczos->explicitmatrix) {
234:       MatDestroy(&lanczos->Z);
235:       MatZCreateContext(svd,&zdata);
236:       MatCreateShell(PetscObjectComm((PetscObject)svd),m+p,n,PETSC_DECIDE,PETSC_DECIDE,zdata,&lanczos->Z);
237:       MatShellSetOperation(lanczos->Z,MATOP_MULT,(void(*)(void))MatMult_Z);
238: #if defined(PETSC_USE_COMPLEX)
239:       MatShellSetOperation(lanczos->Z,MATOP_MULT_HERMITIAN_TRANSPOSE,(void(*)(void))MatMultTranspose_Z);
240: #else
241:       MatShellSetOperation(lanczos->Z,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Z);
242: #endif
243:       MatShellSetOperation(lanczos->Z,MATOP_CREATE_VECS,(void(*)(void))MatCreateVecs_Z);
244:       MatShellSetOperation(lanczos->Z,MATOP_DESTROY,(void(*)(void))MatDestroy_Z);
245:     }
246:     /* Explicit matrix is created here, when updating the scale */
247:     MatZUpdateScale(svd);
248:   }
249:   DSSetType(svd->ds,dstype);
250:   DSSetCompact(svd->ds,PETSC_TRUE);
251:   DSSetExtraRow(svd->ds,PETSC_TRUE);
252:   DSAllocate(svd->ds,svd->ncv+1);
253:   return 0;
254: }

256: static PetscErrorCode SVDOneSideTRLanczosMGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
257: {
258:   PetscReal      a,b;
259:   PetscInt       i,k=nconv+l;
260:   Vec            ui,ui1,vi;

262:   BVGetColumn(V,k,&vi);
263:   BVGetColumn(U,k,&ui);
264:   MatMult(svd->A,vi,ui);
265:   BVRestoreColumn(V,k,&vi);
266:   BVRestoreColumn(U,k,&ui);
267:   if (l>0) {
268:     BVSetActiveColumns(U,nconv,n);
269:     for (i=0;i<l;i++) work[i]=beta[i+nconv];
270:     BVMultColumn(U,-1.0,1.0,k,work);
271:   }
272:   BVNormColumn(U,k,NORM_2,&a);
273:   BVScaleColumn(U,k,1.0/a);
274:   alpha[k] = a;

276:   for (i=k+1;i<n;i++) {
277:     BVGetColumn(V,i,&vi);
278:     BVGetColumn(U,i-1,&ui1);
279:     MatMult(svd->AT,ui1,vi);
280:     BVRestoreColumn(V,i,&vi);
281:     BVRestoreColumn(U,i-1,&ui1);
282:     BVOrthonormalizeColumn(V,i,PETSC_FALSE,&b,NULL);
283:     beta[i-1] = b;

285:     BVGetColumn(V,i,&vi);
286:     BVGetColumn(U,i,&ui);
287:     MatMult(svd->A,vi,ui);
288:     BVRestoreColumn(V,i,&vi);
289:     BVGetColumn(U,i-1,&ui1);
290:     VecAXPY(ui,-b,ui1);
291:     BVRestoreColumn(U,i-1,&ui1);
292:     BVRestoreColumn(U,i,&ui);
293:     BVNormColumn(U,i,NORM_2,&a);
294:     BVScaleColumn(U,i,1.0/a);
295:     alpha[i] = a;
296:   }

298:   BVGetColumn(V,n,&vi);
299:   BVGetColumn(U,n-1,&ui1);
300:   MatMult(svd->AT,ui1,vi);
301:   BVRestoreColumn(V,n,&vi);
302:   BVRestoreColumn(U,n-1,&ui1);
303:   BVOrthogonalizeColumn(V,n,NULL,&b,NULL);
304:   beta[n-1] = b;
305:   return 0;
306: }

308: /*
309:   Custom CGS orthogonalization, preprocess after first orthogonalization
310: */
311: static PetscErrorCode SVDOrthogonalizeCGS(BV V,PetscInt i,PetscScalar* h,PetscReal a,BVOrthogRefineType refine,PetscReal eta,PetscReal *norm)
312: {
313:   PetscReal      sum,onorm;
314:   PetscScalar    dot;
315:   PetscInt       j;

317:   switch (refine) {
318:   case BV_ORTHOG_REFINE_NEVER:
319:     BVNormColumn(V,i,NORM_2,norm);
320:     break;
321:   case BV_ORTHOG_REFINE_ALWAYS:
322:     BVSetActiveColumns(V,0,i);
323:     BVDotColumn(V,i,h);
324:     BVMultColumn(V,-1.0,1.0,i,h);
325:     BVNormColumn(V,i,NORM_2,norm);
326:     break;
327:   case BV_ORTHOG_REFINE_IFNEEDED:
328:     dot = h[i];
329:     onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
330:     sum = 0.0;
331:     for (j=0;j<i;j++) {
332:       sum += PetscRealPart(h[j] * PetscConj(h[j]));
333:     }
334:     *norm = PetscRealPart(dot)/(a*a) - sum;
335:     if (*norm>0.0) *norm = PetscSqrtReal(*norm);
336:     else BVNormColumn(V,i,NORM_2,norm);
337:     if (*norm < eta*onorm) {
338:       BVSetActiveColumns(V,0,i);
339:       BVDotColumn(V,i,h);
340:       BVMultColumn(V,-1.0,1.0,i,h);
341:       BVNormColumn(V,i,NORM_2,norm);
342:     }
343:     break;
344:   }
345:   return 0;
346: }

348: static PetscErrorCode SVDOneSideTRLanczosCGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
349: {
350:   PetscReal          a,b,eta;
351:   PetscInt           i,j,k=nconv+l;
352:   Vec                ui,ui1,vi;
353:   BVOrthogRefineType refine;

355:   BVGetColumn(V,k,&vi);
356:   BVGetColumn(U,k,&ui);
357:   MatMult(svd->A,vi,ui);
358:   BVRestoreColumn(V,k,&vi);
359:   BVRestoreColumn(U,k,&ui);
360:   if (l>0) {
361:     BVSetActiveColumns(U,nconv,n);
362:     for (i=0;i<l;i++) work[i]=beta[i+nconv];
363:     BVMultColumn(U,-1.0,1.0,k,work);
364:   }
365:   BVGetOrthogonalization(V,NULL,&refine,&eta,NULL);

367:   for (i=k+1;i<n;i++) {
368:     BVGetColumn(V,i,&vi);
369:     BVGetColumn(U,i-1,&ui1);
370:     MatMult(svd->AT,ui1,vi);
371:     BVRestoreColumn(V,i,&vi);
372:     BVRestoreColumn(U,i-1,&ui1);
373:     BVNormColumnBegin(U,i-1,NORM_2,&a);
374:     if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
375:       BVSetActiveColumns(V,0,i+1);
376:       BVGetColumn(V,i,&vi);
377:       BVDotVecBegin(V,vi,work);
378:     } else {
379:       BVSetActiveColumns(V,0,i);
380:       BVDotColumnBegin(V,i,work);
381:     }
382:     BVNormColumnEnd(U,i-1,NORM_2,&a);
383:     if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
384:       BVDotVecEnd(V,vi,work);
385:       BVRestoreColumn(V,i,&vi);
386:       BVSetActiveColumns(V,0,i);
387:     } else BVDotColumnEnd(V,i,work);

389:     BVScaleColumn(U,i-1,1.0/a);
390:     for (j=0;j<i;j++) work[j] = work[j] / a;
391:     BVMultColumn(V,-1.0,1.0/a,i,work);
392:     SVDOrthogonalizeCGS(V,i,work,a,refine,eta,&b);
393:     BVScaleColumn(V,i,1.0/b);

396:     BVGetColumn(V,i,&vi);
397:     BVGetColumn(U,i,&ui);
398:     BVGetColumn(U,i-1,&ui1);
399:     MatMult(svd->A,vi,ui);
400:     VecAXPY(ui,-b,ui1);
401:     BVRestoreColumn(V,i,&vi);
402:     BVRestoreColumn(U,i,&ui);
403:     BVRestoreColumn(U,i-1,&ui1);

405:     alpha[i-1] = a;
406:     beta[i-1] = b;
407:   }

409:   BVGetColumn(V,n,&vi);
410:   BVGetColumn(U,n-1,&ui1);
411:   MatMult(svd->AT,ui1,vi);
412:   BVRestoreColumn(V,n,&vi);
413:   BVRestoreColumn(U,n-1,&ui1);

415:   BVNormColumnBegin(svd->U,n-1,NORM_2,&a);
416:   if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
417:     BVSetActiveColumns(V,0,n+1);
418:     BVGetColumn(V,n,&vi);
419:     BVDotVecBegin(V,vi,work);
420:   } else {
421:     BVSetActiveColumns(V,0,n);
422:     BVDotColumnBegin(V,n,work);
423:   }
424:   BVNormColumnEnd(svd->U,n-1,NORM_2,&a);
425:   if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
426:     BVDotVecEnd(V,vi,work);
427:     BVRestoreColumn(V,n,&vi);
428:   } else BVDotColumnEnd(V,n,work);

430:   BVScaleColumn(U,n-1,1.0/a);
431:   for (j=0;j<n;j++) work[j] = work[j] / a;
432:   BVMultColumn(V,-1.0,1.0/a,n,work);
433:   SVDOrthogonalizeCGS(V,n,work,a,refine,eta,&b);
434:   BVSetActiveColumns(V,nconv,n);
435:   alpha[n-1] = a;
436:   beta[n-1] = b;
437:   return 0;
438: }

440: PetscErrorCode SVDSolve_TRLanczos(SVD svd)
441: {
442:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
443:   PetscReal      *alpha,*beta;
444:   PetscScalar    *swork=NULL,*w;
445:   PetscInt       i,k,l,nv,ld;
446:   Mat            U,V;
447:   PetscBool      breakdown=PETSC_FALSE;
448:   BVOrthogType   orthog;

450:   PetscCitationsRegister(citation,&cited);
451:   /* allocate working space */
452:   DSGetLeadingDimension(svd->ds,&ld);
453:   BVGetOrthogonalization(svd->V,&orthog,NULL,NULL,NULL);
454:   PetscMalloc1(ld,&w);
455:   if (lanczos->oneside) PetscMalloc1(svd->ncv+1,&swork);

457:   /* normalize start vector */
458:   if (!svd->nini) {
459:     BVSetRandomColumn(svd->V,0);
460:     BVOrthonormalizeColumn(svd->V,0,PETSC_TRUE,NULL,NULL);
461:   }

463:   l = 0;
464:   while (svd->reason == SVD_CONVERGED_ITERATING) {
465:     svd->its++;

467:     /* inner loop */
468:     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
469:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
470:     beta = alpha + ld;
471:     if (lanczos->oneside) {
472:       if (orthog == BV_ORTHOG_MGS) SVDOneSideTRLanczosMGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);
473:       else SVDOneSideTRLanczosCGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);
474:     } else SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv+l,&nv,&breakdown);
475:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
476:     BVScaleColumn(svd->V,nv,1.0/beta[nv-1]);
477:     BVSetActiveColumns(svd->V,svd->nconv,nv);
478:     BVSetActiveColumns(svd->U,svd->nconv,nv);

480:     /* solve projected problem */
481:     DSSetDimensions(svd->ds,nv,svd->nconv,svd->nconv+l);
482:     DSSVDSetDimensions(svd->ds,nv);
483:     DSSetState(svd->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE);
484:     DSSolve(svd->ds,w,NULL);
485:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
486:     DSUpdateExtraRow(svd->ds);
487:     DSSynchronize(svd->ds,w,NULL);
488:     for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);

490:     /* check convergence */
491:     SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,1.0,&k);
492:     (*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx);

494:     /* update l */
495:     if (svd->reason != SVD_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
496:     else l = PetscMax(1,(PetscInt)((nv-k)*lanczos->keep));
497:     if (!lanczos->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged triplets */
498:     if (l) PetscInfo(svd,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l);

500:     if (svd->reason == SVD_CONVERGED_ITERATING) {
501:       if (PetscUnlikely(breakdown || k==nv)) {
502:         /* Start a new bidiagonalization */
503:         PetscInfo(svd,"Breakdown in bidiagonalization (it=%" PetscInt_FMT ")\n",svd->its);
504:         if (k<svd->nsv) {
505:           BVSetRandomColumn(svd->V,k);
506:           BVOrthonormalizeColumn(svd->V,k,PETSC_FALSE,NULL,&breakdown);
507:           if (breakdown) {
508:             svd->reason = SVD_DIVERGED_BREAKDOWN;
509:             PetscInfo(svd,"Unable to generate more start vectors\n");
510:           }
511:         }
512:       } else DSTruncate(svd->ds,k+l,PETSC_FALSE);
513:     }

515:     /* compute converged singular vectors and restart vectors */
516:     DSGetMat(svd->ds,DS_MAT_V,&V);
517:     BVMultInPlace(svd->V,V,svd->nconv,k+l);
518:     DSRestoreMat(svd->ds,DS_MAT_V,&V);
519:     DSGetMat(svd->ds,DS_MAT_U,&U);
520:     BVMultInPlace(svd->U,U,svd->nconv,k+l);
521:     DSRestoreMat(svd->ds,DS_MAT_U,&U);

523:     /* copy the last vector to be the next initial vector */
524:     if (svd->reason == SVD_CONVERGED_ITERATING && !breakdown) BVCopyColumn(svd->V,nv,k+l);

526:     svd->nconv = k;
527:     SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
528:   }

530:   /* orthonormalize U columns in one side method */
531:   if (lanczos->oneside) {
532:     for (i=0;i<svd->nconv;i++) BVOrthonormalizeColumn(svd->U,i,PETSC_FALSE,NULL,NULL);
533:   }

535:   /* free working space */
536:   PetscFree(w);
537:   if (swork) PetscFree(swork);
538:   DSTruncate(svd->ds,svd->nconv,PETSC_TRUE);
539:   return 0;
540: }

542: /* Given n computed generalized singular values in sigmain, backtransform them
543:    in sigmaout by undoing scaling and reciprocating if swapped=true. Also updates vectors V
544:    if given. If sigmaout=NULL then the result overwrites sigmain. */
545: static PetscErrorCode SVDLanczosBackTransform(SVD svd,PetscInt n,PetscReal *sigmain,PetscReal *sigmaout,BV V)
546: {
547:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
548:   PetscInt      i;
549:   PetscReal     c,s,r,f,scalef;

551:   scalef = svd->swapped? 1.0/lanczos->scalef: lanczos->scalef;
552:   for (i=0;i<n;i++) {
553:     if (V && scalef != 1.0) {
554:       s = 1.0/PetscSqrtReal(1.0+sigmain[i]*sigmain[i]);
555:       c = sigmain[i]*s;
556:       r = s/scalef;
557:       f = 1.0/PetscSqrtReal(c*c+r*r);
558:       BVScaleColumn(V,i,f);
559:     }
560:     if (sigmaout) {
561:       if (svd->swapped) sigmaout[i] = 1.0/(sigmain[i]*scalef);
562:       else sigmaout[i] = sigmain[i]*scalef;
563:     } else {
564:       sigmain[i] *= scalef;
565:       if (svd->swapped) sigmain[i] = 1.0/sigmain[i];
566:     }
567:   }
568:   return 0;
569: }

571: static PetscErrorCode SVDLanczosGSingle(SVD svd,PetscReal *alpha,PetscReal *beta,Mat Z,BV V,BV U,KSP ksp,PetscInt k,PetscInt *n,PetscBool *breakdown)
572: {
573:   SVD_TRLANCZOS     *lanczos = (SVD_TRLANCZOS*)svd->data;
574:   PetscInt          i,j,m;
575:   const PetscScalar *carr;
576:   PetscScalar       *arr;
577:   Vec               u,v,ut=svd->workl[0],x=svd->workr[0],v1,u1,u2;
578:   PetscBool         lindep=PETSC_FALSE;

580:   MatCreateVecsEmpty(svd->A,NULL,&v1);
581:   BVGetColumn(V,k,&v);
582:   BVGetColumn(U,k,&u);

584:   /* Form ut=[u;0] */
585:   VecZeroEntries(ut);
586:   VecGetLocalSize(u,&m);
587:   VecGetArrayRead(u,&carr);
588:   VecGetArray(ut,&arr);
589:   for (j=0; j<m; j++) arr[j] = carr[j];
590:   VecRestoreArrayRead(u,&carr);
591:   VecRestoreArray(ut,&arr);

593:   /* Solve least squares problem */
594:   KSPSolve(ksp,ut,x);

596:   MatMult(Z,x,v);

598:   BVRestoreColumn(U,k,&u);
599:   BVRestoreColumn(V,k,&v);
600:   BVOrthonormalizeColumn(V,k,PETSC_FALSE,alpha+k,&lindep);
601:   if (PetscUnlikely(lindep)) {
602:     *n = k;
603:     if (breakdown) *breakdown = lindep;
604:     return 0;
605:   }

607:   for (i=k+1; i<*n; i++) {

609:     /* Compute vector i of BV U */
610:     BVGetColumn(V,i-1,&v);
611:     VecGetArray(v,&arr);
612:     VecPlaceArray(v1,arr);
613:     VecRestoreArray(v,&arr);
614:     BVRestoreColumn(V,i-1,&v);
615:     BVInsertVec(U,i,v1);
616:     VecResetArray(v1);
617:     BVOrthonormalizeColumn(U,i,PETSC_FALSE,beta+i-1,&lindep);
618:     if (PetscUnlikely(lindep)) {
619:       *n = i;
620:       break;
621:     }

623:     /* Compute vector i of BV V */

625:     BVGetColumn(V,i,&v);
626:     BVGetColumn(U,i,&u);

628:     /* Form ut=[u;0] */
629:     VecGetArrayRead(u,&carr);
630:     VecGetArray(ut,&arr);
631:     for (j=0; j<m; j++) arr[j] = carr[j];
632:     VecRestoreArrayRead(u,&carr);
633:     VecRestoreArray(ut,&arr);

635:     /* Solve least squares problem */
636:     KSPSolve(ksp,ut,x);

638:     MatMult(Z,x,v);

640:     BVRestoreColumn(U,i,&u);
641:     BVRestoreColumn(V,i,&v);
642:     if (!lanczos->oneside || i==k+1) BVOrthonormalizeColumn(V,i,PETSC_FALSE,alpha+i,&lindep);
643:     else {  /* cheap computation of V[i], if restart (i==k+1) do a full reorthogonalization */
644:       BVGetColumn(V,i,&u2);
645:       BVGetColumn(V,i-1,&u1);
646:       VecAXPY(u2,-beta[i-1],u1);
647:       BVRestoreColumn(V,i-1,&u1);
648:       VecNorm(u2,NORM_2,&alpha[i]);
649:       if (alpha[i]==0.0) lindep = PETSC_TRUE;
650:       else VecScale(u2,1.0/alpha[i]);
651:       BVRestoreColumn(V,i,&u2);
652:     }
653:     if (PetscUnlikely(lindep)) {
654:       *n = i;
655:       break;
656:     }
657:   }

659:   /* Compute vector n of BV U */
660:   if (!lindep) {
661:     BVGetColumn(V,*n-1,&v);
662:     VecGetArray(v,&arr);
663:     VecPlaceArray(v1,arr);
664:     VecRestoreArray(v,&arr);
665:     BVRestoreColumn(V,*n-1,&v);
666:     BVInsertVec(U,*n,v1);
667:     VecResetArray(v1);
668:     BVOrthonormalizeColumn(U,*n,PETSC_FALSE,beta+*n-1,&lindep);
669:   }
670:   if (breakdown) *breakdown = lindep;
671:   VecDestroy(&v1);
672:   return 0;
673: }

675: /* solve generalized problem with single bidiagonalization of Q_A */
676: PetscErrorCode SVDSolve_TRLanczosGSingle(SVD svd,BV U1,BV V)
677: {
678:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
679:   PetscReal      *alpha,*beta,normr,scaleth,sigma0,*sigma;
680:   PetscScalar    *w;
681:   PetscInt       i,k,l,nv,ld;
682:   Mat            U,VV;
683:   PetscBool      breakdown=PETSC_FALSE;

685:   DSGetLeadingDimension(svd->ds,&ld);
686:   PetscMalloc2(ld,&w,ld,&sigma);
687:   normr = (svd->conv==SVD_CONV_ABS)? PetscMax(svd->nrma,svd->nrmb*lanczos->scalef): 1.0;
688:   /* Convert scale threshold th=c/s to the corresponding c */
689:   scaleth = (lanczos->scaleth!=0)? lanczos->scaleth/PetscSqrtReal(lanczos->scaleth*lanczos->scaleth+1): 0.0;

691:   /* normalize start vector */
692:   if (!svd->ninil) {
693:     BVSetRandomColumn(U1,0);
694:     BVOrthonormalizeColumn(U1,0,PETSC_TRUE,NULL,NULL);
695:   }

697:   l = 0;
698:   while (svd->reason == SVD_CONVERGED_ITERATING) {
699:     svd->its++;

701:     /* inner loop */
702:     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
703:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
704:     beta = alpha + ld;
705:     SVDLanczosGSingle(svd,alpha,beta,lanczos->Z,V,U1,lanczos->ksp,svd->nconv+l,&nv,&breakdown);
706:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
707:     BVSetActiveColumns(V,svd->nconv,nv);
708:     BVSetActiveColumns(U1,svd->nconv,nv);

710:     /* solve projected problem */
711:     DSSetDimensions(svd->ds,nv,svd->nconv,svd->nconv+l);
712:     DSSVDSetDimensions(svd->ds,nv);
713:     DSSetState(svd->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE);
714:     DSSolve(svd->ds,w,NULL);
715:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
716:     DSUpdateExtraRow(svd->ds);
717:     DSSynchronize(svd->ds,w,NULL);
718:     for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);

720:     /* check convergence */
721:     SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,normr,&k);
722:     (*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx);

724:     sigma0 = svd->which==SVD_LARGEST? svd->sigma[0] : 1.0/svd->sigma[0];
725:     if (scaleth!=0 && k==0 && sigma0>scaleth) {

727:       /* Scale and start from scratch */
728:       lanczos->scalef *= svd->sigma[0]/PetscSqrtReal(1-svd->sigma[0]*svd->sigma[0]);
729:       PetscInfo(svd,"Scaling by factor %g and starting from scratch\n",(double)lanczos->scalef);
730:       MatZUpdateScale(svd);
731:       if (svd->conv==SVD_CONV_ABS) normr = PetscMax(svd->nrma,svd->nrmb*lanczos->scalef);
732:       l = 0;

734:     } else {

736:       /* update l */
737:       if (svd->reason != SVD_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
738:       else l = PetscMax(1,(PetscInt)((nv-k)*lanczos->keep));
739:       if (!lanczos->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged triplets */
740:       if (l) PetscInfo(svd,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l);

742:       if (svd->reason == SVD_CONVERGED_ITERATING) {
743:         if (PetscUnlikely(breakdown || k==nv)) {
744:           /* Start a new bidiagonalization */
745:           PetscInfo(svd,"Breakdown in bidiagonalization (it=%" PetscInt_FMT ")\n",svd->its);
746:           if (k<svd->nsv) {
747:             BVSetRandomColumn(U1,k);
748:             BVOrthonormalizeColumn(U1,k,PETSC_FALSE,NULL,&breakdown);
749:             if (breakdown) {
750:               svd->reason = SVD_DIVERGED_BREAKDOWN;
751:               PetscInfo(svd,"Unable to generate more start vectors\n");
752:             }
753:           }
754:         } else DSTruncate(svd->ds,k+l,PETSC_FALSE);
755:       }

757:       /* compute converged singular vectors and restart vectors */
758:       DSGetMat(svd->ds,DS_MAT_U,&U);
759:       BVMultInPlace(V,U,svd->nconv,k+l);
760:       DSRestoreMat(svd->ds,DS_MAT_U,&U);
761:       DSGetMat(svd->ds,DS_MAT_V,&VV);
762:       BVMultInPlace(U1,VV,svd->nconv,k+l);
763:       DSRestoreMat(svd->ds,DS_MAT_V,&VV);

765:       /* copy the last vector to be the next initial vector */
766:       if (svd->reason == SVD_CONVERGED_ITERATING && !breakdown) BVCopyColumn(U1,nv,k+l);
767:     }

769:     svd->nconv = k;
770:     SVDLanczosBackTransform(svd,nv,svd->sigma,sigma,NULL);
771:     SVDMonitor(svd,svd->its,svd->nconv,sigma,svd->errest,nv);
772:   }

774:   PetscFree2(w,sigma);
775:   return 0;
776: }

778: /* Move generalized left singular vectors (0..nconv) from U1 and U2 to its final destination svd->U (single variant) */
779: static inline PetscErrorCode SVDLeftSingularVectors_Single(SVD svd,BV U1,BV U2)
780: {
781:   PetscInt          i,k,m,p;
782:   Vec               u,u1,u2;
783:   PetscScalar       *ua,*u2a;
784:   const PetscScalar *u1a;
785:   PetscReal         s;

787:   MatGetLocalSize(svd->A,&m,NULL);
788:   MatGetLocalSize(svd->B,&p,NULL);
789:   for (i=0;i<svd->nconv;i++) {
790:     BVGetColumn(U1,i,&u1);
791:     BVGetColumn(U2,i,&u2);
792:     BVGetColumn(svd->U,i,&u);
793:     VecGetArrayRead(u1,&u1a);
794:     VecGetArray(u,&ua);
795:     VecGetArray(u2,&u2a);
796:     /* Copy column from U1 to upper part of u */
797:     for (k=0;k<m;k++) ua[k] = u1a[k];
798:     /* Copy column from lower part of U to U2. Orthogonalize column in U2 and copy back to U */
799:     for (k=0;k<p;k++) u2a[k] = ua[m+k];
800:     VecRestoreArray(u2,&u2a);
801:     BVRestoreColumn(U2,i,&u2);
802:     BVOrthonormalizeColumn(U2,i,PETSC_FALSE,&s,NULL);
803:     BVGetColumn(U2,i,&u2);
804:     VecGetArray(u2,&u2a);
805:     for (k=0;k<p;k++) ua[m+k] = u2a[k];
806:     /* Update singular value */
807:     svd->sigma[i] /= s;
808:     VecRestoreArrayRead(u1,&u1a);
809:     VecRestoreArray(u,&ua);
810:     VecRestoreArray(u2,&u2a);
811:     BVRestoreColumn(U1,i,&u1);
812:     BVRestoreColumn(U2,i,&u2);
813:     BVRestoreColumn(svd->U,i,&u);
814:   }
815:   return 0;
816: }

818: static PetscErrorCode SVDLanczosGUpper(SVD svd,PetscReal *alpha,PetscReal *beta,PetscReal *alphah,PetscReal *betah,Mat Z,BV U1,BV U2,BV V,KSP ksp,PetscInt k,PetscInt *n,PetscBool *breakdown)
819: {
820:   SVD_TRLANCZOS     *lanczos = (SVD_TRLANCZOS*)svd->data;
821:   PetscInt          i,j,m,p;
822:   const PetscScalar *carr;
823:   PetscScalar       *arr,*u2arr;
824:   Vec               u,v,ut=svd->workl[0],x=svd->workr[0],v1,u1,u2;
825:   PetscBool         lindep=PETSC_FALSE,lindep1=PETSC_FALSE,lindep2=PETSC_FALSE;

827:   MatCreateVecsEmpty(svd->A,NULL,&v1);
828:   MatGetLocalSize(svd->A,&m,NULL);
829:   MatGetLocalSize(svd->B,&p,NULL);

831:   for (i=k; i<*n; i++) {
832:     /* Compute vector i of BV U1 */
833:     BVGetColumn(V,i,&v);
834:     VecGetArrayRead(v,&carr);
835:     VecPlaceArray(v1,carr);
836:     BVInsertVec(U1,i,v1);
837:     VecResetArray(v1);
838:     if (!lanczos->oneside || i==k) BVOrthonormalizeColumn(U1,i,PETSC_FALSE,alpha+i,&lindep1);
839:     else {  /* cheap computation of U1[i], if restart (i==k) do a full reorthogonalization */
840:       BVGetColumn(U1,i,&u2);
841:       if (i>0) {
842:         BVGetColumn(U1,i-1,&u1);
843:         VecAXPY(u2,-beta[i-1],u1);
844:         BVRestoreColumn(U1,i-1,&u1);
845:       }
846:       VecNorm(u2,NORM_2,&alpha[i]);
847:       if (alpha[i]==0.0) lindep = PETSC_TRUE;
848:       else VecScale(u2,1.0/alpha[i]);
849:       BVRestoreColumn(U1,i,&u2);
850:     }

852:     /* Compute vector i of BV U2 */
853:     BVGetColumn(U2,i,&u2);
854:     VecGetArray(u2,&u2arr);
855:     if (i%2) {
856:       for (j=0; j<p; j++) u2arr[j] = -carr[m+j];
857:     } else {
858:       for (j=0; j<p; j++) u2arr[j] = carr[m+j];
859:     }
860:     VecRestoreArray(u2,&u2arr);
861:     VecRestoreArrayRead(v,&carr);
862:     BVRestoreColumn(V,i,&v);
863:     if (lanczos->oneside && i>k) {  /* cheap computation of U2[i], if restart (i==k) do a full reorthogonalization */
864:       if (i>0) {
865:         BVGetColumn(U2,i-1,&u1);
866:         VecAXPY(u2,(i%2)?betah[i-1]:-betah[i-1],u1);
867:         BVRestoreColumn(U2,i-1,&u1);
868:       }
869:       VecNorm(u2,NORM_2,&alphah[i]);
870:       if (alphah[i]==0.0) lindep = PETSC_TRUE;
871:       else VecScale(u2,1.0/alphah[i]);
872:     }
873:     BVRestoreColumn(U2,i,&u2);
874:     if (!lanczos->oneside || i==k) BVOrthonormalizeColumn(U2,i,PETSC_FALSE,alphah+i,&lindep2);
875:     if (i%2) alphah[i] = -alphah[i];
876:     if (PetscUnlikely(lindep1 || lindep2)) {
877:       lindep = PETSC_TRUE;
878:       *n = i;
879:       break;
880:     }

882:     /* Compute vector i+1 of BV V */
883:     BVGetColumn(V,i+1,&v);
884:     /* Form ut=[u;0] */
885:     BVGetColumn(U1,i,&u);
886:     VecZeroEntries(ut);
887:     VecGetArrayRead(u,&carr);
888:     VecGetArray(ut,&arr);
889:     for (j=0; j<m; j++) arr[j] = carr[j];
890:     VecRestoreArrayRead(u,&carr);
891:     VecRestoreArray(ut,&arr);
892:     /* Solve least squares problem */
893:     KSPSolve(ksp,ut,x);
894:     MatMult(Z,x,v);
895:     BVRestoreColumn(U1,i,&u);
896:     BVRestoreColumn(V,i+1,&v);
897:     BVOrthonormalizeColumn(V,i+1,PETSC_FALSE,beta+i,&lindep);
898:     betah[i] = -alpha[i]*beta[i]/alphah[i];
899:     if (PetscUnlikely(lindep)) {
900:       *n = i;
901:       break;
902:     }
903:   }
904:   if (breakdown) *breakdown = lindep;
905:   VecDestroy(&v1);
906:   return 0;
907: }

909: /* generate random initial vector in column k for joint upper-upper bidiagonalization */
910: static inline PetscErrorCode SVDInitialVectorGUpper(SVD svd,BV V,BV U1,PetscInt k,PetscBool *breakdown)
911: {
912:   SVD_TRLANCZOS     *lanczos = (SVD_TRLANCZOS*)svd->data;
913:   Vec               u,v,ut=svd->workl[0],x=svd->workr[0];
914:   PetscInt          m,j;
915:   PetscScalar       *arr;
916:   const PetscScalar *carr;

918:   /* Form ut=[u;0] where u is the k-th column of U1 */
919:   VecZeroEntries(ut);
920:   BVGetColumn(U1,k,&u);
921:   VecGetLocalSize(u,&m);
922:   VecGetArrayRead(u,&carr);
923:   VecGetArray(ut,&arr);
924:   for (j=0; j<m; j++) arr[j] = carr[j];
925:   VecRestoreArrayRead(u,&carr);
926:   VecRestoreArray(ut,&arr);
927:   BVRestoreColumn(U1,k,&u);
928:   /* Solve least squares problem Z*x=ut for x. Then set v=Z*x */
929:   KSPSolve(lanczos->ksp,ut,x);
930:   BVGetColumn(V,k,&v);
931:   MatMult(lanczos->Z,x,v);
932:   BVRestoreColumn(V,k,&v);
933:   if (breakdown) BVOrthonormalizeColumn(V,k,PETSC_FALSE,NULL,breakdown);
934:   else BVOrthonormalizeColumn(V,k,PETSC_TRUE,NULL,NULL);
935:   return 0;
936: }

938: /* solve generalized problem with joint upper-upper bidiagonalization */
939: PetscErrorCode SVDSolve_TRLanczosGUpper(SVD svd,BV U1,BV U2,BV V)
940: {
941:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
942:   PetscReal      *alpha,*beta,*alphah,*betah,normr,sigma0,*sigma;
943:   PetscScalar    *w;
944:   PetscInt       i,k,l,nv,ld;
945:   Mat            U,Vmat,X;
946:   PetscBool      breakdown=PETSC_FALSE;

948:   DSGetLeadingDimension(svd->ds,&ld);
949:   PetscMalloc2(ld,&w,ld,&sigma);
950:   normr = (svd->conv==SVD_CONV_ABS)? PetscMax(svd->nrma,svd->nrmb*lanczos->scalef): 1.0;

952:   /* normalize start vector */
953:   if (!svd->ninil) BVSetRandomColumn(U1,0);
954:   SVDInitialVectorGUpper(svd,V,U1,0,NULL);

956:   l = 0;
957:   while (svd->reason == SVD_CONVERGED_ITERATING) {
958:     svd->its++;

960:     /* inner loop */
961:     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
962:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
963:     DSGetArrayReal(svd->ds,DS_MAT_D,&alphah);
964:     beta = alpha + ld;
965:     betah = alpha + 2*ld;
966:     SVDLanczosGUpper(svd,alpha,beta,alphah,betah,lanczos->Z,U1,U2,V,lanczos->ksp,svd->nconv+l,&nv,&breakdown);
967:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
968:     DSRestoreArrayReal(svd->ds,DS_MAT_D,&alphah);
969:     BVSetActiveColumns(V,svd->nconv,nv);
970:     BVSetActiveColumns(U1,svd->nconv,nv);
971:     BVSetActiveColumns(U2,svd->nconv,nv);

973:     /* solve projected problem */
974:     DSSetDimensions(svd->ds,nv,svd->nconv,svd->nconv+l);
975:     DSGSVDSetDimensions(svd->ds,nv,nv);
976:     DSSetState(svd->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE);
977:     DSSolve(svd->ds,w,NULL);
978:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
979:     DSUpdateExtraRow(svd->ds);
980:     DSSynchronize(svd->ds,w,NULL);
981:     for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);

983:     /* check convergence */
984:     SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,normr,&k);
985:     (*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx);

987:     sigma0 = svd->which==SVD_LARGEST? svd->sigma[0] : 1.0/svd->sigma[0];
988:     if (lanczos->scaleth!=0 && k==0 && sigma0>lanczos->scaleth) {

990:       /* Scale and start from scratch */
991:       lanczos->scalef *= svd->sigma[0];
992:       PetscInfo(svd,"Scaling by factor %g and starting from scratch\n",(double)lanczos->scalef);
993:       MatZUpdateScale(svd);
994:       if (svd->conv==SVD_CONV_ABS) normr = PetscMax(svd->nrma,svd->nrmb*lanczos->scalef);
995:       l = 0;
996:       if (!svd->ninil) BVSetRandomColumn(U1,0);
997:       SVDInitialVectorGUpper(svd,V,U1,0,NULL);

999:     } else {

1001:       /* update l */
1002:       if (svd->reason != SVD_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
1003:       else l = PetscMax(1,(PetscInt)((nv-k)*lanczos->keep));
1004:       if (!lanczos->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged triplets */
1005:       if (l) PetscInfo(svd,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l);

1007:       if (svd->reason == SVD_CONVERGED_ITERATING) {
1008:         if (PetscUnlikely(breakdown || k==nv)) {
1009:           /* Start a new bidiagonalization */
1010:           PetscInfo(svd,"Breakdown in bidiagonalization (it=%" PetscInt_FMT ")\n",svd->its);
1011:           if (k<svd->nsv) {
1012:             BVSetRandomColumn(U1,k);
1013:             SVDInitialVectorGUpper(svd,V,U1,k,&breakdown);
1014:             if (breakdown) {
1015:               svd->reason = SVD_DIVERGED_BREAKDOWN;
1016:               PetscInfo(svd,"Unable to generate more start vectors\n");
1017:             }
1018:           }
1019:         } else DSTruncate(svd->ds,k+l,PETSC_FALSE);
1020:       }
1021:       /* compute converged singular vectors and restart vectors */
1022:       DSGetMat(svd->ds,DS_MAT_X,&X);
1023:       BVMultInPlace(V,X,svd->nconv,k+l);
1024:       DSRestoreMat(svd->ds,DS_MAT_X,&X);
1025:       DSGetMat(svd->ds,DS_MAT_U,&U);
1026:       BVMultInPlace(U1,U,svd->nconv,k+l);
1027:       DSRestoreMat(svd->ds,DS_MAT_U,&U);
1028:       DSGetMat(svd->ds,DS_MAT_V,&Vmat);
1029:       BVMultInPlace(U2,Vmat,svd->nconv,k+l);
1030:       DSRestoreMat(svd->ds,DS_MAT_V,&Vmat);

1032:       /* copy the last vector to be the next initial vector */
1033:       if (svd->reason == SVD_CONVERGED_ITERATING && !breakdown) BVCopyColumn(V,nv,k+l);
1034:     }

1036:     svd->nconv = k;
1037:     SVDLanczosBackTransform(svd,nv,svd->sigma,sigma,NULL);
1038:     SVDMonitor(svd,svd->its,svd->nconv,sigma,svd->errest,nv);
1039:   }

1041:   PetscFree2(w,sigma);
1042:   return 0;
1043: }

1045: /* Move generalized left singular vectors (0..nconv) from U1 and U2 to its final destination svd->U (upper and lower variants) */
1046: static inline PetscErrorCode SVDLeftSingularVectors(SVD svd,BV U1,BV U2)
1047: {
1048:   PetscInt          i,k,m,p;
1049:   Vec               u,u1,u2;
1050:   PetscScalar       *ua;
1051:   const PetscScalar *u1a,*u2a;

1053:   BVGetSizes(U1,&m,NULL,NULL);
1054:   BVGetSizes(U2,&p,NULL,NULL);
1055:   for (i=0;i<svd->nconv;i++) {
1056:     BVGetColumn(U1,i,&u1);
1057:     BVGetColumn(U2,i,&u2);
1058:     BVGetColumn(svd->U,i,&u);
1059:     VecGetArrayRead(u1,&u1a);
1060:     VecGetArrayRead(u2,&u2a);
1061:     VecGetArray(u,&ua);
1062:     /* Copy column from u1 to upper part of u */
1063:     for (k=0;k<m;k++) ua[k] = u1a[k];
1064:     /* Copy column from u2 to lower part of u */
1065:     for (k=0;k<p;k++) ua[m+k] = u2a[k];
1066:     VecRestoreArrayRead(u1,&u1a);
1067:     VecRestoreArrayRead(u2,&u2a);
1068:     VecRestoreArray(u,&ua);
1069:     BVRestoreColumn(U1,i,&u1);
1070:     BVRestoreColumn(U2,i,&u2);
1071:     BVRestoreColumn(svd->U,i,&u);
1072:   }
1073:   return 0;
1074: }

1076: static PetscErrorCode SVDLanczosGLower(SVD svd,PetscReal *alpha,PetscReal *beta,PetscReal *alphah,PetscReal *betah,Mat Z,BV U1,BV U2,BV V,KSP ksp,PetscInt k,PetscInt *n,PetscBool *breakdown)
1077: {
1078:   SVD_TRLANCZOS     *lanczos = (SVD_TRLANCZOS*)svd->data;
1079:   PetscInt          i,j,m,p;
1080:   const PetscScalar *carr;
1081:   PetscScalar       *arr,*u2arr;
1082:   Vec               u,v,ut=svd->workl[0],x=svd->workr[0],v1,u1,u2;
1083:   PetscBool         lindep=PETSC_FALSE;

1085:   MatCreateVecsEmpty(svd->A,NULL,&v1);
1086:   MatGetLocalSize(svd->A,&m,NULL);
1087:   MatGetLocalSize(svd->B,&p,NULL);

1089:   for (i=k; i<*n; i++) {
1090:     /* Compute vector i of BV U2 */
1091:     BVGetColumn(V,i,&v);
1092:     VecGetArrayRead(v,&carr);
1093:     BVGetColumn(U2,i,&u2);
1094:     VecGetArray(u2,&u2arr);
1095:     if (i%2) {
1096:       for (j=0; j<p; j++) u2arr[j] = -carr[m+j];
1097:     } else {
1098:       for (j=0; j<p; j++) u2arr[j] = carr[m+j];
1099:     }
1100:     VecRestoreArray(u2,&u2arr);
1101:     if (lanczos->oneside && i>k) {  /* cheap computation of U2[i], if restart (i==k) do a full reorthogonalization */
1102:       if (i>0) {
1103:         BVGetColumn(U2,i-1,&u1);
1104:         VecAXPY(u2,(i%2)?betah[i-1]:-betah[i-1],u1);
1105:         BVRestoreColumn(U2,i-1,&u1);
1106:       }
1107:       VecNorm(u2,NORM_2,&alphah[i]);
1108:       if (alphah[i]==0.0) lindep = PETSC_TRUE;
1109:       else VecScale(u2,1.0/alphah[i]);
1110:     }
1111:     BVRestoreColumn(U2,i,&u2);
1112:     if (!lanczos->oneside || i==k) BVOrthonormalizeColumn(U2,i,PETSC_FALSE,alphah+i,&lindep);
1113:     if (i%2) alphah[i] = -alphah[i];
1114:     if (PetscUnlikely(lindep)) {
1115:       BVRestoreColumn(V,i,&v);
1116:       *n = i;
1117:       break;
1118:     }

1120:     /* Compute vector i+1 of BV U1 */
1121:     VecPlaceArray(v1,carr);
1122:     BVInsertVec(U1,i+1,v1);
1123:     VecResetArray(v1);
1124:     BVOrthonormalizeColumn(U1,i+1,PETSC_FALSE,beta+i,&lindep);
1125:     VecRestoreArrayRead(v,&carr);
1126:     BVRestoreColumn(V,i,&v);
1127:     if (PetscUnlikely(lindep)) {
1128:       *n = i+1;
1129:       break;
1130:     }

1132:     /* Compute vector i+1 of BV V */
1133:     BVGetColumn(V,i+1,&v);
1134:     /* Form ut=[u;0] where u is column i+1 of BV U1 */
1135:     BVGetColumn(U1,i+1,&u);
1136:     VecZeroEntries(ut);
1137:     VecGetArrayRead(u,&carr);
1138:     VecGetArray(ut,&arr);
1139:     for (j=0; j<m; j++) arr[j] = carr[j];
1140:     VecRestoreArrayRead(u,&carr);
1141:     VecRestoreArray(ut,&arr);
1142:     /* Solve least squares problem */
1143:     KSPSolve(ksp,ut,x);
1144:     MatMult(Z,x,v);
1145:     BVRestoreColumn(U1,i+1,&u);
1146:     BVRestoreColumn(V,i+1,&v);
1147:     if (!lanczos->oneside || i==k) BVOrthonormalizeColumn(V,i+1,PETSC_FALSE,alpha+i+1,&lindep);
1148:     else {  /* cheap computation of V[i+1], if restart (i==k) do a full reorthogonalization */
1149:       BVGetColumn(V,i+1,&u2);
1150:       BVGetColumn(V,i,&u1);
1151:       VecAXPY(u2,-beta[i],u1);
1152:       BVRestoreColumn(V,i,&u1);
1153:       VecNorm(u2,NORM_2,&alpha[i+1]);
1154:       if (alpha[i+1]==0.0) lindep = PETSC_TRUE;
1155:       else VecScale(u2,1.0/alpha[i+1]);
1156:       BVRestoreColumn(V,i+1,&u2);
1157:     }
1158:     betah[i] = -alpha[i+1]*beta[i]/alphah[i];
1159:     if (PetscUnlikely(lindep)) {
1160:       *n = i+1;
1161:       break;
1162:     }
1163:   }
1164:   if (breakdown) *breakdown = lindep;
1165:   VecDestroy(&v1);
1166:   return 0;
1167: }

1169: /* generate random initial vector in column k for joint lower-upper bidiagonalization */
1170: static inline PetscErrorCode SVDInitialVectorGLower(SVD svd,BV V,BV U1,BV U2,PetscInt k,PetscBool *breakdown)
1171: {
1172:   SVD_TRLANCZOS     *lanczos = (SVD_TRLANCZOS*)svd->data;
1173:   const PetscScalar *carr;
1174:   PetscScalar       *arr;
1175:   PetscReal         *alpha;
1176:   PetscInt          j,m,p;
1177:   Vec               u,uh,v,ut=svd->workl[0],x=svd->workr[0];

1179:   MatGetLocalSize(svd->A,&m,NULL);
1180:   MatGetLocalSize(svd->B,&p,NULL);
1181:   /* Form ut=[0;uh], where uh is the k-th column of U2 */
1182:   BVGetColumn(U2,k,&uh);
1183:   VecZeroEntries(ut);
1184:   VecGetArrayRead(uh,&carr);
1185:   VecGetArray(ut,&arr);
1186:   for (j=0; j<p; j++) arr[m+j] = carr[j];
1187:   VecRestoreArrayRead(uh,&carr);
1188:   VecRestoreArray(ut,&arr);
1189:   BVRestoreColumn(U2,k,&uh);
1190:   /* Solve least squares problem Z*x=ut for x. Then set ut=Z*x */
1191:   KSPSolve(lanczos->ksp,ut,x);
1192:   MatMult(lanczos->Z,x,ut);
1193:   /* Form u, column k of BV U1, as the upper part of ut and orthonormalize */
1194:   MatCreateVecsEmpty(svd->A,NULL,&u);
1195:   VecGetArrayRead(ut,&carr);
1196:   VecPlaceArray(u,carr);
1197:   BVInsertVec(U1,k,u);
1198:   VecResetArray(u);
1199:   VecRestoreArrayRead(ut,&carr);
1200:   VecDestroy(&u);
1201:   if (breakdown) BVOrthonormalizeColumn(U1,k,PETSC_FALSE,NULL,breakdown);
1202:   else BVOrthonormalizeColumn(U1,k,PETSC_TRUE,NULL,NULL);

1204:   if (!breakdown || !*breakdown) {
1205:     MatGetLocalSize(svd->A,&m,NULL);
1206:     /* Compute k-th vector of BV V */
1207:     BVGetColumn(V,k,&v);
1208:     /* Form ut=[u;0] where u is the 1st column of U1 */
1209:     BVGetColumn(U1,k,&u);
1210:     VecZeroEntries(ut);
1211:     VecGetArrayRead(u,&carr);
1212:     VecGetArray(ut,&arr);
1213:     for (j=0; j<m; j++) arr[j] = carr[j];
1214:     VecRestoreArrayRead(u,&carr);
1215:     VecRestoreArray(ut,&arr);
1216:     /* Solve least squares problem */
1217:     KSPSolve(lanczos->ksp,ut,x);
1218:     MatMult(lanczos->Z,x,v);
1219:     BVRestoreColumn(U1,k,&u);
1220:     BVRestoreColumn(V,k,&v);
1221:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
1222:     if (breakdown) BVOrthonormalizeColumn(V,k,PETSC_FALSE,alpha+k,breakdown);
1223:     else BVOrthonormalizeColumn(V,k,PETSC_TRUE,alpha+k,NULL);
1224:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
1225:   }
1226:   return 0;
1227: }

1229: /* solve generalized problem with joint lower-upper bidiagonalization */
1230: PetscErrorCode SVDSolve_TRLanczosGLower(SVD svd,BV U1,BV U2,BV V)
1231: {
1232:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
1233:   PetscReal      *alpha,*beta,*alphah,*betah,normr,scalef,*sigma,sigma0;
1234:   PetscScalar    *w;
1235:   PetscInt       i,k,l,nv,ld;
1236:   Mat            U,Vmat,X;
1237:   PetscBool      breakdown=PETSC_FALSE,inverted;

1239:   DSGetLeadingDimension(svd->ds,&ld);
1240:   PetscMalloc2(ld,&w,ld,&sigma);
1241:   inverted = ((svd->which==SVD_LARGEST && svd->swapped) || (svd->which==SVD_SMALLEST && !svd->swapped))? PETSC_TRUE: PETSC_FALSE;
1242:   scalef = svd->swapped? 1.0/lanczos->scalef : lanczos->scalef;
1243:   normr = (svd->conv==SVD_CONV_ABS)? PetscMax(svd->nrma,svd->nrmb*scalef): 1.0;

1245:   /* normalize start vector */
1246:   if (!svd->ninil) BVSetRandomColumn(U2,0);
1247:   SVDInitialVectorGLower(svd,V,U1,U2,0,NULL);

1249:   l = 0;
1250:   while (svd->reason == SVD_CONVERGED_ITERATING) {
1251:     svd->its++;

1253:     /* inner loop */
1254:     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
1255:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
1256:     DSGetArrayReal(svd->ds,DS_MAT_D,&alphah);
1257:     beta = alpha + ld;
1258:     betah = alpha + 2*ld;
1259:     SVDLanczosGLower(svd,alpha,beta,alphah,betah,lanczos->Z,U1,U2,V,lanczos->ksp,svd->nconv+l,&nv,&breakdown);
1260:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
1261:     DSRestoreArrayReal(svd->ds,DS_MAT_D,&alphah);
1262:     BVSetActiveColumns(V,svd->nconv,nv);
1263:     BVSetActiveColumns(U1,svd->nconv,nv+1);
1264:     BVSetActiveColumns(U2,svd->nconv,nv);

1266:     /* solve projected problem */
1267:     DSSetDimensions(svd->ds,nv+1,svd->nconv,svd->nconv+l);
1268:     DSGSVDSetDimensions(svd->ds,nv,nv);
1269:     DSSetState(svd->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE);
1270:     DSSolve(svd->ds,w,NULL);
1271:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
1272:     DSUpdateExtraRow(svd->ds);
1273:     DSSynchronize(svd->ds,w,NULL);
1274:     for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);

1276:     /* check convergence */
1277:     SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,normr,&k);
1278:     (*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx);

1280:     sigma0 = inverted? 1.0/svd->sigma[0] : svd->sigma[0];
1281:     if (lanczos->scaleth!=0 && k==0 && sigma0>lanczos->scaleth) {

1283:       /* Scale and start from scratch */
1284:       lanczos->scalef *= svd->swapped? 1.0/svd->sigma[0] : svd->sigma[0];
1285:       PetscInfo(svd,"Scaling by factor %g and starting from scratch\n",(double)lanczos->scalef);
1286:       MatZUpdateScale(svd);
1287:       scalef = svd->swapped? 1.0/lanczos->scalef : lanczos->scalef;
1288:       if (svd->conv==SVD_CONV_ABS) normr = PetscMax(svd->nrma,svd->nrmb*scalef);
1289:       l = 0;
1290:       if (!svd->ninil) BVSetRandomColumn(U2,0);
1291:       SVDInitialVectorGLower(svd,V,U1,U2,0,NULL);

1293:     } else {

1295:       /* update l */
1296:       if (svd->reason != SVD_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
1297:       else l = PetscMax(1,(PetscInt)((nv-k)*lanczos->keep));
1298:       if (!lanczos->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged triplets */
1299:       if (l) PetscInfo(svd,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l);

1301:       if (svd->reason == SVD_CONVERGED_ITERATING) {
1302:         if (PetscUnlikely(breakdown || k==nv)) {
1303:           /* Start a new bidiagonalization */
1304:           PetscInfo(svd,"Breakdown in bidiagonalization (it=%" PetscInt_FMT ")\n",svd->its);
1305:           if (k<svd->nsv) {
1306:             BVSetRandomColumn(U2,k);
1307:             SVDInitialVectorGLower(svd,V,U1,U2,k,&breakdown);
1308:             if (breakdown) {
1309:               svd->reason = SVD_DIVERGED_BREAKDOWN;
1310:               PetscInfo(svd,"Unable to generate more start vectors\n");
1311:             }
1312:           }
1313:         } else DSTruncate(svd->ds,k+l,PETSC_FALSE);
1314:       }

1316:       /* compute converged singular vectors and restart vectors */
1317:       DSGetMat(svd->ds,DS_MAT_X,&X);
1318:       BVMultInPlace(V,X,svd->nconv,k+l);
1319:       DSRestoreMat(svd->ds,DS_MAT_X,&X);
1320:       DSGetMat(svd->ds,DS_MAT_U,&U);
1321:       BVMultInPlace(U1,U,svd->nconv,k+l+1);
1322:       DSRestoreMat(svd->ds,DS_MAT_U,&U);
1323:       DSGetMat(svd->ds,DS_MAT_V,&Vmat);
1324:       BVMultInPlace(U2,Vmat,svd->nconv,k+l);
1325:       DSRestoreMat(svd->ds,DS_MAT_V,&Vmat);

1327:       /* copy the last vector to be the next initial vector */
1328:       if (svd->reason == SVD_CONVERGED_ITERATING && !breakdown) BVCopyColumn(V,nv,k+l);
1329:     }

1331:     svd->nconv = k;
1332:     SVDLanczosBackTransform(svd,nv,svd->sigma,sigma,NULL);
1333:     SVDMonitor(svd,svd->its,svd->nconv,sigma,svd->errest,nv);
1334:   }

1336:   PetscFree2(w,sigma);
1337:   return 0;
1338: }

1340: PetscErrorCode SVDSolve_TRLanczos_GSVD(SVD svd)
1341: {
1342:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
1343:   PetscInt       k,m,p;
1344:   PetscBool      convchg=PETSC_FALSE;
1345:   BV             U1,U2,UU;
1346:   BVType         type;
1347:   Mat            U,V;
1348:   SlepcSC        sc;

1350:   PetscCitationsRegister(citationg,&citedg);

1352:   if (svd->swapped) {
1353:     DSGetSlepcSC(svd->ds,&sc);
1354:     if (svd->which==SVD_LARGEST) sc->comparison = SlepcCompareSmallestReal;
1355:     else                         sc->comparison = SlepcCompareLargestReal;
1356:   }
1357:   if (svd->converged==SVDConvergedNorm) {  /* override temporarily since computed residual is already relative to the norms */
1358:     svd->converged = SVDConvergedAbsolute;
1359:     convchg = PETSC_TRUE;
1360:   }
1361:   MatGetLocalSize(svd->A,&m,NULL);
1362:   MatGetLocalSize(svd->B,&p,NULL);

1364:   /* Create BV for U1 */
1365:   BVCreate(PetscObjectComm((PetscObject)svd),&U1);
1366:   BVGetType(svd->U,&type);
1367:   BVSetType(U1,type);
1368:   BVGetSizes(svd->U,NULL,NULL,&k);
1369:   BVSetSizes(U1,m,PETSC_DECIDE,k);

1371:   /* Create BV for U2 */
1372:   BVCreate(PetscObjectComm((PetscObject)svd),&U2);
1373:   BVSetType(U2,type);
1374:   BVSetSizes(U2,p,PETSC_DECIDE,k);

1376:   /* Copy initial vectors from svd->U to U1 and U2 */
1377:   if (svd->ninil) {
1378:     Vec u, uh, nest, aux[2];
1379:     BVGetColumn(U1,0,&u);
1380:     BVGetColumn(U2,0,&uh);
1381:     aux[0] = u;
1382:     aux[1] = uh;
1383:     VecCreateNest(PetscObjectComm((PetscObject)svd),2,NULL,aux,&nest);
1384:     BVCopyVec(svd->U,0,nest);
1385:     BVRestoreColumn(U1,0,&u);
1386:     BVRestoreColumn(U2,0,&uh);
1387:     VecDestroy(&nest);
1388:   }

1390:   switch (lanczos->bidiag) {
1391:     case SVD_TRLANCZOS_GBIDIAG_SINGLE:
1392:       SVDSolve_TRLanczosGSingle(svd,U1,svd->U);
1393:       break;
1394:     case SVD_TRLANCZOS_GBIDIAG_UPPER:
1395:       SVDSolve_TRLanczosGUpper(svd,U1,U2,svd->U);
1396:       break;
1397:     case SVD_TRLANCZOS_GBIDIAG_LOWER:
1398:       SVDSolve_TRLanczosGLower(svd,U1,U2,svd->U);
1399:       break;
1400:   }

1402:   /* Compute converged right singular vectors */
1403:   BVSetActiveColumns(svd->U,0,svd->nconv);
1404:   BVSetActiveColumns(svd->V,0,svd->nconv);
1405:   BVGetMat(svd->U,&U);
1406:   BVGetMat(svd->V,&V);
1407:   KSPMatSolve(lanczos->ksp,U,V);
1408:   BVRestoreMat(svd->U,&U);
1409:   BVRestoreMat(svd->V,&V);

1411:   /* Finish computing left singular vectors and move them to its place */
1412:   if (svd->swapped) SWAP(U1,U2,UU);
1413:   switch (lanczos->bidiag) {
1414:     case SVD_TRLANCZOS_GBIDIAG_SINGLE:
1415:       SVDLeftSingularVectors_Single(svd,U1,U2);
1416:       break;
1417:     case SVD_TRLANCZOS_GBIDIAG_UPPER:
1418:     case SVD_TRLANCZOS_GBIDIAG_LOWER:
1419:       SVDLeftSingularVectors(svd,U1,U2);
1420:       break;
1421:   }

1423:   /* undo scaling and compute the reciprocals of sigma if matrices were swapped */
1424:   SVDLanczosBackTransform(svd,svd->nconv,svd->sigma,NULL,svd->V);

1426:   BVDestroy(&U1);
1427:   BVDestroy(&U2);
1428:   DSTruncate(svd->ds,svd->nconv,PETSC_TRUE);
1429:   if (convchg) svd->converged = SVDConvergedNorm;
1430:   return 0;
1431: }

1433: PetscErrorCode SVDSetFromOptions_TRLanczos(SVD svd,PetscOptionItems *PetscOptionsObject)
1434: {
1435:   SVD_TRLANCZOS       *lanczos = (SVD_TRLANCZOS*)svd->data;
1436:   PetscBool           flg,val,lock;
1437:   PetscReal           keep,scale;
1438:   SVDTRLanczosGBidiag bidiag;

1440:   PetscOptionsHeadBegin(PetscOptionsObject,"SVD TRLanczos Options");

1442:     PetscOptionsBool("-svd_trlanczos_oneside","Use one-side reorthogonalization","SVDTRLanczosSetOneSide",lanczos->oneside,&val,&flg);
1443:     if (flg) SVDTRLanczosSetOneSide(svd,val);

1445:     PetscOptionsReal("-svd_trlanczos_restart","Proportion of vectors kept after restart","SVDTRLanczosSetRestart",0.5,&keep,&flg);
1446:     if (flg) SVDTRLanczosSetRestart(svd,keep);

1448:     PetscOptionsBool("-svd_trlanczos_locking","Choose between locking and non-locking variants","SVDTRLanczosSetLocking",PETSC_TRUE,&lock,&flg);
1449:     if (flg) SVDTRLanczosSetLocking(svd,lock);

1451:     PetscOptionsEnum("-svd_trlanczos_gbidiag","Bidiagonalization choice for Generalized Problem","SVDTRLanczosSetGBidiag",SVDTRLanczosGBidiags,(PetscEnum)lanczos->bidiag,(PetscEnum*)&bidiag,&flg);
1452:     if (flg) SVDTRLanczosSetGBidiag(svd,bidiag);

1454:     PetscOptionsBool("-svd_trlanczos_explicitmatrix","Build explicit matrix for KSP solver","SVDTRLanczosSetExplicitMatrix",lanczos->explicitmatrix,&val,&flg);
1455:     if (flg) SVDTRLanczosSetExplicitMatrix(svd,val);

1457:     SVDTRLanczosGetScale(svd,&scale);
1458:     PetscOptionsReal("-svd_trlanczos_scale","Scale parameter for matrix B","SVDTRLanczosSetScale",scale,&scale,&flg);
1459:     if (flg) SVDTRLanczosSetScale(svd,scale);

1461:   PetscOptionsHeadEnd();

1463:   if (svd->OPb) {
1464:     if (!lanczos->ksp) SVDTRLanczosGetKSP(svd,&lanczos->ksp);
1465:     KSPSetFromOptions(lanczos->ksp);
1466:   }
1467:   return 0;
1468: }

1470: static PetscErrorCode SVDTRLanczosSetOneSide_TRLanczos(SVD svd,PetscBool oneside)
1471: {
1472:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

1474:   if (lanczos->oneside != oneside) {
1475:     lanczos->oneside = oneside;
1476:     svd->state = SVD_STATE_INITIAL;
1477:   }
1478:   return 0;
1479: }

1481: /*@
1482:    SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method
1483:    to be used is one-sided or two-sided.

1485:    Logically Collective on svd

1487:    Input Parameters:
1488: +  svd     - singular value solver
1489: -  oneside - boolean flag indicating if the method is one-sided or not

1491:    Options Database Key:
1492: .  -svd_trlanczos_oneside <boolean> - Indicates the boolean flag

1494:    Notes:
1495:    By default, a two-sided variant is selected, which is sometimes slightly
1496:    more robust. However, the one-sided variant is faster because it avoids
1497:    the orthogonalization associated to left singular vectors.

1499:    One-sided orthogonalization is also available for the GSVD, in which case
1500:    two orthogonalizations out of three are avoided.

1502:    Level: advanced

1504: .seealso: SVDLanczosSetOneSide()
1505: @*/
1506: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscBool oneside)
1507: {
1510:   PetscTryMethod(svd,"SVDTRLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
1511:   return 0;
1512: }

1514: static PetscErrorCode SVDTRLanczosGetOneSide_TRLanczos(SVD svd,PetscBool *oneside)
1515: {
1516:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

1518:   *oneside = lanczos->oneside;
1519:   return 0;
1520: }

1522: /*@
1523:    SVDTRLanczosGetOneSide - Gets if the variant of the Lanczos method
1524:    to be used is one-sided or two-sided.

1526:    Not Collective

1528:    Input Parameters:
1529: .  svd     - singular value solver

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

1534:    Level: advanced

1536: .seealso: SVDTRLanczosSetOneSide()
1537: @*/
1538: PetscErrorCode SVDTRLanczosGetOneSide(SVD svd,PetscBool *oneside)
1539: {
1542:   PetscUseMethod(svd,"SVDTRLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
1543:   return 0;
1544: }

1546: static PetscErrorCode SVDTRLanczosSetGBidiag_TRLanczos(SVD svd,SVDTRLanczosGBidiag bidiag)
1547: {
1548:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

1550:   switch (bidiag) {
1551:     case SVD_TRLANCZOS_GBIDIAG_SINGLE:
1552:     case SVD_TRLANCZOS_GBIDIAG_UPPER:
1553:     case SVD_TRLANCZOS_GBIDIAG_LOWER:
1554:       if (lanczos->bidiag != bidiag) {
1555:         lanczos->bidiag = bidiag;
1556:         svd->state = SVD_STATE_INITIAL;
1557:       }
1558:       break;
1559:     default:
1560:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid bidiagonalization choice");
1561:   }
1562:   return 0;
1563: }

1565: /*@
1566:    SVDTRLanczosSetGBidiag - Sets the bidiagonalization choice to use in
1567:    the GSVD TRLanczos solver.

1569:    Logically Collective on svd

1571:    Input Parameters:
1572: +  svd - the singular value solver
1573: -  bidiag - the bidiagonalization choice

1575:    Options Database Key:
1576: .  -svd_trlanczos_gbidiag - Sets the bidiagonalization choice (either 's' or 'juu'
1577:    or 'jlu')

1579:    Level: advanced

1581: .seealso: SVDTRLanczosGetGBidiag(), SVDTRLanczosGBidiag
1582: @*/
1583: PetscErrorCode SVDTRLanczosSetGBidiag(SVD svd,SVDTRLanczosGBidiag bidiag)
1584: {
1587:   PetscTryMethod(svd,"SVDTRLanczosSetGBidiag_C",(SVD,SVDTRLanczosGBidiag),(svd,bidiag));
1588:   return 0;
1589: }

1591: static PetscErrorCode SVDTRLanczosGetGBidiag_TRLanczos(SVD svd,SVDTRLanczosGBidiag *bidiag)
1592: {
1593:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

1595:   *bidiag = lanczos->bidiag;
1596:   return 0;
1597: }

1599: /*@
1600:    SVDTRLanczosGetGBidiag - Gets the bidiagonalization choice used in the GSVD
1601:    TRLanczos solver.

1603:    Not Collective

1605:    Input Parameter:
1606: .  svd - the singular value solver

1608:    Output Parameter:
1609: .  bidiag - the bidiagonalization choice

1611:    Level: advanced

1613: .seealso: SVDTRLanczosSetGBidiag(), SVDTRLanczosGBidiag
1614: @*/
1615: PetscErrorCode SVDTRLanczosGetGBidiag(SVD svd,SVDTRLanczosGBidiag *bidiag)
1616: {
1619:   PetscUseMethod(svd,"SVDTRLanczosGetGBidiag_C",(SVD,SVDTRLanczosGBidiag*),(svd,bidiag));
1620:   return 0;
1621: }

1623: static PetscErrorCode SVDTRLanczosSetKSP_TRLanczos(SVD svd,KSP ksp)
1624: {
1625:   SVD_TRLANCZOS  *ctx = (SVD_TRLANCZOS*)svd->data;

1627:   PetscObjectReference((PetscObject)ksp);
1628:   KSPDestroy(&ctx->ksp);
1629:   ctx->ksp   = ksp;
1630:   svd->state = SVD_STATE_INITIAL;
1631:   return 0;
1632: }

1634: /*@
1635:    SVDTRLanczosSetKSP - Associate a linear solver object (KSP) to the SVD solver.

1637:    Collective on svd

1639:    Input Parameters:
1640: +  svd - SVD solver
1641: -  ksp - the linear solver object

1643:    Note:
1644:    Only used for the GSVD problem.

1646:    Level: advanced

1648: .seealso: SVDTRLanczosGetKSP()
1649: @*/
1650: PetscErrorCode SVDTRLanczosSetKSP(SVD svd,KSP ksp)
1651: {
1655:   PetscTryMethod(svd,"SVDTRLanczosSetKSP_C",(SVD,KSP),(svd,ksp));
1656:   return 0;
1657: }

1659: static PetscErrorCode SVDTRLanczosGetKSP_TRLanczos(SVD svd,KSP *ksp)
1660: {
1661:   SVD_TRLANCZOS  *ctx = (SVD_TRLANCZOS*)svd->data;
1662:   PC             pc;

1664:   if (!ctx->ksp) {
1665:     /* Create linear solver */
1666:     KSPCreate(PetscObjectComm((PetscObject)svd),&ctx->ksp);
1667:     PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)svd,1);
1668:     KSPSetOptionsPrefix(ctx->ksp,((PetscObject)svd)->prefix);
1669:     KSPAppendOptionsPrefix(ctx->ksp,"svd_trlanczos_");
1670:     PetscObjectSetOptions((PetscObject)ctx->ksp,((PetscObject)svd)->options);
1671:     KSPSetType(ctx->ksp,KSPLSQR);
1672:     KSPGetPC(ctx->ksp,&pc);
1673:     PCSetType(pc,PCNONE);
1674:     KSPSetErrorIfNotConverged(ctx->ksp,PETSC_TRUE);
1675:     KSPSetTolerances(ctx->ksp,SlepcDefaultTol(svd->tol)/10.0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
1676:   }
1677:   *ksp = ctx->ksp;
1678:   return 0;
1679: }

1681: /*@
1682:    SVDTRLanczosGetKSP - Retrieve the linear solver object (KSP) associated with
1683:    the SVD solver.

1685:    Not Collective

1687:    Input Parameter:
1688: .  svd - SVD solver

1690:    Output Parameter:
1691: .  ksp - the linear solver object

1693:    Level: advanced

1695: .seealso: SVDTRLanczosSetKSP()
1696: @*/
1697: PetscErrorCode SVDTRLanczosGetKSP(SVD svd,KSP *ksp)
1698: {
1701:   PetscUseMethod(svd,"SVDTRLanczosGetKSP_C",(SVD,KSP*),(svd,ksp));
1702:   return 0;
1703: }

1705: static PetscErrorCode SVDTRLanczosSetRestart_TRLanczos(SVD svd,PetscReal keep)
1706: {
1707:   SVD_TRLANCZOS *ctx = (SVD_TRLANCZOS*)svd->data;

1709:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
1710:   else {
1712:     ctx->keep = keep;
1713:   }
1714:   return 0;
1715: }

1717: /*@
1718:    SVDTRLanczosSetRestart - Sets the restart parameter for the thick-restart
1719:    Lanczos method, in particular the proportion of basis vectors that must be
1720:    kept after restart.

1722:    Logically Collective on svd

1724:    Input Parameters:
1725: +  svd  - the singular value solver
1726: -  keep - the number of vectors to be kept at restart

1728:    Options Database Key:
1729: .  -svd_trlanczos_restart - Sets the restart parameter

1731:    Notes:
1732:    Allowed values are in the range [0.1,0.9]. The default is 0.5.

1734:    Level: advanced

1736: .seealso: SVDTRLanczosGetRestart()
1737: @*/
1738: PetscErrorCode SVDTRLanczosSetRestart(SVD svd,PetscReal keep)
1739: {
1742:   PetscTryMethod(svd,"SVDTRLanczosSetRestart_C",(SVD,PetscReal),(svd,keep));
1743:   return 0;
1744: }

1746: static PetscErrorCode SVDTRLanczosGetRestart_TRLanczos(SVD svd,PetscReal *keep)
1747: {
1748:   SVD_TRLANCZOS *ctx = (SVD_TRLANCZOS*)svd->data;

1750:   *keep = ctx->keep;
1751:   return 0;
1752: }

1754: /*@
1755:    SVDTRLanczosGetRestart - Gets the restart parameter used in the thick-restart
1756:    Lanczos method.

1758:    Not Collective

1760:    Input Parameter:
1761: .  svd - the singular value solver

1763:    Output Parameter:
1764: .  keep - the restart parameter

1766:    Level: advanced

1768: .seealso: SVDTRLanczosSetRestart()
1769: @*/
1770: PetscErrorCode SVDTRLanczosGetRestart(SVD svd,PetscReal *keep)
1771: {
1774:   PetscUseMethod(svd,"SVDTRLanczosGetRestart_C",(SVD,PetscReal*),(svd,keep));
1775:   return 0;
1776: }

1778: static PetscErrorCode SVDTRLanczosSetLocking_TRLanczos(SVD svd,PetscBool lock)
1779: {
1780:   SVD_TRLANCZOS *ctx = (SVD_TRLANCZOS*)svd->data;

1782:   ctx->lock = lock;
1783:   return 0;
1784: }

1786: /*@
1787:    SVDTRLanczosSetLocking - Choose between locking and non-locking variants of
1788:    the thick-restart Lanczos method.

1790:    Logically Collective on svd

1792:    Input Parameters:
1793: +  svd  - the singular value solver
1794: -  lock - true if the locking variant must be selected

1796:    Options Database Key:
1797: .  -svd_trlanczos_locking - Sets the locking flag

1799:    Notes:
1800:    The default is to lock converged singular triplets when the method restarts.
1801:    This behaviour can be changed so that all directions are kept in the
1802:    working subspace even if already converged to working accuracy (the
1803:    non-locking variant).

1805:    Level: advanced

1807: .seealso: SVDTRLanczosGetLocking()
1808: @*/
1809: PetscErrorCode SVDTRLanczosSetLocking(SVD svd,PetscBool lock)
1810: {
1813:   PetscTryMethod(svd,"SVDTRLanczosSetLocking_C",(SVD,PetscBool),(svd,lock));
1814:   return 0;
1815: }

1817: static PetscErrorCode SVDTRLanczosGetLocking_TRLanczos(SVD svd,PetscBool *lock)
1818: {
1819:   SVD_TRLANCZOS *ctx = (SVD_TRLANCZOS*)svd->data;

1821:   *lock = ctx->lock;
1822:   return 0;
1823: }

1825: /*@
1826:    SVDTRLanczosGetLocking - Gets the locking flag used in the thick-restart
1827:    Lanczos method.

1829:    Not Collective

1831:    Input Parameter:
1832: .  svd - the singular value solver

1834:    Output Parameter:
1835: .  lock - the locking flag

1837:    Level: advanced

1839: .seealso: SVDTRLanczosSetLocking()
1840: @*/
1841: PetscErrorCode SVDTRLanczosGetLocking(SVD svd,PetscBool *lock)
1842: {
1845:   PetscUseMethod(svd,"SVDTRLanczosGetLocking_C",(SVD,PetscBool*),(svd,lock));
1846:   return 0;
1847: }

1849: static PetscErrorCode SVDTRLanczosSetExplicitMatrix_TRLanczos(SVD svd,PetscBool explicitmat)
1850: {
1851:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

1853:   if (lanczos->explicitmatrix != explicitmat) {
1854:     lanczos->explicitmatrix = explicitmat;
1855:     svd->state = SVD_STATE_INITIAL;
1856:   }
1857:   return 0;
1858: }

1860: /*@
1861:    SVDTRLanczosSetExplicitMatrix - Indicate if the matrix Z=[A;B] must
1862:    be built explicitly.

1864:    Logically Collective on svd

1866:    Input Parameters:
1867: +  svd         - singular value solver
1868: -  explicitmat - Boolean flag indicating if Z=[A;B] is built explicitly

1870:    Options Database Key:
1871: .  -svd_trlanczos_explicitmatrix <boolean> - Indicates the boolean flag

1873:    Notes:
1874:    This option is relevant for the GSVD case only.
1875:    Z is the coefficient matrix of the KSP solver used internally.

1877:    Level: advanced

1879: .seealso: SVDTRLanczosGetExplicitMatrix()
1880: @*/
1881: PetscErrorCode SVDTRLanczosSetExplicitMatrix(SVD svd,PetscBool explicitmat)
1882: {
1885:   PetscTryMethod(svd,"SVDTRLanczosSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmat));
1886:   return 0;
1887: }

1889: static PetscErrorCode SVDTRLanczosGetExplicitMatrix_TRLanczos(SVD svd,PetscBool *explicitmat)
1890: {
1891:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

1893:   *explicitmat = lanczos->explicitmatrix;
1894:   return 0;
1895: }

1897: /*@
1898:    SVDTRLanczosGetExplicitMatrix - Returns the flag indicating if Z=[A;B] is built explicitly.

1900:    Not Collective

1902:    Input Parameter:
1903: .  svd  - singular value solver

1905:    Output Parameter:
1906: .  explicitmat - the mode flag

1908:    Level: advanced

1910: .seealso: SVDTRLanczosSetExplicitMatrix()
1911: @*/
1912: PetscErrorCode SVDTRLanczosGetExplicitMatrix(SVD svd,PetscBool *explicitmat)
1913: {
1916:   PetscUseMethod(svd,"SVDTRLanczosGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmat));
1917:   return 0;
1918: }

1920: static PetscErrorCode SVDTRLanczosSetScale_TRLanczos(SVD svd,PetscReal scale)
1921: {
1922:   SVD_TRLANCZOS *ctx = (SVD_TRLANCZOS*)svd->data;

1924:   if (scale<0) {
1925:     ctx->scalef  = 1.0;
1926:     ctx->scaleth = -scale;
1927:   } else {
1928:     ctx->scalef  = scale;
1929:     ctx->scaleth = 0.0;
1930:   }
1931:   return 0;
1932: }

1934: /*@
1935:    SVDTRLanczosSetScale - Sets the scale parameter for the GSVD.

1937:    Logically Collective on svd

1939:    Input Parameters:
1940: +  svd   - singular value solver
1941: -  scale - scale parameter

1943:    Options Database Key:
1944: .  -svd_trlanczos_scale <real> - scale factor/threshold

1946:    Notes:
1947:    This parameter is relevant for the GSVD case only. If the parameter is
1948:    positive, it indicates the scale factor for B in matrix Z=[A;B]. If
1949:    negative, its absolute value is the threshold for automatic scaling.
1950:    In automatic scaling, whenever the largest approximate generalized singular
1951:    value (or the inverse of the smallest value, if SVD_SMALLEST is used)
1952:    exceeds the threshold, the computation is restarted with matrix B
1953:    scaled by that value.

1955:    Level: advanced

1957: .seealso: SVDTRLanczosGetScale()
1958: @*/
1959: PetscErrorCode SVDTRLanczosSetScale(SVD svd,PetscReal scale)
1960: {
1963:   PetscTryMethod(svd,"SVDTRLanczosSetScale_C",(SVD,PetscReal),(svd,scale));
1964:   return 0;
1965: }

1967: static PetscErrorCode SVDTRLanczosGetScale_TRLanczos(SVD svd,PetscReal *scale)
1968: {
1969:   SVD_TRLANCZOS *ctx = (SVD_TRLANCZOS*)svd->data;

1971:   if (ctx->scaleth==0) *scale = ctx->scalef;
1972:   else                 *scale = -ctx->scaleth;
1973:   return 0;
1974: }

1976: /*@
1977:    SVDTRLanczosGetScale - Gets the scale parameter for the GSVD.

1979:    Not Collective

1981:    Input Parameter:
1982: .  svd - the singular value solver

1984:    Output Parameter:
1985: .  scale - the scale parameter

1987:    Notes:
1988:    This parameter is relevant for the GSVD case only. If the parameter is
1989:    positive, it indicates the scale factor for B in matrix Z=[A;B]. If
1990:    negative, its absolute value is the threshold for automatic scaling.

1992:    Level: advanced

1994: .seealso: SVDTRLanczosSetScale()
1995: @*/
1996: PetscErrorCode SVDTRLanczosGetScale(SVD svd,PetscReal *scale)
1997: {
2000:   PetscUseMethod(svd,"SVDTRLanczosGetScale_C",(SVD,PetscReal*),(svd,scale));
2001:   return 0;
2002: }

2004: PetscErrorCode SVDReset_TRLanczos(SVD svd)
2005: {
2006:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;

2008:   if (svd->isgeneralized || (!svd->problem_type && svd->OPb)) {
2009:     KSPReset(lanczos->ksp);
2010:     MatDestroy(&lanczos->Z);
2011:   }
2012:   return 0;
2013: }

2015: PetscErrorCode SVDDestroy_TRLanczos(SVD svd)
2016: {
2017:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;

2019:   if (svd->isgeneralized || (!svd->problem_type && svd->OPb)) KSPDestroy(&lanczos->ksp);
2020:   PetscFree(svd->data);
2021:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",NULL);
2022:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",NULL);
2023:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetGBidiag_C",NULL);
2024:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetGBidiag_C",NULL);
2025:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetKSP_C",NULL);
2026:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetKSP_C",NULL);
2027:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetRestart_C",NULL);
2028:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetRestart_C",NULL);
2029:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetLocking_C",NULL);
2030:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetLocking_C",NULL);
2031:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetExplicitMatrix_C",NULL);
2032:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetExplicitMatrix_C",NULL);
2033:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetScale_C",NULL);
2034:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetScale_C",NULL);
2035:   return 0;
2036: }

2038: PetscErrorCode SVDView_TRLanczos(SVD svd,PetscViewer viewer)
2039: {
2040:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
2041:   PetscBool      isascii;

2043:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
2044:   if (isascii) {
2045:     PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*lanczos->keep));
2046:     PetscViewerASCIIPrintf(viewer,"  using the %slocking variant\n",lanczos->lock?"":"non-");
2047:     if (svd->isgeneralized) {
2048:       const char *bidiag="";

2050:       switch (lanczos->bidiag) {
2051:         case SVD_TRLANCZOS_GBIDIAG_SINGLE: bidiag = "single"; break;
2052:         case SVD_TRLANCZOS_GBIDIAG_UPPER:  bidiag = "joint upper-upper"; break;
2053:         case SVD_TRLANCZOS_GBIDIAG_LOWER:  bidiag = "joint lower-upper"; break;
2054:       }
2055:       PetscViewerASCIIPrintf(viewer,"  bidiagonalization choice: %s\n",bidiag);
2056:       PetscViewerASCIIPrintf(viewer,"  %s matrix\n",lanczos->explicitmatrix?"explicit":"implicit");
2057:       if (lanczos->scaleth==0) PetscViewerASCIIPrintf(viewer,"  scale factor for matrix B: %g\n",(double)lanczos->scalef);
2058:       else PetscViewerASCIIPrintf(viewer,"  automatic scaling for matrix B with threshold: %g\n",(double)lanczos->scaleth);
2059:       if (!lanczos->ksp) SVDTRLanczosGetKSP(svd,&lanczos->ksp);
2060:       PetscViewerASCIIPushTab(viewer);
2061:       KSPView(lanczos->ksp,viewer);
2062:       PetscViewerASCIIPopTab(viewer);
2063:     } else PetscViewerASCIIPrintf(viewer,"  %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
2064:   }
2065:   return 0;
2066: }

2068: SLEPC_EXTERN PetscErrorCode SVDCreate_TRLanczos(SVD svd)
2069: {
2070:   SVD_TRLANCZOS  *ctx;

2072:   PetscNew(&ctx);
2073:   svd->data = (void*)ctx;

2075:   ctx->lock    = PETSC_TRUE;
2076:   ctx->bidiag  = SVD_TRLANCZOS_GBIDIAG_LOWER;
2077:   ctx->scalef  = 1.0;
2078:   ctx->scaleth = 0.0;

2080:   svd->ops->setup          = SVDSetUp_TRLanczos;
2081:   svd->ops->solve          = SVDSolve_TRLanczos;
2082:   svd->ops->solveg         = SVDSolve_TRLanczos_GSVD;
2083:   svd->ops->destroy        = SVDDestroy_TRLanczos;
2084:   svd->ops->reset          = SVDReset_TRLanczos;
2085:   svd->ops->setfromoptions = SVDSetFromOptions_TRLanczos;
2086:   svd->ops->view           = SVDView_TRLanczos;
2087:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",SVDTRLanczosSetOneSide_TRLanczos);
2088:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",SVDTRLanczosGetOneSide_TRLanczos);
2089:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetGBidiag_C",SVDTRLanczosSetGBidiag_TRLanczos);
2090:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetGBidiag_C",SVDTRLanczosGetGBidiag_TRLanczos);
2091:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetKSP_C",SVDTRLanczosSetKSP_TRLanczos);
2092:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetKSP_C",SVDTRLanczosGetKSP_TRLanczos);
2093:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetRestart_C",SVDTRLanczosSetRestart_TRLanczos);
2094:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetRestart_C",SVDTRLanczosGetRestart_TRLanczos);
2095:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetLocking_C",SVDTRLanczosSetLocking_TRLanczos);
2096:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetLocking_C",SVDTRLanczosGetLocking_TRLanczos);
2097:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetExplicitMatrix_C",SVDTRLanczosSetExplicitMatrix_TRLanczos);
2098:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetExplicitMatrix_C",SVDTRLanczosGetExplicitMatrix_TRLanczos);
2099:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetScale_C",SVDTRLanczosSetScale_TRLanczos);
2100:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetScale_C",SVDTRLanczosGetScale_TRLanczos);
2101:   return 0;
2102: }