Actual source code: arnoldi.c
1: /*
3: SLEPc eigensolver: "arnoldi"
5: Method: Explicitly Restarted Arnoldi
7: Algorithm:
9: Arnoldi method with explicit restart and deflation.
11: References:
13: [1] "Arnoldi Methods in SLEPc", SLEPc Technical Report STR-4,
14: available at http://www.grycap.upv.es/slepc.
16: Last update: Feb 2009
18: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
19: SLEPc - Scalable Library for Eigenvalue Problem Computations
20: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
22: This file is part of SLEPc.
24: SLEPc is free software: you can redistribute it and/or modify it under the
25: terms of version 3 of the GNU Lesser General Public License as published by
26: the Free Software Foundation.
28: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
29: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
30: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
31: more details.
33: You should have received a copy of the GNU Lesser General Public License
34: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
35: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: */
38: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
39: #include <slepcblaslapack.h>
41: PetscErrorCode EPSSolve_Arnoldi(EPS);
43: typedef struct {
44: PetscBool delayed;
45: } EPS_ARNOLDI;
49: PetscErrorCode EPSSetUp_Arnoldi(EPS eps)
50: {
54: if (eps->ncv) { /* ncv set */
55: if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
56: } else if (eps->mpd) { /* mpd set */
57: eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
58: } else { /* neither set: defaults depend on nev being small or large */
59: if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
60: else {
61: eps->mpd = 500;
62: eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
63: }
64: }
65: if (!eps->mpd) eps->mpd = eps->ncv;
66: if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
67: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
68: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
69: if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
71: if (!eps->extraction) {
72: EPSSetExtraction(eps,EPS_RITZ);
73: }
74: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
76: EPSAllocateSolution(eps);
77: DSSetType(eps->ds,DSNHEP);
78: if (eps->extraction==EPS_REFINED || eps->extraction==EPS_REFINED_HARMONIC) {
79: DSSetRefined(eps->ds,PETSC_TRUE);
80: }
81: DSSetExtraRow(eps->ds,PETSC_TRUE);
82: DSAllocate(eps->ds,eps->ncv+1);
83: EPSSetWorkVecs(eps,1);
85: /* dispatch solve method */
86: if (eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Left vectors not supported in this solver");
87: if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems");
88: eps->ops->solve = EPSSolve_Arnoldi;
89: return(0);
90: }
94: /*
95: EPSDelayedArnoldi - This function is equivalent to EPSBasicArnoldi but
96: performs the computation in a different way. The main idea is that
97: reorthogonalization is delayed to the next Arnoldi step. This version is
98: more scalable but in some cases convergence may stagnate.
99: */
100: PetscErrorCode EPSDelayedArnoldi(EPS eps,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscBool *breakdown)
101: {
103: PetscInt i,j,m=*M;
104: Vec u,t;
105: PetscScalar shh[100],*lhh,dot,dot2;
106: PetscReal norm1=0.0,norm2;
109: if (m<=100) lhh = shh;
110: else {
111: PetscMalloc(m*sizeof(PetscScalar),&lhh);
112: }
113: VecDuplicate(f,&u);
114: VecDuplicate(f,&t);
116: for (j=k;j<m;j++) {
117: STApply(eps->st,V[j],f);
118: IPOrthogonalize(eps->ip,0,NULL,eps->nds,NULL,eps->defl,f,NULL,NULL,NULL);
120: IPMInnerProductBegin(eps->ip,f,j+1,V,H+ldh*j);
121: if (j>k) {
122: IPMInnerProductBegin(eps->ip,V[j],j,V,lhh);
123: IPInnerProductBegin(eps->ip,V[j],V[j],&dot);
124: }
125: if (j>k+1) {
126: IPNormBegin(eps->ip,u,&norm2);
127: VecDotBegin(u,V[j-2],&dot2);
128: }
130: IPMInnerProductEnd(eps->ip,f,j+1,V,H+ldh*j);
131: if (j>k) {
132: IPMInnerProductEnd(eps->ip,V[j],j,V,lhh);
133: IPInnerProductEnd(eps->ip,V[j],V[j],&dot);
134: }
135: if (j>k+1) {
136: IPNormEnd(eps->ip,u,&norm2);
137: VecDotEnd(u,V[j-2],&dot2);
138: if (PetscAbsScalar(dot2/norm2) > PETSC_MACHINE_EPSILON) {
139: *breakdown = PETSC_TRUE;
140: *M = j-1;
141: *beta = norm2;
143: if (m>100) { PetscFree(lhh); }
144: VecDestroy(&u);
145: VecDestroy(&t);
146: return(0);
147: }
148: }
150: if (j>k) {
151: norm1 = PetscSqrtReal(PetscRealPart(dot));
152: for (i=0;i<j;i++)
153: H[ldh*j+i] = H[ldh*j+i]/norm1;
154: H[ldh*j+j] = H[ldh*j+j]/dot;
156: VecCopy(V[j],t);
157: VecScale(V[j],1.0/norm1);
158: VecScale(f,1.0/norm1);
159: }
161: SlepcVecMAXPBY(f,1.0,-1.0,j+1,H+ldh*j,V);
163: if (j>k) {
164: SlepcVecMAXPBY(t,1.0,-1.0,j,lhh,V);
165: for (i=0;i<j;i++)
166: H[ldh*(j-1)+i] += lhh[i];
167: }
169: if (j>k+1) {
170: VecCopy(u,V[j-1]);
171: VecScale(V[j-1],1.0/norm2);
172: H[ldh*(j-2)+j-1] = norm2;
173: }
175: if (j<m-1) {
176: VecCopy(f,V[j+1]);
177: VecCopy(t,u);
178: }
179: }
181: IPNorm(eps->ip,t,&norm2);
182: VecScale(t,1.0/norm2);
183: VecCopy(t,V[m-1]);
184: H[ldh*(m-2)+m-1] = norm2;
186: IPMInnerProduct(eps->ip,f,m,V,lhh);
188: SlepcVecMAXPBY(f,1.0,-1.0,m,lhh,V);
189: for (i=0;i<m;i++)
190: H[ldh*(m-1)+i] += lhh[i];
192: IPNorm(eps->ip,f,beta);
193: VecScale(f,1.0 / *beta);
194: *breakdown = PETSC_FALSE;
196: if (m>100) { PetscFree(lhh); }
197: VecDestroy(&u);
198: VecDestroy(&t);
199: return(0);
200: }
204: /*
205: EPSDelayedArnoldi1 - This function is similar to EPSDelayedArnoldi1,
206: but without reorthogonalization (only delayed normalization).
207: */
208: PetscErrorCode EPSDelayedArnoldi1(EPS eps,PetscScalar *H,PetscInt ldh,Vec *V,PetscInt k,PetscInt *M,Vec f,PetscReal *beta,PetscBool *breakdown)
209: {
211: PetscInt i,j,m=*M;
212: PetscScalar dot;
213: PetscReal norm=0.0;
216: for (j=k;j<m;j++) {
217: STApply(eps->st,V[j],f);
218: IPOrthogonalize(eps->ip,0,NULL,eps->nds,NULL,eps->defl,f,NULL,NULL,NULL);
220: IPMInnerProductBegin(eps->ip,f,j+1,V,H+ldh*j);
221: if (j>k) {
222: IPInnerProductBegin(eps->ip,V[j],V[j],&dot);
223: }
225: IPMInnerProductEnd(eps->ip,f,j+1,V,H+ldh*j);
226: if (j>k) {
227: IPInnerProductEnd(eps->ip,V[j],V[j],&dot);
228: }
230: if (j>k) {
231: norm = PetscSqrtReal(PetscRealPart(dot));
232: VecScale(V[j],1.0/norm);
233: H[ldh*(j-1)+j] = norm;
235: for (i=0;i<j;i++)
236: H[ldh*j+i] = H[ldh*j+i]/norm;
237: H[ldh*j+j] = H[ldh*j+j]/dot;
238: VecScale(f,1.0/norm);
239: }
241: SlepcVecMAXPBY(f,1.0,-1.0,j+1,H+ldh*j,V);
243: if (j<m-1) {
244: VecCopy(f,V[j+1]);
245: }
246: }
248: IPNorm(eps->ip,f,beta);
249: VecScale(f,1.0 / *beta);
250: *breakdown = PETSC_FALSE;
251: return(0);
252: }
256: PetscErrorCode EPSSolve_Arnoldi(EPS eps)
257: {
258: PetscErrorCode ierr;
259: PetscInt k,nv,ld;
260: Vec f=eps->work[0];
261: PetscScalar *H,*U,*X;
262: PetscReal beta,gamma=1.0;
263: PetscBool breakdown,harmonic,refined;
264: IPOrthogRefineType orthog_ref;
265: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
268: DSGetLeadingDimension(eps->ds,&ld);
269: DSGetRefined(eps->ds,&refined);
270: harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
271: IPGetOrthogonalization(eps->ip,NULL,&orthog_ref,NULL);
273: /* Get the starting Arnoldi vector */
274: EPSGetStartVector(eps,0,eps->V[0],NULL);
276: /* Restart loop */
277: while (eps->reason == EPS_CONVERGED_ITERATING) {
278: eps->its++;
280: /* Compute an nv-step Arnoldi factorization */
281: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
282: DSSetDimensions(eps->ds,nv,0,eps->nconv,0);
283: DSGetArray(eps->ds,DS_MAT_A,&H);
284: if (!arnoldi->delayed) {
285: EPSBasicArnoldi(eps,PETSC_FALSE,H,ld,eps->V,eps->nconv,&nv,f,&beta,&breakdown);
286: H[nv+(nv-1)*ld] = beta;
287: } else if (orthog_ref == IP_ORTHOG_REFINE_NEVER) {
288: EPSDelayedArnoldi1(eps,H,ld,eps->V,eps->nconv,&nv,f,&beta,&breakdown);
289: } else {
290: EPSDelayedArnoldi(eps,H,ld,eps->V,eps->nconv,&nv,f,&beta,&breakdown);
291: }
292: DSRestoreArray(eps->ds,DS_MAT_A,&H);
293: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
295: /* Compute translation of Krylov decomposition if harmonic extraction used */
296: if (harmonic) {
297: DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,NULL,&gamma);
298: }
300: /* Solve projected problem */
301: DSSolve(eps->ds,eps->eigr,eps->eigi);
302: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
303: DSUpdateExtraRow(eps->ds);
305: /* Check convergence */
306: EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,eps->V,nv,beta,gamma,&k);
307: if (refined) {
308: DSGetArray(eps->ds,DS_MAT_X,&X);
309: SlepcVecMAXPBY(eps->V[k],0.0,1.0,nv,X+k*ld,eps->V);
310: DSRestoreArray(eps->ds,DS_MAT_X,&X);
311: DSGetArray(eps->ds,DS_MAT_Q,&U);
312: SlepcUpdateVectors(nv,eps->V,eps->nconv,PetscMin(k,nv),U,ld,PETSC_FALSE);
313: DSRestoreArray(eps->ds,DS_MAT_Q,&U);
314: IPOrthogonalize(eps->ip,0,NULL,k,NULL,eps->V,eps->V[k],NULL,NULL,NULL);
315: } else {
316: DSGetArray(eps->ds,DS_MAT_Q,&U);
317: SlepcUpdateVectors(nv,eps->V,eps->nconv,PetscMin(k+1,nv),U,ld,PETSC_FALSE);
318: DSRestoreArray(eps->ds,DS_MAT_Q,&U);
319: }
320: eps->nconv = k;
322: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
323: if (breakdown) {
324: PetscInfo2(eps,"Breakdown in Arnoldi method (it=%D norm=%G)\n",eps->its,beta);
325: EPSGetStartVector(eps,k,eps->V[k],&breakdown);
326: if (breakdown) {
327: eps->reason = EPS_DIVERGED_BREAKDOWN;
328: PetscInfo(eps,"Unable to generate more start vectors\n");
329: }
330: }
331: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
332: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
333: }
335: /* truncate Schur decomposition and change the state to raw so that
336: PSVectors() computes eigenvectors from scratch */
337: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
338: DSSetState(eps->ds,DS_STATE_RAW);
339: return(0);
340: }
344: PetscErrorCode EPSSetFromOptions_Arnoldi(EPS eps)
345: {
347: PetscBool set,val;
348: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
351: PetscOptionsHead("EPS Arnoldi Options");
352: PetscOptionsBool("-eps_arnoldi_delayed","Arnoldi with delayed reorthogonalization","EPSArnoldiSetDelayed",arnoldi->delayed,&val,&set);
353: if (set) {
354: EPSArnoldiSetDelayed(eps,val);
355: }
356: PetscOptionsTail();
357: return(0);
358: }
362: static PetscErrorCode EPSArnoldiSetDelayed_Arnoldi(EPS eps,PetscBool delayed)
363: {
364: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
367: arnoldi->delayed = delayed;
368: return(0);
369: }
373: /*@
374: EPSArnoldiSetDelayed - Activates or deactivates delayed reorthogonalization
375: in the Arnoldi iteration.
377: Logically Collective on EPS
379: Input Parameters:
380: + eps - the eigenproblem solver context
381: - delayed - boolean flag
383: Options Database Key:
384: . -eps_arnoldi_delayed - Activates delayed reorthogonalization in Arnoldi
386: Note:
387: Delayed reorthogonalization is an aggressive optimization for the Arnoldi
388: eigensolver than may provide better scalability, but sometimes makes the
389: solver converge less than the default algorithm.
391: Level: advanced
393: .seealso: EPSArnoldiGetDelayed()
394: @*/
395: PetscErrorCode EPSArnoldiSetDelayed(EPS eps,PetscBool delayed)
396: {
402: PetscTryMethod(eps,"EPSArnoldiSetDelayed_C",(EPS,PetscBool),(eps,delayed));
403: return(0);
404: }
408: static PetscErrorCode EPSArnoldiGetDelayed_Arnoldi(EPS eps,PetscBool *delayed)
409: {
410: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
413: *delayed = arnoldi->delayed;
414: return(0);
415: }
419: /*@C
420: EPSArnoldiGetDelayed - Gets the type of reorthogonalization used during the Arnoldi
421: iteration.
423: Not Collective
425: Input Parameter:
426: . eps - the eigenproblem solver context
428: Input Parameter:
429: . delayed - boolean flag indicating if delayed reorthogonalization has been enabled
431: Level: advanced
433: .seealso: EPSArnoldiSetDelayed()
434: @*/
435: PetscErrorCode EPSArnoldiGetDelayed(EPS eps,PetscBool *delayed)
436: {
442: PetscTryMethod(eps,"EPSArnoldiGetDelayed_C",(EPS,PetscBool*),(eps,delayed));
443: return(0);
444: }
448: PetscErrorCode EPSDestroy_Arnoldi(EPS eps)
449: {
453: PetscFree(eps->data);
454: PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiSetDelayed_C",NULL);
455: PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiGetDelayed_C",NULL);
456: return(0);
457: }
461: PetscErrorCode EPSView_Arnoldi(EPS eps,PetscViewer viewer)
462: {
464: PetscBool isascii;
465: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
468: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
469: if (isascii) {
470: if (arnoldi->delayed) {
471: PetscViewerASCIIPrintf(viewer," Arnoldi: using delayed reorthogonalization\n");
472: }
473: }
474: return(0);
475: }
479: PETSC_EXTERN PetscErrorCode EPSCreate_Arnoldi(EPS eps)
480: {
484: PetscNewLog(eps,EPS_ARNOLDI,&eps->data);
485: eps->ops->setup = EPSSetUp_Arnoldi;
486: eps->ops->setfromoptions = EPSSetFromOptions_Arnoldi;
487: eps->ops->destroy = EPSDestroy_Arnoldi;
488: eps->ops->reset = EPSReset_Default;
489: eps->ops->view = EPSView_Arnoldi;
490: eps->ops->backtransform = EPSBackTransform_Default;
491: eps->ops->computevectors = EPSComputeVectors_Schur;
492: PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiSetDelayed_C",EPSArnoldiSetDelayed_Arnoldi);
493: PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiGetDelayed_C",EPSArnoldiGetDelayed_Arnoldi);
494: return(0);
495: }