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-2010, Universidad 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: PetscTruth dvd_isrestarting_fullV(dvdDashboard *d);
29: PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d);
30: PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d);
31: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d);
32: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d);
33: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d);
34: PetscErrorCode dvd_updateV_conv_finish(dvdDashboard *d);
35: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d, PetscInt s, PetscInt pre,
36: PetscInt e, Vec *auxV, PetscScalar *auxS,
37: PetscInt *nConv);
38: PetscErrorCode dvd_updateV_restartV_aux(Vec *V, PetscInt size_V,
39: PetscScalar *U, PetscInt ldU,
40: PetscScalar *pX, PetscInt ldpX,
41: PetscInt cpX, PetscScalar *oldpX,
42: PetscInt ldoldpX, PetscInt roldpX,
43: PetscInt coldpX, PetscScalar *auxS,
44: PetscInt size_auxS,
45: PetscInt *new_size_V);
46: PetscErrorCode dvd_updateV_YtWx(PetscScalar *S, PetscInt ldS,
47: Vec *Y, PetscInt cY, Vec *y, PetscInt cy,
48: Vec *W, PetscInt cW, PetscScalar *x,
49: PetscInt ldx,
50: PetscInt rx, PetscInt cx, Vec *auxV);
51: typedef struct {
52: PetscInt bs, /* common number of approximated eigenpairs obtained */
53: real_max_size_V,
54: /* real max size of V */
55: min_size_V, /* restart with this number of eigenvectors */
56: plusk, /* when restart, save plusk vectors from last iteration */
57: mpd; /* max size of the searching subspace */
58: Vec *real_V, /* real start vectors V */
59: *new_cY; /* new left converged eigenvectors from the last iter */
60: void
61: *old_updateV_data;
62: /* old updateV data */
63: isRestarting_type
64: old_isRestarting;
65: /* old isRestarting */
66: PetscScalar
67: *oldU, /* previous projected right igenvectors */
68: *oldV; /* previous projected left eigenvectors */
69: PetscInt
70: ldoldU, /* leading dimension of oldU */
71: size_oldU, /* size of oldU */
72: size_new_cY; /* size of new_cY */
73: PetscTruth
74: allResiduals; /* if computing all the residuals */
75: } dvdManagV_basic;
79: PetscErrorCode dvd_managementV_basic(dvdDashboard *d, dvdBlackboard *b,
80: PetscInt bs, PetscInt max_size_V,
81: PetscInt mpd, PetscInt min_size_V,
82: PetscInt plusk, PetscTruth harm,
83: PetscTruth allResiduals)
84: {
85: PetscErrorCode ierr;
86: dvdManagV_basic *data;
87: PetscInt i;
91: /* Setting configuration constrains */
92: mpd = PetscMin(mpd, max_size_V);
93: min_size_V = PetscMin(min_size_V, mpd-bs);
94: b->max_size_X = PetscMax(b->max_size_X, PetscMax(bs, min_size_V));
95: b->max_size_auxV = PetscMax(PetscMax(b->max_size_auxV,
96: b->max_size_X /* updateV_conv_gen */ ),
97: 2 /* testConv */ );
98: b->max_size_auxS = PetscMax(PetscMax(PetscMax(b->max_size_auxS,
99: max_size_V*b->max_size_X /* YtWx */ ),
100: max_size_V*2 /* SlepcDenseOrth */ ),
101: max_size_V*b->max_size_X /* testConv:res_0 */ );
102: b->max_size_V = mpd;
103: b->size_V = max_size_V;
104: b->own_vecs+= max_size_V*(harm?2:1); /* V, W? */
105: b->own_scalars+= max_size_V*2 /* eigr, eigr */ +
106: max_size_V /* nR */ +
107: max_size_V /* nX */ +
108: max_size_V /* errest */ +
109: 2*b->max_size_V*b->max_size_V*(harm?2:1)
110: /* MTX,MTY?,oldU,oldV? */;
111: // b->max_size_oldX = plusk;
113: /* Setup the step */
114: if (b->state >= DVD_STATE_CONF) {
115: PetscMalloc(sizeof(dvdManagV_basic), &data);
116: data->real_max_size_V = max_size_V;
117: data->min_size_V = min_size_V;
118: data->real_V = b->free_vecs; b->free_vecs+= max_size_V;
119: data->bs = bs;
120: data->plusk = plusk;
121: data->new_cY = PETSC_NULL;
122: data->size_new_cY = 0;
123: data->mpd = mpd;
124: data->allResiduals = allResiduals;
126: d->V = data->real_V;
127: d->max_size_V = data->real_max_size_V;
128: d->cX = data->real_V;
129: d->eigr = b->free_scalars; b->free_scalars+= max_size_V;
130: d->eigi = b->free_scalars; b->free_scalars+= max_size_V;
131: #ifdef PETSC_USE_COMPLEX
132: for(i=0; i<max_size_V; i++) d->eigi[i] = 0.0;
133: #endif
134: d->nR = (PetscReal*)b->free_scalars;
135: b->free_scalars = (PetscScalar*)(d->nR + max_size_V);
136: for(i=0; i<max_size_V; i++) d->nR[i] = PETSC_MAX;
137: d->nX = (PetscReal*)b->free_scalars;
138: b->free_scalars = (PetscScalar*)(d->nX + max_size_V);
139: d->errest = (PetscReal*)b->free_scalars;
140: b->free_scalars = (PetscScalar*)(d->errest + max_size_V);
141: d->ceigr = d->eigr;
142: d->ceigi = d->eigi;
143: d->MTX = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
144: data->oldU = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
145: data->ldoldU = 0;
146: data->oldV = PETSC_NULL;
147: d->W = PETSC_NULL;
148: d->MTY = PETSC_NULL;
149: d->ldMTY = 0;
150: if (harm) {
151: d->W = b->free_vecs; b->free_vecs+= max_size_V;
152: d->MTY = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
153: data->oldV = b->free_scalars; b->free_scalars+= b->max_size_V*b->max_size_V;
154: }
156: data->size_oldU = 0;
157: data->old_updateV_data = d->updateV_data;
158: d->updateV_data = data;
159: data->old_isRestarting = d->isRestarting;
160: d->isRestarting = dvd_isrestarting_fullV;
161: d->updateV = dvd_updateV_extrapol;
162: d->preTestConv = dvd_updateV_testConv;
163: DVD_FL_ADD(d->endList, dvd_updateV_conv_finish);
164: DVD_FL_ADD(d->destroyList, dvd_managementV_basic_d);
165: }
167: return(0);
168: }
172: PetscTruth dvd_isrestarting_fullV(dvdDashboard *d)
173: {
174: PetscTruth restart;
175: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
179: /* Take into account the possibility of conjugate eigenpairs */
180: #if defined(PETSC_USE_COMPLEX)
181: #define ONE 0
182: #else
183: #define ONE 1
184: #endif
186: restart = (d->size_V + data->bs + ONE > PetscMin(data->mpd,d->max_size_V))?
187: PETSC_TRUE:PETSC_FALSE;
189: #undef ONE
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: /// Temporal! Copy the converged left eigenvectors to cY
227: if (data->size_new_cY > 0) {
228: for (i=0; i<data->size_new_cY; i++) {
229: VecCopy(data->new_cY[i], d->cY[d->size_cY+i]);
230: }
231: d->size_cY+= data->size_new_cY;
232: data->size_new_cY = 0;
233: }
235: d->calcpairs_selectPairs(d, data->min_size_V);
237: /* If the subspaces doesn't need restart, add new vector */
238: if (!d->isRestarting(d)) {
239: i = d->size_V;
240: dvd_updateV_update_gen(d);
242: /* If some vector were add, exit */
243: if (i < d->size_V) { return(0); }
244: }
246: /* If some eigenpairs were converged, lock them */
247: if (d->npreconv > 0) {
248: i = d->npreconv;
249: dvd_updateV_conv_gen(d);
251: /* If some eigenpair was locked, exit */
252: if (i > d->npreconv) { return(0); }
253: }
255: /* Else, a restarting is performed */
256: dvd_updateV_restart_gen(d);
258: return(0);
259: }
263: PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
264: {
265: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
266: PetscInt i, j, npreconv, ldpX, inc_V, cMT, size_cX, size_cy;
267: PetscErrorCode ierr;
268: PetscScalar *pX;
269: PetscReal norm;
270: Vec *new_cY=0, *cX, *cy;
271: PetscTruth lindep;
275: /* If the left subspace is present, constrains the converged pairs to the
276: number of free vectors in V */
277: if (d->cY && d->pY)
278: npreconv = PetscMin(d->max_size_V-d->size_V, d->npreconv);
279: else
280: npreconv = d->npreconv;
282: /* Constrains the converged pairs to nev */
283: #ifndef PETSC_USE_COMPLEX
284: /* Tries to maintain together conjugate eigenpairs */
285: for(i = 0;
286: (i + (d->eigi[i]!=0.0?1:0) < npreconv) && (d->nconv + i < d->nev);
287: i+= (d->eigi[i]!=0.0?2:1));
288: npreconv = i;
289: #else
290: npreconv = PetscMax(PetscMin(d->nev - d->nconv, npreconv), 0);
291: #endif
293: if (d->npreconv == 0) { return(0); }
295: /* f.cY <- [f.cY f.V*f.pY(0:npreconv-1)] */
296: if (d->cY && d->pY) {
297: new_cY = &d->V[d->size_V];
298: SlepcUpdateVectorsZ(new_cY, 0.0, 1.0, d->W?d->W:d->V, d->size_V,
299: d->pY, d->ldpY, d->size_H, npreconv);
300:
302: /// Temporal! Postpone the copy of new_cY to cY
303: data->new_cY = new_cY;
304: data->size_new_cY = npreconv;
305: }
307: #if !defined(PETSC_USE_COMPLEX)
308: /* Correct the order of the conjugate eigenpairs */
309: if (d->T) for (i=0; i<npreconv; i++) if (d->eigi[i] != 0.0) {
310: if (d->eigi[i] < 0.0) {
311: d->eigi[i]*= -1.0;
312: d->eigi[i+1]*= -1.0;
313: for (j=0; j<d->size_H; j++) d->pX[j+(i+1)*d->ldpX]*= -1.0;
314: for (j=0; j<d->size_H; j++) d->S[j+(i+1)*d->ldS]*= -1.0;
315: for (j=0; j<d->size_H; j++) d->T[j+(i+1)*d->ldT]*= -1.0;
316: }
317: i++;
318: }
319: #endif
321: /* BcX <- orth(BcX, B*cX),
322: auxV = B*X, being X the last converged eigenvectors */
323: if (d->BcX) for(i=0; i<npreconv; i++) {
324: /* BcX <- [BcX auxV(i)] */
325: VecCopy(d->auxV[i], d->BcX[d->size_cX+i]);
326: IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_cX+i, PETSC_NULL,
327: d->BcX, d->BcX[d->size_cX+i], PETSC_NULL,
328: &norm, &lindep);
329: if(lindep) SETERRQ(1, "Error during orth(BcX, B*cX(new))!");
330: VecScale(d->BcX[d->size_cX+i], 1.0/norm);
331: }
333: /* Harmonics restarts wiht right eigenvectors, and other with
334: the left ones */
335: pX = (d->W||!d->cY||d->BcX)?d->pX:d->pY;
336: ldpX = (d->W||!d->cY||d->BcX)?d->ldpX:d->ldpY;
338: /* If BcX, f.V <- orth(BcX, f.V) */
339: if (d->BcX) cMT = 0;
340: else cMT = d->size_H - npreconv;
342: /* f.MTX <- pY(npreconv:size_H-1), f.MTY <- f.pY(npreconv:size_H-1) */
343: d->ldMTX = d->ldMTY = d->size_H;
344: d->size_MT = d->size_H;
345: d->MT_type = DVD_MT_ORTHO;
346: SlepcDenseCopy(d->MTX, d->ldMTX, &pX[ldpX*npreconv], ldpX,
347: d->size_H, cMT);
348: if (d->W && d->pY) {
349: SlepcDenseCopy(d->MTY, d->ldMTY, &d->pY[d->ldpY*npreconv], d->ldpY,
350: d->size_H, cMT);
351: }
353: /* [f.cX(f.nconv) f.V] <- f.V*[f.pX(0:npreconv-1) f.pY(npreconv:f.size_V-1)] */
354: if (&d->cX[d->nconv] == d->V) { /* cX and V are contiguous */
355: SlepcDenseCopy(pX, ldpX, d->pX, d->ldpX, d->size_H, npreconv);
356:
357: SlepcUpdateVectorsZ(d->V, 0.0, 1.0, d->V, d->size_V, pX,
358: ldpX, d->size_H, d->size_H);
359: d->V+= npreconv;
360: inc_V = npreconv;
361: d->max_size_V-= npreconv;
362: } else {
363: SETERRQ(1, "Untested case!");
364: /*SlepcUpdateVectorsZ(&d->cX[d->nconv], 0.0, 1.0, d->V, d->size_V,
365: d->pX, d->ldpX, d->size_H, npreconv);
366:
367: SlepcUpdateVectorsZ(d->V, 0.0, 1.0, d->V, d->size_V,
368: &pX[ldpX*npreconv], ldpX,
369: d->size_H, cMT);
370: inc_V = 0;*/
371: }
372: d->size_cX+= npreconv;
373: d->size_V -= npreconv;
375: /* Udpate cS and cT, if needed */
376: if (d->cS) {
377: PetscInt size_cS = d->size_cX-npreconv;
378: cX = d->cY?d->cY:d->cX; size_cX = d->cY?d->size_cY:d->size_cX;
379: cy = d->cY?new_cY:0; size_cy = d->cY?npreconv:0;
380:
381: dvd_updateV_YtWx(&d->cS[d->ldcS*size_cS], d->ldcS, cX, size_cX, cy,
382: size_cy, d->AV, d->size_AV,
383: d->pX, d->ldpX,
384: d->size_H, npreconv, d->auxV);
386: if (DVD_ISNOT(d->sEP, DVD_EP_STD)) {if (d->BV) {
387:
388: dvd_updateV_YtWx(&d->cT[d->ldcT*size_cS], d->ldcT, cX, size_cX, cy,
389: size_cy, d->BV, d->size_AV,
390: d->pX, d->ldpX,
391: d->size_H, npreconv, d->auxV);
392: } else if (!d->B) {
393: VecsMultIa(&d->cT[d->ldcT*size_cS], 0, d->ldcT, cX, 0, size_cX,
394: &d->cX[size_cS], 0, npreconv);
395: VecsMultIa(&d->cT[d->ldcT*size_cS+size_cX], 0, d->ldcT, cy, 0,
396: size_cy, &d->cX[size_cS], 0, npreconv);
397: } else {
398: /* TODO: Only for nprecond==1 */
399: MatMult(d->B, d->cX[d->size_cX-1], d->auxV[0]);
400: VecsMultIa(&d->cT[d->ldcT*size_cS], 0, d->ldcT, cX, 0, size_cX,
401: &d->auxV[0], 0, npreconv);
402: VecsMultIa(&d->cT[d->ldcT*size_cS+size_cX], 0, d->ldcT, cy, 0,
403: size_cy, &d->auxV[0], 0, npreconv);
404: }}
405: }
407: /* f.W <- f.W * f.MTY */
408: if (d->W) {
409: SlepcUpdateVectorsZ(d->W, 0.0, 1.0, d->W, d->size_V+npreconv,
410: d->MTY, d->ldMTY, d->size_H, cMT);
411:
412: }
414: /* Lock the converged pairs */
415: d->eigr+= npreconv;
416: #ifndef PETSC_USE_COMPLEX
417: if (d->eigi) d->eigi+= npreconv;
418: #endif
419: d->nconv+= npreconv;
420: d->errest+= npreconv;
422: /* Notify the changes in V and update the other subspaces */
423: d->V_imm_s = inc_V; d->V_imm_e = inc_V;
424: d->V_tra_s = 0; d->V_tra_e = cMT;
425: d->V_new_s = d->V_tra_e; d->V_new_e = d->size_V;
427: /* Remove oldU */
428: data->size_oldU = 0;
430: /* f.pX <- I, f.pY <- I */
431: d->ldpX = d->size_H;
432: if (d->pY) d->ldpY = d->size_H;
433: for(i=0; i<d->size_H; i++) {
434: for(j=0; j<d->size_H; j++) {
435: d->pX[d->ldpX*i+j] = 0.0;
436: if (d->pY) d->pY[d->ldpY*i+j] = 0.0;
437: }
438: d->pX[d->ldpX*i+i] = 1.0;
439: if (d->pY) d->pY[d->ldpY*i+i] = 1.0;
440: }
442: d->npreconv-= npreconv;
444: return(0);
445: }
449: PetscErrorCode dvd_updateV_conv_finish(dvdDashboard *d)
450: {
451: PetscErrorCode ierr;
452: #if defined(PETSC_USE_COMPLEX)
453: PetscInt i, j;
454: PetscScalar s;
455: #endif
459: /* Finish cS and cT */
460: VecsMultIb(d->cS, 0, d->ldcS, d->nconv, d->nconv, d->auxS, d->V[0]);
461:
462: if (d->cT) {
463: VecsMultIb(d->cT, 0, d->ldcT, d->nconv, d->nconv, d->auxS, d->V[0]);
464:
465: }
467: /* Some functions need the diagonal elements in cT be real */
468: #if defined(PETSC_USE_COMPLEX)
469: if (d->cT) for(i=0; i<d->nconv; i++) {
470: s = PetscConj(d->cT[d->ldcT*i+i])/PetscAbsScalar(d->cT[d->ldcT*i+i]);
471: for(j=0; j<=i; j++)
472: d->cT[d->ldcT*i+j] = PetscRealPart(d->cT[d->ldcT*i+j]*s),
473: d->cS[d->ldcS*i+j]*= s;
474: VecScale(d->cX[i], s);
475: }
476: #endif
478: return(0);
479: }
480:
483: PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
484: {
485: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
486: PetscInt size_plusk, size_X, new_size_X, i, j;
487: PetscErrorCode ierr;
491: /* Select size_X desired pairs from V */
492: size_X = PetscMin(PetscMin(data->min_size_V,
493: d->size_V ),
494: d->max_size_V );
496: /* Add plusk eigenvectors from the previous iteration */
497: size_plusk = PetscMax(0, PetscMin(PetscMin(data->plusk,
498: data->size_oldU ),
499: d->max_size_V - size_X ));
501: /* Modify the subspaces */
502: d->ldMTX = d->size_MT = d->size_V;
503: dvd_updateV_restartV_aux(d->V, d->size_V, d->MTX, d->ldMTX, d->pX,
504: d->ldpX, size_X, data->oldU, data->ldoldU,
505: data->size_oldU, size_plusk, d->auxS,
506: d->size_auxS, &new_size_X);
507: if (d->W && d->pY) {
508: PetscInt new_size_Y;
509: d->ldMTY = d->size_V;
510: dvd_updateV_restartV_aux(d->W, d->size_V, d->MTY, d->ldMTY, d->pY,
511: d->ldpY, size_X, data->oldV, data->ldoldU,
512: data->size_oldU, new_size_X-size_X,
513: d->auxS, d->size_auxS, &new_size_Y);
514:
515: new_size_X = PetscMin(new_size_X, new_size_Y);
516: }
518: /* Notify the changes in V and update the other subspaces */
519: d->size_V = new_size_X;
520: d->MT_type = DVD_MT_ORTHO;
521: d->V_imm_s = 0; d->V_imm_e = 0;
522: d->V_tra_s = 0; d->V_tra_e = new_size_X;
523: d->V_new_s = d->V_tra_e; d->V_new_e = d->V_tra_e;
525: /* Remove oldU */
526: data->size_oldU = 0;
528: /* Remove npreconv */
529: d->npreconv = 0;
530:
531: /* f.pX <- I, f.pY <- I */
532: d->ldpX = d->size_H;
533: if (d->pY) d->ldpY = d->size_H;
534: for(i=0; i<d->size_H; i++) {
535: for(j=0; j<d->size_H; j++) {
536: d->pX[d->ldpX*i+j] = 0.0;
537: if (d->pY) d->pY[d->ldpY*i+j] = 0.0;
538: }
539: d->pX[d->ldpX*i+i] = 1.0;
540: if (d->pY) d->pY[d->ldpY*i+i] = 1.0;
541: }
543: return(0);
544: }
549: PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
550: {
551: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
552: PetscInt size_D;
553: PetscErrorCode ierr;
557: /* Select the desired pairs */
558: size_D = PetscMin(PetscMin(PetscMin(data->bs,
559: d->size_V ),
560: d->max_size_V-d->size_V ),
561: d->size_H );
562: if (size_D == 0) {
563: PetscPrintf(PETSC_COMM_WORLD, "MON: D:%d H:%d\n", size_D, d->size_H);
564: d->initV(d);
565: d->calcPairs(d);
566: //SETERRQ(1, "D == 0!\n");
567: //PetscFunctionReturn(1);
568: }
570: // PetscPrintf(PETSC_COMM_WORLD, "EIGS: ");
571: // for(i=0; i<d->size_H; i++) PetscPrintf(PETSC_COMM_WORLD, "%d:%g ", i, d->eigr[i]);
572: // PetscPrintf(PETSC_COMM_WORLD, "\n");
574: /* Fill V with D */
575: d->improveX(d, d->V+d->size_V, d->max_size_V-d->size_V, 0, size_D,
576: &size_D);
578: /* If D is empty, exit */
579: if (size_D == 0) { return(0); }
581: /* Get the converged pairs */
582: dvd_updateV_testConv(d, 0, size_D,
583: data->allResiduals?d->size_V:size_D, d->auxV, d->auxS,
584: &d->npreconv);
586: /* Notify the changes in V */
587: d->size_V+= size_D;
588: d->size_D = size_D;
589: d->V_imm_s = 0; d->V_imm_e = d->size_V-size_D;
590: d->V_tra_s = 0; d->V_tra_e = 0;
591: d->V_new_s = d->size_V-size_D; d->V_new_e = d->size_V;
593: /* Save the projected eigenvectors */
594: if (data->plusk > 0) {
595: data->ldoldU = data->size_oldU = d->size_H;
596: SlepcDenseCopy(data->oldU, data->ldoldU, d->pX, d->ldpX, d->size_H,
597: d->size_H);
598: if (d->pY && d->W) {
599: SlepcDenseCopy(data->oldV, data->ldoldU, d->pY, d->ldpY, d->size_H,
600: d->size_H);
601: }
602: }
604: return(0);
605: }
610: PetscErrorCode dvd_updateV_testConv(dvdDashboard *d, PetscInt s, PetscInt pre,
611: PetscInt e, Vec *auxV, PetscScalar *auxS,
612: PetscInt *nConv)
613: {
614: PetscInt i;
615: #ifndef PETSC_USE_COMPLEX
616: PetscInt j;
617: #endif
618: PetscReal norm;
619: PetscErrorCode ierr;
620: PetscTruth conv, c;
621: dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
624:
625: *nConv = s;
626: for(i=s, conv=PETSC_TRUE;
627: (conv || data->allResiduals) && (i < e);
628: i++) {
629: if (i >= pre) {
630: d->calcpairs_X(d, i, i+1, &auxV[1]);
631: if ((d->B) && DVD_IS(d->sEP, DVD_EP_STD)) {
632: MatMult(d->B, auxV[1], auxV[0]);
633: }
634: d->calcpairs_residual(d, i, i+1, auxV, auxS, auxV[1]);
635:
636: }
637: norm = d->nR[i]/d->nX[i];
638: c = d->testConv(d, d->eigr[i], d->eigi[i], norm, &d->errest[i]);
639: if (conv && c) *nConv = i+1;
640: else conv = PETSC_FALSE;
641: }
642: pre = PetscMax(pre, i);
644: #ifndef PETSC_USE_COMPLEX
645: /* Enforce converged conjugate conjugate complex eigenpairs */
646: for(j=0; j<*nConv; j++) if(d->eigi[j] != 0.0) j++;
647: if(j > *nConv) (*nConv)--;
648: #endif
649: for(i=pre; i<e; i++) d->errest[i] = d->nR[i] = PETSC_MAX;
650:
651: return(0);
652: }
655: /*
656: U <- [pX(0:size_X-1) gs(pX(0:size_X-1), oldpX(0:size_plusk-1))]
657: V <- V * U,
658: where
659: new_size_V, return the new size of V
660: auxS, auxiliar vector of size 2*ldpX, at least
661: size_auxS, the size of auxS
662: */
665: PetscErrorCode dvd_updateV_restartV_aux(Vec *V, PetscInt size_V,
666: PetscScalar *U, PetscInt ldU,
667: PetscScalar *pX, PetscInt ldpX,
668: PetscInt cpX, PetscScalar *oldpX,
669: PetscInt ldoldpX, PetscInt roldpX,
670: PetscInt coldpX, PetscScalar *auxS,
671: PetscInt size_auxS,
672: PetscInt *new_size_V)
673: {
674: PetscInt i, j;
675: PetscErrorCode ierr;
679: /* Add the size_X best eigenpairs */
680: SlepcDenseCopy(U, ldU, pX, ldpX, size_V, cpX);
682: /* Add plusk eigenvectors from the previous iteration */
683: SlepcDenseCopy(&U[cpX*ldU], ldU, oldpX, ldoldpX, roldpX, coldpX);
684:
685: for(i=cpX; i<cpX+coldpX; i++)
686: for(j=roldpX; j<size_V; j++)
687: U[i*ldU+j] = 0.0;
689: /* U <- orth(U) */
690: /// Temporal! Correct sentence: U <- orth(U(0:size_X-1), U(size_X:size_X+size_plusk))
691: if (coldpX > 0) {
692: SlepcDenseOrth(U, ldU, size_V, cpX+coldpX, auxS, size_auxS,
693: new_size_V);
694: } else
695: *new_size_V = cpX;
697: /* V <- V * U */
698: SlepcUpdateVectorsZ(V, 0.0, 1.0, V, size_V, U, ldU, size_V,
699: *new_size_V);
701: return(0);
702: }
704: /*
705: Compute S = [ Y' * W * x
706: y' * W * x ]
707: where
708: ldS, the leading dimension of S,
709: cY, number of vectors of Y,
710: ldy, the leading dimension of y,
711: ry,cy, rows and columns of y,
712: cW, number of vectors of W,
713: ldH, the leading dimension of H,
714: rH,cH, rows and columns of H,
715: ldx, the leading dimension of y,
716: rx,cx, rows and columns of x,
717: r, a reduction,
718: sr, a permanent space,
719: auxV, array of auxiliar vectors of size cx (at the end, auxV <- W*x),
720: auxS, auxiliar scalar vector of size rH*cx.
721: */
724: PetscErrorCode dvd_updateV_YtWx(PetscScalar *S, PetscInt ldS,
725: Vec *Y, PetscInt cY, Vec *y, PetscInt cy,
726: Vec *W, PetscInt cW, PetscScalar *x,
727: PetscInt ldx,
728: PetscInt rx, PetscInt cx, Vec *auxV)
729: {
730: PetscErrorCode ierr;
734: /* auxV <- W * x */
735: SlepcUpdateVectorsZ(auxV, 0.0, 1.0, W, cW, x, ldx, rx, cx);
736:
738: /* S(0:cY-1, 0:cx-1) <- Y' * auxV */
739: VecsMultIa(S, 0, ldS, Y, 0, cY, auxV, 0, cx);
741: /* S(cY:cY+cy-1, 0:cx-1) <- y' * auxV */
742: VecsMultIa(&S[cY], 0, ldS, y, 0, cy, auxV, 0, cx);
744: return(0);
745: }