Actual source code: ks-hamilt.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

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

 13:    Method: Hamiltonian Krylov-Schur

 15:    References:

 17:        [1] P. Benner et al., "A Hamiltonian Krylov-Schur-type method based on the
 18:            symplectic Lanczos process", Linear Algebra Appl. 435(3), 2011.

 20:        [2] D. Watkins, "The Matrix Eigenvalue Problem", SIAM, 2007.

 22: */
 23: #include <slepc/private/epsimpl.h>
 24: #include "krylovschur.h"
 25: #include <slepcblaslapack.h>

 27: /* J-orthogonalize vector [x1;x2] against first j vectors in U and V (full orthogonalization)
 28:    Coeffs c are updated, coeffs d are computed from scratch.
 29: */
 30: static PetscErrorCode Orthog_Hamilt(Vec x1,Vec x2,BV U1,BV U2,BV V1,BV V2,PetscInt j,PetscScalar *c,PetscScalar *d,PetscScalar *w,PetscBool *breakdown)
 31: {
 32:   PetscInt i;
 33:   Vec      Jx1,Jx2;

 35:   PetscFunctionBegin;
 36:   PetscCall(BVSetActiveColumns(U1,0,j));
 37:   PetscCall(BVSetActiveColumns(U2,0,j));
 38:   PetscCall(BVSetActiveColumns(V1,0,j));
 39:   PetscCall(BVSetActiveColumns(V2,0,j));
 40:   PetscCall(VecDuplicate(x1,&Jx1));
 41:   PetscCall(VecDuplicate(x2,&Jx2));
 42:   PetscCall(VecCopy(x2,Jx1));
 43:   PetscCall(VecCopy(x1,Jx2));
 44:   PetscCall(VecScale(Jx2,-1.0));
 45:   PetscCall(VecConjugate(Jx1));
 46:   PetscCall(VecConjugate(Jx2));
 47:   for (i=0;i<j;i++) w[j+i] = c[i]; /* Copy initial values of c */
 48:   /* c = -V.'*J*x  computed as -conj(c) = V'*conj(J*x) */
 49:   PetscCall(BVDotVecBegin(V1,Jx1,c));
 50:   PetscCall(BVDotVecBegin(V2,Jx2,w));
 51:   PetscCall(BVDotVecEnd(V1,Jx1,c));
 52:   PetscCall(BVDotVecEnd(V2,Jx2,w));
 53:   for (i=0;i<j;i++) c[i] = -PetscConj(c[i]+w[i]);
 54:   /* e = U.'*J*x  computed as conj(c) = U'*conj(J*x) */
 55:   PetscCall(BVDotVecBegin(U1,Jx1,d));
 56:   PetscCall(BVDotVecBegin(U2,Jx2,w));
 57:   PetscCall(BVDotVecEnd(U1,Jx1,d));
 58:   PetscCall(BVDotVecEnd(U2,Jx2,w));
 59:   for (i=0;i<j;i++) d[i] = PetscConj(d[i]+w[i]);
 60:   /* x = x-U*c-V*d */
 61:   PetscCall(BVMultVec(U1,-1.0,1.0,x1,c));
 62:   PetscCall(BVMultVec(U2,-1.0,1.0,x2,c));
 63:   PetscCall(BVMultVec(V1,-1.0,1.0,x1,d));
 64:   PetscCall(BVMultVec(V2,-1.0,1.0,x2,d));
 65:   for (i=0;i<j;i++) c[i] += w[j+i]; /* Add initial values of c */
 66:   PetscCall(VecDestroy(&Jx1));
 67:   PetscCall(VecDestroy(&Jx2));
 68:   PetscFunctionReturn(PETSC_SUCCESS);
 69: }

 71: /* J-orthogonalize vector [x1;x2] against first j vectors in U and V (local+full orthogonalization)*/
 72: static PetscErrorCode OrthogonalizeVector_Hamilt(Vec x1,Vec x2,BV U1,BV U2,BV V1,BV V2,PetscInt j,PetscReal *alpha,PetscReal *beta,PetscInt k,PetscScalar *h,PetscBool *breakdown)
 73: {
 74:   PetscScalar p,p1,p2;
 75:   Vec         v1,v2;
 76:   PetscInt    i,l;

 78:   PetscFunctionBegin;
 79:   PetscCall(PetscArrayzero(h,4*j));

 81:   /* Local orthogonalization */
 82:   l = j==k+1?0:j-2;  /* 1st column to orthogonalize against */
 83:   /* compute p=-V(:,j).'*J*x */
 84:   PetscCall(BVGetColumn(V1,j-1,&v1));
 85:   PetscCall(BVGetColumn(V2,j-1,&v2));
 86:   PetscCall(VecTDotBegin(x1,v2,&p1));
 87:   PetscCall(VecTDotBegin(x2,v1,&p2));
 88:   PetscCall(VecTDotEnd(x1,v2,&p1));
 89:   PetscCall(VecTDotEnd(x2,v1,&p2));
 90:   PetscCall(BVRestoreColumn(V1,j-1,&v1));
 91:   PetscCall(BVRestoreColumn(V2,j-1,&v2));
 92:   p = p1-p2;
 93:   for (i=l; i<j-1; i++) h[i] = beta[i];
 94:   h[j-1] = p;
 95:   /* x = x - U(:,l:j)*h(l:j) */
 96:   PetscCall(BVSetActiveColumns(U1,l,j));
 97:   PetscCall(BVSetActiveColumns(U2,l,j));
 98:   PetscCall(BVMultVec(U1,-1.0,1.0,x1,h+l));
 99:   PetscCall(BVMultVec(U2,-1.0,1.0,x2,h+l));

101:   /* J-orthogonalize x (full orthogonalization) */
102:   PetscCall(Orthog_Hamilt(x1,x2,U1,U2,V1,V2,j,h,h+j,h+2*j,breakdown));
103:   alpha[j-1] = PetscRealPart(h[j-1]);
104:   PetscFunctionReturn(PETSC_SUCCESS);
105: }

