Actual source code: dsnhepts.c
slepc-3.18.1 2022-11-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: */
11: #include <slepc/private/dsimpl.h>
12: #include <slepcblaslapack.h>
14: typedef struct {
15: PetscScalar *wr,*wi; /* eigenvalues of B */
16: } DS_NHEPTS;
18: PetscErrorCode DSAllocate_NHEPTS(DS ds,PetscInt ld)
19: {
20: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
22: DSAllocateMat_Private(ds,DS_MAT_A);
23: DSAllocateMat_Private(ds,DS_MAT_B);
24: DSAllocateMat_Private(ds,DS_MAT_Q);
25: DSAllocateMat_Private(ds,DS_MAT_Z);
26: PetscFree(ds->perm);
27: PetscMalloc1(ld,&ds->perm);
28: PetscMalloc1(ld,&ctx->wr);
29: #if !defined(PETSC_USE_COMPLEX)
30: PetscMalloc1(ld,&ctx->wi);
31: #endif
32: return 0;
33: }
35: PetscErrorCode DSView_NHEPTS(DS ds,PetscViewer viewer)
36: {
37: PetscViewerFormat format;
39: PetscViewerGetFormat(viewer,&format);
40: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return 0;
41: DSViewMat(ds,viewer,DS_MAT_A);
42: DSViewMat(ds,viewer,DS_MAT_B);
43: if (ds->state>DS_STATE_INTERMEDIATE) {
44: DSViewMat(ds,viewer,DS_MAT_Q);
45: DSViewMat(ds,viewer,DS_MAT_Z);
46: }
47: if (ds->omat[DS_MAT_X]) DSViewMat(ds,viewer,DS_MAT_X);
48: if (ds->omat[DS_MAT_Y]) DSViewMat(ds,viewer,DS_MAT_Y);
49: return 0;
50: }
52: static PetscErrorCode DSVectors_NHEPTS_Eigen_Some(DS ds,PetscInt *k,PetscReal *rnorm,PetscBool left)
53: {
54: PetscInt i;
55: PetscBLASInt mm=1,mout,info,ld,n,*select,inc=1,cols=1,zero=0;
56: PetscScalar sone=1.0,szero=0.0;
57: PetscReal norm,done=1.0;
58: PetscBool iscomplex = PETSC_FALSE;
59: PetscScalar *X,*Y;
60: const PetscScalar *A,*Q;
62: PetscBLASIntCast(ds->n,&n);
63: PetscBLASIntCast(ds->ld,&ld);
64: DSAllocateWork_Private(ds,0,0,ld);
65: select = ds->iwork;
66: for (i=0;i<n;i++) select[i] = (PetscBLASInt)PETSC_FALSE;
68: /* compute k-th eigenvector Y of A */
69: MatDenseGetArrayRead(ds->omat[left?DS_MAT_B:DS_MAT_A],&A);
70: MatDenseGetArray(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X);
71: Y = X+(*k)*ld;
72: select[*k] = (PetscBLASInt)PETSC_TRUE;
73: #if !defined(PETSC_USE_COMPLEX)
74: if ((*k)<n-1 && A[(*k)+1+(*k)*ld]!=0.0) iscomplex = PETSC_TRUE;
75: mm = iscomplex? 2: 1;
76: if (iscomplex) select[(*k)+1] = (PetscBLASInt)PETSC_TRUE;
77: DSAllocateWork_Private(ds,3*ld,0,0);
78: PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,(PetscScalar*)A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,&info));
79: #else
80: DSAllocateWork_Private(ds,2*ld,ld,0);
81: PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","S",select,&n,(PetscScalar*)A,&ld,Y,&ld,Y,&ld,&mm,&mout,ds->work,ds->rwork,&info));
82: #endif
83: SlepcCheckLapackInfo("trevc",info);
85: MatDenseRestoreArrayRead(ds->omat[left?DS_MAT_B:DS_MAT_A],&A);
87: /* accumulate and normalize eigenvectors */
88: if (ds->state>=DS_STATE_CONDENSED) {
89: MatDenseGetArrayRead(ds->omat[left?DS_MAT_Z:DS_MAT_Q],&Q);
90: PetscArraycpy(ds->work,Y,mout*ld);
91: PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work,&inc,&szero,Y,&inc));
92: #if !defined(PETSC_USE_COMPLEX)
93: if (iscomplex) PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,Q,&ld,ds->work+ld,&inc,&szero,Y+ld,&inc));
94: #endif
95: MatDenseRestoreArrayRead(ds->omat[left?DS_MAT_Z:DS_MAT_Q],&Q);
96: cols = 1;
97: norm = BLASnrm2_(&n,Y,&inc);
98: #if !defined(PETSC_USE_COMPLEX)
99: if (iscomplex) {
100: norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Y+ld,&inc));
101: cols = 2;
102: }
103: #endif
104: PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Y,&ld,&info));
105: SlepcCheckLapackInfo("lascl",info);
106: }
108: /* set output arguments */
109: if (iscomplex) (*k)++;
110: if (rnorm) {
111: if (iscomplex) *rnorm = SlepcAbsEigenvalue(Y[n-1],Y[n-1+ld]);
112: else *rnorm = PetscAbsScalar(Y[n-1]);
113: }
114: MatDenseRestoreArray(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X);
115: return 0;
116: }
118: static PetscErrorCode DSVectors_NHEPTS_Eigen_All(DS ds,PetscBool left)
119: {
120: PetscInt i;
121: PetscBLASInt n,ld,mout,info,inc=1,cols,zero=0;
122: PetscBool iscomplex;
123: PetscScalar *X;
124: const PetscScalar *A;
125: PetscReal norm,done=1.0;
126: const char *back;
128: PetscBLASIntCast(ds->n,&n);
129: PetscBLASIntCast(ds->ld,&ld);
130: MatDenseGetArrayRead(ds->omat[left?DS_MAT_B:DS_MAT_A],&A);
131: MatDenseGetArrayWrite(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X);
132: if (ds->state>=DS_STATE_CONDENSED) {
133: /* DSSolve() has been called, backtransform with matrix Q */
134: back = "B";
135: MatCopy(ds->omat[left?DS_MAT_Z:DS_MAT_Q],ds->omat[left?DS_MAT_Y:DS_MAT_X],SAME_NONZERO_PATTERN);
136: } else back = "A";
137: #if !defined(PETSC_USE_COMPLEX)
138: DSAllocateWork_Private(ds,3*ld,0,0);
139: PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_("R",back,NULL,&n,(PetscScalar*)A,&ld,X,&ld,X,&ld,&n,&mout,ds->work,&info));
140: #else
141: DSAllocateWork_Private(ds,2*ld,ld,0);
142: PetscCallBLAS("LAPACKtrevc",LAPACKtrevc_("R",back,NULL,&n,(PetscScalar*)A,&ld,X,&ld,X,&ld,&n,&mout,ds->work,ds->rwork,&info));
143: #endif
144: SlepcCheckLapackInfo("trevc",info);
146: /* normalize eigenvectors */
147: for (i=0;i<n;i++) {
148: iscomplex = (i<n-1 && A[i+1+i*ld]!=0.0)? PETSC_TRUE: PETSC_FALSE;
149: cols = 1;
150: norm = BLASnrm2_(&n,X+i*ld,&inc);
151: #if !defined(PETSC_USE_COMPLEX)
152: if (iscomplex) {
153: norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,X+(i+1)*ld,&inc));
154: cols = 2;
155: }
156: #endif
157: PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,X+i*ld,&ld,&info));
158: SlepcCheckLapackInfo("lascl",info);
159: if (iscomplex) i++;
160: }
161: MatDenseRestoreArrayRead(ds->omat[left?DS_MAT_B:DS_MAT_A],&A);
162: MatDenseRestoreArrayWrite(ds->omat[left?DS_MAT_Y:DS_MAT_X],&X);
163: return 0;
164: }
166: PetscErrorCode DSVectors_NHEPTS(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
167: {
168: switch (mat) {
169: case DS_MAT_X:
171: if (j) DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_FALSE);
172: else DSVectors_NHEPTS_Eigen_All(ds,PETSC_FALSE);
173: break;
174: case DS_MAT_Y:
176: if (j) DSVectors_NHEPTS_Eigen_Some(ds,j,rnorm,PETSC_TRUE);
177: else DSVectors_NHEPTS_Eigen_All(ds,PETSC_TRUE);
178: break;
179: case DS_MAT_U:
180: case DS_MAT_V:
181: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
182: default:
183: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
184: }
185: return 0;
186: }
188: PetscErrorCode DSSort_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
189: {
190: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
191: PetscInt i,j,cont,id=0,*p,*idx,*idx2;
192: PetscReal s,t;
193: #if defined(PETSC_USE_COMPLEX)
194: Mat A,U;
195: #endif
198: PetscMalloc3(ds->ld,&idx,ds->ld,&idx2,ds->ld,&p);
199: DSSort_NHEP_Total(ds,DS_MAT_A,DS_MAT_Q,wr,wi);
200: #if defined(PETSC_USE_COMPLEX)
201: DSGetMat(ds,DS_MAT_B,&A);
202: MatConjugate(A);
203: DSRestoreMat(ds,DS_MAT_B,&A);
204: DSGetMat(ds,DS_MAT_Z,&U);
205: MatConjugate(U);
206: DSRestoreMat(ds,DS_MAT_Z,&U);
207: for (i=0;i<ds->n;i++) ctx->wr[i] = PetscConj(ctx->wr[i]);
208: #endif
209: DSSort_NHEP_Total(ds,DS_MAT_B,DS_MAT_Z,ctx->wr,ctx->wi);
210: /* check correct eigenvalue correspondence */
211: cont = 0;
212: for (i=0;i<ds->n;i++) {
213: if (SlepcAbsEigenvalue(ctx->wr[i]-wr[i],ctx->wi[i]-wi[i])>PETSC_SQRT_MACHINE_EPSILON) {idx2[cont] = i; idx[cont++] = i;}
214: p[i] = -1;
215: }
216: if (cont) {
217: for (i=0;i<cont;i++) {
218: t = PETSC_MAX_REAL;
219: for (j=0;j<cont;j++) if (idx2[j]!=-1 && (s=SlepcAbsEigenvalue(ctx->wr[idx[j]]-wr[idx[i]],ctx->wi[idx[j]]-wi[idx[i]]))<t) { id = j; t = s; }
220: p[idx[i]] = idx[id];
221: idx2[id] = -1;
222: }
223: for (i=0;i<ds->n;i++) if (p[i]==-1) p[i] = i;
224: DSSortWithPermutation_NHEP_Private(ds,p,DS_MAT_B,DS_MAT_Z,ctx->wr,ctx->wi);
225: }
226: #if defined(PETSC_USE_COMPLEX)
227: DSGetMat(ds,DS_MAT_B,&A);
228: MatConjugate(A);
229: DSRestoreMat(ds,DS_MAT_B,&A);
230: DSGetMat(ds,DS_MAT_Z,&U);
231: MatConjugate(U);
232: DSRestoreMat(ds,DS_MAT_Z,&U);
233: #endif
234: PetscFree3(idx,idx2,p);
235: return 0;
236: }
238: PetscErrorCode DSUpdateExtraRow_NHEPTS(DS ds)
239: {
240: PetscInt i;
241: PetscBLASInt n,ld,incx=1;
242: PetscScalar *A,*x,*y,one=1.0,zero=0.0;
243: const PetscScalar *Q;
245: PetscBLASIntCast(ds->n,&n);
246: PetscBLASIntCast(ds->ld,&ld);
247: DSAllocateWork_Private(ds,2*ld,0,0);
248: x = ds->work;
249: y = ds->work+ld;
250: MatDenseGetArray(ds->omat[DS_MAT_A],&A);
251: MatDenseGetArrayRead(ds->omat[DS_MAT_Q],&Q);
252: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
253: PetscCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
254: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
255: MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
256: MatDenseRestoreArrayRead(ds->omat[DS_MAT_Q],&Q);
257: MatDenseGetArray(ds->omat[DS_MAT_B],&A);
258: MatDenseGetArrayRead(ds->omat[DS_MAT_Z],&Q);
259: for (i=0;i<n;i++) x[i] = PetscConj(A[n+i*ld]);
260: PetscCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&one,Q,&ld,x,&incx,&zero,y,&incx));
261: for (i=0;i<n;i++) A[n+i*ld] = PetscConj(y[i]);
262: MatDenseRestoreArray(ds->omat[DS_MAT_B],&A);
263: MatDenseRestoreArrayRead(ds->omat[DS_MAT_Z],&Q);
264: ds->k = n;
265: return 0;
266: }
268: PetscErrorCode DSSolve_NHEPTS(DS ds,PetscScalar *wr,PetscScalar *wi)
269: {
270: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
272: #if !defined(PETSC_USE_COMPLEX)
274: #endif
275: DSSolve_NHEP_Private(ds,DS_MAT_A,DS_MAT_Q,wr,wi);
276: DSSolve_NHEP_Private(ds,DS_MAT_B,DS_MAT_Z,ctx->wr,ctx->wi);
277: return 0;
278: }
280: #if !defined(PETSC_HAVE_MPIUNI)
281: PetscErrorCode DSSynchronize_NHEPTS(DS ds,PetscScalar eigr[],PetscScalar eigi[])
282: {
283: PetscInt ld=ds->ld,l=ds->l,k;
284: PetscMPIInt n,rank,off=0,size,ldn;
285: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
286: PetscScalar *A,*B,*Q,*Z;
288: k = 2*(ds->n-l)*ld;
289: if (ds->state>DS_STATE_RAW) k += 2*(ds->n-l)*ld;
290: if (eigr) k += ds->n-l;
291: if (eigi) k += ds->n-l;
292: if (ctx->wr) k += ds->n-l;
293: if (ctx->wi) k += ds->n-l;
294: DSAllocateWork_Private(ds,k,0,0);
295: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
296: PetscMPIIntCast(ds->n-l,&n);
297: PetscMPIIntCast(ld*(ds->n-l),&ldn);
298: MatDenseGetArray(ds->omat[DS_MAT_A],&A);
299: MatDenseGetArray(ds->omat[DS_MAT_B],&B);
300: if (ds->state>DS_STATE_RAW) {
301: MatDenseGetArray(ds->omat[DS_MAT_Q],&Q);
302: MatDenseGetArray(ds->omat[DS_MAT_Z],&Z);
303: }
304: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
305: if (!rank) {
306: MPI_Pack(A+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
307: MPI_Pack(B+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
308: if (ds->state>DS_STATE_RAW) {
309: MPI_Pack(Q+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
310: MPI_Pack(Z+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
311: }
312: if (eigr) MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
313: #if !defined(PETSC_USE_COMPLEX)
314: if (eigi) MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
315: #endif
316: if (ctx->wr) MPI_Pack(ctx->wr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
317: if (ctx->wi) MPI_Pack(ctx->wi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
318: }
319: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
320: if (rank) {
321: MPI_Unpack(ds->work,size,&off,A+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
322: MPI_Unpack(ds->work,size,&off,B+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
323: if (ds->state>DS_STATE_RAW) {
324: MPI_Unpack(ds->work,size,&off,Q+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
325: MPI_Unpack(ds->work,size,&off,Z+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
326: }
327: if (eigr) MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
328: #if !defined(PETSC_USE_COMPLEX)
329: if (eigi) MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
330: #endif
331: if (ctx->wr) MPI_Unpack(ds->work,size,&off,ctx->wr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
332: if (ctx->wi) MPI_Unpack(ds->work,size,&off,ctx->wi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
333: }
334: MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
335: MatDenseRestoreArray(ds->omat[DS_MAT_B],&B);
336: if (ds->state>DS_STATE_RAW) {
337: MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q);
338: MatDenseRestoreArray(ds->omat[DS_MAT_Z],&Z);
339: }
340: return 0;
341: }
342: #endif
344: PetscErrorCode DSGetTruncateSize_NHEPTS(DS ds,PetscInt l,PetscInt n,PetscInt *k)
345: {
346: #if !defined(PETSC_USE_COMPLEX)
347: const PetscScalar *A,*B;
348: #endif
350: #if !defined(PETSC_USE_COMPLEX)
351: MatDenseGetArrayRead(ds->omat[DS_MAT_A],&A);
352: MatDenseGetArrayRead(ds->omat[DS_MAT_B],&B);
353: if (A[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0 || B[l+(*k)+(l+(*k)-1)*ds->ld] != 0.0) {
354: if (l+(*k)<n-1) (*k)++;
355: else (*k)--;
356: }
357: MatDenseRestoreArrayRead(ds->omat[DS_MAT_A],&A);
358: MatDenseRestoreArrayRead(ds->omat[DS_MAT_B],&B);
359: #endif
360: return 0;
361: }
363: PetscErrorCode DSTruncate_NHEPTS(DS ds,PetscInt n,PetscBool trim)
364: {
365: PetscInt i,ld=ds->ld,l=ds->l;
366: PetscScalar *A,*B;
368: MatDenseGetArray(ds->omat[DS_MAT_A],&A);
369: MatDenseGetArray(ds->omat[DS_MAT_B],&B);
370: #if defined(PETSC_USE_DEBUG)
371: /* make sure diagonal 2x2 block is not broken */
373: #endif
374: if (trim) {
375: if (ds->extrarow) { /* clean extra row */
376: for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
377: }
378: ds->l = 0;
379: ds->k = 0;
380: ds->n = n;
381: ds->t = ds->n; /* truncated length equal to the new dimension */
382: } else {
383: if (ds->extrarow && ds->k==ds->n) {
384: /* copy entries of extra row to the new position, then clean last row */
385: for (i=l;i<n;i++) { A[n+i*ld] = A[ds->n+i*ld]; B[n+i*ld] = B[ds->n+i*ld]; }
386: for (i=l;i<ds->n;i++) { A[ds->n+i*ld] = 0.0; B[ds->n+i*ld] = 0.0; }
387: }
388: ds->k = (ds->extrarow)? n: 0;
389: ds->t = ds->n; /* truncated length equal to previous dimension */
390: ds->n = n;
391: }
392: MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
393: MatDenseRestoreArray(ds->omat[DS_MAT_B],&B);
394: return 0;
395: }
397: PetscErrorCode DSDestroy_NHEPTS(DS ds)
398: {
399: DS_NHEPTS *ctx = (DS_NHEPTS*)ds->data;
401: if (ctx->wr) PetscFree(ctx->wr);
402: if (ctx->wi) PetscFree(ctx->wi);
403: PetscFree(ds->data);
404: return 0;
405: }
407: PetscErrorCode DSMatGetSize_NHEPTS(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
408: {
409: *rows = ((t==DS_MAT_A || t==DS_MAT_B) && ds->extrarow)? ds->n+1: ds->n;
410: *cols = ds->n;
411: return 0;
412: }
414: /*MC
415: DSNHEPTS - Dense Non-Hermitian Eigenvalue Problem (special variant intended
416: for two-sided Krylov solvers).
418: Level: beginner
420: Notes:
421: Two related problems are solved, A*X = X*Lambda and B*Y = Y*Lambda', where A and
422: B are supposed to come from the Arnoldi factorizations of a certain matrix and its
423: (conjugate) transpose, respectively. Hence, in exact arithmetic the columns of Y
424: are equal to the left eigenvectors of A. Lambda is a diagonal matrix whose diagonal
425: elements are the arguments of DSSolve(). After solve, A is overwritten with the
426: upper quasi-triangular matrix T of the (real) Schur form, A*Q = Q*T, and similarly
427: another (real) Schur relation is computed, B*Z = Z*S, overwriting B.
429: In the intermediate state A and B are reduced to upper Hessenberg form.
431: When left eigenvectors DS_MAT_Y are requested, right eigenvectors of B are returned,
432: while DS_MAT_X contains right eigenvectors of A.
434: Used DS matrices:
435: + DS_MAT_A - first problem matrix obtained from Arnoldi
436: . DS_MAT_B - second problem matrix obtained from Arnoldi on the transpose
437: . DS_MAT_Q - orthogonal/unitary transformation that reduces A to Hessenberg form
438: (intermediate step) or matrix of orthogonal Schur vectors of A
439: - DS_MAT_Z - orthogonal/unitary transformation that reduces B to Hessenberg form
440: (intermediate step) or matrix of orthogonal Schur vectors of B
442: Implemented methods:
443: . 0 - Implicit QR (_hseqr)
445: .seealso: DSCreate(), DSSetType(), DSType
446: M*/
447: SLEPC_EXTERN PetscErrorCode DSCreate_NHEPTS(DS ds)
448: {
449: DS_NHEPTS *ctx;
451: PetscNew(&ctx);
452: ds->data = (void*)ctx;
454: ds->ops->allocate = DSAllocate_NHEPTS;
455: ds->ops->view = DSView_NHEPTS;
456: ds->ops->vectors = DSVectors_NHEPTS;
457: ds->ops->solve[0] = DSSolve_NHEPTS;
458: ds->ops->sort = DSSort_NHEPTS;
459: #if !defined(PETSC_HAVE_MPIUNI)
460: ds->ops->synchronize = DSSynchronize_NHEPTS;
461: #endif
462: ds->ops->gettruncatesize = DSGetTruncateSize_NHEPTS;
463: ds->ops->truncate = DSTruncate_NHEPTS;
464: ds->ops->update = DSUpdateExtraRow_NHEPTS;
465: ds->ops->destroy = DSDestroy_NHEPTS;
466: ds->ops->matgetsize = DSMatGetSize_NHEPTS;
467: return 0;
468: }