LCOV - code coverage report
Current view: top level - svd/impls/lanczos - gklanczos.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 236 252 93.7 %
Date: 2024-11-21 00:34:55 Functions: 12 13 92.3 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       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: "lanczos"
      12             : 
      13             :    Method: Explicitly restarted Lanczos
      14             : 
      15             :    Algorithm:
      16             : 
      17             :        Golub-Kahan-Lanczos bidiagonalization with explicit restart.
      18             : 
      19             :    References:
      20             : 
      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.
      24             : 
      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             : */
      30             : 
      31             : #include <slepc/private/svdimpl.h>                /*I "slepcsvd.h" I*/
      32             : 
      33             : typedef struct {
      34             :   PetscBool oneside;
      35             : } SVD_LANCZOS;
      36             : 
      37          19 : static PetscErrorCode SVDSetUp_Lanczos(SVD svd)
      38             : {
      39          19 :   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;
      40          19 :   PetscInt       N;
      41             : 
      42          19 :   PetscFunctionBegin;
      43          19 :   SVDCheckStandard(svd);
      44          19 :   SVDCheckDefinite(svd);
      45          19 :   PetscCall(MatGetSize(svd->A,NULL,&N));
      46          19 :   PetscCall(SVDSetDimensions_Default(svd));
      47          19 :   PetscCheck(svd->ncv<=svd->nsv+svd->mpd,PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"The value of ncv must not be larger than nev+mpd");
      48          19 :   if (svd->max_it==PETSC_DETERMINE) svd->max_it = PetscMax(N/svd->ncv,100);
      49          19 :   svd->leftbasis = PetscNot(lanczos->oneside);
      50          19 :   PetscCall(SVDAllocateSolution(svd,1));
      51          19 :   PetscCall(DSSetType(svd->ds,DSSVD));
      52          19 :   PetscCall(DSSetCompact(svd->ds,PETSC_TRUE));
      53          19 :   PetscCall(DSSetExtraRow(svd->ds,PETSC_TRUE));
      54          19 :   PetscCall(DSAllocate(svd->ds,svd->ncv+1));
      55          19 :   PetscFunctionReturn(PETSC_SUCCESS);
      56             : }
      57             : 
      58         966 : PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt *n,PetscBool *breakdown)
      59             : {
      60         966 :   PetscInt       i;
      61         966 :   Vec            u,v;
      62         966 :   PetscBool      lindep=PETSC_FALSE;
      63             : 
      64         966 :   PetscFunctionBegin;
      65         966 :   PetscCall(BVGetColumn(svd->V,k,&v));
      66         966 :   PetscCall(BVGetColumn(svd->U,k,&u));
      67         966 :   PetscCall(MatMult(svd->A,v,u));
      68         966 :   PetscCall(BVRestoreColumn(svd->V,k,&v));
      69         966 :   PetscCall(BVRestoreColumn(svd->U,k,&u));
      70         966 :   PetscCall(BVOrthonormalizeColumn(svd->U,k,PETSC_FALSE,alpha+k,&lindep));
      71         966 :   if (PetscUnlikely(lindep)) {
      72           0 :     *n = k;
      73           0 :     if (breakdown) *breakdown = lindep;
      74           0 :     PetscFunctionReturn(PETSC_SUCCESS);
      75             :   }
      76             : 
      77        7342 :   for (i=k+1;i<*n;i++) {
      78        6376 :     PetscCall(BVGetColumn(svd->V,i,&v));
      79        6376 :     PetscCall(BVGetColumn(svd->U,i-1,&u));
      80        6376 :     PetscCall(MatMult(svd->AT,u,v));
      81        6376 :     PetscCall(BVRestoreColumn(svd->V,i,&v));
      82        6376 :     PetscCall(BVRestoreColumn(svd->U,i-1,&u));
      83        6376 :     PetscCall(BVOrthonormalizeColumn(svd->V,i,PETSC_FALSE,beta+i-1,&lindep));
      84        6376 :     if (PetscUnlikely(lindep)) {
      85           0 :       *n = i;
      86           0 :       break;
      87             :     }
      88        6376 :     PetscCall(BVGetColumn(svd->V,i,&v));
      89        6376 :     PetscCall(BVGetColumn(svd->U,i,&u));
      90        6376 :     PetscCall(MatMult(svd->A,v,u));
      91        6376 :     PetscCall(BVRestoreColumn(svd->V,i,&v));
      92        6376 :     PetscCall(BVRestoreColumn(svd->U,i,&u));
      93        6376 :     PetscCall(BVOrthonormalizeColumn(svd->U,i,PETSC_FALSE,alpha+i,&lindep));
      94        6376 :     if (PetscUnlikely(lindep)) {
      95           0 :       *n = i;
      96           0 :       break;
      97             :     }
      98             :   }
      99             : 
     100         966 :   if (!lindep) {
     101         966 :     PetscCall(BVGetColumn(svd->V,*n,&v));
     102         966 :     PetscCall(BVGetColumn(svd->U,*n-1,&u));
     103         966 :     PetscCall(MatMult(svd->AT,u,v));
     104         966 :     PetscCall(BVRestoreColumn(svd->V,*n,&v));
     105         966 :     PetscCall(BVRestoreColumn(svd->U,*n-1,&u));
     106         966 :     PetscCall(BVOrthogonalizeColumn(svd->V,*n,NULL,beta+*n-1,&lindep));
     107             :   }
     108         966 :   if (breakdown) *breakdown = lindep;
     109         966 :   PetscFunctionReturn(PETSC_SUCCESS);
     110             : }
     111             : 
     112          54 : static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
     113             : {
     114          54 :   PetscInt       i,bvl,bvk;
     115          54 :   PetscReal      a,b;
     116          54 :   Vec            z,temp;
     117             : 
     118          54 :   PetscFunctionBegin;
     119          54 :   PetscCall(BVGetActiveColumns(V,&bvl,&bvk));
     120          54 :   PetscCall(BVGetColumn(V,k,&z));
     121          54 :   PetscCall(MatMult(svd->A,z,u));
     122          54 :   PetscCall(BVRestoreColumn(V,k,&z));
     123             : 
     124         468 :   for (i=k+1;i<n;i++) {
     125         414 :     PetscCall(BVGetColumn(V,i,&z));
     126         414 :     PetscCall(MatMult(svd->AT,u,z));
     127         414 :     PetscCall(BVRestoreColumn(V,i,&z));
     128         414 :     PetscCall(VecNormBegin(u,NORM_2,&a));
     129         414 :     PetscCall(BVSetActiveColumns(V,0,i));
     130         414 :     PetscCall(BVDotColumnBegin(V,i,work));
     131         414 :     PetscCall(VecNormEnd(u,NORM_2,&a));
     132         414 :     PetscCall(BVDotColumnEnd(V,i,work));
     133         414 :     PetscCall(VecScale(u,1.0/a));
     134         414 :     PetscCall(BVMultColumn(V,-1.0/a,1.0/a,i,work));
     135             : 
     136             :     /* h = V^* z, z = z - V h  */
     137         414 :     PetscCall(BVDotColumn(V,i,work));
     138         414 :     PetscCall(BVMultColumn(V,-1.0,1.0,i,work));
     139         414 :     PetscCall(BVNormColumn(V,i,NORM_2,&b));
     140         414 :     PetscCheck(PetscAbsReal(b)>10*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)svd),PETSC_ERR_PLIB,"Recurrence generated a zero vector; use a two-sided variant");
     141         414 :     PetscCall(BVScaleColumn(V,i,1.0/b));
     142             : 
     143         414 :     PetscCall(BVGetColumn(V,i,&z));
     144         414 :     PetscCall(MatMult(svd->A,z,u_1));
     145         414 :     PetscCall(BVRestoreColumn(V,i,&z));
     146         414 :     PetscCall(VecAXPY(u_1,-b,u));
     147         414 :     alpha[i-1] = a;
     148         414 :     beta[i-1] = b;
     149         414 :     temp = u;
     150         414 :     u = u_1;
     151         414 :     u_1 = temp;
     152             :   }
     153             : 
     154          54 :   PetscCall(BVGetColumn(V,n,&z));
     155          54 :   PetscCall(MatMult(svd->AT,u,z));
     156          54 :   PetscCall(BVRestoreColumn(V,n,&z));
     157          54 :   PetscCall(VecNormBegin(u,NORM_2,&a));
     158          54 :   PetscCall(BVDotColumnBegin(V,n,work));
     159          54 :   PetscCall(VecNormEnd(u,NORM_2,&a));
     160          54 :   PetscCall(BVDotColumnEnd(V,n,work));
     161          54 :   PetscCall(VecScale(u,1.0/a));
     162          54 :   PetscCall(BVMultColumn(V,-1.0/a,1.0/a,n,work));
     163             : 
     164             :   /* h = V^* z, z = z - V h  */
     165          54 :   PetscCall(BVDotColumn(V,n,work));
     166          54 :   PetscCall(BVMultColumn(V,-1.0,1.0,n,work));
     167          54 :   PetscCall(BVNormColumn(V,i,NORM_2,&b));
     168             : 
     169          54 :   alpha[n-1] = a;
     170          54 :   beta[n-1] = b;
     171          54 :   PetscCall(BVSetActiveColumns(V,bvl,bvk));
     172          54 :   PetscFunctionReturn(PETSC_SUCCESS);
     173             : }
     174             : 
     175             : /*
     176             :    SVDKrylovConvergence - Implements the loop that checks for convergence
     177             :    in Krylov methods.
     178             : 
     179             :    Input Parameters:
     180             :      svd     - the solver; some error estimates are updated in svd->errest
     181             :      getall  - whether all residuals must be computed
     182             :      kini    - initial value of k (the loop variable)
     183             :      nits    - number of iterations of the loop
     184             :      normr   - norm of triangular factor of qr([A;B]), used only in GSVD
     185             : 
     186             :    Output Parameter:
     187             :      kout  - the first index where the convergence test failed
     188             : */
     189        2271 : PetscErrorCode SVDKrylovConvergence(SVD svd,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal normr,PetscInt *kout)
     190             : {
     191        2271 :   PetscInt       k,marker,ld;
     192        2271 :   PetscReal      *alpha,*beta,*betah,resnorm;
     193        2271 :   PetscBool      extra;
     194             : 
     195        2271 :   PetscFunctionBegin;
     196        2271 :   if (PetscUnlikely(svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it)) *kout = svd->nsv;
     197             :   else {
     198        2271 :     PetscCall(DSGetLeadingDimension(svd->ds,&ld));
     199        2271 :     PetscCall(DSGetExtraRow(svd->ds,&extra));
     200        2271 :     PetscCheck(extra,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Only implemented for DS with extra row");
     201        2271 :     marker = -1;
     202        2271 :     if (svd->trackall) getall = PETSC_TRUE;
     203        2271 :     PetscCall(DSGetArrayReal(svd->ds,DS_MAT_T,&alpha));
     204        2271 :     beta = alpha + ld;
     205        2271 :     betah = alpha + 2*ld;
     206        3004 :     for (k=kini;k<kini+nits;k++) {
     207        2999 :       if (svd->isgeneralized) resnorm = SlepcAbs(beta[k],betah[k])*normr;
     208        1874 :       else resnorm = PetscAbsReal(beta[k]);
     209        2999 :       PetscCall((*svd->converged)(svd,svd->sigma[k],resnorm,&svd->errest[k],svd->convergedctx));
     210        2999 :       if (marker==-1 && svd->errest[k] >= svd->tol) marker = k;
     211        2999 :       if (marker!=-1 && !getall) break;
     212             :     }
     213        2271 :     PetscCall(DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha));
     214        2271 :     if (marker!=-1) k = marker;
     215        2271 :     *kout = k;
     216             :   }
     217        2271 :   PetscFunctionReturn(PETSC_SUCCESS);
     218             : }
     219             : 
     220          19 : static PetscErrorCode SVDSolve_Lanczos(SVD svd)
     221             : {
     222          19 :   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;
     223          19 :   PetscReal      *alpha,*beta;
     224          19 :   PetscScalar    *swork,*w,*P;
     225          19 :   PetscInt       i,k,j,nv,ld;
     226          19 :   Vec            u=NULL,u_1=NULL;
     227          19 :   Mat            U,V;
     228             : 
     229          19 :   PetscFunctionBegin;
     230             :   /* allocate working space */
     231          19 :   PetscCall(DSGetLeadingDimension(svd->ds,&ld));
     232          19 :   PetscCall(PetscMalloc2(ld,&w,svd->ncv,&swork));
     233             : 
     234          19 :   if (lanczos->oneside) {
     235           2 :     PetscCall(MatCreateVecs(svd->A,NULL,&u));
     236           2 :     PetscCall(MatCreateVecs(svd->A,NULL,&u_1));
     237             :   }
     238             : 
     239             :   /* normalize start vector */
     240          19 :   if (!svd->nini) {
     241          15 :     PetscCall(BVSetRandomColumn(svd->V,0));
     242          15 :     PetscCall(BVOrthonormalizeColumn(svd->V,0,PETSC_TRUE,NULL,NULL));
     243             :   }
     244             : 
     245         603 :   while (svd->reason == SVD_CONVERGED_ITERATING) {
     246         584 :     svd->its++;
     247             : 
     248             :     /* inner loop */
     249         584 :     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
     250         584 :     PetscCall(DSGetArrayReal(svd->ds,DS_MAT_T,&alpha));
     251         584 :     beta = alpha + ld;
     252         584 :     if (lanczos->oneside) PetscCall(SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork));
     253             :     else {
     254         530 :       PetscCall(SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,&nv,NULL));
     255         530 :       PetscCall(BVSetActiveColumns(svd->U,svd->nconv,nv));
     256             :     }
     257         584 :     PetscCall(DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha));
     258         584 :     PetscCall(BVSetActiveColumns(svd->V,svd->nconv,nv));
     259             : 
     260             :     /* compute SVD of bidiagonal matrix */
     261         584 :     PetscCall(DSSetDimensions(svd->ds,nv,svd->nconv,0));
     262         584 :     PetscCall(DSSVDSetDimensions(svd->ds,nv));
     263         584 :     PetscCall(DSSetState(svd->ds,DS_STATE_INTERMEDIATE));
     264         584 :     PetscCall(DSSolve(svd->ds,w,NULL));
     265         584 :     PetscCall(DSSort(svd->ds,w,NULL,NULL,NULL,NULL));
     266         584 :     PetscCall(DSUpdateExtraRow(svd->ds));
     267         584 :     PetscCall(DSSynchronize(svd->ds,w,NULL));
     268        6073 :     for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);
     269             : 
     270             :     /* check convergence */
     271         584 :     PetscCall(SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,1.0,&k));
     272         584 :     PetscCall((*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx));
     273             : 
     274             :     /* compute restart vector */
     275         584 :     if (svd->reason == SVD_CONVERGED_ITERATING) {
     276         565 :       if (k<nv) {
     277         565 :         PetscCall(DSGetArray(svd->ds,DS_MAT_V,&P));
     278        5889 :         for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PetscConj(P[j+k*ld]);
     279         565 :         PetscCall(DSRestoreArray(svd->ds,DS_MAT_V,&P));
     280         565 :         PetscCall(BVMultColumn(svd->V,1.0,0.0,nv,swork));
     281             :       } else {
     282             :         /* all approximations have converged, generate a new initial vector */
     283           0 :         PetscCall(BVSetRandomColumn(svd->V,nv));
     284           0 :         PetscCall(BVOrthonormalizeColumn(svd->V,nv,PETSC_FALSE,NULL,NULL));
     285             :       }
     286             :     }
     287             : 
     288             :     /* compute converged singular vectors */
     289         584 :     PetscCall(DSGetMat(svd->ds,DS_MAT_V,&V));
     290         584 :     PetscCall(BVMultInPlace(svd->V,V,svd->nconv,k));
     291         584 :     PetscCall(DSRestoreMat(svd->ds,DS_MAT_V,&V));
     292         584 :     if (!lanczos->oneside) {
     293         530 :       PetscCall(DSGetMat(svd->ds,DS_MAT_U,&U));
     294         530 :       PetscCall(BVMultInPlace(svd->U,U,svd->nconv,k));
     295         530 :       PetscCall(DSRestoreMat(svd->ds,DS_MAT_U,&U));
     296             :     }
     297             : 
     298             :     /* copy restart vector from the last column */
     299         584 :     if (svd->reason == SVD_CONVERGED_ITERATING) PetscCall(BVCopyColumn(svd->V,nv,k));
     300             : 
     301         584 :     svd->nconv = k;
     302         603 :     PetscCall(SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv));
     303             :   }
     304             : 
     305             :   /* free working space */
     306          19 :   PetscCall(VecDestroy(&u));
     307          19 :   PetscCall(VecDestroy(&u_1));
     308          19 :   PetscCall(PetscFree2(w,swork));
     309          19 :   PetscFunctionReturn(PETSC_SUCCESS);
     310             : }
     311             : 
     312          13 : static PetscErrorCode SVDSetFromOptions_Lanczos(SVD svd,PetscOptionItems *PetscOptionsObject)
     313             : {
     314          13 :   PetscBool      set,val;
     315          13 :   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;
     316             : 
     317          13 :   PetscFunctionBegin;
     318          13 :   PetscOptionsHeadBegin(PetscOptionsObject,"SVD Lanczos Options");
     319             : 
     320          13 :     PetscCall(PetscOptionsBool("-svd_lanczos_oneside","Use one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set));
     321          13 :     if (set) PetscCall(SVDLanczosSetOneSide(svd,val));
     322             : 
     323          13 :   PetscOptionsHeadEnd();
     324          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     325             : }
     326             : 
     327           2 : static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
     328             : {
     329           2 :   SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
     330             : 
     331           2 :   PetscFunctionBegin;
     332           2 :   if (lanczos->oneside != oneside) {
     333           2 :     lanczos->oneside = oneside;
     334           2 :     svd->state = SVD_STATE_INITIAL;
     335             :   }
     336           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     337             : }
     338             : 
     339             : /*@
     340             :    SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
     341             :    to be used is one-sided or two-sided.
     342             : 
     343             :    Logically Collective
     344             : 
     345             :    Input Parameters:
     346             : +  svd     - singular value solver
     347             : -  oneside - boolean flag indicating if the method is one-sided or not
     348             : 
     349             :    Options Database Key:
     350             : .  -svd_lanczos_oneside <boolean> - Indicates the boolean flag
     351             : 
     352             :    Note:
     353             :    By default, a two-sided variant is selected, which is sometimes slightly
     354             :    more robust. However, the one-sided variant is faster because it avoids
     355             :    the orthogonalization associated to left singular vectors. It also saves
     356             :    the memory required for storing such vectors.
     357             : 
     358             :    Level: advanced
     359             : 
     360             : .seealso: SVDTRLanczosSetOneSide()
     361             : @*/
     362           2 : PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
     363             : {
     364           2 :   PetscFunctionBegin;
     365           2 :   PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
     366           6 :   PetscValidLogicalCollectiveBool(svd,oneside,2);
     367           2 :   PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
     368           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     369             : }
     370             : 
     371           2 : static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
     372             : {
     373           2 :   SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
     374             : 
     375           2 :   PetscFunctionBegin;
     376           2 :   *oneside = lanczos->oneside;
     377           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     378             : }
     379             : 
     380             : /*@
     381             :    SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
     382             :    to be used is one-sided or two-sided.
     383             : 
     384             :    Not Collective
     385             : 
     386             :    Input Parameters:
     387             : .  svd     - singular value solver
     388             : 
     389             :    Output Parameters:
     390             : .  oneside - boolean flag indicating if the method is one-sided or not
     391             : 
     392             :    Level: advanced
     393             : 
     394             : .seealso: SVDLanczosSetOneSide()
     395             : @*/
     396           2 : PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
     397             : {
     398           2 :   PetscFunctionBegin;
     399           2 :   PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
     400           2 :   PetscAssertPointer(oneside,2);
     401           2 :   PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
     402           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     403             : }
     404             : 
     405          15 : static PetscErrorCode SVDDestroy_Lanczos(SVD svd)
     406             : {
     407          15 :   PetscFunctionBegin;
     408          15 :   PetscCall(PetscFree(svd->data));
     409          15 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL));
     410          15 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL));
     411          15 :   PetscFunctionReturn(PETSC_SUCCESS);
     412             : }
     413             : 
     414           0 : static PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
     415             : {
     416           0 :   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;
     417           0 :   PetscBool      isascii;
     418             : 
     419           0 :   PetscFunctionBegin;
     420           0 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
     421           0 :   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer,"  %s-sided reorthogonalization\n",lanczos->oneside? "one": "two"));
     422           0 :   PetscFunctionReturn(PETSC_SUCCESS);
     423             : }
     424             : 
     425          15 : SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
     426             : {
     427          15 :   SVD_LANCZOS    *ctx;
     428             : 
     429          15 :   PetscFunctionBegin;
     430          15 :   PetscCall(PetscNew(&ctx));
     431          15 :   svd->data = (void*)ctx;
     432             : 
     433          15 :   svd->ops->setup          = SVDSetUp_Lanczos;
     434          15 :   svd->ops->solve          = SVDSolve_Lanczos;
     435          15 :   svd->ops->destroy        = SVDDestroy_Lanczos;
     436          15 :   svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
     437          15 :   svd->ops->view           = SVDView_Lanczos;
     438          15 :   svd->ops->computevectors = SVDComputeVectors_Left;
     439          15 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos));
     440          15 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos));
     441          15 :   PetscFunctionReturn(PETSC_SUCCESS);
     442             : }

Generated by: LCOV version 1.14