107: /* Normalize u,v so that u.'*J*v=1 */
108: static PetscErrorCode NormalizeVector_Hamilt(Vec u1,Vec u2,Vec v1,Vec v2,PetscReal *rr,PetscReal *ss,PetscBool *breakdown)
109: {
110:   PetscScalar p,p1,p2;
111:   PetscReal   r,s;

113:   PetscFunctionBegin;
114:   PetscCall(VecTDotBegin(v1,u2,&p1));
115:   PetscCall(VecTDotBegin(v2,u1,&p2));
116:   PetscCall(VecTDotEnd(v1,u2,&p1));
117:   PetscCall(VecTDotEnd(v2,u1,&p2));
118:   p = p2-p1;
119:   if (breakdown) {
120:     *breakdown = p==0.0? PETSC_TRUE: PETSC_FALSE;
121:     if (*breakdown) PetscFunctionReturn(PETSC_SUCCESS);
122:   }
123:   r = PetscSqrtReal(PetscAbsScalar(p));
124:   s = PetscSign(PetscRealPart(p));
125:   PetscCall(VecScale(u1,1.0/r));
126:   PetscCall(VecScale(u2,1.0/r));
127:   PetscCall(VecScale(v1,s/r));
128:   PetscCall(VecScale(v2,s/r));
129:   if (rr) *rr = r;
130:   if (ss) *ss = s;
131:   PetscFunctionReturn(PETSC_SUCCESS);
132: }

