LCOV - code coverage report
Current view: top level - svd/impls/cross - cross.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 381 386 98.7 %
Date: 2024-11-21 00:34:55 Functions: 22 22 100.0 %
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: "cross"
      12             : 
      13             :    Method: Uses a Hermitian eigensolver for A^T*A
      14             : */
      15             : 
      16             : #include <slepc/private/svdimpl.h>                /*I "slepcsvd.h" I*/
      17             : 
      18             : typedef struct {
      19             :   PetscBool explicitmatrix;
      20             :   EPS       eps;
      21             :   PetscBool usereps;
      22             :   Mat       C,D;
      23             : } SVD_CROSS;
      24             : 
      25             : typedef struct {
      26             :   Mat       A,AT;
      27             :   Vec       w,diag,omega;
      28             :   PetscBool swapped;
      29             : } SVD_CROSS_SHELL;
      30             : 
      31       29431 : static PetscErrorCode MatMult_Cross(Mat B,Vec x,Vec y)
      32             : {
      33       29431 :   SVD_CROSS_SHELL *ctx;
      34             : 
      35       29431 :   PetscFunctionBegin;
      36       29431 :   PetscCall(MatShellGetContext(B,&ctx));
      37       29431 :   PetscCall(MatMult(ctx->A,x,ctx->w));
      38       29431 :   if (ctx->omega && !ctx->swapped) PetscCall(VecPointwiseMult(ctx->w,ctx->w,ctx->omega));
      39       29431 :   PetscCall(MatMult(ctx->AT,ctx->w,y));
      40       29431 :   PetscFunctionReturn(PETSC_SUCCESS);
      41             : }
      42             : 
      43           9 : static PetscErrorCode MatGetDiagonal_Cross(Mat B,Vec d)
      44             : {
      45           9 :   SVD_CROSS_SHELL   *ctx;
      46           9 :   PetscMPIInt       len;
      47           9 :   PetscInt          N,n,i,j,start,end,ncols;
      48           9 :   PetscScalar       *work1,*work2,*diag;
      49           9 :   const PetscInt    *cols;
      50           9 :   const PetscScalar *vals;
      51             : 
      52           9 :   PetscFunctionBegin;
      53           9 :   PetscCall(MatShellGetContext(B,&ctx));
      54           9 :   if (!ctx->diag) {
      55             :     /* compute diagonal from rows and store in ctx->diag */
      56           9 :     PetscCall(VecDuplicate(d,&ctx->diag));
      57           9 :     PetscCall(MatGetSize(ctx->A,NULL,&N));
      58           9 :     PetscCall(MatGetLocalSize(ctx->A,NULL,&n));
      59           9 :     PetscCall(PetscCalloc2(N,&work1,N,&work2));
      60           9 :     if (ctx->swapped) {
      61           2 :       PetscCall(MatGetOwnershipRange(ctx->AT,&start,&end));
      62          42 :       for (i=start;i<end;i++) {
      63          40 :         PetscCall(MatGetRow(ctx->AT,i,&ncols,NULL,&vals));
      64         120 :         for (j=0;j<ncols;j++) work1[i] += vals[j]*vals[j];
      65          40 :         PetscCall(MatRestoreRow(ctx->AT,i,&ncols,NULL,&vals));
      66             :       }
      67             :     } else {
      68           7 :       PetscCall(MatGetOwnershipRange(ctx->A,&start,&end));
      69         338 :       for (i=start;i<end;i++) {
      70         331 :         PetscCall(MatGetRow(ctx->A,i,&ncols,&cols,&vals));
      71        2214 :         for (j=0;j<ncols;j++) work1[cols[j]] += vals[j]*vals[j];
      72         331 :         PetscCall(MatRestoreRow(ctx->A,i,&ncols,&cols,&vals));
      73             :       }
      74             :     }
      75           9 :     PetscCall(PetscMPIIntCast(N,&len));
      76          27 :     PetscCallMPI(MPIU_Allreduce(work1,work2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)B)));
      77           9 :     PetscCall(VecGetOwnershipRange(ctx->diag,&start,&end));
      78           9 :     PetscCall(VecGetArrayWrite(ctx->diag,&diag));
      79         364 :     for (i=start;i<end;i++) diag[i-start] = work2[i];
      80           9 :     PetscCall(VecRestoreArrayWrite(ctx->diag,&diag));
      81           9 :     PetscCall(PetscFree2(work1,work2));
      82             :   }
      83           9 :   PetscCall(VecCopy(ctx->diag,d));
      84           9 :   PetscFunctionReturn(PETSC_SUCCESS);
      85             : }
      86             : 
      87          45 : static PetscErrorCode MatDestroy_Cross(Mat B)
      88             : {
      89          45 :   SVD_CROSS_SHELL *ctx;
      90             : 
      91          45 :   PetscFunctionBegin;
      92          45 :   PetscCall(MatShellGetContext(B,&ctx));
      93          45 :   PetscCall(VecDestroy(&ctx->w));
      94          45 :   PetscCall(VecDestroy(&ctx->diag));
      95          45 :   PetscCall(PetscFree(ctx));
      96          45 :   PetscFunctionReturn(PETSC_SUCCESS);
      97             : }
      98             : 
      99          72 : static PetscErrorCode SVDCrossGetProductMat(SVD svd,Mat A,Mat AT,Mat *C)
     100             : {
     101          72 :   SVD_CROSS       *cross = (SVD_CROSS*)svd->data;
     102          72 :   SVD_CROSS_SHELL *ctx;
     103          72 :   PetscInt        n;
     104          72 :   VecType         vtype;
     105          72 :   Mat             B;
     106             : 
     107          72 :   PetscFunctionBegin;
     108          72 :   if (cross->explicitmatrix) {
     109          50 :     if (!svd->ishyperbolic || svd->swapped) B = (!svd->expltrans && svd->swapped)? AT: A;
     110             :     else {  /* duplicate A and scale by signature */
     111           2 :       PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B));
     112           2 :       PetscCall(MatDiagonalScale(B,svd->omega,NULL));
     113             :     }
     114          27 :     if (svd->expltrans) {  /* explicit transpose */
     115          23 :       PetscCall(MatProductCreate(AT,B,NULL,C));
     116          23 :       PetscCall(MatProductSetType(*C,MATPRODUCT_AB));
     117             :     } else {  /* implicit transpose */
     118             : #if defined(PETSC_USE_COMPLEX)
     119             :       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Must use explicit transpose with complex scalars");
     120             : #else
     121           4 :       if (!svd->swapped) {
     122           2 :         PetscCall(MatProductCreate(A,B,NULL,C));
     123           2 :         PetscCall(MatProductSetType(*C,MATPRODUCT_AtB));
     124             :       } else {
     125           2 :         PetscCall(MatProductCreate(B,AT,NULL,C));
     126           2 :         PetscCall(MatProductSetType(*C,MATPRODUCT_ABt));
     127             :       }
     128             : #endif
     129             :     }
     130          27 :     PetscCall(MatProductSetFromOptions(*C));
     131          27 :     PetscCall(MatProductSymbolic(*C));
     132          27 :     PetscCall(MatProductNumeric(*C));
     133          27 :     if (svd->ishyperbolic && !svd->swapped) PetscCall(MatDestroy(&B));
     134             :   } else {
     135          45 :     PetscCall(PetscNew(&ctx));
     136          45 :     ctx->A       = A;
     137          45 :     ctx->AT      = AT;
     138          45 :     ctx->omega   = svd->omega;
     139          45 :     ctx->swapped = svd->swapped;
     140          45 :     PetscCall(MatCreateVecs(A,NULL,&ctx->w));
     141          45 :     PetscCall(MatGetLocalSize(A,NULL,&n));
     142          45 :     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)svd),n,n,PETSC_DETERMINE,PETSC_DETERMINE,(void*)ctx,C));
     143          45 :     PetscCall(MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cross));
     144          45 :     if (!svd->ishyperbolic || svd->swapped) PetscCall(MatShellSetOperation(*C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cross));
     145          45 :     PetscCall(MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))MatDestroy_Cross));
     146          45 :     PetscCall(MatGetVecType(A,&vtype));
     147          45 :     PetscCall(MatSetVecType(*C,vtype));
     148             :   }
     149          72 :   PetscFunctionReturn(PETSC_SUCCESS);
     150             : }
     151             : 
     152             : /* Convergence test relative to the norm of R (used in GSVD only) */
     153         445 : static PetscErrorCode EPSConv_Cross(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
     154             : {
     155         445 :   SVD svd = (SVD)ctx;
     156             : 
     157         445 :   PetscFunctionBegin;
     158         445 :   *errest = res/PetscMax(svd->nrma,svd->nrmb);
     159         445 :   PetscFunctionReturn(PETSC_SUCCESS);
     160             : }
     161             : 
     162          58 : static PetscErrorCode SVDSetUp_Cross(SVD svd)
     163             : {
     164          58 :   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
     165          58 :   ST             st;
     166          58 :   PetscBool      trackall,issinv,isks;
     167          58 :   EPSProblemType ptype;
     168          58 :   EPSWhich       which;
     169          58 :   Mat            Omega;
     170          58 :   MatType        Atype;
     171          58 :   PetscInt       n,N;
     172             : 
     173          58 :   PetscFunctionBegin;
     174          58 :   if (!cross->eps) PetscCall(SVDCrossGetEPS(svd,&cross->eps));
     175          58 :   PetscCall(MatDestroy(&cross->C));
     176          58 :   PetscCall(MatDestroy(&cross->D));
     177          58 :   PetscCall(SVDCrossGetProductMat(svd,svd->A,svd->AT,&cross->C));
     178          58 :   if (svd->isgeneralized) {
     179          14 :     PetscCall(SVDCrossGetProductMat(svd,svd->B,svd->BT,&cross->D));
     180          14 :     PetscCall(EPSSetOperators(cross->eps,cross->C,cross->D));
     181          14 :     PetscCall(EPSGetProblemType(cross->eps,&ptype));
     182          14 :     if (!ptype) PetscCall(EPSSetProblemType(cross->eps,EPS_GHEP));
     183          44 :   } else if (svd->ishyperbolic && svd->swapped) {
     184           4 :     PetscCall(MatGetType(svd->OP,&Atype));
     185           4 :     PetscCall(MatGetSize(svd->A,NULL,&N));
     186           4 :     PetscCall(MatGetLocalSize(svd->A,NULL,&n));
     187           4 :     PetscCall(MatCreate(PetscObjectComm((PetscObject)svd),&Omega));
     188           4 :     PetscCall(MatSetSizes(Omega,n,n,N,N));
     189           4 :     PetscCall(MatSetType(Omega,Atype));
     190           4 :     PetscCall(MatDiagonalSet(Omega,svd->omega,INSERT_VALUES));
     191           4 :     PetscCall(EPSSetOperators(cross->eps,cross->C,Omega));
     192           4 :     PetscCall(EPSSetProblemType(cross->eps,EPS_GHIEP));
     193           4 :     PetscCall(MatDestroy(&Omega));
     194             :   } else {
     195          40 :     PetscCall(EPSSetOperators(cross->eps,cross->C,NULL));
     196          40 :     PetscCall(EPSSetProblemType(cross->eps,EPS_HEP));
     197             :   }
     198          58 :   if (!cross->usereps) {
     199          56 :     PetscCall(EPSGetST(cross->eps,&st));
     200          56 :     PetscCall(PetscObjectTypeCompare((PetscObject)st,STSINVERT,&issinv));
     201          56 :     PetscCall(PetscObjectTypeCompare((PetscObject)cross->eps,EPSKRYLOVSCHUR,&isks));
     202          56 :     if (svd->isgeneralized && svd->which==SVD_SMALLEST) {
     203           4 :       if (cross->explicitmatrix && isks && !issinv) {  /* default to shift-and-invert */
     204           2 :         PetscCall(STSetType(st,STSINVERT));
     205           2 :         PetscCall(EPSSetTarget(cross->eps,0.0));
     206             :         which = EPS_TARGET_REAL;
     207           2 :       } else which = issinv?EPS_TARGET_REAL:EPS_SMALLEST_REAL;
     208             :     } else {
     209          52 :       if (issinv) which = EPS_TARGET_MAGNITUDE;
     210          51 :       else if (svd->ishyperbolic) which = svd->which==SVD_LARGEST?EPS_LARGEST_MAGNITUDE:EPS_SMALLEST_MAGNITUDE;
     211          39 :       else which = svd->which==SVD_LARGEST?EPS_LARGEST_REAL:EPS_SMALLEST_REAL;
     212             :     }
     213          56 :     PetscCall(EPSSetWhichEigenpairs(cross->eps,which));
     214          56 :     PetscCall(EPSSetDimensions(cross->eps,svd->nsv,svd->ncv,svd->mpd));
     215          74 :     PetscCall(EPSSetTolerances(cross->eps,svd->tol==(PetscReal)PETSC_DETERMINE?SLEPC_DEFAULT_TOL/10.0:svd->tol,svd->max_it));
     216          56 :     switch (svd->conv) {
     217           3 :     case SVD_CONV_ABS:
     218           3 :       PetscCall(EPSSetConvergenceTest(cross->eps,EPS_CONV_ABS));break;
     219          39 :     case SVD_CONV_REL:
     220          39 :       PetscCall(EPSSetConvergenceTest(cross->eps,EPS_CONV_REL));break;
     221          14 :     case SVD_CONV_NORM:
     222          14 :       if (svd->isgeneralized) {
     223          14 :         if (!svd->nrma) PetscCall(MatNorm(svd->OP,NORM_INFINITY,&svd->nrma));
     224          14 :         if (!svd->nrmb) PetscCall(MatNorm(svd->OPb,NORM_INFINITY,&svd->nrmb));
     225          14 :         PetscCall(EPSSetConvergenceTestFunction(cross->eps,EPSConv_Cross,svd,NULL));
     226             :       } else {
     227           0 :         PetscCall(EPSSetConvergenceTest(cross->eps,EPS_CONV_NORM));break;
     228             :       }
     229             :       break;
     230           0 :     case SVD_CONV_MAXIT:
     231           0 :       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Maxit convergence test not supported in this solver");
     232           0 :     case SVD_CONV_USER:
     233           0 :       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
     234             :     }
     235           2 :   }
     236          58 :   SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
     237             :   /* Transfer the trackall option from svd to eps */
     238          58 :   PetscCall(SVDGetTrackAll(svd,&trackall));
     239          58 :   PetscCall(EPSSetTrackAll(cross->eps,trackall));
     240             :   /* Transfer the initial space from svd to eps */
     241          58 :   if (svd->nini<0) {
     242           5 :     PetscCall(EPSSetInitialSpace(cross->eps,-svd->nini,svd->IS));
     243           5 :     PetscCall(SlepcBasisDestroy_Private(&svd->nini,&svd->IS));
     244             :   }
     245          58 :   PetscCall(EPSSetUp(cross->eps));
     246          58 :   PetscCall(EPSGetDimensions(cross->eps,NULL,&svd->ncv,&svd->mpd));
     247          58 :   PetscCall(EPSGetTolerances(cross->eps,NULL,&svd->max_it));
     248          58 :   if (svd->tol==(PetscReal)PETSC_DETERMINE) svd->tol = SLEPC_DEFAULT_TOL;
     249             : 
     250          58 :   svd->leftbasis = PETSC_FALSE;
     251          58 :   PetscCall(SVDAllocateSolution(svd,0));
     252          58 :   PetscFunctionReturn(PETSC_SUCCESS);
     253             : }
     254             : 
     255          58 : static PetscErrorCode SVDSolve_Cross(SVD svd)
     256             : {
     257          58 :   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
     258          58 :   PetscInt       i;
     259          58 :   PetscScalar    lambda;
     260          58 :   PetscReal      sigma;
     261             : 
     262          58 :   PetscFunctionBegin;
     263          58 :   PetscCall(EPSSolve(cross->eps));
     264          58 :   PetscCall(EPSGetConverged(cross->eps,&svd->nconv));
     265          58 :   PetscCall(EPSGetIterationNumber(cross->eps,&svd->its));
     266          58 :   PetscCall(EPSGetConvergedReason(cross->eps,(EPSConvergedReason*)&svd->reason));
     267         860 :   for (i=0;i<svd->nconv;i++) {
     268         802 :     PetscCall(EPSGetEigenvalue(cross->eps,i,&lambda,NULL));
     269         802 :     sigma = PetscRealPart(lambda);
     270         802 :     if (svd->ishyperbolic) svd->sigma[i] = PetscSqrtReal(PetscAbsReal(sigma));
     271             :     else {
     272         340 :       PetscCheck(sigma>-10*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)svd),PETSC_ERR_FP,"Negative eigenvalue computed by EPS: %g",(double)sigma);
     273         340 :       if (sigma<0.0) {
     274           2 :         PetscCall(PetscInfo(svd,"Negative eigenvalue computed by EPS: %g, resetting to 0\n",(double)sigma));
     275             :         sigma = 0.0;
     276             :       }
     277         340 :       svd->sigma[i] = PetscSqrtReal(sigma);
     278             :     }
     279             :   }
     280          58 :   PetscFunctionReturn(PETSC_SUCCESS);
     281             : }
     282             : 
     283          52 : static PetscErrorCode SVDComputeVectors_Cross(SVD svd)
     284             : {
     285          52 :   SVD_CROSS         *cross = (SVD_CROSS*)svd->data;
     286          52 :   PetscInt          i,mloc,ploc;
     287          52 :   Vec               u,v,x,uv,w,omega2=NULL;
     288          52 :   Mat               Omega;
     289          52 :   PetscScalar       *dst,alpha,lambda,*varray;
     290          52 :   const PetscScalar *src;
     291          52 :   PetscReal         nrm;
     292             : 
     293          52 :   PetscFunctionBegin;
     294          52 :   if (svd->isgeneralized) {
     295          14 :     PetscCall(MatCreateVecs(svd->A,NULL,&u));
     296          14 :     PetscCall(VecGetLocalSize(u,&mloc));
     297          14 :     PetscCall(MatCreateVecs(svd->B,NULL,&v));
     298          14 :     PetscCall(VecGetLocalSize(v,&ploc));
     299         124 :     for (i=0;i<svd->nconv;i++) {
     300         110 :       PetscCall(BVGetColumn(svd->V,i,&x));
     301         110 :       PetscCall(EPSGetEigenpair(cross->eps,i,&lambda,NULL,x,NULL));
     302         110 :       PetscCall(MatMult(svd->A,x,u));     /* u_i*c_i/alpha = A*x_i */
     303         110 :       PetscCall(VecNormalize(u,NULL));
     304         110 :       PetscCall(MatMult(svd->B,x,v));     /* v_i*s_i/alpha = B*x_i */
     305         110 :       PetscCall(VecNormalize(v,&nrm));    /* ||v||_2 = s_i/alpha   */
     306         110 :       alpha = 1.0/(PetscSqrtReal(1.0+PetscRealPart(lambda))*nrm);    /* alpha=s_i/||v||_2 */
     307         110 :       PetscCall(VecScale(x,alpha));
     308         110 :       PetscCall(BVRestoreColumn(svd->V,i,&x));
     309             :       /* copy [u;v] to U[i] */
     310         110 :       PetscCall(BVGetColumn(svd->U,i,&uv));
     311         110 :       PetscCall(VecGetArrayWrite(uv,&dst));
     312         110 :       PetscCall(VecGetArrayRead(u,&src));
     313         110 :       PetscCall(PetscArraycpy(dst,src,mloc));
     314         110 :       PetscCall(VecRestoreArrayRead(u,&src));
     315         110 :       PetscCall(VecGetArrayRead(v,&src));
     316         110 :       PetscCall(PetscArraycpy(dst+mloc,src,ploc));
     317         110 :       PetscCall(VecRestoreArrayRead(v,&src));
     318         110 :       PetscCall(VecRestoreArrayWrite(uv,&dst));
     319         110 :       PetscCall(BVRestoreColumn(svd->U,i,&uv));
     320             :     }
     321          14 :     PetscCall(VecDestroy(&v));
     322          14 :     PetscCall(VecDestroy(&u));
     323          38 :   } else if (svd->ishyperbolic && svd->swapped) {  /* was solved as GHIEP, set u=Omega*u and normalize */
     324           4 :     PetscCall(EPSGetOperators(cross->eps,NULL,&Omega));
     325           4 :     PetscCall(MatCreateVecs(Omega,&w,NULL));
     326           4 :     PetscCall(VecCreateSeq(PETSC_COMM_SELF,svd->ncv,&omega2));
     327           4 :     PetscCall(VecGetArrayWrite(omega2,&varray));
     328         120 :     for (i=0;i<svd->nconv;i++) {
     329         116 :       PetscCall(BVGetColumn(svd->V,i,&v));
     330         116 :       PetscCall(EPSGetEigenvector(cross->eps,i,v,NULL));
     331         116 :       PetscCall(MatMult(Omega,v,w));
     332         116 :       PetscCall(VecDot(v,w,&alpha));
     333         116 :       svd->sign[i] = PetscSign(PetscRealPart(alpha));
     334         116 :       varray[i] = svd->sign[i];
     335         116 :       alpha = 1.0/PetscSqrtScalar(PetscAbsScalar(alpha));
     336         116 :       PetscCall(VecScale(w,alpha));
     337         116 :       PetscCall(VecCopy(w,v));
     338         116 :       PetscCall(BVRestoreColumn(svd->V,i,&v));
     339             :     }
     340           4 :     PetscCall(BVSetSignature(svd->V,omega2));
     341           4 :     PetscCall(VecRestoreArrayWrite(omega2,&varray));
     342           4 :     PetscCall(VecDestroy(&omega2));
     343           4 :     PetscCall(VecDestroy(&w));
     344           4 :     PetscCall(SVDComputeVectors_Left(svd));
     345             :   } else {
     346         595 :     for (i=0;i<svd->nconv;i++) {
     347         561 :       PetscCall(BVGetColumn(svd->V,i,&v));
     348         561 :       PetscCall(EPSGetEigenvector(cross->eps,i,v,NULL));
     349         561 :       PetscCall(BVRestoreColumn(svd->V,i,&v));
     350             :     }
     351          34 :     PetscCall(SVDComputeVectors_Left(svd));
     352             :   }
     353          52 :   PetscFunctionReturn(PETSC_SUCCESS);
     354             : }
     355             : 
     356        1430 : static PetscErrorCode EPSMonitor_Cross(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
     357             : {
     358        1430 :   PetscInt       i;
     359        1430 :   SVD            svd = (SVD)ctx;
     360        1430 :   PetscScalar    er,ei;
     361        1430 :   ST             st;
     362             : 
     363        1430 :   PetscFunctionBegin;
     364        1430 :   PetscCall(EPSGetST(eps,&st));
     365       22899 :   for (i=0;i<PetscMin(nest,svd->ncv);i++) {
     366       21469 :     er = eigr[i]; ei = eigi[i];
     367       21469 :     PetscCall(STBackTransform(st,1,&er,&ei));
     368       21469 :     svd->sigma[i] = PetscSqrtReal(PetscAbsReal(PetscRealPart(er)));
     369       21469 :     svd->errest[i] = errest[i];
     370             :   }
     371        1430 :   PetscCall(SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest));
     372        1430 :   PetscFunctionReturn(PETSC_SUCCESS);
     373             : }
     374             : 
     375          45 : static PetscErrorCode SVDSetFromOptions_Cross(SVD svd,PetscOptionItems *PetscOptionsObject)
     376             : {
     377          45 :   PetscBool      set,val;
     378          45 :   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
     379          45 :   ST             st;
     380             : 
     381          45 :   PetscFunctionBegin;
     382          45 :   PetscOptionsHeadBegin(PetscOptionsObject,"SVD Cross Options");
     383             : 
     384          45 :     PetscCall(PetscOptionsBool("-svd_cross_explicitmatrix","Use cross explicit matrix","SVDCrossSetExplicitMatrix",cross->explicitmatrix,&val,&set));
     385          45 :     if (set) PetscCall(SVDCrossSetExplicitMatrix(svd,val));
     386             : 
     387          45 :   PetscOptionsHeadEnd();
     388             : 
     389          45 :   if (!cross->eps) PetscCall(SVDCrossGetEPS(svd,&cross->eps));
     390          45 :   if (!cross->explicitmatrix && !cross->usereps) {
     391             :     /* use as default an ST with shell matrix and Jacobi */
     392          26 :     PetscCall(EPSGetST(cross->eps,&st));
     393          26 :     PetscCall(STSetMatMode(st,ST_MATMODE_SHELL));
     394             :   }
     395          45 :   PetscCall(EPSSetFromOptions(cross->eps));
     396          45 :   PetscFunctionReturn(PETSC_SUCCESS);
     397             : }
     398             : 
     399          26 : static PetscErrorCode SVDCrossSetExplicitMatrix_Cross(SVD svd,PetscBool explicitmatrix)
     400             : {
     401          26 :   SVD_CROSS *cross = (SVD_CROSS*)svd->data;
     402             : 
     403          26 :   PetscFunctionBegin;
     404          26 :   if (cross->explicitmatrix != explicitmatrix) {
     405          17 :     cross->explicitmatrix = explicitmatrix;
     406          17 :     svd->state = SVD_STATE_INITIAL;
     407             :   }
     408          26 :   PetscFunctionReturn(PETSC_SUCCESS);
     409             : }
     410             : 
     411             : /*@
     412             :    SVDCrossSetExplicitMatrix - Indicate if the eigensolver operator A^T*A must
     413             :    be computed explicitly.
     414             : 
     415             :    Logically Collective
     416             : 
     417             :    Input Parameters:
     418             : +  svd         - singular value solver
     419             : -  explicitmat - boolean flag indicating if A^T*A is built explicitly
     420             : 
     421             :    Options Database Key:
     422             : .  -svd_cross_explicitmatrix <boolean> - Indicates the boolean flag
     423             : 
     424             :    Level: advanced
     425             : 
     426             : .seealso: SVDCrossGetExplicitMatrix()
     427             : @*/
     428          26 : PetscErrorCode SVDCrossSetExplicitMatrix(SVD svd,PetscBool explicitmat)
     429             : {
     430          26 :   PetscFunctionBegin;
     431          26 :   PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
     432          78 :   PetscValidLogicalCollectiveBool(svd,explicitmat,2);
     433          26 :   PetscTryMethod(svd,"SVDCrossSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmat));
     434          26 :   PetscFunctionReturn(PETSC_SUCCESS);
     435             : }
     436             : 
     437           2 : static PetscErrorCode SVDCrossGetExplicitMatrix_Cross(SVD svd,PetscBool *explicitmat)
     438             : {
     439           2 :   SVD_CROSS *cross = (SVD_CROSS*)svd->data;
     440             : 
     441           2 :   PetscFunctionBegin;
     442           2 :   *explicitmat = cross->explicitmatrix;
     443           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     444             : }
     445             : 
     446             : /*@
     447             :    SVDCrossGetExplicitMatrix - Returns the flag indicating if A^T*A is built explicitly.
     448             : 
     449             :    Not Collective
     450             : 
     451             :    Input Parameter:
     452             : .  svd  - singular value solver
     453             : 
     454             :    Output Parameter:
     455             : .  explicitmat - the mode flag
     456             : 
     457             :    Level: advanced
     458             : 
     459             : .seealso: SVDCrossSetExplicitMatrix()
     460             : @*/
     461           2 : PetscErrorCode SVDCrossGetExplicitMatrix(SVD svd,PetscBool *explicitmat)
     462             : {
     463           2 :   PetscFunctionBegin;
     464           2 :   PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
     465           2 :   PetscAssertPointer(explicitmat,2);
     466           2 :   PetscUseMethod(svd,"SVDCrossGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmat));
     467           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     468             : }
     469             : 
     470           2 : static PetscErrorCode SVDCrossSetEPS_Cross(SVD svd,EPS eps)
     471             : {
     472           2 :   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
     473             : 
     474           2 :   PetscFunctionBegin;
     475           2 :   PetscCall(PetscObjectReference((PetscObject)eps));
     476           2 :   PetscCall(EPSDestroy(&cross->eps));
     477           2 :   cross->eps     = eps;
     478           2 :   cross->usereps = PETSC_TRUE;
     479           2 :   svd->state     = SVD_STATE_INITIAL;
     480           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     481             : }
     482             : 
     483             : /*@
     484             :    SVDCrossSetEPS - Associate an eigensolver object (EPS) to the
     485             :    singular value solver.
     486             : 
     487             :    Collective
     488             : 
     489             :    Input Parameters:
     490             : +  svd - singular value solver
     491             : -  eps - the eigensolver object
     492             : 
     493             :    Level: advanced
     494             : 
     495             : .seealso: SVDCrossGetEPS()
     496             : @*/
     497           2 : PetscErrorCode SVDCrossSetEPS(SVD svd,EPS eps)
     498             : {
     499           2 :   PetscFunctionBegin;
     500           2 :   PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
     501           2 :   PetscValidHeaderSpecific(eps,EPS_CLASSID,2);
     502           2 :   PetscCheckSameComm(svd,1,eps,2);
     503           2 :   PetscTryMethod(svd,"SVDCrossSetEPS_C",(SVD,EPS),(svd,eps));
     504           2 :   PetscFunctionReturn(PETSC_SUCCESS);
     505             : }
     506             : 
     507          45 : static PetscErrorCode SVDCrossGetEPS_Cross(SVD svd,EPS *eps)
     508             : {
     509          45 :   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
     510             : 
     511          45 :   PetscFunctionBegin;
     512          45 :   if (!cross->eps) {
     513          45 :     PetscCall(EPSCreate(PetscObjectComm((PetscObject)svd),&cross->eps));
     514          45 :     PetscCall(PetscObjectIncrementTabLevel((PetscObject)cross->eps,(PetscObject)svd,1));
     515          45 :     PetscCall(EPSSetOptionsPrefix(cross->eps,((PetscObject)svd)->prefix));
     516          45 :     PetscCall(EPSAppendOptionsPrefix(cross->eps,"svd_cross_"));
     517          45 :     PetscCall(PetscObjectSetOptions((PetscObject)cross->eps,((PetscObject)svd)->options));
     518          45 :     PetscCall(EPSSetWhichEigenpairs(cross->eps,EPS_LARGEST_REAL));
     519          45 :     PetscCall(EPSMonitorSet(cross->eps,EPSMonitor_Cross,svd,NULL));
     520             :   }
     521          45 :   *eps = cross->eps;
     522          45 :   PetscFunctionReturn(PETSC_SUCCESS);
     523             : }
     524             : 
     525             : /*@
     526             :    SVDCrossGetEPS - Retrieve the eigensolver object (EPS) associated
     527             :    to the singular value solver.
     528             : 
     529             :    Collective
     530             : 
     531             :    Input Parameter:
     532             : .  svd - singular value solver
     533             : 
     534             :    Output Parameter:
     535             : .  eps - the eigensolver object
     536             : 
     537             :    Level: advanced
     538             : 
     539             : .seealso: SVDCrossSetEPS()
     540             : @*/
     541          45 : PetscErrorCode SVDCrossGetEPS(SVD svd,EPS *eps)
     542             : {
     543          45 :   PetscFunctionBegin;
     544          45 :   PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
     545          45 :   PetscAssertPointer(eps,2);
     546          45 :   PetscUseMethod(svd,"SVDCrossGetEPS_C",(SVD,EPS*),(svd,eps));
     547          45 :   PetscFunctionReturn(PETSC_SUCCESS);
     548             : }
     549             : 
     550           1 : static PetscErrorCode SVDView_Cross(SVD svd,PetscViewer viewer)
     551             : {
     552           1 :   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
     553           1 :   PetscBool      isascii;
     554             : 
     555           1 :   PetscFunctionBegin;
     556           1 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
     557           1 :   if (isascii) {
     558           1 :     if (!cross->eps) PetscCall(SVDCrossGetEPS(svd,&cross->eps));
     559           2 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  %s matrix\n",cross->explicitmatrix?"explicit":"implicit"));
     560           1 :     PetscCall(PetscViewerASCIIPushTab(viewer));
     561           1 :     PetscCall(EPSView(cross->eps,viewer));
     562           1 :     PetscCall(PetscViewerASCIIPopTab(viewer));
     563             :   }
     564           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     565             : }
     566             : 
     567          48 : static PetscErrorCode SVDReset_Cross(SVD svd)
     568             : {
     569          48 :   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
     570             : 
     571          48 :   PetscFunctionBegin;
     572          48 :   PetscCall(EPSReset(cross->eps));
     573          48 :   PetscCall(MatDestroy(&cross->C));
     574          48 :   PetscCall(MatDestroy(&cross->D));
     575          48 :   PetscFunctionReturn(PETSC_SUCCESS);
     576             : }
     577             : 
     578          47 : static PetscErrorCode SVDDestroy_Cross(SVD svd)
     579             : {
     580          47 :   SVD_CROSS      *cross = (SVD_CROSS*)svd->data;
     581             : 
     582          47 :   PetscFunctionBegin;
     583          47 :   PetscCall(EPSDestroy(&cross->eps));
     584          47 :   PetscCall(PetscFree(svd->data));
     585          47 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",NULL));
     586          47 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",NULL));
     587          47 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",NULL));
     588          47 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",NULL));
     589          47 :   PetscFunctionReturn(PETSC_SUCCESS);
     590             : }
     591             : 
     592          47 : SLEPC_EXTERN PetscErrorCode SVDCreate_Cross(SVD svd)
     593             : {
     594          47 :   SVD_CROSS      *cross;
     595             : 
     596          47 :   PetscFunctionBegin;
     597          47 :   PetscCall(PetscNew(&cross));
     598          47 :   svd->data = (void*)cross;
     599             : 
     600          47 :   svd->ops->solve          = SVDSolve_Cross;
     601          47 :   svd->ops->solveg         = SVDSolve_Cross;
     602          47 :   svd->ops->solveh         = SVDSolve_Cross;
     603          47 :   svd->ops->setup          = SVDSetUp_Cross;
     604          47 :   svd->ops->setfromoptions = SVDSetFromOptions_Cross;
     605          47 :   svd->ops->destroy        = SVDDestroy_Cross;
     606          47 :   svd->ops->reset          = SVDReset_Cross;
     607          47 :   svd->ops->view           = SVDView_Cross;
     608          47 :   svd->ops->computevectors = SVDComputeVectors_Cross;
     609          47 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetEPS_C",SVDCrossSetEPS_Cross));
     610          47 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetEPS_C",SVDCrossGetEPS_Cross));
     611          47 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossSetExplicitMatrix_C",SVDCrossSetExplicitMatrix_Cross));
     612          47 :   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCrossGetExplicitMatrix_C",SVDCrossGetExplicitMatrix_Cross));
     613          47 :   PetscFunctionReturn(PETSC_SUCCESS);
     614             : }

Generated by: LCOV version 1.14