Actual source code: dsnhep.c
slepc-3.17.0 2022-03-31
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: */
11: #include <slepc/private/dsimpl.h>
12: #include <slepcblaslapack.h>
14: PetscErrorCode DSAllocate_NHEP(DS ds,PetscInt ld)
15: {
16: DSAllocateMat_Private(ds,DS_MAT_A);
17: DSAllocateMat_Private(ds,DS_MAT_Q);
18: PetscFree(ds->perm);
19: PetscMalloc1(ld,&ds->perm);
20: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
21: PetscFunctionReturn(0);
22: }
24: PetscErrorCode DSView_NHEP(DS ds,PetscViewer viewer)
25: {
26: PetscViewerFormat format;
28: PetscViewerGetFormat(viewer,&format);
29: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
30: DSViewMat(ds,viewer,DS_MAT_A);
31: if (ds->state>DS_STATE_INTERMEDIATE) DSViewMat(ds,viewer,DS_MAT_Q);
32: if (ds->mat[DS_MAT_X]) DSViewMat(ds,viewer,DS_MAT_X);
33: if (ds->mat[DS_MAT_Y]) DSViewMat(ds,viewer,DS_MAT_Y);
34: PetscFunctionReturn(0);
35: }
37: static PetscErrorCode DSVectors_NHEP_Refined_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
38: {
39: PetscInt i,j;
40: PetscBLASInt info,ld,n,n1,lwork,inc=1;
41: PetscScalar sdummy,done=1.0,zero=0.0;
42: PetscReal *sigma;
43: PetscBool iscomplex = PETSC_FALSE;
44: PetscScalar *A = ds->mat[DS_MAT_A];
45: PetscScalar *Q = ds->mat[DS_MAT_Q];
46: PetscScalar *X = ds->mat[left?DS_MAT_Y:DS_MAT_X];
47: PetscScalar *W;
50: PetscBLASIntCast(ds->n,&n);
51: PetscBLASIntCast(ds->ld,&ld);
52: n1 = n+1;
53: if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
55: DSAllocateWork_Private(ds,5*ld,6*ld,0);
56: DSAllocateMat_Private(ds,DS_MAT_W);
57: W = ds->mat[DS_MAT_W];
58: lwork = 5*ld;
59: sigma = ds->rwork+5*ld;
61: /* build A-w*I in W */
62: for (j=0;j<n;j++)
63: for (i=0;i<=n;i++)
64: W[i+j*ld] = A[i+j*ld];
65: for (i=0;i<n;i++)
66: W[i+i*ld] -= A[(*k)+(*k)*ld];
68: /* compute SVD of W */
69: #if !defined(PETSC_USE_COMPLEX)
70: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&n1,&n,W,&ld,sigma,&sdummy,&ld,&sdummy,&ld,ds->work,&lwork,&info));
71: #else
72: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&n1,&n,W,&ld,sigma,&sdummy,&ld,&sdummy,&ld,ds->work,&lwork,ds->rwork,&info));
73: #endif
74: SlepcCheckLapackInfo("gesvd",info);
76: /* the smallest singular value is the new error estimate */
77: if (rnorm) *rnorm = sigma[n-1];
79: /* update vector with right singular vector associated to smallest singular value,
80: accumulating the transformation matrix Q */
81: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&done,Q,&ld,W+n-1,&ld,&zero,X+(*k)*ld,&inc));
82: PetscFunctionReturn(0);
83: }
85: static PetscErrorCode DSVectors_NHEP_Refined_All(DS ds,PetscBool left)
86: {
87: PetscInt i;
89: for (i=0;i<ds->n;i++) DSVectors_NHEP_Refined_Some(ds,&i,NULL,left);
90: PetscFunctionReturn(0);
91: }
93: static PetscErrorCode DSVectors_NHEP_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
94: {
95: PetscInt i;
96: PetscBLASInt mm=1,mout,info,ld,n,*select,inc=1,cols=1,zero=0;
97: PetscScalar sone=1.0,szero=0.0;
98: PetscReal norm,done=1.0;
99: PetscBool iscomplex = PETSC_FALSE;
100: PetscScalar *A = ds->mat[DS_MAT_A];
101: PetscScalar *Q = ds->mat[DS_MAT_Q];
102: PetscScalar *X = ds->mat[left?DS_MAT_Y:DS_MAT_X];
103: PetscScalar *Y;
105: PetscBLASIntCast(ds->n,&n);
106: PetscBLASIntCast(ds->ld,&ld);
107: DSAllocateWork_Private(ds,0,0,ld);
108: select = ds->iwork;
109: for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;
111: /* compute k-th eigenvector Y of A */
112: Y = X+(*k)*ld;
113: select[*k] = (PetscBLASInt)PETSC_TRUE;
114: #if !defined(PETSC_USE_COMPLEX)
115: if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
116: mm = iscomplex? 2: 1;
117: if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
118: DSAllocateWork_Private(ds,3*ld,0,0);
119: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_(left?"L":"R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,&info));
120: #else
121: DSAllocateWork_Private(ds,2*ld,ld,0);
122: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_(left?"L":"R","S",select,&n,A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,ds->rwork,&info));
123: #endif
124: SlepcCheckLapackInfo("trevc",info);
127: /* accumulate and normalize eigenvectors */
128: if (ds->state>=DS_STATE_CONDENSED) {
129: PetscArraycpy(ds->work,Y,mout*ld);
130: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work,&inc,&szero,Y,&inc));
131: #if !defined(PETSC_USE_COMPLEX)
132: if (iscomplex) PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work+ld,&inc,&szero,Y+ld,&inc));
133: #endif
134: cols = 1;
135: norm = BLASnrm2_(&n,Y,&inc);
136: #if !defined(PETSC_USE_COMPLEX)
137: if (iscomplex) {
138: norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Y+ld,&inc));
139: cols = 2;
140: }
141: #endif
142: PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Y,&ld,&info));
143: SlepcCheckLapackInfo("lascl",info);
144: }
146: /* set output arguments */
147: if (iscomplex) (*k)++;
148: if (rnorm) {
149: if (iscomplex) *rnorm = SlepcAbsEigenvalue(Y[n-1],Y[n-1+ld]);
150: else *rnorm = PetscAbsScalar(Y[n-1]);
151: }
152: PetscFunctionReturn(0);
153: }
155: static PetscErrorCode DSVectors_NHEP_Eigen_All(DS ds,PetscBool left)
156: {
157: PetscInt i;
158: PetscBLASInt n,ld,mout,info,inc=1,cols,zero=0;
159: PetscBool iscomplex;
160: PetscScalar *X,*Y,*Z,*A = ds->mat[DS_MAT_A];
161: PetscReal norm,done=1.0;
162: const char *side,*back;
164: PetscBLASIntCast(ds->n,&n);
165: PetscBLASIntCast(ds->ld,&ld);
166: if (left) {
167: X = NULL;
168: Y = ds->mat[DS_MAT_Y];
169: side = "L";
170: } else {
171: X = ds->mat[DS_MAT_X];
172: Y = NULL;
173: side = "R";
174: }
175: Z = left? Y: X;
176: if (ds->state>=DS_STATE_CONDENSED) {
177: /* DSSolve() has been called, backtransform with matrix Q */
178: back = "B";
179: PetscArraycpy(Z,ds->mat[DS_MAT_Q],ld*ld);
180: } else back = "A";
181: #if !defined(PETSC_USE_COMPLEX)
182: DSAllocateWork_Private(ds,3*ld,0,0);
183: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_(side,back,NULL,&n,A,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,&info));
184: #else
185: DSAllocateWork_Private(ds,2*ld,ld,0);
186: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_(side,back,NULL,&n,A,&ld,Y,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
187: #endif
188: SlepcCheckLapackInfo("trevc",info);
190: /* normalize eigenvectors */
191: for (i=0;i<n;i++) {
192: iscomplex = (i<n-1 && A[i+1+i*ld]!=0.0)? PETSC_TRUE: PETSC_FALSE;
193: cols = 1;
194: norm = BLASnrm2_(&n,Z+i*ld,&inc);
195: #if !defined(PETSC_USE_COMPLEX)
196: if (iscomplex) {
197: norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Z+(i+1)*ld,&inc));
198: cols = 2;
199: }
200: #endif
201: PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Z+i*ld,&ld,&info));
202: SlepcCheckLapackInfo("lascl",info);
203: if (iscomplex) i++;
204: }
205: PetscFunctionReturn(0);
206: }
208: PetscErrorCode DSVectors_NHEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
209: {
210: switch (mat) {
211: case DS_MAT_X:
212: if (ds->refined) {
214: if (j) DSVectors_NHEP_Refined_Some(ds,j,rnorm,PETSC_FALSE);
215: else DSVectors_NHEP_Refined_All(ds,PETSC_FALSE);
216: } else {
217: if (j) DSVectors_NHEP_Eigen_Some(ds,j,rnorm,PETSC_FALSE);
218: else DSVectors_NHEP_Eigen_All(ds,PETSC_FALSE);
219: }
220: break;
221: case DS_MAT_Y:
223: if (j) DSVectors_NHEP_Eigen_Some(ds,j,rnorm,PETSC_TRUE);
224: else DSVectors_NHEP_Eigen_All(ds,PETSC_TRUE);
225: break;
226: case DS_MAT_U:
227: case DS_MAT_V:
228: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
229: default:
230: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
231: }
232: PetscFunctionReturn(0);
233: }
235: static PetscErrorCode DSSort_NHEP_Arbitrary(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
236: {
237: PetscInt i;
238: PetscBLASInt info,n,ld,mout,lwork,*selection;
239: PetscScalar *T = ds->mat[DS_MAT_A],*Q = ds->mat[DS_MAT_Q],*work;
240: PetscReal dummy;
241: #if !defined(PETSC_USE_COMPLEX)
242: PetscBLASInt *iwork,liwork;
243: #endif
246: PetscBLASIntCast(ds->n,&n);
247: PetscBLASIntCast(ds->ld,&ld);
248: #if !defined(PETSC_USE_COMPLEX)
249: lwork = n;
250: liwork = 1;
251: DSAllocateWork_Private(ds,lwork,0,liwork+n);
252: work = ds->work;
253: lwork = ds->lwork;
254: selection = ds->iwork;
255: iwork = ds->iwork + n;
256: liwork = ds->liwork - n;
257: #else
258: lwork = 1;
259: DSAllocateWork_Private(ds,lwork,0,n);
260: work = ds->work;
261: selection = ds->iwork;
262: #endif
263: /* Compute the selected eigenvalue to be in the leading position */
264: DSSortEigenvalues_Private(ds,rr,ri,ds->perm,PETSC_FALSE);
265: PetscArrayzero(selection,n);
266: for (i=0;i<*k;i++) selection[ds->perm[i]] = 1;
267: #if !defined(PETSC_USE_COMPLEX)
268: PetscStackCallBLAS("LAPACKtrsen",LAPACKtrsen_("N","V",selection,&n,T,&ld,Q,&ld,wr,wi,&mout,&dummy,&dummy,work,&lwork,iwork,&liwork,&info));
269: #else
270: PetscStackCallBLAS("LAPACKtrsen",LAPACKtrsen_("N","V",selection,&n,T,&ld,Q,&ld,wr,&mout,&dummy,&dummy,work,&lwork,&info));
271: #endif
272: SlepcCheckLapackInfo("trsen",info);
273: *k = mout;
274: PetscFunctionReturn(0);
275: }
277: PetscErrorCode DSSort_NHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
278: {
279: if (!rr || wr == rr) DSSort_NHEP_Total(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
280: else DSSort_NHEP_Arbitrary(ds,wr,wi,rr,ri,k);
281: PetscFunctionReturn(0);
282: }
284: static PetscErrorCode DSSortWithPermutation_NHEP(DS ds,PetscInt *perm,PetscScalar *wr,PetscScalar *wi)
285: {
286: DSSortWithPermutation_NHEP_Private(ds,perm,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
287: PetscFunctionReturn(0);
288: }
290: PetscErrorCode DSUpdateExtraRow_NHEP(DS ds)
291: {
292: PetscInt i;
293: PetscBLASInt n,ld,incx=1;
294: PetscScalar *A,*Q,*x,*y,one=1.0,zero=0.0;
296: PetscBLASIntCast(ds->n,&n);
297: PetscBLASIntCast(ds->ld,&ld);
298: A = ds->mat[DS_MAT_A];
299: Q = ds->mat[DS_MAT_Q];
300: DSAllocateWork_Private(ds,2*ld,0,0);
301: x = ds->work;
302: y = ds->work+ld;
303: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
304: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
305: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
306: ds->k = n;
307: PetscFunctionReturn(0);
308: }
310: PetscErrorCode DSSolve_NHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
311: {
312: #if !defined(PETSC_USE_COMPLEX)
314: #endif
315: DSSolve_NHEP_Private(ds,ds->mat[DS_MAT_A],ds->mat[DS_MAT_Q],wr,wi);
316: PetscFunctionReturn(0);
317: }
319: PetscErrorCode DSSynchronize_NHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
320: {
321: PetscInt ld=ds->ld,l=ds->l,k;
322: PetscMPIInt n,rank,off=0,size,ldn;
324: k = (ds->n-l)*ld;
325: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
326: if (eigr) k += ds->n-l;
327: if (eigi) k += ds->n-l;
328: DSAllocateWork_Private(ds,k,0,0);
329: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
330: PetscMPIIntCast(ds->n-l,&n);
331: PetscMPIIntCast(ld*(ds->n-l),&ldn);
332: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
333: if (!rank) {
334: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
335: if (ds->state>DS_STATE_RAW) MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
336: if (eigr) MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
337: #if !defined(PETSC_USE_COMPLEX)
338: if (eigi) MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
339: #endif
340: }
341: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
342: if (rank) {
343: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
344: if (ds->state>DS_STATE_RAW) MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
345: if (eigr) MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
346: #if !defined(PETSC_USE_COMPLEX)
347: if (eigi) MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
348: #endif
349: }
350: PetscFunctionReturn(0);
351: }
353: PetscErrorCode DSTruncate_NHEP(DS ds,PetscInt n,PetscBool trim)
354: {
355: PetscInt i,ld=ds->ld,l=ds->l;
356: PetscScalar *A = ds->mat[DS_MAT_A];
358: #if defined(PETSC_USE_DEBUG)
359: /* make sure diagonal 2x2 block is not broken */
361: #endif
362: if (trim) {
363: if (ds->extrarow) { /* clean extra row */
364: for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
365: }
366: ds->l = 0;
367: ds->k = 0;
368: ds->n = n;
369: ds->t = ds->n; /* truncated length equal to the new dimension */
370: } else {
371: if (ds->extrarow && ds->k==ds->n) {
372: /* copy entries of extra row to the new position, then clean last row */
373: for (i=l;i<n;i++) A[n+i*ld] = A[ds->n+i*ld];
374: for (i=l;i<ds->n;i++) A[ds->n+i*ld] = 0.0;
375: }
376: ds->k = (ds->extrarow)? n: 0;
377: ds->t = ds->n; /* truncated length equal to previous dimension */
378: ds->n = n;
379: }
380: PetscFunctionReturn(0);
381: }
383: PetscErrorCode DSCond_NHEP(DS ds,PetscReal *cond)
384: {
385: PetscScalar *work;
386: PetscReal *rwork;
387: PetscBLASInt *ipiv;
388: PetscBLASInt lwork,info,n,ld;
389: PetscReal hn,hin;
390: PetscScalar *A;
392: PetscBLASIntCast(ds->n,&n);
393: PetscBLASIntCast(ds->ld,&ld);
394: lwork = 8*ld;
395: DSAllocateWork_Private(ds,lwork,ld,ld);
396: work = ds->work;
397: rwork = ds->rwork;
398: ipiv = ds->iwork;
400: /* use workspace matrix W to avoid overwriting A */
401: DSAllocateMat_Private(ds,DS_MAT_W);
402: A = ds->mat[DS_MAT_W];
403: PetscArraycpy(A,ds->mat[DS_MAT_A],ds->ld*ds->ld);
405: /* norm of A */
406: if (ds->state<DS_STATE_INTERMEDIATE) hn = LAPACKlange_("I",&n,&n,A,&ld,rwork);
407: else hn = LAPACKlanhs_("I",&n,A,&ld,rwork);
409: /* norm of inv(A) */
410: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,A,&ld,ipiv,&info));
411: SlepcCheckLapackInfo("getrf",info);
412: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,A,&ld,ipiv,work,&lwork,&info));
413: SlepcCheckLapackInfo("getri",info);
414: hin = LAPACKlange_("I",&n,&n,A,&ld,rwork);
416: *cond = hn*hin;
417: PetscFunctionReturn(0);
418: }
420: PetscErrorCode DSTranslateHarmonic_NHEP(DS ds,PetscScalar tau,PetscReal beta,PetscBool recover,PetscScalar *gin,PetscReal *gammaout)
421: {
422: PetscInt i,j;
423: PetscBLASInt *ipiv,info,n,ld,one=1,ncol;
424: PetscScalar *A,*B,*Q,*g=gin,*ghat;
425: PetscScalar done=1.0,dmone=-1.0,dzero=0.0;
426: PetscReal gamma=1.0;
428: PetscBLASIntCast(ds->n,&n);
429: PetscBLASIntCast(ds->ld,&ld);
430: A = ds->mat[DS_MAT_A];
432: if (!recover) {
434: DSAllocateWork_Private(ds,0,0,ld);
435: ipiv = ds->iwork;
436: if (!g) {
437: DSAllocateWork_Private(ds,ld,0,0);
438: g = ds->work;
439: }
440: /* use workspace matrix W to factor A-tau*eye(n) */
441: DSAllocateMat_Private(ds,DS_MAT_W);
442: B = ds->mat[DS_MAT_W];
443: PetscArraycpy(B,A,ld*ld);
445: /* Vector g initially stores b = beta*e_n^T */
446: PetscArrayzero(g,n);
447: g[n-1] = beta;
449: /* g = (A-tau*eye(n))'\b */
450: for (i=0;i<n;i++) B[i+i*ld] -= tau;
451: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,B,&ld,ipiv,&info));
452: SlepcCheckLapackInfo("getrf",info);
453: PetscLogFlops(2.0*n*n*n/3.0);
454: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n,&one,B,&ld,ipiv,g,&ld,&info));
455: SlepcCheckLapackInfo("getrs",info);
456: PetscLogFlops(2.0*n*n-n);
458: /* A = A + g*b' */
459: for (i=0;i<n;i++) A[i+(n-1)*ld] += g[i]*beta;
461: } else { /* recover */
463: DSAllocateWork_Private(ds,ld,0,0);
464: ghat = ds->work;
465: Q = ds->mat[DS_MAT_Q];
467: /* g^ = -Q(:,idx)'*g */
468: PetscBLASIntCast(ds->l+ds->k,&ncol);
469: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&ncol,&dmone,Q,&ld,g,&one,&dzero,ghat,&one));
471: /* A = A + g^*b' */
472: for (i=0;i<ds->l+ds->k;i++)
473: for (j=ds->l;j<ds->l+ds->k;j++)
474: A[i+j*ld] += ghat[i]*Q[n-1+j*ld]*beta;
476: /* g~ = (I-Q(:,idx)*Q(:,idx)')*g = g+Q(:,idx)*g^ */
477: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&ncol,&done,Q,&ld,ghat,&one,&done,g,&one));
478: }
480: /* Compute gamma factor */
481: if (gammaout || (recover && ds->extrarow)) gamma = SlepcAbs(1.0,BLASnrm2_(&n,g,&one));
482: if (gammaout) *gammaout = gamma;
483: if (recover && ds->extrarow) {
484: for (j=ds->l;j<ds->l+ds->k;j++) A[ds->n+j*ld] *= gamma;
485: }
486: PetscFunctionReturn(0);
487: }
489: /*MC
490: DSNHEP - Dense Non-Hermitian Eigenvalue Problem.
492: Level: beginner
494: Notes:
495: The problem is expressed as A*X = X*Lambda, where A is the input matrix.
496: Lambda is a diagonal matrix whose diagonal elements are the arguments of
497: DSSolve(). After solve, A is overwritten with the upper quasi-triangular
498: matrix T of the (real) Schur form, A*Q = Q*T.
500: In the intermediate state A is reduced to upper Hessenberg form.
502: Computation of left eigenvectors is supported, but two-sided Krylov solvers
503: usually rely on the related DSNHEPTS.
505: Used DS matrices:
506: + DS_MAT_A - problem matrix
507: - DS_MAT_Q - orthogonal/unitary transformation that reduces to Hessenberg form
508: (intermediate step) or matrix of orthogonal Schur vectors
510: Implemented methods:
511: . 0 - Implicit QR (_hseqr)
513: .seealso: DSCreate(), DSSetType(), DSType
514: M*/
515: SLEPC_EXTERN PetscErrorCode DSCreate_NHEP(DS ds)
516: {
517: ds->ops->allocate = DSAllocate_NHEP;
518: ds->ops->view = DSView_NHEP;
519: ds->ops->vectors = DSVectors_NHEP;
520: ds->ops->solve[0] = DSSolve_NHEP;
521: ds->ops->sort = DSSort_NHEP;
522: ds->ops->sortperm = DSSortWithPermutation_NHEP;
523: ds->ops->synchronize = DSSynchronize_NHEP;
524: ds->ops->gettruncatesize = DSGetTruncateSize_Default;
525: ds->ops->truncate = DSTruncate_NHEP;
526: ds->ops->update = DSUpdateExtraRow_NHEP;
527: ds->ops->cond = DSCond_NHEP;
528: ds->ops->transharm = DSTranslateHarmonic_NHEP;
529: PetscFunctionReturn(0);
530: }