134: static PetscErrorCode EPSHamiltonianKS(EPS eps,BV U,BV V,PetscReal *alpha,PetscReal *beta,PetscReal *omega,PetscInt k,PetscInt *M,PetscBool *symmlost,PetscBool *breakdown)
135: {
136:   PetscInt       i,j,m = *M,ld,l,work_len=4*m;
137:   Vec            u,v,u1,u2,v1,v2;
138:   Mat            H;
139:   IS             is[2];
140:   BV             U1,U2,V1,V2;
141:   PetscScalar    *hwork,lhwork[100];
142:   PetscReal      sym=0.0,fro=0.0;
143:   PetscBLASInt   j2_,one=1;

145:   PetscFunctionBegin;
146:   PetscCall(DSGetLeadingDimension(eps->ds,&ld));
147:   PetscCall(DSGetDimensions(eps->ds,NULL,&l,NULL,NULL));
148:   if (work_len > 100) PetscCall(PetscMalloc1(work_len,&hwork));
149:   else hwork = lhwork;
150:   PetscCall(STGetMatrix(eps->st,0,&H));
151:   PetscCall(MatNestGetISs(H,is,NULL));
152:   PetscCall(BVGetSplitRows(U,is[0],is[1],&U1,&U2));
153:   PetscCall(BVGetSplitRows(V,is[0],is[1],&V1,&V2));

155:   /* normalize initial vector */
156:   if (k==0) {
157:     if (eps->nini==0) PetscCall(BVSetRandomColumn(eps->V,0));
158:     PetscCall(BVGetColumn(U,0,&u));
159:     PetscCall(BVGetColumn(V,0,&v));
160:     PetscCall(STApply(eps->st,u,v));
161:     PetscCall(BVRestoreColumn(U,0,&u));
162:     PetscCall(BVRestoreColumn(V,0,&v));
163:     PetscCall(BVGetColumn(U1,0,&u1));
164:     PetscCall(BVGetColumn(U2,0,&u2));
165:     PetscCall(BVGetColumn(V1,0,&v1));
166:     PetscCall(BVGetColumn(V2,0,&v2));
167:     PetscCall(NormalizeVector_Hamilt(u1,u2,v1,v2,NULL,&omega[0],breakdown));
168:     PetscCheck(!*breakdown,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Breakdown in Hamiltonian Krylov-Schur");
169:     PetscCall(BVRestoreColumn(U1,0,&u1));
170:     PetscCall(BVRestoreColumn(U2,0,&u2));
171:     PetscCall(BVRestoreColumn(V1,0,&v1));
172:     PetscCall(BVRestoreColumn(V2,0,&v2));
173:   }

175:   for (j=k;j<m;j++) {
176:     PetscCall(BVGetColumn(U,j+1,&u));
177:     PetscCall(BVGetColumn(V,j,&v));
178:     PetscCall(STApply(eps->st,v,u));
179:     PetscCall(BVRestoreColumn(V,j,&v));
180:     PetscCall(BVGetColumn(U1,j+1,&u1));
181:     PetscCall(BVGetColumn(U2,j+1,&u2));
182:     PetscCall(OrthogonalizeVector_Hamilt(u1,u2,U1,U2,V1,V2,j+1,alpha,beta,k,hwork,breakdown));
183:     PetscCall(BVGetColumn(V,j+1,&v));
184:     PetscCall(STApply(eps->st,u,v));
185:     PetscCall(BVRestoreColumn(U,j+1,&u));
186:     PetscCall(BVRestoreColumn(V,j+1,&v));
187:     PetscCall(BVGetColumn(V1,j+1,&v1));
188:     PetscCall(BVGetColumn(V2,j+1,&v2));
189:     PetscCall(NormalizeVector_Hamilt(u1,u2,v1,v2,&beta[j],&omega[j+1],breakdown));

191:     /* Update norm of asymmetry (sym) */
192:     if (j==k) {
193:       PetscReal *f;

195:       PetscCall(DSGetArrayReal(eps->ds,DS_MAT_T,&f));
196:       for (i=0;i<l;i++) hwork[i]  = 0.0;
197:       for (;i<j-1;i++)  hwork[i] -= f[2*ld+i];
198:       PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_T,&f));
199:     }
200:     if (j>0) {
201:       hwork[j-1] -= beta[j-1];
202:       hwork[j] = 0;
203:       PetscCall(PetscBLASIntCast(2*(j+1),&j2_));
204:       sym = SlepcAbs(BLASnrm2_(&j2_,hwork,&one),sym);
205:     }
206:     /* Update norm of tridiagonal elements (fro) */
207:     fro = SlepcAbs(fro,SlepcAbs(alpha[j],beta[j]));
208:     if (j>0) fro = SlepcAbs(fro,beta[j-1]);

210:     PetscCheck(!*breakdown,PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Breakdown in Hamiltonian Krylov-Schur");
211:     PetscCall(BVRestoreColumn(U1,j+1,&u1));
212:     PetscCall(BVRestoreColumn(U2,j+1,&u2));
213:     PetscCall(BVRestoreColumn(V1,j+1,&v1));
214:     PetscCall(BVRestoreColumn(V2,j+1,&v2));

216:     /* Check if relative asymmetry is too large */
217:     if (sym/fro>PetscMax(PETSC_SQRT_MACHINE_EPSILON,100*eps->tol)) {
218:       *symmlost = PETSC_TRUE;
219:       *M=j;
220:       break;
221:     }
222:   }
223:   PetscCall(BVRestoreSplitRows(U,is[0],is[1],&U1,&U2));
224:   PetscCall(BVRestoreSplitRows(V,is[0],is[1],&V1,&V2));
225:   if (work_len > 100) PetscCall(PetscFree(hwork));
226:   PetscFunctionReturn(PETSC_SUCCESS);
227: }

