Actual source code: dvdimprovex.c
slepc-3.15.2 2021-09-20
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: */
10: /*
11: SLEPc eigensolver: "davidson"
13: Step: improve the eigenvectors X
14: */
16: #include "davidson.h"
17: #include <slepcblaslapack.h>
19: /**** JD update step (I - Kfg'/(g'Kf)) K(A - sB) (I - Kfg'/(g'Kf)) t = (I - Kfg'/(g'Kf))r ****/
21: typedef struct {
22: PetscInt size_X;
23: KSP ksp; /* correction equation solver */
24: Vec friends; /* reference vector for composite vectors */
25: PetscScalar theta[4],thetai[2]; /* the shifts used in the correction eq. */
26: PetscInt maxits; /* maximum number of iterations */
27: PetscInt r_s,r_e; /* the selected eigenpairs to improve */
28: PetscInt ksp_max_size; /* the ksp maximum subvectors size */
29: PetscReal tol; /* the maximum solution tolerance */
30: PetscReal lastTol; /* last tol for dynamic stopping criterion */
31: PetscReal fix; /* tolerance for using the approx. eigenvalue */
32: PetscBool dynamic; /* if the dynamic stopping criterion is applied */
33: dvdDashboard *d; /* the currect dvdDashboard reference */
34: PC old_pc; /* old pc in ksp */
35: BV KZ; /* KZ vecs for the projector KZ*inv(X'*KZ)*X' */
36: BV U; /* new X vectors */
37: PetscScalar *XKZ; /* X'*KZ */
38: PetscScalar *iXKZ; /* inverse of XKZ */
39: PetscInt ldXKZ; /* leading dimension of XKZ */
40: PetscInt size_iXKZ; /* size of iXKZ */
41: PetscInt ldiXKZ; /* leading dimension of iXKZ */
42: PetscInt size_cX; /* last value of d->size_cX */
43: PetscInt old_size_X; /* last number of improved vectors */
44: PetscBLASInt *iXKZPivots; /* array of pivots */
45: } dvdImprovex_jd;
47: /*
48: Compute (I - KZ*iXKZ*X')*V where,
49: V, the vectors to apply the projector,
50: cV, the number of vectors in V,
51: */
52: static PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d,Vec *V,PetscInt cV)
53: {
55: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
56: PetscInt i,ldh,k,l;
57: PetscScalar *h;
58: PetscBLASInt cV_,n,info,ld;
59: #if defined(PETSC_USE_COMPLEX)
60: PetscInt j;
61: #endif
64: if (cV > 2) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
66: /* h <- X'*V */
67: PetscMalloc1(data->size_iXKZ*cV,&h);
68: ldh = data->size_iXKZ;
69: BVGetActiveColumns(data->U,&l,&k);
70: if (ldh!=k) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
71: BVSetActiveColumns(data->U,0,k);
72: for (i=0;i<cV;i++) {
73: BVDotVec(data->U,V[i],&h[ldh*i]);
74: #if defined(PETSC_USE_COMPLEX)
75: for (j=0; j<k; j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
76: #endif
77: }
78: BVSetActiveColumns(data->U,l,k);
80: /* h <- iXKZ\h */
81: PetscBLASIntCast(cV,&cV_);
82: PetscBLASIntCast(data->size_iXKZ,&n);
83: PetscBLASIntCast(data->ldiXKZ,&ld);
85: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
86: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n,&cV_,data->iXKZ,&ld,data->iXKZPivots,h,&n,&info));
87: PetscFPTrapPop();
88: SlepcCheckLapackInfo("getrs",info);
90: /* V <- V - KZ*h */
91: BVSetActiveColumns(data->KZ,0,k);
92: for (i=0;i<cV;i++) {
93: BVMultVec(data->KZ,-1.0,1.0,V[i],&h[ldh*i]);
94: }
95: BVSetActiveColumns(data->KZ,l,k);
96: PetscFree(h);
97: return(0);
98: }
100: /*
101: Compute (I - X*iXKZ*KZ')*V where,
102: V, the vectors to apply the projector,
103: cV, the number of vectors in V,
104: */
105: static PetscErrorCode dvd_improvex_applytrans_proj(dvdDashboard *d,Vec *V,PetscInt cV)
106: {
108: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
109: PetscInt i,ldh,k,l;
110: PetscScalar *h;
111: PetscBLASInt cV_, n, info, ld;
112: #if defined(PETSC_USE_COMPLEX)
113: PetscInt j;
114: #endif
117: if (cV > 2) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
119: /* h <- KZ'*V */
120: PetscMalloc1(data->size_iXKZ*cV,&h);
121: ldh = data->size_iXKZ;
122: BVGetActiveColumns(data->U,&l,&k);
123: if (ldh!=k) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
124: BVSetActiveColumns(data->KZ,0,k);
125: for (i=0;i<cV;i++) {
126: BVDotVec(data->KZ,V[i],&h[ldh*i]);
127: #if defined(PETSC_USE_COMPLEX)
128: for (j=0;j<k;j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
129: #endif
130: }
131: BVSetActiveColumns(data->KZ,l,k);
133: /* h <- iXKZ\h */
134: PetscBLASIntCast(cV,&cV_);
135: PetscBLASIntCast(data->size_iXKZ,&n);
136: PetscBLASIntCast(data->ldiXKZ,&ld);
138: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
139: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n,&cV_,data->iXKZ,&ld,data->iXKZPivots,h,&n,&info));
140: PetscFPTrapPop();
141: SlepcCheckLapackInfo("getrs",info);
143: /* V <- V - U*h */
144: BVSetActiveColumns(data->U,0,k);
145: for (i=0;i<cV;i++) {
146: BVMultVec(data->U,-1.0,1.0,V[i],&h[ldh*i]);
147: }
148: BVSetActiveColumns(data->U,l,k);
149: PetscFree(h);
150: return(0);
151: }
153: static PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d)
154: {
156: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
159: VecDestroy(&data->friends);
161: /* Restore the pc of ksp */
162: if (data->old_pc) {
163: KSPSetPC(data->ksp, data->old_pc);
164: PCDestroy(&data->old_pc);
165: }
166: return(0);
167: }
169: static PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d)
170: {
172: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
175: /* Free local data and objects */
176: PetscFree(data->XKZ);
177: PetscFree(data->iXKZ);
178: PetscFree(data->iXKZPivots);
179: BVDestroy(&data->KZ);
180: BVDestroy(&data->U);
181: PetscFree(data);
182: return(0);
183: }
185: /*
186: y <- theta[1]A*x - theta[0]*B*x
187: auxV, two auxiliary vectors
188: */
189: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmult(dvdImprovex_jd *data,const Vec *x,const Vec *y)
190: {
192: PetscInt n,i;
193: const Vec *Bx;
194: Vec *auxV;
197: n = data->r_e - data->r_s;
198: for (i=0;i<n;i++) {
199: MatMult(data->d->A,x[i],y[i]);
200: }
202: SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
203: for (i=0;i<n;i++) {
204: #if !defined(PETSC_USE_COMPLEX)
205: if (data->d->eigi[data->r_s+i] != 0.0) {
206: if (data->d->B) {
207: MatMult(data->d->B,x[i],auxV[0]);
208: MatMult(data->d->B,x[i+1],auxV[1]);
209: Bx = auxV;
210: } else Bx = &x[i];
212: /* y_i <- [ t_2i+1*A*x_i - t_2i*Bx_i + ti_i*Bx_i+1;
213: y_i+1 t_2i+1*A*x_i+1 - ti_i*Bx_i - t_2i*Bx_i+1 ] */
214: VecAXPBYPCZ(y[i],-data->theta[2*i],data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
215: VecAXPBYPCZ(y[i+1],-data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
216: i++;
217: } else
218: #endif
219: {
220: if (data->d->B) {
221: MatMult(data->d->B,x[i],auxV[0]);
222: Bx = auxV;
223: } else Bx = &x[i];
224: VecAXPBY(y[i],-data->theta[i*2],data->theta[i*2+1],Bx[0]);
225: }
226: }
227: SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
228: return(0);
229: }
231: /*
232: y <- theta[1]'*A'*x - theta[0]'*B'*x
233: */
234: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmulttrans(dvdImprovex_jd *data,const Vec *x,const Vec *y)
235: {
237: PetscInt n,i;
238: const Vec *Bx;
239: Vec *auxV;
242: n = data->r_e - data->r_s;
243: for (i=0;i<n;i++) {
244: MatMultTranspose(data->d->A,x[i],y[i]);
245: }
247: SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
248: for (i=0;i<n;i++) {
249: #if !defined(PETSC_USE_COMPLEX)
250: if (data->d->eigi[data->r_s+i] != 0.0) {
251: if (data->d->B) {
252: MatMultTranspose(data->d->B,x[i],auxV[0]);
253: MatMultTranspose(data->d->B,x[i+1],auxV[1]);
254: Bx = auxV;
255: } else Bx = &x[i];
257: /* y_i <- [ t_2i+1*A*x_i - t_2i*Bx_i - ti_i*Bx_i+1;
258: y_i+1 t_2i+1*A*x_i+1 + ti_i*Bx_i - t_2i*Bx_i+1 ] */
259: VecAXPBYPCZ(y[i],-data->theta[2*i],-data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
260: VecAXPBYPCZ(y[i+1],data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
261: i++;
262: } else
263: #endif
264: {
265: if (data->d->B) {
266: MatMultTranspose(data->d->B,x[i],auxV[0]);
267: Bx = auxV;
268: } else Bx = &x[i];
269: VecAXPBY(y[i],PetscConj(-data->theta[i*2]),PetscConj(data->theta[i*2+1]),Bx[0]);
270: }
271: }
272: SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
273: return(0);
274: }
276: static PetscErrorCode PCApplyBA_dvd(PC pc,PCSide side,Vec in,Vec out,Vec w)
277: {
279: dvdImprovex_jd *data;
280: PetscInt n,i;
281: const Vec *inx,*outx,*wx;
282: Vec *auxV;
283: Mat A;
286: PCGetOperators(pc,&A,NULL);
287: MatShellGetContext(A,(void**)&data);
288: VecCompGetSubVecs(in,NULL,&inx);
289: VecCompGetSubVecs(out,NULL,&outx);
290: VecCompGetSubVecs(w,NULL,&wx);
291: n = data->r_e - data->r_s;
292: SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
293: switch (side) {
294: case PC_LEFT:
295: /* aux <- theta[1]A*in - theta[0]*B*in */
296: dvd_aux_matmult(data,inx,auxV);
298: /* out <- K * aux */
299: for (i=0;i<n;i++) {
300: data->d->improvex_precond(data->d,data->r_s+i,auxV[i],outx[i]);
301: }
302: break;
303: case PC_RIGHT:
304: /* aux <- K * in */
305: for (i=0;i<n;i++) {
306: data->d->improvex_precond(data->d,data->r_s+i,inx[i],auxV[i]);
307: }
309: /* out <- theta[1]A*auxV - theta[0]*B*auxV */
310: dvd_aux_matmult(data,auxV,outx);
311: break;
312: case PC_SYMMETRIC:
313: /* aux <- K^{1/2} * in */
314: for (i=0;i<n;i++) {
315: PCApplySymmetricRight(data->old_pc,inx[i],auxV[i]);
316: }
318: /* wx <- theta[1]A*auxV - theta[0]*B*auxV */
319: dvd_aux_matmult(data,auxV,wx);
321: /* aux <- K^{1/2} * in */
322: for (i=0;i<n;i++) {
323: PCApplySymmetricLeft(data->old_pc,wx[i],outx[i]);
324: }
325: break;
326: default:
327: SETERRQ(PETSC_COMM_SELF,1,"Unsupported KSP side");
328: }
329: /* out <- out - v*(u'*out) */
330: dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
331: SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
332: return(0);
333: }
335: static PetscErrorCode PCApply_dvd(PC pc,Vec in,Vec out)
336: {
338: dvdImprovex_jd *data;
339: PetscInt n,i;
340: const Vec *inx, *outx;
341: Mat A;
344: PCGetOperators(pc,&A,NULL);
345: MatShellGetContext(A,(void**)&data);
346: VecCompGetSubVecs(in,NULL,&inx);
347: VecCompGetSubVecs(out,NULL,&outx);
348: n = data->r_e - data->r_s;
349: /* out <- K * in */
350: for (i=0;i<n;i++) {
351: data->d->improvex_precond(data->d,data->r_s+i,inx[i],outx[i]);
352: }
353: /* out <- out - v*(u'*out) */
354: dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
355: return(0);
356: }
358: static PetscErrorCode PCApplyTranspose_dvd(PC pc,Vec in,Vec out)
359: {
361: dvdImprovex_jd *data;
362: PetscInt n,i;
363: const Vec *inx, *outx;
364: Vec *auxV;
365: Mat A;
368: PCGetOperators(pc,&A,NULL);
369: MatShellGetContext(A,(void**)&data);
370: VecCompGetSubVecs(in,NULL,&inx);
371: VecCompGetSubVecs(out,NULL,&outx);
372: n = data->r_e - data->r_s;
373: SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
374: /* auxV <- in */
375: for (i=0;i<n;i++) {
376: VecCopy(inx[i],auxV[i]);
377: }
378: /* auxV <- auxV - u*(v'*auxV) */
379: dvd_improvex_applytrans_proj(data->d,auxV,n);
380: /* out <- K' * aux */
381: for (i=0;i<n;i++) {
382: PCApplyTranspose(data->old_pc,auxV[i],outx[i]);
383: }
384: SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
385: return(0);
386: }
388: static PetscErrorCode MatMult_dvd_jd(Mat A,Vec in,Vec out)
389: {
391: dvdImprovex_jd *data;
392: PetscInt n;
393: const Vec *inx, *outx;
394: PCSide side;
397: MatShellGetContext(A,(void**)&data);
398: VecCompGetSubVecs(in,NULL,&inx);
399: VecCompGetSubVecs(out,NULL,&outx);
400: n = data->r_e - data->r_s;
401: /* out <- theta[1]A*in - theta[0]*B*in */
402: dvd_aux_matmult(data,inx,outx);
403: KSPGetPCSide(data->ksp,&side);
404: if (side == PC_RIGHT) {
405: /* out <- out - v*(u'*out) */
406: dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
407: }
408: return(0);
409: }
411: static PetscErrorCode MatMultTranspose_dvd_jd(Mat A,Vec in,Vec out)
412: {
414: dvdImprovex_jd *data;
415: PetscInt n,i;
416: const Vec *inx,*outx,*r;
417: Vec *auxV;
418: PCSide side;
421: MatShellGetContext(A,(void**)&data);
422: VecCompGetSubVecs(in,NULL,&inx);
423: VecCompGetSubVecs(out,NULL,&outx);
424: n = data->r_e - data->r_s;
425: KSPGetPCSide(data->ksp,&side);
426: if (side == PC_RIGHT) {
427: /* auxV <- in */
428: SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
429: for (i=0;i<n;i++) {
430: VecCopy(inx[i],auxV[i]);
431: }
432: /* auxV <- auxV - v*(u'*auxV) */
433: dvd_improvex_applytrans_proj(data->d,auxV,n);
434: r = auxV;
435: } else r = inx;
436: /* out <- theta[1]A*r - theta[0]*B*r */
437: dvd_aux_matmulttrans(data,r,outx);
438: if (side == PC_RIGHT) {
439: SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
440: }
441: return(0);
442: }
444: static PetscErrorCode MatCreateVecs_dvd_jd(Mat A,Vec *right,Vec *left)
445: {
447: Vec *r,*l;
448: dvdImprovex_jd *data;
449: PetscInt n,i;
452: MatShellGetContext(A,(void**)&data);
453: n = data->ksp_max_size;
454: if (right) {
455: PetscMalloc1(n,&r);
456: }
457: if (left) {
458: PetscMalloc1(n,&l);
459: }
460: for (i=0;i<n;i++) {
461: MatCreateVecs(data->d->A,right?&r[i]:NULL,left?&l[i]:NULL);
462: }
463: if (right) {
464: VecCreateCompWithVecs(r,n,data->friends,right);
465: for (i=0;i<n;i++) {
466: VecDestroy(&r[i]);
467: }
468: }
469: if (left) {
470: VecCreateCompWithVecs(l,n,data->friends,left);
471: for (i=0;i<n;i++) {
472: VecDestroy(&l[i]);
473: }
474: }
476: if (right) {
477: PetscFree(r);
478: }
479: if (left) {
480: PetscFree(l);
481: }
482: return(0);
483: }
485: static PetscErrorCode dvd_improvex_jd_start(dvdDashboard *d)
486: {
488: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
489: PetscInt rA, cA, rlA, clA;
490: Mat A;
491: PetscBool t;
492: PC pc;
493: Vec v0[2];
496: data->size_cX = data->old_size_X = 0;
497: data->lastTol = data->dynamic?0.5:0.0;
499: /* Setup the ksp */
500: if (data->ksp) {
501: /* Create the reference vector */
502: BVGetColumn(d->eps->V,0,&v0[0]);
503: v0[1] = v0[0];
504: VecCreateCompWithVecs(v0,data->ksp_max_size,NULL,&data->friends);
505: BVRestoreColumn(d->eps->V,0,&v0[0]);
506: PetscLogObjectParent((PetscObject)d->eps,(PetscObject)data->friends);
508: /* Save the current pc and set a PCNONE */
509: KSPGetPC(data->ksp, &data->old_pc);
510: PetscObjectTypeCompare((PetscObject)data->old_pc,PCNONE,&t);
511: data->lastTol = 0.5;
512: if (t) data->old_pc = 0;
513: else {
514: PetscObjectReference((PetscObject)data->old_pc);
515: PCCreate(PetscObjectComm((PetscObject)d->eps),&pc);
516: PCSetType(pc,PCSHELL);
517: PCSetOperators(pc,d->A,d->A);
518: PCSetReusePreconditioner(pc,PETSC_TRUE);
519: PCShellSetApply(pc,PCApply_dvd);
520: PCShellSetApplyBA(pc,PCApplyBA_dvd);
521: PCShellSetApplyTranspose(pc,PCApplyTranspose_dvd);
522: KSPSetPC(data->ksp,pc);
523: PCDestroy(&pc);
524: }
526: /* Create the (I-v*u')*K*(A-s*B) matrix */
527: MatGetSize(d->A,&rA,&cA);
528: MatGetLocalSize(d->A,&rlA,&clA);
529: MatCreateShell(PetscObjectComm((PetscObject)d->A),rlA*data->ksp_max_size,clA*data->ksp_max_size,rA*data->ksp_max_size,cA*data->ksp_max_size,data,&A);
530: MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_dvd_jd);
531: MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_dvd_jd);
532: MatShellSetOperation(A,MATOP_CREATE_VECS,(void(*)(void))MatCreateVecs_dvd_jd);
534: /* Try to avoid KSPReset */
535: KSPGetOperatorsSet(data->ksp,&t,NULL);
536: if (t) {
537: Mat M;
538: PetscInt rM;
539: KSPGetOperators(data->ksp,&M,NULL);
540: MatGetSize(M,&rM,NULL);
541: if (rM != rA*data->ksp_max_size) { KSPReset(data->ksp); }
542: }
543: KSPSetOperators(data->ksp,A,A);
544: KSPSetReusePreconditioner(data->ksp,PETSC_TRUE);
545: KSPSetUp(data->ksp);
546: MatDestroy(&A);
547: } else {
548: data->old_pc = 0;
549: data->friends = NULL;
550: }
551: BVSetActiveColumns(data->KZ,0,0);
552: BVSetActiveColumns(data->U,0,0);
553: return(0);
554: }
556: /*
557: Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
558: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
559: where
560: pX,pY, the right and left eigenvectors of the projected system
561: ld, the leading dimension of pX and pY
562: */
563: static PetscErrorCode dvd_improvex_jd_proj_cuv(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *kr,PetscScalar *theta,PetscScalar *thetai,PetscScalar *pX,PetscScalar *pY,PetscInt ld)
564: {
565: PetscErrorCode ierr;
566: PetscInt n=i_e-i_s,size_KZ,V_new,rm,i,lv,kv,lKZ,kKZ;
567: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
568: const PetscScalar *array;
569: Mat M;
570: Vec u[2],v[2];
571: PetscBLASInt s,ldXKZ,info;
574: /* Check consistency */
575: BVGetActiveColumns(d->eps->V,&lv,&kv);
576: V_new = lv - data->size_cX;
577: if (V_new > data->old_size_X) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
578: data->old_size_X = n;
579: data->size_cX = lv;
581: /* KZ <- KZ(rm:rm+max_cX-1) */
582: BVGetActiveColumns(data->KZ,&lKZ,&kKZ);
583: rm = PetscMax(V_new+lKZ,0);
584: if (rm > 0) {
585: for (i=0;i<lKZ;i++) {
586: BVCopyColumn(data->KZ,i+rm,i);
587: BVCopyColumn(data->U,i+rm,i);
588: }
589: }
591: /* XKZ <- XKZ(rm:rm+max_cX-1,rm:rm+max_cX-1) */
592: if (rm > 0) {
593: for (i=0;i<lKZ;i++) {
594: PetscArraycpy(&data->XKZ[i*data->ldXKZ+i],&data->XKZ[(i+rm)*data->ldXKZ+i+rm],lKZ);
595: }
596: }
597: lKZ = PetscMin(0,lKZ+V_new);
598: BVSetActiveColumns(data->KZ,lKZ,lKZ+n);
599: BVSetActiveColumns(data->U,lKZ,lKZ+n);
601: /* Compute X, KZ and KR */
602: BVGetColumn(data->U,lKZ,u);
603: if (n>1) { BVGetColumn(data->U,lKZ+1,&u[1]); }
604: BVGetColumn(data->KZ,lKZ,v);
605: if (n>1) { BVGetColumn(data->KZ,lKZ+1,&v[1]); }
606: d->improvex_jd_proj_uv(d,i_s,i_e,u,v,kr,theta,thetai,pX,pY,ld);
607: BVRestoreColumn(data->U,lKZ,u);
608: if (n>1) { BVRestoreColumn(data->U,lKZ+1,&u[1]); }
609: BVRestoreColumn(data->KZ,lKZ,v);
610: if (n>1) { BVRestoreColumn(data->KZ,lKZ+1,&v[1]); }
612: /* XKZ <- U'*KZ */
613: MatCreateSeqDense(PETSC_COMM_SELF,lKZ+n,lKZ+n,NULL,&M);
614: BVMatProject(data->KZ,NULL,data->U,M);
615: MatDenseGetArrayRead(M,&array);
616: for (i=lKZ;i<lKZ+n;i++) { /* upper part */
617: PetscArraycpy(&data->XKZ[data->ldXKZ*i],&array[i*(lKZ+n)],lKZ);
618: }
619: for (i=0;i<lKZ+n;i++) { /* lower part */
620: PetscArraycpy(&data->XKZ[data->ldXKZ*i+lKZ],&array[i*(lKZ+n)+lKZ],n);
621: }
622: MatDenseRestoreArrayRead(M,&array);
623: MatDestroy(&M);
625: /* iXKZ <- inv(XKZ) */
626: size_KZ = lKZ+n;
627: PetscBLASIntCast(lKZ+n,&s);
628: data->ldiXKZ = data->size_iXKZ = size_KZ;
629: for (i=0;i<size_KZ;i++) {
630: PetscArraycpy(&data->iXKZ[data->ldiXKZ*i],&data->XKZ[data->ldXKZ*i],size_KZ);
631: }
632: PetscBLASIntCast(data->ldiXKZ,&ldXKZ);
633: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
634: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&s,&s,data->iXKZ,&ldXKZ,data->iXKZPivots,&info));
635: PetscFPTrapPop();
636: SlepcCheckLapackInfo("getrf",info);
637: return(0);
638: }
640: static PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
641: {
642: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
644: PetscInt i,j,n,maxits,maxits0,lits,s,ld,k,max_size_D,lV,kV;
645: PetscScalar *pX,*pY;
646: PetscReal tol,tol0;
647: Vec *kr,kr_comp,D_comp,D[2],kr0[2];
648: PetscBool odd_situation = PETSC_FALSE;
651: BVGetActiveColumns(d->eps->V,&lV,&kV);
652: max_size_D = d->eps->ncv-kV;
653: /* Quick exit */
654: if ((max_size_D == 0) || r_e-r_s <= 0) {
655: *size_D = 0;
656: return(0);
657: }
659: n = PetscMin(PetscMin(data->size_X, max_size_D), r_e-r_s);
660: if (n == 0) SETERRQ(PETSC_COMM_SELF,1,"n == 0");
661: if (data->size_X < r_e-r_s) SETERRQ(PETSC_COMM_SELF,1,"size_X < r_e-r_s");
663: DSGetLeadingDimension(d->eps->ds,&ld);
665: /* Restart lastTol if a new pair converged */
666: if (data->dynamic && data->size_cX < lV)
667: data->lastTol = 0.5;
669: for (i=0,s=0;i<n;i+=s) {
670: /* If the selected eigenvalue is complex, but the arithmetic is real... */
671: #if !defined(PETSC_USE_COMPLEX)
672: if (d->eigi[r_s+i] != 0.0) {
673: if (i+2 <= max_size_D) s=2;
674: else break;
675: } else
676: #endif
677: s=1;
679: data->r_s = r_s+i;
680: data->r_e = r_s+i+s;
681: SlepcVecPoolGetVecs(d->auxV,s,&kr);
683: /* Compute theta, maximum iterations and tolerance */
684: maxits = 0;
685: tol = 1;
686: for (j=0;j<s;j++) {
687: d->improvex_jd_lit(d,r_s+i+j,&data->theta[2*j],&data->thetai[j],&maxits0,&tol0);
688: maxits += maxits0;
689: tol *= tol0;
690: }
691: maxits/= s;
692: tol = data->dynamic?data->lastTol:PetscExpReal(PetscLogReal(tol)/s);
694: /* Compute u, v and kr */
695: k = r_s+i;
696: DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
697: k = r_s+i;
698: DSVectors(d->eps->ds,DS_MAT_Y,&k,NULL);
699: DSGetArray(d->eps->ds,DS_MAT_X,&pX);
700: DSGetArray(d->eps->ds,DS_MAT_Y,&pY);
701: dvd_improvex_jd_proj_cuv(d,r_s+i,r_s+i+s,kr,data->theta,data->thetai,pX,pY,ld);
702: DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);
703: DSRestoreArray(d->eps->ds,DS_MAT_Y,&pY);
705: /* Check if the first eigenpairs are converged */
706: if (i == 0) {
707: PetscInt oldnpreconv = d->npreconv;
708: d->preTestConv(d,0,r_s+s,r_s+s,&d->npreconv);
709: if (d->npreconv > oldnpreconv) break;
710: }
712: /* Test the odd situation of solving Ax=b with A=I */
713: #if !defined(PETSC_USE_COMPLEX)
714: odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && data->thetai[0] == 0. && d->B == NULL)? PETSC_TRUE: PETSC_FALSE;
715: #else
716: odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && d->B == NULL)? PETSC_TRUE: PETSC_FALSE;
717: #endif
718: /* If JD */
719: if (data->ksp && !odd_situation) {
720: /* kr <- -kr */
721: for (j=0;j<s;j++) {
722: VecScale(kr[j],-1.0);
723: }
725: /* Compose kr and D */
726: kr0[0] = kr[0];
727: kr0[1] = (s==2 ? kr[1] : NULL);
728: VecCreateCompWithVecs(kr0,data->ksp_max_size,data->friends,&kr_comp);
729: BVGetColumn(d->eps->V,kV+i,&D[0]);
730: if (s==2) { BVGetColumn(d->eps->V,kV+i+1,&D[1]); }
731: else D[1] = NULL;
732: VecCreateCompWithVecs(D,data->ksp_max_size,data->friends,&D_comp);
733: VecCompSetSubVecs(data->friends,s,NULL);
735: /* Solve the correction equation */
736: KSPSetTolerances(data->ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,maxits);
737: KSPSolve(data->ksp,kr_comp,D_comp);
738: KSPGetIterationNumber(data->ksp,&lits);
740: /* Destroy the composed ks and D */
741: VecDestroy(&kr_comp);
742: VecDestroy(&D_comp);
743: BVRestoreColumn(d->eps->V,kV+i,&D[0]);
744: if (s==2) { BVRestoreColumn(d->eps->V,kV+i+1,&D[1]); }
746: /* If GD */
747: } else {
748: BVGetColumn(d->eps->V,kV+i,&D[0]);
749: if (s==2) { BVGetColumn(d->eps->V,kV+i+1,&D[1]); }
750: for (j=0;j<s;j++) {
751: d->improvex_precond(d,r_s+i+j,kr[j],D[j]);
752: }
753: dvd_improvex_apply_proj(d,D,s);
754: BVRestoreColumn(d->eps->V,kV+i,&D[0]);
755: if (s==2) { BVRestoreColumn(d->eps->V,kV+i+1,&D[1]); }
756: }
757: /* Prevent that short vectors are discarded in the orthogonalization */
758: if (i == 0 && d->eps->errest[d->nconv+r_s] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s] < PETSC_MAX_REAL) {
759: for (j=0;j<s;j++) {
760: BVScaleColumn(d->eps->V,kV+i+j,1.0/d->eps->errest[d->nconv+r_s]);
761: }
762: }
763: SlepcVecPoolRestoreVecs(d->auxV,s,&kr);
764: }
765: *size_D = i;
766: if (data->dynamic) data->lastTol = PetscMax(data->lastTol/2.0,PETSC_MACHINE_EPSILON*10.0);
767: return(0);
768: }
770: PetscErrorCode dvd_improvex_jd(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs,PetscBool dynamic)
771: {
773: dvdImprovex_jd *data;
774: PetscBool useGD;
775: PC pc;
776: PetscInt size_P;
779: /* Setting configuration constrains */
780: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&useGD);
782: /* If the arithmetic is real and the problem is not Hermitian, then
783: the block size is incremented in one */
784: #if !defined(PETSC_USE_COMPLEX)
785: if (!DVD_IS(d->sEP,DVD_EP_HERMITIAN)) {
786: max_bs++;
787: b->max_size_P = PetscMax(b->max_size_P,2);
788: } else
789: #endif
790: {
791: b->max_size_P = PetscMax(b->max_size_P,1);
792: }
793: b->max_size_X = PetscMax(b->max_size_X,max_bs);
794: size_P = b->max_size_P;
796: /* Setup the preconditioner */
797: if (ksp) {
798: KSPGetPC(ksp,&pc);
799: dvd_static_precond_PC(d,b,pc);
800: } else {
801: dvd_static_precond_PC(d,b,0);
802: }
804: /* Setup the step */
805: if (b->state >= DVD_STATE_CONF) {
806: PetscNewLog(d->eps,&data);
807: data->dynamic = dynamic;
808: PetscMalloc1(size_P*size_P,&data->XKZ);
809: PetscMalloc1(size_P*size_P,&data->iXKZ);
810: PetscMalloc1(size_P,&data->iXKZPivots);
811: data->ldXKZ = size_P;
812: data->size_X = b->max_size_X;
813: d->improveX_data = data;
814: data->ksp = useGD? NULL: ksp;
815: data->d = d;
816: d->improveX = dvd_improvex_jd_gen;
817: #if !defined(PETSC_USE_COMPLEX)
818: if (!DVD_IS(d->sEP,DVD_EP_HERMITIAN)) data->ksp_max_size = 2;
819: else
820: #endif
821: data->ksp_max_size = 1;
822: /* Create various vector basis */
823: BVDuplicateResize(d->eps->V,size_P,&data->KZ);
824: BVSetMatrix(data->KZ,NULL,PETSC_FALSE);
825: BVDuplicate(data->KZ,&data->U);
827: EPSDavidsonFLAdd(&d->startList,dvd_improvex_jd_start);
828: EPSDavidsonFLAdd(&d->endList,dvd_improvex_jd_end);
829: EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_jd_d);
830: }
831: return(0);
832: }
834: #if !defined(PETSC_USE_COMPLEX)
835: PETSC_STATIC_INLINE PetscErrorCode dvd_complex_rayleigh_quotient(Vec ur,Vec ui,Vec Axr,Vec Axi,Vec Bxr,Vec Bxi,PetscScalar *eigr,PetscScalar *eigi)
836: {
838: PetscScalar rAr,iAr,rAi,iAi,rBr,iBr,rBi,iBi,b0,b2,b4,b6,b7;
841: /* eigr = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k
842: eigi = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k
843: k = (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr) */
844: VecDotBegin(Axr,ur,&rAr); /* r*A*r */
845: VecDotBegin(Axr,ui,&iAr); /* i*A*r */
846: VecDotBegin(Axi,ur,&rAi); /* r*A*i */
847: VecDotBegin(Axi,ui,&iAi); /* i*A*i */
848: VecDotBegin(Bxr,ur,&rBr); /* r*B*r */
849: VecDotBegin(Bxr,ui,&iBr); /* i*B*r */
850: VecDotBegin(Bxi,ur,&rBi); /* r*B*i */
851: VecDotBegin(Bxi,ui,&iBi); /* i*B*i */
852: VecDotEnd(Axr,ur,&rAr); /* r*A*r */
853: VecDotEnd(Axr,ui,&iAr); /* i*A*r */
854: VecDotEnd(Axi,ur,&rAi); /* r*A*i */
855: VecDotEnd(Axi,ui,&iAi); /* i*A*i */
856: VecDotEnd(Bxr,ur,&rBr); /* r*B*r */
857: VecDotEnd(Bxr,ui,&iBr); /* i*B*r */
858: VecDotEnd(Bxi,ur,&rBi); /* r*B*i */
859: VecDotEnd(Bxi,ui,&iBi); /* i*B*i */
860: b0 = rAr+iAi; /* rAr+iAi */
861: b2 = rAi-iAr; /* rAi-iAr */
862: b4 = rBr+iBi; /* rBr+iBi */
863: b6 = rBi-iBr; /* rBi-iBr */
864: b7 = b4*b4 + b6*b6; /* k */
865: *eigr = (b0*b4 + b2*b6) / b7; /* eig_r */
866: *eigi = (b2*b4 - b0*b6) / b7; /* eig_i */
867: return(0);
868: }
869: #endif
871: PETSC_STATIC_INLINE PetscErrorCode dvd_compute_n_rr(PetscInt i_s,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,Vec *u,Vec *Ax,Vec *Bx)
872: {
874: PetscInt i;
875: PetscScalar b0,b1;
878: for (i=0; i<n; i++) {
879: #if !defined(PETSC_USE_COMPLEX)
880: if (eigi[i_s+i] != 0.0) {
881: PetscScalar eigr0=0.0,eigi0=0.0;
882: dvd_complex_rayleigh_quotient(u[i],u[i+1],Ax[i],Ax[i+1],Bx[i],Bx[i+1],&eigr0,&eigi0);
883: if (PetscAbsScalar(eigr[i_s+i]-eigr0)/PetscAbsScalar(eigr[i_s+i]) > 1e-10 || PetscAbsScalar(eigi[i_s+i]-eigi0)/PetscAbsScalar(eigi[i_s+i]) > 1e-10) {
884: PetscInfo4(u[0],"The eigenvalue %g%+gi is far from its Rayleigh quotient value %g%+gi\n",(double)eigr[i_s+i],(double)eigi[i_s+i],(double)eigr0,(double)eigi0);
885: }
886: i++;
887: } else
888: #endif
889: {
890: VecDotBegin(Ax[i],u[i],&b0);
891: VecDotBegin(Bx[i],u[i],&b1);
892: VecDotEnd(Ax[i],u[i],&b0);
893: VecDotEnd(Bx[i],u[i],&b1);
894: b0 = b0/b1;
895: if (PetscAbsScalar(eigr[i_s+i]-b0)/PetscAbsScalar(eigr[i_s+i]) > 1e-10) {
896: PetscInfo4(u[0],"The eigenvalue %g+%g is far from its Rayleigh quotient value %g+%g\n",(double)PetscRealPart(eigr[i_s+i]),(double)PetscImaginaryPart(eigr[i_s+i]),(double)PetscRealPart(b0),(double)PetscImaginaryPart(b0));
897: }
898: }
899: }
900: return(0);
901: }
903: /*
904: Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
905: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
906: where
907: pX,pY, the right and left eigenvectors of the projected system
908: ld, the leading dimension of pX and pY
909: */
910: PetscErrorCode dvd_improvex_jd_proj_uv_KZX(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *u,Vec *v,Vec *kr,PetscScalar *theta,PetscScalar *thetai,PetscScalar *pX,PetscScalar *pY,PetscInt ld)
911: {
913: PetscInt n = i_e-i_s,i;
914: PetscScalar *b;
915: Vec *Ax,*Bx,*r;
916: Mat M;
917: BV X;
920: BVDuplicateResize(d->eps->V,4,&X);
921: MatCreateSeqDense(PETSC_COMM_SELF,4,4,NULL,&M);
922: /* u <- X(i) */
923: dvd_improvex_compute_X(d,i_s,i_e,u,pX,ld);
925: /* v <- theta[0]A*u + theta[1]*B*u */
927: /* Bx <- B*X(i) */
928: Bx = kr;
929: if (d->BX) {
930: for (i=i_s; i<i_e; ++i) {
931: BVMultVec(d->BX,1.0,0.0,Bx[i-i_s],&pX[ld*i]);
932: }
933: } else {
934: for (i=0;i<n;i++) {
935: if (d->B) {
936: MatMult(d->B, u[i], Bx[i]);
937: } else {
938: VecCopy(u[i], Bx[i]);
939: }
940: }
941: }
943: /* Ax <- A*X(i) */
944: SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&r);
945: Ax = r;
946: for (i=i_s; i<i_e; ++i) {
947: BVMultVec(d->AX,1.0,0.0,Ax[i-i_s],&pX[ld*i]);
948: }
950: /* v <- Y(i) */
951: for (i=i_s; i<i_e; ++i) {
952: BVMultVec(d->W?d->W:d->eps->V,1.0,0.0,v[i-i_s],&pY[ld*i]);
953: }
955: /* Recompute the eigenvalue */
956: dvd_compute_n_rr(i_s,n,d->eigr,d->eigi,v,Ax,Bx);
958: for (i=0;i<n;i++) {
959: #if !defined(PETSC_USE_COMPLEX)
960: if (d->eigi[i_s+i] != 0.0) {
961: /* [r_i r_i+1 kr_i kr_i+1]*= [ theta_2i' 0 1 0
962: 0 theta_2i' 0 1
963: theta_2i+1 -thetai_i -eigr_i -eigi_i
964: thetai_i theta_2i+1 eigi_i -eigr_i ] */
965: MatDenseGetArrayWrite(M,&b);
966: b[0] = b[5] = PetscConj(theta[2*i]);
967: b[2] = b[7] = -theta[2*i+1];
968: b[6] = -(b[3] = thetai[i]);
969: b[1] = b[4] = 0.0;
970: b[8] = b[13] = 1.0/d->nX[i_s+i];
971: b[10] = b[15] = -d->eigr[i_s+i]/d->nX[i_s+i];
972: b[14] = -(b[11] = d->eigi[i_s+i]/d->nX[i_s+i]);
973: b[9] = b[12] = 0.0;
974: MatDenseRestoreArrayWrite(M,&b);
975: BVInsertVec(X,0,Ax[i]);
976: BVInsertVec(X,1,Ax[i+1]);
977: BVInsertVec(X,2,Bx[i]);
978: BVInsertVec(X,3,Bx[i+1]);
979: BVSetActiveColumns(X,0,4);
980: BVMultInPlace(X,M,0,4);
981: BVCopyVec(X,0,Ax[i]);
982: BVCopyVec(X,1,Ax[i+1]);
983: BVCopyVec(X,2,Bx[i]);
984: BVCopyVec(X,3,Bx[i+1]);
985: i++;
986: } else
987: #endif
988: {
989: /* [Ax_i Bx_i]*= [ theta_2i' 1/nX_i
990: theta_2i+1 -eig_i/nX_i ] */
991: MatDenseGetArrayWrite(M,&b);
992: b[0] = PetscConj(theta[i*2]);
993: b[1] = theta[i*2+1];
994: b[4] = 1.0/d->nX[i_s+i];
995: b[5] = -d->eigr[i_s+i]/d->nX[i_s+i];
996: MatDenseRestoreArrayWrite(M,&b);
997: BVInsertVec(X,0,Ax[i]);
998: BVInsertVec(X,1,Bx[i]);
999: BVSetActiveColumns(X,0,2);
1000: BVMultInPlace(X,M,0,2);
1001: BVCopyVec(X,0,Ax[i]);
1002: BVCopyVec(X,1,Bx[i]);
1003: }
1004: }
1005: for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;
1007: /* v <- K^{-1} r = K^{-1}(theta_2i'*Ax + theta_2i+1*Bx) */
1008: for (i=0;i<n;i++) {
1009: d->improvex_precond(d,i_s+i,r[i],v[i]);
1010: }
1011: SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&r);
1013: /* kr <- P*(Ax - eig_i*Bx) */
1014: d->calcpairs_proj_res(d,i_s,i_e,kr);
1015: BVDestroy(&X);
1016: MatDestroy(&M);
1017: return(0);
1018: }
1020: static PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d,PetscInt i,PetscScalar* theta,PetscScalar* thetai,PetscInt *maxits,PetscReal *tol)
1021: {
1022: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
1023: PetscReal a;
1026: a = SlepcAbsEigenvalue(d->eigr[i],d->eigi[i]);
1028: if (d->nR[i] < data->fix*a) {
1029: theta[0] = d->eigr[i];
1030: theta[1] = 1.0;
1031: #if !defined(PETSC_USE_COMPLEX)
1032: *thetai = d->eigi[i];
1033: #endif
1034: } else {
1035: theta[0] = d->target[0];
1036: theta[1] = d->target[1];
1037: #if !defined(PETSC_USE_COMPLEX)
1038: *thetai = 0.0;
1039: #endif
1040: }
1042: #if defined(PETSC_USE_COMPLEX)
1043: if (thetai) *thetai = 0.0;
1044: #endif
1045: *maxits = data->maxits;
1046: *tol = data->tol;
1047: return(0);
1048: }
1050: PetscErrorCode dvd_improvex_jd_lit_const(dvdDashboard *d,dvdBlackboard *b,PetscInt maxits,PetscReal tol,PetscReal fix)
1051: {
1052: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
1055: /* Setup the step */
1056: if (b->state >= DVD_STATE_CONF) {
1057: data->maxits = maxits;
1058: data->tol = tol;
1059: data->fix = fix;
1060: d->improvex_jd_lit = dvd_improvex_jd_lit_const_0;
1061: }
1062: return(0);
1063: }
1065: PetscErrorCode dvd_improvex_jd_proj_uv(dvdDashboard *d,dvdBlackboard *b)
1066: {
1068: /* Setup the step */
1069: if (b->state >= DVD_STATE_CONF) {
1070: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KZX;
1071: }
1072: return(0);
1073: }
1075: PetscErrorCode dvd_improvex_compute_X(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *u_,PetscScalar *pX,PetscInt ld)
1076: {
1078: PetscInt n = i_e - i_s,i;
1079: Vec *u;
1082: if (u_) u = u_;
1083: else if (d->correctXnorm) {
1084: SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&u);
1085: }
1086: if (u_ || d->correctXnorm) {
1087: for (i=0; i<n; i++) {
1088: BVMultVec(d->eps->V,1.0,0.0,u[i],&pX[ld*(i+i_s)]);
1089: }
1090: }
1091: /* nX(i) <- ||X(i)|| */
1092: if (d->correctXnorm) {
1093: for (i=0; i<n; i++) {
1094: VecNormBegin(u[i],NORM_2,&d->nX[i_s+i]);
1095: }
1096: for (i=0; i<n; i++) {
1097: VecNormEnd(u[i],NORM_2,&d->nX[i_s+i]);
1098: }
1099: #if !defined(PETSC_USE_COMPLEX)
1100: for (i=0;i<n;i++) {
1101: if (d->eigi[i_s+i] != 0.0) {
1102: d->nX[i_s+i] = d->nX[i_s+i+1] = PetscSqrtScalar(d->nX[i_s+i]*d->nX[i_s+i]+d->nX[i_s+i+1]*d->nX[i_s+i+1]);
1103: i++;
1104: }
1105: }
1106: #endif
1107: } else {
1108: for (i=0;i<n;i++) d->nX[i_s+i] = 1.0;
1109: }
1110: if (d->correctXnorm && !u_) {
1111: SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&u);
1112: }
1113: return(0);
1114: }