Actual source code: rqcg.c
slepc-3.5.2 2014-10-10
1: /*
3: SLEPc eigensolver: "rqcg"
5: Method: Rayleigh Quotient Conjugate Gradient
7: Algorithm:
9: Conjugate Gradient minimization of the Rayleigh quotient with
10: periodic Rayleigh-Ritz acceleration.
12: References:
14: [1] L. Bergamaschi et al., "Parallel preconditioned conjugate gradient
15: optimization of the Rayleigh quotient for the solution of sparse
16: eigenproblems", Appl. Math. Comput. 175(2):1694-1715, 2006.
18: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
19: SLEPc - Scalable Library for Eigenvalue Problem Computations
20: Copyright (c) 2002-2014, 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*/
40: PetscErrorCode EPSSolve_RQCG(EPS);
42: typedef struct {
43: PetscInt nrest;
44: BV AV,W,P,G;
45: } EPS_RQCG;
49: PetscErrorCode EPSSetUp_RQCG(EPS eps)
50: {
52: PetscBool precond;
53: PetscInt nmat;
54: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
57: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"RQCG only works for Hermitian problems");
58: EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
59: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
60: if (!eps->which) eps->which = EPS_SMALLEST_REAL;
61: if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
62: if (!eps->extraction) {
63: EPSSetExtraction(eps,EPS_RITZ);
64: } else if (eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
65: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
66: /* Set STPrecond as the default ST */
67: if (!((PetscObject)eps->st)->type_name) {
68: STSetType(eps->st,STPRECOND);
69: }
70: PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&precond);
71: if (!precond) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"RQCG only works with precond ST");
73: if (!ctx->nrest) ctx->nrest = 20;
75: EPSAllocateSolution(eps,0);
76: EPS_SetInnerProduct(eps);
77: BVDuplicateResize(eps->V,eps->mpd,&ctx->AV);
78: PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->AV);
79: STGetNumMatrices(eps->st,&nmat);
80: if (nmat>1) {
81: BVDuplicate(ctx->AV,&ctx->W);
82: PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->W);
83: }
84: BVDuplicate(ctx->AV,&ctx->P);
85: PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->P);
86: BVDuplicate(ctx->AV,&ctx->G);
87: PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->G);
88: DSSetType(eps->ds,DSHEP);
89: DSAllocate(eps->ds,eps->ncv);
90: EPSSetWorkVecs(eps,1);
91: return(0);
92: }
96: /*
97: ExtractSubmatrix - Returns B = A(k+1:end,k+1:end).
98: */
99: static PetscErrorCode ExtractSubmatrix(Mat A,PetscInt k,Mat *B)
100: {
102: PetscInt j,m,n;
103: PetscScalar *pA,*pB;
106: MatGetSize(A,&m,&n);
107: MatCreateSeqDense(PETSC_COMM_SELF,m-k,n-k,NULL,B);
108: MatDenseGetArray(A,&pA);
109: MatDenseGetArray(*B,&pB);
110: for (j=k;j<n;j++) {
111: PetscMemcpy(pB+(j-k)*(m-k),pA+j*m+k,(m-k)*sizeof(PetscScalar));
112: }
113: MatDenseRestoreArray(A,&pA);
114: MatDenseRestoreArray(*B,&pB);
115: return(0);
116: }
120: PetscErrorCode EPSSolve_RQCG(EPS eps)
121: {
123: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
124: PetscInt i,j,k,ld,nv,ncv = eps->ncv,kini,nmat;
125: PetscScalar *C,*gamma,g,pap,pbp,pbx,pax,nu,mu,alpha,beta;
126: PetscReal resnorm,norm,a,b,c,disc,t;
127: PetscBool reset,breakdown;
128: Mat A,B,Q,Q1;
129: Vec v,av,bv,p,w=eps->work[0];
132: DSGetLeadingDimension(eps->ds,&ld);
133: STGetNumMatrices(eps->st,&nmat);
134: STGetOperators(eps->st,0,&A);
135: if (nmat>1) { STGetOperators(eps->st,1,&B); }
136: else B = NULL;
137: PetscMalloc1(eps->mpd,&gamma);
139: kini = eps->nini;
140: while (eps->reason == EPS_CONVERGED_ITERATING) {
141: eps->its++;
142: nv = PetscMin(eps->nconv+eps->mpd,ncv);
143: DSSetDimensions(eps->ds,nv,0,eps->nconv,0);
144: /* Generate more initial vectors if necessary */
145: while (kini<nv) {
146: BVSetRandomColumn(eps->V,kini,eps->rand);
147: BVOrthogonalizeColumn(eps->V,kini,NULL,&norm,&breakdown);
148: if (norm>0.0 && !breakdown) {
149: BVScaleColumn(eps->V,kini,1.0/norm);
150: kini++;
151: }
152: }
153: reset = (eps->its>1 && (eps->its-1)%ctx->nrest==0)? PETSC_TRUE: PETSC_FALSE;
155: if (reset) {
156: /* Prevent BVDotVec below to use B-product, restored a the end */
157: BVSetMatrix(eps->V,NULL,PETSC_FALSE);
159: /* Compute Rayleigh quotient */
160: BVSetActiveColumns(eps->V,eps->nconv,nv);
161: BVSetActiveColumns(ctx->AV,0,nv-eps->nconv);
162: BVMatMult(eps->V,A,ctx->AV);
163: DSGetArray(eps->ds,DS_MAT_A,&C);
164: for (i=eps->nconv;i<nv;i++) {
165: BVSetActiveColumns(eps->V,eps->nconv,i+1);
166: BVGetColumn(ctx->AV,i-eps->nconv,&av);
167: BVDotVec(eps->V,av,C+eps->nconv+i*ld);
168: BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
169: for (j=eps->nconv;j<i-1;j++) C[i+j*ld] = C[j+i*ld];
170: }
171: DSRestoreArray(eps->ds,DS_MAT_A,&C);
172: DSSetState(eps->ds,DS_STATE_RAW);
174: /* Solve projected problem */
175: DSSolve(eps->ds,eps->eigr,eps->eigi);
176: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
178: /* Update vectors V(:,idx) = V * Y(:,idx) */
179: DSGetMat(eps->ds,DS_MAT_Q,&Q);
180: BVMultInPlace(eps->V,Q,eps->nconv,nv);
181: ExtractSubmatrix(Q,eps->nconv,&Q1);
182: BVMultInPlace(ctx->AV,Q1,0,nv-eps->nconv);
183: MatDestroy(&Q);
184: MatDestroy(&Q1);
185: if (B) { BVSetMatrix(eps->V,B,PETSC_FALSE); }
186: } else {
187: /* No need to do Rayleigh-Ritz, just take diag(V'*A*V) */
188: for (i=eps->nconv;i<nv;i++) {
189: BVGetColumn(eps->V,i,&v);
190: BVGetColumn(ctx->AV,i-eps->nconv,&av);
191: MatMult(A,v,av);
192: VecDot(av,v,eps->eigr+i);
193: BVRestoreColumn(eps->V,i,&v);
194: BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
195: }
196: }
198: /* Compute gradient and check convergence */
199: k = -1;
200: for (i=eps->nconv;i<nv;i++) {
201: BVGetColumn(eps->V,i,&v);
202: BVGetColumn(ctx->AV,i-eps->nconv,&av);
203: BVGetColumn(ctx->G,i-eps->nconv,&p);
204: if (B) {
205: BVGetColumn(ctx->W,i-eps->nconv,&bv);
206: MatMult(B,v,bv);
207: VecWAXPY(p,-eps->eigr[i],bv,av);
208: BVRestoreColumn(ctx->W,i-eps->nconv,&bv);
209: } else {
210: VecWAXPY(p,-eps->eigr[i],v,av);
211: }
212: BVRestoreColumn(eps->V,i,&v);
213: BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
214: VecNorm(p,NORM_2,&resnorm);
215: BVRestoreColumn(ctx->G,i-eps->nconv,&p);
216: (*eps->converged)(eps,eps->eigr[i],0.0,resnorm,&eps->errest[i],eps->convergedctx);
217: if (k==-1 && eps->errest[i] >= eps->tol) k = i;
218: }
219: if (k==-1) k = nv;
220: if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
221: if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
223: /* The next lines are necessary to avoid DS zeroing eigr */
224: DSGetArray(eps->ds,DS_MAT_A,&C);
225: for (i=eps->nconv;i<k;i++) C[i+i*ld] = eps->eigr[i];
226: DSRestoreArray(eps->ds,DS_MAT_A,&C);
228: if (eps->reason == EPS_CONVERGED_ITERATING) {
230: /* Search direction */
231: for (i=0;i<nv-eps->nconv;i++) {
232: BVGetColumn(ctx->G,i,&v);
233: STMatSolve(eps->st,v,w);
234: VecDot(v,w,&g);
235: BVRestoreColumn(ctx->G,i,&v);
236: beta = (!reset && eps->its>1)? g/gamma[i]: 0.0;
237: gamma[i] = g;
238: BVGetColumn(ctx->P,i,&v);
239: VecAXPBY(v,1.0,beta,w);
240: if (i+eps->nconv>0) {
241: BVSetActiveColumns(eps->V,0,i+eps->nconv);
242: BVOrthogonalizeVec(eps->V,v,NULL,NULL,NULL);
243: }
244: BVRestoreColumn(ctx->P,i,&v);
245: }
247: /* Minimization problem */
248: for (i=eps->nconv;i<nv;i++) {
249: BVGetColumn(eps->V,i,&v);
250: BVGetColumn(ctx->AV,i-eps->nconv,&av);
251: BVGetColumn(ctx->P,i-eps->nconv,&p);
252: VecDot(v,av,&nu);
253: VecDot(p,av,&pax);
254: MatMult(A,p,w);
255: VecDot(p,w,&pap);
256: if (B) {
257: BVGetColumn(ctx->W,i-eps->nconv,&bv);
258: VecDot(v,bv,&mu);
259: VecDot(p,bv,&pbx);
260: BVRestoreColumn(ctx->W,i-eps->nconv,&bv);
261: MatMult(B,p,w);
262: VecDot(p,w,&pbp);
263: } else {
264: VecDot(v,v,&mu);
265: VecDot(p,v,&pbx);
266: VecDot(p,p,&pbp);
267: }
268: BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
269: a = PetscRealPart(pap*pbx-pax*pbp);
270: b = PetscRealPart(nu*pbp-mu*pap);
271: c = PetscRealPart(mu*pax-nu*pbx);
272: t = PetscMax(PetscMax(PetscAbsReal(a),PetscAbsReal(b)),PetscAbsReal(c));
273: if (t!=0.0) { a /= t; b /= t; c /= t; }
274: disc = PetscSqrtReal(PetscAbsReal(b*b-4.0*a*c));
275: if (b>=0.0 && a!=0.0) alpha = (b+disc)/(2.0*a);
276: else if (b!=disc) alpha = 2.0*c/(b-disc);
277: else alpha = 0;
278: /* Next iterate */
279: if (alpha!=0.0) {
280: VecAXPY(v,alpha,p);
281: }
282: BVRestoreColumn(eps->V,i,&v);
283: BVRestoreColumn(ctx->P,i-eps->nconv,&p);
284: BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&breakdown);
285: if (!breakdown && norm!=0.0) {
286: BVScaleColumn(eps->V,i,1.0/norm);
287: }
288: }
289: }
291: EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv);
292: eps->nconv = k;
293: }
295: PetscFree(gamma);
296: /* truncate Schur decomposition and change the state to raw so that
297: PSVectors() computes eigenvectors from scratch */
298: DSSetDimensions(eps->ds,eps->nconv,0,0,0);
299: DSSetState(eps->ds,DS_STATE_RAW);
300: return(0);
301: }
305: static PetscErrorCode EPSRQCGSetReset_RQCG(EPS eps,PetscInt nrest)
306: {
307: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
310: ctx->nrest = nrest;
311: return(0);
312: }
316: /*@
317: EPSRQCGSetReset - Sets the reset parameter of the RQCG iteration. Every
318: nrest iterations, the solver performs a Rayleigh-Ritz projection step.
320: Logically Collective on EPS
322: Input Parameters:
323: + eps - the eigenproblem solver context
324: - nrest - the number of iterations between resets
326: Options Database Key:
327: . -eps_rqcg_reset - Sets the reset parameter
329: Level: advanced
331: .seealso: EPSRQCGGetReset()
332: @*/
333: PetscErrorCode EPSRQCGSetReset(EPS eps,PetscInt nrest)
334: {
340: PetscTryMethod(eps,"EPSRQCGSetReset_C",(EPS,PetscInt),(eps,nrest));
341: return(0);
342: }
346: static PetscErrorCode EPSRQCGGetReset_RQCG(EPS eps,PetscInt *nrest)
347: {
348: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
351: *nrest = ctx->nrest;
352: return(0);
353: }
357: /*@
358: EPSRQCGGetReset - Gets the reset parameter used in the RQCG method.
360: Not Collective
362: Input Parameter:
363: . eps - the eigenproblem solver context
365: Output Parameter:
366: . nrest - the reset parameter
368: Level: advanced
370: .seealso: EPSRQCGSetReset()
371: @*/
372: PetscErrorCode EPSRQCGGetReset(EPS eps,PetscInt *nrest)
373: {
379: PetscTryMethod(eps,"EPSRQCGGetReset_C",(EPS,PetscInt*),(eps,nrest));
380: return(0);
381: }
385: PetscErrorCode EPSReset_RQCG(EPS eps)
386: {
388: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
391: BVDestroy(&ctx->AV);
392: BVDestroy(&ctx->W);
393: BVDestroy(&ctx->P);
394: BVDestroy(&ctx->G);
395: ctx->nrest = 0;
396: return(0);
397: }
401: PetscErrorCode EPSSetFromOptions_RQCG(EPS eps)
402: {
404: PetscBool flg;
405: PetscInt nrest;
408: PetscOptionsHead("EPS RQCG Options");
409: PetscOptionsInt("-eps_rqcg_reset","RQCG reset parameter","EPSRQCGSetReset",20,&nrest,&flg);
410: if (flg) {
411: EPSRQCGSetReset(eps,nrest);
412: }
413: PetscOptionsTail();
414: return(0);
415: }
419: PetscErrorCode EPSDestroy_RQCG(EPS eps)
420: {
424: PetscFree(eps->data);
425: PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGSetReset_C",NULL);
426: PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGGetReset_C",NULL);
427: return(0);
428: }
432: PetscErrorCode EPSView_RQCG(EPS eps,PetscViewer viewer)
433: {
435: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
436: PetscBool isascii;
439: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
440: if (isascii) {
441: PetscViewerASCIIPrintf(viewer," RQCG: reset every %D iterations\n",ctx->nrest);
442: }
443: return(0);
444: }
448: PETSC_EXTERN PetscErrorCode EPSCreate_RQCG(EPS eps)
449: {
450: EPS_RQCG *rqcg;
454: PetscNewLog(eps,&rqcg);
455: eps->data = (void*)rqcg;
457: eps->ops->setup = EPSSetUp_RQCG;
458: eps->ops->solve = EPSSolve_RQCG;
459: eps->ops->setfromoptions = EPSSetFromOptions_RQCG;
460: eps->ops->destroy = EPSDestroy_RQCG;
461: eps->ops->reset = EPSReset_RQCG;
462: eps->ops->view = EPSView_RQCG;
463: eps->ops->backtransform = EPSBackTransform_Default;
464: STSetType(eps->st,STPRECOND);
465: STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
466: PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGSetReset_C",EPSRQCGSetReset_RQCG);
467: PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGGetReset_C",EPSRQCGGetReset_RQCG);
468: return(0);
469: }