Actual source code: dvdupdatev.c
slepc-3.17.1 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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;
35: for (i=0;i<d->eps->ncv;i++) d->eigi[i] = 0.0;
36: d->nR = d->real_nR;
37: for (i=0;i<d->eps->ncv;i++) d->nR[i] = 1.0;
38: d->nX = d->real_nX;
39: for (i=0;i<d->eps->ncv;i++) d->errest[i] = 1.0;
40: data->size_oldU = 0;
41: d->nconv = 0;
42: d->npreconv = 0;
43: d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
44: d->size_D = 0;
45: PetscFunctionReturn(0);
46: }
48: static PetscErrorCode dvd_isrestarting_fullV(dvdDashboard *d,PetscBool *r)
49: {
50: PetscInt l,k;
51: PetscBool restart;
52: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
54: BVGetActiveColumns(d->eps->V,&l,&k);
55: restart = (k+2 > d->eps->ncv)? PETSC_TRUE: PETSC_FALSE;
57: /* Check old isRestarting function */
58: if (PetscUnlikely(!restart && data->old_isRestarting)) data->old_isRestarting(d,&restart);
59: *r = restart;
60: PetscFunctionReturn(0);
61: }
63: static PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
64: {
65: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
67: /* Restore changes in dvdDashboard */
68: d->updateV_data = data->old_updateV_data;
70: /* Free local data */
71: MatDestroy(&data->oldU);
72: MatDestroy(&data->oldV);
73: PetscFree(d->real_nR);
74: PetscFree(d->real_nX);
75: PetscFree(data);
76: PetscFunctionReturn(0);
77: }
79: static PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
80: {
81: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
82: PetscInt npreconv,cMT,cMTX,lV,kV,nV;
83: Mat Z;
84: PetscBool t;
85: #if !defined(PETSC_USE_COMPLEX)
86: PetscInt i;
87: #endif
89: npreconv = d->npreconv;
90: /* Constrains the converged pairs to nev */
91: #if !defined(PETSC_USE_COMPLEX)
92: /* Tries to maintain together conjugate eigenpairs */
93: 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));
94: npreconv = i;
95: #else
96: npreconv = PetscMax(PetscMin(d->nev-d->nconv,npreconv),0);
97: #endif
98: /* For GHEP without B-ortho, converge all of the requested pairs at once */
99: PetscObjectTypeCompare((PetscObject)d->eps->ds,DSGHEP,&t);
100: if (t && d->nconv+npreconv<d->nev) npreconv = 0;
101: /* Quick exit */
102: if (npreconv == 0) PetscFunctionReturn(0);
104: BVGetActiveColumns(d->eps->V,&lV,&kV);
105: nV = kV - lV;
106: cMT = nV - npreconv;
107: /* Harmonics restarts with right eigenvectors, and other with the left ones.
108: If the problem is standard or hermitian, left and right vectors are the same */
109: if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
110: /* ps.Q <- [ps.Q(0:npreconv-1) ps.Z(npreconv:size_H-1)] */
111: DSGetMat(d->eps->ds,DS_MAT_Z,&Z);
112: DSCopyMat(d->eps->ds,DS_MAT_Q,0,npreconv,Z,0,npreconv,nV,cMT,PETSC_FALSE);
113: MatDestroy(&Z);
114: }
115: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,nV,d->nBds,&cMTX,d->nBds);
116: else DSOrthogonalize(d->eps->ds,DS_MAT_Q,nV,&cMTX);
117: cMT = cMTX - npreconv;
119: if (d->W) {
120: DSOrthogonalize(d->eps->ds,DS_MAT_Z,nV,&cMTX);
121: cMT = PetscMin(cMT,cMTX - npreconv);
122: }
124: /* Lock the converged pairs */
125: d->eigr+= npreconv;
126: #if !defined(PETSC_USE_COMPLEX)
127: if (d->eigi) d->eigi+= npreconv;
128: #endif
129: d->nconv+= npreconv;
130: d->errest+= npreconv;
131: /* Notify the changes in V and update the other subspaces */
132: d->V_tra_s = npreconv; d->V_tra_e = nV;
133: d->V_new_s = cMT; d->V_new_e = d->V_new_s;
134: /* Remove oldU */
135: data->size_oldU = 0;
137: d->npreconv-= npreconv;
138: PetscFunctionReturn(0);
139: }
141: static PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
142: {
143: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
144: PetscInt lV,kV,nV,size_plusk,size_X,cMTX,cMTY,max_restart_size;
145: Mat Z;
147: /* Select size_X desired pairs from V */
148: /* The restarted basis should:
149: - have at least one spot to add a new direction;
150: - keep converged vectors, npreconv;
151: - keep at least 1 oldU direction if possible.
152: */
153: BVGetActiveColumns(d->eps->V,&lV,&kV);
154: nV = kV - lV;
155: max_restart_size = PetscMax(0,PetscMin(d->eps->mpd - 1,d->eps->ncv - lV - 2));
156: size_X = PetscMin(PetscMin(data->min_size_V+d->npreconv,max_restart_size - (max_restart_size - d->npreconv > 1 && data->plusk > 0 && data->size_oldU > 0 ? 1 : 0)), nV);
158: /* Add plusk eigenvectors from the previous iteration */
159: size_plusk = PetscMax(0,PetscMin(PetscMin(PetscMin(data->plusk,data->size_oldU),max_restart_size - size_X),nV - size_X));
161: d->size_MT = nV;
162: /* ps.Q <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
163: /* Harmonics restarts with right eigenvectors, and other with the left ones.
164: If the problem is standard or hermitian, left and right vectors are the same */
165: if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
166: DSGetMat(d->eps->ds,DS_MAT_Z,&Z);
167: DSCopyMat(d->eps->ds,DS_MAT_Q,0,0,Z,0,0,nV,size_X,PETSC_FALSE);
168: MatDestroy(&Z);
169: }
171: if (size_plusk > 0) DSCopyMat(d->eps->ds,DS_MAT_Q,0,size_X,data->oldU,0,0,nV,size_plusk,PETSC_FALSE);
172: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,size_X,d->nBds,&cMTX,d->nBds);
173: else DSOrthogonalize(d->eps->ds,DS_MAT_Q,size_X+size_plusk,&cMTX);
175: if (d->W && size_plusk > 0) {
176: /* ps.Z <- orth([ps.Z(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
177: DSCopyMat(d->eps->ds,DS_MAT_Z,0,size_X,data->oldV,0,0,nV,size_plusk,PETSC_FALSE);
178: DSOrthogonalize(d->eps->ds,DS_MAT_Z,size_X+size_plusk,&cMTY);
179: cMTX = PetscMin(cMTX, cMTY);
180: }
181: PetscAssert(cMTX<=size_X+size_plusk,PETSC_COMM_SELF,PETSC_ERR_SUP,"Invalid number of columns to restart");
183: /* Notify the changes in V and update the other subspaces */
184: d->V_tra_s = 0; d->V_tra_e = cMTX;
185: d->V_new_s = d->V_tra_e; d->V_new_e = d->V_new_s;
187: /* Remove oldU */
188: data->size_oldU = 0;
190: /* Remove npreconv */
191: d->npreconv = 0;
192: PetscFunctionReturn(0);
193: }
195: static PetscErrorCode dvd_updateV_testConv(dvdDashboard *d,PetscInt s,PetscInt pre,PetscInt e,PetscInt *nConv)
196: {
197: PetscInt i,j,b;
198: PetscReal norm;
199: PetscBool conv, c;
200: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
202: if (nConv) *nConv = s;
203: for (i=s,conv=PETSC_TRUE;(conv || data->allResiduals) && (i < e);i+=b) {
204: #if !defined(PETSC_USE_COMPLEX)
205: b = d->eigi[i]!=0.0?2:1;
206: #else
207: b = 1;
208: #endif
209: if (i+b-1 >= pre) d->calcpairs_residual(d,i,i+b);
210: /* Test the Schur vector */
211: for (j=0,c=PETSC_TRUE;j<b && c;j++) {
212: norm = d->nR[i+j]/d->nX[i+j];
213: c = d->testConv(d,d->eigr[i+j],d->eigi[i+j],norm,&d->errest[i+j]);
214: }
215: if (conv && c) { if (nConv) *nConv = i+b; }
216: else conv = PETSC_FALSE;
217: }
218: pre = PetscMax(pre,i);
220: #if !defined(PETSC_USE_COMPLEX)
221: /* Enforce converged conjugate complex eigenpairs */
222: if (nConv) {
223: for (j=0;j<*nConv;j++) if (d->eigi[j] != 0.0) j++;
224: if (j>*nConv) (*nConv)--;
225: }
226: #endif
227: for (i=pre;i<e;i++) d->errest[i] = d->nR[i] = 1.0;
228: PetscFunctionReturn(0);
229: }
231: static PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
232: {
233: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
234: PetscInt size_D,s,lV,kV,nV;
236: /* Select the desired pairs */
237: BVGetActiveColumns(d->eps->V,&lV,&kV);
238: nV = kV - lV;
239: size_D = PetscMin(PetscMin(PetscMin(d->bs,nV),d->eps->ncv-nV),nV);
240: if (size_D == 0) PetscFunctionReturn(0);
242: /* Fill V with D */
243: d->improveX(d,d->npreconv,d->npreconv+size_D,&size_D);
245: /* If D is empty, exit */
246: d->size_D = size_D;
247: if (size_D == 0) PetscFunctionReturn(0);
249: /* Get the residual of all pairs */
250: #if !defined(PETSC_USE_COMPLEX)
251: s = (d->eigi[0]!=0.0)? 2: 1;
252: #else
253: s = 1;
254: #endif
255: BVGetActiveColumns(d->eps->V,&lV,&kV);
256: nV = kV - lV;
257: dvd_updateV_testConv(d,s,s,data->allResiduals?nV:size_D,NULL);
259: /* Notify the changes in V */
260: d->V_tra_s = 0; d->V_tra_e = 0;
261: d->V_new_s = nV; d->V_new_e = nV+size_D;
263: /* Save the projected eigenvectors */
264: if (data->plusk > 0) {
265: MatZeroEntries(data->oldU);
266: data->size_oldU = nV;
267: DSCopyMat(d->eps->ds,DS_MAT_Q,0,0,data->oldU,0,0,nV,nV,PETSC_TRUE);
268: if (d->W) {
269: MatZeroEntries(data->oldV);
270: DSCopyMat(d->eps->ds,DS_MAT_Z,0,0,data->oldV,0,0,nV,nV,PETSC_TRUE);
271: }
272: }
273: PetscFunctionReturn(0);
274: }
276: static PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
277: {
278: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
279: PetscInt i;
280: PetscBool restart,t;
282: /* TODO: restrict select pairs to each case */
283: d->calcpairs_selectPairs(d, data->min_size_V+d->npreconv);
285: /* If the subspaces doesn't need restart, add new vector */
286: d->isRestarting(d,&restart);
287: if (!restart) {
288: d->size_D = 0;
289: dvd_updateV_update_gen(d);
291: /* If no vector were converged, exit */
292: /* For GHEP without B-ortho, converge all of the requested pairs at once */
293: PetscObjectTypeCompare((PetscObject)d->eps->ds,DSGHEP,&t);
294: if (d->nconv+d->npreconv < d->nev && (t || d->npreconv == 0)) PetscFunctionReturn(0);
295: }
297: /* If some eigenpairs were converged, lock them */
298: if (d->npreconv > 0) {
299: i = d->npreconv;
300: dvd_updateV_conv_gen(d);
302: /* If some eigenpair was locked, exit */
303: if (i > d->npreconv) PetscFunctionReturn(0);
304: }
306: /* Else, a restarting is performed */
307: dvd_updateV_restart_gen(d);
308: PetscFunctionReturn(0);
309: }
311: PetscErrorCode dvd_managementV_basic(dvdDashboard *d,dvdBlackboard *b,PetscInt bs,PetscInt mpd,PetscInt min_size_V,PetscInt plusk,PetscBool harm,PetscBool allResiduals)
312: {
313: dvdManagV_basic *data;
314: #if !defined(PETSC_USE_COMPLEX)
315: PetscBool her_probl,std_probl;
316: #endif
318: /* Setting configuration constrains */
319: #if !defined(PETSC_USE_COMPLEX)
320: /* if the last converged eigenvalue is complex its conjugate pair is also
321: converged */
322: her_probl = DVD_IS(d->sEP,DVD_EP_HERMITIAN)? PETSC_TRUE: PETSC_FALSE;
323: std_probl = DVD_IS(d->sEP,DVD_EP_STD)? PETSC_TRUE: PETSC_FALSE;
324: b->max_size_X = PetscMax(b->max_size_X,bs+((her_probl && std_probl)?0:1));
325: #else
326: b->max_size_X = PetscMax(b->max_size_X,bs);
327: #endif
329: b->max_size_V = PetscMax(b->max_size_V,mpd);
330: min_size_V = PetscMin(min_size_V,mpd-bs);
331: b->size_V = PetscMax(b->size_V,b->max_size_V+b->max_size_P+b->max_nev);
332: b->max_size_oldX = plusk;
334: /* Setup the step */
335: if (b->state >= DVD_STATE_CONF) {
336: PetscNewLog(d->eps,&data);
337: data->mpd = b->max_size_V;
338: data->min_size_V = min_size_V;
339: d->bs = bs;
340: data->plusk = plusk;
341: data->allResiduals = allResiduals;
343: d->eigr = d->eps->eigr;
344: d->eigi = d->eps->eigi;
345: d->errest = d->eps->errest;
346: PetscMalloc1(d->eps->ncv,&d->real_nR);
347: PetscMalloc1(d->eps->ncv,&d->real_nX);
348: if (plusk > 0) MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldU);
349: else data->oldU = NULL;
350: if (harm && plusk>0) MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldV);
351: else data->oldV = NULL;
353: data->old_updateV_data = d->updateV_data;
354: d->updateV_data = data;
355: data->old_isRestarting = d->isRestarting;
356: d->isRestarting = dvd_isrestarting_fullV;
357: d->updateV = dvd_updateV_extrapol;
358: d->preTestConv = dvd_updateV_testConv;
359: EPSDavidsonFLAdd(&d->startList,dvd_updateV_start);
360: EPSDavidsonFLAdd(&d->destroyList,dvd_managementV_basic_d);
361: }
362: PetscFunctionReturn(0);
363: }