229: PetscErrorCode EPSSolve_KrylovSchur_Hamilt(EPS eps)
230: {
231:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
232:   PetscInt        i,k,l,ld,nv,t,nconv=0,nevsave;
233:   Mat             Q,W,D;
234:   Vec             vomega,vomegaold;
235:   BV              U,V;
236:   PetscReal       *a,*a2,*b,*omega,beta,u_norm;
237:   PetscBool       breakdown=PETSC_FALSE,symmlost=PETSC_FALSE;
238:   PetscComplex    eig;

240:   PetscFunctionBegin;
241:   PetscCall(DSGetLeadingDimension(eps->ds,&ld));

243:   /* Get the split bases */
244:   PetscCall(BVSetActiveColumns(eps->V,eps->ncv/2+1,eps->ncv+2));
245:   PetscCall(BVGetSplit(eps->V,&U,&V));

247:   nevsave  = eps->nev;
248:   eps->nev = (eps->nev+1)/2;
249:   l = 0;

251:   /* Restart loop */
252:   while (eps->reason == EPS_CONVERGED_ITERATING) {
253:     eps->its++;

255:     /* Compute an nv-step Krylov factorization */
256:     nv = PetscMin(eps->nconv+eps->mpd/2,eps->ncv/2);
257:     PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l));
258:     PetscCall(DSGetArrayReal(eps->ds,DS_MAT_T,&a));
259:     /* FIXME: low level access to array a internals. This is obscure */
260:     a2 = a + ld;
261:     b = a + 2*ld;
262:     PetscCall(DSGetArrayReal(eps->ds,DS_MAT_D,&omega));
263:     PetscCall(EPSHamiltonianKS(eps,U,V,a,b,omega,eps->nconv+l,&nv,&symmlost,&breakdown));
264:     for (i=eps->nconv+l; i<nv; i++)
265:       a2[i] = b[i];
266:     beta = b[nv-1];
267:     PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_T,&a));
268:     PetscCall(DSRestoreArrayReal(eps->ds,DS_MAT_D,&omega));
269:     if (symmlost) {
270:       eps->reason = EPS_DIVERGED_SYMMETRY_LOST;
271:       if (nv==eps->nconv+l) { eps->nconv = nconv; break; }
272:     }
273:     PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,eps->nconv+l));
274:     PetscCall(DSSetState(eps->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE));
275:     PetscCall(BVSetActiveColumns(U,eps->nconv,nv));
276:     PetscCall(BVSetActiveColumns(V,eps->nconv,nv));

278:     /* Save a copy of the signature */
279:     PetscCall(DSGetMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));
280:     PetscCall(VecDuplicate(vomega,&vomegaold));
281:     PetscCall(VecCopy(vomega,vomegaold));
282:     PetscCall(DSRestoreMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));

284:     /* Solve projected problem */
285:     PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
286:     PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL));
287:     PetscCall(DSUpdateExtraRow(eps->ds));
288:     PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));

