LCOV - code coverage report
Current view: top level - pep/impls/ciss - pciss.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 580 640 90.6 %
Date: 2024-12-18 00:51:33 Functions: 28 29 96.6 %
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 eigensolver: "ciss"
      12             : 
      13             :    Method: Contour Integral Spectral Slicing
      14             : 
      15             :    Algorithm:
      16             : 
      17             :        Contour integral based on Sakurai-Sugiura method to construct a
      18             :        subspace, with various eigenpair extractions (Rayleigh-Ritz,
      19             :        explicit moment).
      20             : 
      21             :    Based on code contributed by Y. Maeda, T. Sakurai.
      22             : 
      23             :    References:
      24             : 
      25             :        [1] J. Asakura, T. Sakurai, H. Tadano, T. Ikegami, K. Kimura, "A
      26             :            numerical method for polynomial eigenvalue problems using contour
      27             :            integral", Japan J. Indust. Appl. Math. 27:73-90, 2010.
      28             : */
      29             : 
      30             : #include <slepc/private/pepimpl.h>         /*I "slepcpep.h" I*/
      31             : #include <slepc/private/slepccontour.h>
      32             : 
      33             : typedef struct {
      34             :   /* parameters */
      35             :   PetscInt          N;             /* number of integration points (32) */
      36             :   PetscInt          L;             /* block size (16) */
      37             :   PetscInt          M;             /* moment degree (N/4 = 4) */
      38             :   PetscReal         delta;         /* threshold of singular value (1e-12) */
      39             :   PetscInt          L_max;         /* maximum number of columns of the source matrix V */
      40             :   PetscReal         spurious_threshold; /* discard spurious eigenpairs */
      41             :   PetscBool         isreal;        /* T(z) is real for real z */
      42             :   PetscInt          npart;         /* number of partitions */
      43             :   PetscInt          refine_inner;
      44             :   PetscInt          refine_blocksize;
      45             :   PEPCISSExtraction extraction;
      46             :   /* private data */
      47             :   SlepcContourData  contour;
      48             :   PetscReal         *sigma;        /* threshold for numerical rank */
      49             :   PetscScalar       *weight;
      50             :   PetscScalar       *omega;
      51             :   PetscScalar       *pp;
      52             :   BV                V;
      53             :   BV                S;
      54             :   BV                Y;
      55             :   PetscBool         useconj;
      56             :   Mat               J,*Psplit;     /* auxiliary matrices */
      57             :   BV                pV;
      58             :   PetscObjectId     rgid;
      59             :   PetscObjectState  rgstate;
      60             : } PEP_CISS;
      61             : 
      62         872 : static PetscErrorCode PEPComputeFunction(PEP pep,PetscScalar lambda,Mat T,Mat P,PetscBool deriv)
      63             : {
      64         872 :   PetscInt         i;
      65         872 :   PetscScalar      *coeff;
      66         872 :   Mat              *A,*K;
      67         872 :   MatStructure     str,strp;
      68         872 :   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
      69         872 :   SlepcContourData contour = ctx->contour;
      70             : 
      71         872 :   PetscFunctionBegin;
      72         872 :   A = (contour->pA)?contour->pA:pep->A;
      73         872 :   K = (contour->pP)?contour->pP:ctx->Psplit;
      74         872 :   PetscCall(PetscMalloc1(pep->nmat,&coeff));
      75         872 :   if (deriv) PetscCall(PEPEvaluateBasisDerivative(pep,lambda,0,coeff,NULL));
      76         372 :   else PetscCall(PEPEvaluateBasis(pep,lambda,0,coeff,NULL));
      77         872 :   PetscCall(STGetMatStructure(pep->st,&str));
      78         872 :   PetscCall(MatZeroEntries(T));
      79         872 :   if (!deriv && T != P) {
      80          64 :     PetscCall(STGetSplitPreconditionerInfo(pep->st,NULL,&strp));
      81          64 :     PetscCall(MatZeroEntries(P));
      82             :   }
      83         872 :   i = deriv?1:0;
      84        4252 :   for (;i<pep->nmat;i++) {
      85        3380 :     PetscCall(MatAXPY(T,coeff[i],A[i],str));
      86        3380 :     if (!deriv && T != P) PetscCall(MatAXPY(P,coeff[i],K[i],strp));
      87             :   }
      88         872 :   PetscCall(PetscFree(coeff));
      89         872 :   PetscFunctionReturn(PETSC_SUCCESS);
      90             : }
      91             : 
      92             : /*
      93             :   Set up KSP solvers for every integration point
      94             : */
      95          13 : static PetscErrorCode PEPCISSSetUp(PEP pep,Mat T,Mat P)
      96             : {
      97          13 :   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
      98          13 :   SlepcContourData contour;
      99          13 :   PetscInt         i,p_id;
     100          13 :   Mat              Amat,Pmat;
     101             : 
     102          13 :   PetscFunctionBegin;
     103          13 :   if (!ctx->contour || !ctx->contour->ksp) PetscCall(PEPCISSGetKSPs(pep,NULL,NULL));
     104          13 :   contour = ctx->contour;
     105          13 :   PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Something went wrong with PEPCISSGetKSPs()");
     106         385 :   for (i=0;i<contour->npoints;i++) {
     107         372 :     p_id = i*contour->subcomm->n + contour->subcomm->color;
     108         372 :     PetscCall(MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&Amat));
     109         372 :     if (T != P) PetscCall(MatDuplicate(P,MAT_DO_NOT_COPY_VALUES,&Pmat)); else Pmat = Amat;
     110         372 :     PetscCall(PEPComputeFunction(pep,ctx->omega[p_id],Amat,Pmat,PETSC_FALSE));
     111         372 :     PetscCall(PEP_KSPSetOperators(contour->ksp[i],Amat,Pmat));
     112         372 :     PetscCall(MatDestroy(&Amat));
     113         372 :     if (T != P) PetscCall(MatDestroy(&Pmat));
     114             :   }
     115          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     116             : }
     117             : 
     118             : /*
     119             :   Y_i = F(z_i)^{-1}Fp(z_i)V for every integration point, Y=[Y_i] is in the context
     120             : */
     121          17 : static PetscErrorCode PEPCISSSolve(PEP pep,Mat dT,BV V,PetscInt L_start,PetscInt L_end)
     122             : {
     123          17 :   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
     124          17 :   SlepcContourData contour;
     125          17 :   PetscInt         i,p_id;
     126          17 :   Mat              MV,BMV=NULL,MC;
     127             : 
     128          17 :   PetscFunctionBegin;
     129          17 :   contour = ctx->contour;
     130          17 :   PetscCall(BVSetActiveColumns(V,L_start,L_end));
     131          17 :   PetscCall(BVGetMat(V,&MV));
     132         517 :   for (i=0;i<contour->npoints;i++) {
     133         500 :     p_id = i*contour->subcomm->n + contour->subcomm->color;
     134         500 :     PetscCall(PEPComputeFunction(pep,ctx->omega[p_id],dT,NULL,PETSC_TRUE));
     135         500 :     PetscCall(BVSetActiveColumns(ctx->Y,i*ctx->L+L_start,i*ctx->L+L_end));
     136         500 :     PetscCall(BVGetMat(ctx->Y,&MC));
     137         500 :     if (!i) {
     138          17 :       PetscCall(MatProductCreate(dT,MV,NULL,&BMV));
     139          17 :       PetscCall(MatProductSetType(BMV,MATPRODUCT_AB));
     140          17 :       PetscCall(MatProductSetFromOptions(BMV));
     141          17 :       PetscCall(MatProductSymbolic(BMV));
     142             :     }
     143         500 :     PetscCall(MatProductNumeric(BMV));
     144         500 :     PetscCall(KSPMatSolve(contour->ksp[i],BMV,MC));
     145         500 :     PetscCall(BVRestoreMat(ctx->Y,&MC));
     146             :   }
     147          17 :   PetscCall(MatDestroy(&BMV));
     148          17 :   PetscCall(BVRestoreMat(V,&MV));
     149          17 :   PetscFunctionReturn(PETSC_SUCCESS);
     150             : }
     151             : 
     152          13 : static PetscErrorCode PEPSetUp_CISS(PEP pep)
     153             : {
     154          13 :   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
     155          13 :   SlepcContourData contour;
     156          13 :   PetscInt         i,nwork,nsplit;
     157          13 :   PetscBool        istrivial,isellipse,flg;
     158          13 :   PetscObjectId    id;
     159          13 :   PetscObjectState state;
     160          13 :   Vec              v0;
     161             : 
     162          13 :   PetscFunctionBegin;
     163          13 :   if (pep->ncv==PETSC_DETERMINE) pep->ncv = ctx->L_max*ctx->M;
     164             :   else {
     165           0 :     ctx->L_max = pep->ncv/ctx->M;
     166           0 :     if (!ctx->L_max) {
     167           0 :       ctx->L_max = 1;
     168           0 :       pep->ncv = ctx->L_max*ctx->M;
     169             :     }
     170             :   }
     171          13 :   ctx->L = PetscMin(ctx->L,ctx->L_max);
     172          13 :   if (pep->max_it==PETSC_DETERMINE) pep->max_it = 5;
     173          13 :   if (pep->mpd==PETSC_DETERMINE) pep->mpd = pep->ncv;
     174          13 :   if (!pep->which) pep->which = PEP_ALL;
     175          13 :   PetscCheck(pep->which==PEP_ALL,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"This solver supports only computing all eigenvalues");
     176          13 :   PEPCheckUnsupported(pep,PEP_FEATURE_STOPPING);
     177          13 :   PEPCheckIgnored(pep,PEP_FEATURE_SCALE);
     178             : 
     179             :   /* check region */
     180          13 :   PetscCall(RGIsTrivial(pep->rg,&istrivial));
     181          13 :   PetscCheck(!istrivial,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"CISS requires a nontrivial region, e.g. -rg_type ellipse ...");
     182          13 :   PetscCall(RGGetComplement(pep->rg,&flg));
     183          13 :   PetscCheck(!flg,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"A region with complement flag set is not allowed");
     184          13 :   PetscCall(PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse));
     185          13 :   PetscCheck(isellipse,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Currently only implemented for elliptic regions");
     186             : 
     187             :   /* if the region has changed, then reset contour data */
     188          13 :   PetscCall(PetscObjectGetId((PetscObject)pep->rg,&id));
     189          13 :   PetscCall(PetscObjectStateGet((PetscObject)pep->rg,&state));
     190          13 :   if (ctx->rgid && (id != ctx->rgid || state != ctx->rgstate)) {
     191           0 :     PetscCall(SlepcContourDataDestroy(&ctx->contour));
     192           0 :     PetscCall(PetscInfo(pep,"Resetting the contour data structure due to a change of region\n"));
     193           0 :     ctx->rgid = id; ctx->rgstate = state;
     194             :   }
     195             : 
     196             :   /* create contour data structure */
     197          13 :   if (!ctx->contour) {
     198           0 :     PetscCall(RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj));
     199           0 :     PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour));
     200             :   }
     201             : 
     202          13 :   PetscCall(PEPAllocateSolution(pep,0));
     203          13 :   if (ctx->weight) PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
     204          13 :   PetscCall(PetscMalloc4(ctx->N,&ctx->weight,ctx->N,&ctx->omega,ctx->N,&ctx->pp,ctx->L_max*ctx->M,&ctx->sigma));
     205             : 
     206             :   /* allocate basis vectors */
     207          13 :   PetscCall(BVDestroy(&ctx->S));
     208          13 :   PetscCall(BVDuplicateResize(pep->V,ctx->L*ctx->M,&ctx->S));
     209          13 :   PetscCall(BVDestroy(&ctx->V));
     210          13 :   PetscCall(BVDuplicateResize(pep->V,ctx->L,&ctx->V));
     211             : 
     212             :   /* check if a user-defined split preconditioner has been set */
     213          13 :   PetscCall(STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL));
     214          13 :   if (nsplit) {
     215           3 :     PetscCall(PetscFree(ctx->Psplit));
     216           3 :     PetscCall(PetscMalloc1(nsplit,&ctx->Psplit));
     217          12 :     for (i=0;i<nsplit;i++) PetscCall(STGetSplitPreconditionerTerm(pep->st,i,&ctx->Psplit[i]));
     218             :   }
     219             : 
     220          13 :   contour = ctx->contour;
     221          13 :   PetscCall(SlepcContourRedundantMat(contour,pep->nmat,pep->A,ctx->Psplit));
     222          13 :   if (!ctx->J) PetscCall(MatDuplicate(contour->pA?contour->pA[0]:pep->A[0],MAT_DO_NOT_COPY_VALUES,&ctx->J));
     223          13 :   if (contour->pA) {
     224           4 :     PetscCall(BVGetColumn(ctx->V,0,&v0));
     225           4 :     PetscCall(SlepcContourScatterCreate(contour,v0));
     226           4 :     PetscCall(BVRestoreColumn(ctx->V,0,&v0));
     227           4 :     PetscCall(BVDestroy(&ctx->pV));
     228           4 :     PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->pV));
     229           4 :     PetscCall(BVSetSizesFromVec(ctx->pV,contour->xsub,pep->n));
     230           4 :     PetscCall(BVSetFromOptions(ctx->pV));
     231           4 :     PetscCall(BVResize(ctx->pV,ctx->L,PETSC_FALSE));
     232             :   }
     233             : 
     234          13 :   PetscCall(BVDestroy(&ctx->Y));
     235          13 :   if (contour->pA) {
     236           4 :     PetscCall(BVCreate(PetscObjectComm((PetscObject)contour->xsub),&ctx->Y));
     237           4 :     PetscCall(BVSetSizesFromVec(ctx->Y,contour->xsub,pep->n));
     238           4 :     PetscCall(BVSetFromOptions(ctx->Y));
     239           4 :     PetscCall(BVResize(ctx->Y,contour->npoints*ctx->L,PETSC_FALSE));
     240           9 :   } else PetscCall(BVDuplicateResize(pep->V,contour->npoints*ctx->L,&ctx->Y));
     241             : 
     242          13 :   if (ctx->extraction == PEP_CISS_EXTRACTION_RITZ) {
     243          10 :     PetscCall(DSPEPSetDegree(pep->ds,pep->nmat-1));
     244          10 :     PetscCall(DSPEPSetCoefficients(pep->ds,pep->pbc));
     245             :   }
     246          13 :   PetscCall(DSAllocate(pep->ds,pep->ncv));
     247          13 :   nwork = 2;
     248          13 :   PetscCall(PEPSetWorkVecs(pep,nwork));
     249          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     250             : }
     251             : 
     252          13 : static PetscErrorCode PEPSolve_CISS(PEP pep)
     253             : {
     254          13 :   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
     255          13 :   SlepcContourData contour = ctx->contour;
     256          13 :   Mat              X,M,E,T,P;
     257          13 :   PetscInt         i,j,ld,L_add=0,nv=0,L_base=ctx->L,inner,*inside,nsplit;
     258          13 :   PetscScalar      *Mu,*H0,*H1,*rr,*temp,center;
     259          13 :   PetscReal        error,max_error,radius,rgscale,est_eig,eta;
     260          13 :   PetscBool        isellipse,*fl1;
     261          13 :   Vec              si;
     262          13 :   SlepcSC          sc;
     263          13 :   PetscRandom      rand;
     264             : 
     265          13 :   PetscFunctionBegin;
     266          13 :   PetscCall(DSGetSlepcSC(pep->ds,&sc));
     267          13 :   sc->comparison    = SlepcCompareLargestMagnitude;
     268          13 :   sc->comparisonctx = NULL;
     269          13 :   sc->map           = NULL;
     270          13 :   sc->mapobj        = NULL;
     271          13 :   PetscCall(DSGetLeadingDimension(pep->ds,&ld));
     272          13 :   PetscCall(RGComputeQuadrature(pep->rg,RG_QUADRULE_TRAPEZOIDAL,ctx->N,ctx->omega,ctx->pp,ctx->weight));
     273          13 :   PetscCall(STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL));
     274          13 :   if (contour->pA) {
     275           4 :     T = contour->pA[0];
     276           4 :     P = nsplit? contour->pP[0]: T;
     277             :   } else {
     278           9 :     T = pep->A[0];
     279           9 :     P = nsplit? ctx->Psplit[0]: T;
     280             :   }
     281          13 :   PetscCall(PEPCISSSetUp(pep,T,P));
     282          13 :   PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
     283          13 :   PetscCall(BVSetRandomSign(ctx->V));
     284          13 :   PetscCall(BVGetRandomContext(ctx->V,&rand));
     285          13 :   if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
     286          13 :   PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L));
     287          13 :   PetscCall(PetscObjectTypeCompare((PetscObject)pep->rg,RGELLIPSE,&isellipse));
     288          13 :   if (isellipse) {
     289          13 :     PetscCall(BVTraceQuadrature(ctx->Y,ctx->V,ctx->L,ctx->L,ctx->weight,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj,&est_eig));
     290          13 :     PetscCall(PetscInfo(pep,"Estimated eigenvalue count: %f\n",(double)est_eig));
     291          13 :     eta = PetscPowReal(10.0,-PetscLog10Real(pep->tol)/ctx->N);
     292          13 :     L_add = PetscMax(0,(PetscInt)PetscCeilReal((est_eig*eta)/ctx->M)-ctx->L);
     293          13 :     if (L_add>ctx->L_max-ctx->L) {
     294           0 :       PetscCall(PetscInfo(pep,"Number of eigenvalues inside the contour path may be too large\n"));
     295           0 :       L_add = ctx->L_max-ctx->L;
     296             :     }
     297             :   }
     298             :   /* Updates L after estimate the number of eigenvalue */
     299          13 :   if (L_add>0) {
     300           2 :     PetscCall(PetscInfo(pep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by Estimate #Eig\n",ctx->L,ctx->L+L_add));
     301           2 :     PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
     302           2 :     PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
     303           2 :     PetscCall(BVSetRandomSign(ctx->V));
     304           2 :     if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
     305           2 :     ctx->L += L_add;
     306           2 :     PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L-L_add,ctx->L));
     307             :   }
     308             : 
     309          13 :   PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
     310          15 :   for (i=0;i<ctx->refine_blocksize;i++) {
     311           5 :     PetscCall(BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj));
     312           5 :     PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
     313           5 :     PetscCall(PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0));
     314           5 :     PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
     315           5 :     PetscCall(PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0));
     316           5 :     if (ctx->sigma[0]<=ctx->delta || nv < ctx->L*ctx->M || ctx->L == ctx->L_max) break;
     317           2 :     L_add = L_base;
     318           2 :     if (ctx->L+L_add>ctx->L_max) L_add = ctx->L_max-ctx->L;
     319           2 :     PetscCall(PetscInfo(pep,"Changing L %" PetscInt_FMT " -> %" PetscInt_FMT " by SVD(H0)\n",ctx->L,ctx->L+L_add));
     320           2 :     PetscCall(BVCISSResizeBases(ctx->S,contour->pA?ctx->pV:ctx->V,ctx->Y,ctx->L,ctx->L+L_add,ctx->M,contour->npoints));
     321           2 :     PetscCall(BVSetActiveColumns(ctx->V,ctx->L,ctx->L+L_add));
     322           2 :     PetscCall(BVSetRandomSign(ctx->V));
     323           2 :     if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
     324           2 :     ctx->L += L_add;
     325           2 :     PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,ctx->L-L_add,ctx->L));
     326           2 :     if (L_add) {
     327           2 :       PetscCall(PetscFree2(Mu,H0));
     328           2 :       PetscCall(PetscMalloc2(ctx->L*ctx->L*ctx->M*2,&Mu,ctx->L*ctx->M*ctx->L*ctx->M,&H0));
     329             :     }
     330             :   }
     331             : 
     332          13 :   PetscCall(RGGetScale(pep->rg,&rgscale));
     333          13 :   PetscCall(RGEllipseGetParameters(pep->rg,&center,&radius,NULL));
     334             : 
     335          13 :   if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) PetscCall(PetscMalloc1(ctx->L*ctx->M*ctx->L*ctx->M,&H1));
     336             : 
     337          26 :   while (pep->reason == PEP_CONVERGED_ITERATING) {
     338          13 :     pep->its++;
     339          13 :     for (inner=0;inner<=ctx->refine_inner;inner++) {
     340          13 :       if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
     341           2 :         PetscCall(BVDotQuadrature(ctx->Y,(contour->pA)?ctx->pV:ctx->V,Mu,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->subcomm,contour->npoints,ctx->useconj));
     342           2 :         PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
     343           2 :         PetscCall(PetscLogEventBegin(PEP_CISS_SVD,pep,0,0,0));
     344           2 :         PetscCall(SlepcCISS_BH_SVD(H0,ctx->L*ctx->M,ctx->delta,ctx->sigma,&nv));
     345           2 :         PetscCall(PetscLogEventEnd(PEP_CISS_SVD,pep,0,0,0));
     346             :       } else {
     347          11 :         PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
     348             :         /* compute SVD of S */
     349          21 :         PetscCall(BVSVDAndRank(ctx->S,ctx->M,ctx->L,ctx->delta,(ctx->extraction==PEP_CISS_EXTRACTION_CAA)?BV_SVD_METHOD_QR_CAA:BV_SVD_METHOD_QR,H0,ctx->sigma,&nv));
     350             :       }
     351          13 :       PetscCall(PetscInfo(pep,"Estimated rank: nv = %" PetscInt_FMT "\n",nv));
     352          13 :       if (ctx->sigma[0]>ctx->delta && nv==ctx->L*ctx->M && inner!=ctx->refine_inner) {
     353           0 :         PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
     354           0 :         PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
     355           0 :         PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
     356           0 :         PetscCall(BVCopy(ctx->S,ctx->V));
     357           0 :         if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
     358           0 :         PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L));
     359             :       } else break;
     360             :     }
     361          13 :     pep->nconv = 0;
     362          13 :     if (nv == 0) { pep->reason = PEP_CONVERGED_TOL; break; }
     363             :     else {
     364             :       /* Extracting eigenpairs */
     365          13 :       PetscCall(DSSetDimensions(pep->ds,nv,0,0));
     366          13 :       PetscCall(DSSetState(pep->ds,DS_STATE_RAW));
     367          13 :       if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
     368           2 :         PetscCall(CISS_BlockHankel(Mu,0,ctx->L,ctx->M,H0));
     369           2 :         PetscCall(CISS_BlockHankel(Mu,1,ctx->L,ctx->M,H1));
     370           2 :         PetscCall(DSGetArray(pep->ds,DS_MAT_A,&temp));
     371          75 :         for (j=0;j<nv;j++)
     372        3750 :           for (i=0;i<nv;i++)
     373        3677 :             temp[i+j*ld] = H1[i+j*ctx->L*ctx->M];
     374           2 :         PetscCall(DSRestoreArray(pep->ds,DS_MAT_A,&temp));
     375           2 :         PetscCall(DSGetArray(pep->ds,DS_MAT_B,&temp));
     376          75 :         for (j=0;j<nv;j++)
     377        3750 :           for (i=0;i<nv;i++)
     378        3677 :             temp[i+j*ld] = H0[i+j*ctx->L*ctx->M];
     379           2 :         PetscCall(DSRestoreArray(pep->ds,DS_MAT_B,&temp));
     380          11 :       } else if (ctx->extraction == PEP_CISS_EXTRACTION_CAA) {
     381           1 :         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
     382           1 :         PetscCall(DSGetArray(pep->ds,DS_MAT_A,&temp));
     383          23 :         for (i=0;i<nv;i++) PetscCall(PetscArraycpy(temp+i*ld,H0+i*nv,nv));
     384           1 :         PetscCall(DSRestoreArray(pep->ds,DS_MAT_A,&temp));
     385             :       } else {
     386          10 :         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
     387          50 :         for (i=0;i<pep->nmat;i++) {
     388          40 :           PetscCall(DSGetMat(pep->ds,DSMatExtra[i],&E));
     389          40 :           PetscCall(BVMatProject(ctx->S,pep->A[i],ctx->S,E));
     390          40 :           PetscCall(DSRestoreMat(pep->ds,DSMatExtra[i],&E));
     391             :         }
     392          10 :         nv = (pep->nmat-1)*nv;
     393             :       }
     394          13 :       PetscCall(DSSolve(pep->ds,pep->eigr,pep->eigi));
     395          13 :       PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi));
     396          13 :       if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
     397          98 :         for (i=0;i<nv;i++) {
     398          95 :           pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
     399             :         }
     400             :       }
     401          13 :       PetscCall(PetscMalloc3(nv,&fl1,nv,&inside,nv,&rr));
     402          13 :       PetscCall(DSVectors(pep->ds,DS_MAT_X,NULL,NULL));
     403          13 :       PetscCall(DSGetMat(pep->ds,DS_MAT_X,&X));
     404          13 :       PetscCall(SlepcCISS_isGhost(X,nv,ctx->sigma,ctx->spurious_threshold,fl1));
     405          13 :       PetscCall(DSRestoreMat(pep->ds,DS_MAT_X,&X));
     406          13 :       PetscCall(RGCheckInside(pep->rg,nv,pep->eigr,pep->eigi,inside));
     407        1406 :       for (i=0;i<nv;i++) {
     408        1393 :         if (fl1[i] && inside[i]>=0) {
     409         101 :           rr[i] = 1.0;
     410         101 :           pep->nconv++;
     411        1292 :         } else rr[i] = 0.0;
     412             :       }
     413          13 :       PetscCall(DSSort(pep->ds,pep->eigr,pep->eigi,rr,NULL,&pep->nconv));
     414          13 :       PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi));
     415          13 :       if (ctx->extraction == PEP_CISS_EXTRACTION_CAA || ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
     416          98 :         for (i=0;i<nv;i++) pep->eigr[i] = (pep->eigr[i]*radius+center)*rgscale;
     417             :       }
     418          13 :       PetscCall(PetscFree3(fl1,inside,rr));
     419          13 :       PetscCall(BVSetActiveColumns(pep->V,0,nv));
     420          13 :       PetscCall(DSVectors(pep->ds,DS_MAT_X,NULL,NULL));
     421          13 :       if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) {
     422           2 :         PetscCall(BVSumQuadrature(ctx->S,ctx->Y,ctx->M,ctx->L,ctx->L,ctx->weight,ctx->pp,contour->scatterin,contour->subcomm,contour->npoints,ctx->useconj));
     423           2 :         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
     424           2 :         PetscCall(BVCopy(ctx->S,pep->V));
     425           2 :         PetscCall(DSGetMat(pep->ds,DS_MAT_X,&X));
     426           2 :         PetscCall(BVMultInPlace(ctx->S,X,0,pep->nconv));
     427           2 :         PetscCall(BVMultInPlace(pep->V,X,0,pep->nconv));
     428           2 :         PetscCall(DSRestoreMat(pep->ds,DS_MAT_X,&X));
     429             :       } else {
     430          11 :         PetscCall(DSGetMat(pep->ds,DS_MAT_X,&X));
     431          11 :         PetscCall(BVMultInPlace(ctx->S,X,0,pep->nconv));
     432          11 :         PetscCall(DSRestoreMat(pep->ds,DS_MAT_X,&X));
     433          11 :         PetscCall(BVSetActiveColumns(ctx->S,0,pep->nconv));
     434          11 :         PetscCall(BVCopy(ctx->S,pep->V));
     435             :       }
     436             :       max_error = 0.0;
     437         114 :       for (i=0;i<pep->nconv;i++) {
     438         101 :         PetscCall(BVGetColumn(pep->V,i,&si));
     439         101 :         PetscCall(VecNormalize(si,NULL));
     440         101 :         PetscCall(PEPComputeResidualNorm_Private(pep,pep->eigr[i],0,si,NULL,pep->work,&error));
     441         101 :         PetscCall((*pep->converged)(pep,pep->eigr[i],0,error,&error,pep->convergedctx));
     442         101 :         PetscCall(BVRestoreColumn(pep->V,i,&si));
     443         132 :         max_error = PetscMax(max_error,error);
     444             :       }
     445          13 :       if (max_error <= pep->tol) pep->reason = PEP_CONVERGED_TOL;
     446           0 :       else if (pep->its > pep->max_it) pep->reason = PEP_DIVERGED_ITS;
     447             :       else {
     448           0 :         if (pep->nconv > ctx->L) nv = pep->nconv;
     449           0 :         else if (ctx->L > nv) nv = ctx->L;
     450           0 :         nv = PetscMin(nv,ctx->L*ctx->M);
     451           0 :         PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,nv,ctx->L,NULL,&M));
     452           0 :         PetscCall(MatSetRandom(M,rand));
     453           0 :         PetscCall(BVSetActiveColumns(ctx->S,0,nv));
     454           0 :         PetscCall(BVMultInPlace(ctx->S,M,0,ctx->L));
     455           0 :         PetscCall(MatDestroy(&M));
     456           0 :         PetscCall(BVSetActiveColumns(ctx->S,0,ctx->L));
     457           0 :         PetscCall(BVSetActiveColumns(ctx->V,0,ctx->L));
     458           0 :         PetscCall(BVCopy(ctx->S,ctx->V));
     459           0 :         if (contour->pA) PetscCall(BVScatter(ctx->V,ctx->pV,contour->scatterin,contour->xdup));
     460          26 :         PetscCall(PEPCISSSolve(pep,ctx->J,(contour->pA)?ctx->pV:ctx->V,0,ctx->L));
     461             :       }
     462             :     }
     463             :   }
     464          13 :   PetscCall(PetscFree2(Mu,H0));
     465          13 :   if (ctx->extraction == PEP_CISS_EXTRACTION_HANKEL) PetscCall(PetscFree(H1));
     466          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     467             : }
     468             : 
     469           9 : static PetscErrorCode PEPCISSSetSizes_CISS(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
     470             : {
     471           9 :   PEP_CISS       *ctx = (PEP_CISS*)pep->data;
     472           9 :   PetscInt       oN,oL,oM,oLmax,onpart;
     473           9 :   PetscMPIInt    size;
     474             : 
     475           9 :   PetscFunctionBegin;
     476           9 :   oN = ctx->N;
     477           9 :   if (ip == PETSC_DETERMINE) {
     478           0 :     if (ctx->N!=32) { ctx->N =32; ctx->M = ctx->N/4; }
     479           9 :   } else if (ip != PETSC_CURRENT) {
     480           9 :     PetscCheck(ip>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be > 0");
     481           9 :     PetscCheck(ip%2==0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ip argument must be an even number");
     482           9 :     if (ctx->N!=ip) { ctx->N = ip; ctx->M = ctx->N/4; }
     483             :   }
     484           9 :   oL = ctx->L;
     485           9 :   if (bs == PETSC_DETERMINE) {
     486           1 :     ctx->L = 16;
     487           8 :   } else if (bs != PETSC_CURRENT) {
     488           8 :     PetscCheck(bs>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The bs argument must be > 0");
     489           8 :     ctx->L = bs;
     490             :   }
     491           9 :   oM = ctx->M;
     492           9 :   if (ms == PETSC_DETERMINE) {
     493           1 :     ctx->M = ctx->N/4;
     494           8 :   } else if (ms != PETSC_CURRENT) {
     495           8 :     PetscCheck(ms>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be > 0");
     496           8 :     PetscCheck(ms<=ctx->N,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The ms argument must be less than or equal to the number of integration points");
     497           8 :     ctx->M = ms;
     498             :   }
     499           9 :   onpart = ctx->npart;
     500           9 :   if (npart == PETSC_DETERMINE) {
     501           0 :     ctx->npart = 1;
     502           9 :   } else if (npart != PETSC_CURRENT) {
     503           9 :     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size));
     504           9 :     PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
     505           9 :     ctx->npart = npart;
     506             :   }
     507           9 :   oLmax = ctx->L_max;
     508           9 :   if (bsmax == PETSC_DETERMINE) {
     509           1 :     ctx->L_max = 64;
     510           8 :   } else if (bsmax != PETSC_CURRENT) {
     511           8 :     PetscCheck(bsmax>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The bsmax argument must be > 0");
     512           8 :     ctx->L_max = PetscMax(bsmax,ctx->L);
     513             :   }
     514           9 :   if (onpart != ctx->npart || oN != ctx->N || realmats != ctx->isreal) {
     515           6 :     PetscCall(SlepcContourDataDestroy(&ctx->contour));
     516           6 :     PetscCall(PetscInfo(pep,"Resetting the contour data structure due to a change of parameters\n"));
     517           6 :     pep->state = PEP_STATE_INITIAL;
     518             :   }
     519           9 :   ctx->isreal = realmats;
     520           9 :   if (oL != ctx->L || oM != ctx->M || oLmax != ctx->L_max) pep->state = PEP_STATE_INITIAL;
     521           9 :   PetscFunctionReturn(PETSC_SUCCESS);
     522             : }
     523             : 
     524             : /*@
     525             :    PEPCISSSetSizes - Sets the values of various size parameters in the CISS solver.
     526             : 
     527             :    Logically Collective
     528             : 
     529             :    Input Parameters:
     530             : +  pep   - the polynomial eigensolver context
     531             : .  ip    - number of integration points
     532             : .  bs    - block size
     533             : .  ms    - moment size
     534             : .  npart - number of partitions when splitting the communicator
     535             : .  bsmax - max block size
     536             : -  realmats - all coefficient matrices of P(.) are real
     537             : 
     538             :    Options Database Keys:
     539             : +  -pep_ciss_integration_points - Sets the number of integration points
     540             : .  -pep_ciss_blocksize - Sets the block size
     541             : .  -pep_ciss_moments - Sets the moment size
     542             : .  -pep_ciss_partitions - Sets the number of partitions
     543             : .  -pep_ciss_maxblocksize - Sets the maximum block size
     544             : -  -pep_ciss_realmats - all coefficient matrices of P(.) are real
     545             : 
     546             :    Notes:
     547             :    For all integer arguments, you can use PETSC_CURRENT to keep the current value, and
     548             :    PETSC_DETERMINE to set them to a default value.
     549             : 
     550             :    The default number of partitions is 1. This means the internal KSP object is shared
     551             :    among all processes of the PEP communicator. Otherwise, the communicator is split
     552             :    into npart communicators, so that npart KSP solves proceed simultaneously.
     553             : 
     554             :    Level: advanced
     555             : 
     556             : .seealso: PEPCISSGetSizes()
     557             : @*/
     558           9 : PetscErrorCode PEPCISSSetSizes(PEP pep,PetscInt ip,PetscInt bs,PetscInt ms,PetscInt npart,PetscInt bsmax,PetscBool realmats)
     559             : {
     560           9 :   PetscFunctionBegin;
     561           9 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     562          27 :   PetscValidLogicalCollectiveInt(pep,ip,2);
     563          27 :   PetscValidLogicalCollectiveInt(pep,bs,3);
     564          27 :   PetscValidLogicalCollectiveInt(pep,ms,4);
     565          27 :   PetscValidLogicalCollectiveInt(pep,npart,5);
     566          27 :   PetscValidLogicalCollectiveInt(pep,bsmax,6);
     567          27 :   PetscValidLogicalCollectiveBool(pep,realmats,7);
     568           9 :   PetscTryMethod(pep,"PEPCISSSetSizes_C",(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool),(pep,ip,bs,ms,npart,bsmax,realmats));
     569           9 :   PetscFunctionReturn(PETSC_SUCCESS);
     570             : }
     571             : 
     572          13 : static PetscErrorCode PEPCISSGetSizes_CISS(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
     573             : {
     574          13 :   PEP_CISS *ctx = (PEP_CISS*)pep->data;
     575             : 
     576          13 :   PetscFunctionBegin;
     577          13 :   if (ip) *ip = ctx->N;
     578          13 :   if (bs) *bs = ctx->L;
     579          13 :   if (ms) *ms = ctx->M;
     580          13 :   if (npart) *npart = ctx->npart;
     581          13 :   if (bsmax) *bsmax = ctx->L_max;
     582          13 :   if (realmats) *realmats = ctx->isreal;
     583          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     584             : }
     585             : 
     586             : /*@
     587             :    PEPCISSGetSizes - Gets the values of various size parameters in the CISS solver.
     588             : 
     589             :    Not Collective
     590             : 
     591             :    Input Parameter:
     592             : .  pep - the polynomial eigensolver context
     593             : 
     594             :    Output Parameters:
     595             : +  ip    - number of integration points
     596             : .  bs    - block size
     597             : .  ms    - moment size
     598             : .  npart - number of partitions when splitting the communicator
     599             : .  bsmax - max block size
     600             : -  realmats - all coefficient matrices of P(.) are real
     601             : 
     602             :    Level: advanced
     603             : 
     604             : .seealso: PEPCISSSetSizes()
     605             : @*/
     606          13 : PetscErrorCode PEPCISSGetSizes(PEP pep,PetscInt *ip,PetscInt *bs,PetscInt *ms,PetscInt *npart,PetscInt *bsmax,PetscBool *realmats)
     607             : {
     608          13 :   PetscFunctionBegin;
     609          13 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     610          13 :   PetscUseMethod(pep,"PEPCISSGetSizes_C",(PEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*),(pep,ip,bs,ms,npart,bsmax,realmats));
     611          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     612             : }
     613             : 
     614           1 : static PetscErrorCode PEPCISSSetThreshold_CISS(PEP pep,PetscReal delta,PetscReal spur)
     615             : {
     616           1 :   PEP_CISS *ctx = (PEP_CISS*)pep->data;
     617             : 
     618           1 :   PetscFunctionBegin;
     619           1 :   if (delta == (PetscReal)PETSC_DETERMINE) {
     620           0 :     ctx->delta = SLEPC_DEFAULT_TOL*1e-4;
     621           1 :   } else if (delta != (PetscReal)PETSC_CURRENT) {
     622           1 :     PetscCheck(delta>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The delta argument must be > 0.0");
     623           1 :     ctx->delta = delta;
     624             :   }
     625           1 :   if (spur == (PetscReal)PETSC_DETERMINE) {
     626           0 :     ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
     627           1 :   } else if (spur != (PetscReal)PETSC_CURRENT) {
     628           1 :     PetscCheck(spur>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The spurious threshold argument must be > 0.0");
     629           1 :     ctx->spurious_threshold = spur;
     630             :   }
     631           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     632             : }
     633             : 
     634             : /*@
     635             :    PEPCISSSetThreshold - Sets the values of various threshold parameters in
     636             :    the CISS solver.
     637             : 
     638             :    Logically Collective
     639             : 
     640             :    Input Parameters:
     641             : +  pep   - the polynomial eigensolver context
     642             : .  delta - threshold for numerical rank
     643             : -  spur  - spurious threshold (to discard spurious eigenpairs)
     644             : 
     645             :    Options Database Keys:
     646             : +  -pep_ciss_delta - Sets the delta
     647             : -  -pep_ciss_spurious_threshold - Sets the spurious threshold
     648             : 
     649             :    Note:
     650             :    PETSC_CURRENT can be used to preserve the current value of any of the
     651             :    arguments, and PETSC_DETERMINE to set them to a default value.
     652             : 
     653             :    Level: advanced
     654             : 
     655             : .seealso: PEPCISSGetThreshold()
     656             : @*/
     657           1 : PetscErrorCode PEPCISSSetThreshold(PEP pep,PetscReal delta,PetscReal spur)
     658             : {
     659           1 :   PetscFunctionBegin;
     660           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     661           3 :   PetscValidLogicalCollectiveReal(pep,delta,2);
     662           3 :   PetscValidLogicalCollectiveReal(pep,spur,3);
     663           1 :   PetscTryMethod(pep,"PEPCISSSetThreshold_C",(PEP,PetscReal,PetscReal),(pep,delta,spur));
     664           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     665             : }
     666             : 
     667          13 : static PetscErrorCode PEPCISSGetThreshold_CISS(PEP pep,PetscReal *delta,PetscReal *spur)
     668             : {
     669          13 :   PEP_CISS *ctx = (PEP_CISS*)pep->data;
     670             : 
     671          13 :   PetscFunctionBegin;
     672          13 :   if (delta) *delta = ctx->delta;
     673          13 :   if (spur)  *spur = ctx->spurious_threshold;
     674          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     675             : }
     676             : 
     677             : /*@
     678             :    PEPCISSGetThreshold - Gets the values of various threshold parameters in
     679             :    the CISS solver.
     680             : 
     681             :    Not Collective
     682             : 
     683             :    Input Parameter:
     684             : .  pep - the polynomial eigensolver context
     685             : 
     686             :    Output Parameters:
     687             : +  delta - threshold for numerical rank
     688             : -  spur  - spurious threshold (to discard spurious eigenpairs)
     689             : 
     690             :    Level: advanced
     691             : 
     692             : .seealso: PEPCISSSetThreshold()
     693             : @*/
     694          13 : PetscErrorCode PEPCISSGetThreshold(PEP pep,PetscReal *delta,PetscReal *spur)
     695             : {
     696          13 :   PetscFunctionBegin;
     697          13 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     698          13 :   PetscUseMethod(pep,"PEPCISSGetThreshold_C",(PEP,PetscReal*,PetscReal*),(pep,delta,spur));
     699          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     700             : }
     701             : 
     702           3 : static PetscErrorCode PEPCISSSetRefinement_CISS(PEP pep,PetscInt inner,PetscInt blsize)
     703             : {
     704           3 :   PEP_CISS *ctx = (PEP_CISS*)pep->data;
     705             : 
     706           3 :   PetscFunctionBegin;
     707           3 :   if (inner == PETSC_DETERMINE) {
     708           0 :     ctx->refine_inner = 0;
     709           3 :   } else if (inner != PETSC_CURRENT) {
     710           3 :     PetscCheck(inner>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The refine inner argument must be >= 0");
     711           3 :     ctx->refine_inner = inner;
     712             :   }
     713           3 :   if (blsize == PETSC_DETERMINE) {
     714           0 :     ctx->refine_blocksize = 0;
     715           3 :   } else if (blsize != PETSC_CURRENT) {
     716           3 :     PetscCheck(blsize>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The refine blocksize argument must be >= 0");
     717           3 :     ctx->refine_blocksize = blsize;
     718             :   }
     719           3 :   PetscFunctionReturn(PETSC_SUCCESS);
     720             : }
     721             : 
     722             : /*@
     723             :    PEPCISSSetRefinement - Sets the values of various refinement parameters
     724             :    in the CISS solver.
     725             : 
     726             :    Logically Collective
     727             : 
     728             :    Input Parameters:
     729             : +  pep    - the polynomial eigensolver context
     730             : .  inner  - number of iterative refinement iterations (inner loop)
     731             : -  blsize - number of iterative refinement iterations (blocksize loop)
     732             : 
     733             :    Options Database Keys:
     734             : +  -pep_ciss_refine_inner - Sets number of inner iterations
     735             : -  -pep_ciss_refine_blocksize - Sets number of blocksize iterations
     736             : 
     737             :    Note:
     738             :    PETSC_CURRENT can be used to preserve the current value of any of the
     739             :    arguments, and PETSC_DETERMINE to set them to a default value.
     740             : 
     741             :    Level: advanced
     742             : 
     743             : .seealso: PEPCISSGetRefinement()
     744             : @*/
     745           3 : PetscErrorCode PEPCISSSetRefinement(PEP pep,PetscInt inner,PetscInt blsize)
     746             : {
     747           3 :   PetscFunctionBegin;
     748           3 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     749           9 :   PetscValidLogicalCollectiveInt(pep,inner,2);
     750           9 :   PetscValidLogicalCollectiveInt(pep,blsize,3);
     751           3 :   PetscTryMethod(pep,"PEPCISSSetRefinement_C",(PEP,PetscInt,PetscInt),(pep,inner,blsize));
     752           3 :   PetscFunctionReturn(PETSC_SUCCESS);
     753             : }
     754             : 
     755          13 : static PetscErrorCode PEPCISSGetRefinement_CISS(PEP pep,PetscInt *inner,PetscInt *blsize)
     756             : {
     757          13 :   PEP_CISS *ctx = (PEP_CISS*)pep->data;
     758             : 
     759          13 :   PetscFunctionBegin;
     760          13 :   if (inner)  *inner = ctx->refine_inner;
     761          13 :   if (blsize) *blsize = ctx->refine_blocksize;
     762          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     763             : }
     764             : 
     765             : /*@
     766             :    PEPCISSGetRefinement - Gets the values of various refinement parameters
     767             :    in the CISS solver.
     768             : 
     769             :    Not Collective
     770             : 
     771             :    Input Parameter:
     772             : .  pep - the polynomial eigensolver context
     773             : 
     774             :    Output Parameters:
     775             : +  inner  - number of iterative refinement iterations (inner loop)
     776             : -  blsize - number of iterative refinement iterations (blocksize loop)
     777             : 
     778             :    Level: advanced
     779             : 
     780             : .seealso: PEPCISSSetRefinement()
     781             : @*/
     782          13 : PetscErrorCode PEPCISSGetRefinement(PEP pep, PetscInt *inner, PetscInt *blsize)
     783             : {
     784          13 :   PetscFunctionBegin;
     785          13 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     786          13 :   PetscUseMethod(pep,"PEPCISSGetRefinement_C",(PEP,PetscInt*,PetscInt*),(pep,inner,blsize));
     787          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     788             : }
     789             : 
     790           4 : static PetscErrorCode PEPCISSSetExtraction_CISS(PEP pep,PEPCISSExtraction extraction)
     791             : {
     792           4 :   PEP_CISS *ctx = (PEP_CISS*)pep->data;
     793             : 
     794           4 :   PetscFunctionBegin;
     795           4 :   if (ctx->extraction != extraction) {
     796           3 :     ctx->extraction = extraction;
     797           3 :     pep->state      = PEP_STATE_INITIAL;
     798             :   }
     799           4 :   PetscFunctionReturn(PETSC_SUCCESS);
     800             : }
     801             : 
     802             : /*@
     803             :    PEPCISSSetExtraction - Sets the extraction technique used in the CISS solver.
     804             : 
     805             :    Logically Collective
     806             : 
     807             :    Input Parameters:
     808             : +  pep        - the polynomial eigensolver context
     809             : -  extraction - the extraction technique
     810             : 
     811             :    Options Database Key:
     812             : .  -pep_ciss_extraction - Sets the extraction technique (either 'ritz', 'hankel' or 'caa')
     813             : 
     814             :    Notes:
     815             :    By default, the Rayleigh-Ritz extraction is used (PEP_CISS_EXTRACTION_RITZ).
     816             : 
     817             :    If the 'hankel' or the 'caa' option is specified (PEP_CISS_EXTRACTION_HANKEL or
     818             :    PEP_CISS_EXTRACTION_CAA), then the Block Hankel method, or the Communication-avoiding
     819             :    Arnoldi method, respectively, is used for extracting eigenpairs.
     820             : 
     821             :    Level: advanced
     822             : 
     823             : .seealso: PEPCISSGetExtraction(), PEPCISSExtraction
     824             : @*/
     825           4 : PetscErrorCode PEPCISSSetExtraction(PEP pep,PEPCISSExtraction extraction)
     826             : {
     827           4 :   PetscFunctionBegin;
     828           4 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     829          12 :   PetscValidLogicalCollectiveEnum(pep,extraction,2);
     830           4 :   PetscTryMethod(pep,"PEPCISSSetExtraction_C",(PEP,PEPCISSExtraction),(pep,extraction));
     831           4 :   PetscFunctionReturn(PETSC_SUCCESS);
     832             : }
     833             : 
     834           1 : static PetscErrorCode PEPCISSGetExtraction_CISS(PEP pep,PEPCISSExtraction *extraction)
     835             : {
     836           1 :   PEP_CISS *ctx = (PEP_CISS*)pep->data;
     837             : 
     838           1 :   PetscFunctionBegin;
     839           1 :   *extraction = ctx->extraction;
     840           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     841             : }
     842             : 
     843             : /*@
     844             :    PEPCISSGetExtraction - Gets the extraction technique used in the CISS solver.
     845             : 
     846             :    Not Collective
     847             : 
     848             :    Input Parameter:
     849             : .  pep - the polynomial eigensolver context
     850             : 
     851             :    Output Parameters:
     852             : .  extraction - extraction technique
     853             : 
     854             :    Level: advanced
     855             : 
     856             : .seealso: PEPCISSSetExtraction() PEPCISSExtraction
     857             : @*/
     858           1 : PetscErrorCode PEPCISSGetExtraction(PEP pep,PEPCISSExtraction *extraction)
     859             : {
     860           1 :   PetscFunctionBegin;
     861           1 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     862           1 :   PetscAssertPointer(extraction,2);
     863           1 :   PetscUseMethod(pep,"PEPCISSGetExtraction_C",(PEP,PEPCISSExtraction*),(pep,extraction));
     864           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     865             : }
     866             : 
     867          13 : static PetscErrorCode PEPCISSGetKSPs_CISS(PEP pep,PetscInt *nsolve,KSP **ksp)
     868             : {
     869          13 :   PEP_CISS         *ctx = (PEP_CISS*)pep->data;
     870          13 :   SlepcContourData contour;
     871          13 :   PetscInt         i,nsplit;
     872          13 :   PC               pc;
     873          13 :   MPI_Comm         child;
     874             : 
     875          13 :   PetscFunctionBegin;
     876          13 :   if (!ctx->contour) {  /* initialize contour data structure first */
     877          13 :     PetscCall(RGCanUseConjugates(pep->rg,ctx->isreal,&ctx->useconj));
     878          13 :     PetscCall(SlepcContourDataCreate(ctx->useconj?ctx->N/2:ctx->N,ctx->npart,(PetscObject)pep,&ctx->contour));
     879             :   }
     880          13 :   contour = ctx->contour;
     881          13 :   if (!contour->ksp) {
     882          13 :     PetscCall(PetscMalloc1(contour->npoints,&contour->ksp));
     883          13 :     PetscCall(PEPGetST(pep,&pep->st));
     884          13 :     PetscCall(STGetSplitPreconditionerInfo(pep->st,&nsplit,NULL));
     885          13 :     PetscCall(PetscSubcommGetChild(contour->subcomm,&child));
     886         385 :     for (i=0;i<contour->npoints;i++) {
     887         372 :       PetscCall(KSPCreate(child,&contour->ksp[i]));
     888         372 :       PetscCall(PetscObjectIncrementTabLevel((PetscObject)contour->ksp[i],(PetscObject)pep,1));
     889         372 :       PetscCall(KSPSetOptionsPrefix(contour->ksp[i],((PetscObject)pep)->prefix));
     890         372 :       PetscCall(KSPAppendOptionsPrefix(contour->ksp[i],"pep_ciss_"));
     891         372 :       PetscCall(PetscObjectSetOptions((PetscObject)contour->ksp[i],((PetscObject)pep)->options));
     892         372 :       PetscCall(KSPSetErrorIfNotConverged(contour->ksp[i],PETSC_TRUE));
     893         712 :       PetscCall(KSPSetTolerances(contour->ksp[i],SlepcDefaultTol(pep->tol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
     894         372 :       PetscCall(KSPGetPC(contour->ksp[i],&pc));
     895         372 :       if (nsplit) {
     896          64 :         PetscCall(KSPSetType(contour->ksp[i],KSPBCGS));
     897          64 :         PetscCall(PCSetType(pc,PCBJACOBI));
     898             :       } else {
     899         308 :         PetscCall(KSPSetType(contour->ksp[i],KSPPREONLY));
     900         372 :         PetscCall(PCSetType(pc,PCLU));
     901             :       }
     902             :     }
     903             :   }
     904          13 :   if (nsolve) *nsolve = contour->npoints;
     905          13 :   if (ksp)    *ksp    = contour->ksp;
     906          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     907             : }
     908             : 
     909             : /*@C
     910             :    PEPCISSGetKSPs - Retrieve the array of linear solver objects associated with
     911             :    the CISS solver.
     912             : 
     913             :    Collective
     914             : 
     915             :    Input Parameter:
     916             : .  pep - polynomial eigenvalue solver
     917             : 
     918             :    Output Parameters:
     919             : +  nsolve - number of solver objects
     920             : -  ksp - array of linear solver object
     921             : 
     922             :    Notes:
     923             :    The number of KSP solvers is equal to the number of integration points divided by
     924             :    the number of partitions. This value is halved in the case of real matrices with
     925             :    a region centered at the real axis.
     926             : 
     927             :    Level: advanced
     928             : 
     929             : .seealso: PEPCISSSetSizes()
     930             : @*/
     931          13 : PetscErrorCode PEPCISSGetKSPs(PEP pep,PetscInt *nsolve,KSP **ksp)
     932             : {
     933          13 :   PetscFunctionBegin;
     934          13 :   PetscValidHeaderSpecific(pep,PEP_CLASSID,1);
     935          13 :   PetscUseMethod(pep,"PEPCISSGetKSPs_C",(PEP,PetscInt*,KSP**),(pep,nsolve,ksp));
     936          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     937             : }
     938             : 
     939          13 : static PetscErrorCode PEPReset_CISS(PEP pep)
     940             : {
     941          13 :   PEP_CISS       *ctx = (PEP_CISS*)pep->data;
     942             : 
     943          13 :   PetscFunctionBegin;
     944          13 :   PetscCall(BVDestroy(&ctx->S));
     945          13 :   PetscCall(BVDestroy(&ctx->V));
     946          13 :   PetscCall(BVDestroy(&ctx->Y));
     947          13 :   PetscCall(SlepcContourDataReset(ctx->contour));
     948          13 :   PetscCall(MatDestroy(&ctx->J));
     949          13 :   PetscCall(BVDestroy(&ctx->pV));
     950          13 :   PetscCall(PetscFree(ctx->Psplit));
     951          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     952             : }
     953             : 
     954          13 : static PetscErrorCode PEPSetFromOptions_CISS(PEP pep,PetscOptionItems *PetscOptionsObject)
     955             : {
     956          13 :   PEP_CISS          *ctx = (PEP_CISS*)pep->data;
     957          13 :   PetscReal         r1,r2;
     958          13 :   PetscInt          i,i1,i2,i3,i4,i5,i6,i7;
     959          13 :   PetscBool         b1,flg,flg2,flg3,flg4,flg5,flg6;
     960          13 :   PEPCISSExtraction extraction;
     961             : 
     962          13 :   PetscFunctionBegin;
     963          13 :   PetscOptionsHeadBegin(PetscOptionsObject,"PEP CISS Options");
     964             : 
     965          13 :     PetscCall(PEPCISSGetSizes(pep,&i1,&i2,&i3,&i4,&i5,&b1));
     966          13 :     PetscCall(PetscOptionsInt("-pep_ciss_integration_points","Number of integration points","PEPCISSSetSizes",i1,&i1,&flg));
     967          13 :     PetscCall(PetscOptionsInt("-pep_ciss_blocksize","Block size","PEPCISSSetSizes",i2,&i2,&flg2));
     968          13 :     PetscCall(PetscOptionsInt("-pep_ciss_moments","Moment size","PEPCISSSetSizes",i3,&i3,&flg3));
     969          13 :     PetscCall(PetscOptionsInt("-pep_ciss_partitions","Number of partitions","PEPCISSSetSizes",i4,&i4,&flg4));
     970          13 :     PetscCall(PetscOptionsInt("-pep_ciss_maxblocksize","Maximum block size","PEPCISSSetSizes",i5,&i5,&flg5));
     971          13 :     PetscCall(PetscOptionsBool("-pep_ciss_realmats","True if all coefficient matrices of P(.) are real","PEPCISSSetSizes",b1,&b1,&flg6));
     972          13 :     if (flg || flg2 || flg3 || flg4 || flg5 || flg6) PetscCall(PEPCISSSetSizes(pep,i1,i2,i3,i4,i5,b1));
     973             : 
     974          13 :     PetscCall(PEPCISSGetThreshold(pep,&r1,&r2));
     975          13 :     PetscCall(PetscOptionsReal("-pep_ciss_delta","Threshold for numerical rank","PEPCISSSetThreshold",r1,&r1,&flg));
     976          13 :     PetscCall(PetscOptionsReal("-pep_ciss_spurious_threshold","Threshold for the spurious eigenpairs","PEPCISSSetThreshold",r2,&r2,&flg2));
     977          13 :     if (flg || flg2) PetscCall(PEPCISSSetThreshold(pep,r1,r2));
     978             : 
     979          13 :     PetscCall(PEPCISSGetRefinement(pep,&i6,&i7));
     980          13 :     PetscCall(PetscOptionsInt("-pep_ciss_refine_inner","Number of inner iterative refinement iterations","PEPCISSSetRefinement",i6,&i6,&flg));
     981          13 :     PetscCall(PetscOptionsInt("-pep_ciss_refine_blocksize","Number of blocksize iterative refinement iterations","PEPCISSSetRefinement",i7,&i7,&flg2));
     982          13 :     if (flg || flg2) PetscCall(PEPCISSSetRefinement(pep,i6,i7));
     983             : 
     984          13 :     PetscCall(PetscOptionsEnum("-pep_ciss_extraction","Extraction technique","PEPCISSSetExtraction",PEPCISSExtractions,(PetscEnum)ctx->extraction,(PetscEnum*)&extraction,&flg));
     985          13 :     if (flg) PetscCall(PEPCISSSetExtraction(pep,extraction));
     986             : 
     987          13 :   PetscOptionsHeadEnd();
     988             : 
     989          13 :   if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
     990          13 :   PetscCall(RGSetFromOptions(pep->rg)); /* this is necessary here to set useconj */
     991          13 :   if (!ctx->contour || !ctx->contour->ksp) PetscCall(PEPCISSGetKSPs(pep,NULL,NULL));
     992          13 :   PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Something went wrong with PEPCISSGetKSPs()");
     993         385 :   for (i=0;i<ctx->contour->npoints;i++) PetscCall(KSPSetFromOptions(ctx->contour->ksp[i]));
     994          13 :   PetscCall(PetscSubcommSetFromOptions(ctx->contour->subcomm));
     995          13 :   PetscFunctionReturn(PETSC_SUCCESS);
     996             : }
     997             : 
     998          13 : static PetscErrorCode PEPDestroy_CISS(PEP pep)
     999             : {
    1000          13 :   PEP_CISS       *ctx = (PEP_CISS*)pep->data;
    1001             : 
    1002          13 :   PetscFunctionBegin;
    1003          13 :   PetscCall(SlepcContourDataDestroy(&ctx->contour));
    1004          13 :   PetscCall(PetscFree4(ctx->weight,ctx->omega,ctx->pp,ctx->sigma));
    1005          13 :   PetscCall(PetscFree(pep->data));
    1006          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",NULL));
    1007          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",NULL));
    1008          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",NULL));
    1009          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",NULL));
    1010          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",NULL));
    1011          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",NULL));
    1012          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",NULL));
    1013          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",NULL));
    1014          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",NULL));
    1015          13 :   PetscFunctionReturn(PETSC_SUCCESS);
    1016             : }
    1017             : 
    1018           0 : static PetscErrorCode PEPView_CISS(PEP pep,PetscViewer viewer)
    1019             : {
    1020           0 :   PEP_CISS       *ctx = (PEP_CISS*)pep->data;
    1021           0 :   PetscBool      isascii;
    1022           0 :   PetscViewer    sviewer;
    1023             : 
    1024           0 :   PetscFunctionBegin;
    1025           0 :   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
    1026           0 :   if (isascii) {
    1027           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  sizes { integration points: %" PetscInt_FMT ", block size: %" PetscInt_FMT ", moment size: %" PetscInt_FMT ", partitions: %" PetscInt_FMT ", maximum block size: %" PetscInt_FMT " }\n",ctx->N,ctx->L,ctx->M,ctx->npart,ctx->L_max));
    1028           0 :     if (ctx->isreal) PetscCall(PetscViewerASCIIPrintf(viewer,"  exploiting symmetry of integration points\n"));
    1029           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  threshold { delta: %g, spurious threshold: %g }\n",(double)ctx->delta,(double)ctx->spurious_threshold));
    1030           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  iterative refinement  { inner: %" PetscInt_FMT ", blocksize: %" PetscInt_FMT " }\n",ctx->refine_inner, ctx->refine_blocksize));
    1031           0 :     PetscCall(PetscViewerASCIIPrintf(viewer,"  extraction: %s\n",PEPCISSExtractions[ctx->extraction]));
    1032           0 :     if (!ctx->contour || !ctx->contour->ksp) PetscCall(PEPCISSGetKSPs(pep,NULL,NULL));
    1033           0 :     PetscAssert(ctx->contour && ctx->contour->ksp,PetscObjectComm((PetscObject)pep),PETSC_ERR_PLIB,"Something went wrong with PEPCISSGetKSPs()");
    1034           0 :     PetscCall(PetscViewerASCIIPushTab(viewer));
    1035           0 :     if (ctx->npart>1 && ctx->contour->subcomm) {
    1036           0 :       PetscCall(PetscViewerGetSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
    1037           0 :       if (!ctx->contour->subcomm->color) PetscCall(KSPView(ctx->contour->ksp[0],sviewer));
    1038           0 :       PetscCall(PetscViewerFlush(sviewer));
    1039           0 :       PetscCall(PetscViewerRestoreSubViewer(viewer,ctx->contour->subcomm->child,&sviewer));
    1040             :       /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
    1041           0 :       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
    1042           0 :     } else PetscCall(KSPView(ctx->contour->ksp[0],viewer));
    1043           0 :     PetscCall(PetscViewerASCIIPopTab(viewer));
    1044             :   }
    1045           0 :   PetscFunctionReturn(PETSC_SUCCESS);
    1046             : }
    1047             : 
    1048          26 : static PetscErrorCode PEPSetDSType_CISS(PEP pep)
    1049             : {
    1050          26 :   PEP_CISS       *ctx = (PEP_CISS*)pep->data;
    1051             : 
    1052          26 :   PetscFunctionBegin;
    1053          26 :   switch (ctx->extraction) {
    1054          20 :     case PEP_CISS_EXTRACTION_RITZ:   PetscCall(DSSetType(pep->ds,DSPEP)); break;
    1055           4 :     case PEP_CISS_EXTRACTION_HANKEL: PetscCall(DSSetType(pep->ds,DSGNHEP)); break;
    1056           2 :     case PEP_CISS_EXTRACTION_CAA:    PetscCall(DSSetType(pep->ds,DSNHEP)); break;
    1057             :   }
    1058          26 :   PetscFunctionReturn(PETSC_SUCCESS);
    1059             : }
    1060             : 
    1061          13 : SLEPC_EXTERN PetscErrorCode PEPCreate_CISS(PEP pep)
    1062             : {
    1063          13 :   PEP_CISS       *ctx = (PEP_CISS*)pep->data;
    1064             : 
    1065          13 :   PetscFunctionBegin;
    1066          13 :   PetscCall(PetscNew(&ctx));
    1067          13 :   pep->data = ctx;
    1068             :   /* set default values of parameters */
    1069          13 :   ctx->N                  = 32;
    1070          13 :   ctx->L                  = 16;
    1071          13 :   ctx->M                  = ctx->N/4;
    1072          13 :   ctx->delta              = SLEPC_DEFAULT_TOL*1e-4;
    1073          13 :   ctx->L_max              = 64;
    1074          13 :   ctx->spurious_threshold = PetscSqrtReal(SLEPC_DEFAULT_TOL);
    1075          13 :   ctx->isreal             = PETSC_FALSE;
    1076          13 :   ctx->npart              = 1;
    1077             : 
    1078          13 :   pep->ops->solve          = PEPSolve_CISS;
    1079          13 :   pep->ops->setup          = PEPSetUp_CISS;
    1080          13 :   pep->ops->setfromoptions = PEPSetFromOptions_CISS;
    1081          13 :   pep->ops->reset          = PEPReset_CISS;
    1082          13 :   pep->ops->destroy        = PEPDestroy_CISS;
    1083          13 :   pep->ops->view           = PEPView_CISS;
    1084          13 :   pep->ops->setdstype      = PEPSetDSType_CISS;
    1085             : 
    1086          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetSizes_C",PEPCISSSetSizes_CISS));
    1087          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetSizes_C",PEPCISSGetSizes_CISS));
    1088          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetThreshold_C",PEPCISSSetThreshold_CISS));
    1089          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetThreshold_C",PEPCISSGetThreshold_CISS));
    1090          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetRefinement_C",PEPCISSSetRefinement_CISS));
    1091          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetRefinement_C",PEPCISSGetRefinement_CISS));
    1092          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSSetExtraction_C",PEPCISSSetExtraction_CISS));
    1093          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetExtraction_C",PEPCISSGetExtraction_CISS));
    1094          13 :   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPCISSGetKSPs_C",PEPCISSGetKSPs_CISS));
    1095          13 :   PetscFunctionReturn(PETSC_SUCCESS);
    1096             : }

Generated by: LCOV version 1.14