Actual source code: dvdupdatev.c
slepc-3.7.3 2016-09-29
1: /*
2: SLEPc eigensolver: "davidson"
4: Step: test for restarting, updateV, restartV
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 min_size_V; /* restart with this number of eigenvectors */
30: PetscInt plusk; /* at restart, save plusk vectors from last iteration */
31: PetscInt mpd; /* max size of the searching subspace */
32: void *old_updateV_data; /* old updateV data */
33: PetscErrorCode (*old_isRestarting)(dvdDashboard*,PetscBool*); /* old isRestarting */
34: Mat oldU; /* previous projected right igenvectors */
35: Mat oldV; /* previous projected left eigenvectors */
36: PetscInt size_oldU; /* size of oldU */
37: PetscBool allResiduals; /* if computing all the residuals */
38: } dvdManagV_basic;
42: static PetscErrorCode dvd_updateV_start(dvdDashboard *d)
43: {
44: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
45: PetscInt i;
48: for (i=0;i<d->eps->ncv;i++) d->eigi[i] = 0.0;
49: d->nR = d->real_nR;
50: for (i=0;i<d->eps->ncv;i++) d->nR[i] = PETSC_MAX_REAL;
51: d->nX = d->real_nX;
52: for (i=0;i<d->eps->ncv;i++) d->errest[i] = PETSC_MAX_REAL;
53: data->size_oldU = 0;
54: d->nconv = 0;
55: d->npreconv = 0;
56: d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
57: d->size_D = 0;
58: return(0);
59: }
63: static PetscErrorCode dvd_isrestarting_fullV(dvdDashboard *d,PetscBool *r)
64: {
65: PetscErrorCode ierr;
66: PetscInt l,k;
67: PetscBool restart;
68: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
71: BVGetActiveColumns(d->eps->V,&l,&k);
72: restart = (k+2 > d->eps->ncv)? PETSC_TRUE: PETSC_FALSE;
74: /* Check old isRestarting function */
75: if (!restart && data->old_isRestarting) {
76: data->old_isRestarting(d,&restart);
77: }
78: *r = restart;
79: return(0);
80: }
84: static PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
85: {
86: PetscErrorCode ierr;
87: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
90: /* Restore changes in dvdDashboard */
91: d->updateV_data = data->old_updateV_data;
93: /* Free local data */
94: if (data->oldU) { MatDestroy(&data->oldU); }
95: if (data->oldV) { MatDestroy(&data->oldV); }
96: PetscFree(d->real_nR);
97: PetscFree(d->real_nX);
98: PetscFree(data);
99: return(0);
100: }
104: static PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
105: {
106: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
107: PetscInt npreconv,cMT,cMTX,lV,kV,nV;
108: PetscErrorCode ierr;
109: Mat Q;
110: #if !defined(PETSC_USE_COMPLEX)
111: PetscInt i;
112: #endif
115: npreconv = d->npreconv;
116: /* Constrains the converged pairs to nev */
117: #if !defined(PETSC_USE_COMPLEX)
118: /* Tries to maintain together conjugate eigenpairs */
119: for (i=0; (i + (d->eigi[i]!=0.0?1:0) < npreconv) && (d->nconv + i < d->nev); i+= (d->eigi[i]!=0.0?2:1));
120: npreconv = i;
121: #else
122: npreconv = PetscMax(PetscMin(d->nev-d->nconv,npreconv),0);
123: #endif
124: /* Quick exit */
125: if (npreconv == 0) return(0);
127: BVGetActiveColumns(d->eps->V,&lV,&kV);
128: nV = kV - lV;
129: cMT = nV - npreconv;
130: /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
131: If the problem is standard or hermitian, left and right vectors are the same */
132: if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
133: /* ps.Q <- [ps.Q(0:npreconv-1) ps.Z(npreconv:size_H-1)] */
134: DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
135: DSCopyMat(d->eps->ds,DS_MAT_Z,0,npreconv,Q,0,npreconv,nV,cMT,PETSC_TRUE);
136: DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q);
137: }
138: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
139: DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,nV,d->nBds,&cMTX,d->nBds);
140: } else {
141: DSOrthogonalize(d->eps->ds,DS_MAT_Q,nV,&cMTX);
142: }
143: cMT = cMTX - npreconv;
145: if (d->W) {
146: DSOrthogonalize(d->eps->ds,DS_MAT_Z,nV,&cMTX);
147: cMT = PetscMin(cMT,cMTX - npreconv);
148: }
150: /* Lock the converged pairs */
151: d->eigr+= npreconv;
152: #if !defined(PETSC_USE_COMPLEX)
153: if (d->eigi) d->eigi+= npreconv;
154: #endif
155: d->nconv+= npreconv;
156: d->errest+= npreconv;
157: /* Notify the changes in V and update the other subspaces */
158: d->V_tra_s = npreconv; d->V_tra_e = nV;
159: d->V_new_s = cMT; d->V_new_e = d->V_new_s;
160: /* Remove oldU */
161: data->size_oldU = 0;
163: d->npreconv-= npreconv;
164: return(0);
165: }
169: static PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
170: {
171: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
172: PetscInt lV,kV,nV,size_plusk,size_X,cMTX,cMTY;
173: Mat Q;
174: PetscErrorCode ierr;
177: /* Select size_X desired pairs from V */
178: BVGetActiveColumns(d->eps->V,&lV,&kV);
179: nV = kV - lV;
180: size_X = PetscMin(data->min_size_V,nV);
182: /* Add plusk eigenvectors from the previous iteration */
183: size_plusk = PetscMax(0,PetscMin(PetscMin(data->plusk,data->size_oldU),d->eps->ncv - size_X));
185: d->size_MT = nV;
186: /* ps.Q <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
187: /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
188: If the problem is standard or hermitian, left and right vectors are the same */
189: if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
190: DSGetMat(d->eps->ds,DS_MAT_Q,&Q);
191: DSCopyMat(d->eps->ds,DS_MAT_Z,0,0,Q,0,0,nV,size_X,PETSC_TRUE);
192: DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q);
193: }
194: if (size_plusk > 0 && DVD_IS(d->sEP,DVD_EP_INDEFINITE)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported plusk>0 in indefinite eigenvalue problems");
195: if (size_plusk > 0) {
196: DSCopyMat(d->eps->ds,DS_MAT_Q,0,size_X,data->oldU,0,0,nV,size_plusk,PETSC_FALSE);
197: }
198: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
199: DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,size_X,d->nBds,&cMTX,d->nBds);
200: } else {
201: DSOrthogonalize(d->eps->ds,DS_MAT_Q,size_X+size_plusk,&cMTX);
202: }
204: if (d->W && size_plusk > 0) {
205: /* ps.Z <- orth([ps.Z(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
206: DSCopyMat(d->eps->ds,DS_MAT_Z,0,size_X,data->oldV,0,0,nV,size_plusk,PETSC_FALSE);
207: DSOrthogonalize(d->eps->ds,DS_MAT_Z,size_X+size_plusk,&cMTY);
208: cMTX = PetscMin(cMTX, cMTY);
209: }
211: /* Notify the changes in V and update the other subspaces */
212: d->V_tra_s = 0; d->V_tra_e = cMTX;
213: d->V_new_s = d->V_tra_e; d->V_new_e = d->V_new_s;
215: /* Remove oldU */
216: data->size_oldU = 0;
218: /* Remove npreconv */
219: d->npreconv = 0;
220: return(0);
221: }
225: static PetscErrorCode dvd_updateV_testConv(dvdDashboard *d,PetscInt s,PetscInt pre,PetscInt e,PetscInt *nConv)
226: {
227: PetscInt i,j,b;
228: PetscReal norm;
229: PetscErrorCode ierr;
230: PetscBool conv, c;
231: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
234: if (nConv) *nConv = s;
235: for (i=s,conv=PETSC_TRUE;(conv || data->allResiduals) && (i < e);i+=b) {
236: #if !defined(PETSC_USE_COMPLEX)
237: b = d->eigi[i]!=0.0?2:1;
238: #else
239: b = 1;
240: #endif
241: if (i+b-1 >= pre) {
242: d->calcpairs_residual(d,i,i+b);
243: }
244: /* Test the Schur vector */
245: for (j=0,c=PETSC_TRUE;j<b && c;j++) {
246: norm = d->nR[i+j]/d->nX[i+j];
247: c = d->testConv(d,d->eigr[i+j],d->eigi[i+j],norm,&d->errest[i+j]);
248: }
249: if (conv && c) { if (nConv) *nConv = i+b; }
250: else conv = PETSC_FALSE;
251: }
252: pre = PetscMax(pre,i);
254: #if !defined(PETSC_USE_COMPLEX)
255: /* Enforce converged conjugate complex eigenpairs */
256: if (nConv) {
257: for (j=0;j<*nConv;j++) if (d->eigi[j] != 0.0) j++;
258: if (j>*nConv) (*nConv)--;
259: }
260: #endif
261: for (i=pre;i<e;i++) d->errest[i] = d->nR[i] = PETSC_MAX_REAL;
262: return(0);
263: }
267: static PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
268: {
269: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
270: PetscInt size_D,s,lV,kV,nV;
271: PetscErrorCode ierr;
274: /* Select the desired pairs */
275: BVGetActiveColumns(d->eps->V,&lV,&kV);
276: nV = kV - lV;
277: size_D = PetscMin(PetscMin(PetscMin(d->bs,nV),d->eps->ncv-nV),nV);
278: if (size_D == 0) return(0);
280: /* Fill V with D */
281: d->improveX(d,0,size_D,&size_D);
283: /* If D is empty, exit */
284: d->size_D = size_D;
285: if (size_D == 0) return(0);
287: /* Get the residual of all pairs */
288: #if !defined(PETSC_USE_COMPLEX)
289: s = (d->eigi[0]!=0.0)? 2: 1;
290: #else
291: s = 1;
292: #endif
293: BVGetActiveColumns(d->eps->V,&lV,&kV);
294: nV = kV - lV;
295: dvd_updateV_testConv(d,s,s,data->allResiduals?nV:size_D,NULL);
297: /* Notify the changes in V */
298: d->V_tra_s = 0; d->V_tra_e = 0;
299: d->V_new_s = nV; d->V_new_e = nV+size_D;
301: /* Save the projected eigenvectors */
302: if (data->plusk > 0) {
303: MatZeroEntries(data->oldU);
304: data->size_oldU = nV;
305: DSCopyMat(d->eps->ds,DS_MAT_Q,0,0,data->oldU,0,0,nV,nV,PETSC_TRUE);
306: if (d->W) {
307: MatZeroEntries(data->oldV);
308: DSCopyMat(d->eps->ds,DS_MAT_Z,0,0,data->oldV,0,0,nV,nV,PETSC_TRUE);
309: }
310: }
311: return(0);
312: }
316: static PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
317: {
318: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
319: PetscInt i;
320: PetscBool restart;
321: PetscErrorCode ierr;
324: /* TODO: restrict select pairs to each case */
325: d->calcpairs_selectPairs(d, data->min_size_V);
327: /* If the subspaces doesn't need restart, add new vector */
328: d->isRestarting(d,&restart);
329: if (!restart) {
330: d->size_D = 0;
331: dvd_updateV_update_gen(d);
333: /* If some vector were add, exit */
334: if (d->size_D > 0) return(0);
335: }
337: /* If some eigenpairs were converged, lock them */
338: if (d->npreconv > 0) {
339: i = d->npreconv;
340: dvd_updateV_conv_gen(d);
342: /* If some eigenpair was locked, exit */
343: if (i > d->npreconv) return(0);
344: }
346: /* Else, a restarting is performed */
347: dvd_updateV_restart_gen(d);
348: return(0);
349: }
353: PetscErrorCode dvd_managementV_basic(dvdDashboard *d,dvdBlackboard *b,PetscInt bs,PetscInt mpd,PetscInt min_size_V,PetscInt plusk,PetscBool harm,PetscBool allResiduals)
354: {
355: PetscErrorCode ierr;
356: dvdManagV_basic *data;
357: #if !defined(PETSC_USE_COMPLEX)
358: PetscBool her_probl,std_probl;
359: #endif
362: /* Setting configuration constrains */
363: #if !defined(PETSC_USE_COMPLEX)
364: /* if the last converged eigenvalue is complex its conjugate pair is also
365: converged */
366: her_probl = DVD_IS(d->sEP,DVD_EP_HERMITIAN)? PETSC_TRUE: PETSC_FALSE;
367: std_probl = DVD_IS(d->sEP,DVD_EP_STD)? PETSC_TRUE: PETSC_FALSE;
368: b->max_size_X = PetscMax(b->max_size_X,bs+((her_probl && std_probl)?0:1));
369: #else
370: b->max_size_X = PetscMax(b->max_size_X,bs);
371: #endif
373: b->max_size_V = PetscMax(b->max_size_V,mpd);
374: min_size_V = PetscMin(min_size_V,mpd-bs);
375: b->size_V = PetscMax(b->size_V,b->max_size_V+b->max_size_P+b->max_nev);
376: b->max_size_oldX = plusk;
378: /* Setup the step */
379: if (b->state >= DVD_STATE_CONF) {
380: PetscNewLog(d->eps,&data);
381: data->mpd = b->max_size_V;
382: data->min_size_V = min_size_V;
383: d->bs = bs;
384: data->plusk = plusk;
385: data->allResiduals = allResiduals;
387: d->eigr = d->eps->eigr;
388: d->eigi = d->eps->eigi;
389: d->errest = d->eps->errest;
390: PetscMalloc1(d->eps->ncv,&d->real_nR);
391: PetscMalloc1(d->eps->ncv,&d->real_nX);
392: if (plusk > 0) { MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldU); }
393: else data->oldU = NULL;
394: if (harm && plusk>0) { MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldV); }
395: else data->oldV = NULL;
397: data->old_updateV_data = d->updateV_data;
398: d->updateV_data = data;
399: data->old_isRestarting = d->isRestarting;
400: d->isRestarting = dvd_isrestarting_fullV;
401: d->updateV = dvd_updateV_extrapol;
402: d->preTestConv = dvd_updateV_testConv;
403: EPSDavidsonFLAdd(&d->startList,dvd_updateV_start);
404: EPSDavidsonFLAdd(&d->destroyList,dvd_managementV_basic_d);
405: }
406: return(0);
407: }