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;
61: #define _Ceil(A,B) ((A)/(B)+((A)%(B)==0?0:1))
62: #define FromRealToScalar(S) ((PetscInt)_Ceil((S)*sizeof(PetscReal),sizeof(PetscScalar)))
66: PetscErrorCode dvd_managementV_basic(dvdDashboard *d, dvdBlackboard *b,
67: PetscInt bs, PetscInt mpd,
68: PetscInt min_size_V,
69: PetscInt plusk, PetscBool harm,
70: PetscBool allResiduals)
71: {
72: PetscErrorCode ierr;
73: dvdManagV_basic *data;
74: #if !defined(PETSC_USE_COMPLEX)
75: PetscBool her_probl, std_probl;
76: #endif
79: /* Setting configuration constrains */
80: #if !defined(PETSC_USE_COMPLEX)
81: /* if the last converged eigenvalue is complex its conjugate pair is also
82: converged */
83: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
84: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
85: b->max_size_X = PetscMax(b->max_size_X, bs+(her_probl && std_probl)?0:1);
86: #else
87: b->max_size_X = PetscMax(b->max_size_X, bs);
88: #endif
90: b->max_size_V = PetscMax(b->max_size_V, mpd);
91: min_size_V = PetscMin(min_size_V, mpd-bs);
92: b->max_size_auxV = PetscMax(b->max_size_auxV, 1); /* dvd_updateV_testConv */
93: b->max_size_auxS = PetscMax(b->max_size_auxS, b->max_size_V*2 /* SlepcDenseOrth */ );
94: b->size_V = PetscMax(b->size_V, b->max_size_V + b->max_size_P + b->max_nev);
95: b->own_scalars+= b->size_V*2 /* eigr, eigr */ +
96: b->size_V /* nR */ +
97: b->size_V /* nX */ +
98: b->size_V /* errest */ +
99: b->max_size_V*b->max_size_V*(harm?2:1)*(plusk>0?2:1)
100: /* MTX,MTY?,oldU,oldV? */;
101: b->max_size_oldX = plusk;
103: /* Setup the step */
104: if (b->state >= DVD_STATE_CONF) {
105: PetscMalloc(sizeof(dvdManagV_basic), &data);
106: data->mpd = b->max_size_V;
107: data->min_size_V = min_size_V;
108: d->bs = bs;
109: d->max_size_X = b->max_size_X;
110: data->plusk = plusk;
111: data->allResiduals = allResiduals;
113: d->size_real_eigr = b->size_V;
114: d->real_eigr = b->free_scalars; b->free_scalars+= b->size_V;
115: d->real_eigi = b->free_scalars; b->free_scalars+= b->size_V;
116: d->real_nR = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(b->size_V);
117: d->real_nX = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(b->size_V);
118: d->real_errest = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(b->size_V);
119: d->MTX = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
120: if (plusk > 0) {
121: data->oldU = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
122: }
123: if (harm) {
124: d->MTY = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
125: if (plusk > 0) {
126: data->oldV = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
127: }
128: } else {
129: d->MTY = PETSC_NULL;
130: data->oldV = PETSC_NULL;
131: }
133: data->old_updateV_data = d->updateV_data;
134: d->updateV_data = data;
135: data->old_isRestarting = d->isRestarting;
136: d->isRestarting = dvd_isrestarting_fullV;
137: d->updateV = dvd_updateV_extrapol;
138: d->preTestConv = dvd_updateV_testConv;
139: DVD_FL_ADD(d->startList, dvd_updateV_start);
140: DVD_FL_ADD(d->endList, dvd_updateV_conv_finish);
141: DVD_FL_ADD(d->destroyList, dvd_managementV_basic_d);
142: }
143: return(0);
144: }
148: PetscErrorCode dvd_updateV_start(dvdDashboard *d)
149: {
150: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
151: PetscInt i;
155: d->size_cX = 0;
156: d->eigr = d->ceigr = d->real_eigr;
157: d->eigi = d->ceigi = d->real_eigi;
158: #if defined(PETSC_USE_COMPLEX)
159: for(i=0; i<d->size_real_V; i++) d->eigi[i] = 0.0;
160: #endif
161: d->nR = d->real_nR;
162: for(i=0; i<d->size_real_V; i++) d->nR[i] = PETSC_MAX_REAL;
163: d->nX = d->real_nX;
164: d->errest = d->real_errest;
165: for(i=0; i<d->size_real_V; i++) d->errest[i] = PETSC_MAX_REAL;
166: data->ldoldU = 0;
167: data->oldV = PETSC_NULL;
168: d->ldMTY = 0;
169: data->size_oldU = 0;
170: d->nconv = 0;
171: d->npreconv = 0;
172: d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
173: d->size_D = 0;
175: return(0);
176: }
181: PetscBool dvd_isrestarting_fullV(dvdDashboard *d)
182: {
183: PetscBool restart;
184: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
188: restart = (d->size_V + d->max_size_X > PetscMin(data->mpd,d->max_size_V))?
189: PETSC_TRUE:PETSC_FALSE;
191: /* Check old isRestarting function */
192: if (!restart && data->old_isRestarting)
193: restart = data->old_isRestarting(d);
195: PetscFunctionReturn(restart);
196: }
200: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
201: {
202: PetscErrorCode ierr;
203: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
207: /* Restore changes in dvdDashboard */
208: d->updateV_data = data->old_updateV_data;
209:
210: /* Free local data */
211: PetscFree(data);
213: return(0);
214: }
218: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
219: {
220: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
221: PetscInt i;
222: PetscErrorCode ierr;
226: d->calcpairs_selectPairs(d, data->min_size_V);
228: /* If the subspaces doesn't need restart, add new vector */
229: if (!d->isRestarting(d)) {
230: d->size_D = 0;
231: dvd_updateV_update_gen(d);
233: /* If some vector were add, exit */
234: if (d->size_D > 0) { return(0); }
235: }
237: /* If some eigenpairs were converged, lock them */
238: if (d->npreconv > 0) {
239: i = d->npreconv;
240: dvd_updateV_conv_gen(d);
242: /* If some eigenpair was locked, exit */
243: if (i > d->npreconv) { return(0); }
244: }
246: /* Else, a restarting is performed */
247: dvd_updateV_restart_gen(d);
249: return(0);
250: }
254: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
255: {
256: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
257: PetscInt npreconv, ldpX, cMT;
258: PetscErrorCode ierr;
259: PetscScalar *pX;
260: #if !defined(PETSC_USE_COMPLEX)
261: PetscInt i, j;
262: #endif
266: npreconv = d->npreconv;
267: /* Constrains the converged pairs to nev */
268: #if !defined(PETSC_USE_COMPLEX)
269: /* Tries to maintain together conjugate eigenpairs */
270: for(i = 0;
271: (i + (d->eigi[i]!=0.0?1:0) < npreconv) && (d->nconv + i < d->nev);
272: i+= (d->eigi[i]!=0.0?2:1));
273: npreconv = i;
274: #else
275: npreconv = PetscMax(PetscMin(d->nev - d->nconv, npreconv), 0);
276: #endif
278: if (npreconv == 0) { return(0); }
279: npreconv+= d->cX_in_H;
281: #if !defined(PETSC_USE_COMPLEX)
282: /* Correct the order of the conjugate eigenpairs */
283: if (d->T) for (i=0; i<npreconv-d->cX_in_H; i++) if (d->eigi[i] != 0.0) {
284: if (d->eigi[i] < 0.0) {
285: d->eigi[i]*= -1.0;
286: d->eigi[i+1]*= -1.0;
287: for (j=0; j<d->size_H; j++) d->pX[j+(d->cX_in_H+i+1)*d->ldpX]*= -1.0;
288: for (j=0; j<d->size_H; j++) d->S[j+(d->cX_in_H+i+1)*d->ldS]*= -1.0;
289: for (j=0; j<d->size_H; j++) d->T[j+(d->cX_in_H+i+1)*d->ldT]*= -1.0;
290: }
291: i++;
292: }
293: #endif
295: /* Harmonics restarts wiht right eigenvectors, and other with
296: the left ones */
297: pX = (d->W||!d->cY||d->BcX)?d->pX:d->pY;
298: ldpX = (d->W||!d->cY||d->BcX)?d->ldpX:d->ldpY;
300: /* MTX <- [d.pX(0:npreconv-1) pX(npreconv:size_H-1)] */
301: d->ldMTX = d->ldMTY = d->size_H;
302: d->size_MT = d->size_H;
303: cMT = d->size_H - npreconv;
304: SlepcDenseCopy(d->MTX, d->ldMTX, d->pX, d->ldpX, d->size_H, npreconv);
305:
306: SlepcDenseCopy(&d->MTX[d->ldMTX*npreconv], d->ldMTX,
307: &pX[ldpX*npreconv], ldpX, d->size_H, cMT);
309: if (d->pY && d->MTY) {
310: SlepcDenseCopy(d->MTY, d->ldMTY, d->pY, d->ldpY, d->size_H,
311: d->size_H);
312: } else
313: d->MTY = PETSC_NULL;
315: /* Lock the converged pairs */
316: d->eigr+= npreconv-d->cX_in_H;
317: #if !defined(PETSC_USE_COMPLEX)
318: if (d->eigi) d->eigi+= npreconv-d->cX_in_H;
319: #endif
320: d->nconv+= npreconv-d->cX_in_H;
321: d->errest+= npreconv-d->cX_in_H;
323: /* Notify the changes in V and update the other subspaces */
324: d->V_tra_s = npreconv; d->V_tra_e = d->size_H;
325: d->V_new_s = cMT; d->V_new_e = d->V_new_s;
327: /* Remove oldU */
328: data->size_oldU = 0;
330: d->npreconv-= npreconv-d->cX_in_H;
332: return(0);
333: }
337: PetscErrorCode dvd_updateV_conv_finish(dvdDashboard *d)
338: {
339: PetscErrorCode ierr;
340: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
341: #if defined(PETSC_USE_COMPLEX)
342: PetscInt i, j;
343: PetscScalar s;
344: #endif
348: /* Some functions need the diagonal elements in cT be real */
349: #if defined(PETSC_USE_COMPLEX)
350: if (d->cT) for(i=0; i<d->nconv; i++) {
351: s = PetscConj(d->cT[d->ldcT*i+i])/PetscAbsScalar(d->cT[d->ldcT*i+i]);
352: for(j=0; j<=i; j++)
353: d->cT[d->ldcT*i+j] = PetscRealPart(d->cT[d->ldcT*i+j]*s),
354: d->cS[d->ldcS*i+j]*= s;
355: VecScale(d->cX[i], s);
356: }
357: #endif
358: d->calcpairs_selectPairs(d, data->min_size_V);
360: return(0);
361: }
362:
365: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
366: {
367: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
368: PetscInt size_plusk, size_X, i, j, ldpX, cMTX, cMTY;
369: PetscScalar *pX;
370: PetscErrorCode ierr;
374: /* Select size_X desired pairs from V */
375: size_X = PetscMin(PetscMin(data->min_size_V,
376: d->size_V ),
377: d->max_size_V );
379: /* Add plusk eigenvectors from the previous iteration */
380: size_plusk = PetscMax(0, PetscMin(PetscMin(data->plusk,
381: data->size_oldU ),
382: d->max_size_V - size_X ));
384: /* Harmonics restarts wiht right eigenvectors, and other with
385: the left ones */
386: pX = (d->W||!d->cY||d->BcX)?d->pX:d->pY;
387: ldpX = (d->W||!d->cY||d->BcX)?d->ldpX:d->ldpY;
389: /* MTX <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
390: d->ldMTX = d->size_MT = d->size_H;
391: SlepcDenseCopy(d->MTX, d->ldMTX, pX, ldpX, d->size_H, size_X);
392:
393: if (size_plusk > 0) {
394: SlepcDenseCopy(&d->MTX[d->ldMTX*size_X], d->ldMTX, data->oldU,
395: data->ldoldU, data->size_oldU, size_plusk);
396: for(i=size_X; i<size_X+size_plusk; i++)
397: for(j=data->size_oldU; j<d->size_H; j++)
398: d->MTX[j*d->ldMTX+i] = 0.0;
399: SlepcDenseOrth(d->MTX, d->ldMTX, d->size_V, size_X+size_plusk,
400: d->auxS, d->size_auxS, &cMTX);
401: } else
402: cMTX = size_X;
404: if (d->pY && d->MTY) {
405: /* MTY <- orth([pY(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
406: d->ldMTY = d->ldMTX;
407: SlepcDenseCopy(d->MTY, d->ldMTY, d->pY, d->ldpY, d->size_H, size_X);
408:
409: if (size_plusk > 0) {
410: SlepcDenseCopy(&d->MTY[d->ldMTY*size_X], d->ldMTY, data->oldV,
411: data->ldoldU, data->size_oldU, size_plusk);
412: for(i=size_X; i<size_X+size_plusk; i++)
413: for(j=data->size_oldU; j<d->size_H; j++)
414: d->MTY[j*d->ldMTY+i] = 0.0;
415: SlepcDenseOrth(d->MTY, d->ldMTY, d->size_V, size_X+size_plusk,
416: d->auxS, d->size_auxS, &cMTY);
417: cMTX = PetscMin(cMTX, cMTY);
418: }
419: }
421: /* Notify the changes in V and update the other subspaces */
422: d->V_tra_s = d->cX_in_H; d->V_tra_e = cMTX;
423: d->V_new_s = d->V_tra_e-d->cX_in_H; d->V_new_e = d->V_new_s;
425: /* Remove oldU */
426: data->size_oldU = 0;
428: /* Remove npreconv */
429: d->npreconv = 0;
430:
431: return(0);
432: }
437: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
438: {
439: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
440: PetscInt size_D;
441: PetscErrorCode ierr;
445: /* Select the desired pairs */
446: size_D = PetscMin(PetscMin(PetscMin(d->bs,
447: d->size_V ),
448: d->max_size_V-d->size_V ),
449: d->size_H );
450: if (size_D == 0) {
451: PetscPrintf(PETSC_COMM_WORLD, "MON: D:%D H:%D\n", size_D, d->size_H);
452: d->initV(d);
453: d->calcPairs(d);
454: }
456: /* Fill V with D */
457: d->improveX(d, d->V+d->size_V, d->max_size_V-d->size_V, 0, size_D, &size_D);
459: /* If D is empty, exit */
460: d->size_D = size_D;
461: if (size_D == 0) { return(0); }
463: /* Get the residual of all pairs */
464: dvd_updateV_testConv(d,size_D,size_D,data->allResiduals?d->size_V:size_D,d->auxV,d->auxS,PETSC_NULL);
466: /* Notify the changes in V */
467: d->V_tra_s = 0; d->V_tra_e = 0;
468: d->V_new_s = d->size_V; d->V_new_e = d->size_V+size_D;
470: /* Save the projected eigenvectors */
471: if (data->plusk > 0) {
472: data->ldoldU = data->size_oldU = d->size_H;
473: SlepcDenseCopy(data->oldU, data->ldoldU, d->pX, d->ldpX, d->size_H,
474: d->size_H);
475: if (d->cY) {
476: SlepcDenseCopy(data->oldV, data->ldoldU, d->pY, d->ldpY, d->size_H,
477: d->size_H);
478: }
479: }
481: return(0);
482: }
487: /* auxV: (by calcpairs_residual_eig) */
488: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d, PetscInt s, PetscInt pre,
489: PetscInt e, Vec *auxV, PetscScalar *auxS,
490: PetscInt *nConv)
491: {
492: PetscInt i,j,b;
493: PetscReal norm;
494: PetscErrorCode ierr;
495: PetscBool conv, c;
496: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
499:
500: if (nConv) *nConv = s;
501: for(i=s, conv=PETSC_TRUE;
502: (conv || data->allResiduals) && (i < e);
503: i+=b) {
504: #if !defined(PETSC_USE_COMPLEX)
505: b = d->eigi[i]!=0.0?2:1;
506: #else
507: b = 1;
508: #endif
509: if (i+b-1 >= pre) {
510: d->calcpairs_residual(d, i, i+b, auxV);
511:
512: }
513: /* Test the Schur vector */
514: for (j=0,c=PETSC_TRUE; j<b && c; j++) {
515: norm = d->nR[i+j]/d->nX[i+j];
516: c = d->testConv(d, d->eigr[i+j], d->eigi[i+j], norm, &d->errest[i+j]);
517: }
518: /* Test the eigenvector */
519: if (d->eps->trueres && conv && c) {
520: d->calcpairs_residual_eig(d,i,i+b,auxV);
521: for (j=0,c=PETSC_TRUE; j<b && c; j++) {
522: norm = d->nR[i+j]/d->nX[i+j];
523: c = d->testConv(d, d->eigr[i+j], d->eigi[i+j], norm, &d->errest[i+j]);
524: }
525: }
526: if (conv && c) { if (nConv) *nConv = i+b; }
527: else conv = PETSC_FALSE;
528: }
529: pre = PetscMax(pre, i);
531: #if !defined(PETSC_USE_COMPLEX)
532: /* Enforce converged conjugate complex eigenpairs */
533: if (nConv) {
534: for(j=0; j<*nConv; j++) if(d->eigi[j] != 0.0) j++;
535: if(j > *nConv) (*nConv)--;
536: }
537: #endif
538: for(i=pre; i<e; i++) d->errest[i] = d->nR[i] = PETSC_MAX_REAL;
539:
540: return(0);
541: }