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