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-2013, 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
27: #include <slepc-private/dsimpl.h> /*I "slepcds.h" I*/
29: PetscErrorCode dvd_updateV_start(dvdDashboard *d);
30: PetscBool dvd_isrestarting_fullV(dvdDashboard *d);
31: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d);
32: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d);
33: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d);
34: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d);
35: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d);
36: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d,PetscInt s,PetscInt pre,PetscInt e,Vec *auxV,PetscScalar *auxS,PetscInt *nConv);
38: typedef struct {
39: PetscInt
40: min_size_V, /* restart with this number of eigenvectors */
41: plusk, /* when restart, save plusk vectors from last iteration */
42: mpd; /* max size of the searching subspace */
43: void
44: *old_updateV_data;
45: /* old updateV data */
46: isRestarting_type
47: old_isRestarting;
48: /* old isRestarting */
49: PetscScalar
50: *oldU, /* previous projected right igenvectors */
51: *oldV; /* previous projected left eigenvectors */
52: PetscInt
53: ldoldU, /* leading dimension of oldU */
54: size_oldU; /* size of oldU */
55: PetscBool
56: allResiduals; /* if computing all the residuals */
57: } dvdManagV_basic;
62: PetscErrorCode dvd_managementV_basic(dvdDashboard *d,dvdBlackboard *b,PetscInt bs,PetscInt mpd,PetscInt min_size_V,PetscInt plusk,PetscBool harm,PetscBool allResiduals)
63: {
64: PetscErrorCode ierr;
65: dvdManagV_basic *data;
66: #if !defined(PETSC_USE_COMPLEX)
67: PetscBool her_probl, std_probl;
68: #endif
71: /* Setting configuration constrains */
72: #if !defined(PETSC_USE_COMPLEX)
73: /* if the last converged eigenvalue is complex its conjugate pair is also
74: converged */
75: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
76: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
77: b->max_size_X = PetscMax(b->max_size_X, bs+(her_probl && std_probl)?0:1);
78: #else
79: b->max_size_X = PetscMax(b->max_size_X, bs);
80: #endif
82: b->max_size_V = PetscMax(b->max_size_V, mpd);
83: min_size_V = PetscMin(min_size_V, mpd-bs);
84: b->max_size_auxV = PetscMax(b->max_size_auxV, 1); /* dvd_updateV_testConv */
85: b->size_V = PetscMax(b->size_V, b->max_size_V + b->max_size_P + b->max_nev);
86: b->own_scalars+= b->size_V*2 /* eigr, eigr */ +
87: b->size_V /* nR */ +
88: b->size_V /* nX */ +
89: b->size_V /* errest */ +
90: b->max_size_V*b->max_size_V*(harm?2:1)*(plusk>0?1:0)
91: /* oldU,oldV? */;
92: b->max_size_oldX = plusk;
94: /* Setup the step */
95: if (b->state >= DVD_STATE_CONF) {
96: PetscMalloc(sizeof(dvdManagV_basic),&data);
97: PetscLogObjectMemory(d->eps,sizeof(dvdManagV_basic));
98: data->mpd = b->max_size_V;
99: data->min_size_V = min_size_V;
100: d->bs = bs;
101: d->max_size_X = b->max_size_X;
102: data->plusk = plusk;
103: data->allResiduals = allResiduals;
105: d->size_real_eigr = b->size_V;
106: d->real_eigr = b->free_scalars; b->free_scalars+= b->size_V;
107: d->real_eigi = b->free_scalars; b->free_scalars+= b->size_V;
108: d->real_nR = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(b->size_V);
109: d->real_nX = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(b->size_V);
110: d->real_errest = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(b->size_V);
111: if (plusk > 0) {
112: data->oldU = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
113: }
114: if (harm) {
115: if (plusk > 0) {
116: data->oldV = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
117: }
118: } else {
119: data->oldV = NULL;
120: }
122: data->old_updateV_data = d->updateV_data;
123: d->updateV_data = data;
124: data->old_isRestarting = d->isRestarting;
125: d->isRestarting = dvd_isrestarting_fullV;
126: d->updateV = dvd_updateV_extrapol;
127: d->preTestConv = dvd_updateV_testConv;
128: DVD_FL_ADD(d->startList, dvd_updateV_start);
129: DVD_FL_ADD(d->destroyList, dvd_managementV_basic_d);
130: }
131: return(0);
132: }
136: PetscErrorCode dvd_updateV_start(dvdDashboard *d)
137: {
138: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
139: PetscInt i;
142: d->size_cX = 0;
143: d->eigr = d->ceigr = d->real_eigr;
144: d->eigi = d->ceigi = d->real_eigi;
145: for (i=0;i<d->size_real_V;i++) d->eigi[i] = 0.0;
146: d->nR = d->real_nR;
147: for (i=0;i<d->size_real_V;i++) d->nR[i] = PETSC_MAX_REAL;
148: d->nX = d->real_nX;
149: d->errest = d->real_errest;
150: for (i=0;i<d->size_real_V;i++) d->errest[i] = PETSC_MAX_REAL;
151: data->ldoldU = 0;
152: data->oldV = NULL;
153: data->size_oldU = 0;
154: d->nconv = 0;
155: d->npreconv = 0;
156: d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
157: d->size_D = 0;
158: return(0);
159: }
163: PetscBool dvd_isrestarting_fullV(dvdDashboard *d)
164: {
165: PetscBool restart;
166: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
169: restart = (d->size_V + d->max_size_X > PetscMin(data->mpd,d->max_size_V))?
170: PETSC_TRUE:PETSC_FALSE;
172: /* Check old isRestarting function */
173: if (!restart && data->old_isRestarting)
174: restart = data->old_isRestarting(d);
175: PetscFunctionReturn(restart);
176: }
180: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
181: {
182: PetscErrorCode ierr;
183: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
186: /* Restore changes in dvdDashboard */
187: d->updateV_data = data->old_updateV_data;
189: /* Free local data */
190: PetscFree(data);
191: return(0);
192: }
196: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
197: {
198: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
199: PetscInt i;
200: PetscErrorCode ierr;
203: d->calcpairs_selectPairs(d, data->min_size_V);
205: /* If the subspaces doesn't need restart, add new vector */
206: if (!d->isRestarting(d)) {
207: d->size_D = 0;
208: dvd_updateV_update_gen(d);
210: /* If some vector were add, exit */
211: if (d->size_D > 0) return(0);
212: }
214: /* If some eigenpairs were converged, lock them */
215: if (d->npreconv > 0) {
216: i = d->npreconv;
217: dvd_updateV_conv_gen(d);
219: /* If some eigenpair was locked, exit */
220: if (i > d->npreconv) return(0);
221: }
223: /* Else, a restarting is performed */
224: dvd_updateV_restart_gen(d);
225: return(0);
226: }
230: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
231: {
232: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
233: PetscInt npreconv,ld,cMT,cMTX;
234: PetscErrorCode ierr;
235: PetscScalar *pQ,*pZ;
236: #if !defined(PETSC_USE_COMPLEX)
237: PetscInt i;
238: #endif
241: npreconv = d->npreconv;
242: /* Constrains the converged pairs to nev */
243: #if !defined(PETSC_USE_COMPLEX)
244: /* Tries to maintain together conjugate eigenpairs */
245: 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));
246: npreconv = i;
247: #else
248: npreconv = PetscMax(PetscMin(d->nev - d->nconv, npreconv), 0);
249: #endif
250: /* Quick exit */
251: if (npreconv == 0) return(0);
253: npreconv+= d->cX_in_H;
254: DSGetLeadingDimension(d->ps,&ld);
255: d->size_MT = d->size_H;
256: cMT = d->size_H - npreconv;
257: /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
258: If the problem is standard or hermitian, left and right vectors are the same */
259: if (!(d->W||!d->cY||d->BcX||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
260: /* ps.Q <- [ps.Q(0:npreconv-1) ps.Z(npreconv:size_H-1)] */
261: DSGetArray(d->ps,DS_MAT_Q,&pQ);
262: DSGetArray(d->ps,DS_MAT_Z,&pZ);
263: SlepcDenseCopy(&pQ[ld*npreconv],ld,&pZ[ld*npreconv],ld,d->size_H,cMT);
264: DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
265: DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
266: }
267: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
268: DSPseudoOrthogonalize(d->ps,DS_MAT_Q,d->size_H,d->nBV-d->cX_in_H,&cMTX,d->nBpX);
269: } else {
270: DSOrthogonalize(d->ps,DS_MAT_Q,d->size_H,&cMTX);
271: }
272: cMT = cMTX - npreconv;
274: if (d->W) {
275: DSOrthogonalize(d->ps,DS_MAT_Z,d->size_H,&cMTX);
276: cMT = PetscMin(cMT,cMTX - npreconv);
277: }
279: /* Lock the converged pairs */
280: d->eigr+= npreconv-d->cX_in_H;
281: #if !defined(PETSC_USE_COMPLEX)
282: if (d->eigi) d->eigi+= npreconv-d->cX_in_H;
283: #endif
284: d->nconv+= npreconv-d->cX_in_H;
285: d->errest+= npreconv-d->cX_in_H;
286: /* Notify the changes in V and update the other subspaces */
287: d->V_tra_s = npreconv; d->V_tra_e = d->size_H;
288: d->V_new_s = cMT; d->V_new_e = d->V_new_s;
289: /* Remove oldU */
290: data->size_oldU = 0;
292: d->npreconv-= npreconv-d->cX_in_H;
293: return(0);
294: }
298: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
299: {
300: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
301: PetscInt size_plusk,size_X,i,j,ld,cMTX,cMTY;
302: PetscScalar *pQ,*pZ;
303: PetscErrorCode ierr;
306: /* Select size_X desired pairs from V */
307: size_X = PetscMin(PetscMin(data->min_size_V,
308: d->size_V),
309: d->max_size_V);
311: /* Add plusk eigenvectors from the previous iteration */
312: size_plusk = PetscMax(0, PetscMin(PetscMin(data->plusk,
313: data->size_oldU),
314: d->max_size_V - size_X));
316: DSGetLeadingDimension(d->ps,&ld);
317: d->size_MT = d->size_H;
318: /* ps.Q <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
319: /* Harmonics restarts wiht right eigenvectors, and other with the left ones.
320: If the problem is standard or hermitian, left and right vectors are the same */
321: DSGetArray(d->ps,DS_MAT_Q,&pQ);
322: if (!(d->W||!d->cY||d->BcX||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
323: DSGetArray(d->ps,DS_MAT_Z,&pZ);
324: SlepcDenseCopy(pQ,ld,pZ,ld,d->size_H,size_X);
325: DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
326: }
327: 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");
328: if (size_plusk > 0) {
329: SlepcDenseCopy(&pQ[ld*size_X],ld,data->oldU,data->ldoldU,data->size_oldU,size_plusk);
330: for (i=size_X;i<size_X+size_plusk;i++) {
331: for (j=data->size_oldU;j<d->size_H;j++) {
332: pQ[j*ld+i] = 0.0;
333: }
334: }
335: }
336: DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
337: if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) {
338: DSPseudoOrthogonalize(d->ps,DS_MAT_Q,size_X,d->nBV-d->cX_in_H,&cMTX,d->nBpX);
339: } else {
340: DSOrthogonalize(d->ps,DS_MAT_Q,size_X+size_plusk,&cMTX);
341: }
343: if (d->W && size_plusk > 0) {
344: /* ps.Z <- orth([ps.Z(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
345: DSGetArray(d->ps,DS_MAT_Z,&pZ);
346: SlepcDenseCopy(&pZ[ld*size_X],ld,data->oldV,data->ldoldU,data->size_oldU,size_plusk);
347: for(i=size_X; i<size_X+size_plusk; i++) {
348: for(j=data->size_oldU; j<d->size_H; j++) {
349: pZ[j*ld+i] = 0.0;
350: }
351: }
352: DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
353: DSOrthogonalize(d->ps,DS_MAT_Z,size_X+size_plusk,&cMTY);
354: cMTX = PetscMin(cMTX, cMTY);
355: }
357: /* Notify the changes in V and update the other subspaces */
358: d->V_tra_s = d->cX_in_H; d->V_tra_e = cMTX;
359: d->V_new_s = d->V_tra_e-d->cX_in_H; d->V_new_e = d->V_new_s;
361: /* Remove oldU */
362: data->size_oldU = 0;
364: /* Remove npreconv */
365: d->npreconv = 0;
366: return(0);
367: }
371: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
372: {
373: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
374: PetscInt size_D,ld,s;
375: PetscScalar *pQ,*pZ;
376: PetscErrorCode ierr;
379: /* Select the desired pairs */
380: size_D = PetscMin(PetscMin(PetscMin(d->bs,
381: d->size_V),
382: d->max_size_V-d->size_V),
383: d->size_H);
384: if (size_D == 0) {
385: PetscInfo2(d->eps, "MON: D:%D H:%D\n", size_D, d->size_H);
386: d->initV(d);
387: d->calcPairs(d);
388: }
390: /* Fill V with D */
391: d->improveX(d, d->V+d->size_V, d->max_size_V-d->size_V, 0, size_D, &size_D);
393: /* If D is empty, exit */
394: d->size_D = size_D;
395: if (size_D == 0) return(0);
397: /* Get the residual of all pairs */
398: #if !defined(PETSC_USE_COMPLEX)
399: s = d->eigi[0]!=0.0?2:1;
400: #else
401: s = 1;
402: #endif
403: dvd_updateV_testConv(d,s,s,data->allResiduals?d->size_V:size_D,d->auxV,d->auxS,NULL);
405: /* Notify the changes in V */
406: d->V_tra_s = 0; d->V_tra_e = 0;
407: d->V_new_s = d->size_V; d->V_new_e = d->size_V+size_D;
409: /* Save the projected eigenvectors */
410: if (data->plusk > 0) {
411: data->ldoldU = data->size_oldU = d->size_H;
412: DSGetLeadingDimension(d->ps,&ld);
413: DSGetArray(d->ps,DS_MAT_Q,&pQ);
414: SlepcDenseCopy(data->oldU,data->ldoldU,pQ,ld,d->size_H,d->size_H);
415: DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
416: if (d->cY) {
417: DSGetArray(d->ps,DS_MAT_Z,&pZ);
418: SlepcDenseCopy(data->oldV,data->ldoldU,pZ,ld,d->size_H,d->size_H);
419: DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
420: }
421: }
422: return(0);
423: }
427: /* auxV: (by calcpairs_residual_eig) */
428: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d,PetscInt s,PetscInt pre,PetscInt e,Vec *auxV,PetscScalar *auxS,PetscInt *nConv)
429: {
430: PetscInt i,j,b;
431: PetscReal norm;
432: PetscErrorCode ierr;
433: PetscBool conv, c;
434: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
437: if (nConv) *nConv = s;
438: for (i=s, conv=PETSC_TRUE;
439: (conv || data->allResiduals) && (i < e);
440: i+=b) {
441: #if !defined(PETSC_USE_COMPLEX)
442: b = d->eigi[i]!=0.0?2:1;
443: #else
444: b = 1;
445: #endif
446: if (i+b-1 >= pre) {
447: d->calcpairs_residual(d, i, i+b, auxV);
448: }
449: /* Test the Schur vector */
450: for (j=0,c=PETSC_TRUE; j<b && c; j++) {
451: norm = d->nR[i+j]/d->nX[i+j];
452: c = d->testConv(d, d->eigr[i+j], d->eigi[i+j], norm, &d->errest[i+j]);
453: }
454: /* Test the eigenvector */
455: if (d->eps->trueres && conv && c) {
456: d->calcpairs_residual_eig(d,i,i+b,auxV);
457: for (j=0,c=PETSC_TRUE; j<b && c; j++) {
458: norm = d->nR[i+j]/d->nX[i+j];
459: c = d->testConv(d, d->eigr[i+j], d->eigi[i+j], norm, &d->errest[i+j]);
460: }
461: }
462: if (conv && c) { if (nConv) *nConv = i+b; }
463: else conv = PETSC_FALSE;
464: }
465: pre = PetscMax(pre, i);
467: #if !defined(PETSC_USE_COMPLEX)
468: /* Enforce converged conjugate complex eigenpairs */
469: if (nConv) {
470: for (j=0;j<*nConv;j++) if (d->eigi[j] != 0.0) j++;
471: if (j>*nConv) (*nConv)--;
472: }
473: #endif
474: for (i=pre;i<e;i++) d->errest[i] = d->nR[i] = PETSC_MAX_REAL;
475: return(0);
476: }