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 == PETSC_FALSE)
133: SETERRQ1(PETSC_ERR_SUP, "%s only works with precond spectral transformation",
134: ((PetscObject)eps)->type_name);
136: /* Extract pc from st->ksp */
137: if (data->pc) { PCDestroy(data->pc); data->pc = 0; }
138: STGetKSP(eps->OP, &ksp);
139: KSPGetPC(ksp, &pc);
140: PetscTypeCompare((PetscObject)pc, PCNONE, &t);
141: if (t == PETSC_TRUE) {
142: pc = 0;
143: } else {
144: PetscObjectReference((PetscObject)pc);
145: data->pc = pc;
146: PCCreate(((PetscObject)eps)->comm, &pc2);
147: PCSetType(pc2, PCNONE);
148: KSPSetPC(ksp, pc2);
149: PCDestroy(pc2);
150: }
152: /* Setup problem specification in dvd */
153: STGetOperators(eps->OP, &A, &B);
154: PetscMemzero(dvd, sizeof(dvdDashboard));
155: dvd->A = A; dvd->B = (eps->isgeneralized==PETSC_TRUE) ? B : PETSC_NULL;
156: ispositive = eps->ispositive;
157: dvd->sA = DVD_MAT_IMPLICIT |
158: (eps->ishermitian == PETSC_TRUE ? DVD_MAT_HERMITIAN : 0) |
159: (((ispositive == PETSC_TRUE) &&
160: (eps->isgeneralized == PETSC_FALSE)) ? DVD_MAT_POS_DEF : 0);
161: /* Asume -eps_hermitian means hermitian-definite in generalized problems */
162: if ((ispositive == PETSC_FALSE) &&
163: (eps->isgeneralized == PETSC_FALSE) &&
164: (eps->ishermitian == PETSC_TRUE)) ispositive = PETSC_TRUE;
165: if (eps->isgeneralized == PETSC_FALSE)
166: dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY |
167: DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
168: else
169: dvd->sB = DVD_MAT_IMPLICIT |
170: (eps->ishermitian == PETSC_TRUE ? DVD_MAT_HERMITIAN : 0) |
171: (ispositive == PETSC_TRUE ? DVD_MAT_POS_DEF : 0);
172: ipB = DVD_IS(dvd->sB, DVD_MAT_POS_DEF)?PETSC_TRUE:PETSC_FALSE;
173: dvd->sEP = ((eps->isgeneralized == PETSC_FALSE) ||
174: ( (eps->isgeneralized == PETSC_TRUE) &&
175: (ipB == PETSC_TRUE) ) ? DVD_EP_STD : 0) |
176: (ispositive == PETSC_TRUE ? DVD_EP_HERMITIAN : 0);
177: dvd->nev = eps->nev;
178: dvd->which = eps->which;
179: switch(eps->which) {
180: case EPS_TARGET_MAGNITUDE:
181: case EPS_TARGET_REAL:
182: case EPS_TARGET_IMAGINARY:
183: dvd->withTarget = PETSC_TRUE;
184: dvd->target[0] = eps->target; dvd->target[1] = 1.0;
185: break;
187: case EPS_LARGEST_MAGNITUDE:
188: case EPS_LARGEST_REAL:
189: case EPS_LARGEST_IMAGINARY: //TODO: think about this case
190: default:
191: dvd->withTarget = PETSC_TRUE;
192: dvd->target[0] = 1.0; dvd->target[1] = 0.0;
193: break;
194:
195: case EPS_SMALLEST_MAGNITUDE:
196: case EPS_SMALLEST_REAL:
197: case EPS_SMALLEST_IMAGINARY: //TODO: think about this case
198: dvd->withTarget = PETSC_TRUE;
199: dvd->target[0] = 0.0; dvd->target[1] = 1.0;
200: }
201: dvd->tol = eps->tol;
202: dvd->eps = eps;
204: /* Setup the extraction technique */
205: switch(eps->extraction) {
206: case 0:
207: case EPS_RITZ: harm = DVD_HARM_NONE; break;
208: case EPS_HARMONIC: harm = DVD_HARM_RR; break;
209: case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
210: case EPS_HARMONIC_RIGHT: harm = DVD_HARM_REIGS; break;
211: case EPS_HARMONIC_LARGEST: harm = DVD_HARM_LEIGS; break;
212: default: SETERRQ(PETSC_ERR_SUP,"Unsupported extraction type");
213: }
215: /* Setup the type of starting subspace */
216: EPSDAVIDSONGetKrylovStart_DAVIDSON(eps, &t);
217: init = t==PETSC_FALSE ? DVD_INITV_CLASSIC : DVD_INITV_KRYLOV;
219: /* Setup IP */
220: if ((ipB == PETSC_TRUE) && (dvd->B)) {
221: IPSetBilinearForm(eps->ip, dvd->B, IP_INNER_HERMITIAN);
222: } else {
223: IPSetBilinearForm(eps->ip, 0, IP_INNER_HERMITIAN);
224: }
226: /* Get the fix parameter */
227: EPSDAVIDSONGetFix_DAVIDSON(eps, &fix);
229: /* Orthonormalize the DS */
230: dvd_orthV(eps->ip, PETSC_NULL, 0, PETSC_NULL, 0, eps->DS, 0, eps->nds,
231: PETSC_NULL, 0, eps->rand);
233: /* Preconfigure dvd */
234: dvd_schm_basic_preconf(dvd, &b, eps->ncv, eps->mpd, min_size_V, bs,
235: initv, eps->IS,
236: eps->nini,
237: plusk, pc, harm,
238: PETSC_NULL, init, eps->trackall);
239:
241: /* Reserve memory */
242: nvecs = b.max_size_auxV + b.own_vecs;
243: nscalars = b.own_scalars + b.max_size_auxS;
244: PetscMalloc((nvecs*eps->nloc+nscalars)*sizeof(PetscScalar), &data->wS);
245:
246: PetscMalloc(nvecs*sizeof(Vec), &data->wV);
247: data->size_wV = nvecs;
248: for (i=0; i<nvecs; i++) {
249: VecCreateMPIWithArray(((PetscObject)eps)->comm, eps->nloc, PETSC_DECIDE,
250: data->wS+i*eps->nloc, &data->wV[i]);
251:
252: }
253: b.free_vecs = data->wV;
254: b.free_scalars = data->wS + nvecs*eps->nloc;
255: dvd->auxV = data->wV + b.own_vecs;
256: dvd->auxS = b.free_scalars + b.own_scalars;
257: dvd->size_auxV = b.max_size_auxV;
258: dvd->size_auxS = b.max_size_auxS;
260: /* Configure dvd for a basic GD */
261: dvd_schm_basic_conf(dvd, &b, eps->ncv, eps->mpd, min_size_V, bs,
262: initv, eps->IS,
263: eps->nini, plusk, pc,
264: eps->ip, harm, dvd->withTarget,
265: eps->target, ksp,
266: fix, init, eps->trackall);
267:
269: /* Associate the eigenvalues to the EPS */
270: eps->eigr = dvd->eigr;
271: eps->eigi = dvd->eigi;
272: eps->errest = dvd->errest;
273: eps->V = dvd->V;
275: return(0);
276: }
281: PetscErrorCode EPSSolve_DAVIDSON(EPS eps) {
282: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
283: dvdDashboard *d = &data->ddb;
284: KSP ksp;
285: PetscErrorCode ierr;
289: /* Call the starting routines */
290: DVD_FL_CALL(d->startList, d);
292: for(eps->its=0; eps->its < eps->max_it; eps->its++) {
293: /* Initialize V, if it is needed */
294: if (d->size_V == 0) { d->initV(d); }
296: /* Find the best approximated eigenpairs in V, X */
297: d->calcPairs(d);
299: /* Expand the subspace */
300: d->updateV(d);
302: /* Monitor progress */
303: eps->nconv = d->nconv;
304: EPSMonitor(eps, eps->its+1, eps->nconv, eps->eigr, eps->eigi, eps->errest, d->size_H+d->nconv);
306: /* Test for convergence */
307: if (eps->nconv >= eps->nev) break;
308: }
310: /* Call the ending routines */
311: DVD_FL_CALL(d->endList, d);
313: if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
314: else eps->reason = EPS_DIVERGED_ITS;
316: /* Merge the pc extracted from st->ksp in EPSSetUp_DAVIDSON */
317: if (data->pc) {
318: STGetKSP(eps->OP, &ksp);
319: KSPSetPC(ksp, data->pc);
320: PCDestroy(data->pc);
321: data->pc = 0;
322: }
324: return(0);
325: }
330: PetscErrorCode EPSDestroy_DAVIDSON(EPS eps) {
331: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
332: dvdDashboard *dvd = &data->ddb;
333: PetscErrorCode ierr;
334: PetscInt i;
338: /* Call step destructors and destroys the list */
339: DVD_FL_CALL(dvd->destroyList, dvd);
340: DVD_FL_DEL(dvd->destroyList);
341: DVD_FL_DEL(dvd->startList);
342: DVD_FL_DEL(dvd->endList);
344: for(i=0; i<data->size_wV; i++) {
345: VecDestroy(data->wV[i]);
346: }
347: PetscFree(data->wV);
348: PetscFree(data->wS);
349: PetscFree(data);
351: return(0);
352: }
357: PetscErrorCode EPSView_DAVIDSON(EPS eps,PetscViewer viewer)
358: {
359: PetscErrorCode ierr;
360: PetscTruth isascii;
361: PetscInt opi, opi0;
362: PetscTruth opb;
363: const char* name;
366:
367: name = ((PetscObject)eps)->type_name;
368: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
369: if (!isascii) {
370: SETERRQ2(1,"Viewer type %s not supported for %s",((PetscObject)viewer)->type_name,name);
371: }
372:
373: EPSDAVIDSONGetBlockSize_DAVIDSON(eps, &opi);
374: PetscViewerASCIIPrintf(viewer,"block size: %d\n", opi);
375: EPSDAVIDSONGetKrylovStart_DAVIDSON(eps, &opb);
376: if(opb == PETSC_FALSE) {
377: PetscViewerASCIIPrintf(viewer,"type of the initial subspace: non-Krylov\n");
378: } else {
379: PetscViewerASCIIPrintf(viewer,"type of the initial subspace: Krylov\n");
380: }
381: EPSDAVIDSONGetRestart_DAVIDSON(eps, &opi, &opi0);
382: PetscViewerASCIIPrintf(viewer,"size of the subspace after restarting: %d\n", opi);
383: PetscViewerASCIIPrintf(viewer,"number of vectors after restarting from the previous iteration: %d\n", opi0);
385: return(0);
386: }
391: PetscErrorCode SLEPcNotImplemented() {
392: SETERRQ(1, "Not call this function!");
393: }
398: PetscErrorCode EPSDAVIDSONSetKrylovStart_DAVIDSON(EPS eps,PetscTruth krylovstart)
399: {
400: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
405: data->krylovstart = krylovstart;
407: return(0);
408: }
412: PetscErrorCode EPSDAVIDSONGetKrylovStart_DAVIDSON(EPS eps,PetscTruth *krylovstart)
413: {
414: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
419: *krylovstart = data->krylovstart;
421: return(0);
422: }
426: PetscErrorCode EPSDAVIDSONSetBlockSize_DAVIDSON(EPS eps,PetscInt blocksize)
427: {
428: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
433: if(blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
434: if(blocksize <= 0)
435: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
436: data->blocksize = blocksize;
438: return(0);
439: }
443: PetscErrorCode EPSDAVIDSONGetBlockSize_DAVIDSON(EPS eps,PetscInt *blocksize)
444: {
445: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
450: *blocksize = data->blocksize;
452: return(0);
453: }
457: PetscErrorCode EPSDAVIDSONSetRestart_DAVIDSON(EPS eps,PetscInt minv,PetscInt plusk)
458: {
459: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
464: if(minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
465: if(minv <= 0)
466: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
467: if(plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = 5;
468: if(plusk < 0)
469: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
470: data->minv = minv;
471: data->plusk = plusk;
473: return(0);
474: }
478: PetscErrorCode EPSDAVIDSONGetRestart_DAVIDSON(EPS eps,PetscInt *minv,PetscInt *plusk)
479: {
480: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
485: *minv = data->minv;
486: *plusk = data->plusk;
488: return(0);
489: }
493: PetscErrorCode EPSDAVIDSONGetInitialSize_DAVIDSON(EPS eps,PetscInt *initialsize)
494: {
495: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
500: *initialsize = data->initialsize;
502: return(0);
503: }
507: PetscErrorCode EPSDAVIDSONSetInitialSize_DAVIDSON(EPS eps,PetscInt initialsize)
508: {
509: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
514: if(initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
515: if(initialsize <= 0)
516: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
517: data->initialsize = initialsize;
519: return(0);
520: }
524: PetscErrorCode EPSDAVIDSONGetFix_DAVIDSON(EPS eps,PetscReal *fix)
525: {
526: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
531: *fix = data->fix;
533: return(0);
534: }
538: PetscErrorCode EPSDAVIDSONSetFix_DAVIDSON(EPS eps,PetscReal fix)
539: {
540: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
545: if(fix == PETSC_DEFAULT || fix == PETSC_DECIDE) fix = 0.01;
546: if(fix < 0.0)
547: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
548: data->fix = fix;
550: return(0);
551: }
556: /*
557: EPSComputeVectors_QZ - Compute eigenvectors from the vectors
558: provided by the eigensolver. This version is intended for solvers
559: that provide Schur vectors from the QZ decompositon. Given the partial
560: Schur decomposition OP*V=V*T, the following steps are performed:
561: 1) compute eigenvectors of (S,T): S*Z=T*Z*D
562: 2) compute eigenvectors of OP: X=V*Z
563: If left eigenvectors are required then also do Z'*Tl=D*Z', Y=W*Z
564: */
565: PetscErrorCode EPSComputeVectors_QZ(EPS eps)
566: {
567: PetscErrorCode ierr;
568: EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;
569: dvdDashboard *d = &data->ddb;
570: PetscScalar *pX, *auxS;
571: PetscInt size_auxS;
575: /* Compute the eigenvectors associated to (cS, cT) */
576: PetscMalloc(sizeof(PetscScalar)*d->nconv*d->nconv, &pX);
577: size_auxS = 11*d->nconv + 4*d->nconv*d->nconv;
578: PetscMalloc(sizeof(PetscScalar)*size_auxS, &auxS);
579: dvd_compute_eigenvectors(d->nconv, d->cS, d->ldcS, d->cT, d->ldcT,
580: pX, d->nconv, PETSC_NULL, 0, auxS,
581: size_auxS, PETSC_FALSE);
583: /* pX[i] <- pX[i] / ||pX[i]|| */
584: SlepcDenseNorm(pX, d->nconv, d->nconv, d->nconv, d->ceigi);
585:
587: /* V <- cX * pX */
588: SlepcUpdateVectorsZ(eps->V, 0.0, 1.0, d->cX, d->size_cX, pX,
589: d->nconv, d->nconv, d->nconv);
591: PetscFree(pX);
592: PetscFree(auxS);
594: eps->evecsavailable = PETSC_TRUE;
596: return(0);
597: }