Actual source code: rqcg.c
1: /*
2:
3: SLEPc eigensolver: "rqcg"
4:
5: Method: Rayleigh Quotient Conjugate Gradient
6:
7: Algorithm:
8:
9: Conjugate Gradient minimization of the Rayleigh quotient with
10: periodic Rayleigh-Ritz acceleration.
11:
12: References:
13:
14: [1] "Parallel preconditioned conjugate gradient optimization of the
15: Rayleigh quotient for the solution of sparse eigenproblems",
16: Appl. Math. Comput. 175, pp. 1694-1715 (2006).
17:
18: Last update: Jul 2012
19:
20: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
21: SLEPc - Scalable Library for Eigenvalue Problem Computations
22: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
23:
24: This file is part of SLEPc.
25:
26: SLEPc is free software: you can redistribute it and/or modify it under the
27: terms of version 3 of the GNU Lesser General Public License as published by
28: the Free Software Foundation.
29:
30: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
31: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
32: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
33: more details.
34:
35: You should have received a copy of the GNU Lesser General Public License
36: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: */
39:
40: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
41: #include <slepcblaslapack.h>
42:
43: PetscErrorCode EPSSolve_RQCG(EPS);
44:
45: typedef struct {
46: PetscInt nrest;
47: Vec *AV,*BV,*P,*G;
48: } EPS_RQCG;
49:
52: PetscErrorCode EPSSetUp_RQCG(EPS eps)
53: {
55: PetscBool precond;
56: Mat B;
57: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
58:
60: if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"RQCG only works for Hermitian problems");
61: if (eps->ncv) { /* ncv set */
62: if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
63: }
64: else if (eps->mpd) { /* mpd set */
65: eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
66: }
67: else { /* neither set: defaults depend on nev being small or large */
68: if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
69: else { eps->mpd = 500; eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd); }
70: }
71: if (!eps->mpd) eps->mpd = eps->ncv;
72: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
73: if (!eps->which) eps->which = EPS_SMALLEST_REAL;
74: if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
75: if (!eps->extraction) {
76: EPSSetExtraction(eps,EPS_RITZ);
77: } else if (eps->extraction!=EPS_RITZ) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Unsupported extraction type");
78: if (eps->arbit_func) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
79: /* Set STPrecond as the default ST */
80: if (!((PetscObject)eps->OP)->type_name) {
81: STSetType(eps->OP,STPRECOND);
82: }
83: PetscObjectTypeCompare((PetscObject)eps->OP,STPRECOND,&precond);
84: if (!precond) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"RQCG only works with precond ST");
85:
86: if (!ctx->nrest) ctx->nrest = 20;
87:
88: EPSAllocateSolution(eps);
89: VecDuplicateVecs(eps->t,eps->mpd,&ctx->AV);
90: STGetOperators(eps->OP,PETSC_NULL,&B);
91: if (B) {
92: VecDuplicateVecs(eps->t,eps->mpd,&ctx->BV);
93: }
94: VecDuplicateVecs(eps->t,eps->mpd,&ctx->P);
95: VecDuplicateVecs(eps->t,eps->mpd,&ctx->G);
96: DSSetType(eps->ds,DSHEP);
97: DSAllocate(eps->ds,eps->ncv);
98: EPSDefaultGetWork(eps,1);
99:
100: /* dispatch solve method */
101: if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
102: eps->ops->solve = EPSSolve_RQCG;
103: return(0);
104: }
105:
108: PetscErrorCode EPSSolve_RQCG(EPS eps)
109: {
111: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
112: PetscInt i,j,k,ld,off,nv,ncv = eps->ncv,kini;
113: PetscScalar *C,*Y,*gamma,g,pap,pbp,pbx,pax,nu,mu,alpha,beta;
114: PetscReal resnorm,norm,a,b,c,disc,t;
115: PetscBool reset,breakdown;
116: Mat A,B;
117: Vec w=eps->work[0];
118:
120: DSGetLeadingDimension(eps->ds,&ld);
121: STGetOperators(eps->OP,&A,&B);
122: PetscMalloc(eps->mpd*sizeof(PetscScalar),&gamma);
123:
124: kini = eps->nini;
125: while (eps->reason == EPS_CONVERGED_ITERATING) {
126: eps->its++;
127: nv = PetscMin(eps->nconv+eps->mpd,ncv);
128: DSSetDimensions(eps->ds,nv,PETSC_IGNORE,eps->nconv,0);
129: /* Generate more initial vectors if necessary */
130: while (kini<nv) {
131: SlepcVecSetRandom(eps->V[kini],eps->rand);
132: IPOrthogonalize(eps->ip,eps->nds,eps->defl,kini,PETSC_NULL,eps->V,eps->V[kini],PETSC_NULL,&norm,&breakdown);
133: if (norm>0.0 && !breakdown) {
134: VecScale(eps->V[kini],1.0/norm);
135: kini++;
136: }
137: }
138: reset = (eps->its>1 && (eps->its-1)%ctx->nrest==0)? PETSC_TRUE: PETSC_FALSE;
139:
140: if (reset) {
141: /* Compute Rayleigh quotient */
142: DSGetArray(eps->ds,DS_MAT_A,&C);
143: for (i=eps->nconv;i<nv;i++) {
144: MatMult(A,eps->V[i],ctx->AV[i-eps->nconv]);
145: VecMDot(ctx->AV[i-eps->nconv],i-eps->nconv+1,eps->V+eps->nconv,C+eps->nconv+i*ld);
146: for (j=eps->nconv;j<i-1;j++) C[i+j*ld] = C[j+i*ld];
147: }
148: DSRestoreArray(eps->ds,DS_MAT_A,&C);
149: DSSetState(eps->ds,DS_STATE_RAW);
150:
151: /* Solve projected problem */
152: DSSolve(eps->ds,eps->eigr,eps->eigi);
153: DSSort(eps->ds,eps->eigr,eps->eigi,PETSC_NULL,PETSC_NULL,PETSC_NULL);
154:
155: /* Update vectors V(:,idx) = V * Y(:,idx) */
156: DSGetArray(eps->ds,DS_MAT_Q,&Y);
157: off = eps->nconv+eps->nconv*ld;
158: SlepcUpdateVectors(nv-eps->nconv,ctx->AV,0,nv-eps->nconv,Y+off,ld,PETSC_FALSE);
159: SlepcUpdateVectors(nv-eps->nconv,eps->V+eps->nconv,0,nv-eps->nconv,Y+off,ld,PETSC_FALSE);
160: DSRestoreArray(eps->ds,DS_MAT_Q,&Y);
161:
162: } else {
163: /* No need to do Rayleigh-Ritz, just take diag(V'*A*V) */
164: for (i=eps->nconv;i<nv;i++) {
165: MatMult(A,eps->V[i],ctx->AV[i-eps->nconv]);
166: VecDot(ctx->AV[i-eps->nconv],eps->V[i],eps->eigr+i);
167: }
168: }
169:
170: /* Compute gradient and check convergence */
171: k = -1;
172: for (i=eps->nconv;i<nv;i++) {
173: if (B) {
174: MatMult(B,eps->V[i],ctx->BV[i-eps->nconv]);
175: VecWAXPY(ctx->G[i-eps->nconv],-eps->eigr[i],ctx->BV[i-eps->nconv],ctx->AV[i-eps->nconv]);
176: } else {
177: VecWAXPY(ctx->G[i-eps->nconv],-eps->eigr[i],eps->V[i],ctx->AV[i-eps->nconv]);
178: }
179: VecNorm(ctx->G[i-eps->nconv],NORM_2,&resnorm);
180: (*eps->conv_func)(eps,eps->eigr[i],0.0,resnorm,&eps->errest[i],eps->conv_ctx);
181: if (k==-1 && eps->errest[i] >= eps->tol) k = i;
182: }
183: if (k==-1) k = nv;
184: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
185: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
186:
187: /* The next lines are necessary to avoid DS zeroing eigr */
188: DSGetArray(eps->ds,DS_MAT_A,&C);
189: for (i=eps->nconv;i<k;i++) C[i+i*ld] = eps->eigr[i];
190: DSRestoreArray(eps->ds,DS_MAT_A,&C);
191:
192: if (eps->reason == EPS_CONVERGED_ITERATING) {
193:
194: /* Search direction */
195: for (i=0;i<nv-eps->nconv;i++) {
196: STAssociatedKSPSolve(eps->OP,ctx->G[i],w);
197: VecDot(ctx->G[i],w,&g);
198: beta = (!reset && eps->its>1)? g/gamma[i]: 0.0;
199: gamma[i] = g;
200: VecAXPBY(ctx->P[i],1.0,beta,w);
201: IPOrthogonalize(eps->ip,eps->nds,eps->defl,i+eps->nconv,PETSC_NULL,eps->V,ctx->P[i],PETSC_NULL,&resnorm,&breakdown);
202: }
203:
204: /* Minimization problem */
205: for (i=eps->nconv;i<nv;i++) {
206: VecDot(eps->V[i],ctx->AV[i-eps->nconv],&nu);
207: VecDot(ctx->P[i-eps->nconv],ctx->AV[i-eps->nconv],&pax);
208: MatMult(A,ctx->P[i-eps->nconv],w);
209: VecDot(ctx->P[i-eps->nconv],w,&pap);
210: if (B) {
211: VecDot(eps->V[i],ctx->BV[i-eps->nconv],&mu);
212: VecDot(ctx->P[i-eps->nconv],ctx->BV[i-eps->nconv],&pbx);
213: MatMult(B,ctx->P[i-eps->nconv],w);
214: VecDot(ctx->P[i-eps->nconv],w,&pbp);
215: } else {
216: VecDot(eps->V[i],eps->V[i],&mu);
217: VecDot(ctx->P[i-eps->nconv],eps->V[i],&pbx);
218: VecDot(ctx->P[i-eps->nconv],ctx->P[i-eps->nconv],&pbp);
219: }
220: a = PetscRealPart(pap*pbx-pax*pbp);
221: b = PetscRealPart(nu*pbp-mu*pap);
222: c = PetscRealPart(mu*pax-nu*pbx);
223: t = PetscMax(PetscMax(PetscAbsReal(a),PetscAbsReal(b)),PetscAbsReal(c));
224: if (t!=0.0) { a /= t; b /= t; c /= t; }
225: disc = PetscSqrtReal(PetscAbsReal(b*b-4.0*a*c));
226: if (b>=0.0 && a!=0.0) alpha = (b+disc)/(2.0*a);
227: else if(b!=disc) alpha = 2.0*c/(b-disc);
228: else alpha = 0;
229: /* Next iterate */
230: if (alpha!=0.0){
231: VecAXPY(eps->V[i],alpha,ctx->P[i-eps->nconv]);
232: }
233: IPOrthogonalize(eps->ip,eps->nds,eps->defl,i,PETSC_NULL,eps->V,eps->V[i],PETSC_NULL,&norm,&breakdown);
234: if (!breakdown && norm!=0.0) {
235: VecScale(eps->V[i],1.0/norm);
236: }
237: }
238: }
239:
240: EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv);
241: eps->nconv = k;
242: }
243:
244: PetscFree(gamma);
245: /* truncate Schur decomposition and change the state to raw so that
246: PSVectors() computes eigenvectors from scratch */
247: DSSetDimensions(eps->ds,eps->nconv,PETSC_IGNORE,0,0);
248: DSSetState(eps->ds,DS_STATE_RAW);
249: return(0);
250: }
251:
252: EXTERN_C_BEGIN
255: PetscErrorCode EPSRQCGSetReset_RQCG(EPS eps,PetscInt nrest)
256: {
257: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
258:
260: ctx->nrest = nrest;
261: return(0);
262: }
263: EXTERN_C_END
264:
267: /*@
268: EPSRQCGSetReset - Sets the reset parameter of the RQCG iteration. Every
269: nrest iterations, the solver performs a Rayleigh-Ritz projection step.
270:
271: Logically Collective on EPS
272:
273: Input Parameters:
274: + eps - the eigenproblem solver context
275: - nrest - the number of iterations between resets
276:
277: Options Database Key:
278: . -eps_rqcg_reset - Sets the reset parameter
279:
280: Level: advanced
281:
282: .seealso: EPSRQCGGetReset()
283: @*/
284: PetscErrorCode EPSRQCGSetReset(EPS eps,PetscInt nrest)
285: {
287:
291: PetscTryMethod(eps,"EPSRQCGSetReset_C",(EPS,PetscInt),(eps,nrest));
292: return(0);
293: }
294:
295: EXTERN_C_BEGIN
298: PetscErrorCode EPSRQCGGetReset_RQCG(EPS eps,PetscInt *nrest)
299: {
300: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
301:
303: *nrest = ctx->nrest;
304: return(0);
305: }
306: EXTERN_C_END
307:
310: /*@
311: EPSRQCGGetReset - Gets the reset parameter used in the RQCG method.
312:
313: Not Collective
314:
315: Input Parameter:
316: . eps - the eigenproblem solver context
317:
318: Output Parameter:
319: . nrest - the reset parameter
320:
321: Level: advanced
322:
323: .seealso: EPSRQCGSetReset()
324: @*/
325: PetscErrorCode EPSRQCGGetReset(EPS eps,PetscInt *nrest)
326: {
328:
332: PetscTryMethod(eps,"EPSRQCGGetReset_C",(EPS,PetscInt*),(eps,nrest));
333: return(0);
334: }
335:
338: PetscErrorCode EPSReset_RQCG(EPS eps)
339: {
341: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
342:
344: VecDestroyVecs(eps->mpd,&ctx->AV);
345: VecDestroyVecs(eps->mpd,&ctx->BV);
346: VecDestroyVecs(eps->mpd,&ctx->P);
347: VecDestroyVecs(eps->mpd,&ctx->G);
348: ctx->nrest = 0;
349: EPSReset_Default(eps);
350: return(0);
351: }
352:
355: PetscErrorCode EPSSetFromOptions_RQCG(EPS eps)
356: {
358: PetscBool flg;
359: PetscInt nrest;
360:
362: PetscOptionsHead("EPS RQCG Options");
363: PetscOptionsInt("-eps_rqcg_reset","RQCG reset parameter","EPSRQCGSetReset",20,&nrest,&flg);
364: if (flg) { EPSRQCGSetReset(eps,nrest); }
365: PetscOptionsTail();
366: return(0);
367: }
368:
371: PetscErrorCode EPSDestroy_RQCG(EPS eps)
372: {
374:
376: PetscFree(eps->data);
377: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSRQCGSetReset_C","",PETSC_NULL);
378: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSRQCGGetReset_C","",PETSC_NULL);
379: return(0);
380: }
381:
384: PetscErrorCode EPSView_RQCG(EPS eps,PetscViewer viewer)
385: {
387: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
388: PetscBool isascii;
389:
391: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
392: if (!isascii) SETERRQ1(((PetscObject)eps)->comm,1,"Viewer type %s not supported for EPS RQCG",((PetscObject)viewer)->type_name);
393: PetscViewerASCIIPrintf(viewer," RQCG: reset every %d iterations\n",ctx->nrest);
394: return(0);
395: }
396:
397: EXTERN_C_BEGIN
400: PetscErrorCode EPSCreate_RQCG(EPS eps)
401: {
403:
405: PetscNewLog(eps,EPS_RQCG,&eps->data);
406: eps->ops->setup = EPSSetUp_RQCG;
407: eps->ops->setfromoptions = EPSSetFromOptions_RQCG;
408: eps->ops->destroy = EPSDestroy_RQCG;
409: eps->ops->reset = EPSReset_RQCG;
410: eps->ops->view = EPSView_RQCG;
411: eps->ops->backtransform = EPSBackTransform_Default;
412: eps->ops->computevectors = EPSComputeVectors_Default;
413: STSetType(eps->OP,STPRECOND);
414: STPrecondSetKSPHasMat(eps->OP,PETSC_TRUE);
415: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSRQCGSetReset_C","EPSRQCGSetReset_RQCG",EPSRQCGSetReset_RQCG);
416: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSRQCGGetReset_C","EPSRQCGGetReset_RQCG",EPSRQCGGetReset_RQCG);
417: return(0);
418: }
419: EXTERN_C_END
420: