Actual source code: dvd_updatev.c
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-2011, Universitat Politecnica de Valencia, Spain
10: This file is part of SLEPc.
11:
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: PetscErrorCode dvd_updateV_start(dvdDashboard *d);
29: PetscBool dvd_isrestarting_fullV(dvdDashboard *d);
30: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d);
31: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d);
32: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d);
33: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d);
34: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d);
35: PetscErrorCode dvd_updateV_conv_finish(dvdDashboard *d);
36: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d, PetscInt s, PetscInt pre,
37: PetscInt e, Vec *auxV, PetscScalar *auxS,
38: PetscInt *nConv);
40: typedef struct {
41: PetscInt
42: min_size_V, /* restart with this number of eigenvectors */
43: plusk, /* when restart, save plusk vectors from last iteration */
44: mpd; /* max size of the searching subspace */
45: void
46: *old_updateV_data;
47: /* old updateV data */
48: isRestarting_type
49: old_isRestarting;
50: /* old isRestarting */
51: PetscScalar
52: *oldU, /* previous projected right igenvectors */
53: *oldV; /* previous projected left eigenvectors */
54: PetscInt
55: ldoldU, /* leading dimension of oldU */
56: size_oldU; /* size of oldU */
57: PetscBool
58: allResiduals; /* if computing all the residuals */
59: } dvdManagV_basic;
63: PetscErrorCode dvd_managementV_basic(dvdDashboard *d, dvdBlackboard *b,
64: PetscInt bs, PetscInt mpd,
65: PetscInt min_size_V,
66: PetscInt plusk, PetscBool harm,
67: PetscBool allResiduals)
68: {
69: PetscErrorCode ierr;
70: dvdManagV_basic *data;
71: PetscBool her_probl, std_probl;
74: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
75: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
77: /* Setting configuration constrains */
78: #if !defined(PETSC_USE_COMPLEX)
79: /* if the last converged eigenvalue is complex its conjugate pair is also
80: converged */
81: b->max_size_X = PetscMax(b->max_size_X, bs+(her_probl && std_probl)?0:1);
82: #else
83: b->max_size_X = PetscMax(b->max_size_X, bs);
84: #endif
86: b->max_size_V = PetscMax(b->max_size_V, mpd);
87: min_size_V = PetscMin(min_size_V, mpd-bs);
88: b->max_size_auxV = PetscMax(b->max_size_auxV, 2 /* testConv */ );
89: b->max_size_auxS = PetscMax(PetscMax(b->max_size_auxS,
90: b->max_size_V*2 /* SlepcDenseOrth */ ),
91: b->max_size_V /* testConv */ );
92: b->size_V = PetscMax(b->size_V, b->max_size_V + b->max_size_P + b->max_nev);
93: b->own_scalars+= b->size_V*2 /* eigr, eigr */ +
94: b->size_V /* nR */ +
95: b->size_V /* nX */ +
96: b->size_V /* errest */ +
97: b->max_size_V*b->max_size_V*(harm?2:1)*(plusk>0?2:1)
98: /* MTX,MTY?,oldU,oldV? */;
99: b->max_size_oldX = plusk;
101: /* Setup the step */
102: if (b->state >= DVD_STATE_CONF) {
103: PetscMalloc(sizeof(dvdManagV_basic), &data);
104: data->mpd = b->max_size_V;
105: data->min_size_V = min_size_V;
106: d->bs = bs;
107: d->max_size_X = b->max_size_X;
108: data->plusk = plusk;
109: data->allResiduals = allResiduals;
111: d->size_real_eigr = b->size_V;
112: d->real_eigr = b->free_scalars; b->free_scalars+= b->size_V;
113: d->real_eigi = b->free_scalars; b->free_scalars+= b->size_V;
114: d->real_nR = (PetscReal*)b->free_scalars;
115: b->free_scalars = (PetscScalar*)(d->real_nR + b->size_V);
116: d->real_nX = (PetscReal*)b->free_scalars;
117: b->free_scalars = (PetscScalar*)(d->real_nX + b->size_V);
118: d->real_errest = (PetscReal*)b->free_scalars;
119: b->free_scalars = (PetscScalar*)(d->real_errest + b->size_V);
120: d->MTX = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
121: if (plusk > 0) {
122: data->oldU = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
123: }
124: if (harm) {
125: d->MTY = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
126: if (plusk > 0) {
127: data->oldV = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
128: }
129: } else {
130: d->MTY = PETSC_NULL;
131: data->oldV = PETSC_NULL;
132: }
134: data->old_updateV_data = d->updateV_data;
135: d->updateV_data = data;
136: data->old_isRestarting = d->isRestarting;
137: d->isRestarting = dvd_isrestarting_fullV;
138: d->updateV = dvd_updateV_extrapol;
139: d->preTestConv = dvd_updateV_testConv;
140: DVD_FL_ADD(d->startList, dvd_updateV_start);
141: DVD_FL_ADD(d->endList, dvd_updateV_conv_finish);
142: DVD_FL_ADD(d->destroyList, dvd_managementV_basic_d);
143: }
145: return(0);
146: }
150: PetscErrorCode dvd_updateV_start(dvdDashboard *d)
151: {
152: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
153: PetscInt i;
157: d->size_cX = 0;
158: d->eigr = d->ceigr = d->real_eigr;
159: d->eigi = d->ceigi = d->real_eigi;
160: #if defined(PETSC_USE_COMPLEX)
161: for(i=0; i<d->size_real_V; i++) d->eigi[i] = 0.0;
162: #endif
163: d->nR = d->real_nR;
164: for(i=0; i<d->size_real_V; i++) d->nR[i] = PETSC_MAX_REAL;
165: d->nX = d->real_nX;
166: d->errest = d->real_errest;
167: for(i=0; i<d->size_real_V; i++) d->errest[i] = PETSC_MAX_REAL;
168: data->ldoldU = 0;
169: data->oldV = PETSC_NULL;
170: d->ldMTY = 0;
171: data->size_oldU = 0;
172: d->nconv = 0;
173: d->npreconv = 0;
174: d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
175: d->size_D = 0;
177: return(0);
178: }
183: PetscBool dvd_isrestarting_fullV(dvdDashboard *d)
184: {
185: PetscBool restart;
186: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
190: restart = (d->size_V + d->max_size_X > PetscMin(data->mpd,d->max_size_V))?
191: PETSC_TRUE:PETSC_FALSE;
193: /* Check old isRestarting function */
194: if (!restart && data->old_isRestarting)
195: restart = data->old_isRestarting(d);
197: PetscFunctionReturn(restart);
198: }
202: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
203: {
204: PetscErrorCode ierr;
205: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
209: /* Restore changes in dvdDashboard */
210: d->updateV_data = data->old_updateV_data;
211:
212: /* Free local data */
213: PetscFree(data);
215: return(0);
216: }
220: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
221: {
222: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
223: PetscInt i;
224: PetscErrorCode ierr;
228: d->calcpairs_selectPairs(d, data->min_size_V);
230: /* If the subspaces doesn't need restart, add new vector */
231: if (!d->isRestarting(d)) {
232: d->size_D = 0;
233: dvd_updateV_update_gen(d);
235: /* If some vector were add, exit */
236: if (d->size_D > 0) { return(0); }
237: }
239: /* If some eigenpairs were converged, lock them */
240: if (d->npreconv > 0) {
241: i = d->npreconv;
242: dvd_updateV_conv_gen(d);
244: /* If some eigenpair was locked, exit */
245: if (i > d->npreconv) { return(0); }
246: }
248: /* Else, a restarting is performed */
249: dvd_updateV_restart_gen(d);
251: return(0);
252: }
256: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
257: {
258: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
259: PetscInt npreconv, ldpX, cMT;
260: PetscErrorCode ierr;
261: PetscScalar *pX;
262: #if !defined(PETSC_USE_COMPLEX)
263: PetscInt i, j;
264: #endif
268: npreconv = d->npreconv;
269: /* Constrains the converged pairs to nev */
270: #if !defined(PETSC_USE_COMPLEX)
271: /* Tries to maintain together conjugate eigenpairs */
272: for(i = 0;
273: (i + (d->eigi[i]!=0.0?1:0) < npreconv) && (d->nconv + i < d->nev);
274: i+= (d->eigi[i]!=0.0?2:1));
275: npreconv = i;
276: #else
277: npreconv = PetscMax(PetscMin(d->nev - d->nconv, npreconv), 0);
278: #endif
280: if (npreconv == 0) { return(0); }
281: npreconv+= d->cX_in_H;
283: #if !defined(PETSC_USE_COMPLEX)
284: /* Correct the order of the conjugate eigenpairs */
285: if (d->T) for (i=0; i<npreconv-d->cX_in_H; i++) if (d->eigi[i] != 0.0) {
286: if (d->eigi[i] < 0.0) {
287: d->eigi[i]*= -1.0;
288: d->eigi[i+1]*= -1.0;
289: for (j=0; j<d->size_H; j++) d->pX[j+(d->cX_in_H+i+1)*d->ldpX]*= -1.0;
290: for (j=0; j<d->size_H; j++) d->S[j+(d->cX_in_H+i+1)*d->ldS]*= -1.0;
291: for (j=0; j<d->size_H; j++) d->T[j+(d->cX_in_H+i+1)*d->ldT]*= -1.0;
292: }
293: i++;
294: }
295: #endif
297: /* Harmonics restarts wiht right eigenvectors, and other with
298: the left ones */
299: pX = (d->W||!d->cY||d->BcX)?d->pX:d->pY;
300: ldpX = (d->W||!d->cY||d->BcX)?d->ldpX:d->ldpY;
302: /* MTX <- [d.pX(0:npreconv-1) pX(npreconv:size_H-1)] */
303: d->ldMTX = d->ldMTY = d->size_H;
304: d->size_MT = d->size_H;
305: cMT = d->size_H - npreconv;
306: SlepcDenseCopy(d->MTX, d->ldMTX, d->pX, d->ldpX, d->size_H, npreconv);
307:
308: SlepcDenseCopy(&d->MTX[d->ldMTX*npreconv], d->ldMTX,
309: &pX[ldpX*npreconv], ldpX, d->size_H, cMT);
311: if (d->pY && d->MTY) {
312: SlepcDenseCopy(d->MTY, d->ldMTY, d->pY, d->ldpY, d->size_H,
313: d->size_H);
314: } else
315: d->MTY = PETSC_NULL;
317: /* Lock the converged pairs */
318: d->eigr+= npreconv-d->cX_in_H;
319: #if !defined(PETSC_USE_COMPLEX)
320: if (d->eigi) d->eigi+= npreconv-d->cX_in_H;
321: #endif
322: d->nconv+= npreconv-d->cX_in_H;
323: d->errest+= npreconv-d->cX_in_H;
325: /* Notify the changes in V and update the other subspaces */
326: d->V_tra_s = npreconv; d->V_tra_e = d->size_H;
327: d->V_new_s = cMT; d->V_new_e = d->V_new_s;
329: /* Remove oldU */
330: data->size_oldU = 0;
332: d->npreconv-= npreconv-d->cX_in_H;
334: return(0);
335: }
339: PetscErrorCode dvd_updateV_conv_finish(dvdDashboard *d)
340: {
341: PetscErrorCode ierr;
342: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
343: #if defined(PETSC_USE_COMPLEX)
344: PetscInt i, j;
345: PetscScalar s;
346: #endif
350: /* Some functions need the diagonal elements in cT be real */
351: #if defined(PETSC_USE_COMPLEX)
352: if (d->cT) for(i=0; i<d->nconv; i++) {
353: s = PetscConj(d->cT[d->ldcT*i+i])/PetscAbsScalar(d->cT[d->ldcT*i+i]);
354: for(j=0; j<=i; j++)
355: d->cT[d->ldcT*i+j] = PetscRealPart(d->cT[d->ldcT*i+j]*s),
356: d->cS[d->ldcS*i+j]*= s;
357: VecScale(d->cX[i], s);
358: }
359: #endif
360: d->calcpairs_selectPairs(d, data->min_size_V);
362: return(0);
363: }
364:
367: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
368: {
369: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
370: PetscInt size_plusk, size_X, i, j, ldpX, cMTX, cMTY;
371: PetscScalar *pX;
372: PetscErrorCode ierr;
376: /* Select size_X desired pairs from V */
377: size_X = PetscMin(PetscMin(data->min_size_V,
378: d->size_V ),
379: d->max_size_V );
381: /* Add plusk eigenvectors from the previous iteration */
382: size_plusk = PetscMax(0, PetscMin(PetscMin(data->plusk,
383: data->size_oldU ),
384: d->max_size_V - size_X ));
386: /* Harmonics restarts wiht right eigenvectors, and other with
387: the left ones */
388: pX = (d->W||!d->cY||d->BcX)?d->pX:d->pY;
389: ldpX = (d->W||!d->cY||d->BcX)?d->ldpX:d->ldpY;
391: /* MTX <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
392: d->ldMTX = d->size_MT = d->size_H;
393: SlepcDenseCopy(d->MTX, d->ldMTX, pX, ldpX, d->size_H, size_X);
394:
395: if (size_plusk > 0) {
396: SlepcDenseCopy(&d->MTX[d->ldMTX*size_X], d->ldMTX, data->oldU,
397: data->ldoldU, data->size_oldU, size_plusk);
398: for(i=size_X; i<size_X+size_plusk; i++)
399: for(j=data->size_oldU; j<d->size_H; j++)
400: d->MTX[j*d->ldMTX+i] = 0.0;
401: SlepcDenseOrth(d->MTX, d->ldMTX, d->size_V, size_X+size_plusk,
402: d->auxS, d->size_auxS, &cMTX);
403: } else
404: cMTX = size_X;
406: if (d->pY && d->MTY) {
407: /* MTY <- orth([pY(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
408: d->ldMTY = d->ldMTX;
409: SlepcDenseCopy(d->MTY, d->ldMTY, d->pY, d->ldpY, d->size_H, size_X);
410:
411: if (size_plusk > 0) {
412: SlepcDenseCopy(&d->MTY[d->ldMTY*size_X], d->ldMTY, data->oldV,
413: data->ldoldU, data->size_oldU, size_plusk);
414: for(i=size_X; i<size_X+size_plusk; i++)
415: for(j=data->size_oldU; j<d->size_H; j++)
416: d->MTY[j*d->ldMTY+i] = 0.0;
417: SlepcDenseOrth(d->MTY, d->ldMTY, d->size_V, size_X+size_plusk,
418: d->auxS, d->size_auxS, &cMTY);
419: cMTX = PetscMin(cMTX, cMTY);
420: }
421: }
423: /* Notify the changes in V and update the other subspaces */
424: d->V_tra_s = d->cX_in_H; d->V_tra_e = cMTX;
425: d->V_new_s = d->V_tra_e-d->cX_in_H; d->V_new_e = d->V_tra_e;
427: /* Remove oldU */
428: data->size_oldU = 0;
430: /* Remove npreconv */
431: d->npreconv = 0;
432:
433: return(0);
434: }
439: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
440: {
441: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
442: PetscInt size_D;
443: PetscErrorCode ierr;
447: /* Select the desired pairs */
448: size_D = PetscMin(PetscMin(PetscMin(d->bs,
449: d->size_V ),
450: d->max_size_V-d->size_V ),
451: d->size_H );
452: if (size_D == 0) {
453: PetscPrintf(PETSC_COMM_WORLD, "MON: D:%D H:%D\n", size_D, d->size_H);
454: d->initV(d);
455: d->calcPairs(d);
456: }
458: /* Fill V with D */
459: d->improveX(d, d->V+d->size_V, d->max_size_V-d->size_V, 0, size_D, &size_D);
461: /* If D is empty, exit */
462: d->size_D = size_D;
463: if (size_D == 0) { return(0); }
465: /* Get the converged pairs */
466: dvd_updateV_testConv(d, 0, size_D,
467: data->allResiduals?d->size_V:size_D, d->auxV, d->auxS,
468: &d->npreconv);
470: /* Notify the changes in V */
471: d->V_tra_s = 0; d->V_tra_e = 0;
472: d->V_new_s = d->size_V; d->V_new_e = d->size_V+size_D;
474: /* Save the projected eigenvectors */
475: if (data->plusk > 0) {
476: data->ldoldU = data->size_oldU = d->size_H;
477: SlepcDenseCopy(data->oldU, data->ldoldU, d->pX, d->ldpX, d->size_H,
478: d->size_H);
479: if (d->cY) {
480: SlepcDenseCopy(data->oldV, data->ldoldU, d->pY, d->ldpY, d->size_H,
481: d->size_H);
482: }
483: }
485: return(0);
486: }
491: /* auxS: max_size_V */
492: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d, PetscInt s, PetscInt pre,
493: PetscInt e, Vec *auxV, PetscScalar *auxS,
494: PetscInt *nConv)
495: {
496: PetscInt i;
497: #if !defined(PETSC_USE_COMPLEX)
498: PetscInt j;
499: #endif
500: PetscReal norm;
501: PetscErrorCode ierr;
502: PetscBool conv, c;
503: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
506:
507: *nConv = s;
508: for(i=s, conv=PETSC_TRUE;
509: (conv || data->allResiduals) && (i < e);
510: i++) {
511: if (i >= pre) {
512: d->calcpairs_residual(d, i, i+1, auxV);
513:
514: }
515: norm = d->nR[i]/d->nX[i];
516: c = d->testConv(d, d->eigr[i], d->eigi[i], norm, &d->errest[i]);
517: if (conv && c) *nConv = i+1;
518: else conv = PETSC_FALSE;
519: }
520: pre = PetscMax(pre, i);
522: #if !defined(PETSC_USE_COMPLEX)
523: /* Enforce converged conjugate conjugate complex eigenpairs */
524: for(j=0; j<*nConv; j++) if(d->eigi[j] != 0.0) j++;
525: if(j > *nConv) (*nConv)--;
526: #endif
527: for(i=pre; i<e; i++) d->errest[i] = d->nR[i] = PETSC_MAX_REAL;
528:
529: return(0);
530: }