Actual source code: dvd_calcpairs.c
slepc-3.5.2 2014-10-10
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-2014, 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>
36: static PetscErrorCode dvd_calcpairs_proj(dvdDashboard *d);
37: static PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d);
38: static PetscErrorCode dvd_calcpairs_qz_d(dvdDashboard *d);
39: static PetscErrorCode dvd_calcpairs_projeig_solve(dvdDashboard *d);
40: static PetscErrorCode dvd_calcpairs_selectPairs(dvdDashboard *d,PetscInt n);
41: static PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d,PetscInt r_s,PetscInt r_e);
42: static PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d,PetscInt r_s,PetscInt r_e,Vec *R);
43: static PetscErrorCode dvd_calcpairs_updateproj(dvdDashboard *d);
44: static PetscErrorCode dvd_calcpairs_apply_arbitrary(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscScalar *rr,PetscScalar *ri);
45: static PetscErrorCode EPSXDUpdateProj(Mat Q,Mat Z,PetscInt l,Mat A,PetscInt lA,PetscInt kA,Mat aux);
46: PETSC_STATIC_INLINE PetscErrorCode dvd_calcpairs_updateBV0_gen(dvdDashboard *d,BV bv,DSMatType MT);
47: static PetscErrorCode EPSXDComputeDSConv(dvdDashboard *d);
49: /**** Control routines ********************************************************/
52: PetscErrorCode dvd_calcpairs_qz(dvdDashboard *d,dvdBlackboard *b,PetscBool borth,PetscInt cX_proj,PetscBool harm)
53: {
55: PetscBool std_probl,her_probl,ind_probl,her_ind_probl;
56: DSType dstype;
57: Vec v1;
60: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
61: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
62: ind_probl = DVD_IS(d->sEP, DVD_EP_INDEFINITE)?PETSC_TRUE:PETSC_FALSE;
63: her_ind_probl = (her_probl || ind_probl)? PETSC_TRUE:PETSC_FALSE;
65: /* Setting configuration constrains */
66: b->max_size_proj = PetscMax(b->max_size_proj, b->max_size_V+cX_proj);
67: d->W_shift = d->B?PETSC_TRUE:PETSC_FALSE;
68: if (d->B && her_ind_probl && !borth) d->BV_shift = PETSC_TRUE;
69: else d->BV_shift = PETSC_FALSE;
71: /* Setup the step */
72: if (b->state >= DVD_STATE_CONF) {
73: d->max_cX_in_proj = cX_proj;
74: d->max_size_P = b->max_size_P;
75: d->max_size_proj = b->max_size_proj;
76: /* Create a DS if the method works with Schur decompositions */
77: d->calcPairs = dvd_calcpairs_proj;
78: d->calcpairs_residual = dvd_calcpairs_res_0;
79: d->calcpairs_proj_res = dvd_calcpairs_proj_res;
80: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs;
81: /* Create and configure a DS for solving the projected problems */
82: if (d->W) { /* If we use harmonics */
83: dstype = DSGNHEP;
84: } else {
85: if (ind_probl) {
86: dstype = DSGHIEP;
87: } else if (std_probl) {
88: dstype = her_probl ? DSHEP : DSNHEP;
89: } else {
90: dstype = her_probl ? DSGHEP : DSGNHEP;
91: }
92: }
93: DSSetType(d->eps->ds,dstype);
94: DSAllocate(d->eps->ds,d->eps->ncv);
95: /* Create various vector basis */
96: if (harm) {
97: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->W);
98: BVSetMatrix(d->W,NULL,PETSC_FALSE);
99: } else d->W = NULL;
100: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->AX);
101: BVSetMatrix(d->AX,NULL,PETSC_FALSE);
102: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->auxBV);
103: BVSetMatrix(d->auxBV,NULL,PETSC_FALSE);
104: if (d->B) {
105: BVDuplicateResize(d->eps->V,d->eps->ncv,&d->BX);
106: BVSetMatrix(d->BX,NULL,PETSC_FALSE);
107: } else d->BX = NULL;
108: MatGetVecs(d->A,&v1,NULL);
109: SlepcVecPoolCreate(v1,0,&d->auxV);
110: VecDestroy(&v1);
111: /* Create projected problem matrices */
112: MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&d->H);
113: if (!std_probl) {
114: MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&d->G);
115: } else d->G = NULL;
116: if (her_probl) {
117: MatSetOption(d->H,MAT_HERMITIAN,PETSC_TRUE);
118: if (d->G) {MatSetOption(d->G,MAT_HERMITIAN,PETSC_TRUE);}
119: }
121: if (ind_probl) {
122: PetscMalloc1(d->eps->ncv,&d->nBds);
123: } else d->nBds = NULL;
124: MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&d->auxM);
126: EPSDavidsonFLAdd(&d->startList,dvd_calcpairs_qz_start);
127: EPSDavidsonFLAdd(&d->endList,EPSXDComputeDSConv);
128: EPSDavidsonFLAdd(&d->destroyList,dvd_calcpairs_qz_d);
129: }
130: return(0);
131: }
135: static PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d)
136: {
140: BVSetActiveColumns(d->eps->V,0,0);
141: if (d->W) { BVSetActiveColumns(d->W,0,0); }
142: BVSetActiveColumns(d->AX,0,0);
143: if (d->BX) { BVSetActiveColumns(d->BX,0,0); }
144: return(0);
145: }
149: static PetscErrorCode dvd_calcpairs_qz_d(dvdDashboard *d)
150: {
151: PetscErrorCode ierr;
154: BVDestroy(&d->W);
155: BVDestroy(&d->AX);
156: BVDestroy(&d->BX);
157: BVDestroy(&d->auxBV);
158: MatDestroy(&d->H);
159: if (d->G) {MatDestroy(&d->G);}
160: MatDestroy(&d->auxM);
161: SlepcVecPoolDestroy(&d->auxV);
162: PetscFree(d->nBds);
163: return(0);
164: }
168: static PetscErrorCode dvd_calcpairs_proj(dvdDashboard *d)
169: {
171: PetscInt i,l,k;
172: Vec v1,v2;
173: PetscScalar *pv;
176: BVGetActiveColumns(d->eps->V,&l,&k);
177: /* Update AV, BV, W and the projected matrices */
178: /* 1. S <- S*MT */
179: if (d->V_tra_s != d->V_tra_e) {
180: dvd_calcpairs_updateBV0_gen(d,d->eps->V,DS_MAT_Q);
181: if (d->W) {dvd_calcpairs_updateBV0_gen(d,d->W,DS_MAT_Z);}
182: dvd_calcpairs_updateBV0_gen(d,d->AX,DS_MAT_Q);
183: if (d->BX) {dvd_calcpairs_updateBV0_gen(d,d->BX,DS_MAT_Q);}
184: dvd_calcpairs_updateproj(d);
185: /* Update signature */
186: if (d->nBds) {
187: VecCreateSeq(PETSC_COMM_SELF,l+d->V_tra_e,&v1);
188: BVSetActiveColumns(d->eps->V,0,l+d->V_tra_e);
189: BVGetSignature(d->eps->V,v1);
190: VecGetArray(v1,&pv);
191: for (i=0; i<d->V_tra_e; ++i) pv[l+i] = d->nBds[i];
192: VecRestoreArray(v1,&pv);
193: BVSetSignature(d->eps->V,v1);
194: BVSetActiveColumns(d->eps->V,l,k);
195: VecDestroy(&v1);
196: }
197: k = l+d->V_tra_e;
198: l+= d->V_tra_s;
199: } else {
200: /* 2. V <- orth(V, V_new) */
201: dvd_orthV(d->eps->V,l+d->V_new_s,l+d->V_new_e,d->eps->rand);
202: /* 3. AV <- [AV A * V(V_new_s:V_new_e-1)] */
203: /* Check consistency */
204: if (k-l != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
205: for (i=l+d->V_new_s; i<l+d->V_new_e; i++) {
206: BVGetColumn(d->eps->V,i,&v1);
207: BVGetColumn(d->AX,i,&v2);
208: MatMult(d->A,v1,v2);
209: BVRestoreColumn(d->eps->V,i,&v1);
210: BVRestoreColumn(d->AX,i,&v2);
211: }
212: /* 4. BV <- [BV B * V(V_new_s:V_new_e-1)] */
213: if (d->BX) {
214: /* Check consistency */
215: if (k-l != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
216: for (i=l+d->V_new_s; i<l+d->V_new_e; i++) {
217: BVGetColumn(d->eps->V,i,&v1);
218: BVGetColumn(d->BX,i,&v2);
219: MatMult(d->B,v1,v2);
220: BVRestoreColumn(d->eps->V,i,&v1);
221: BVRestoreColumn(d->BX,i,&v2);
222: }
223: }
224: /* 5. W <- [W f(AV,BV)] */
225: if (d->W) {
226: d->calcpairs_W(d);
227: dvd_orthV(d->W,l+d->V_new_s,l+d->V_new_e,d->eps->rand);
228: }
229: /* 6. H <- W' * AX; G <- W' * BX */
230: BVSetActiveColumns(d->eps->V,l+d->V_new_s,l+d->V_new_e);
231: BVSetActiveColumns(d->AX,l+d->V_new_s,l+d->V_new_e);
232: if (d->BX) {BVSetActiveColumns(d->BX,l+d->V_new_s,l+d->V_new_e);}
233: if (d->W) {BVSetActiveColumns(d->W,l+d->V_new_s,l+d->V_new_e);}
234: BVMatProject(d->AX,NULL,d->W?d->W:d->eps->V,d->H);
235: if (d->G) {BVMatProject(d->BX?d->BX:d->eps->V,NULL,d->W?d->W:d->eps->V,d->G);}
236: BVSetActiveColumns(d->eps->V,l,k);
237: BVSetActiveColumns(d->AX,l,k);
238: if (d->BX) {BVSetActiveColumns(d->BX,l,k);}
239: if (d->W) {BVSetActiveColumns(d->W,l,k);}
241: /* Perform the transformation on the projected problem */
242: if (d->W) {
243: d->calcpairs_proj_trans(d);
244: }
245: k = l+d->V_new_e;
246: }
247: BVSetActiveColumns(d->eps->V,l,k);
248: BVSetActiveColumns(d->AX,l,k);
249: if (d->BX) {BVSetActiveColumns(d->BX,l,k);}
250: if (d->W) {BVSetActiveColumns(d->W,l,k);}
253: /* Solve the projected problem */
254: dvd_calcpairs_projeig_solve(d);
256: d->V_tra_s = d->V_tra_e = 0;
257: d->V_new_s = d->V_new_e;
258: return(0);
259: }
261: /**** Basic routines **********************************************************/
265: static PetscErrorCode dvd_calcpairs_updateproj(dvdDashboard *d)
266: {
267: PetscErrorCode ierr;
268: Mat Q,Z;
269: PetscInt lV,kV;
270: PetscBool symm;
273: DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
274: if (d->W) {DSGetMat(d->eps->ds,DS_MAT_Z,&Z);}
275: else Z = Q;
276: BVGetActiveColumns(d->eps->V,&lV,&kV);
277: EPSXDUpdateProj(Q,Z,0,d->H,lV,lV+d->V_tra_e,d->auxM);
278: if (d->G) {EPSXDUpdateProj(Q,Z,0,d->G,lV,lV+d->V_tra_e,d->auxM);}
279: DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q);
280: if (d->W) {DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z);}
282: PetscObjectTypeCompareAny((PetscObject)d->eps->ds,&symm,DSHEP,DSGHIEP,"");
283: if (d->V_tra_s==0 || symm) return(0);
284: /* Compute upper part of H (and G): H(0:l-1,l:k-1) <- W(0:l-1)' * AV(l:k-1), where
285: k=l+d->V_tra_s */
286: BVSetActiveColumns(d->W?d->W:d->eps->V,0,lV);
287: BVSetActiveColumns(d->AX,lV,lV+d->V_tra_s);
288: BVDot(d->AX,d->W?d->W:d->eps->V,d->H);
289: if (d->G) {
290: BVSetActiveColumns(d->BX?d->BX:d->eps->V,lV,lV+d->V_tra_s);
291: BVDot(d->BX?d->BX:d->eps->V,d->W?d->W:d->eps->V,d->G);
292: }
293: PetscObjectTypeCompareAny((PetscObject)d->eps->ds,&symm,DSGHEP,"");
294: if (!symm) {
295: BVSetActiveColumns(d->W?d->W:d->eps->V,lV,lV+d->V_tra_s);
296: BVSetActiveColumns(d->AX,0,lV);
297: BVDot(d->AX,d->W?d->W:d->eps->V,d->H);
298: if (d->G) {
299: BVSetActiveColumns(d->BX?d->BX:d->eps->V,0,lV);
300: BVDot(d->BX?d->BX:d->eps->V,d->W?d->W:d->eps->V,d->G);
301: }
302: }
303: BVSetActiveColumns(d->eps->V,lV,kV);
304: BVSetActiveColumns(d->AX,lV,kV);
305: if (d->BX) {BVSetActiveColumns(d->BX,lV,kV);}
306: if (d->W) {BVSetActiveColumns(d->W,lV,kV);}
307: if (d->W) {dvd_harm_updateproj(d);}
308: return(0);
309: }
311: /* in complex, d->size_H real auxiliar values are needed */
314: static PetscErrorCode dvd_calcpairs_projeig_solve(dvdDashboard *d)
315: {
316: PetscErrorCode ierr;
317: Mat A;
318: Vec v;
319: PetscScalar *pA,*pv;
320: PetscInt i,lV,kV,n,ld;
323: BVGetActiveColumns(d->eps->V,&lV,&kV);
324: n = kV-lV;
325: DSSetDimensions(d->eps->ds,n,0,0,0);
326: DSGetMat(d->eps->ds,DS_MAT_A,&A);
327: SlepcMatDenseCopy(d->H,lV,lV,A,0,0,n,n);
328: DSRestoreMat(d->eps->ds,DS_MAT_A,&A);
329: if (d->G) {
330: DSGetMat(d->eps->ds,DS_MAT_B,&A);
331: SlepcMatDenseCopy(d->G,lV,lV,A,0,0,n,n);
332: DSRestoreMat(d->eps->ds,DS_MAT_B,&A);
333: }
334: /* Set the signature on projected matrix B */
335: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
336: DSGetLeadingDimension(d->eps->ds,&ld);
337: DSGetArray(d->eps->ds,DS_MAT_B,&pA);
338: PetscMemzero(pA,sizeof(PetscScalar)*n*ld);
339: VecCreateSeq(PETSC_COMM_SELF,kV,&v);
340: BVGetSignature(d->eps->V,v);
341: VecGetArray(v,&pv);
342: for (i=0; i<n; i++) {
343: pA[i+ld*i] = d->nBds[i] = PetscRealPart(pv[lV+i]);
344: }
345: VecRestoreArray(v,&pv);
346: VecDestroy(&v);
347: DSRestoreArray(d->eps->ds,DS_MAT_B,&pA);
348: }
349: DSSetState(d->eps->ds,DS_STATE_RAW);
350: DSSolve(d->eps->ds,d->eigr,d->eigi);
351: return(0);
352: }
356: static PetscErrorCode dvd_calcpairs_apply_arbitrary(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscScalar *rr,PetscScalar *ri)
357: {
358: PetscInt i,k,ld;
359: PetscScalar *pX;
360: Vec *X,xr,xi;
361: PetscErrorCode ierr;
362: #if defined(PETSC_USE_COMPLEX)
363: PetscInt N=1;
364: #else
365: PetscInt N=2,j;
366: #endif
369: /* Quick exit without neither arbitrary selection nor harmonic extraction */
370: if (!d->eps->arbitrary && !d->calcpairs_eig_backtrans) return(0);
372: /* Quick exit without arbitrary selection, but with harmonic extraction */
373: if (d->calcpairs_eig_backtrans) {
374: for (i=r_s; i<r_e; i++) {
375: d->calcpairs_eig_backtrans(d,d->eigr[i],d->eigi[i],&rr[i-r_s],&ri[i-r_s]);
376: }
377: }
378: if (!d->eps->arbitrary) return(0);
380: SlepcVecPoolGetVecs(d->auxV,N,&X);
381: DSGetLeadingDimension(d->eps->ds,&ld);
382: for (i=r_s; i<r_e; i++) {
383: k = i;
384: DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
385: DSNormalize(d->eps->ds,DS_MAT_X,i);
386: DSGetArray(d->eps->ds,DS_MAT_X,&pX);
387: dvd_improvex_compute_X(d,i,k+1,X,pX,ld);
388: DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);
389: #if !defined(PETSC_USE_COMPLEX)
390: if (d->nX[i] != 1.0) {
391: for (j=i; j<k+1; j++) {
392: VecScale(X[j-i],1.0/d->nX[i]);
393: }
394: }
395: xr = X[0];
396: xi = X[1];
397: if (i == k) {
398: VecZeroEntries(xi);
399: }
400: #else
401: xr = X[0];
402: xi = NULL;
403: if (d->nX[i] != 1.0) {
404: VecScale(xr,1.0/d->nX[i]);
405: }
406: #endif
407: (d->eps->arbitrary)(rr[i-r_s],ri[i-r_s],xr,xi,&rr[i-r_s],&ri[i-r_s],d->eps->arbitraryctx);
408: #if !defined(PETSC_USE_COMPLEX)
409: if (i != k) {
410: rr[i+1-r_s] = rr[i-r_s];
411: ri[i+1-r_s] = ri[i-r_s];
412: i++;
413: }
414: #endif
415: }
416: SlepcVecPoolRestoreVecs(d->auxV,N,&X);
417: return(0);
418: }
422: static PetscErrorCode dvd_calcpairs_selectPairs(dvdDashboard *d,PetscInt n)
423: {
424: PetscInt k,lV,kV,nV;
425: PetscScalar *rr,*ri;
426: PetscErrorCode ierr;
429: BVGetActiveColumns(d->eps->V,&lV,&kV);
430: nV = kV - lV;
431: n = PetscMin(n,nV);
432: /* Put the best n pairs at the beginning. Useful for restarting */
433: if (d->eps->arbitrary || d->calcpairs_eig_backtrans) {
434: PetscMalloc1(nV,&rr);
435: PetscMalloc1(nV,&ri);
436: dvd_calcpairs_apply_arbitrary(d,0,nV,rr,ri);
437: } else {
438: rr = d->eigr;
439: ri = d->eigi;
440: }
441: k = n;
442: DSSort(d->eps->ds,d->eigr,d->eigi,rr,ri,&k);
443: /* Put the best pair at the beginning. Useful to check its residual */
444: #if !defined(PETSC_USE_COMPLEX)
445: if (n != 1 && (n != 2 || d->eigi[0] == 0.0))
446: #else
447: if (n != 1)
448: #endif
449: {
450: dvd_calcpairs_apply_arbitrary(d,0,nV,rr,ri);
451: k = 1;
452: DSSort(d->eps->ds,d->eigr,d->eigi,rr,ri,&k);
453: }
454: if (d->calcpairs_eigs_trans) {
455: d->calcpairs_eigs_trans(d);
456: }
457: if (d->eps->arbitrary || d->calcpairs_eig_backtrans) {
458: PetscFree(rr);
459: PetscFree(ri);
460: }
461: return(0);
462: }
466: static PetscErrorCode EPSXDComputeDSConv(dvdDashboard *d)
467: {
469: PetscInt i,ld;
470: Mat A;
471: Vec v;
472: PetscScalar *pA,*pv;
473: PetscBool symm;
476: BVSetActiveColumns(d->eps->V,0,d->eps->nconv);
477: PetscObjectTypeCompareAny((PetscObject)d->eps->ds,&symm,DSHEP,"");
478: if (symm) return(0);
479: DSSetDimensions(d->eps->ds,d->eps->nconv,0,0,0);
480: DSGetMat(d->eps->ds,DS_MAT_A,&A);
481: SlepcMatDenseCopy(d->H,0,0,A,0,0,d->eps->nconv,d->eps->nconv);
482: DSRestoreMat(d->eps->ds,DS_MAT_A,&A);
483: if (d->G) {
484: DSGetMat(d->eps->ds,DS_MAT_B,&A);
485: SlepcMatDenseCopy(d->G,0,0,A,0,0,d->eps->nconv,d->eps->nconv);
486: DSRestoreMat(d->eps->ds,DS_MAT_B,&A);
487: }
488: /* Set the signature on projected matrix B */
489: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
490: DSGetLeadingDimension(d->eps->ds,&ld);
491: DSGetArray(d->eps->ds,DS_MAT_B,&pA);
492: PetscMemzero(pA,sizeof(PetscScalar)*d->eps->nconv*ld);
493: VecCreateSeq(PETSC_COMM_SELF,d->eps->nconv,&v);
494: BVGetSignature(d->eps->V,v);
495: VecGetArray(v,&pv);
496: for (i=0; i<d->eps->nconv; i++) {
497: pA[i+ld*i] = pv[i];
498: }
499: VecRestoreArray(v,&pv);
500: VecDestroy(&v);
501: DSRestoreArray(d->eps->ds,DS_MAT_B,&pA);
502: }
503: DSSetState(d->eps->ds,DS_STATE_RAW);
504: DSSolve(d->eps->ds,d->eps->eigr,d->eps->eigi);
505: if (d->W) {
506: for (i=0; i<d->eps->nconv; i++) {
507: d->calcpairs_eig_backtrans(d,d->eps->eigr[i],d->eps->eigi[i],&d->eps->eigr[i],&d->eps->eigi[i]);
508: }
509: }
510: return(0);
511: }
515: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
516: the norm associated to the Schur pair, where i = r_s..r_e
517: */
518: static PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d,PetscInt r_s,PetscInt r_e)
519: {
520: PetscInt i,ldpX;
521: PetscScalar *pX;
522: PetscErrorCode ierr;
523: BV BX = d->BX?d->BX:d->eps->V;
524: Vec *R;
527: DSGetLeadingDimension(d->eps->ds,&ldpX);
528: DSGetArray(d->eps->ds,DS_MAT_Q,&pX);
529: /* nX(i) <- ||X(i)|| */
530: dvd_improvex_compute_X(d,r_s,r_e,NULL,pX,ldpX);
531: SlepcVecPoolGetVecs(d->auxV,r_e-r_s,&R);
532: for (i=r_s; i<r_e; i++) {
533: /* R(i-r_s) <- AV*pX(i) */
534: BVMultVec(d->AX,1.0,0.0,R[i-r_s],&pX[ldpX*i]);
535: /* R(i-r_s) <- R(i-r_s) - eigr(i)*BX*pX(i) */
536: BVMultVec(BX,-d->eigr[i],1.0,R[i-r_s],&pX[ldpX*i]);
537: }
538: DSRestoreArray(d->eps->ds,DS_MAT_Q,&pX);
539: d->calcpairs_proj_res(d,r_s,r_e,R);
540: SlepcVecPoolRestoreVecs(d->auxV,r_e-r_s,&R);
541: return(0);
542: }
546: static PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d,PetscInt r_s,PetscInt r_e,Vec *R)
547: {
548: PetscInt i,l,k;
549: PetscErrorCode ierr;
550: PetscBool lindep;
551: BV cX;
554: /* If exists left subspace, R <- orth(cY, R), nR[i] <- ||R[i]|| */
555: if (d->W) cX = d->W;
557: /* If not HEP, R <- orth(cX, R), nR[i] <- ||R[i]|| */
558: else if (!(DVD_IS(d->sEP, DVD_EP_STD) && DVD_IS(d->sEP, DVD_EP_HERMITIAN))) cX = d->eps->V;
560: /* Otherwise, nR[i] <- ||R[i]|| */
561: else cX = NULL;
563: if (cX) {
564: BVGetActiveColumns(cX,&l,&k);
565: BVSetActiveColumns(cX,0,l);
566: for (i=0; i<r_e-r_s; i++) {
567: BVOrthogonalizeVec(cX,R[i],NULL,&d->nR[r_s+i],&lindep);
568: }
569: BVSetActiveColumns(cX,l,k);
570: if (lindep || (PetscAbs(d->nR[r_s+i]) < PETSC_MACHINE_EPSILON)) {
571: PetscInfo2(d->eps,"The computed eigenvector residual %D is too low, %g!\n",r_s+i,(double)(d->nR[r_s+i]));
572: }
573: } else {
574: for (i=0;i<r_e-r_s;i++) {
575: VecNormBegin(R[i],NORM_2,&d->nR[r_s+i]);
576: }
577: for (i=0;i<r_e-r_s;i++) {
578: VecNormEnd(R[i],NORM_2,&d->nR[r_s+i]);
579: }
580: }
581: return(0);
582: }
584: /**** Pattern routines ********************************************************/
586: /* BV <- BV*MT */
589: PETSC_STATIC_INLINE PetscErrorCode dvd_calcpairs_updateBV0_gen(dvdDashboard *d,BV bv,DSMatType mat)
590: {
591: PetscErrorCode ierr;
592: PetscInt l,k,n;
593: Mat MT,auxM;
596: BVGetActiveColumns(d->eps->V,&l,&k);
597: DSGetMat(d->eps->ds,mat,&MT);
598: MatCreateSeqDense(PETSC_COMM_SELF,k,k,NULL,&auxM);
599: MatZeroEntries(auxM);
600: MatGetSize(MT,&n,0);
601: if (k-l!=n) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
602: SlepcMatDenseCopy(MT,0,0,auxM,l,l,n,d->V_tra_e);
603: DSRestoreMat(d->eps->ds,mat,&MT);
604: BVMultInPlace(bv,auxM,l,l+d->V_tra_e);
605: MatDestroy(&auxM);
606: return(0);
607: }
611: /* 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 */
612: static PetscErrorCode EPSXDUpdateProj(Mat Q,Mat Z,PetscInt l,Mat A,PetscInt lA,PetscInt kA,Mat aux)
613: {
615: PetscScalar one=1.0,zero=0.0;
616: PetscInt dA_=kA-lA,m0,n0,ldA_,ldQ_,ldZ_,nQ_;
617: PetscBLASInt dA,nQ,ldA,ldQ,ldZ;
618: PetscScalar *pA,*pQ,*pZ,*pW;
619: PetscBool symm=PETSC_FALSE,set,flg;
622: MatGetSize(A,&m0,&n0); ldA_=m0;
623: if (m0!=n0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"A should be square");
624: if (lA<0 || lA>m0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial row, column in A");
625: if (kA<0 || kA<lA || kA>m0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid final row, column in A");
626: MatIsHermitianKnown(A,&set,&flg);
627: symm = set? flg: PETSC_FALSE;
628: MatGetSize(Q,&m0,&n0); ldQ_=nQ_=m0;
629: if (l<0 || l>n0 || l+dA_>n0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial column in Q");
630: MatGetSize(Z,&m0,&n0); ldZ_=m0;
631: if (l<0 || l>n0 || l+dA_>n0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial column in Z");
632: MatGetSize(aux,&m0,&n0);
633: if (m0*n0<nQ_*dA_) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"aux should be larger");
634: PetscBLASIntCast(dA_,&dA);
635: PetscBLASIntCast(nQ_,&nQ);
636: PetscBLASIntCast(ldA_,&ldA);
637: PetscBLASIntCast(ldQ_,&ldQ);
638: PetscBLASIntCast(ldZ_,&ldZ);
639: MatDenseGetArray(A,&pA);
640: MatDenseGetArray(Q,&pQ);
641: if (Q!=Z) {MatDenseGetArray(Z,&pZ);}
642: else pZ = pQ;
643: #if PETSC_USE_DEBUG
644: /* Avoid valgrind warning in xgemm and xsymm */
645: MatZeroEntries(aux);
646: #endif
647: MatDenseGetArray(aux,&pW);
648: /* W = A*Q */
649: if (symm) {
650: PetscStackCallBLAS("BLASsymm",BLASsymm_("L","U",&nQ,&dA,&one,&pA[ldA*lA+lA],&ldA,&pQ[ldQ*l+l],&ldQ,&zero,pW,&nQ));
651: } else {
652: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&nQ,&dA,&nQ,&one,&pA[ldA*lA+lA],&ldA,&pQ[ldQ*l+l],&ldQ,&zero,pW,&nQ));
653: }
654: /* A = Q'*W */
655: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&dA,&dA,&nQ,&one,&pZ[ldZ*l+l],&ldZ,pW,&nQ,&zero,&pA[ldA*lA+lA],&ldA));
656: MatDenseGetArray(A,&pA);
657: MatDenseGetArray(Q,&pQ);
658: if (Q!=Z) {MatDenseGetArray(Z,&pZ);}
659: else pZ = pQ;
660: MatDenseGetArray(aux,&pW);
661: return(0);
662: }