Actual source code: dvd_calcpairs.c
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-2010, Universidad Politecnica de Valencia, Spain
17: This file is part of SLEPc.
18:
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 slepc.h
34: #include private/epsimpl.h
35: #include davidson.h
36: #include slepcblaslapack.h
38: PetscErrorCode dvd_calcpairs_proj_qz(dvdDashboard *d);
39: PetscErrorCode dvd_calcpairs_proj_qz_harm(dvdDashboard *d);
40: PetscErrorCode dvd_calcpairs_updateV(dvdDashboard *d);
41: PetscErrorCode dvd_calcpairs_updateW(dvdDashboard *d);
42: PetscErrorCode dvd_calcpairs_updateAV(dvdDashboard *d);
43: PetscErrorCode dvd_calcpairs_updateBV(dvdDashboard *d);
44: PetscErrorCode dvd_calcpairs_VtAV_gen(dvdDashboard *d, DvdReduction *r,
45: DvdMult_copy_func *sr);
46: PetscErrorCode dvd_calcpairs_VtBV_gen(dvdDashboard *d, DvdReduction *r,
47: DvdMult_copy_func *sr);
48: PetscErrorCode dvd_calcpairs_projeig_qz_std(dvdDashboard *d);
49: PetscErrorCode dvd_calcpairs_projeig_qz_gen(dvdDashboard *d);
50: PetscErrorCode dvd_calcpairs_selectPairs_qz(dvdDashboard *d, PetscInt n);
51: PetscErrorCode dvd_calcpairs_X(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
52: Vec *X);
53: PetscErrorCode dvd_calcpairs_Y(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
54: Vec *Y);
55: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
56: Vec *R, PetscScalar *auxS, Vec auxV);
57: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
58: PetscInt r_e, Vec *R);
59: PetscErrorCode dvd_calcpairs_updateMatV(Mat A, Vec **AV, PetscInt *size_AV,
60: dvdDashboard *d);
61: PetscErrorCode dvd_calcpairs_WtMatV_gen(PetscScalar **H, MatType_t sH,
62: PetscInt ldH, PetscInt *size_H, PetscScalar *MTY, PetscInt ldMTY,
63: PetscScalar *MTX, PetscInt ldMTX, PetscInt rMT, PetscInt cMT, Vec *W,
64: Vec *V, PetscInt size_V, PetscScalar *auxS, DvdReduction *r,
65: DvdMult_copy_func *sr, dvdDashboard *d);
67: /**** Control routines ********************************************************/
70: PetscErrorCode dvd_calcpairs_qz(dvdDashboard *d, dvdBlackboard *b, IP ipI)
71: {
72: PetscTruth std_probl, her_probl;
73: PetscInt i;
77: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
78: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
80: /* Setting configuration constrains */
81: #ifndef PETSC_USE_COMPLEX
82: /* if the last converged eigenvalue is complex its conjugate pair is also
83: converged */
84: b->max_nev = PetscMax(b->max_nev, d->nev+1);
85: #else
86: b->max_nev = PetscMax(b->max_nev, d->nev);
87: #endif
88: b->own_vecs+= b->size_V*(d->B?2:1) /* AV, BV? */;
89: b->own_scalars+= b->max_size_V*b->max_size_V*2*(std_probl?1:2);
90: /* H, G?, S, T? */
91: b->own_scalars+= b->max_size_V*b->max_size_V*(std_probl?1:2);
92: /* pX, pY? */
93: b->own_scalars+= b->max_nev*b->max_nev*(std_probl?1:2); /* cS, cT? */
94: b->max_size_auxS = PetscMax(PetscMax(
95: b->max_size_auxS,
96: b->max_size_V*b->max_size_V*4
97: /* SlepcReduction */ ),
98: std_probl?0:(b->max_size_V*11+16) /* projeig */);
100: /* Setup the step */
101: if (b->state >= DVD_STATE_CONF) {
102: d->size_AV = 0;
103: d->real_AV = d->AV = b->free_vecs; b->free_vecs+= b->size_V;
104: d->max_size_AV = b->size_V;
105: d->size_H = 0;
106: d->H = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
107: d->real_H = d->H;
108: d->ldH = b->max_size_V;
109: d->pX = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
110: d->S = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
111: d->cS = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
112: for (i=0; i<b->max_nev*b->max_nev; i++) d->cS[i] = 0.0;
113: d->ldcS = b->max_nev;
114: d->ipV = ipI;
115: d->ipW = ipI;
116: d->size_cX = d->size_cY = 0;
117: d->cY = PETSC_NULL;
118: d->pY = PETSC_NULL;
119: d->T = PETSC_NULL;
120: d->ldcT = PETSC_NULL;
121: d->cT = 0;
122: if (d->B) {
123: d->size_BV = 0;
124: d->real_BV = d->BV = b->free_vecs; b->free_vecs+= b->size_V;
125: d->max_size_BV = b->size_V;
126: }
127: if (!std_probl) {
128: d->size_G = 0;
129: d->G = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
130: d->real_G = d->G;
131: d->T = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
132: d->cT = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
133: for (i=0; i<b->max_nev*b->max_nev; i++) d->cT[i] = 0.0;
134: d->ldcT = b->max_nev;
135: d->pY = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
136: /* If the problem is GHEP without B-orthonormalization, active BcX */
137: if(her_probl) d->BcX = d->AV;
139: /* Else, active the left and right converged invariant subspaces */
140: else d->cY = d->AV;
141: }
143: d->calcPairs = d->W?dvd_calcpairs_proj_qz_harm:dvd_calcpairs_proj_qz;
144: d->calcpairs_residual = dvd_calcpairs_res_0;
145: d->calcpairs_proj_res = dvd_calcpairs_proj_res;
146: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_qz;
147: d->calcpairs_X = dvd_calcpairs_X;
148: d->calcpairs_Y = dvd_calcpairs_Y;
149: d->ipI = ipI;
150: }
152: return(0);
153: }
158: PetscErrorCode dvd_calcpairs_proj_qz(dvdDashboard *d)
159: {
160: PetscErrorCode ierr;
161: DvdReduction r;
162: DvdReductionChunk
163: ops[2];
164: DvdMult_copy_func
165: sr[2];
166: PetscInt size_in = 2*d->size_V*d->size_V;
167: PetscScalar *in = d->auxS, *out = in+size_in;
171: /* Prepare reductions */
172: SlepcAllReduceSumBegin(ops, 2, in, out, size_in, &r,
173: ((PetscObject)d->V[0])->comm);
175: /* Update AV, BV and the projected matrices */
176: dvd_calcpairs_updateV(d);
177: dvd_calcpairs_updateAV(d);
178: dvd_calcpairs_VtAV_gen(d, &r, &sr[0]);
179: if (d->BV) { dvd_calcpairs_updateBV(d); }
180: if (DVD_ISNOT(d->sEP, DVD_EP_STD)) {
181: dvd_calcpairs_VtBV_gen(d, &r, &sr[1]);
182: }
184: /* Do reductions */
185: SlepcAllReduceSumEnd(&r);
187: if (d->MT_type != DVD_MT_IDENTITY) {
188: d->MT_type = DVD_MT_IDENTITY;
189: // d->pX_type|= DVD_MAT_IDENTITY;
190: d->V_tra_s = d->V_tra_e = 0;
191: }
192: d->pX_type = 0;
193: //TODO: uncomment this condition
194: // if(d->V_new_e - d->V_new_s > 0) {
195: if (DVD_IS(d->sEP, DVD_EP_STD)) {
196: dvd_calcpairs_projeig_qz_std(d);
197: } else {
198: dvd_calcpairs_projeig_qz_gen(d);
199: }
200: // }
201: d->V_new_s = d->V_new_e;
203: /* Check consistency */
204: if ((d->size_V != d->V_new_e) || (d->size_V != d->size_H) ||
205: (d->size_V != d->size_AV) || (DVD_ISNOT(d->sEP, DVD_EP_STD) && (
206: (d->size_V != d->size_G) || (d->size_V != d->size_BV) ))) {
207: SETERRQ(1, "Consistency broken!");
208: }
210: return(0);
211: }
215: PetscErrorCode dvd_calcpairs_proj_qz_harm(dvdDashboard *d)
216: {
217: PetscErrorCode ierr;
218: DvdReduction r;
219: DvdReductionChunk
220: ops[2];
221: DvdMult_copy_func
222: sr[2];
223: PetscInt size_in = 2*d->size_V*d->size_V;
224: PetscScalar *in = d->auxS, *out = in+size_in;
228: /* Prepare reductions */
229: SlepcAllReduceSumBegin(ops, 2, in, out, size_in, &r,
230: ((PetscObject)d->V[0])->comm);
232: /* Update AV, BV and the projected matrices */
233: dvd_calcpairs_updateV(d);
234: dvd_calcpairs_updateAV(d);
235: if (d->BV) { dvd_calcpairs_updateBV(d); }
236: dvd_calcpairs_updateW(d);
237: dvd_calcpairs_VtAV_gen(d, &r, &sr[0]);
238: if (DVD_ISNOT(d->sEP, DVD_EP_STD)) {
239: dvd_calcpairs_VtBV_gen(d, &r, &sr[1]);
240: }
242: /* Do reductions */
243: SlepcAllReduceSumEnd(&r);
245: /* Perform the transformation on the projected problem */
246: d->calcpairs_proj_trans(d);
248: if (d->MT_type != DVD_MT_IDENTITY) {
249: d->MT_type = DVD_MT_IDENTITY;
250: // d->pX_type|= DVD_MAT_IDENTITY;
251: d->V_tra_s = d->V_tra_e = 0;
252: }
253: d->pX_type = 0;
254: //TODO: uncomment this condition
255: // if(d->V_new_e - d->V_new_s > 0) {
256: if (DVD_IS(d->sEP, DVD_EP_STD)) {
257: dvd_calcpairs_projeig_qz_std(d);
258: } else {
259: dvd_calcpairs_projeig_qz_gen(d);
260: }
261: // }
262: d->V_new_s = d->V_new_e;
264: return(0);
265: }
267: /**** Basic routines **********************************************************/
271: PetscErrorCode dvd_calcpairs_updateV(dvdDashboard *d)
272: {
273: PetscErrorCode ierr;
274: Vec *cX = d->BcX? d->BcX : ( (d->cY && !d->W)? d->cY : d->cX );
278: /* V <- gs([cX f.V(0:f.V_new_s-1)], f.V(V_new_s:V_new_e-1)) */
279: dvd_orthV(d->ipV, d->eps->DS, d->eps->nds, cX, d->size_cX, d->V,
280: d->V_new_s, d->V_new_e, d->auxS, d->auxV[0], d->eps->rand);
281:
283: return(0);
284: }
288: PetscErrorCode dvd_calcpairs_updateW(dvdDashboard *d)
289: {
290: PetscErrorCode ierr;
294: /* Update W */
295: d->calcpairs_W(d);
297: /* W <- gs([cY f.W(0:f.V_new_s-1)], f.W(V_new_s:V_new_e-1)) */
298: dvd_orthV(d->ipW, PETSC_NULL, 0, d->cY, d->size_cY, d->W, d->V_new_s,
299: d->V_new_e, d->auxS, d->auxV[0], d->eps->rand);
300:
302: return(0);
303: }
308: PetscErrorCode dvd_calcpairs_updateAV(dvdDashboard *d)
309: {
310: PetscErrorCode ierr;
314: /* f.AV(f.V_tra) = f.AV * f.MT; f.AV(f.V_new) = A*f.V(f.V_new) */
315: dvd_calcpairs_updateMatV(d->A, &d->AV, &d->size_AV, d);
317: return(0);
318: }
322: PetscErrorCode dvd_calcpairs_updateBV(dvdDashboard *d)
323: {
324: PetscErrorCode ierr;
328: /* f.BV(f.V_tra) = f.BV * f.MT; f.BV(f.V_new) = B*f.V(f.V_new) */
329: dvd_calcpairs_updateMatV(d->B, &d->BV, &d->size_BV, d);
331: return(0);
332: }
336: PetscErrorCode dvd_calcpairs_VtAV_gen(dvdDashboard *d, DvdReduction *r,
337: DvdMult_copy_func *sr)
338: {
339: PetscInt ldMTY = d->MTY?d->ldMTY:d->ldMTX;
340: /* WARNING: auxS uses space assigned to r */
341: PetscScalar *auxS = r->out,
342: *MTY = d->MTY?d->MTY:d->MTX;
343: Vec *W = d->W?d->W:d->V;
344: PetscErrorCode ierr;
348: /* f.H = [f.H(f.V_imm,f.V_imm) f.V(f.V_imm)'*f.AV(f.V_new);
349: f.V(f.V_new)'*f.AV(f.V_imm) f.V(f.V_new)'*f.AV(f.V_new) ] */
350: if (DVD_IS(d->sA,DVD_MAT_HERMITIAN))
351: d->sH = DVD_MAT_HERMITIAN | DVD_MAT_IMPLICIT | DVD_MAT_UTRIANG;
352: if ((d->V_imm_e - d->V_imm_s == 0) && (d->V_tra_e - d->V_tra_s == 0))
353: d->size_H = 0;
354: dvd_calcpairs_WtMatV_gen(&d->H, d->sH, d->ldH, &d->size_H,
355: &MTY[ldMTY*d->V_tra_s], ldMTY,
356: &d->MTX[d->ldMTX*d->V_tra_s], d->ldMTX,
357: d->size_MT, d->V_tra_e-d->V_tra_s,
358: W, d->AV, d->size_V,
359: auxS, r, sr, d);
361: return(0);
362: }
367: PetscErrorCode dvd_calcpairs_VtBV_gen(dvdDashboard *d, DvdReduction *r,
368: DvdMult_copy_func *sr)
369: {
370: PetscErrorCode ierr;
371: PetscInt ldMTY = d->MTY?d->ldMTY:d->ldMTX;
372: /* WARNING: auxS uses space assigned to r */
373: PetscScalar *auxS = r->out,
374: *MTY = d->MTY?d->MTY:d->MTX;
375: Vec *W = d->W?d->W:d->V;
379: /* f.G = [f.G(f.V_imm,f.V_imm) f.V(f.V_imm)'*f.BV(f.V_new);
380: f.V(f.V_new)'*f.BV(f.V_imm) f.V(f.V_new)'*f.BV(f.V_new) ] */
381: if (DVD_IS(d->sB,DVD_MAT_HERMITIAN))
382: d->sG = DVD_MAT_HERMITIAN | DVD_MAT_IMPLICIT | DVD_MAT_UTRIANG;
383: if ((d->V_imm_e - d->V_imm_s == 0) && (d->V_tra_e - d->V_tra_s == 0))
384: d->size_G = 0;
385: dvd_calcpairs_WtMatV_gen(&d->G, d->sG, d->ldH, &d->size_G,
386: &MTY[ldMTY*d->V_tra_s], ldMTY,
387: &d->MTX[d->ldMTX*d->V_tra_s], d->ldMTX,
388: d->size_MT, d->V_tra_e-d->V_tra_s,
389: W, d->BV?d->BV:d->V, d->size_V,
390: auxS, r, sr, d);
392: return(0);
393: }
398: PetscErrorCode dvd_calcpairs_projeig_qz_std(dvdDashboard *d)
399: {
400: PetscErrorCode ierr;
404: /* S <- H */
405: d->ldS = d->ldpX = d->size_H;
406: SlepcDenseCopyTriang(d->S, 0, d->size_H, d->H, d->sH, d->ldH,
407: d->size_H, d->size_H);
409: /* S = pX' * H * pX */
410: EPSDenseHessenberg(d->size_H, 0, d->S, d->ldS, d->pX);
411: EPSDenseSchur(d->size_H, 0, d->S, d->ldS, d->pX, d->eigr, d->eigi);
412:
414: d->pX_type = (d->pX_type & !DVD_MAT_IDENTITY) | DVD_MAT_UNITARY;
416: return(0);
417: }
419: /*
420: auxS(dgges) = size_H (beta) + 8*size_H+16 (work)
421: auxS(zgges) = size_H (beta) + 1+2*size_H (work) + 8*size_H (rwork)
422: */
425: PetscErrorCode dvd_calcpairs_projeig_qz_gen(dvdDashboard *d)
426: {
427: PetscErrorCode ierr;
428: PetscScalar *beta = d->auxS;
429: #if !defined(PETSC_USE_COMPLEX)
430: PetscScalar *auxS = beta + d->size_H;
431: PetscBLASInt n_auxS = d->size_auxS - d->size_H;
432: #else
433: PetscReal *auxR = (PetscReal*)(beta + d->size_H);
434: PetscScalar *auxS = (PetscScalar*)(auxR+8*d->size_H);
435: PetscBLASInt n_auxS = d->size_auxS - 9*d->size_H;
436: #endif
437: PetscInt i;
438: PetscBLASInt info,n,a;
442: /* S <- H, T <- G */
443: d->ldS = d->ldT = d->ldpX = d->ldpY = d->size_H;
444: SlepcDenseCopyTriang(d->S, 0, d->size_H, d->H, d->sH, d->ldH,
445: d->size_H, d->size_H);
446: SlepcDenseCopyTriang(d->T, 0, d->size_H, d->G, d->sG, d->ldH,
447: d->size_H, d->size_H);
449: /* S = Z'*H*Q, T = Z'*G*Q */
450: n = d->size_H;
451: #if !defined(PETSC_USE_COMPLEX)
452: LAPACKgges_(d->pY?"V":"N", "V", "N", PETSC_NULL, &n, d->S, &n, d->T, &n,
453: &a, d->eigr, d->eigi, beta, d->pY, &n, d->pX, &n,
454: auxS, &n_auxS, PETSC_NULL, &info);
455: #else
456: LAPACKgges_(d->pY?"V":"N", "V", "N", PETSC_NULL, &n, d->S, &n, d->T, &n,
457: &a, d->eigr, beta, d->pY, &n, d->pX, &n,
458: auxS, &n_auxS, auxR, PETSC_NULL, &info);
459: #endif
460: if (info) SETERRQ1(PETSC_ERR_LIB, "Error in Lapack GGES %d", info);
462: /* eigr[i] <- eigr[i] / beta[i] */
463: for (i=0; i<d->size_H; i++)
464: d->eigr[i] /= beta[i],
465: d->eigi[i] /= beta[i];
467: d->pX_type = (d->pX_type & !DVD_MAT_IDENTITY) | DVD_MAT_UNITARY;
468: d->pY_type = (d->pY_type & !DVD_MAT_IDENTITY) | DVD_MAT_UNITARY;
470: return(0);
471: }
476: PetscErrorCode dvd_calcpairs_selectPairs_qz(dvdDashboard *d, PetscInt n)
477: {
478: PetscErrorCode ierr;
479: #if defined(PETSC_USE_COMPLEX)
480: PetscScalar s;
481: PetscInt i, j;
482: #endif
485: if ((d->ldpX != d->size_H) ||
486: ( d->T &&
487: ((d->ldS != d->ldT) || (d->ldpX != d->ldpY) ||
488: (d->ldpX != d->size_H)) ) ) {
489: SETERRQ(1, "Error before ordering eigenpairs!");
490: }
492: if (d->T) {
493: EPSSortDenseSchurGeneralized(d->eps, d->size_H, 0, n, d->S, d->T,
494: d->ldS, d->pY, d->pX, d->eigr,
495: d->eigi);
496: } else {
497: EPSSortDenseSchur(d->eps, d->size_H, 0, d->S, d->ldS, d->pX,
498: d->eigr, d->eigi);
499: }
501: if (d->calcpairs_eigs_trans) {
502: d->calcpairs_eigs_trans(d);
503: }
505: /* Some functions need the diagonal elements in T be real */
506: #if defined(PETSC_USE_COMPLEX)
507: if (d->T) for(i=0; i<d->size_H; i++)
508: if (PetscImaginaryPart(d->T[d->ldT*i+i]) != 0.0) {
509: s = PetscConj(d->T[d->ldT*i+i])/PetscAbsScalar(d->T[d->ldT*i+i]);
510: for(j=0; j<=i; j++)
511: d->T[d->ldT*i+j] = PetscRealPart(d->T[d->ldT*i+j]*s),
512: d->S[d->ldS*i+j]*= s;
513: for(j=0; j<d->size_H; j++) d->pX[d->ldpX*i+j]*= s;
514: }
515: #endif
517: return(0);
518: }
523: PetscErrorCode dvd_calcpairs_X(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
524: Vec *X)
525: {
526: PetscInt i;
527: PetscErrorCode ierr;
531: /* X = V * U(0:n-1) */
532: if (DVD_IS(d->pX_type,DVD_MAT_IDENTITY)) {
533: if (d->V != X) for(i=r_s; i<r_e; i++) {
534: VecCopy(d->V[i], X[i]);
535: }
536: } else {
537: SlepcUpdateVectorsZ(X, 0.0, 1.0, d->V, d->size_H, &d->pX[d->ldpX*r_s],
538: d->ldpX, d->size_H, r_e-r_s);
539: }
541: /* nX[i] <- ||X[i]|| */
542: if (d->correctXnorm) for(i=0; i<r_e-r_s; i++) {
543: VecNorm(X[i], NORM_2, &d->nX[r_s+i]);
544: } else for(i=0; i<r_e-r_s; i++) {
545: d->nX[r_s+i] = 1.0;
546: }
548: return(0);
549: }
554: PetscErrorCode dvd_calcpairs_Y(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
555: Vec *Y)
556: {
557: PetscInt i, ldpX = d->pY?d->ldpY:d->ldpX;
558: PetscErrorCode ierr;
559: Vec *V = d->W?d->W:d->V;
560: PetscScalar *pX = d->pY?d->pY:d->pX;
564: /* Y = V * pX(0:n-1) */
565: if (DVD_IS(d->pX_type,DVD_MAT_IDENTITY)) {
566: if (V != Y) for(i=r_s; i<r_e; i++) {
567: VecCopy(V[i], Y[i]);
568: }
569: } else {
570: SlepcUpdateVectorsZ(Y, 0.0, 1.0, V, d->size_H, &pX[ldpX*r_s], ldpX,
571: d->size_H, r_e-r_s);
572: }
574: return(0);
575: }
577: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
578: the norm, where
579: i <- r_s..r_e,
580: UL, auxiliar scalar matrix of size size_H*(r_e-r_s),
581: auxV, auxiliar global vector.
582: */
585: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
586: Vec *R, PetscScalar *UL, Vec auxV)
587: {
588: PetscInt i, j;
589: PetscErrorCode ierr;
593: /* If the eigenproblem is not reduced to standard */
594: if ((d->B == PETSC_NULL) || DVD_ISNOT(d->sEP, DVD_EP_STD)) {
595: /* UL = f.U(0:n-1) * diag(f.pL(0:n-1)) */
596: for(i=r_s; i<r_e; i++) for(j=0; j<d->size_H; j++)
597: UL[d->size_H*(i-r_s)+j] = d->pX[d->ldpX*i+j]*d->eigr[i];
599: if (d->B == PETSC_NULL) {
600: /* R <- V * UL */
601: SlepcUpdateVectorsZ(R, 0.0, 1.0, d->V, d->size_V, UL, d->size_H,
602: d->size_H, r_e-r_s);
603: } else {
604: /* R <- BV * UL */
605: SlepcUpdateVectorsZ(R, 0.0, 1.0, d->BV, d->size_BV, UL,
606: d->size_H, d->size_H, r_e-r_s);
607:
608: }
609: /* R <- AV*U - R */
610: SlepcUpdateVectorsZ(R, -1.0, 1.0, d->AV, d->size_AV,
611: &d->pX[d->ldpX*r_s], d->ldpX, d->size_H, r_e-r_s);
612:
614: /* If the problem was reduced to standard, R[i] = B*X[i] */
615: } else {
616: /* R[i] <- R[i] * eigr[i] */
617: for(i=r_s; i<r_e; i++) {
618: VecScale(R[i-r_s], d->eigr[i]);
619: }
620:
621: /* R <- AV*U - R */
622: SlepcUpdateVectorsZ(R, -1.0, 1.0, d->AV, d->size_AV,
623: &d->pX[d->ldpX*r_s], d->ldpX, d->size_H, r_e-r_s);
624:
625: }
627: d->calcpairs_proj_res(d, r_s, r_e, R);
629: return(0);
630: }
634: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
635: PetscInt r_e, Vec *R)
636: {
637: PetscInt i;
638: PetscErrorCode ierr;
639: PetscTruth lindep;
640: Vec *cX;
644: /* If exists the BcX, R <- orth(BcX, R), nR[i] <- ||R[i]|| */
645: if (d->BcX)
646: cX = d->BcX;
648: /* If exists left subspace, R <- orth(cY, R), nR[i] <- ||R[i]|| */
649: else if (d->cY)
650: cX = d->cY;
652: /* If fany configurations, R <- orth(cX, R), nR[i] <- ||R[i]|| */
653: else if (!(DVD_IS(d->sEP, DVD_EP_STD) && DVD_IS(d->sEP, DVD_EP_HERMITIAN)))
654: cX = d->cX;
656: /* Otherwise, nR[i] <- ||R[i]|| */
657: else
658: cX = PETSC_NULL;
660: if (cX) for (i=0; i<r_e-r_s; i++) {
661: IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_cX, PETSC_NULL,
662: cX, R[i], PETSC_NULL, &d->nR[r_s+i], &lindep);
663:
664: if(lindep || (d->nR[r_s+i] < PETSC_MACHINE_EPSILON)) {
665: SETERRQ(1, "Error during the residual computation of the eigenvectors!");
666: }
668: } else for(i=0; i<r_e-r_s; i++) {
669: VecNorm(R[i], NORM_2, &d->nR[r_s+i]);
670: }
672: return(0);
673: }
675: /**** Patterns implementation *************************************************/
678: PetscErrorCode dvd_calcpairs_updateMatV(Mat A, Vec **AV, PetscInt *size_AV,
679: dvdDashboard *d)
680: {
681: PetscInt i;
682: PetscErrorCode ierr;
686: /* f.AV((0:f.V_tra.size)+f.imm.s) = f.AV * f.U(f.V_tra) */
687: if (d->MT_type == DVD_MT_pX) {
688: SlepcUpdateVectorsZ(*AV+d->V_imm_e, 0.0, 1.0, *AV, *size_AV,
689: &d->pX[d->ldpX*d->V_tra_s], d->ldpX,
690: *size_AV, d->V_tra_e-d->V_tra_s);
691: } else if (d->MT_type == DVD_MT_ORTHO) {
692: SlepcUpdateVectorsZ(*AV+d->V_imm_e, 0.0, 1.0, *AV, *size_AV,
693: &d->MTX[d->ldMTX*d->V_tra_s], d->ldMTX,
694: *size_AV, d->V_tra_e-d->V_tra_s);
695: }
696: *AV = *AV+d->V_imm_s;
698: /* f.AV(f.V_new) = A*f.V(f.V_new) */
699: if (d->V_imm_e-d->V_imm_s + d->V_tra_e-d->V_tra_s != d->V_new_s) {
700: SETERRQ(1, "d->V_imm_e-d->V_imm_s + d->V_tra_e-d->V_tra_s != d->V_new_s !");
701: }
703: for (i=d->V_new_s; i<d->V_new_e; i++) {
704: MatMult(A, d->V[i], (*AV)[i]);
705: }
706: *size_AV = d->V_new_e;
708: return(0);
709: }
711: /*
712: Compute f.H = [MTY'*H*MTX W(tra)'*V(new);
713: W(new)'*V(tra) W(new)'*V(new) ]
714: where
715: tra = 0:cMT-1,
716: new = cMT:size_V-1,
717: ldH, the leading dimension of H,
718: auxS, auxiliary scalar vector of size ldH*max(tra,size_V),
719: */
722: PetscErrorCode dvd_calcpairs_WtMatV_gen(PetscScalar **H, MatType_t sH,
723: PetscInt ldH, PetscInt *size_H, PetscScalar *MTY, PetscInt ldMTY,
724: PetscScalar *MTX, PetscInt ldMTX, PetscInt rMT, PetscInt cMT, Vec *W,
725: Vec *V, PetscInt size_V, PetscScalar *auxS, DvdReduction *r,
726: DvdMult_copy_func *sr, dvdDashboard *d)
727: {
728: PetscErrorCode ierr;
732: /* H <- MTY^T * (H * MTX) */
733: if (cMT > 0) {
734: SlepcDenseMatProdTriang(auxS, 0, ldH,
735: *H, sH, ldH, *size_H, *size_H, PETSC_FALSE,
736: MTX, 0, ldMTX, rMT, cMT, PETSC_FALSE);
737:
738: SlepcDenseMatProdTriang(*H, sH, ldH,
739: MTY, 0, ldMTY, rMT, cMT, PETSC_TRUE,
740: auxS, 0, ldH, *size_H, cMT, PETSC_FALSE);
741:
742: *size_H = cMT;
743: }
745: /* H = [H W(tra)'*W(new);
746: W(new)'*V(tra) W(new)'*V(new) ] */
747: VecsMultS(*H, sH, ldH, W, *size_H, size_V, V, *size_H, size_V, r, sr);
748:
749: *size_H = size_V;
751: return(0);
752: }