Actual source code: dvdupdatev.c
slepc-3.8.0 2017-10-20
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, 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: test for restarting, updateV, restartV
14: */
16: #include "davidson.h"
18: typedef struct {
19: PetscInt min_size_V; /* restart with this number of eigenvectors */
20: PetscInt plusk; /* at restart, save plusk vectors from last iteration */
21: PetscInt mpd; /* max size of the searching subspace */
22: void *old_updateV_data; /* old updateV data */
23: PetscErrorCode (*old_isRestarting)(dvdDashboard*,PetscBool*); /* old isRestarting */
24: Mat oldU; /* previous projected right igenvectors */
25: Mat oldV; /* previous projected left eigenvectors */
26: PetscInt size_oldU; /* size of oldU */
27: PetscBool allResiduals; /* if computing all the residuals */
28: } dvdManagV_basic;
30: static PetscErrorCode dvd_updateV_start(dvdDashboard *d)
31: {
32: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
33: PetscInt i;
36: for (i=0;i<d->eps->ncv;i++) d->eigi[i] = 0.0;
37: d->nR = d->real_nR;
38: for (i=0;i<d->eps->ncv;i++) d->nR[i] = 1.0;
39: d->nX = d->real_nX;
40: for (i=0;i<d->eps->ncv;i++) d->errest[i] = 1.0;
41: data->size_oldU = 0;
42: d->nconv = 0;
43: d->npreconv = 0;
44: d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
45: d->size_D = 0;
46: return(0);
47: }
49: static PetscErrorCode dvd_isrestarting_fullV(dvdDashboard *d,PetscBool *r)
50: {
51: PetscErrorCode ierr;
52: PetscInt l,k;
53: PetscBool restart;
54: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
57: BVGetActiveColumns(d->eps->V,&l,&k);
58: restart = (k+2 > d->eps->ncv)? PETSC_TRUE: PETSC_FALSE;
60: /* Check old isRestarting function */
61: if (!restart && data->old_isRestarting) {
62: data->old_isRestarting(d,&restart);
63: }
64: *r = restart;
65: return(0);
66: }
68: static PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
69: {
70: PetscErrorCode ierr;
71: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
74: /* Restore changes in dvdDashboard */
75: d->updateV_data = data->old_updateV_data;
77: /* Free local data */
78: MatDestroy(&data->oldU);
79: MatDestroy(&data->oldV);
80: PetscFree(d->real_nR);
81: PetscFree(d->real_nX);
82: PetscFree(data);
83: return(0);
84: }
86: static PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
87: {
88: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
89: PetscInt npreconv,cMT,cMTX,lV,kV,nV;
90: PetscErrorCode ierr;
91: Mat Z;
92: PetscBool t;
93: #if !defined(PETSC_USE_COMPLEX)
94: PetscInt i;
95: #endif
98: npreconv = d->npreconv;
99: /* Constrains the converged pairs to nev */
100: #if !defined(PETSC_USE_COMPLEX)
101: /* Tries to maintain together conjugate eigenpairs */
102: 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));
103: npreconv = i;
104: #else
105: npreconv = PetscMax(PetscMin(d->nev-d->nconv,npreconv),0);
106: #endif
107: /* For GHEP without B-ortho, converge all of the requested pairs at once */
108: PetscObjectTypeCompare((PetscObject)d->eps->ds,DSGHEP,&t);
109: if (t && d->nconv+npreconv<d->nev) npreconv = 0;
110: /* Quick exit */
111: if (npreconv == 0) return(0);
113: BVGetActiveColumns(d->eps->V,&lV,&kV);
114: nV = kV - lV;
115: cMT = nV - npreconv;
116: /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
117: If the problem is standard or hermitian, left and right vectors are the same */
118: if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
119: /* ps.Q <- [ps.Q(0:npreconv-1) ps.Z(npreconv:size_H-1)] */
120: DSGetMat(d->eps->ds,DS_MAT_Z,&Z);
121: DSCopyMat(d->eps->ds,DS_MAT_Q,0,npreconv,Z,0,npreconv,nV,cMT,PETSC_FALSE);
122: MatDestroy(&Z);
123: }
124: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
125: DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,nV,d->nBds,&cMTX,d->nBds);
126: } else {
127: DSOrthogonalize(d->eps->ds,DS_MAT_Q,nV,&cMTX);
128: }
129: cMT = cMTX - npreconv;
131: if (d->W) {
132: DSOrthogonalize(d->eps->ds,DS_MAT_Z,nV,&cMTX);
133: cMT = PetscMin(cMT,cMTX - npreconv);
134: }
136: /* Lock the converged pairs */
137: d->eigr+= npreconv;
138: #if !defined(PETSC_USE_COMPLEX)
139: if (d->eigi) d->eigi+= npreconv;
140: #endif
141: d->nconv+= npreconv;
142: d->errest+= npreconv;
143: /* Notify the changes in V and update the other subspaces */
144: d->V_tra_s = npreconv; d->V_tra_e = nV;
145: d->V_new_s = cMT; d->V_new_e = d->V_new_s;
146: /* Remove oldU */
147: data->size_oldU = 0;
149: d->npreconv-= npreconv;
150: return(0);
151: }
153: static PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
154: {
155: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
156: PetscInt lV,kV,nV,size_plusk,size_X,cMTX,cMTY;
157: Mat Z;
158: PetscErrorCode ierr;
161: /* Select size_X desired pairs from V */
162: BVGetActiveColumns(d->eps->V,&lV,&kV);
163: nV = kV - lV;
164: size_X = PetscMin(data->min_size_V+d->npreconv,nV);
166: /* Add plusk eigenvectors from the previous iteration */
167: size_plusk = PetscMax(0,PetscMin(PetscMin(data->plusk,data->size_oldU),d->eps->ncv - size_X));
169: d->size_MT = nV;
170: /* ps.Q <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
171: /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
172: If the problem is standard or hermitian, left and right vectors are the same */
173: if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
174: DSGetMat(d->eps->ds,DS_MAT_Z,&Z);
175: DSCopyMat(d->eps->ds,DS_MAT_Q,0,0,Z,0,0,nV,size_X,PETSC_FALSE);
176: MatDestroy(&Z);
177: }
178: 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");
179: if (size_plusk > 0) {
180: DSCopyMat(d->eps->ds,DS_MAT_Q,0,size_X,data->oldU,0,0,nV,size_plusk,PETSC_FALSE);
181: }
182: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
183: DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,size_X,d->nBds,&cMTX,d->nBds);
184: } else {
185: DSOrthogonalize(d->eps->ds,DS_MAT_Q,size_X+size_plusk,&cMTX);
186: }
188: if (d->W && size_plusk > 0) {
189: /* ps.Z <- orth([ps.Z(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
190: DSCopyMat(d->eps->ds,DS_MAT_Z,0,size_X,data->oldV,0,0,nV,size_plusk,PETSC_FALSE);
191: DSOrthogonalize(d->eps->ds,DS_MAT_Z,size_X+size_plusk,&cMTY);
192: cMTX = PetscMin(cMTX, cMTY);
193: }
195: /* Notify the changes in V and update the other subspaces */
196: d->V_tra_s = 0; d->V_tra_e = cMTX;
197: d->V_new_s = d->V_tra_e; d->V_new_e = d->V_new_s;
199: /* Remove oldU */
200: data->size_oldU = 0;
202: /* Remove npreconv */
203: d->npreconv = 0;
204: return(0);
205: }
207: static PetscErrorCode dvd_updateV_testConv(dvdDashboard *d,PetscInt s,PetscInt pre,PetscInt e,PetscInt *nConv)
208: {
209: PetscInt i,j,b;
210: PetscReal norm;
211: PetscErrorCode ierr;
212: PetscBool conv, c;
213: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
216: if (nConv) *nConv = s;
217: for (i=s,conv=PETSC_TRUE;(conv || data->allResiduals) && (i < e);i+=b) {
218: #if !defined(PETSC_USE_COMPLEX)
219: b = d->eigi[i]!=0.0?2:1;
220: #else
221: b = 1;
222: #endif
223: if (i+b-1 >= pre) {
224: d->calcpairs_residual(d,i,i+b);
225: }
226: /* Test the Schur vector */
227: for (j=0,c=PETSC_TRUE;j<b && c;j++) {
228: norm = d->nR[i+j]/d->nX[i+j];
229: c = d->testConv(d,d->eigr[i+j],d->eigi[i+j],norm,&d->errest[i+j]);
230: }
231: if (conv && c) { if (nConv) *nConv = i+b; }
232: else conv = PETSC_FALSE;
233: }
234: pre = PetscMax(pre,i);
236: #if !defined(PETSC_USE_COMPLEX)
237: /* Enforce converged conjugate complex eigenpairs */
238: if (nConv) {
239: for (j=0;j<*nConv;j++) if (d->eigi[j] != 0.0) j++;
240: if (j>*nConv) (*nConv)--;
241: }
242: #endif
243: for (i=pre;i<e;i++) d->errest[i] = d->nR[i] = 1.0;
244: return(0);
245: }
247: static PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
248: {
249: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
250: PetscInt size_D,s,lV,kV,nV;
251: PetscErrorCode ierr;
254: /* Select the desired pairs */
255: BVGetActiveColumns(d->eps->V,&lV,&kV);
256: nV = kV - lV;
257: size_D = PetscMin(PetscMin(PetscMin(d->bs,nV),d->eps->ncv-nV),nV);
258: if (size_D == 0) return(0);
260: /* Fill V with D */
261: d->improveX(d,d->npreconv,d->npreconv+size_D,&size_D);
263: /* If D is empty, exit */
264: d->size_D = size_D;
265: if (size_D == 0) return(0);
267: /* Get the residual of all pairs */
268: #if !defined(PETSC_USE_COMPLEX)
269: s = (d->eigi[0]!=0.0)? 2: 1;
270: #else
271: s = 1;
272: #endif
273: BVGetActiveColumns(d->eps->V,&lV,&kV);
274: nV = kV - lV;
275: dvd_updateV_testConv(d,s,s,data->allResiduals?nV:size_D,NULL);
277: /* Notify the changes in V */
278: d->V_tra_s = 0; d->V_tra_e = 0;
279: d->V_new_s = nV; d->V_new_e = nV+size_D;
281: /* Save the projected eigenvectors */
282: if (data->plusk > 0) {
283: MatZeroEntries(data->oldU);
284: data->size_oldU = nV;
285: DSCopyMat(d->eps->ds,DS_MAT_Q,0,0,data->oldU,0,0,nV,nV,PETSC_TRUE);
286: if (d->W) {
287: MatZeroEntries(data->oldV);
288: DSCopyMat(d->eps->ds,DS_MAT_Z,0,0,data->oldV,0,0,nV,nV,PETSC_TRUE);
289: }
290: }
291: return(0);
292: }
294: static PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
295: {
296: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
297: PetscInt i;
298: PetscBool restart,t;
299: PetscErrorCode ierr;
302: /* TODO: restrict select pairs to each case */
303: d->calcpairs_selectPairs(d, data->min_size_V+d->npreconv);
305: /* If the subspaces doesn't need restart, add new vector */
306: d->isRestarting(d,&restart);
307: if (!restart) {
308: d->size_D = 0;
309: dvd_updateV_update_gen(d);
311: /* If no vector were converged, exit */
312: /* For GHEP without B-ortho, converge all of the requested pairs at once */
313: PetscObjectTypeCompareAny((PetscObject)d->eps->ds,&t,DSGHEP,"");
314: if (d->nconv+d->npreconv < d->nev && (t || d->npreconv == 0)) return(0);
315: }
317: /* If some eigenpairs were converged, lock them */
318: if (d->npreconv > 0) {
319: i = d->npreconv;
320: dvd_updateV_conv_gen(d);
322: /* If some eigenpair was locked, exit */
323: if (i > d->npreconv) return(0);
324: }
326: /* Else, a restarting is performed */
327: dvd_updateV_restart_gen(d);
328: return(0);
329: }
331: PetscErrorCode dvd_managementV_basic(dvdDashboard *d,dvdBlackboard *b,PetscInt bs,PetscInt mpd,PetscInt min_size_V,PetscInt plusk,PetscBool harm,PetscBool allResiduals)
332: {
333: PetscErrorCode ierr;
334: dvdManagV_basic *data;
335: #if !defined(PETSC_USE_COMPLEX)
336: PetscBool her_probl,std_probl;
337: #endif
340: /* Setting configuration constrains */
341: #if !defined(PETSC_USE_COMPLEX)
342: /* if the last converged eigenvalue is complex its conjugate pair is also
343: converged */
344: her_probl = DVD_IS(d->sEP,DVD_EP_HERMITIAN)? PETSC_TRUE: PETSC_FALSE;
345: std_probl = DVD_IS(d->sEP,DVD_EP_STD)? PETSC_TRUE: PETSC_FALSE;
346: b->max_size_X = PetscMax(b->max_size_X,bs+((her_probl && std_probl)?0:1));
347: #else
348: b->max_size_X = PetscMax(b->max_size_X,bs);
349: #endif
351: b->max_size_V = PetscMax(b->max_size_V,mpd);
352: min_size_V = PetscMin(min_size_V,mpd-bs);
353: b->size_V = PetscMax(b->size_V,b->max_size_V+b->max_size_P+b->max_nev);
354: b->max_size_oldX = plusk;
356: /* Setup the step */
357: if (b->state >= DVD_STATE_CONF) {
358: PetscNewLog(d->eps,&data);
359: data->mpd = b->max_size_V;
360: data->min_size_V = min_size_V;
361: d->bs = bs;
362: data->plusk = plusk;
363: data->allResiduals = allResiduals;
365: d->eigr = d->eps->eigr;
366: d->eigi = d->eps->eigi;
367: d->errest = d->eps->errest;
368: PetscMalloc1(d->eps->ncv,&d->real_nR);
369: PetscMalloc1(d->eps->ncv,&d->real_nX);
370: if (plusk > 0) { MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldU); }
371: else data->oldU = NULL;
372: if (harm && plusk>0) { MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldV); }
373: else data->oldV = NULL;
375: data->old_updateV_data = d->updateV_data;
376: d->updateV_data = data;
377: data->old_isRestarting = d->isRestarting;
378: d->isRestarting = dvd_isrestarting_fullV;
379: d->updateV = dvd_updateV_extrapol;
380: d->preTestConv = dvd_updateV_testConv;
381: EPSDavidsonFLAdd(&d->startList,dvd_updateV_start);
382: EPSDavidsonFLAdd(&d->destroyList,dvd_managementV_basic_d);
383: }
384: return(0);
385: }