Actual source code: dvdcalcpairs.c
slepc-3.7.1 2016-05-27
1: /*
2: SLEPc eigensolver: "davidson"
4: Step: calc the best eigenpairs in the subspace V.
6: For that, performs these steps:
7: 1) Update W <- A * V
8: 2) Update H <- V' * W
9: 3) Obtain eigenpairs of H
10: 4) Select some eigenpairs
11: 5) Compute the Ritz pairs of the selected ones
13: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
14: SLEPc - Scalable Library for Eigenvalue Problem Computations
15: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
17: This file is part of SLEPc.
19: SLEPc is free software: you can redistribute it and/or modify it under the
20: terms of version 3 of the GNU Lesser General Public License as published by
21: the Free Software Foundation.
23: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
24: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
25: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
26: more details.
28: You should have received a copy of the GNU Lesser General Public License
29: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
30: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31: */
33: #include davidson.h
34: #include <slepcblaslapack.h>
38: static PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d)
39: {
43: BVSetActiveColumns(d->eps->V,0,0);
44: if (d->W) { BVSetActiveColumns(d->W,0,0); }
45: BVSetActiveColumns(d->AX,0,0);
46: if (d->BX) { BVSetActiveColumns(d->BX,0,0); }
47: return(0);
48: }
52: static PetscErrorCode dvd_calcpairs_qz_d(dvdDashboard *d)
53: {
54: PetscErrorCode ierr;
57: BVDestroy(&d->W);
58: BVDestroy(&d->AX);
59: BVDestroy(&d->BX);
60: BVDestroy(&d->auxBV);
61: MatDestroy(&d->H);
62: if (d->G) { MatDestroy(&d->G); }
63: MatDestroy(&d->auxM);
64: SlepcVecPoolDestroy(&d->auxV);
65: PetscFree(d->nBds);
66: return(0);
67: }
69: /* in complex, d->size_H real auxiliar values are needed */
72: static PetscErrorCode dvd_calcpairs_projeig_solve(dvdDashboard *d)
73: {
74: PetscErrorCode ierr;
75: Vec v;
76: PetscScalar *pA;
77: const PetscScalar *pv;
78: PetscInt i,lV,kV,n,ld;
81: BVGetActiveColumns(d->eps->V,&lV,&kV);
82: n = kV-lV;
83: DSSetDimensions(d->eps->ds,n,0,0,0);
84: DSCopyMat(d->eps->ds,DS_MAT_A,0,0,d->H,lV,lV,n,n,PETSC_FALSE);
85: if (d->G) {
86: DSCopyMat(d->eps->ds,DS_MAT_B,0,0,d->G,lV,lV,n,n,PETSC_FALSE);
87: }
88: /* Set the signature on projected matrix B */
89: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
90: DSGetLeadingDimension(d->eps->ds,&ld);
91: DSGetArray(d->eps->ds,DS_MAT_B,&pA);
92: PetscMemzero(pA,sizeof(PetscScalar)*n*ld);
93: VecCreateSeq(PETSC_COMM_SELF,kV,&v);
94: BVGetSignature(d->eps->V,v);
95: VecGetArrayRead(v,&pv);
96: for (i=0;i<n;i++) {
97: pA[i+ld*i] = d->nBds[i] = PetscRealPart(pv[lV+i]);
98: }
99: VecRestoreArrayRead(v,&pv);
100: VecDestroy(&v);
101: DSRestoreArray(d->eps->ds,DS_MAT_B,&pA);
102: }
103: DSSetState(d->eps->ds,DS_STATE_RAW);
104: DSSolve(d->eps->ds,d->eigr,d->eigi);
105: return(0);
106: }
110: /*
111: 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
112: */
113: static PetscErrorCode EPSXDUpdateProj(Mat Q,Mat Z,PetscInt l,Mat A,PetscInt lA,PetscInt kA,Mat aux)
114: {
116: PetscScalar one=1.0,zero=0.0;
117: PetscInt i,j,dA_=kA-lA,m0,n0,ldA_,ldQ_,ldZ_,nQ_;
118: PetscBLASInt dA,nQ,ldA,ldQ,ldZ;
119: PetscScalar *pA,*pQ,*pZ,*pW;
120: PetscBool symm=PETSC_FALSE,set,flg;
123: MatGetSize(A,&m0,&n0); ldA_=m0;
124: if (m0!=n0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"A should be square");
125: if (lA<0 || lA>m0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial row, column in A");
126: if (kA<0 || kA<lA || kA>m0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid final row, column in A");
127: MatIsHermitianKnown(A,&set,&flg);
128: symm = set? flg: PETSC_FALSE;
129: MatGetSize(Q,&m0,&n0); ldQ_=nQ_=m0;
130: if (l<0 || l>n0 || l+dA_>n0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial column in Q");
131: MatGetSize(Z,&m0,&n0); ldZ_=m0;
132: if (l<0 || l>n0 || l+dA_>n0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial column in Z");
133: MatGetSize(aux,&m0,&n0);
134: if (m0*n0<nQ_*dA_) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"aux should be larger");
135: PetscBLASIntCast(dA_,&dA);
136: PetscBLASIntCast(nQ_,&nQ);
137: PetscBLASIntCast(ldA_,&ldA);
138: PetscBLASIntCast(ldQ_,&ldQ);
139: PetscBLASIntCast(ldZ_,&ldZ);
140: MatDenseGetArray(A,&pA);
141: MatDenseGetArray(Q,&pQ);
142: if (Q!=Z) { MatDenseGetArray(Z,&pZ); }
143: else pZ = pQ;
144: #if PETSC_USE_DEBUG
145: /* Avoid valgrind warning in xgemm and xsymm */
146: MatZeroEntries(aux);
147: #endif
148: MatDenseGetArray(aux,&pW);
149: /* W = A*Q */
150: if (symm) {
151: /* symmetrize before multiplying */
152: for (i=lA+1;i<lA+nQ;i++) {
153: for (j=lA;j<i;j++) pA[i+j*ldA] = PetscConj(pA[j+i*ldA]);
154: }
155: }
156: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&nQ,&dA,&nQ,&one,&pA[ldA*lA+lA],&ldA,&pQ[ldQ*l+l],&ldQ,&zero,pW,&nQ));
157: /* A = Q'*W */
158: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&dA,&dA,&nQ,&one,&pZ[ldZ*l+l],&ldZ,pW,&nQ,&zero,&pA[ldA*lA+lA],&ldA));
159: MatDenseGetArray(A,&pA);
160: MatDenseGetArray(Q,&pQ);
161: if (Q!=Z) { MatDenseGetArray(Z,&pZ); }
162: else pZ = pQ;
163: MatDenseGetArray(aux,&pW);
164: return(0);
165: }
169: static PetscErrorCode dvd_calcpairs_updateproj(dvdDashboard *d)
170: {
172: Mat Q,Z;
173: PetscInt lV,kV;
174: PetscBool symm;
177: DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
178: if (d->W) { DSGetMat(d->eps->ds,DS_MAT_Z,&Z); }
179: else Z = Q;
180: BVGetActiveColumns(d->eps->V,&lV,&kV);
181: EPSXDUpdateProj(Q,Z,0,d->H,lV,lV+d->V_tra_e,d->auxM);
182: if (d->G) { EPSXDUpdateProj(Q,Z,0,d->G,lV,lV+d->V_tra_e,d->auxM); }
183: DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q);
184: if (d->W) { DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z); }
186: PetscObjectTypeCompareAny((PetscObject)d->eps->ds,&symm,DSHEP,DSGHIEP,"");
187: if (d->V_tra_s==0 || symm) return(0);
188: /* Compute upper part of H (and G): H(0:l-1,l:k-1) <- W(0:l-1)' * AV(l:k-1), where
189: k=l+d->V_tra_s */
190: BVSetActiveColumns(d->W?d->W:d->eps->V,0,lV);
191: BVSetActiveColumns(d->AX,lV,lV+d->V_tra_s);
192: BVDot(d->AX,d->W?d->W:d->eps->V,d->H);
193: if (d->G) {
194: BVSetActiveColumns(d->BX?d->BX:d->eps->V,lV,lV+d->V_tra_s);
195: BVDot(d->BX?d->BX:d->eps->V,d->W?d->W:d->eps->V,d->G);
196: }
197: PetscObjectTypeCompareAny((PetscObject)d->eps->ds,&symm,DSGHEP,"");
198: if (!symm) {
199: BVSetActiveColumns(d->W?d->W:d->eps->V,lV,lV+d->V_tra_s);
200: BVSetActiveColumns(d->AX,0,lV);
201: BVDot(d->AX,d->W?d->W:d->eps->V,d->H);
202: if (d->G) {
203: BVSetActiveColumns(d->BX?d->BX:d->eps->V,0,lV);
204: BVDot(d->BX?d->BX:d->eps->V,d->W?d->W:d->eps->V,d->G);
205: }
206: }
207: BVSetActiveColumns(d->eps->V,lV,kV);
208: BVSetActiveColumns(d->AX,lV,kV);
209: if (d->BX) { BVSetActiveColumns(d->BX,lV,kV); }
210: if (d->W) { BVSetActiveColumns(d->W,lV,kV); }
211: if (d->W) { dvd_harm_updateproj(d); }
212: return(0);
213: }
217: /*
218: BV <- BV*MT
219: */
220: PETSC_STATIC_INLINE PetscErrorCode dvd_calcpairs_updateBV0_gen(dvdDashboard *d,BV bv,DSMatType mat)
221: {
223: PetscInt l,k,n;
224: Mat auxM;
227: BVGetActiveColumns(d->eps->V,&l,&k);
228: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&auxM);
229: MatZeroEntries(auxM);
230: DSGetDimensions(d->eps->ds,&n,NULL,NULL,NULL,NULL);
231: if (k-l!=n) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
232: DSCopyMat(d->eps->ds,mat,0,0,auxM,l,l,n,d->V_tra_e,PETSC_TRUE);
233: BVMultInPlace(bv,auxM,l,l+d->V_tra_e);
234: MatDestroy(&auxM);
235: return(0);
236: }
240: static PetscErrorCode dvd_calcpairs_proj(dvdDashboard *d)
241: {
243: PetscInt i,l,k;
244: Vec v1,v2;
245: PetscScalar *pv;
248: BVGetActiveColumns(d->eps->V,&l,&k);
249: /* Update AV, BV, W and the projected matrices */
250: /* 1. S <- S*MT */
251: if (d->V_tra_s != d->V_tra_e) {
252: dvd_calcpairs_updateBV0_gen(d,d->eps->V,DS_MAT_Q);
253: if (d->W) { dvd_calcpairs_updateBV0_gen(d,d->W,DS_MAT_Z); }
254: dvd_calcpairs_updateBV0_gen(d,d->AX,DS_MAT_Q);
255: if (d->BX) { dvd_calcpairs_updateBV0_gen(d,d->BX,DS_MAT_Q); }
256: dvd_calcpairs_updateproj(d);
257: /* Update signature */
258: if (d->nBds) {
259: VecCreateSeq(PETSC_COMM_SELF,l+d->V_tra_e,&v1);
260: BVSetActiveColumns(d->eps->V,0,l+d->V_tra_e);
261: BVGetSignature(d->eps->V,v1);
262: VecGetArray(v1,&pv);
263: for (i=0;i<d->V_tra_e;i++) pv[l+i] = d->nBds[i];
264: VecRestoreArray(v1,&pv);
265: BVSetSignature(d->eps->V,v1);
266: BVSetActiveColumns(d->eps->V,l,k);
267: VecDestroy(&v1);
268: }
269: k = l+d->V_tra_e;
270: l+= d->V_tra_s;
271: } else {
272: /* 2. V <- orth(V, V_new) */
273: dvd_orthV(d->eps->V,l+d->V_new_s,l+d->V_new_e);
274: /* 3. AV <- [AV A * V(V_new_s:V_new_e-1)] */
275: /* Check consistency */
276: if (k-l != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
277: for (i=l+d->V_new_s;i<l+d->V_new_e;i++) {
278: BVGetColumn(d->eps->V,i,&v1);
279: BVGetColumn(d->AX,i,&v2);
280: MatMult(d->A,v1,v2);
281: BVRestoreColumn(d->eps->V,i,&v1);
282: BVRestoreColumn(d->AX,i,&v2);
283: }
284: /* 4. BV <- [BV B * V(V_new_s:V_new_e-1)] */
285: if (d->BX) {
286: /* Check consistency */
287: if (k-l != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
288: for (i=l+d->V_new_s;i<l+d->V_new_e;i++) {
289: BVGetColumn(d->eps->V,i,&v1);
290: BVGetColumn(d->BX,i,&v2);
291: MatMult(d->B,v1,v2);
292: BVRestoreColumn(d->eps->V,i,&v1);
293: BVRestoreColumn(d->BX,i,&v2);
294: }
295: }
296: /* 5. W <- [W f(AV,BV)] */
297: if (d->W) {
298: d->calcpairs_W(d);
299: dvd_orthV(d->W,l+d->V_new_s,l+d->V_new_e);
300: }
301: /* 6. H <- W' * AX; G <- W' * BX */
302: BVSetActiveColumns(d->eps->V,l+d->V_new_s,l+d->V_new_e);
303: BVSetActiveColumns(d->AX,l+d->V_new_s,l+d->V_new_e);
304: if (d->BX) { BVSetActiveColumns(d->BX,l+d->V_new_s,l+d->V_new_e); }
305: if (d->W) { BVSetActiveColumns(d->W,l+d->V_new_s,l+d->V_new_e); }
306: BVMatProject(d->AX,NULL,d->W?d->W:d->eps->V,d->H);
307: if (d->G) { BVMatProject(d->BX?d->BX:d->eps->V,NULL,d->W?d->W:d->eps->V,d->G); }
308: BVSetActiveColumns(d->eps->V,l,k);
309: BVSetActiveColumns(d->AX,l,k);
310: if (d->BX) { BVSetActiveColumns(d->BX,l,k); }
311: if (d->W) { BVSetActiveColumns(d->W,l,k); }
313: /* Perform the transformation on the projected problem */
314: if (d->W) {
315: d->calcpairs_proj_trans(d);
316: }
317: k = l+d->V_new_e;
318: }
319: BVSetActiveColumns(d->eps->V,l,k);
320: BVSetActiveColumns(d->AX,l,k);
321: if (d->BX) { BVSetActiveColumns(d->BX,l,k); }
322: if (d->W) { BVSetActiveColumns(d->W,l,k); }
324: /* Solve the projected problem */
325: dvd_calcpairs_projeig_solve(d);
327: d->V_tra_s = d->V_tra_e = 0;
328: d->V_new_s = d->V_new_e;
329: return(0);
330: }
334: static PetscErrorCode dvd_calcpairs_apply_arbitrary(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscScalar *rr,PetscScalar *ri)
335: {
336: PetscInt i,k,ld;
337: PetscScalar *pX;
338: Vec *X,xr,xi;
340: #if defined(PETSC_USE_COMPLEX)
341: PetscInt N=1;
342: #else
343: PetscInt N=2,j;
344: #endif
347: /* Quick exit without neither arbitrary selection nor harmonic extraction */
348: if (!d->eps->arbitrary && !d->calcpairs_eig_backtrans) return(0);
350: /* Quick exit without arbitrary selection, but with harmonic extraction */
351: if (d->calcpairs_eig_backtrans) {
352: for (i=r_s; i<r_e; i++) {
353: d->calcpairs_eig_backtrans(d,d->eigr[i],d->eigi[i],&rr[i-r_s],&ri[i-r_s]);
354: }
355: }
356: if (!d->eps->arbitrary) return(0);
358: SlepcVecPoolGetVecs(d->auxV,N,&X);
359: DSGetLeadingDimension(d->eps->ds,&ld);
360: for (i=r_s;i<r_e;i++) {
361: k = i;
362: DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
363: DSGetArray(d->eps->ds,DS_MAT_X,&pX);
364: dvd_improvex_compute_X(d,i,k+1,X,pX,ld);
365: DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);
366: #if !defined(PETSC_USE_COMPLEX)
367: if (d->nX[i] != 1.0) {
368: for (j=i;j<k+1;j++) {
369: VecScale(X[j-i],1.0/d->nX[i]);
370: }
371: }
372: xr = X[0];
373: xi = X[1];
374: if (i == k) {
375: VecSet(xi,0.0);
376: }
377: #else
378: xr = X[0];
379: xi = NULL;
380: if (d->nX[i] != 1.0) {
381: VecScale(xr,1.0/d->nX[i]);
382: }
383: #endif
384: (d->eps->arbitrary)(rr[i-r_s],ri[i-r_s],xr,xi,&rr[i-r_s],&ri[i-r_s],d->eps->arbitraryctx);
385: #if !defined(PETSC_USE_COMPLEX)
386: if (i != k) {
387: rr[i+1-r_s] = rr[i-r_s];
388: ri[i+1-r_s] = ri[i-r_s];
389: i++;
390: }
391: #endif
392: }
393: SlepcVecPoolRestoreVecs(d->auxV,N,&X);
394: return(0);
395: }
399: static PetscErrorCode dvd_calcpairs_selectPairs(dvdDashboard *d,PetscInt n)
400: {
401: PetscInt k,lV,kV,nV;
402: PetscScalar *rr,*ri;
406: BVGetActiveColumns(d->eps->V,&lV,&kV);
407: nV = kV - lV;
408: n = PetscMin(n,nV);
409: /* Put the best n pairs at the beginning. Useful for restarting */
410: if (d->eps->arbitrary || d->calcpairs_eig_backtrans) {
411: PetscMalloc1(nV,&rr);
412: PetscMalloc1(nV,&ri);
413: dvd_calcpairs_apply_arbitrary(d,0,nV,rr,ri);
414: } else {
415: rr = d->eigr;
416: ri = d->eigi;
417: }
418: k = n;
419: DSSort(d->eps->ds,d->eigr,d->eigi,rr,ri,&k);
420: /* Put the best pair at the beginning. Useful to check its residual */
421: #if !defined(PETSC_USE_COMPLEX)
422: if (n != 1 && (n != 2 || d->eigi[0] == 0.0))
423: #else
424: if (n != 1)
425: #endif
426: {
427: dvd_calcpairs_apply_arbitrary(d,0,nV,rr,ri);
428: k = 1;
429: DSSort(d->eps->ds,d->eigr,d->eigi,rr,ri,&k);
430: }
431: if (d->calcpairs_eigs_trans) {
432: d->calcpairs_eigs_trans(d);
433: }
434: if (d->eps->arbitrary || d->calcpairs_eig_backtrans) {
435: PetscFree(rr);
436: PetscFree(ri);
437: }
438: return(0);
439: }
443: static PetscErrorCode EPSXDComputeDSConv(dvdDashboard *d)
444: {
445: PetscErrorCode ierr;
446: PetscInt i,ld;
447: Vec v;
448: PetscScalar *pA;
449: const PetscScalar *pv;
450: PetscBool symm;
453: BVSetActiveColumns(d->eps->V,0,d->eps->nconv);
454: PetscObjectTypeCompareAny((PetscObject)d->eps->ds,&symm,DSHEP,"");
455: if (symm) return(0);
456: DSSetDimensions(d->eps->ds,d->eps->nconv,0,0,0);
457: DSCopyMat(d->eps->ds,DS_MAT_A,0,0,d->H,0,0,d->eps->nconv,d->eps->nconv,PETSC_FALSE);
458: if (d->G) {
459: DSCopyMat(d->eps->ds,DS_MAT_B,0,0,d->G,0,0,d->eps->nconv,d->eps->nconv,PETSC_FALSE);
460: }
461: /* Set the signature on projected matrix B */
462: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
463: DSGetLeadingDimension(d->eps->ds,&ld);
464: DSGetArray(d->eps->ds,DS_MAT_B,&pA);
465: PetscMemzero(pA,sizeof(PetscScalar)*d->eps->nconv*ld);
466: VecCreateSeq(PETSC_COMM_SELF,d->eps->nconv,&v);
467: BVGetSignature(d->eps->V,v);
468: VecGetArrayRead(v,&pv);
469: for (i=0;i<d->eps->nconv;i++) pA[i+ld*i] = pv[i];
470: VecRestoreArrayRead(v,&pv);
471: VecDestroy(&v);
472: DSRestoreArray(d->eps->ds,DS_MAT_B,&pA);
473: }
474: DSSetState(d->eps->ds,DS_STATE_RAW);
475: DSSolve(d->eps->ds,d->eps->eigr,d->eps->eigi);
476: if (d->W) {
477: for (i=0; i<d->eps->nconv; i++) {
478: d->calcpairs_eig_backtrans(d,d->eps->eigr[i],d->eps->eigi[i],&d->eps->eigr[i],&d->eps->eigi[i]);
479: }
480: }
481: return(0);
482: }
486: /*
487: Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
488: the norm associated to the Schur pair, where i = r_s..r_e
489: */
490: static PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d,PetscInt r_s,PetscInt r_e)
491: {
492: PetscInt i,ldpX;
493: PetscScalar *pX;
495: BV BX = d->BX?d->BX:d->eps->V;
496: Vec *R;
499: DSGetLeadingDimension(d->eps->ds,&ldpX);
500: DSGetArray(d->eps->ds,DS_MAT_Q,&pX);
501: /* nX(i) <- ||X(i)|| */
502: dvd_improvex_compute_X(d,r_s,r_e,NULL,pX,ldpX);
503: SlepcVecPoolGetVecs(d->auxV,r_e-r_s,&R);
504: for (i=r_s;i<r_e;i++) {
505: /* R(i-r_s) <- AV*pX(i) */
506: BVMultVec(d->AX,1.0,0.0,R[i-r_s],&pX[ldpX*i]);
507: /* R(i-r_s) <- R(i-r_s) - eigr(i)*BX*pX(i) */
508: BVMultVec(BX,-d->eigr[i],1.0,R[i-r_s],&pX[ldpX*i]);
509: }
510: DSRestoreArray(d->eps->ds,DS_MAT_Q,&pX);
511: d->calcpairs_proj_res(d,r_s,r_e,R);
512: SlepcVecPoolRestoreVecs(d->auxV,r_e-r_s,&R);
513: return(0);
514: }
518: static PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d,PetscInt r_s,PetscInt r_e,Vec *R)
519: {
520: PetscInt i,l,k;
522: PetscBool lindep=PETSC_FALSE;
523: BV cX;
526: if (d->W) cX = d->W; /* If left subspace exists, R <- orth(cY, R), nR[i] <- ||R[i]|| */
527: 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]|| */
528: else cX = NULL; /* Otherwise, nR[i] <- ||R[i]|| */
530: if (cX) {
531: BVGetActiveColumns(cX,&l,&k);
532: BVSetActiveColumns(cX,0,l);
533: for (i=0;i<r_e-r_s;i++) {
534: BVOrthogonalizeVec(cX,R[i],NULL,&d->nR[r_s+i],&lindep);
535: }
536: BVSetActiveColumns(cX,l,k);
537: if (lindep || (PetscAbs(d->nR[r_s+i]) < PETSC_MACHINE_EPSILON)) {
538: PetscInfo2(d->eps,"The computed eigenvector residual %D is too low, %g!\n",r_s+i,(double)(d->nR[r_s+i]));
539: }
540: } else {
541: for (i=0;i<r_e-r_s;i++) {
542: VecNormBegin(R[i],NORM_2,&d->nR[r_s+i]);
543: }
544: for (i=0;i<r_e-r_s;i++) {
545: VecNormEnd(R[i],NORM_2,&d->nR[r_s+i]);
546: }
547: }
548: return(0);
549: }
553: PetscErrorCode dvd_calcpairs_qz(dvdDashboard *d,dvdBlackboard *b,PetscBool borth,PetscInt cX_proj,PetscBool harm)
554: {
556: PetscBool std_probl,her_probl,ind_probl,her_ind_probl;
557: DSType dstype;
558: Vec v1;
561: std_probl = DVD_IS(d->sEP,DVD_EP_STD)? PETSC_TRUE: PETSC_FALSE;
562: her_probl = DVD_IS(d->sEP,DVD_EP_HERMITIAN)? PETSC_TRUE: PETSC_FALSE;
563: ind_probl = DVD_IS(d->sEP,DVD_EP_INDEFINITE)? PETSC_TRUE: PETSC_FALSE;
564: her_ind_probl = (her_probl||ind_probl)? PETSC_TRUE: PETSC_FALSE;
566: /* Setting configuration constrains */
567: b->max_size_proj = PetscMax(b->max_size_proj,b->max_size_V+cX_proj);
568: d->W_shift = d->B? PETSC_TRUE: PETSC_FALSE;
569: if (d->B && her_ind_probl && !borth) d->BV_shift = PETSC_TRUE;
570: else d->BV_shift = PETSC_FALSE;
572: /* Setup the step */
573: if (b->state >= DVD_STATE_CONF) {
574: d->max_cX_in_proj = cX_proj;
575: d->max_size_P = b->max_size_P;
576: d->max_size_proj = b->max_size_proj;
577: /* Create a DS if the method works with Schur decompositions */
578: d->calcPairs = dvd_calcpairs_proj;
579: d->calcpairs_residual = dvd_calcpairs_res_0;
580: d->calcpairs_proj_res = dvd_calcpairs_proj_res;
581: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs;
582: /* Create and configure a DS for solving the projected problems */
583: if (d->W) dstype = DSGNHEP; /* If we use harmonics */
584: else {
585: if (ind_probl) dstype = DSGHIEP;
586: else if (std_probl) dstype = her_probl? DSHEP : DSNHEP;
587: else dstype = her_probl? DSGHEP : DSGNHEP;
588: }
589: DSSetType(d->eps->ds,dstype);
590: DSAllocate(d->eps->ds,d->eps->ncv);
591: /* Create various vector basis */
592: if (harm) {
593: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->W);
594: BVSetMatrix(d->W,NULL,PETSC_FALSE);
595: } else d->W = NULL;
596: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->AX);
597: BVSetMatrix(d->AX,NULL,PETSC_FALSE);
598: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->auxBV);
599: BVSetMatrix(d->auxBV,NULL,PETSC_FALSE);
600: if (d->B) {
601: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->BX);
602: BVSetMatrix(d->BX,NULL,PETSC_FALSE);
603: } else d->BX = NULL;
604: MatCreateVecs(d->A,&v1,NULL);
605: SlepcVecPoolCreate(v1,0,&d->auxV);
606: VecDestroy(&v1);
607: /* Create projected problem matrices */
608: MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&d->H);
609: if (!std_probl) {
610: MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&d->G);
611: } else d->G = NULL;
612: if (her_probl) {
613: MatSetOption(d->H,MAT_HERMITIAN,PETSC_TRUE);
614: if (d->G) { MatSetOption(d->G,MAT_HERMITIAN,PETSC_TRUE); }
615: }
617: if (ind_probl) {
618: PetscMalloc1(d->eps->ncv,&d->nBds);
619: } else d->nBds = NULL;
620: MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&d->auxM);
622: EPSDavidsonFLAdd(&d->startList,dvd_calcpairs_qz_start);
623: EPSDavidsonFLAdd(&d->endList,EPSXDComputeDSConv);
624: EPSDavidsonFLAdd(&d->destroyList,dvd_calcpairs_qz_d);
625: }
626: return(0);
627: }