Actual source code: ptoar.c
 
   slepc-3.19.1 2023-05-15
   
  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 polynomial eigensolver: "toar"
 13:    Method: TOAR
 15:    Algorithm:
 17:        Two-Level Orthogonal Arnoldi.
 19:    References:
 21:        [1] Y. Su, J. Zhang and Z. Bai, "A compact Arnoldi algorithm for
 22:            polynomial eigenvalue problems", talk presented at RANMEP 2008.
 24:        [2] C. Campos and J.E. Roman, "Parallel Krylov solvers for the
 25:            polynomial eigenvalue problem in SLEPc", SIAM J. Sci. Comput.
 26:            38(5):S385-S411, 2016.
 28:        [3] D. Lu, Y. Su and Z. Bai, "Stability analysis of the two-level
 29:            orthogonal Arnoldi procedure", SIAM J. Matrix Anal. App.
 30:            37(1):195-214, 2016.
 31: */
 33: #include <slepc/private/pepimpl.h>
 34: #include "../src/pep/impls/krylov/pepkrylov.h"
 35: #include <slepcblaslapack.h>
 37: static PetscBool  cited = PETSC_FALSE;
 38: static const char citation[] =
 39:   "@Article{slepc-pep,\n"
 40:   "   author = \"C. Campos and J. E. Roman\",\n"
 41:   "   title = \"Parallel {Krylov} solvers for the polynomial eigenvalue problem in {SLEPc}\",\n"
 42:   "   journal = \"{SIAM} J. Sci. Comput.\",\n"
 43:   "   volume = \"38\",\n"
 44:   "   number = \"5\",\n"
 45:   "   pages = \"S385--S411\",\n"
 46:   "   year = \"2016,\"\n"
 47:   "   doi = \"https://doi.org/10.1137/15M1022458\"\n"
 48:   "}\n";
 50: PetscErrorCode PEPSetUp_TOAR(PEP pep)
 51: {
 52:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
 53:   PetscBool      sinv,flg;
 54:   PetscInt       i;
 56:   PetscFunctionBegin;
 57:   PEPCheckShiftSinvert(pep);
 58:   PetscCall(PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd));
 59:   PetscCheck(ctx->lock || pep->mpd>=pep->ncv,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Should not use mpd parameter in non-locking variant");
 60:   if (pep->max_it==PETSC_DEFAULT) pep->max_it = PetscMax(100,2*(pep->nmat-1)*pep->n/pep->ncv);
 61:   if (!pep->which) PetscCall(PEPSetWhichEigenpairs_Default(pep));
 62:   PetscCheck(pep->which!=PEP_ALL,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues");
 63:   if (pep->problem_type!=PEP_GENERAL) PetscCall(PetscInfo(pep,"Problem type ignored, performing a non-symmetric linearization\n"));
 65:   if (!ctx->keep) ctx->keep = 0.5;
 67:   PetscCall(PEPAllocateSolution(pep,pep->nmat-1));
 68:   PetscCall(PEPSetWorkVecs(pep,3));
 69:   PetscCall(DSSetType(pep->ds,DSNHEP));
 70:   PetscCall(DSSetExtraRow(pep->ds,PETSC_TRUE));
 71:   PetscCall(DSAllocate(pep->ds,pep->ncv+1));
 73:   PetscCall(PEPBasisCoefficients(pep,pep->pbc));
 74:   PetscCall(STGetTransform(pep->st,&flg));
 75:   if (!flg) {
 76:     PetscCall(PetscFree(pep->solvematcoeffs));
 77:     PetscCall(PetscMalloc1(pep->nmat,&pep->solvematcoeffs));
 78:     PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv));
 79:     if (sinv) PetscCall(PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL));
 80:     else {
 81:       for (i=0;i<pep->nmat-1;i++) pep->solvematcoeffs[i] = 0.0;
 82:       pep->solvematcoeffs[pep->nmat-1] = 1.0;
 83:     }
 84:   }
 85:   PetscCall(BVDestroy(&ctx->V));
 86:   PetscCall(BVCreateTensor(pep->V,pep->nmat-1,&ctx->V));
 87:   PetscFunctionReturn(PETSC_SUCCESS);
 88: }
 90: /*
 91:   Extend the TOAR basis by applying the matrix operator
 92:   over a vector which is decomposed in the TOAR way
 93:   Input:
 94:     - pbc: array containing the polynomial basis coefficients
 95:     - S,V: define the latest Arnoldi vector (nv vectors in V)
 96:   Output:
 97:     - t: new vector extending the TOAR basis
 98:     - r: temporary coefficients to compute the TOAR coefficients
 99:          for the new Arnoldi vector
100:   Workspace: t_ (two vectors)
101: */
102: static PetscErrorCode PEPTOARExtendBasis(PEP pep,PetscBool sinvert,PetscScalar sigma,PetscScalar *S,PetscInt ls,PetscInt nv,BV V,Vec t,PetscScalar *r,PetscInt lr,Vec *t_)
103: {
104:   PetscInt       nmat=pep->nmat,deg=nmat-1,k,j,off=0,lss;
105:   Vec            v=t_[0],ve=t_[1],q=t_[2];
106:   PetscScalar    alpha=1.0,*ss,a;
107:   PetscReal      *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat;
108:   PetscBool      flg;
110:   PetscFunctionBegin;
111:   PetscCall(BVSetActiveColumns(pep->V,0,nv));
112:   PetscCall(STGetTransform(pep->st,&flg));
113:   if (sinvert) {
114:     for (j=0;j<nv;j++) {
115:       if (deg>1) r[lr+j] = S[j]/ca[0];
116:       if (deg>2) r[2*lr+j] = (S[ls+j]+(sigma-cb[1])*r[lr+j])/ca[1];
117:     }
118:     for (k=2;k<deg-1;k++) {
119:       for (j=0;j<nv;j++) r[(k+1)*lr+j] = (S[k*ls+j]+(sigma-cb[k])*r[k*lr+j]-cg[k]*r[(k-1)*lr+j])/ca[k];
120:     }
121:     k = deg-1;
122:     for (j=0;j<nv;j++) r[j] = (S[k*ls+j]+(sigma-cb[k])*r[k*lr+j]-cg[k]*r[(k-1)*lr+j])/ca[k];
123:     ss = r; lss = lr; off = 1; alpha = -1.0; a = pep->sfactor;
124:   } else {
125:     ss = S; lss = ls; off = 0; alpha = -ca[deg-1]; a = 1.0;
126:   }
127:   PetscCall(BVMultVec(V,1.0,0.0,v,ss+off*lss));
128:   if (PetscUnlikely(pep->Dr)) { /* balancing */
129:     PetscCall(VecPointwiseMult(v,v,pep->Dr));
130:   }
131:   PetscCall(STMatMult(pep->st,off,v,q));
132:   PetscCall(VecScale(q,a));
133:   for (j=1+off;j<deg+off-1;j++) {
134:     PetscCall(BVMultVec(V,1.0,0.0,v,ss+j*lss));
135:     if (PetscUnlikely(pep->Dr)) PetscCall(VecPointwiseMult(v,v,pep->Dr));
136:     PetscCall(STMatMult(pep->st,j,v,t));
137:     a *= pep->sfactor;
138:     PetscCall(VecAXPY(q,a,t));
139:   }
140:   if (sinvert) {
141:     PetscCall(BVMultVec(V,1.0,0.0,v,ss));
142:     if (PetscUnlikely(pep->Dr)) PetscCall(VecPointwiseMult(v,v,pep->Dr));
143:     PetscCall(STMatMult(pep->st,deg,v,t));
144:     a *= pep->sfactor;
145:     PetscCall(VecAXPY(q,a,t));
146:   } else {
147:     PetscCall(BVMultVec(V,1.0,0.0,ve,ss+(deg-1)*lss));
148:     if (PetscUnlikely(pep->Dr)) PetscCall(VecPointwiseMult(ve,ve,pep->Dr));
149:     a *= pep->sfactor;
150:     PetscCall(STMatMult(pep->st,deg-1,ve,t));
151:     PetscCall(VecAXPY(q,a,t));
152:     a *= pep->sfactor;
153:   }
154:   if (flg || !sinvert) alpha /= a;
155:   PetscCall(STMatSolve(pep->st,q,t));
156:   PetscCall(VecScale(t,alpha));
157:   if (!sinvert) {
158:     PetscCall(VecAXPY(t,cg[deg-1],v));
159:     PetscCall(VecAXPY(t,cb[deg-1],ve));
160:   }
161:   if (PetscUnlikely(pep->Dr)) PetscCall(VecPointwiseDivide(t,t,pep->Dr));
162:   PetscFunctionReturn(PETSC_SUCCESS);
163: }
165: /*
166:   Compute TOAR coefficients of the blocks of the new Arnoldi vector computed
167: */
168: static PetscErrorCode PEPTOARCoefficients(PEP pep,PetscBool sinvert,PetscScalar sigma,PetscInt nv,PetscScalar *S,PetscInt ls,PetscScalar *r,PetscInt lr,PetscScalar *x)
169: {
170:   PetscInt    k,j,nmat=pep->nmat,d=nmat-1;
171:   PetscReal   *ca=pep->pbc,*cb=pep->pbc+nmat,*cg=pep->pbc+2*nmat;
172:   PetscScalar t=1.0,tp=0.0,tt;
174:   PetscFunctionBegin;
175:   if (sinvert) {
176:     for (k=1;k<d;k++) {
177:       tt = t;
178:       t = ((sigma-cb[k-1])*t-cg[k-1]*tp)/ca[k-1]; /* k-th basis polynomial */
179:       tp = tt;
180:       for (j=0;j<=nv;j++) r[k*lr+j] += t*x[j];
181:     }
182:   } else {
183:     for (j=0;j<=nv;j++) r[j] = (cb[0]-sigma)*S[j]+ca[0]*S[ls+j];
184:     for (k=1;k<d-1;k++) {
185:       for (j=0;j<=nv;j++) r[k*lr+j] = (cb[k]-sigma)*S[k*ls+j]+ca[k]*S[(k+1)*ls+j]+cg[k]*S[(k-1)*ls+j];
186:     }
187:     if (sigma!=0.0) for (j=0;j<=nv;j++) r[(d-1)*lr+j] -= sigma*S[(d-1)*ls+j];
188:   }
189:   PetscFunctionReturn(PETSC_SUCCESS);
190: }
192: /*
193:   Compute a run of Arnoldi iterations dim(work)=ld
194: */
195: static PetscErrorCode PEPTOARrun(PEP pep,PetscScalar sigma,Mat A,PetscInt k,PetscInt *M,PetscReal *beta,PetscBool *breakdown,Vec *t_)
196: {
197:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
198:   PetscInt       j,m=*M,deg=pep->nmat-1,ld;
199:   PetscInt       ldh,lds,nqt,l;
200:   Vec            t;
201:   PetscReal      norm=0.0;
202:   PetscBool      flg,sinvert=PETSC_FALSE,lindep;
203:   PetscScalar    *H,*x,*S;
204:   Mat            MS;
206:   PetscFunctionBegin;
207:   *beta = 0.0;
208:   PetscCall(MatDenseGetArray(A,&H));
209:   PetscCall(MatDenseGetLDA(A,&ldh));
210:   PetscCall(BVTensorGetFactors(ctx->V,NULL,&MS));
211:   PetscCall(MatDenseGetArray(MS,&S));
212:   PetscCall(BVGetSizes(pep->V,NULL,NULL,&ld));
213:   lds = ld*deg;
214:   PetscCall(BVGetActiveColumns(pep->V,&l,&nqt));
215:   PetscCall(STGetTransform(pep->st,&flg));
216:   if (!flg) {
217:     /* spectral transformation handled by the solver */
218:     PetscCall(PetscObjectTypeCompareAny((PetscObject)pep->st,&flg,STSINVERT,STSHIFT,""));
219:     PetscCheck(flg,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"ST type not supported for TOAR without transforming matrices");
220:     PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinvert));
221:   }
222:   PetscCall(BVSetActiveColumns(ctx->V,0,m));
223:   for (j=k;j<m;j++) {
224:     /* apply operator */
225:     PetscCall(BVGetColumn(pep->V,nqt,&t));
226:     PetscCall(PEPTOARExtendBasis(pep,sinvert,sigma,S+j*lds,ld,nqt,pep->V,t,S+(j+1)*lds,ld,t_));
227:     PetscCall(BVRestoreColumn(pep->V,nqt,&t));
229:     /* orthogonalize */
230:     if (sinvert) x = S+(j+1)*lds;
231:     else x = S+(deg-1)*ld+(j+1)*lds;
232:     PetscCall(BVOrthogonalizeColumn(pep->V,nqt,x,&norm,&lindep));
233:     if (!lindep) {
234:       x[nqt] = norm;
235:       PetscCall(BVScaleColumn(pep->V,nqt,1.0/norm));
236:       nqt++;
237:     }
239:     PetscCall(PEPTOARCoefficients(pep,sinvert,sigma,nqt-1,S+j*lds,ld,S+(j+1)*lds,ld,x));
241:     /* level-2 orthogonalization */
242:     PetscCall(BVOrthogonalizeColumn(ctx->V,j+1,H+j*ldh,&norm,breakdown));
243:     H[j+1+ldh*j] = norm;
244:     if (PetscUnlikely(*breakdown)) {
245:       *M = j+1;
246:       break;
247:     }
248:     PetscCall(BVScaleColumn(ctx->V,j+1,1.0/norm));
249:     PetscCall(BVSetActiveColumns(pep->V,l,nqt));
250:   }
251:   *beta = norm;
252:   PetscCall(BVSetActiveColumns(ctx->V,0,*M));
253:   PetscCall(MatDenseRestoreArray(MS,&S));
254:   PetscCall(BVTensorRestoreFactors(ctx->V,NULL,&MS));
255:   PetscCall(MatDenseRestoreArray(A,&H));
256:   PetscFunctionReturn(PETSC_SUCCESS);
257: }
259: /*
260:   Computes T_j = phi_idx(T). In T_j and T_p are phi_{idx-1}(T)
261:    and phi_{idx-2}(T) respectively or null if idx=0,1.
262:    Tp and Tj are input/output arguments
263: */
264: static PetscErrorCode PEPEvaluateBasisM(PEP pep,PetscInt k,PetscScalar *T,PetscInt ldt,PetscInt idx,PetscScalar **Tp,PetscScalar **Tj)
265: {
266:   PetscInt       i;
267:   PetscReal      *ca,*cb,*cg;
268:   PetscScalar    *pt,g,a;
269:   PetscBLASInt   k_,ldt_;
271:   PetscFunctionBegin;
272:   if (idx==0) {
273:     PetscCall(PetscArrayzero(*Tj,k*k));
274:     PetscCall(PetscArrayzero(*Tp,k*k));
275:     for (i=0;i<k;i++) (*Tj)[i+i*k] = 1.0;
276:   } else {
277:     PetscCall(PetscBLASIntCast(ldt,&ldt_));
278:     PetscCall(PetscBLASIntCast(k,&k_));
279:     ca = pep->pbc; cb = pep->pbc+pep->nmat; cg = pep->pbc+2*pep->nmat;
280:     for (i=0;i<k;i++) T[i*ldt+i] -= cb[idx-1];
281:     a = 1/ca[idx-1];
282:     g = (idx==1)?0.0:-cg[idx-1]/ca[idx-1];
283:     PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&a,T,&ldt_,*Tj,&k_,&g,*Tp,&k_));
284:     pt = *Tj; *Tj = *Tp; *Tp = pt;
285:     for (i=0;i<k;i++) T[i*ldt+i] += cb[idx-1];
286:   }
287:   PetscFunctionReturn(PETSC_SUCCESS);
288: }
290: static PetscErrorCode PEPExtractInvariantPair(PEP pep,PetscScalar sigma,PetscInt sr,PetscInt k,PetscScalar *S,PetscInt ld,PetscInt deg,Mat HH)
291: {
292:   PetscInt       i,j,jj,ldh,lds,ldt,d=pep->nmat-1,idxcpy=0;
293:   PetscScalar    *H,*At,*Bt,*Hj,*Hp,*T,sone=1.0,g,a,*pM,*work;
294:   PetscBLASInt   k_,sr_,lds_,ldh_,info,*p,lwork,ldt_;
295:   PetscBool      transf=PETSC_FALSE,flg;
296:   PetscReal      norm,maxnrm,*rwork;
297:   BV             *R,Y;
298:   Mat            M,*A;
300:   PetscFunctionBegin;
301:   if (k==0) PetscFunctionReturn(PETSC_SUCCESS);
302:   PetscCall(MatDenseGetArray(HH,&H));
303:   PetscCall(MatDenseGetLDA(HH,&ldh));
304:   lds = deg*ld;
305:   PetscCall(PetscCalloc6(k,&p,sr*k,&At,k*k,&Bt,k*k,&Hj,k*k,&Hp,sr*k,&work));
306:   PetscCall(PetscBLASIntCast(sr,&sr_));
307:   PetscCall(PetscBLASIntCast(k,&k_));
308:   PetscCall(PetscBLASIntCast(lds,&lds_));
309:   PetscCall(PetscBLASIntCast(ldh,&ldh_));
310:   PetscCall(STGetTransform(pep->st,&flg));
311:   if (!flg) {
312:     PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&flg));
313:     if (flg || sigma!=0.0) transf=PETSC_TRUE;
314:   }
315:   if (transf) {
316:     PetscCall(PetscMalloc1(k*k,&T));
317:     ldt = k;
318:     for (i=0;i<k;i++) PetscCall(PetscArraycpy(T+k*i,H+i*ldh,k));
319:     if (flg) {
320:       PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
321:       PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&k_,&k_,T,&k_,p,&info));
322:       SlepcCheckLapackInfo("getrf",info);
323:       PetscCall(PetscBLASIntCast(sr*k,&lwork));
324:       PetscCallBLAS("LAPACKgetri",LAPACKgetri_(&k_,T,&k_,p,work,&lwork,&info));
325:       SlepcCheckLapackInfo("getri",info);
326:       PetscCall(PetscFPTrapPop());
327:     }
328:     if (sigma!=0.0) for (i=0;i<k;i++) T[i+k*i] += sigma;
329:   } else {
330:     T = H; ldt = ldh;
331:   }
332:   PetscCall(PetscBLASIntCast(ldt,&ldt_));
333:   switch (pep->extract) {
334:   case PEP_EXTRACT_NONE:
335:     break;
336:   case PEP_EXTRACT_NORM:
337:     if (pep->basis == PEP_BASIS_MONOMIAL) {
338:       PetscCall(PetscBLASIntCast(ldt,&ldt_));
339:       PetscCall(PetscMalloc1(k,&rwork));
340:       norm = LAPACKlange_("F",&k_,&k_,T,&ldt_,rwork);
341:       PetscCall(PetscFree(rwork));
342:       if (norm>1.0) idxcpy = d-1;
343:     } else {
344:       PetscCall(PetscBLASIntCast(ldt,&ldt_));
345:       PetscCall(PetscMalloc1(k,&rwork));
346:       maxnrm = 0.0;
347:       for (i=0;i<pep->nmat-1;i++) {
348:         PetscCall(PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj));
349:         norm = LAPACKlange_("F",&k_,&k_,Hj,&k_,rwork);
350:         if (norm > maxnrm) {
351:           idxcpy = i;
352:           maxnrm = norm;
353:         }
354:       }
355:       PetscCall(PetscFree(rwork));
356:     }
357:     if (idxcpy>0) {
358:       /* copy block idxcpy of S to the first one */
359:       for (j=0;j<k;j++) PetscCall(PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr));
360:     }
361:     break;
362:   case PEP_EXTRACT_RESIDUAL:
363:     PetscCall(STGetTransform(pep->st,&flg));
364:     if (flg) {
365:       PetscCall(PetscMalloc1(pep->nmat,&A));
366:       for (i=0;i<pep->nmat;i++) PetscCall(STGetMatrixTransformed(pep->st,i,A+i));
367:     } else A = pep->A;
368:     PetscCall(PetscMalloc1(pep->nmat-1,&R));
369:     for (i=0;i<pep->nmat-1;i++) PetscCall(BVDuplicateResize(pep->V,k,R+i));
370:     PetscCall(BVDuplicateResize(pep->V,sr,&Y));
371:     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,sr,k,NULL,&M));
372:     g = 0.0; a = 1.0;
373:     PetscCall(BVSetActiveColumns(pep->V,0,sr));
374:     for (j=0;j<pep->nmat;j++) {
375:       PetscCall(BVMatMult(pep->V,A[j],Y));
376:       PetscCall(PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj));
377:       for (i=0;i<pep->nmat-1;i++) {
378:         PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&sr_,&k_,&k_,&a,S+i*ld,&lds_,Hj,&k_,&g,At,&sr_));
379:         PetscCall(MatDenseGetArray(M,&pM));
380:         for (jj=0;jj<k;jj++) PetscCall(PetscArraycpy(pM+jj*sr,At+jj*sr,sr));
381:         PetscCall(MatDenseRestoreArray(M,&pM));
382:         PetscCall(BVMult(R[i],1.0,(i==0)?0.0:1.0,Y,M));
383:       }
384:     }
386:     /* frobenius norm */
387:     maxnrm = 0.0;
388:     for (i=0;i<pep->nmat-1;i++) {
389:       PetscCall(BVNorm(R[i],NORM_FROBENIUS,&norm));
390:       if (maxnrm > norm) {
391:         maxnrm = norm;
392:         idxcpy = i;
393:       }
394:     }
395:     if (idxcpy>0) {
396:       /* copy block idxcpy of S to the first one */
397:       for (j=0;j<k;j++) PetscCall(PetscArraycpy(S+j*lds,S+idxcpy*ld+j*lds,sr));
398:     }
399:     if (flg) PetscCall(PetscFree(A));
400:     for (i=0;i<pep->nmat-1;i++) PetscCall(BVDestroy(&R[i]));
401:     PetscCall(PetscFree(R));
402:     PetscCall(BVDestroy(&Y));
403:     PetscCall(MatDestroy(&M));
404:     break;
405:   case PEP_EXTRACT_STRUCTURED:
406:     for (j=0;j<k;j++) Bt[j+j*k] = 1.0;
407:     for (j=0;j<sr;j++) {
408:       for (i=0;i<k;i++) At[j*k+i] = PetscConj(S[i*lds+j]);
409:     }
410:     PetscCall(PEPEvaluateBasisM(pep,k,T,ldt,0,&Hp,&Hj));
411:     for (i=1;i<deg;i++) {
412:       PetscCall(PEPEvaluateBasisM(pep,k,T,ldt,i,&Hp,&Hj));
413:       PetscCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&sr_,&k_,&sone,Hj,&k_,S+i*ld,&lds_,&sone,At,&k_));
414:       PetscCallBLAS("BLASgemm",BLASgemm_("N","C",&k_,&k_,&k_,&sone,Hj,&k_,Hj,&k_,&sone,Bt,&k_));
415:     }
416:     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
417:     PetscCallBLAS("LAPACKgesv",LAPACKgesv_(&k_,&sr_,Bt,&k_,p,At,&k_,&info));
418:     PetscCall(PetscFPTrapPop());
419:     SlepcCheckLapackInfo("gesv",info);
420:     for (j=0;j<sr;j++) {
421:       for (i=0;i<k;i++) S[i*lds+j] = PetscConj(At[j*k+i]);
422:     }
423:     break;
424:   }
425:   if (transf) PetscCall(PetscFree(T));
426:   PetscCall(PetscFree6(p,At,Bt,Hj,Hp,work));
427:   PetscCall(MatDenseRestoreArray(HH,&H));
428:   PetscFunctionReturn(PETSC_SUCCESS);
429: }
431: PetscErrorCode PEPSolve_TOAR(PEP pep)
432: {
433:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
434:   PetscInt       i,j,k,l,nv=0,ld,lds,nq=0,nconv=0;
435:   PetscInt       nmat=pep->nmat,deg=nmat-1;
436:   PetscScalar    *S,sigma;
437:   PetscReal      beta;
438:   PetscBool      breakdown=PETSC_FALSE,flg,falselock=PETSC_FALSE,sinv=PETSC_FALSE;
439:   Mat            H,MS,MQ;
441:   PetscFunctionBegin;
442:   PetscCall(PetscCitationsRegister(citation,&cited));
443:   if (ctx->lock) {
444:     /* undocumented option to use a cheaper locking instead of the true locking */
445:     PetscCall(PetscOptionsGetBool(NULL,NULL,"-pep_toar_falselocking",&falselock,NULL));
446:   }
447:   PetscCall(STGetShift(pep->st,&sigma));
449:   /* update polynomial basis coefficients */
450:   PetscCall(STGetTransform(pep->st,&flg));
451:   if (pep->sfactor!=1.0) {
452:     for (i=0;i<nmat;i++) {
453:       pep->pbc[nmat+i] /= pep->sfactor;
454:       pep->pbc[2*nmat+i] /= pep->sfactor*pep->sfactor;
455:     }
456:     if (!flg) {
457:       pep->target /= pep->sfactor;
458:       PetscCall(RGPushScale(pep->rg,1.0/pep->sfactor));
459:       PetscCall(STScaleShift(pep->st,1.0/pep->sfactor));
460:       sigma /= pep->sfactor;
461:     } else {
462:       PetscCall(PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv));
463:       pep->target = sinv?pep->target*pep->sfactor:pep->target/pep->sfactor;
464:       PetscCall(RGPushScale(pep->rg,sinv?pep->sfactor:1.0/pep->sfactor));
465:       PetscCall(STScaleShift(pep->st,sinv?pep->sfactor:1.0/pep->sfactor));
466:     }
467:   }
469:   if (flg) sigma = 0.0;
471:   /* clean projected matrix (including the extra-arrow) */
472:   PetscCall(DSSetDimensions(pep->ds,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
473:   PetscCall(DSGetMat(pep->ds,DS_MAT_A,&H));
474:   PetscCall(MatZeroEntries(H));
475:   PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&H));
477:   /* Get the starting Arnoldi vector */
478:   PetscCall(BVTensorBuildFirstColumn(ctx->V,pep->nini));
480:   /* restart loop */
481:   l = 0;
482:   while (pep->reason == PEP_CONVERGED_ITERATING) {
483:     pep->its++;
485:     /* compute an nv-step Lanczos factorization */
486:     nv = PetscMax(PetscMin(nconv+pep->mpd,pep->ncv),nv);
487:     PetscCall(DSGetMat(pep->ds,DS_MAT_A,&H));
488:     PetscCall(PEPTOARrun(pep,sigma,H,pep->nconv+l,&nv,&beta,&breakdown,pep->work));
489:     PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&H));
490:     PetscCall(DSSetDimensions(pep->ds,nv,pep->nconv,pep->nconv+l));
491:     PetscCall(DSSetState(pep->ds,l?DS_STATE_RAW:DS_STATE_INTERMEDIATE));
492:     PetscCall(BVSetActiveColumns(ctx->V,pep->nconv,nv));
494:     /* solve projected problem */
495:     PetscCall(DSSolve(pep->ds,pep->eigr,pep->eigi));
496:     PetscCall(DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL));
497:     PetscCall(DSUpdateExtraRow(pep->ds));
498:     PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi));
500:     /* check convergence */
501:     PetscCall(PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k));
502:     PetscCall((*pep->stopping)(pep,pep->its,pep->max_it,k,pep->nev,&pep->reason,pep->stoppingctx));
504:     /* update l */
505:     if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
506:     else {
507:       l = (nv==k)?0:PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
508:       PetscCall(DSGetTruncateSize(pep->ds,k,nv,&l));
509:       if (!breakdown) {
510:         /* prepare the Rayleigh quotient for restart */
511:         PetscCall(DSTruncate(pep->ds,k+l,PETSC_FALSE));
512:       }
513:     }
514:     nconv = k;
515:     if (!ctx->lock && pep->reason == PEP_CONVERGED_ITERATING && !breakdown) { l += k; k = 0; } /* non-locking variant: reset no. of converged pairs */
516:     if (l) PetscCall(PetscInfo(pep,"Preparing to restart keeping l=%" PetscInt_FMT " vectors\n",l));
518:     /* update S */
519:     PetscCall(DSGetMat(pep->ds,DS_MAT_Q,&MQ));
520:     PetscCall(BVMultInPlace(ctx->V,MQ,pep->nconv,k+l));
521:     PetscCall(DSRestoreMat(pep->ds,DS_MAT_Q,&MQ));
523:     /* copy last column of S */
524:     PetscCall(BVCopyColumn(ctx->V,nv,k+l));
526:     if (PetscUnlikely(breakdown && pep->reason == PEP_CONVERGED_ITERATING)) {
527:       /* stop if breakdown */
528:       PetscCall(PetscInfo(pep,"Breakdown TOAR method (it=%" PetscInt_FMT " norm=%g)\n",pep->its,(double)beta));
529:       pep->reason = PEP_DIVERGED_BREAKDOWN;
530:     }
531:     if (pep->reason != PEP_CONVERGED_ITERATING) l--;
532:     /* truncate S */
533:     PetscCall(BVGetActiveColumns(pep->V,NULL,&nq));
534:     if (k+l+deg<=nq) {
535:       PetscCall(BVSetActiveColumns(ctx->V,pep->nconv,k+l+1));
536:       if (!falselock && ctx->lock) PetscCall(BVTensorCompress(ctx->V,k-pep->nconv));
537:       else PetscCall(BVTensorCompress(ctx->V,0));
538:     }
539:     pep->nconv = k;
540:     PetscCall(PEPMonitor(pep,pep->its,nconv,pep->eigr,pep->eigi,pep->errest,nv));
541:   }
542:   if (pep->nconv>0) {
543:     /* {V*S_nconv^i}_{i=0}^{d-1} has rank nconv instead of nconv+d-1. Force zeros in each S_nconv^i block */
544:     PetscCall(BVSetActiveColumns(ctx->V,0,pep->nconv));
545:     PetscCall(BVGetActiveColumns(pep->V,NULL,&nq));
546:     PetscCall(BVSetActiveColumns(pep->V,0,nq));
547:     if (nq>pep->nconv) {
548:       PetscCall(BVTensorCompress(ctx->V,pep->nconv));
549:       PetscCall(BVSetActiveColumns(pep->V,0,pep->nconv));
550:       nq = pep->nconv;
551:     }
553:     /* perform Newton refinement if required */
554:     if (pep->refine==PEP_REFINE_MULTIPLE && pep->rits>0) {
555:       /* extract invariant pair */
556:       PetscCall(BVTensorGetFactors(ctx->V,NULL,&MS));
557:       PetscCall(MatDenseGetArray(MS,&S));
558:       PetscCall(DSGetMat(pep->ds,DS_MAT_A,&H));
559:       PetscCall(BVGetSizes(pep->V,NULL,NULL,&ld));
560:       lds = deg*ld;
561:       PetscCall(PEPExtractInvariantPair(pep,sigma,nq,pep->nconv,S,ld,deg,H));
562:       PetscCall(DSRestoreMat(pep->ds,DS_MAT_A,&H));
563:       PetscCall(DSSetDimensions(pep->ds,pep->nconv,0,0));
564:       PetscCall(DSSetState(pep->ds,DS_STATE_RAW));
565:       PetscCall(PEPNewtonRefinement_TOAR(pep,sigma,&pep->rits,NULL,pep->nconv,S,lds));
566:       PetscCall(DSSolve(pep->ds,pep->eigr,pep->eigi));
567:       PetscCall(DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL));
568:       PetscCall(DSSynchronize(pep->ds,pep->eigr,pep->eigi));
569:       PetscCall(DSGetMat(pep->ds,DS_MAT_Q,&MQ));
570:       PetscCall(BVMultInPlace(ctx->V,MQ,0,pep->nconv));
571:       PetscCall(DSRestoreMat(pep->ds,DS_MAT_Q,&MQ));
572:       PetscCall(MatDenseRestoreArray(MS,&S));
573:       PetscCall(BVTensorRestoreFactors(ctx->V,NULL,&MS));
574:     }
575:   }
576:   PetscCall(STGetTransform(pep->st,&flg));
577:   if (pep->refine!=PEP_REFINE_MULTIPLE || pep->rits==0) {
578:     if (!flg) PetscTryTypeMethod(pep,backtransform);
579:     if (pep->sfactor!=1.0) {
580:       for (j=0;j<pep->nconv;j++) {
581:         pep->eigr[j] *= pep->sfactor;
582:         pep->eigi[j] *= pep->sfactor;
583:       }
584:       /* restore original values */
585:       for (i=0;i<pep->nmat;i++) {
586:         pep->pbc[pep->nmat+i] *= pep->sfactor;
587:         pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
588:       }
589:     }
590:   }
591:   /* restore original values */
592:   if (!flg) {
593:     pep->target *= pep->sfactor;
594:     PetscCall(STScaleShift(pep->st,pep->sfactor));
595:   } else {
596:     PetscCall(STScaleShift(pep->st,sinv?1.0/pep->sfactor:pep->sfactor));
597:     pep->target = (sinv)?pep->target/pep->sfactor:pep->target*pep->sfactor;
598:   }
599:   if (pep->sfactor!=1.0) PetscCall(RGPopScale(pep->rg));
601:   /* change the state to raw so that DSVectors() computes eigenvectors from scratch */
602:   PetscCall(DSSetDimensions(pep->ds,pep->nconv,0,0));
603:   PetscCall(DSSetState(pep->ds,DS_STATE_RAW));
604:   PetscFunctionReturn(PETSC_SUCCESS);
605: }
607: static PetscErrorCode PEPTOARSetRestart_TOAR(PEP pep,PetscReal keep)
608: {
609:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
611:   PetscFunctionBegin;
612:   if (keep==(PetscReal)PETSC_DEFAULT) ctx->keep = 0.5;
613:   else {
614:     PetscCheck(keep>=0.1 && keep<=0.9,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
615:     ctx->keep = keep;
616:   }
617:   PetscFunctionReturn(PETSC_SUCCESS);
618: }
620: /*@
621:    PEPTOARSetRestart - Sets the restart parameter for the TOAR
622:    method, in particular the proportion of basis vectors that must be kept
623:    after restart.
625:    Logically Collective
627:    Input Parameters:
628: +  pep  - the eigenproblem solver context
629: -  keep - the number of vectors to be kept at restart
631:    Options Database Key:
632: .  -pep_toar_restart - Sets the restart parameter
634:    Notes:
635:    Allowed values are in the range [0.1,0.9]. The default is 0.5.
637:    Level: advanced
639: .seealso: PEPTOARGetRestart()
640: @*/
641: PetscErrorCode PEPTOARSetRestart(PEP pep,PetscReal keep)
642: {
643:   PetscFunctionBegin;
646:   PetscTryMethod(pep,"PEPTOARSetRestart_C",(PEP,PetscReal),(pep,keep));
647:   PetscFunctionReturn(PETSC_SUCCESS);
648: }
650: static PetscErrorCode PEPTOARGetRestart_TOAR(PEP pep,PetscReal *keep)
651: {
652:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
654:   PetscFunctionBegin;
655:   *keep = ctx->keep;
656:   PetscFunctionReturn(PETSC_SUCCESS);
657: }
659: /*@
660:    PEPTOARGetRestart - Gets the restart parameter used in the TOAR method.
662:    Not Collective
664:    Input Parameter:
665: .  pep - the eigenproblem solver context
667:    Output Parameter:
668: .  keep - the restart parameter
670:    Level: advanced
672: .seealso: PEPTOARSetRestart()
673: @*/
674: PetscErrorCode PEPTOARGetRestart(PEP pep,PetscReal *keep)
675: {
676:   PetscFunctionBegin;
679:   PetscUseMethod(pep,"PEPTOARGetRestart_C",(PEP,PetscReal*),(pep,keep));
680:   PetscFunctionReturn(PETSC_SUCCESS);
681: }
683: static PetscErrorCode PEPTOARSetLocking_TOAR(PEP pep,PetscBool lock)
684: {
685:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
687:   PetscFunctionBegin;
688:   ctx->lock = lock;
689:   PetscFunctionReturn(PETSC_SUCCESS);
690: }
692: /*@
693:    PEPTOARSetLocking - Choose between locking and non-locking variants of
694:    the TOAR method.
696:    Logically Collective
698:    Input Parameters:
699: +  pep  - the eigenproblem solver context
700: -  lock - true if the locking variant must be selected
702:    Options Database Key:
703: .  -pep_toar_locking - Sets the locking flag
705:    Notes:
706:    The default is to lock converged eigenpairs when the method restarts.
707:    This behaviour can be changed so that all directions are kept in the
708:    working subspace even if already converged to working accuracy (the
709:    non-locking variant).
711:    Level: advanced
713: .seealso: PEPTOARGetLocking()
714: @*/
715: PetscErrorCode PEPTOARSetLocking(PEP pep,PetscBool lock)
716: {
717:   PetscFunctionBegin;
720:   PetscTryMethod(pep,"PEPTOARSetLocking_C",(PEP,PetscBool),(pep,lock));
721:   PetscFunctionReturn(PETSC_SUCCESS);
722: }
724: static PetscErrorCode PEPTOARGetLocking_TOAR(PEP pep,PetscBool *lock)
725: {
726:   PEP_TOAR *ctx = (PEP_TOAR*)pep->data;
728:   PetscFunctionBegin;
729:   *lock = ctx->lock;
730:   PetscFunctionReturn(PETSC_SUCCESS);
731: }
733: /*@
734:    PEPTOARGetLocking - Gets the locking flag used in the TOAR method.
736:    Not Collective
738:    Input Parameter:
739: .  pep - the eigenproblem solver context
741:    Output Parameter:
742: .  lock - the locking flag
744:    Level: advanced
746: .seealso: PEPTOARSetLocking()
747: @*/
748: PetscErrorCode PEPTOARGetLocking(PEP pep,PetscBool *lock)
749: {
750:   PetscFunctionBegin;
753:   PetscUseMethod(pep,"PEPTOARGetLocking_C",(PEP,PetscBool*),(pep,lock));
754:   PetscFunctionReturn(PETSC_SUCCESS);
755: }
757: PetscErrorCode PEPSetFromOptions_TOAR(PEP pep,PetscOptionItems *PetscOptionsObject)
758: {
759:   PetscBool      flg,lock;
760:   PetscReal      keep;
762:   PetscFunctionBegin;
763:   PetscOptionsHeadBegin(PetscOptionsObject,"PEP TOAR Options");
765:     PetscCall(PetscOptionsReal("-pep_toar_restart","Proportion of vectors kept after restart","PEPTOARSetRestart",0.5,&keep,&flg));
766:     if (flg) PetscCall(PEPTOARSetRestart(pep,keep));
768:     PetscCall(PetscOptionsBool("-pep_toar_locking","Choose between locking and non-locking variants","PEPTOARSetLocking",PETSC_FALSE,&lock,&flg));
769:     if (flg) PetscCall(PEPTOARSetLocking(pep,lock));
771:   PetscOptionsHeadEnd();
772:   PetscFunctionReturn(PETSC_SUCCESS);
773: }
775: PetscErrorCode PEPView_TOAR(PEP pep,PetscViewer viewer)
776: {
777:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
778:   PetscBool      isascii;
780:   PetscFunctionBegin;
781:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
782:   if (isascii) {
783:     PetscCall(PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep)));
784:     PetscCall(PetscViewerASCIIPrintf(viewer,"  using the %slocking variant\n",ctx->lock?"":"non-"));
785:   }
786:   PetscFunctionReturn(PETSC_SUCCESS);
787: }
789: PetscErrorCode PEPDestroy_TOAR(PEP pep)
790: {
791:   PEP_TOAR       *ctx = (PEP_TOAR*)pep->data;
793:   PetscFunctionBegin;
794:   PetscCall(BVDestroy(&ctx->V));
795:   PetscCall(PetscFree(pep->data));
796:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",NULL));
797:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",NULL));
798:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",NULL));
799:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",NULL));
800:   PetscFunctionReturn(PETSC_SUCCESS);
801: }
803: SLEPC_EXTERN PetscErrorCode PEPCreate_TOAR(PEP pep)
804: {
805:   PEP_TOAR       *ctx;
807:   PetscFunctionBegin;
808:   PetscCall(PetscNew(&ctx));
809:   pep->data = (void*)ctx;
811:   pep->lineariz = PETSC_TRUE;
812:   ctx->lock     = PETSC_TRUE;
814:   pep->ops->solve          = PEPSolve_TOAR;
815:   pep->ops->setup          = PEPSetUp_TOAR;
816:   pep->ops->setfromoptions = PEPSetFromOptions_TOAR;
817:   pep->ops->destroy        = PEPDestroy_TOAR;
818:   pep->ops->view           = PEPView_TOAR;
819:   pep->ops->backtransform  = PEPBackTransform_Default;
820:   pep->ops->computevectors = PEPComputeVectors_Default;
821:   pep->ops->extractvectors = PEPExtractVectors_TOAR;
823:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetRestart_C",PEPTOARSetRestart_TOAR));
824:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetRestart_C",PEPTOARGetRestart_TOAR));
825:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARSetLocking_C",PEPTOARSetLocking_TOAR));
826:   PetscCall(PetscObjectComposeFunction((PetscObject)pep,"PEPTOARGetLocking_C",PEPTOARGetLocking_TOAR));
827:   PetscFunctionReturn(PETSC_SUCCESS);
828: }