290:     /* Check convergence */
291:     PetscCall(DSGetDimensions(eps->ds,NULL,NULL,NULL,&t));
292:     PetscCall(BVNormColumn(U,nv,NORM_2,&u_norm));
293:     PetscCall(EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,0.0,u_norm,&k));
294:     EPSSetCtxThreshold(eps,eps->eigr,eps->eigi,k);
295:     if (!symmlost) PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx));
296:     nconv = k;

298:     /* Update l */
299:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown || k==nv) l = 0;
300:     else {
301:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
302:       l = PetscMin(l,t);
303:       PetscCall(DSGetTruncateSize(eps->ds,k,t,&l));
304:     }
305:     if (!ctx->lock && l>0) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
306:     if (l) PetscCall(PetscInfo(eps,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l));

308:     if (eps->reason == EPS_CONVERGED_ITERATING) {
309:       PetscCheck(!breakdown,PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Breakdown in Hamiltonian Krylov-Schur (beta=%g)",(double)beta);
310:       /* Prepare the Rayleigh quotient for restart */
311:       PetscCall(DSTruncate(eps->ds,k+l,PETSC_FALSE));
312:     }
313:     /* Update the corresponding vectors
314:        U(:,idx) = U*Omega*Q(:,idx)*Omega2
315:        V(:,idx) = V*Q(:,idx) */
316:     PetscCall(DSGetMat(eps->ds,DS_MAT_Q,&Q));
317:     PetscCall(MatDuplicate(Q,MAT_COPY_VALUES,&W));
318:     PetscCall(DSGetMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));
319:     PetscCall(MatDiagonalScale(W,vomegaold,vomega));
320:     PetscCall(DSRestoreMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));
321:     PetscCall(VecDestroy(&vomegaold));
322:     PetscCall(BVMultInPlace(U,W,eps->nconv,k+l));
323:     PetscCall(MatDestroy(&W));
324:     PetscCall(BVMultInPlace(V,Q,eps->nconv,k+l));
325:     PetscCall(DSRestoreMat(eps->ds,DS_MAT_Q,&Q));

327:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
328:       PetscCall(BVCopyColumn(U,nv,k+l));
329:       PetscCall(BVCopyColumn(V,nv,k+l));
330:       if (eps->stop==EPS_STOP_THRESHOLD && nv-k<5) {  /* reallocate */
331:         eps->ncv = eps->mpd+k;
332:         PetscCall(BVRestoreSplit(eps->V,&U,&V));
333:         PetscCall(EPSReallocateSolution(eps,eps->ncv+2));
334:         PetscCall(BVSetActiveColumns(eps->V,eps->ncv/2+1,eps->ncv+2));
335:         PetscCall(BVGetSplit(eps->V,&U,&V));
336:         for (i=nv;i<eps->ncv;i++) eps->perm[i] = i;
337:         PetscCall(DSReallocate(eps->ds,eps->ncv/2+1));
338:         PetscCall(DSGetLeadingDimension(eps->ds,&ld));
339:       }
340:     }
341:     eps->nconv = k;
342:     for (i=0;i<nv;i++) {
343: #if defined(PETSC_USE_COMPLEX)
344:       eig = PetscSqrtScalar(eps->eigr[i]);
345: #else
346:       eig = PetscSqrtComplex(PetscCMPLX(eps->eigr[i],eps->eigi[i]));
347: #endif
348:       eps->eigr[i] = PetscRealPartComplex(eig);
349:       eps->eigi[i] = PetscImaginaryPartComplex(eig);
350:     }
351:     PetscCall(EPSMonitor(eps,eps->its,nconv,eps->eigr,eps->eigi,eps->errest,nv));
352:   }

354:   eps->nev = nevsave;

356:   PetscCall(DSTruncate(eps->ds,eps->nconv,PETSC_TRUE));
357:   PetscCall(BVRestoreSplit(eps->V,&U,&V));
358:   PetscFunctionReturn(PETSC_SUCCESS);
359: }

