Actual source code: arnoldi.c
slepc-3.6.1 2015-09-03
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://slepc.upv.es.
16: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
17: SLEPc - Scalable Library for Eigenvalue Problem Computations
18: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
20: This file is part of SLEPc.
22: SLEPc is free software: you can redistribute it and/or modify it under the
23: terms of version 3 of the GNU Lesser General Public License as published by
24: the Free Software Foundation.
26: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
27: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
28: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
29: more details.
31: You should have received a copy of the GNU Lesser General Public License
32: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
33: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34: */
36: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
38: PetscErrorCode EPSSolve_Arnoldi(EPS);
40: typedef struct {
41: PetscBool delayed;
42: } EPS_ARNOLDI;
46: PetscErrorCode EPSSetUp_Arnoldi(EPS eps)
47: {
51: EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
52: if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
53: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
54: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
55: if (eps->ishermitian && eps->ispositive && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
57: if (!eps->extraction) {
58: EPSSetExtraction(eps,EPS_RITZ);
59: }
60: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
62: EPSAllocateSolution(eps,1);
63: EPS_SetInnerProduct(eps);
64: DSSetType(eps->ds,DSNHEP);
65: if (eps->extraction==EPS_REFINED || eps->extraction==EPS_REFINED_HARMONIC) {
66: DSSetRefined(eps->ds,PETSC_TRUE);
67: }
68: DSSetExtraRow(eps->ds,PETSC_TRUE);
69: DSAllocate(eps->ds,eps->ncv+1);
71: /* dispatch solve method */
72: if (eps->isgeneralized && eps->ishermitian && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method does not work for indefinite problems");
73: eps->ops->solve = EPSSolve_Arnoldi;
74: return(0);
75: }
79: PetscErrorCode EPSSolve_Arnoldi(EPS eps)
80: {
81: PetscErrorCode ierr;
82: PetscInt k,nv,ld;
83: Mat U;
84: PetscScalar *H,*X;
85: PetscReal beta,gamma=1.0;
86: PetscBool breakdown,harmonic,refined;
87: BVOrthogRefineType orthog_ref;
88: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
91: DSGetLeadingDimension(eps->ds,&ld);
92: DSGetRefined(eps->ds,&refined);
93: harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
94: BVGetOrthogonalization(eps->V,NULL,&orthog_ref,NULL,NULL);
96: /* Get the starting Arnoldi vector */
97: EPSGetStartVector(eps,0,NULL);
99: /* Restart loop */
100: while (eps->reason == EPS_CONVERGED_ITERATING) {
101: eps->its++;
103: /* Compute an nv-step Arnoldi factorization */
104: nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
105: DSSetDimensions(eps->ds,nv,0,eps->nconv,0);
106: DSGetArray(eps->ds,DS_MAT_A,&H);
107: if (!arnoldi->delayed) {
108: EPSBasicArnoldi(eps,PETSC_FALSE,H,ld,eps->nconv,&nv,&beta,&breakdown);
109: } else if (orthog_ref == BV_ORTHOG_REFINE_NEVER) {
110: EPSDelayedArnoldi1(eps,H,ld,eps->nconv,&nv,&beta,&breakdown);
111: } else {
112: EPSDelayedArnoldi(eps,H,ld,eps->nconv,&nv,&beta,&breakdown);
113: }
114: DSRestoreArray(eps->ds,DS_MAT_A,&H);
115: DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
116: BVSetActiveColumns(eps->V,eps->nconv,nv);
118: /* Compute translation of Krylov decomposition if harmonic extraction used */
119: if (harmonic) {
120: DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,NULL,&gamma);
121: }
123: /* Solve projected problem */
124: DSSolve(eps->ds,eps->eigr,eps->eigi);
125: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
126: DSUpdateExtraRow(eps->ds);
128: /* Check convergence */
129: EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,gamma,&k);
130: if (refined) {
131: DSGetArray(eps->ds,DS_MAT_X,&X);
132: BVMultColumn(eps->V,1.0,0.0,k,X+k*ld);
133: DSRestoreArray(eps->ds,DS_MAT_X,&X);
134: DSGetMat(eps->ds,DS_MAT_Q,&U);
135: BVMultInPlace(eps->V,U,eps->nconv,nv);
136: MatDestroy(&U);
137: BVOrthogonalizeColumn(eps->V,k,NULL,NULL,NULL);
138: } else {
139: DSGetMat(eps->ds,DS_MAT_Q,&U);
140: BVMultInPlace(eps->V,U,eps->nconv,nv);
141: MatDestroy(&U);
142: }
143: eps->nconv = k;
145: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
146: if (breakdown && k<eps->nev) {
147: PetscInfo2(eps,"Breakdown in Arnoldi method (it=%D norm=%g)\n",eps->its,(double)beta);
148: EPSGetStartVector(eps,k,&breakdown);
149: if (breakdown) {
150: eps->reason = EPS_DIVERGED_BREAKDOWN;
151: PetscInfo(eps,"Unable to generate more start vectors\n");
152: }
153: }
154: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
155: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
156: }
158: /* truncate Schur decomposition and change the state to raw so that
159: DSVectors() computes eigenvectors from scratch */
160: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
161: DSSetState(eps->ds,DS_STATE_RAW);
162: return(0);
163: }
167: PetscErrorCode EPSSetFromOptions_Arnoldi(PetscOptions *PetscOptionsObject,EPS eps)
168: {
170: PetscBool set,val;
171: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
174: PetscOptionsHead(PetscOptionsObject,"EPS Arnoldi Options");
175: PetscOptionsBool("-eps_arnoldi_delayed","Arnoldi with delayed reorthogonalization","EPSArnoldiSetDelayed",arnoldi->delayed,&val,&set);
176: if (set) {
177: EPSArnoldiSetDelayed(eps,val);
178: }
179: PetscOptionsTail();
180: return(0);
181: }
185: static PetscErrorCode EPSArnoldiSetDelayed_Arnoldi(EPS eps,PetscBool delayed)
186: {
187: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
190: arnoldi->delayed = delayed;
191: return(0);
192: }
196: /*@
197: EPSArnoldiSetDelayed - Activates or deactivates delayed reorthogonalization
198: in the Arnoldi iteration.
200: Logically Collective on EPS
202: Input Parameters:
203: + eps - the eigenproblem solver context
204: - delayed - boolean flag
206: Options Database Key:
207: . -eps_arnoldi_delayed - Activates delayed reorthogonalization in Arnoldi
209: Note:
210: Delayed reorthogonalization is an aggressive optimization for the Arnoldi
211: eigensolver than may provide better scalability, but sometimes makes the
212: solver converge less than the default algorithm.
214: Level: advanced
216: .seealso: EPSArnoldiGetDelayed()
217: @*/
218: PetscErrorCode EPSArnoldiSetDelayed(EPS eps,PetscBool delayed)
219: {
225: PetscTryMethod(eps,"EPSArnoldiSetDelayed_C",(EPS,PetscBool),(eps,delayed));
226: return(0);
227: }
231: static PetscErrorCode EPSArnoldiGetDelayed_Arnoldi(EPS eps,PetscBool *delayed)
232: {
233: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
236: *delayed = arnoldi->delayed;
237: return(0);
238: }
242: /*@
243: EPSArnoldiGetDelayed - Gets the type of reorthogonalization used during the Arnoldi
244: iteration.
246: Not Collective
248: Input Parameter:
249: . eps - the eigenproblem solver context
251: Input Parameter:
252: . delayed - boolean flag indicating if delayed reorthogonalization has been enabled
254: Level: advanced
256: .seealso: EPSArnoldiSetDelayed()
257: @*/
258: PetscErrorCode EPSArnoldiGetDelayed(EPS eps,PetscBool *delayed)
259: {
265: PetscTryMethod(eps,"EPSArnoldiGetDelayed_C",(EPS,PetscBool*),(eps,delayed));
266: return(0);
267: }
271: PetscErrorCode EPSDestroy_Arnoldi(EPS eps)
272: {
276: PetscFree(eps->data);
277: PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiSetDelayed_C",NULL);
278: PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiGetDelayed_C",NULL);
279: return(0);
280: }
284: PetscErrorCode EPSView_Arnoldi(EPS eps,PetscViewer viewer)
285: {
287: PetscBool isascii;
288: EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;
291: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
292: if (isascii) {
293: if (arnoldi->delayed) {
294: PetscViewerASCIIPrintf(viewer," Arnoldi: using delayed reorthogonalization\n");
295: }
296: }
297: return(0);
298: }
302: PETSC_EXTERN PetscErrorCode EPSCreate_Arnoldi(EPS eps)
303: {
304: EPS_ARNOLDI *ctx;
308: PetscNewLog(eps,&ctx);
309: eps->data = (void*)ctx;
311: eps->ops->setup = EPSSetUp_Arnoldi;
312: eps->ops->setfromoptions = EPSSetFromOptions_Arnoldi;
313: eps->ops->destroy = EPSDestroy_Arnoldi;
314: eps->ops->view = EPSView_Arnoldi;
315: eps->ops->backtransform = EPSBackTransform_Default;
316: eps->ops->computevectors = EPSComputeVectors_Schur;
317: PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiSetDelayed_C",EPSArnoldiSetDelayed_Arnoldi);
318: PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiGetDelayed_C",EPSArnoldiGetDelayed_Arnoldi);
319: return(0);
320: }