Actual source code: invit.c
slepc-3.16.2 2022-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: struct HRtr
15: {
16: PetscScalar *data;
17: PetscInt m;
18: PetscInt idx[2];
19: PetscInt n[2];
20: PetscScalar tau[2];
21: PetscReal alpha;
22: PetscReal cs;
23: PetscReal sn;
24: PetscInt type;
25: };
27: /*
28: Generates a hyperbolic rotation
29: if x1*x1 - x2*x2 != 0
30: r = sqrt(|x1*x1 - x2*x2|)
31: c = x1/r s = x2/r
33: | c -s||x1| |d*r|
34: |-s c||x2| = | 0 |
35: where d = 1 for type==1 and -1 for type==2
36: Returns the condition number of the reduction
37: */
38: static PetscErrorCode HRGen(PetscReal x1,PetscReal x2,PetscInt *type,PetscReal *c,PetscReal *s,PetscReal *r,PetscReal *cond)
39: {
40: PetscReal t,n2,xa,xb;
41: PetscInt type_;
44: if (x2==0.0) {
45: *r = PetscAbsReal(x1); *c = (x1>=0.0)?1.0:-1.0; *s = 0.0;
46: if (type) *type = 1;
47: return(0);
48: }
49: if (PetscAbsReal(x1) == PetscAbsReal(x2)) {
50: /* hyperbolic rotation doesn't exist */
51: *c = *s = *r = 0.0;
52: if (type) *type = 0;
53: *cond = PETSC_MAX_REAL;
54: return(0);
55: }
57: if (PetscAbsReal(x1)>PetscAbsReal(x2)) {
58: xa = x1; xb = x2; type_ = 1;
59: } else {
60: xa = x2; xb = x1; type_ = 2;
61: }
62: t = xb/xa;
63: n2 = PetscAbsReal(1 - t*t);
64: *r = PetscSqrtReal(n2)*PetscAbsReal(xa);
65: *c = x1/(*r);
66: *s = x2/(*r);
67: if (type_ == 2) *r *= -1;
68: if (type) *type = type_;
69: if (cond) *cond = (PetscAbsReal(*c) + PetscAbsReal(*s))/PetscAbsReal(PetscAbsReal(*c) - PetscAbsReal(*s));
70: return(0);
71: }
73: /*
74: |c s|
75: Applies an hyperbolic rotator |s c|
76: |c s|
77: [x1 x2]|s c|
78: */
79: static PetscErrorCode HRApply(PetscInt n,PetscScalar *x1,PetscInt inc1,PetscScalar *x2,PetscInt inc2,PetscReal c,PetscReal s)
80: {
81: PetscInt i;
82: PetscReal t;
83: PetscScalar tmp;
86: if (PetscAbsReal(c)>PetscAbsReal(s)) { /* Type I */
87: t = s/c;
88: for (i=0;i<n;i++) {
89: x1[i*inc1] = c*x1[i*inc1] + s*x2[i*inc2];
90: x2[i*inc2] = t*x1[i*inc1] + x2[i*inc2]/c;
91: }
92: } else { /* Type II */
93: t = c/s;
94: for (i=0;i<n;i++) {
95: tmp = x1[i*inc1];
96: x1[i*inc1] = c*x1[i*inc1] + s*x2[i*inc2];
97: x2[i*inc2] = t*x1[i*inc1] + tmp/s;
98: }
99: }
100: return(0);
101: }
103: /*
104: Reduction to tridiagonal-diagonal form (see F. Tisseur, SIMAX 26(1), 2004).
106: Input:
107: A symmetric (only lower triangular part is referred)
108: s vector +1 and -1 (signature matrix)
109: Output:
110: d,e
111: s
112: Q s-orthogonal matrix with Q^T*A*Q = T (symmetric tridiagonal matrix)
113: */
114: static PetscErrorCode TridiagDiag_HHR(PetscInt n,PetscScalar *A,PetscInt lda,PetscReal *s,PetscScalar* Q,PetscInt ldq,PetscBool flip,PetscReal *d,PetscReal *e,PetscInt *perm_,PetscScalar *work,PetscReal *rwork,PetscBLASInt *iwork)
115: {
117: PetscInt i,j,k,*ii,*jj,i0=0,ik=0,tmp,type;
118: PetscInt nwu=0;
119: PetscReal *ss,cond=1.0,cs,sn,r;
120: PetscScalar tau,t,*AA;
121: PetscBLASInt n0,n1,ni,inc=1,m,n_,lda_,ldq_,*perm;
122: PetscBool breakdown = PETSC_TRUE;
125: if (n<3) {
126: if (n==1) Q[0]=1;
127: if (n==2) {
128: Q[0] = Q[1+ldq] = 1;
129: Q[1] = Q[ldq] = 0;
130: }
131: return(0);
132: }
133: PetscBLASIntCast(lda,&lda_);
134: PetscBLASIntCast(n,&n_);
135: PetscBLASIntCast(ldq,&ldq_);
136: ss = rwork;
137: perm = iwork;
138: AA = work;
139: for (i=0;i<n;i++) {
140: PetscArraycpy(AA+i*n,A+i*lda,n);
141: }
142: nwu += n*n;
143: k=0;
144: while (breakdown && k<n) {
145: breakdown = PETSC_FALSE;
146: /* Classify (and flip) A and s according to sign */
147: if (flip) {
148: for (i=0;i<n;i++) {
149: perm[i] = n-1-perm_[i];
150: if (perm[i]==0) i0 = i;
151: if (perm[i]==k) ik = i;
152: }
153: } else {
154: for (i=0;i<n;i++) {
155: perm[i] = perm_[i];
156: if (perm[i]==0) i0 = i;
157: if (perm[i]==k) ik = i;
158: }
159: }
160: perm[ik] = 0;
161: perm[i0] = k;
162: i=1;
163: while (i<n-1 && s[perm[i-1]]==s[perm[0]]) {
164: if (s[perm[i]]!=s[perm[0]]) {
165: j=i+1;
166: while (j<n-1 && s[perm[j]]!=s[perm[0]])j++;
167: tmp = perm[i]; perm[i] = perm[j]; perm[j] = tmp;
168: }
169: i++;
170: }
171: for (i=0;i<n;i++) {
172: ss[i] = s[perm[i]];
173: }
174: if (flip) {
175: ii = &j;
176: jj = &i;
177: } else {
178: ii = &i;
179: jj = &j;
180: }
181: for (i=0;i<n;i++)
182: for (j=0;j<n;j++)
183: A[i+j*lda] = AA[perm[*ii]+perm[*jj]*n];
184: /* Initialize Q */
185: for (i=0;i<n;i++) {
186: PetscArrayzero(Q+i*ldq,n);
187: Q[perm[i]+i*ldq] = 1.0;
188: }
189: for (ni=1;ni<n && ss[ni]==ss[0]; ni++);
190: n0 = ni-1;
191: n1 = n_-ni;
192: for (j=0;j<n-2;j++) {
193: PetscBLASIntCast(n-j-1,&m);
194: /* Forming and applying reflectors */
195: if (n0 > 1) {
196: PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n0,A+ni-n0+j*lda,A+ni-n0+j*lda+1,&inc,&tau));
197: /* Apply reflector */
198: if (PetscAbsScalar(tau) != 0.0) {
199: t=*(A+ni-n0+j*lda); *(A+ni-n0+j*lda)=1.0;
200: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&m,&n0,A+ni-n0+j*lda,&inc,&tau,A+j+1+(j+1)*lda,&lda_,work+nwu));
201: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0,&m,A+ni-n0+j*lda,&inc,&tau,A+j+1+(j+1)*lda,&lda_,work+nwu));
202: /* Update Q */
203: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n0,A+ni-n0+j*lda,&inc,&tau,Q+(j+1)*ldq,&ldq_,work+nwu));
204: *(A+ni-n0+j*lda) = t;
205: for (i=1;i<n0;i++) {
206: *(A+ni-n0+j*lda+i) = 0.0; *(A+j+(ni-n0+i)*lda) = 0.0;
207: }
208: *(A+j+(ni-n0)*lda) = *(A+ni-n0+j*lda);
209: }
210: }
211: if (n1 > 1) {
212: PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n1,A+n-n1+j*lda,A+n-n1+j*lda+1,&inc,&tau));
213: /* Apply reflector */
214: if (PetscAbsScalar(tau) != 0.0) {
215: t=*(A+n-n1+j*lda); *(A+n-n1+j*lda)=1.0;
216: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&m,&n1,A+n-n1+j*lda,&inc,&tau,A+j+1+(n-n1)*lda,&lda_,work+nwu));
217: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n1,&m,A+n-n1+j*lda,&inc,&tau,A+n-n1+(j+1)*lda,&lda_,work+nwu));
218: /* Update Q */
219: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n1,A+n-n1+j*lda,&inc,&tau,Q+(n-n1)*ldq,&ldq_,work+nwu));
220: *(A+n-n1+j*lda) = t;
221: for (i=1;i<n1;i++) {
222: *(A+n-n1+i+j*lda) = 0.0; *(A+j+(n-n1+i)*lda) = 0.0;
223: }
224: *(A+j+(n-n1)*lda) = *(A+n-n1+j*lda);
225: }
226: }
227: /* Hyperbolic rotation */
228: if (n0 > 0 && n1 > 0) {
229: HRGen(PetscRealPart(A[ni-n0+j*lda]),PetscRealPart(A[n-n1+j*lda]),&type,&cs,&sn,&r,&cond);
230: /* Check condition number */
231: if (cond > 1.0/(10*PETSC_SQRT_MACHINE_EPSILON)) {
232: breakdown = PETSC_TRUE;
233: k++;
234: if (k==n || flip) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Breakdown in construction of hyperbolic transformation");
235: break;
236: }
237: A[ni-n0+j*lda] = r; A[n-n1+j*lda] = 0.0;
238: A[j+(ni-n0)*lda] = r; A[j+(n-n1)*lda] = 0.0;
239: /* Apply to A */
240: HRApply(m,A+j+1+(ni-n0)*lda,1,A+j+1+(n-n1)*lda,1,cs,-sn);
241: HRApply(m,A+ni-n0+(j+1)*lda,lda,A+n-n1+(j+1)*lda,lda,cs,-sn);
243: /* Update Q */
244: HRApply(n,Q+(ni-n0)*ldq,1,Q+(n-n1)*ldq,1,cs,-sn);
245: if (type==2) {
246: ss[ni-n0] = -ss[ni-n0]; ss[n-n1] = -ss[n-n1];
247: n0++;ni++;n1--;
248: }
249: }
250: if (n0>0) n0--;
251: else n1--;
252: }
253: }
255: /* flip matrices */
256: if (flip) {
257: for (i=0;i<n-1;i++) {
258: d[i] = PetscRealPart(A[n-i-1+(n-i-1)*lda]);
259: e[i] = PetscRealPart(A[n-i-1+(n-i-2)*lda]);
260: s[i] = ss[n-i-1];
261: }
262: s[n-1] = ss[0];
263: d[n-1] = PetscRealPart(A[0]);
264: for (i=0;i<n;i++) {
265: PetscArraycpy(work+i*n,Q+i*ldq,n);
266: }
267: for (i=0;i<n;i++)
268: for (j=0;j<n;j++)
269: Q[i+j*ldq] = work[i+(n-j-1)*n];
270: } else {
271: for (i=0;i<n-1;i++) {
272: d[i] = PetscRealPart(A[i+i*lda]);
273: e[i] = PetscRealPart(A[i+1+i*lda]);
274: s[i] = ss[i];
275: }
276: s[n-1] = ss[n-1];
277: d[n-1] = PetscRealPart(A[n-1 + (n-1)*lda]);
278: }
279: return(0);
280: }
282: static PetscErrorCode MadeHRtr(PetscInt sz,PetscInt n,PetscInt idx0,PetscInt n0,PetscInt idx1,PetscInt n1,struct HRtr *tr1,struct HRtr *tr2,PetscReal *ncond,PetscScalar *work)
283: {
285: PetscScalar *x,*y;
286: PetscReal ncond2=1.0;
287: PetscBLASInt n0_,n1_,inc=1;
290: /* Hyperbolic transformation to make zeros in x */
291: x = tr1->data;
292: tr1->n[0] = n0;
293: tr1->n[1] = n1;
294: tr1->idx[0] = idx0;
295: tr1->idx[1] = idx1;
296: PetscBLASIntCast(tr1->n[0],&n0_);
297: PetscBLASIntCast(tr1->n[1],&n1_);
298: if (tr1->n[0] > 1) {
299: PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n0_,x+tr1->idx[0],x+tr1->idx[0]+1,&inc,tr1->tau));
300: }
301: if (tr1->n[1]> 1) {
302: PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n1_,x+tr1->idx[1],x+tr1->idx[1]+1,&inc,tr1->tau+1));
303: }
304: if (tr1->idx[0]<tr1->idx[1]) {
305: HRGen(PetscRealPart(x[tr1->idx[0]]),PetscRealPart(x[tr1->idx[1]]),&(tr1->type),&(tr1->cs),&(tr1->sn),&(tr1->alpha),ncond);
306: } else {
307: tr1->alpha = PetscRealPart(x[tr1->idx[0]]);
308: *ncond = 1.0;
309: }
310: if (sz==2) {
311: y = tr2->data;
312: /* Apply first transformation to second column */
313: if (tr1->n[0] > 1 && PetscAbsScalar(tr1->tau[0])!=0.0) {
314: x[tr1->idx[0]] = 1.0;
315: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0_,&inc,x+tr1->idx[0],&inc,tr1->tau,y+tr1->idx[0],&n0_,work));
316: }
317: if (tr1->n[1] > 1 && PetscAbsScalar(tr1->tau[1])!=0.0) {
318: x[tr1->idx[1]] = 1.0;
319: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n1_,&inc,x+tr1->idx[1],&inc,tr1->tau+1,y+tr1->idx[1],&n1_,work));
320: }
321: if (tr1->idx[0]<tr1->idx[1]) {
322: HRApply(1,y+tr1->idx[0],1,y+tr1->idx[1],1,tr1->cs,-tr1->sn);
323: }
324: tr2->n[0] = tr1->n[0];
325: tr2->n[1] = tr1->n[1];
326: tr2->idx[0] = tr1->idx[0];
327: tr2->idx[1] = tr1->idx[1];
328: if (tr1->idx[0]<tr1->idx[1] && tr1->type==2) {
329: tr2->idx[1]++; tr2->n[1]--; tr2->n[0]++;
330: }
331: if (tr2->n[0]>0) {
332: tr2->n[0]--; tr2->idx[0]++;
333: if (tr2->n[1]==0) tr2->idx[1] = tr2->idx[0];
334: } else {
335: tr2->n[1]--; tr2->idx[1]++; tr2->idx[0] = tr2->idx[1];
336: }
337: /* Hyperbolic transformation to make zeros in y */
338: PetscBLASIntCast(tr2->n[0],&n0_);
339: PetscBLASIntCast(tr2->n[1],&n1_);
340: if (tr2->n[0] > 1) {
341: PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n0_,y+tr2->idx[0],y+tr2->idx[0]+1,&inc,tr2->tau));
342: }
343: if (tr2->n[1]> 1) {
344: PetscStackCallBLAS("LAPACKlarfg",LAPACKlarfg_(&n1_,y+tr2->idx[1],y+tr2->idx[1]+1,&inc,tr2->tau+1));
345: }
346: if (tr2->idx[0]<tr2->idx[1]) {
347: HRGen(PetscRealPart(y[tr2->idx[0]]),PetscRealPart(y[tr2->idx[1]]),&(tr2->type),&(tr2->cs),&(tr2->sn),&(tr2->alpha),&ncond2);
348: } else {
349: tr2->alpha = PetscRealPart(y[tr2->idx[0]]);
350: ncond2 = 1.0;
351: }
352: if (ncond2>*ncond) *ncond = ncond2;
353: }
354: return(0);
355: }
357: /*
358: Auxiliary function to try perform one iteration of hr routine,
359: checking condition number. If it is < tolD, apply the
360: transformation to H and R, if not, ok=false and it do nothing
361: tolE, tolerance to exchange complex pairs to improve conditioning
362: */
363: static PetscErrorCode TryHRIt(PetscInt n,PetscInt j,PetscInt sz,PetscScalar *H,PetscInt ldh,PetscScalar *R,PetscInt ldr,PetscReal *s,PetscBool *exg,PetscBool *ok,PetscInt *n0,PetscInt *n1,PetscInt *idx0,PetscInt *idx1,PetscReal *cond,PetscScalar *work)
364: {
366: struct HRtr *tr1,*tr2,tr1_t,tr2_t,tr1_te,tr2_te;
367: PetscScalar *x,*y;
368: PetscReal ncond,ncond_e;
369: PetscInt nwu=0,i,d=1;
370: PetscBLASInt n0_,n1_,inc=1,mh,mr,n_,ldr_,ldh_;
371: PetscReal tolD = 1e+5;
374: if (cond) *cond = 1.0;
375: PetscBLASIntCast(n,&n_);
376: PetscBLASIntCast(ldr,&ldr_);
377: PetscBLASIntCast(ldh,&ldh_);
378: x = work+nwu;
379: nwu += n;
380: PetscArraycpy(x,R+j*ldr,n);
381: *exg = PETSC_FALSE;
382: *ok = PETSC_TRUE;
383: tr1_t.data = x;
384: if (sz==1) {
385: /* Hyperbolic transformation to make zeros in x */
386: MadeHRtr(sz,n,*idx0,*n0,*idx1,*n1,&tr1_t,NULL,&ncond,work+nwu);
387: /* Check condition number to single column*/
388: if (ncond>tolD) *ok = PETSC_FALSE;
389: tr1 = &tr1_t;
390: tr2 = &tr2_t;
391: } else {
392: y = work+nwu;
393: nwu += n;
394: PetscArraycpy(y,R+(j+1)*ldr,n);
395: tr2_t.data = y;
396: MadeHRtr(sz,n,*idx0,*n0,*idx1,*n1,&tr1_t,&tr2_t,&ncond,work+nwu);
397: /* Computing hyperbolic transformations also for exchanged vectors */
398: tr1_te.data = work+nwu;
399: nwu += n;
400: PetscArraycpy(tr1_te.data,R+(j+1)*ldr,n);
401: tr2_te.data = work+nwu;
402: nwu += n;
403: PetscArraycpy(tr2_te.data,R+j*ldr,n);
404: MadeHRtr(sz,n,*idx0,*n0,*idx1,*n1,&tr1_te,&tr2_te,&ncond_e,work+nwu);
405: if (ncond > d*ncond_e) {
406: *exg = PETSC_TRUE;
407: tr1 = &tr1_te;
408: tr2 = &tr2_te;
409: ncond = ncond_e;
410: } else {
411: tr1 = &tr1_t;
412: tr2 = &tr2_t;
413: }
414: if (ncond>tolD) *ok = PETSC_FALSE;
415: }
416: if (*ok) {
417: /* Everything is OK, apply transformations to R and H */
418: /* First column */
419: if (cond && *cond<ncond) *cond = ncond;
420: x = tr1->data;
421: PetscBLASIntCast(tr1->n[0],&n0_);
422: PetscBLASIntCast(tr1->n[1],&n1_);
423: PetscBLASIntCast(n-j-sz,&mr);
424: if (tr1->n[0] > 1 && PetscAbsScalar(tr1->tau[0])!=0.0) {
425: x[tr1->idx[0]] = 1.0;
426: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0_,&mr,x+tr1->idx[0],&inc,tr1->tau,R+(j+sz)*ldr+tr1->idx[0],&ldr_,work+nwu));
427: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n0_,x+tr1->idx[0],&inc,tr1->tau,H+(tr1->idx[0])*ldh,&ldh_,work+nwu));
428: }
429: if (tr1->n[1] > 1 && PetscAbsScalar(tr1->tau[1])!=0.0) {
430: x[tr1->idx[1]] = 1.0;
431: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n1_,&mr,x+tr1->idx[1],&inc,tr1->tau+1,R+(j+sz)*ldr+tr1->idx[1],&ldr_,work+nwu));
432: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n1_,x+tr1->idx[1],&inc,tr1->tau+1,H+(tr1->idx[1])*ldh,&ldh_,work+nwu));
433: }
434: if (tr1->idx[0]<tr1->idx[1]) {
435: HRApply(mr,R+(j+sz)*ldr+tr1->idx[0],ldr,R+(j+sz)*ldr+tr1->idx[1],ldr,tr1->cs,-tr1->sn);
436: if (tr1->type==1) {
437: HRApply(n,H+(tr1->idx[0])*ldh,1,H+(tr1->idx[1])*ldh,1,tr1->cs,tr1->sn);
438: } else {
439: HRApply(n,H+(tr1->idx[0])*ldh,1,H+(tr1->idx[1])*ldh,1,-tr1->cs,-tr1->sn);
440: s[tr1->idx[0]] = -s[tr1->idx[0]];
441: s[tr1->idx[1]] = -s[tr1->idx[1]];
442: }
443: }
444: for (i=0;i<tr1->idx[0];i++) *(R+j*ldr+i) = x[i];
445: for (i=tr1->idx[0]+1;i<n;i++) *(R+j*ldr+i) = 0.0;
446: *(R+j*ldr+tr1->idx[0]) = tr1->alpha;
447: if (sz==2) {
448: y = tr2->data;
449: /* Second column */
450: PetscBLASIntCast(tr2->n[0],&n0_);
451: PetscBLASIntCast(tr2->n[1],&n1_);
452: PetscBLASIntCast(n-j-sz,&mr);
453: PetscBLASIntCast(n-tr2->idx[0],&mh);
454: if (tr2->n[0] > 1 && PetscAbsScalar(tr2->tau[0])!=0.0) {
455: y[tr2->idx[0]] = 1.0;
456: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n0_,&mr,y+tr2->idx[0],&inc,tr2->tau,R+(j+2)*ldr+tr2->idx[0],&ldr_,work+nwu));
457: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n0_,y+tr2->idx[0],&inc,tr2->tau,H+(tr2->idx[0])*ldh,&ldh_,work+nwu));
458: }
459: if (tr2->n[1] > 1 && PetscAbsScalar(tr2->tau[1])!=0.0) {
460: y[tr2->idx[1]] = 1.0;
461: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("L",&n1_,&mr,y+tr2->idx[1],&inc,tr2->tau+1,R+(j+2)*ldr+tr2->idx[1],&ldr_,work+nwu));
462: PetscStackCallBLAS("LAPACKlarf",LAPACKlarf_("R",&n_,&n1_,y+tr2->idx[1],&inc,tr2->tau+1,H+(tr2->idx[1])*ldh,&ldh_,work+nwu));
463: }
464: if (tr2->idx[0]<tr2->idx[1]) {
465: HRApply(mr,R+(j+2)*ldr+tr2->idx[0],ldr,R+(j+2)*ldr+tr2->idx[1],ldr,tr2->cs,-tr2->sn);
466: if (tr2->type==1) {
467: HRApply(n,H+(tr2->idx[0])*ldh,1,H+(tr2->idx[1])*ldh,1,tr2->cs,tr2->sn);
468: } else {
469: HRApply(n,H+(tr2->idx[0])*ldh,1,H+(tr2->idx[1])*ldh,1,-tr2->cs,-tr2->sn);
470: s[tr2->idx[0]] = -s[tr2->idx[0]];
471: s[tr2->idx[1]] = -s[tr2->idx[1]];
472: }
473: }
474: for (i=0;i<tr2->idx[0]-1;i++) *(R+(j+1)*ldr+i) = y[i];
475: *(R+(j+1)*ldr+tr2->idx[0]-1) = y[tr2->idx[0]-1];
476: for (i=tr2->idx[0]+1;i<n;i++) *(R+(j+1)*ldr+i) = 0.0;
477: *(R+(j+1)*ldr+tr2->idx[0]) = tr2->alpha;
478: *n0 = tr2->n[0];
479: *n1 = tr2->n[1];
480: *idx0 = tr2->idx[0];
481: *idx1 = tr2->idx[1];
482: if (tr2->idx[0]<tr2->idx[1] && tr2->type==2) {
483: (*idx1)++; (*n1)--; (*n0)++;
484: }
485: } else {
486: *n0 = tr1->n[0];
487: *n1 = tr1->n[1];
488: *idx0 = tr1->idx[0];
489: *idx1 = tr1->idx[1];
490: if (tr1->idx[0]<tr1->idx[1] && tr1->type==2) {
491: (*idx1)++; (*n1)--; (*n0)++;
492: }
493: }
494: if (*n0>0) {
495: (*n0)--; (*idx0)++;
496: if (*n1==0) *idx1 = *idx0;
497: } else {
498: (*n1)--; (*idx1)++; *idx0 = *idx1;
499: }
500: }
501: return(0);
502: }
504: /*
505: compute V = HR whit H s-orthogonal and R upper triangular
506: */
507: static PetscErrorCode PseudoOrthog_HR(PetscInt *nv,PetscScalar *V,PetscInt ldv,PetscReal *s,PetscScalar *R,PetscInt ldr,PetscBLASInt *perm,PetscBLASInt *cmplxEig,PetscBool *breakdown,PetscScalar *work)
508: {
510: PetscInt i,j,n,n0,n1,np,idx0,idx1,sz=1,k=0,t1,t2,nwu=0;
511: PetscScalar *col1,*col2;
512: PetscBool exg=PETSC_FALSE,ok=PETSC_FALSE;
515: n = *nv;
516: col1 = work+nwu;
517: nwu += n;
518: col2 = work+nwu;
519: nwu += n;
520: /* Sort R and s according to sing(s) */
521: np = 0;
522: for (i=0;i<n;i++) if (s[i]>0) np++;
523: if (s[0]>0) n1 = np;
524: else n1 = n-np;
525: n0 = 0;
526: for (i=0;i<n;i++) {
527: if (s[i]==s[0]) {
528: s[n0] = s[0];
529: perm[n0++] = i;
530: } else perm[n1++] = i;
531: }
532: for (i=n0;i<n;i++) s[i] = -s[0];
533: n1 -= n0;
534: idx0 = 0;
535: idx1 = n0;
536: if (idx1==n) idx1=idx0;
537: for (i=0;i<n;i++) {
538: for (j=0;j<n;j++) R[j*ldr+i] = V[j*ldv+perm[i]];
539: }
540: /* Initialize H */
541: for (i=0;i<n;i++) {
542: PetscArrayzero(V+i*ldv,n);
543: V[perm[i]+i*ldv] = 1.0;
544: }
545: for (i=0;i<n;i++) perm[i] = i;
546: j = 0;
547: while (j<n-k) {
548: if (cmplxEig[j]==0) sz=1;
549: else sz=2;
550: TryHRIt(n,j,sz,V,ldv,R,ldr,s,&exg,&ok,&n0,&n1,&idx0,&idx1,NULL,work+nwu);
551: if (ok) {
552: if (exg) cmplxEig[j] = -cmplxEig[j];
553: j = j+sz;
554: } else { /* to be discarded */
555: k = k+1;
556: if (cmplxEig[j]==0) {
557: if (j<n) {
558: t1 = perm[j];
559: for (i=j;i<n-1;i++) perm[i] = perm[i+1];
560: perm[n-1] = t1;
561: t1 = cmplxEig[j];
562: for (i=j;i<n-1;i++) cmplxEig[i] = cmplxEig[i+1];
563: cmplxEig[n-1] = t1;
564: PetscArraycpy(col1,R+j*ldr,n);
565: for (i=j;i<n-1;i++) {
566: PetscArraycpy(R+i*ldr,R+(i+1)*ldr,n);
567: }
568: PetscArraycpy(R+(n-1)*ldr,col1,n);
569: }
570: } else {
571: k = k+1;
572: if (j<n-1) {
573: t1 = perm[j]; t2 = perm[j+1];
574: for (i=j;i<n-2;i++) perm[i] = perm[i+2];
575: perm[n-2] = t1; perm[n-1] = t2;
576: t1 = cmplxEig[j]; t2 = cmplxEig[j+1];
577: for (i=j;i<n-2;i++) cmplxEig[i] = cmplxEig[i+2];
578: cmplxEig[n-2] = t1; cmplxEig[n-1] = t2;
579: PetscArraycpy(col1,R+j*ldr,n);
580: PetscArraycpy(col2,R+(j+1)*ldr,n);
581: for (i=j;i<n-2;i++) {
582: PetscArraycpy(R+i*ldr,R+(i+2)*ldr,n);
583: }
584: PetscArraycpy(R+(n-2)*ldr,col1,n);
585: PetscArraycpy(R+(n-1)*ldr,col2,n);
586: }
587: }
588: }
589: }
590: if (k!=0) {
591: if (breakdown) *breakdown = PETSC_TRUE;
592: *nv = n-k;
593: }
594: return(0);
595: }
597: PetscErrorCode DSGHIEPOrthogEigenv(DS ds,DSMatType mat,PetscScalar *wr,PetscScalar *wi,PetscBool accum)
598: {
600: PetscInt lws,nwus=0,nwui=0,lwi;
601: PetscInt off,n,nv,ld,i,ldr,l;
602: PetscScalar *W,*X,*R,*ts,zeroS=0.0,oneS=1.0;
603: PetscReal *s,vi,vr,tr,*d,*e;
604: PetscBLASInt ld_,n_,nv_,*perm,*cmplxEig;
607: l = ds->l;
608: n = ds->n-l;
609: PetscBLASIntCast(n,&n_);
610: ld = ds->ld;
611: PetscBLASIntCast(ld,&ld_);
612: off = l*ld+l;
613: s = ds->rmat[DS_MAT_D];
614: if (!ds->compact) {
615: for (i=l;i<ds->n;i++) s[i] = PetscRealPart(*(ds->mat[DS_MAT_B]+i*ld+i));
616: }
617: lws = n*n+7*n;
618: lwi = 2*n;
619: DSAllocateWork_Private(ds,lws,0,lwi);
620: R = ds->work+nwus;
621: nwus += n*n;
622: ldr = n;
623: perm = ds->iwork + nwui;
624: nwui += n;
625: cmplxEig = ds->iwork+nwui;
626: X = ds->mat[mat];
627: for (i=0;i<n;i++) {
628: #if defined(PETSC_USE_COMPLEX)
629: vi = PetscImaginaryPart(wr[l+i]);
630: #else
631: vi = PetscRealPart(wi[l+i]);
632: #endif
633: if (vi!=0) {
634: cmplxEig[i] = 1;
635: cmplxEig[i+1] = 2;
636: i++;
637: } else cmplxEig[i] = 0;
638: }
639: nv = n;
641: /* Perform HR decomposition */
642: /* Hyperbolic rotators */
643: PseudoOrthog_HR(&nv,X+off,ld,s+l,R,ldr,perm,cmplxEig,NULL,ds->work+nwus);
644: /* Sort wr,wi perm */
645: ts = ds->work+nwus;
646: PetscArraycpy(ts,wr+l,n);
647: for (i=0;i<n;i++) wr[i+l] = ts[perm[i]];
648: #if !defined(PETSC_USE_COMPLEX)
649: PetscArraycpy(ts,wi+l,n);
650: for (i=0;i<n;i++) wi[i+l] = ts[perm[i]];
651: #endif
652: /* Projected Matrix */
653: PetscArrayzero(ds->rmat[DS_MAT_T]+2*ld,ld);
654: d = ds->rmat[DS_MAT_T];
655: e = d+ld;
656: d[l+nv-1] = PetscRealPart(wr[l+nv-1]*s[l+nv-1]);
657: for (i=0;i<nv-1;i++) {
658: if (cmplxEig[i]==0) { /* Real */
659: d[l+i] = PetscRealPart(wr[l+i]*s[l+i]);
660: e[l+i] = 0.0;
661: } else {
662: vr = PetscRealPart(wr[l+i]);
663: #if defined(PETSC_USE_COMPLEX)
664: vi = PetscImaginaryPart(wr[l+i]);
665: #else
666: vi = PetscRealPart(wi[l+i]);
667: #endif
668: if (cmplxEig[i]==-1) vi = -vi;
669: tr = PetscRealPart((R[i+(i+1)*ldr]/R[i+i*ldr]))*vi;
670: d[l+i] = (vr-tr)*s[l+i];
671: d[l+i+1] = (vr+tr)*s[l+i+1];
672: e[l+i] = PetscRealPart(s[l+i]*(R[(i+1)+(i+1)*ldr]/R[i+i*ldr])*vi);
673: if (i<nv-2) e[l+i+1] = 0.0;
674: i++;
675: }
676: }
677: /* accumulate previous Q */
678: if (accum) {
679: PetscBLASIntCast(nv,&nv_);
680: DSAllocateMat_Private(ds,DS_MAT_W);
681: W = ds->mat[DS_MAT_W];
682: DSCopyMatrix_Private(ds,DS_MAT_W,DS_MAT_Q);
683: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&nv_,&n_,&oneS,W+off,&ld_,X+off,&ld_,&zeroS,ds->mat[DS_MAT_Q]+off,&ld_));
684: } else {
685: PetscArrayzero(ds->mat[DS_MAT_Q],ld*ld);
686: for (i=0;i<ds->l;i++) *(ds->mat[DS_MAT_Q]+i+i*ld) = 1.0;
687: for (i=0;i<n;i++) { PetscArraycpy(ds->mat[DS_MAT_Q]+off+i*ld,X+off+i*ld,n); }
688: }
689: ds->t = nv+l;
690: if (!ds->compact) { DSSwitchFormat_GHIEP(ds,PETSC_FALSE); }
691: return(0);
692: }
694: /*
695: Reduce to tridiagonal-diagonal pair by means of TridiagDiag_HHR.
696: */
697: PetscErrorCode DSIntermediate_GHIEP(DS ds)
698: {
700: PetscInt i,ld,off;
701: PetscInt nwall,nwallr,nwalli;
702: PetscScalar *A,*B,*Q;
703: PetscReal *d,*e,*s;
706: ld = ds->ld;
707: A = ds->mat[DS_MAT_A];
708: B = ds->mat[DS_MAT_B];
709: Q = ds->mat[DS_MAT_Q];
710: d = ds->rmat[DS_MAT_T];
711: e = ds->rmat[DS_MAT_T]+ld;
712: s = ds->rmat[DS_MAT_D];
713: off = ds->l+ds->l*ld;
714: PetscArrayzero(Q,ld*ld);
715: nwall = ld*ld+ld;
716: nwallr = ld;
717: nwalli = ld;
718: DSAllocateWork_Private(ds,nwall,nwallr,nwalli);
719: for (i=0;i<ds->n;i++) Q[i+i*ld]=1.0;
720: for (i=0;i<ds->n-ds->l;i++) *(ds->perm+i)=i;
721: if (ds->compact) {
722: if (ds->state < DS_STATE_INTERMEDIATE) {
723: DSSwitchFormat_GHIEP(ds,PETSC_FALSE);
724: TridiagDiag_HHR(ds->k-ds->l+1,A+off,ld,s+ds->l,Q+off,ld,PETSC_TRUE,d+ds->l,e+ds->l,ds->perm,ds->work,ds->rwork,ds->iwork);
725: ds->k = ds->l;
726: PetscArrayzero(d+2*ld+ds->l,ds->n-ds->l);
727: }
728: } else {
729: if (ds->state < DS_STATE_INTERMEDIATE) {
730: for (i=0;i<ds->n;i++) s[i] = PetscRealPart(B[i+i*ld]);
731: TridiagDiag_HHR(ds->n-ds->l,A+off,ld,s+ds->l,Q+off,ld,PETSC_FALSE,d+ds->l,e+ds->l,ds->perm,ds->work,ds->rwork,ds->iwork);
732: PetscArrayzero(d+2*ld,ds->n);
733: ds->k = ds->l;
734: DSSwitchFormat_GHIEP(ds,PETSC_FALSE);
735: } else {
736: DSSwitchFormat_GHIEP(ds,PETSC_TRUE);
737: }
738: }
739: return(0);
740: }