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