Actual source code: davidson.c
1: /*
2: Method: General Davidson Method (includes GD and JD)
4: References:
5: - Ernest R. Davidson. Super-matrix methods. Computer Physics Communications,
6: 53:49–60, May 1989.
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: SLEPc - Scalable Library for Eigenvalue Problem Computations
10: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
12: This file is part of SLEPc.
13:
14: SLEPc is free software: you can redistribute it and/or modify it under the
15: terms of version 3 of the GNU Lesser General Public License as published by
16: the Free Software Foundation.
18: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
19: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
21: more details.
23: You should have received a copy of the GNU Lesser General Public License
24: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
25: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26: */
28: #include private/epsimpl.h
29: #include private/stimpl.h
30: #include davidson.h
31: #include slepcblaslapack.h
33: PetscErrorCode EPSView_DAVIDSON(EPS eps,PetscViewer viewer);
37: PetscErrorCode EPSCreate_DAVIDSON(EPS eps) {
38: PetscErrorCode ierr;
39: EPS_DAVIDSON *data;
43: STSetType(eps->OP, STPRECOND);
44: STPrecondSetKSPHasMat(eps->OP, PETSC_FALSE);
46: eps->OP->ops->getbilinearform = STGetBilinearForm_Default;
47: eps->ops->solve = EPSSolve_DAVIDSON;
48: eps->ops->setup = EPSSetUp_DAVIDSON;
49: eps->ops->destroy = EPSDestroy_DAVIDSON;
50: eps->ops->backtransform = EPSBackTransform_Default;
51: eps->ops->computevectors = EPSComputeVectors_QZ;
52: eps->ops->view = EPSView_DAVIDSON;
54: PetscMalloc(sizeof(EPS_DAVIDSON), &data);
55: eps->data = data;
56: data->pc = 0;
58: /* Set default values */
59: EPSDAVIDSONSetKrylovStart_DAVIDSON(eps, PETSC_FALSE);
60: EPSDAVIDSONSetBlockSize_DAVIDSON(eps, 1);
61: EPSDAVIDSONSetRestart_DAVIDSON(eps, 6, 0);
62: EPSDAVIDSONSetInitialSize_DAVIDSON(eps, 5);
63: EPSDAVIDSONSetFix_DAVIDSON(eps, 0.01);
65: dvd_prof_init();
67: return(0);
68: }
73: PetscErrorCode EPSSetUp_DAVIDSON(EPS eps) {
74: PetscErrorCode ierr;
75: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
76: dvdDashboard *dvd = &data->ddb;
77: dvdBlackboard b;
78: PetscInt i,nvecs,nscalars,min_size_V,plusk,bs,initv;
79: Mat A,B;
80: KSP ksp;
81: PC pc, pc2;
82: PetscTruth t,ipB,ispositive;
83: HarmType_t harm;
84: InitType_t init;
85: PetscReal fix;
89: /* Setup EPS options and get the problem specification */
90: EPSDAVIDSONGetBlockSize_DAVIDSON(eps, &bs);
91: if (bs <= 0) bs = 1;
92: if(eps->ncv) {
93: if (eps->ncv<eps->nev) SETERRQ(PETSC_ERR_SUP,"The value of ncv must be at least nev");
94: } else if (eps->mpd) eps->ncv = eps->mpd + eps->nev + bs;
95: else if (eps->nev<500)
96: eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15))+bs;
97: else
98: eps->ncv = PetscMin(eps->n,eps->nev+500)+bs;
99: if (!eps->mpd) eps->mpd = eps->ncv;
100: if (eps->mpd > eps->ncv)
101: SETERRQ(PETSC_ERR_SUP,"The mpd has to be less or equal than ncv");
102: if (eps->mpd < 2)
103: SETERRQ(PETSC_ERR_SUP,"The mpd has to be greater than 2");
104: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
105: if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
106: if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY))
107: SETERRQ(PETSC_ERR_SUP,"Wrong value of eps->which");
108: if (!(eps->nev + bs <= eps->ncv))
109: SETERRQ(PETSC_ERR_SUP, "The ncv has to be greater than nev plus blocksize!");
111: EPSDAVIDSONGetRestart_DAVIDSON(eps, &min_size_V, &plusk);
112:
113: if (!min_size_V) min_size_V = PetscMin(PetscMax(bs,5), eps->mpd/2);
114: if (!(min_size_V+bs <= eps->mpd))
115: SETERRQ(PETSC_ERR_SUP, "The value of minv must be less than mpd minus blocksize");
116: EPSDAVIDSONGetInitialSize_DAVIDSON(eps, &initv);
117: if (eps->mpd < initv)
118: SETERRQ(PETSC_ERR_SUP,"The initv has to be less or equal than mpd");
120: /* Davidson solvers do not support left eigenvectors */
121: if (eps->leftvecs) SETERRQ(PETSC_ERR_SUP,"Left vectors not supported in this solver");
123: /* Change the default sigma to inf if necessary */
124: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
125: eps->which == EPS_LARGEST_IMAGINARY) {
126: STSetDefaultShift(eps->OP, 3e300);
127: }
128:
129: /* Davidson solvers only support STPRECOND */
130: STSetUp(eps->OP);
131: PetscTypeCompare((PetscObject)eps->OP, STPRECOND, &t);
132: if (!t) SETERRQ1(PETSC_ERR_SUP, "%s only works with precond spectral transformation",
133: ((PetscObject)eps)->type_name);
135: /* Extract pc from st->ksp */
136: if (data->pc) { PCDestroy(data->pc); data->pc = 0; }
137: STGetKSP(eps->OP, &ksp);
138: KSPGetPC(ksp, &pc);
139: PetscTypeCompare((PetscObject)pc, PCNONE, &t);
140: if (t) {
141: pc = 0;
142: } else {
143: PetscObjectReference((PetscObject)pc);
144: data->pc = pc;
145: PCCreate(((PetscObject)eps)->comm, &pc2);
146: PCSetType(pc2, PCNONE);
147: KSPSetPC(ksp, pc2);
148: PCDestroy(pc2);
149: }
151: /* Setup problem specification in dvd */
152: STGetOperators(eps->OP, &A, &B);
153: PetscMemzero(dvd, sizeof(dvdDashboard));
154: dvd->A = A; dvd->B = eps->isgeneralized? B : PETSC_NULL;
155: ispositive = eps->ispositive;
156: dvd->sA = DVD_MAT_IMPLICIT |
157: (eps->ishermitian? DVD_MAT_HERMITIAN : 0) |
158: ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF : 0);
159: /* Asume -eps_hermitian means hermitian-definite in generalized problems */
160: if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
161: if (!eps->isgeneralized)
162: dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY |
163: DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
164: else
165: dvd->sB = DVD_MAT_IMPLICIT |
166: (eps->ishermitian? DVD_MAT_HERMITIAN : 0) |
167: (ispositive? DVD_MAT_POS_DEF : 0);
168: ipB = DVD_IS(dvd->sB, DVD_MAT_POS_DEF)?PETSC_TRUE:PETSC_FALSE;
169: dvd->sEP = ((!eps->isgeneralized || (eps->isgeneralized && ipB))? DVD_EP_STD : 0) |
170: (ispositive? DVD_EP_HERMITIAN : 0);
171: dvd->nev = eps->nev;
172: dvd->which = eps->which;
173: switch(eps->which) {
174: case EPS_TARGET_MAGNITUDE:
175: case EPS_TARGET_REAL:
176: case EPS_TARGET_IMAGINARY:
177: dvd->withTarget = PETSC_TRUE;
178: dvd->target[0] = eps->target; dvd->target[1] = 1.0;
179: break;
181: case EPS_LARGEST_MAGNITUDE:
182: case EPS_LARGEST_REAL:
183: case EPS_LARGEST_IMAGINARY: //TODO: think about this case
184: default:
185: dvd->withTarget = PETSC_TRUE;
186: dvd->target[0] = 1.0; dvd->target[1] = 0.0;
187: break;
188:
189: case EPS_SMALLEST_MAGNITUDE:
190: case EPS_SMALLEST_REAL:
191: case EPS_SMALLEST_IMAGINARY: //TODO: think about this case
192: dvd->withTarget = PETSC_TRUE;
193: dvd->target[0] = 0.0; dvd->target[1] = 1.0;
194: }
195: dvd->tol = eps->tol;
196: dvd->eps = eps;
198: /* Setup the extraction technique */
199: switch(eps->extraction) {
200: case 0:
201: case EPS_RITZ: harm = DVD_HARM_NONE; break;
202: case EPS_HARMONIC: harm = DVD_HARM_RR; break;
203: case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
204: case EPS_HARMONIC_RIGHT: harm = DVD_HARM_REIGS; break;
205: case EPS_HARMONIC_LARGEST: harm = DVD_HARM_LEIGS; break;
206: default: SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type");
207: }
209: /* Setup the type of starting subspace */
210: EPSDAVIDSONGetKrylovStart_DAVIDSON(eps, &t);
211: init = (!t)? DVD_INITV_CLASSIC : DVD_INITV_KRYLOV;
213: /* Setup IP */
214: if (ipB && dvd->B) {
215: IPSetBilinearForm(eps->ip, dvd->B, IP_INNER_HERMITIAN);
216: } else {
217: IPSetBilinearForm(eps->ip, 0, IP_INNER_HERMITIAN);
218: }
220: /* Get the fix parameter */
221: EPSDAVIDSONGetFix_DAVIDSON(eps, &fix);
223: /* Orthonormalize the DS */
224: dvd_orthV(eps->ip, PETSC_NULL, 0, PETSC_NULL, 0, eps->DS, 0,
225: PetscAbs(eps->nds), PETSC_NULL, 0, eps->rand);
227: /* Preconfigure dvd */
228: dvd_schm_basic_preconf(dvd, &b, eps->ncv, eps->mpd, min_size_V, bs,
229: initv,
230: PetscAbs(eps->nini),
231: plusk, pc, harm,
232: PETSC_NULL, init, eps->trackall);
233:
235: /* Reserve memory */
236: nvecs = b.max_size_auxV + b.own_vecs;
237: nscalars = b.own_scalars + b.max_size_auxS;
238: PetscMalloc((nvecs*eps->nloc+nscalars)*sizeof(PetscScalar), &data->wS);
239:
240: PetscMalloc(nvecs*sizeof(Vec), &data->wV);
241: data->size_wV = nvecs;
242: for (i=0; i<nvecs; i++) {
243: VecCreateMPIWithArray(((PetscObject)eps)->comm, eps->nloc, PETSC_DECIDE,
244: data->wS+i*eps->nloc, &data->wV[i]);
245:
246: }
247: b.free_vecs = data->wV;
248: b.free_scalars = data->wS + nvecs*eps->nloc;
249: dvd->auxV = data->wV + b.own_vecs;
250: dvd->auxS = b.free_scalars + b.own_scalars;
251: dvd->size_auxV = b.max_size_auxV;
252: dvd->size_auxS = b.max_size_auxS;
254: /* Configure dvd for a basic GD */
255: dvd_schm_basic_conf(dvd, &b, eps->ncv, eps->mpd, min_size_V, bs,
256: initv,
257: PetscAbs(eps->nini), plusk, pc,
258: eps->ip, harm, dvd->withTarget,
259: eps->target, ksp,
260: fix, init, eps->trackall);
261:
263: /* Associate the eigenvalues to the EPS */
264: eps->eigr = dvd->eigr;
265: eps->eigi = dvd->eigi;
266: eps->errest = dvd->errest;
267: eps->V = dvd->V;
269: return(0);
270: }
275: PetscErrorCode EPSSolve_DAVIDSON(EPS eps) {
276: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
277: dvdDashboard *d = &data->ddb;
278: KSP ksp;
279: PetscErrorCode ierr;
283: /* Call the starting routines */
284: DVD_FL_CALL(d->startList, d);
286: for(eps->its=0; eps->its < eps->max_it; eps->its++) {
287: /* Initialize V, if it is needed */
288: if (d->size_V == 0) { d->initV(d); }
290: /* Find the best approximated eigenpairs in V, X */
291: d->calcPairs(d);
293: /* Expand the subspace */
294: d->updateV(d);
296: /* Monitor progress */
297: eps->nconv = d->nconv;
298: EPSMonitor(eps, eps->its+1, eps->nconv, eps->eigr, eps->eigi, eps->errest, d->size_H+d->nconv);
300: /* Test for convergence */
301: if (eps->nconv >= eps->nev) break;
302: }
304: /* Call the ending routines */
305: DVD_FL_CALL(d->endList, d);
307: if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
308: else eps->reason = EPS_DIVERGED_ITS;
310: /* Merge the pc extracted from st->ksp in EPSSetUp_DAVIDSON */
311: if (data->pc) {
312: STGetKSP(eps->OP, &ksp);
313: KSPSetPC(ksp, data->pc);
314: PCDestroy(data->pc);
315: data->pc = 0;
316: }
318: return(0);
319: }
324: PetscErrorCode EPSDestroy_DAVIDSON(EPS eps) {
325: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
326: dvdDashboard *dvd = &data->ddb;
327: PetscErrorCode ierr;
328: PetscInt i;
332: /* Call step destructors and destroys the list */
333: DVD_FL_CALL(dvd->destroyList, dvd);
334: DVD_FL_DEL(dvd->destroyList);
335: DVD_FL_DEL(dvd->startList);
336: DVD_FL_DEL(dvd->endList);
338: for(i=0; i<data->size_wV; i++) {
339: VecDestroy(data->wV[i]);
340: }
341: PetscFree(data->wV);
342: PetscFree(data->wS);
343: PetscFree(data);
345: return(0);
346: }
351: PetscErrorCode EPSView_DAVIDSON(EPS eps,PetscViewer viewer)
352: {
353: PetscErrorCode ierr;
354: PetscTruth isascii;
355: PetscInt opi, opi0;
356: PetscTruth opb;
357: const char* name;
360:
361: name = ((PetscObject)eps)->type_name;
362: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
363: if (!isascii) {
364: SETERRQ2(1,"Viewer type %s not supported for %s",((PetscObject)viewer)->type_name,name);
365: }
366:
367: EPSDAVIDSONGetBlockSize_DAVIDSON(eps, &opi);
368: PetscViewerASCIIPrintf(viewer,"block size: %d\n", opi);
369: EPSDAVIDSONGetKrylovStart_DAVIDSON(eps, &opb);
370: if(!opb) {
371: PetscViewerASCIIPrintf(viewer,"type of the initial subspace: non-Krylov\n");
372: } else {
373: PetscViewerASCIIPrintf(viewer,"type of the initial subspace: Krylov\n");
374: }
375: EPSDAVIDSONGetRestart_DAVIDSON(eps, &opi, &opi0);
376: PetscViewerASCIIPrintf(viewer,"size of the subspace after restarting: %d\n", opi);
377: PetscViewerASCIIPrintf(viewer,"number of vectors after restarting from the previous iteration: %d\n", opi0);
379: return(0);
380: }
385: PetscErrorCode SLEPcNotImplemented() {
386: SETERRQ(1, "Not call this function!");
387: }
392: PetscErrorCode EPSDAVIDSONSetKrylovStart_DAVIDSON(EPS eps,PetscTruth krylovstart)
393: {
394: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
399: data->krylovstart = krylovstart;
401: return(0);
402: }
406: PetscErrorCode EPSDAVIDSONGetKrylovStart_DAVIDSON(EPS eps,PetscTruth *krylovstart)
407: {
408: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
413: *krylovstart = data->krylovstart;
415: return(0);
416: }
420: PetscErrorCode EPSDAVIDSONSetBlockSize_DAVIDSON(EPS eps,PetscInt blocksize)
421: {
422: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
427: if(blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
428: if(blocksize <= 0)
429: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
430: data->blocksize = blocksize;
432: return(0);
433: }
437: PetscErrorCode EPSDAVIDSONGetBlockSize_DAVIDSON(EPS eps,PetscInt *blocksize)
438: {
439: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
444: *blocksize = data->blocksize;
446: return(0);
447: }
451: PetscErrorCode EPSDAVIDSONSetRestart_DAVIDSON(EPS eps,PetscInt minv,PetscInt plusk)
452: {
453: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
458: if(minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
459: if(minv <= 0)
460: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
461: if(plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = 5;
462: if(plusk < 0)
463: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
464: data->minv = minv;
465: data->plusk = plusk;
467: return(0);
468: }
472: PetscErrorCode EPSDAVIDSONGetRestart_DAVIDSON(EPS eps,PetscInt *minv,PetscInt *plusk)
473: {
474: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
479: *minv = data->minv;
480: *plusk = data->plusk;
482: return(0);
483: }
487: PetscErrorCode EPSDAVIDSONGetInitialSize_DAVIDSON(EPS eps,PetscInt *initialsize)
488: {
489: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
494: *initialsize = data->initialsize;
496: return(0);
497: }
501: PetscErrorCode EPSDAVIDSONSetInitialSize_DAVIDSON(EPS eps,PetscInt initialsize)
502: {
503: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
508: if(initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
509: if(initialsize <= 0)
510: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
511: data->initialsize = initialsize;
513: return(0);
514: }
518: PetscErrorCode EPSDAVIDSONGetFix_DAVIDSON(EPS eps,PetscReal *fix)
519: {
520: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
525: *fix = data->fix;
527: return(0);
528: }
532: PetscErrorCode EPSDAVIDSONSetFix_DAVIDSON(EPS eps,PetscReal fix)
533: {
534: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
539: if(fix == PETSC_DEFAULT || fix == PETSC_DECIDE) fix = 0.01;
540: if(fix < 0.0)
541: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
542: data->fix = fix;
544: return(0);
545: }
550: /*
551: EPSComputeVectors_QZ - Compute eigenvectors from the vectors
552: provided by the eigensolver. This version is intended for solvers
553: that provide Schur vectors from the QZ decompositon. Given the partial
554: Schur decomposition OP*V=V*T, the following steps are performed:
555: 1) compute eigenvectors of (S,T): S*Z=T*Z*D
556: 2) compute eigenvectors of OP: X=V*Z
557: If left eigenvectors are required then also do Z'*Tl=D*Z', Y=W*Z
558: */
559: PetscErrorCode EPSComputeVectors_QZ(EPS eps)
560: {
561: PetscErrorCode ierr;
562: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
563: dvdDashboard *d = &data->ddb;
564: PetscScalar *pX, *auxS;
565: PetscInt size_auxS;
569: /* Compute the eigenvectors associated to (cS, cT) */
570: PetscMalloc(sizeof(PetscScalar)*d->nconv*d->nconv, &pX);
571: size_auxS = 11*d->nconv + 4*d->nconv*d->nconv;
572: PetscMalloc(sizeof(PetscScalar)*size_auxS, &auxS);
573: dvd_compute_eigenvectors(d->nconv, d->cS, d->ldcS, d->cT, d->ldcT,
574: pX, d->nconv, PETSC_NULL, 0, auxS,
575: size_auxS, PETSC_FALSE);
577: /* pX[i] <- pX[i] / ||pX[i]|| */
578: SlepcDenseNorm(pX, d->nconv, d->nconv, d->nconv, d->ceigi);
579:
581: /* V <- cX * pX */
582: SlepcUpdateVectorsZ(eps->V, 0.0, 1.0, d->cX, d->size_cX, pX,
583: d->nconv, d->nconv, d->nconv);
585: PetscFree(pX);
586: PetscFree(auxS);
588: eps->evecsavailable = PETSC_TRUE;
590: return(0);
591: }