Actual source code: dvdcalcpairs.c
slepc-3.15.0 2021-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: "davidson"
13: Step: calculate the best eigenpairs in the subspace V
15: For that, performs these steps:
16: 1) Update W <- A * V
17: 2) Update H <- V' * W
18: 3) Obtain eigenpairs of H
19: 4) Select some eigenpairs
20: 5) Compute the Ritz pairs of the selected ones
21: */
23: #include "davidson.h"
24: #include <slepcblaslapack.h>
26: static PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d)
27: {
31: BVSetActiveColumns(d->eps->V,0,0);
32: if (d->W) { BVSetActiveColumns(d->W,0,0); }
33: BVSetActiveColumns(d->AX,0,0);
34: if (d->BX) { BVSetActiveColumns(d->BX,0,0); }
35: return(0);
36: }
38: static PetscErrorCode dvd_calcpairs_qz_d(dvdDashboard *d)
39: {
40: PetscErrorCode ierr;
43: BVDestroy(&d->W);
44: BVDestroy(&d->AX);
45: BVDestroy(&d->BX);
46: BVDestroy(&d->auxBV);
47: MatDestroy(&d->H);
48: MatDestroy(&d->G);
49: MatDestroy(&d->auxM);
50: SlepcVecPoolDestroy(&d->auxV);
51: PetscFree(d->nBds);
52: return(0);
53: }
55: /* in complex, d->size_H real auxiliary values are needed */
56: static PetscErrorCode dvd_calcpairs_projeig_solve(dvdDashboard *d)
57: {
58: PetscErrorCode ierr;
59: Vec v;
60: PetscScalar *pA;
61: const PetscScalar *pv;
62: PetscInt i,lV,kV,n,ld;
65: BVGetActiveColumns(d->eps->V,&lV,&kV);
66: n = kV-lV;
67: DSSetDimensions(d->eps->ds,n,0,0,0);
68: DSCopyMat(d->eps->ds,DS_MAT_A,0,0,d->H,lV,lV,n,n,PETSC_FALSE);
69: if (d->G) {
70: DSCopyMat(d->eps->ds,DS_MAT_B,0,0,d->G,lV,lV,n,n,PETSC_FALSE);
71: }
72: /* Set the signature on projected matrix B */
73: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
74: DSGetLeadingDimension(d->eps->ds,&ld);
75: DSGetArray(d->eps->ds,DS_MAT_B,&pA);
76: PetscArrayzero(pA,n*ld);
77: VecCreateSeq(PETSC_COMM_SELF,kV,&v);
78: BVGetSignature(d->eps->V,v);
79: VecGetArrayRead(v,&pv);
80: for (i=0;i<n;i++) {
81: pA[i+ld*i] = d->nBds[i] = PetscRealPart(pv[lV+i]);
82: }
83: VecRestoreArrayRead(v,&pv);
84: VecDestroy(&v);
85: DSRestoreArray(d->eps->ds,DS_MAT_B,&pA);
86: }
87: DSSetState(d->eps->ds,DS_STATE_RAW);
88: DSSolve(d->eps->ds,d->eigr,d->eigi);
89: return(0);
90: }
92: /*
93: A(lA:kA-1,lA:kA-1) <- Z(l:k-1)'*A(l:k-1,l:k-1)*Q(l,k-1), where k=l+kA-lA
94: */
95: static PetscErrorCode EPSXDUpdateProj(Mat Q,Mat Z,PetscInt l,Mat A,PetscInt lA,PetscInt kA,Mat aux)
96: {
97: PetscErrorCode ierr;
98: PetscScalar one=1.0,zero=0.0;
99: PetscInt i,j,dA_=kA-lA,m0,n0,ldA_,ldQ_,ldZ_,nQ_;
100: PetscBLASInt dA,nQ,ldA,ldQ,ldZ;
101: PetscScalar *pA,*pW;
102: const PetscScalar *pQ,*pZ;
103: PetscBool symm=PETSC_FALSE,set,flg;
106: MatGetSize(A,&m0,&n0); ldA_=m0;
107: if (m0!=n0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"A should be square");
108: if (lA<0 || lA>m0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial row, column in A");
109: if (kA<0 || kA<lA || kA>m0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid final row, column in A");
110: MatIsHermitianKnown(A,&set,&flg);
111: symm = set? flg: PETSC_FALSE;
112: MatGetSize(Q,&m0,&n0); ldQ_=nQ_=m0;
113: if (l<0 || l>n0 || l+dA_>n0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial column in Q");
114: MatGetSize(Z,&m0,&n0); ldZ_=m0;
115: if (l<0 || l>n0 || l+dA_>n0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial column in Z");
116: MatGetSize(aux,&m0,&n0);
117: if (m0*n0<nQ_*dA_) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"aux should be larger");
118: PetscBLASIntCast(dA_,&dA);
119: PetscBLASIntCast(nQ_,&nQ);
120: PetscBLASIntCast(ldA_,&ldA);
121: PetscBLASIntCast(ldQ_,&ldQ);
122: PetscBLASIntCast(ldZ_,&ldZ);
123: MatDenseGetArray(A,&pA);
124: MatDenseGetArrayRead(Q,&pQ);
125: if (Q!=Z) { MatDenseGetArrayRead(Z,&pZ); }
126: else pZ = pQ;
127: MatDenseGetArrayWrite(aux,&pW);
128: /* W = A*Q */
129: if (symm) {
130: /* symmetrize before multiplying */
131: for (i=lA+1;i<lA+nQ;i++) {
132: for (j=lA;j<i;j++) pA[i+j*ldA] = PetscConj(pA[j+i*ldA]);
133: }
134: }
135: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&nQ,&dA,&nQ,&one,&pA[ldA*lA+lA],&ldA,&pQ[ldQ*l+l],&ldQ,&zero,pW,&nQ));
136: /* A = Q'*W */
137: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&dA,&dA,&nQ,&one,&pZ[ldZ*l+l],&ldZ,pW,&nQ,&zero,&pA[ldA*lA+lA],&ldA));
138: MatDenseRestoreArray(A,&pA);
139: MatDenseRestoreArrayRead(Q,&pQ);
140: if (Q!=Z) { MatDenseRestoreArrayRead(Z,&pZ); }
141: MatDenseRestoreArrayWrite(aux,&pW);
142: return(0);
143: }
145: static PetscErrorCode dvd_calcpairs_updateproj(dvdDashboard *d)
146: {
148: Mat Q,Z;
149: PetscInt lV,kV;
150: PetscBool symm;
153: DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
154: if (d->W) { DSGetMat(d->eps->ds,DS_MAT_Z,&Z); }
155: else Z = Q;
156: BVGetActiveColumns(d->eps->V,&lV,&kV);
157: EPSXDUpdateProj(Q,Z,0,d->H,lV,lV+d->V_tra_e,d->auxM);
158: if (d->G) { EPSXDUpdateProj(Q,Z,0,d->G,lV,lV+d->V_tra_e,d->auxM); }
159: MatDestroy(&Q);
160: if (d->W) { MatDestroy(&Z); }
162: PetscObjectTypeCompareAny((PetscObject)d->eps->ds,&symm,DSHEP,DSGHIEP,DSGHEP,"");
163: if (d->V_tra_s==0 || symm) return(0);
164: /* Compute upper part of H (and G): H(0:l-1,l:k-1) <- W(0:l-1)' * AV(l:k-1), where
165: k=l+d->V_tra_s */
166: BVSetActiveColumns(d->W?d->W:d->eps->V,0,lV);
167: BVSetActiveColumns(d->AX,lV,lV+d->V_tra_s);
168: BVDot(d->AX,d->W?d->W:d->eps->V,d->H);
169: if (d->G) {
170: BVSetActiveColumns(d->BX?d->BX:d->eps->V,lV,lV+d->V_tra_s);
171: BVDot(d->BX?d->BX:d->eps->V,d->W?d->W:d->eps->V,d->G);
172: }
173: PetscObjectTypeCompare((PetscObject)d->eps->ds,DSGHEP,&symm);
174: if (!symm) {
175: BVSetActiveColumns(d->W?d->W:d->eps->V,lV,lV+d->V_tra_s);
176: BVSetActiveColumns(d->AX,0,lV);
177: BVDot(d->AX,d->W?d->W:d->eps->V,d->H);
178: if (d->G) {
179: BVSetActiveColumns(d->BX?d->BX:d->eps->V,0,lV);
180: BVDot(d->BX?d->BX:d->eps->V,d->W?d->W:d->eps->V,d->G);
181: }
182: }
183: BVSetActiveColumns(d->eps->V,lV,kV);
184: BVSetActiveColumns(d->AX,lV,kV);
185: if (d->BX) { BVSetActiveColumns(d->BX,lV,kV); }
186: if (d->W) { BVSetActiveColumns(d->W,lV,kV); }
187: if (d->W) { dvd_harm_updateproj(d); }
188: return(0);
189: }
191: /*
192: BV <- BV*MT
193: */
194: PETSC_STATIC_INLINE PetscErrorCode dvd_calcpairs_updateBV0_gen(dvdDashboard *d,BV bv,DSMatType mat)
195: {
197: PetscInt l,k,n;
198: Mat auxM;
201: BVGetActiveColumns(d->eps->V,&l,&k);
202: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&auxM);
203: MatZeroEntries(auxM);
204: DSGetDimensions(d->eps->ds,&n,NULL,NULL,NULL,NULL);
205: if (k-l!=n) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
206: DSCopyMat(d->eps->ds,mat,0,0,auxM,l,l,n,d->V_tra_e,PETSC_TRUE);
207: BVMultInPlace(bv,auxM,l,l+d->V_tra_e);
208: MatDestroy(&auxM);
209: return(0);
210: }
212: static PetscErrorCode dvd_calcpairs_proj(dvdDashboard *d)
213: {
215: PetscInt i,l,k;
216: Vec v1,v2;
217: PetscScalar *pv;
220: BVGetActiveColumns(d->eps->V,&l,&k);
221: /* Update AV, BV, W and the projected matrices */
222: /* 1. S <- S*MT */
223: if (d->V_tra_s != d->V_tra_e || d->V_tra_e > 0) {
224: dvd_calcpairs_updateBV0_gen(d,d->eps->V,DS_MAT_Q);
225: if (d->W) { dvd_calcpairs_updateBV0_gen(d,d->W,DS_MAT_Z); }
226: dvd_calcpairs_updateBV0_gen(d,d->AX,DS_MAT_Q);
227: if (d->BX) { dvd_calcpairs_updateBV0_gen(d,d->BX,DS_MAT_Q); }
228: dvd_calcpairs_updateproj(d);
229: /* Update signature */
230: if (d->nBds) {
231: VecCreateSeq(PETSC_COMM_SELF,l+d->V_tra_e,&v1);
232: BVSetActiveColumns(d->eps->V,0,l+d->V_tra_e);
233: BVGetSignature(d->eps->V,v1);
234: VecGetArray(v1,&pv);
235: for (i=0;i<d->V_tra_e;i++) pv[l+i] = d->nBds[i];
236: VecRestoreArray(v1,&pv);
237: BVSetSignature(d->eps->V,v1);
238: BVSetActiveColumns(d->eps->V,l,k);
239: VecDestroy(&v1);
240: }
241: k = l+d->V_tra_e;
242: l+= d->V_tra_s;
243: } else {
244: /* 2. V <- orth(V, V_new) */
245: dvd_orthV(d->eps->V,l+d->V_new_s,l+d->V_new_e);
246: /* 3. AV <- [AV A * V(V_new_s:V_new_e-1)] */
247: /* Check consistency */
248: if (k-l != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
249: for (i=l+d->V_new_s;i<l+d->V_new_e;i++) {
250: BVGetColumn(d->eps->V,i,&v1);
251: BVGetColumn(d->AX,i,&v2);
252: MatMult(d->A,v1,v2);
253: BVRestoreColumn(d->eps->V,i,&v1);
254: BVRestoreColumn(d->AX,i,&v2);
255: }
256: /* 4. BV <- [BV B * V(V_new_s:V_new_e-1)] */
257: if (d->BX) {
258: /* Check consistency */
259: if (k-l != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
260: for (i=l+d->V_new_s;i<l+d->V_new_e;i++) {
261: BVGetColumn(d->eps->V,i,&v1);
262: BVGetColumn(d->BX,i,&v2);
263: MatMult(d->B,v1,v2);
264: BVRestoreColumn(d->eps->V,i,&v1);
265: BVRestoreColumn(d->BX,i,&v2);
266: }
267: }
268: /* 5. W <- [W f(AV,BV)] */
269: if (d->W) {
270: d->calcpairs_W(d);
271: dvd_orthV(d->W,l+d->V_new_s,l+d->V_new_e);
272: }
273: /* 6. H <- W' * AX; G <- W' * BX */
274: BVSetActiveColumns(d->eps->V,l+d->V_new_s,l+d->V_new_e);
275: BVSetActiveColumns(d->AX,l+d->V_new_s,l+d->V_new_e);
276: if (d->BX) { BVSetActiveColumns(d->BX,l+d->V_new_s,l+d->V_new_e); }
277: if (d->W) { BVSetActiveColumns(d->W,l+d->V_new_s,l+d->V_new_e); }
278: BVMatProject(d->AX,NULL,d->W?d->W:d->eps->V,d->H);
279: if (d->G) { BVMatProject(d->BX?d->BX:d->eps->V,NULL,d->W?d->W:d->eps->V,d->G); }
280: BVSetActiveColumns(d->eps->V,l,k);
281: BVSetActiveColumns(d->AX,l,k);
282: if (d->BX) { BVSetActiveColumns(d->BX,l,k); }
283: if (d->W) { BVSetActiveColumns(d->W,l,k); }
285: /* Perform the transformation on the projected problem */
286: if (d->W) {
287: d->calcpairs_proj_trans(d);
288: }
289: k = l+d->V_new_e;
290: }
291: BVSetActiveColumns(d->eps->V,l,k);
292: BVSetActiveColumns(d->AX,l,k);
293: if (d->BX) { BVSetActiveColumns(d->BX,l,k); }
294: if (d->W) { BVSetActiveColumns(d->W,l,k); }
296: /* Solve the projected problem */
297: dvd_calcpairs_projeig_solve(d);
299: d->V_tra_s = d->V_tra_e = 0;
300: d->V_new_s = d->V_new_e;
301: return(0);
302: }
304: static PetscErrorCode dvd_calcpairs_apply_arbitrary(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscScalar *rr,PetscScalar *ri)
305: {
306: PetscInt i,k,ld;
307: PetscScalar *pX;
308: Vec *X,xr,xi;
310: #if defined(PETSC_USE_COMPLEX)
311: PetscInt N=1;
312: #else
313: PetscInt N=2,j;
314: #endif
317: /* Quick exit without neither arbitrary selection nor harmonic extraction */
318: if (!d->eps->arbitrary && !d->calcpairs_eig_backtrans) return(0);
320: /* Quick exit without arbitrary selection, but with harmonic extraction */
321: if (d->calcpairs_eig_backtrans) {
322: for (i=r_s; i<r_e; i++) {
323: d->calcpairs_eig_backtrans(d,d->eigr[i],d->eigi[i],&rr[i-r_s],&ri[i-r_s]);
324: }
325: }
326: if (!d->eps->arbitrary) return(0);
328: SlepcVecPoolGetVecs(d->auxV,N,&X);
329: DSGetLeadingDimension(d->eps->ds,&ld);
330: for (i=r_s;i<r_e;i++) {
331: k = i;
332: DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
333: DSGetArray(d->eps->ds,DS_MAT_X,&pX);
334: dvd_improvex_compute_X(d,i,k+1,X,pX,ld);
335: DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);
336: #if !defined(PETSC_USE_COMPLEX)
337: if (d->nX[i] != 1.0) {
338: for (j=i;j<k+1;j++) {
339: VecScale(X[j-i],1.0/d->nX[i]);
340: }
341: }
342: xr = X[0];
343: xi = X[1];
344: if (i == k) {
345: VecSet(xi,0.0);
346: }
347: #else
348: xr = X[0];
349: xi = NULL;
350: if (d->nX[i] != 1.0) {
351: VecScale(xr,1.0/d->nX[i]);
352: }
353: #endif
354: (d->eps->arbitrary)(rr[i-r_s],ri[i-r_s],xr,xi,&rr[i-r_s],&ri[i-r_s],d->eps->arbitraryctx);
355: #if !defined(PETSC_USE_COMPLEX)
356: if (i != k) {
357: rr[i+1-r_s] = rr[i-r_s];
358: ri[i+1-r_s] = ri[i-r_s];
359: i++;
360: }
361: #endif
362: }
363: SlepcVecPoolRestoreVecs(d->auxV,N,&X);
364: return(0);
365: }
367: static PetscErrorCode dvd_calcpairs_selectPairs(dvdDashboard *d,PetscInt n)
368: {
369: PetscInt k,lV,kV,nV;
370: PetscScalar *rr,*ri;
374: BVGetActiveColumns(d->eps->V,&lV,&kV);
375: nV = kV - lV;
376: n = PetscMin(n,nV);
377: if (n <= 0) return(0);
378: /* Put the best n pairs at the beginning. Useful for restarting */
379: if (d->eps->arbitrary || d->calcpairs_eig_backtrans) {
380: PetscMalloc1(nV,&rr);
381: PetscMalloc1(nV,&ri);
382: dvd_calcpairs_apply_arbitrary(d,0,nV,rr,ri);
383: } else {
384: rr = d->eigr;
385: ri = d->eigi;
386: }
387: k = n;
388: DSSort(d->eps->ds,d->eigr,d->eigi,rr,ri,&k);
389: /* Put the best pair at the beginning. Useful to check its residual */
390: #if !defined(PETSC_USE_COMPLEX)
391: if (n != 1 && (n != 2 || d->eigi[0] == 0.0))
392: #else
393: if (n != 1)
394: #endif
395: {
396: dvd_calcpairs_apply_arbitrary(d,0,nV,rr,ri);
397: k = 1;
398: DSSort(d->eps->ds,d->eigr,d->eigi,rr,ri,&k);
399: }
400: DSSynchronize(d->eps->ds,d->eigr,d->eigi);
402: if (d->calcpairs_eigs_trans) {
403: d->calcpairs_eigs_trans(d);
404: }
405: if (d->eps->arbitrary || d->calcpairs_eig_backtrans) {
406: PetscFree(rr);
407: PetscFree(ri);
408: }
409: return(0);
410: }
412: static PetscErrorCode EPSXDComputeDSConv(dvdDashboard *d)
413: {
414: PetscErrorCode ierr;
415: PetscInt i,ld;
416: Vec v;
417: PetscScalar *pA;
418: const PetscScalar *pv;
419: PetscBool symm;
422: BVSetActiveColumns(d->eps->V,0,d->eps->nconv);
423: PetscObjectTypeCompare((PetscObject)d->eps->ds,DSHEP,&symm);
424: if (symm) return(0);
425: DSSetDimensions(d->eps->ds,d->eps->nconv,0,0,0);
426: DSCopyMat(d->eps->ds,DS_MAT_A,0,0,d->H,0,0,d->eps->nconv,d->eps->nconv,PETSC_FALSE);
427: if (d->G) {
428: DSCopyMat(d->eps->ds,DS_MAT_B,0,0,d->G,0,0,d->eps->nconv,d->eps->nconv,PETSC_FALSE);
429: }
430: /* Set the signature on projected matrix B */
431: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
432: DSGetLeadingDimension(d->eps->ds,&ld);
433: DSGetArray(d->eps->ds,DS_MAT_B,&pA);
434: PetscArrayzero(pA,d->eps->nconv*ld);
435: VecCreateSeq(PETSC_COMM_SELF,d->eps->nconv,&v);
436: BVGetSignature(d->eps->V,v);
437: VecGetArrayRead(v,&pv);
438: for (i=0;i<d->eps->nconv;i++) pA[i+ld*i] = pv[i];
439: VecRestoreArrayRead(v,&pv);
440: VecDestroy(&v);
441: DSRestoreArray(d->eps->ds,DS_MAT_B,&pA);
442: }
443: DSSetState(d->eps->ds,DS_STATE_RAW);
444: DSSolve(d->eps->ds,d->eps->eigr,d->eps->eigi);
445: DSSynchronize(d->eps->ds,d->eps->eigr,d->eps->eigi);
446: if (d->W) {
447: for (i=0; i<d->eps->nconv; i++) {
448: d->calcpairs_eig_backtrans(d,d->eps->eigr[i],d->eps->eigi[i],&d->eps->eigr[i],&d->eps->eigi[i]);
449: }
450: }
451: return(0);
452: }
454: /*
455: Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
456: the norm associated to the Schur pair, where i = r_s..r_e
457: */
458: static PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d,PetscInt r_s,PetscInt r_e)
459: {
460: PetscInt i,ldpX;
461: PetscScalar *pX;
463: BV BX = d->BX?d->BX:d->eps->V;
464: Vec *R;
467: DSGetLeadingDimension(d->eps->ds,&ldpX);
468: DSGetArray(d->eps->ds,DS_MAT_Q,&pX);
469: /* nX(i) <- ||X(i)|| */
470: dvd_improvex_compute_X(d,r_s,r_e,NULL,pX,ldpX);
471: SlepcVecPoolGetVecs(d->auxV,r_e-r_s,&R);
472: for (i=r_s;i<r_e;i++) {
473: /* R(i-r_s) <- AV*pX(i) */
474: BVMultVec(d->AX,1.0,0.0,R[i-r_s],&pX[ldpX*i]);
475: /* R(i-r_s) <- R(i-r_s) - eigr(i)*BX*pX(i) */
476: BVMultVec(BX,-d->eigr[i],1.0,R[i-r_s],&pX[ldpX*i]);
477: }
478: DSRestoreArray(d->eps->ds,DS_MAT_Q,&pX);
479: d->calcpairs_proj_res(d,r_s,r_e,R);
480: SlepcVecPoolRestoreVecs(d->auxV,r_e-r_s,&R);
481: return(0);
482: }
484: static PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d,PetscInt r_s,PetscInt r_e,Vec *R)
485: {
486: PetscInt i,l,k;
488: PetscBool lindep=PETSC_FALSE;
489: BV cX;
492: if (d->W) cX = d->W; /* If left subspace exists, R <- orth(cY, R), nR[i] <- ||R[i]|| */
493: else if (!(DVD_IS(d->sEP, DVD_EP_STD) && DVD_IS(d->sEP, DVD_EP_HERMITIAN))) cX = d->eps->V; /* If not HEP, R <- orth(cX, R), nR[i] <- ||R[i]|| */
494: else cX = NULL; /* Otherwise, nR[i] <- ||R[i]|| */
496: if (cX) {
497: BVGetActiveColumns(cX,&l,&k);
498: BVSetActiveColumns(cX,0,l);
499: for (i=0;i<r_e-r_s;i++) {
500: BVOrthogonalizeVec(cX,R[i],NULL,&d->nR[r_s+i],&lindep);
501: }
502: BVSetActiveColumns(cX,l,k);
503: if (lindep || (PetscAbs(d->nR[r_s+i]) < PETSC_MACHINE_EPSILON)) {
504: PetscInfo2(d->eps,"The computed eigenvector residual %D is too low, %g!\n",r_s+i,(double)(d->nR[r_s+i]));
505: }
506: } else {
507: for (i=0;i<r_e-r_s;i++) {
508: VecNormBegin(R[i],NORM_2,&d->nR[r_s+i]);
509: }
510: for (i=0;i<r_e-r_s;i++) {
511: VecNormEnd(R[i],NORM_2,&d->nR[r_s+i]);
512: }
513: }
514: return(0);
515: }
517: PetscErrorCode dvd_calcpairs_qz(dvdDashboard *d,dvdBlackboard *b,PetscBool borth,PetscBool harm)
518: {
520: PetscBool std_probl,her_probl,ind_probl;
521: DSType dstype;
522: Vec v1;
525: std_probl = DVD_IS(d->sEP,DVD_EP_STD)? PETSC_TRUE: PETSC_FALSE;
526: her_probl = DVD_IS(d->sEP,DVD_EP_HERMITIAN)? PETSC_TRUE: PETSC_FALSE;
527: ind_probl = DVD_IS(d->sEP,DVD_EP_INDEFINITE)? PETSC_TRUE: PETSC_FALSE;
529: /* Setting configuration constrains */
530: b->max_size_proj = PetscMax(b->max_size_proj,b->max_size_V);
531: d->W_shift = d->B? PETSC_TRUE: PETSC_FALSE;
533: /* Setup the step */
534: if (b->state >= DVD_STATE_CONF) {
535: d->max_size_P = b->max_size_P;
536: d->max_size_proj = b->max_size_proj;
537: /* Create a DS if the method works with Schur decompositions */
538: d->calcPairs = dvd_calcpairs_proj;
539: d->calcpairs_residual = dvd_calcpairs_res_0;
540: d->calcpairs_proj_res = dvd_calcpairs_proj_res;
541: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs;
542: /* Create and configure a DS for solving the projected problems */
543: if (d->W) dstype = DSGNHEP; /* If we use harmonics */
544: else {
545: if (ind_probl) dstype = DSGHIEP;
546: else if (std_probl) dstype = her_probl? DSHEP : DSNHEP;
547: else dstype = her_probl? DSGHEP : DSGNHEP;
548: }
549: DSSetType(d->eps->ds,dstype);
550: DSAllocate(d->eps->ds,d->eps->ncv);
551: /* Create various vector basis */
552: if (harm) {
553: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->W);
554: BVSetMatrix(d->W,NULL,PETSC_FALSE);
555: } else d->W = NULL;
556: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->AX);
557: BVSetMatrix(d->AX,NULL,PETSC_FALSE);
558: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->auxBV);
559: BVSetMatrix(d->auxBV,NULL,PETSC_FALSE);
560: if (d->B) {
561: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->BX);
562: BVSetMatrix(d->BX,NULL,PETSC_FALSE);
563: } else d->BX = NULL;
564: MatCreateVecsEmpty(d->A,&v1,NULL);
565: SlepcVecPoolCreate(v1,0,&d->auxV);
566: VecDestroy(&v1);
567: /* Create projected problem matrices */
568: MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&d->H);
569: if (!std_probl) {
570: MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&d->G);
571: } else d->G = NULL;
572: if (her_probl) {
573: MatSetOption(d->H,MAT_HERMITIAN,PETSC_TRUE);
574: if (d->G) { MatSetOption(d->G,MAT_HERMITIAN,PETSC_TRUE); }
575: }
577: if (ind_probl) {
578: PetscMalloc1(d->eps->ncv,&d->nBds);
579: } else d->nBds = NULL;
580: MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&d->auxM);
582: EPSDavidsonFLAdd(&d->startList,dvd_calcpairs_qz_start);
583: EPSDavidsonFLAdd(&d->endList,EPSXDComputeDSConv);
584: EPSDavidsonFLAdd(&d->destroyList,dvd_calcpairs_qz_d);
585: }
586: return(0);
587: }