Actual source code: dvdgd2.c
slepc-3.7.1 2016-05-27
1: /*
2: SLEPc eigensolver: "davidson"
4: Step: improve the eigenvectors X with GD2
6: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7: SLEPc - Scalable Library for Eigenvalue Problem Computations
8: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
10: This file is part of SLEPc.
12: SLEPc is free software: you can redistribute it and/or modify it under the
13: terms of version 3 of the GNU Lesser General Public License as published by
14: the Free Software Foundation.
16: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
17: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
19: more details.
21: You should have received a copy of the GNU Lesser General Public License
22: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
23: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24: */
26: #include davidson.h
28: typedef struct {
29: PetscInt size_X;
30: } dvdImprovex_gd2;
34: static PetscErrorCode dvd_improvex_gd2_d(dvdDashboard *d)
35: {
36: PetscErrorCode ierr;
37: dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;
40: /* Free local data and objects */
41: PetscFree(data);
42: return(0);
43: }
47: static PetscErrorCode dvd_improvex_gd2_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
48: {
49: dvdImprovex_gd2 *data = (dvdImprovex_gd2*)d->improveX_data;
50: PetscErrorCode ierr;
51: PetscInt i,j,n,s,ld,k,lv,kv,max_size_D;
52: PetscScalar *pX,*b;
53: Vec *Ax,*Bx,v,*x;
54: Mat M;
55: BV X;
58: /* Compute the number of pairs to improve */
59: BVGetActiveColumns(d->eps->V,&lv,&kv);
60: max_size_D = d->eps->ncv-kv;
61: n = PetscMin(PetscMin(data->size_X*2,max_size_D),(r_e-r_s)*2)/2;
62: #if !defined(PETSC_USE_COMPLEX)
63: /* If the last eigenvalue is a complex conjugate pair, n is increased by one */
64: for (i=0; i<n; i++) {
65: if (d->eigi[i] != 0.0) i++;
66: }
67: if (i > n) {
68: n = PetscMin(PetscMin(data->size_X*2,max_size_D),(n+1)*2)/2;
69: if (i > n) n--;
70: }
71: #endif
73: /* Quick exit */
74: if (max_size_D == 0 || r_e-r_s <= 0 || n == 0) {
75: *size_D = 0;
76: return(0);
77: }
79: BVDuplicateResize(d->eps->V,4,&X);
80: MatCreateSeqDense(PETSC_COMM_SELF,4,2,NULL,&M);
82: /* Compute the eigenvectors of the selected pairs */
83: for (i=0;i<n;) {
84: k = r_s+i;
85: DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
86: i = k+1; /* skip complex conjugate pairs */
87: }
88: DSGetArray(d->eps->ds,DS_MAT_X,&pX);
89: DSGetLeadingDimension(d->eps->ds,&ld);
91: SlepcVecPoolGetVecs(d->auxV,n,&Ax);
92: SlepcVecPoolGetVecs(d->auxV,n,&Bx);
94: /* Bx <- B*X(i) */
95: if (d->BX) {
96: /* Compute the norms of the eigenvectors */
97: if (d->correctXnorm) {
98: dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
99: } else {
100: for (i=0;i<n;i++) d->nX[r_s+i] = 1.0;
101: }
102: for (i=0;i<n;i++) {
103: BVMultVec(d->BX,1.0,0.0,Bx[i],&pX[ld*(r_s+i)]);
104: }
105: } else if (d->B) {
106: SlepcVecPoolGetVecs(d->auxV,1,&x);
107: for (i=0;i<n;i++) {
108: /* auxV(0) <- X(i) */
109: dvd_improvex_compute_X(d,r_s+i,r_s+i+1,x,pX,ld);
110: /* Bx(i) <- B*auxV(0) */
111: MatMult(d->B,x[0],Bx[i]);
112: }
113: SlepcVecPoolRestoreVecs(d->auxV,1,&x);
114: } else {
115: /* Bx <- X */
116: dvd_improvex_compute_X(d,r_s,r_s+n,Bx,pX,ld);
117: }
119: /* Ax <- A*X(i) */
120: for (i=0;i<n;i++) {
121: BVMultVec(d->AX,1.0,0.0,Ax[i],&pX[ld*(i+r_s)]);
122: }
124: DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);
126: for (i=0,s=0;i<n;i+=s) {
127: #if !defined(PETSC_USE_COMPLEX)
128: if (d->eigi[r_s+i] != 0.0) {
129: /* [Ax_i Ax_i+1 Bx_i Bx_i+1]*= [ 1 0
130: 0 1
131: -eigr_i -eigi_i
132: eigi_i -eigr_i] */
133: MatDenseGetArray(M,&b);
134: b[0] = b[5] = 1.0/d->nX[r_s+i];
135: b[2] = b[7] = -d->eigr[r_s+i]/d->nX[r_s+i];
136: b[6] = -(b[3] = d->eigi[r_s+i]/d->nX[r_s+i]);
137: b[1] = b[4] = 0.0;
138: MatDenseRestoreArray(M,&b);
139: BVInsertVec(X,0,Ax[i]);
140: BVInsertVec(X,1,Ax[i+1]);
141: BVInsertVec(X,2,Bx[i]);
142: BVInsertVec(X,3,Bx[i+1]);
143: BVSetActiveColumns(X,0,4);
144: BVMultInPlace(X,M,0,2);
145: BVCopyVec(X,0,Ax[i]);
146: BVCopyVec(X,1,Ax[i+1]);
147: s = 2;
148: } else
149: #endif
150: {
151: /* [Ax_i Bx_i]*= [ 1/nX_i conj(eig_i/nX_i)
152: -eig_i/nX_i 1/nX_i ] */
153: MatDenseGetArray(M,&b);
154: b[0] = 1.0/d->nX[r_s+i];
155: b[1] = -d->eigr[r_s+i]/d->nX[r_s+i];
156: b[4] = PetscConj(d->eigr[r_s+i]/d->nX[r_s+i]);
157: b[5] = 1.0/d->nX[r_s+i];
158: MatDenseRestoreArray(M,&b);
159: BVInsertVec(X,0,Ax[i]);
160: BVInsertVec(X,1,Bx[i]);
161: BVSetActiveColumns(X,0,2);
162: BVMultInPlace(X,M,0,2);
163: BVCopyVec(X,0,Ax[i]);
164: BVCopyVec(X,1,Bx[i]);
165: s = 1;
166: }
167: for (j=0;j<s;j++) d->nX[r_s+i+j] = 1.0;
169: /* Ax = R <- P*(Ax - eig_i*Bx) */
170: d->calcpairs_proj_res(d,r_s+i,r_s+i+s,&Ax[i]);
172: /* Check if the first eigenpairs are converged */
173: if (i == 0) {
174: d->preTestConv(d,0,s,s,&d->npreconv);
175: if (d->npreconv > 0) break;
176: }
177: }
179: /* D <- K*[Ax Bx] */
180: if (d->npreconv == 0) {
181: for (i=0;i<n;i++) {
182: BVGetColumn(d->eps->V,kv+i,&v);
183: d->improvex_precond(d,r_s+i,Ax[i],v);
184: BVRestoreColumn(d->eps->V,kv+i,&v);
185: }
186: for (i=n;i<n*2;i++) {
187: BVGetColumn(d->eps->V,kv+i,&v);
188: d->improvex_precond(d,r_s+i-n,Bx[i-n],v);
189: BVRestoreColumn(d->eps->V,kv+i,&v);
190: }
191: *size_D = 2*n;
192: #if !defined(PETSC_USE_COMPLEX)
193: if (d->eigi[r_s] != 0.0) {
194: s = 4;
195: } else
196: #endif
197: {
198: s = 2;
199: }
200: /* Prevent that short vectors are discarded in the orthogonalization */
201: for (i=0; i<s && i<*size_D; i++) {
202: if (d->eps->errest[d->nconv+r_s+i] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s+i] < PETSC_MAX_REAL) {
203: BVScaleColumn(d->eps->V,i+kv,1.0/d->eps->errest[d->nconv+r_s+i]);
204: }
205: }
206: } else *size_D = 0;
208: SlepcVecPoolRestoreVecs(d->auxV,n,&Bx);
209: SlepcVecPoolRestoreVecs(d->auxV,n,&Ax);
210: BVDestroy(&X);
211: MatDestroy(&M);
212: return(0);
213: }
217: PetscErrorCode dvd_improvex_gd2(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs)
218: {
219: PetscErrorCode ierr;
220: dvdImprovex_gd2 *data;
221: PC pc;
224: /* Setting configuration constrains */
225: /* If the arithmetic is real and the problem is not Hermitian, then
226: the block size is incremented in one */
227: #if !defined(PETSC_USE_COMPLEX)
228: if (!DVD_IS(d->sEP, DVD_EP_HERMITIAN)) {
229: max_bs++;
230: b->max_size_P = PetscMax(b->max_size_P,2);
231: } else
232: #endif
233: {
234: b->max_size_P = PetscMax(b->max_size_P,1);
235: }
236: b->max_size_X = PetscMax(b->max_size_X,max_bs);
238: /* Setup the preconditioner */
239: if (ksp) {
240: KSPGetPC(ksp,&pc);
241: dvd_static_precond_PC(d,b,pc);
242: } else {
243: dvd_static_precond_PC(d,b,0);
244: }
246: /* Setup the step */
247: if (b->state >= DVD_STATE_CONF) {
248: PetscNewLog(d->eps,&data);
249: d->improveX_data = data;
250: data->size_X = b->max_size_X;
251: d->improveX = dvd_improvex_gd2_gen;
253: EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_gd2_d);
254: }
255: return(0);
256: }