361: static PetscErrorCode EPSComputeVectors_Hamilt(EPS eps)
362: {
363:   Mat         X,W,D;
364:   Vec         vomega,vomegacopy;
365:   BV          U,V;
366:   PetscInt    n;

368:   PetscFunctionBegin;
369:   /* Get the split bases */
370:   PetscCall(BVSetActiveColumns(eps->V,eps->ncv/2+1,eps->ncv+2));
371:   PetscCall(BVGetSplit(eps->V,&U,&V));
372:   PetscCall(BVSetActiveColumns(eps->V,0,eps->nconv));

374:   /* Update bases again, needed due to 2x2 blocks */
375:   PetscCall(DSGetDimensions(eps->ds,&n,NULL,NULL,NULL));
376:   PetscCall(DSGetMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));
377:   PetscCall(VecDuplicate(vomega,&vomegacopy));
378:   PetscCall(VecCopy(vomega,vomegacopy));
379:   PetscCall(DSRestoreMatAndColumn(eps->ds,DS_MAT_D,0,&D,&vomega));
380:   PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
381:   PetscCall(DSGetMat(eps->ds,DS_MAT_X,&X));
382:   PetscCall(MatDuplicate(X,MAT_COPY_VALUES,&W));
383:   PetscCall(MatDiagonalScale(W,vomegacopy,NULL));
384:   PetscCall(VecDestroy(&vomegacopy));
385:   PetscCall(BVSetActiveColumns(U,0,n));
386:   PetscCall(BVMultInPlace(U,W,0,n));
387:   PetscCall(MatDestroy(&W));
388:   PetscCall(BVSetActiveColumns(V,0,n));
389:   PetscCall(BVMultInPlace(V,X,0,n));
390:   PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&X));

392:   PetscCall(BVRestoreSplit(eps->V,&U,&V));
393:   PetscFunctionReturn(PETSC_SUCCESS);
394: }

396: PetscErrorCode EPSSetUp_KrylovSchur_Hamilt(EPS eps)
397: {
398:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
399:   PetscBool       flg,sinvert;

401:   PetscFunctionBegin;
402:   PetscCheck(eps->problem_type==EPS_HAMILT,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Problem type should be Hamiltonian");
403:   EPSCheckUnsupportedCondition(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_EXTRACTION | EPS_FEATURE_BALANCE,PETSC_TRUE," with Hamiltonian structure");
404:   PetscCall(EPSSetDimensions_Default(eps,&eps->nev,&eps->ncv,&eps->mpd));
405:   PetscCheck(eps->ncv<=eps->nev+eps->mpd,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must not be larger than nev+mpd");
406:   if (eps->max_it==PETSC_DETERMINE) eps->max_it = PetscMax(100,2*eps->n/eps->ncv)*((eps->stop==EPS_STOP_THRESHOLD)?10:1);

408:   PetscCall(PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STSINVERT,STSHIFT,""));
409:   PetscCheck(flg,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hamiltonian Krylov-Schur only supports shift and shift-and-invert ST");
410:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&sinvert));
411:   PetscCheck(!sinvert,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hamiltonian Krylov-Schur does not currently support shift-and-invert");
412:   if (!eps->which) {
413:     if (sinvert) eps->which = EPS_TARGET_MAGNITUDE;
414:     else eps->which = EPS_LARGEST_MAGNITUDE;
415:   }

417:   if (!ctx->keep) ctx->keep = 0.5;
418:   PetscCall(STSetStructured(eps->st,PETSC_FALSE));

420:   PetscCall(EPSAllocateSolution(eps,2));
421:   eps->ops->solve = EPSSolve_KrylovSchur_Hamilt;
422:   eps->ops->computevectors = EPSComputeVectors_Hamilt;
423:   PetscCall(DSSetType(eps->ds,DSGHIEP));
424:   PetscCall(DSSetCompact(eps->ds,PETSC_TRUE));
425:   PetscCall(DSSetExtraRow(eps->ds,PETSC_TRUE));
426:   PetscCall(DSAllocate(eps->ds,eps->ncv/2+1));
427:   PetscFunctionReturn(PETSC_SUCCESS);
428: }