Actual source code: rqcg.c
slepc-main 2024-11-09
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: SLEPc eigensolver: "rqcg"
13: Method: Rayleigh Quotient Conjugate Gradient
15: Algorithm:
17: Conjugate Gradient minimization of the Rayleigh quotient with
18: periodic Rayleigh-Ritz acceleration.
20: References:
22: [1] L. Bergamaschi et al., "Parallel preconditioned conjugate gradient
23: optimization of the Rayleigh quotient for the solution of sparse
24: eigenproblems", Appl. Math. Comput. 175(2):1694-1715, 2006.
25: */
27: #include <slepc/private/epsimpl.h>
29: static PetscErrorCode EPSSolve_RQCG(EPS);
31: typedef struct {
32: PetscInt nrest; /* user-provided reset parameter */
33: PetscInt allocsize; /* number of columns of work BV's allocated at setup */
34: BV AV,W,P,G;
35: } EPS_RQCG;
37: static PetscErrorCode EPSSetUp_RQCG(EPS eps)
38: {
39: PetscInt nmat;
40: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
42: PetscFunctionBegin;
43: EPSCheckHermitianDefinite(eps);
44: EPSCheckNotStructured(eps);
45: PetscCall(EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd));
46: if (eps->max_it==PETSC_DETERMINE) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
47: if (!eps->which) eps->which = EPS_SMALLEST_REAL;
48: PetscCheck(eps->which==EPS_SMALLEST_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only smallest real eigenvalues");
49: EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_EXTRACTION);
50: EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);
52: if (!ctx->nrest) ctx->nrest = 20;
54: PetscCall(EPSAllocateSolution(eps,0));
55: PetscCall(EPS_SetInnerProduct(eps));
57: PetscCall(STGetNumMatrices(eps->st,&nmat));
58: if (!ctx->allocsize) {
59: ctx->allocsize = eps->mpd;
60: PetscCall(BVDuplicateResize(eps->V,eps->mpd,&ctx->AV));
61: if (nmat>1) PetscCall(BVDuplicate(ctx->AV,&ctx->W));
62: PetscCall(BVDuplicate(ctx->AV,&ctx->P));
63: PetscCall(BVDuplicate(ctx->AV,&ctx->G));
64: } else if (ctx->allocsize!=eps->mpd) {
65: ctx->allocsize = eps->mpd;
66: PetscCall(BVResize(ctx->AV,eps->mpd,PETSC_FALSE));
67: if (nmat>1) PetscCall(BVResize(ctx->W,eps->mpd,PETSC_FALSE));
68: PetscCall(BVResize(ctx->P,eps->mpd,PETSC_FALSE));
69: PetscCall(BVResize(ctx->G,eps->mpd,PETSC_FALSE));
70: }
71: PetscCall(DSSetType(eps->ds,DSHEP));
72: PetscCall(DSAllocate(eps->ds,eps->ncv));
73: PetscCall(EPSSetWorkVecs(eps,1));
74: PetscFunctionReturn(PETSC_SUCCESS);
75: }
77: static PetscErrorCode EPSSolve_RQCG(EPS eps)
78: {
79: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
80: PetscInt i,j,k,ld,nv,ncv = eps->ncv,kini,nmat;
81: PetscScalar *C,*gamma,g,pap,pbp,pbx,pax,nu,mu,alpha,beta;
82: PetscReal resnorm,a,b,c,d,disc,t;
83: PetscBool reset;
84: Mat A,B,Q,Q1;
85: Vec v,av,bv,p,w=eps->work[0];
87: PetscFunctionBegin;
88: PetscCall(DSGetLeadingDimension(eps->ds,&ld));
89: PetscCall(STGetNumMatrices(eps->st,&nmat));
90: PetscCall(STGetMatrix(eps->st,0,&A));
91: if (nmat>1) PetscCall(STGetMatrix(eps->st,1,&B));
92: else B = NULL;
93: PetscCall(PetscMalloc1(eps->mpd,&gamma));
95: kini = eps->nini;
96: while (eps->reason == EPS_CONVERGED_ITERATING) {
97: eps->its++;
98: nv = PetscMin(eps->nconv+eps->mpd,ncv);
99: PetscCall(DSSetDimensions(eps->ds,nv,eps->nconv,0));
100: for (;kini<nv;kini++) { /* Generate more initial vectors if necessary */
101: PetscCall(BVSetRandomColumn(eps->V,kini));
102: PetscCall(BVOrthonormalizeColumn(eps->V,kini,PETSC_TRUE,NULL,NULL));
103: }
104: reset = (eps->its>1 && (eps->its-1)%ctx->nrest==0)? PETSC_TRUE: PETSC_FALSE;
106: if (reset) {
107: /* Prevent BVDotVec below to use B-product, restored at the end */
108: PetscCall(BVSetMatrix(eps->V,NULL,PETSC_FALSE));
110: /* Compute Rayleigh quotient */
111: PetscCall(BVSetActiveColumns(eps->V,eps->nconv,nv));
112: PetscCall(BVSetActiveColumns(ctx->AV,0,nv-eps->nconv));
113: PetscCall(BVMatMult(eps->V,A,ctx->AV));
114: PetscCall(DSGetArray(eps->ds,DS_MAT_A,&C));
115: for (i=eps->nconv;i<nv;i++) {
116: PetscCall(BVSetActiveColumns(eps->V,eps->nconv,i+1));
117: PetscCall(BVGetColumn(ctx->AV,i-eps->nconv,&av));
118: PetscCall(BVDotVec(eps->V,av,C+eps->nconv+i*ld));
119: PetscCall(BVRestoreColumn(ctx->AV,i-eps->nconv,&av));
120: for (j=eps->nconv;j<i-1;j++) C[i+j*ld] = PetscConj(C[j+i*ld]);
121: }
122: PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&C));
123: PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
125: /* Solve projected problem */
126: PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
127: PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL));
128: PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));
130: /* Update vectors V(:,idx) = V * Y(:,idx) */
131: PetscCall(DSGetMat(eps->ds,DS_MAT_Q,&Q));
132: PetscCall(BVMultInPlace(eps->V,Q,eps->nconv,nv));
133: PetscCall(MatDenseGetSubMatrix(Q,eps->nconv,PETSC_DECIDE,eps->nconv,PETSC_DECIDE,&Q1));
134: PetscCall(BVMultInPlace(ctx->AV,Q1,0,nv-eps->nconv));
135: PetscCall(MatDenseRestoreSubMatrix(Q,&Q1));
136: PetscCall(DSRestoreMat(eps->ds,DS_MAT_Q,&Q));
137: if (B) PetscCall(BVSetMatrix(eps->V,B,PETSC_FALSE));
138: } else {
139: /* No need to do Rayleigh-Ritz, just take diag(V'*A*V) */
140: for (i=eps->nconv;i<nv;i++) {
141: PetscCall(BVGetColumn(eps->V,i,&v));
142: PetscCall(BVGetColumn(ctx->AV,i-eps->nconv,&av));
143: PetscCall(MatMult(A,v,av));
144: PetscCall(VecDot(av,v,eps->eigr+i));
145: PetscCall(BVRestoreColumn(eps->V,i,&v));
146: PetscCall(BVRestoreColumn(ctx->AV,i-eps->nconv,&av));
147: }
148: }
150: /* Compute gradient and check convergence */
151: k = -1;
152: for (i=eps->nconv;i<nv;i++) {
153: PetscCall(BVGetColumn(eps->V,i,&v));
154: PetscCall(BVGetColumn(ctx->AV,i-eps->nconv,&av));
155: PetscCall(BVGetColumn(ctx->G,i-eps->nconv,&p));
156: if (B) {
157: PetscCall(BVGetColumn(ctx->W,i-eps->nconv,&bv));
158: PetscCall(MatMult(B,v,bv));
159: PetscCall(VecWAXPY(p,-eps->eigr[i],bv,av));
160: PetscCall(BVRestoreColumn(ctx->W,i-eps->nconv,&bv));
161: } else PetscCall(VecWAXPY(p,-eps->eigr[i],v,av));
162: PetscCall(BVRestoreColumn(eps->V,i,&v));
163: PetscCall(BVRestoreColumn(ctx->AV,i-eps->nconv,&av));
164: PetscCall(VecNorm(p,NORM_2,&resnorm));
165: PetscCall(BVRestoreColumn(ctx->G,i-eps->nconv,&p));
166: PetscCall((*eps->converged)(eps,eps->eigr[i],0.0,resnorm,&eps->errest[i],eps->convergedctx));
167: if (k==-1 && eps->errest[i] >= eps->tol) k = i;
168: }
169: if (k==-1) k = nv;
170: PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx));
172: /* The next lines are necessary to avoid DS zeroing eigr */
173: PetscCall(DSGetArray(eps->ds,DS_MAT_A,&C));
174: for (i=eps->nconv;i<k;i++) C[i+i*ld] = eps->eigr[i];
175: PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&C));
177: if (eps->reason == EPS_CONVERGED_ITERATING) {
179: /* Search direction */
180: for (i=0;i<nv-eps->nconv;i++) {
181: PetscCall(BVGetColumn(ctx->G,i,&v));
182: PetscCall(STApply(eps->st,v,w));
183: PetscCall(VecDot(w,v,&g));
184: PetscCall(BVRestoreColumn(ctx->G,i,&v));
185: beta = (!reset && eps->its>1)? g/gamma[i]: 0.0;
186: gamma[i] = g;
187: PetscCall(BVGetColumn(ctx->P,i,&v));
188: PetscCall(VecAXPBY(v,1.0,beta,w));
189: if (i+eps->nconv>0) {
190: PetscCall(BVSetActiveColumns(eps->V,0,i+eps->nconv));
191: PetscCall(BVOrthogonalizeVec(eps->V,v,NULL,NULL,NULL));
192: }
193: PetscCall(BVRestoreColumn(ctx->P,i,&v));
194: }
196: /* Minimization problem */
197: for (i=eps->nconv;i<nv;i++) {
198: PetscCall(BVGetColumn(eps->V,i,&v));
199: PetscCall(BVGetColumn(ctx->AV,i-eps->nconv,&av));
200: PetscCall(BVGetColumn(ctx->P,i-eps->nconv,&p));
201: PetscCall(VecDot(av,v,&nu));
202: PetscCall(VecDot(av,p,&pax));
203: PetscCall(MatMult(A,p,w));
204: PetscCall(VecDot(w,p,&pap));
205: if (B) {
206: PetscCall(BVGetColumn(ctx->W,i-eps->nconv,&bv));
207: PetscCall(VecDot(bv,v,&mu));
208: PetscCall(VecDot(bv,p,&pbx));
209: PetscCall(BVRestoreColumn(ctx->W,i-eps->nconv,&bv));
210: PetscCall(MatMult(B,p,w));
211: PetscCall(VecDot(w,p,&pbp));
212: } else {
213: PetscCall(VecDot(v,v,&mu));
214: PetscCall(VecDot(v,p,&pbx));
215: PetscCall(VecDot(p,p,&pbp));
216: }
217: PetscCall(BVRestoreColumn(ctx->AV,i-eps->nconv,&av));
218: a = PetscRealPart(pap*pbx-pax*pbp);
219: b = PetscRealPart(nu*pbp-mu*pap);
220: c = PetscRealPart(mu*pax-nu*pbx);
221: t = PetscMax(PetscMax(PetscAbsReal(a),PetscAbsReal(b)),PetscAbsReal(c));
222: if (t!=0.0) { a /= t; b /= t; c /= t; }
223: disc = b*b-4.0*a*c;
224: d = PetscSqrtReal(PetscAbsReal(disc));
225: if (b>=0.0 && a!=0.0) alpha = (b+d)/(2.0*a);
226: else if (b!=d) alpha = 2.0*c/(b-d);
227: else alpha = 0;
228: /* Next iterate */
229: if (alpha!=0.0) PetscCall(VecAXPY(v,alpha,p));
230: PetscCall(BVRestoreColumn(eps->V,i,&v));
231: PetscCall(BVRestoreColumn(ctx->P,i-eps->nconv,&p));
232: PetscCall(BVOrthonormalizeColumn(eps->V,i,PETSC_TRUE,NULL,NULL));
233: }
234: }
236: PetscCall(EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv));
237: eps->nconv = k;
238: }
240: PetscCall(PetscFree(gamma));
241: PetscFunctionReturn(PETSC_SUCCESS);
242: }
244: static PetscErrorCode EPSRQCGSetReset_RQCG(EPS eps,PetscInt nrest)
245: {
246: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
248: PetscFunctionBegin;
249: if (nrest==PETSC_DEFAULT || nrest==PETSC_DECIDE) {
250: ctx->nrest = 0;
251: eps->state = EPS_STATE_INITIAL;
252: } else {
253: PetscCheck(nrest>0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reset parameter must be >0");
254: ctx->nrest = nrest;
255: }
256: PetscFunctionReturn(PETSC_SUCCESS);
257: }
259: /*@
260: EPSRQCGSetReset - Sets the reset parameter of the RQCG iteration. Every
261: nrest iterations, the solver performs a Rayleigh-Ritz projection step.
263: Logically Collective
265: Input Parameters:
266: + eps - the eigenproblem solver context
267: - nrest - the number of iterations between resets
269: Options Database Key:
270: . -eps_rqcg_reset - Sets the reset parameter
272: Level: advanced
274: .seealso: EPSRQCGGetReset()
275: @*/
276: PetscErrorCode EPSRQCGSetReset(EPS eps,PetscInt nrest)
277: {
278: PetscFunctionBegin;
281: PetscTryMethod(eps,"EPSRQCGSetReset_C",(EPS,PetscInt),(eps,nrest));
282: PetscFunctionReturn(PETSC_SUCCESS);
283: }
285: static PetscErrorCode EPSRQCGGetReset_RQCG(EPS eps,PetscInt *nrest)
286: {
287: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
289: PetscFunctionBegin;
290: *nrest = ctx->nrest;
291: PetscFunctionReturn(PETSC_SUCCESS);
292: }
294: /*@
295: EPSRQCGGetReset - Gets the reset parameter used in the RQCG method.
297: Not Collective
299: Input Parameter:
300: . eps - the eigenproblem solver context
302: Output Parameter:
303: . nrest - the reset parameter
305: Level: advanced
307: .seealso: EPSRQCGSetReset()
308: @*/
309: PetscErrorCode EPSRQCGGetReset(EPS eps,PetscInt *nrest)
310: {
311: PetscFunctionBegin;
313: PetscAssertPointer(nrest,2);
314: PetscUseMethod(eps,"EPSRQCGGetReset_C",(EPS,PetscInt*),(eps,nrest));
315: PetscFunctionReturn(PETSC_SUCCESS);
316: }
318: static PetscErrorCode EPSReset_RQCG(EPS eps)
319: {
320: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
322: PetscFunctionBegin;
323: PetscCall(BVDestroy(&ctx->AV));
324: PetscCall(BVDestroy(&ctx->W));
325: PetscCall(BVDestroy(&ctx->P));
326: PetscCall(BVDestroy(&ctx->G));
327: ctx->allocsize = 0;
328: PetscFunctionReturn(PETSC_SUCCESS);
329: }
331: static PetscErrorCode EPSSetFromOptions_RQCG(EPS eps,PetscOptionItems *PetscOptionsObject)
332: {
333: PetscBool flg;
334: PetscInt nrest;
336: PetscFunctionBegin;
337: PetscOptionsHeadBegin(PetscOptionsObject,"EPS RQCG Options");
339: PetscCall(PetscOptionsInt("-eps_rqcg_reset","Reset parameter","EPSRQCGSetReset",20,&nrest,&flg));
340: if (flg) PetscCall(EPSRQCGSetReset(eps,nrest));
342: PetscOptionsHeadEnd();
343: PetscFunctionReturn(PETSC_SUCCESS);
344: }
346: static PetscErrorCode EPSDestroy_RQCG(EPS eps)
347: {
348: PetscFunctionBegin;
349: PetscCall(PetscFree(eps->data));
350: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGSetReset_C",NULL));
351: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGGetReset_C",NULL));
352: PetscFunctionReturn(PETSC_SUCCESS);
353: }
355: static PetscErrorCode EPSView_RQCG(EPS eps,PetscViewer viewer)
356: {
357: EPS_RQCG *ctx = (EPS_RQCG*)eps->data;
358: PetscBool isascii;
360: PetscFunctionBegin;
361: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
362: if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer," reset every %" PetscInt_FMT " iterations\n",ctx->nrest));
363: PetscFunctionReturn(PETSC_SUCCESS);
364: }
366: SLEPC_EXTERN PetscErrorCode EPSCreate_RQCG(EPS eps)
367: {
368: EPS_RQCG *rqcg;
370: PetscFunctionBegin;
371: PetscCall(PetscNew(&rqcg));
372: eps->data = (void*)rqcg;
374: eps->useds = PETSC_TRUE;
375: eps->categ = EPS_CATEGORY_PRECOND;
377: eps->ops->solve = EPSSolve_RQCG;
378: eps->ops->setup = EPSSetUp_RQCG;
379: eps->ops->setupsort = EPSSetUpSort_Default;
380: eps->ops->setfromoptions = EPSSetFromOptions_RQCG;
381: eps->ops->destroy = EPSDestroy_RQCG;
382: eps->ops->reset = EPSReset_RQCG;
383: eps->ops->view = EPSView_RQCG;
384: eps->ops->backtransform = EPSBackTransform_Default;
385: eps->ops->setdefaultst = EPSSetDefaultST_GMRES;
387: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGSetReset_C",EPSRQCGSetReset_RQCG));
388: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGGetReset_C",EPSRQCGGetReset_RQCG));
389: PetscFunctionReturn(PETSC_SUCCESS);
390: }