Actual source code: rqcg.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  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: 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: PetscErrorCode EPSSetUp_RQCG(EPS eps)
 38: {
 39:   PetscInt       nmat;
 40:   EPS_RQCG       *ctx = (EPS_RQCG*)eps->data;

 42:   EPSCheckHermitianDefinite(eps);
 43:   EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
 44:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 45:   if (!eps->which) eps->which = EPS_SMALLEST_REAL;
 47:   EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_EXTRACTION);
 48:   EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);

 50:   if (!ctx->nrest) ctx->nrest = 20;

 52:   EPSAllocateSolution(eps,0);
 53:   EPS_SetInnerProduct(eps);

 55:   STGetNumMatrices(eps->st,&nmat);
 56:   if (!ctx->allocsize) {
 57:     ctx->allocsize = eps->mpd;
 58:     BVDuplicateResize(eps->V,eps->mpd,&ctx->AV);
 59:     PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->AV);
 60:     if (nmat>1) {
 61:       BVDuplicate(ctx->AV,&ctx->W);
 62:       PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->W);
 63:     }
 64:     BVDuplicate(ctx->AV,&ctx->P);
 65:     PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->P);
 66:     BVDuplicate(ctx->AV,&ctx->G);
 67:     PetscLogObjectParent((PetscObject)eps,(PetscObject)ctx->G);
 68:   } else if (ctx->allocsize!=eps->mpd) {
 69:     ctx->allocsize = eps->mpd;
 70:     BVResize(ctx->AV,eps->mpd,PETSC_FALSE);
 71:     if (nmat>1) BVResize(ctx->W,eps->mpd,PETSC_FALSE);
 72:     BVResize(ctx->P,eps->mpd,PETSC_FALSE);
 73:     BVResize(ctx->G,eps->mpd,PETSC_FALSE);
 74:   }
 75:   DSSetType(eps->ds,DSHEP);
 76:   DSAllocate(eps->ds,eps->ncv);
 77:   EPSSetWorkVecs(eps,1);
 78:   PetscFunctionReturn(0);
 79: }

 81: /*
 82:    ExtractSubmatrix - Returns B = A(k+1:end,k+1:end).
 83: */
 84: static PetscErrorCode ExtractSubmatrix(Mat A,PetscInt k,Mat *B)
 85: {
 86:   PetscInt          j,m,n;
 87:   PetscScalar       *pB;
 88:   const PetscScalar *pA;

 90:   MatGetSize(A,&m,&n);
 91:   MatCreateSeqDense(PETSC_COMM_SELF,m-k,n-k,NULL,B);
 92:   MatDenseGetArrayRead(A,&pA);
 93:   MatDenseGetArrayWrite(*B,&pB);
 94:   for (j=k;j<n;j++) PetscArraycpy(pB+(j-k)*(m-k),pA+j*m+k,m-k);
 95:   MatDenseRestoreArrayRead(A,&pA);
 96:   MatDenseRestoreArrayWrite(*B,&pB);
 97:   PetscFunctionReturn(0);
 98: }

100: PetscErrorCode EPSSolve_RQCG(EPS eps)
101: {
102:   EPS_RQCG       *ctx = (EPS_RQCG*)eps->data;
103:   PetscInt       i,j,k,ld,nv,ncv = eps->ncv,kini,nmat;
104:   PetscScalar    *C,*gamma,g,pap,pbp,pbx,pax,nu,mu,alpha,beta;
105:   PetscReal      resnorm,a,b,c,d,disc,t;
106:   PetscBool      reset;
107:   Mat            A,B,Q,Q1;
108:   Vec            v,av,bv,p,w=eps->work[0];

110:   DSGetLeadingDimension(eps->ds,&ld);
111:   STGetNumMatrices(eps->st,&nmat);
112:   STGetMatrix(eps->st,0,&A);
113:   if (nmat>1) STGetMatrix(eps->st,1,&B);
114:   else B = NULL;
115:   PetscMalloc1(eps->mpd,&gamma);

117:   kini = eps->nini;
118:   while (eps->reason == EPS_CONVERGED_ITERATING) {
119:     eps->its++;
120:     nv = PetscMin(eps->nconv+eps->mpd,ncv);
121:     DSSetDimensions(eps->ds,nv,eps->nconv,0);
122:     for (;kini<nv;kini++) { /* Generate more initial vectors if necessary */
123:       BVSetRandomColumn(eps->V,kini);
124:       BVOrthonormalizeColumn(eps->V,kini,PETSC_TRUE,NULL,NULL);
125:     }
126:     reset = (eps->its>1 && (eps->its-1)%ctx->nrest==0)? PETSC_TRUE: PETSC_FALSE;

128:     if (reset) {
129:       /* Prevent BVDotVec below to use B-product, restored at the end */
130:       BVSetMatrix(eps->V,NULL,PETSC_FALSE);

132:       /* Compute Rayleigh quotient */
133:       BVSetActiveColumns(eps->V,eps->nconv,nv);
134:       BVSetActiveColumns(ctx->AV,0,nv-eps->nconv);
135:       BVMatMult(eps->V,A,ctx->AV);
136:       DSGetArray(eps->ds,DS_MAT_A,&C);
137:       for (i=eps->nconv;i<nv;i++) {
138:         BVSetActiveColumns(eps->V,eps->nconv,i+1);
139:         BVGetColumn(ctx->AV,i-eps->nconv,&av);
140:         BVDotVec(eps->V,av,C+eps->nconv+i*ld);
141:         BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
142:         for (j=eps->nconv;j<i-1;j++) C[i+j*ld] = PetscConj(C[j+i*ld]);
143:       }
144:       DSRestoreArray(eps->ds,DS_MAT_A,&C);
145:       DSSetState(eps->ds,DS_STATE_RAW);

147:       /* Solve projected problem */
148:       DSSolve(eps->ds,eps->eigr,eps->eigi);
149:       DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
150:       DSSynchronize(eps->ds,eps->eigr,eps->eigi);

152:       /* Update vectors V(:,idx) = V * Y(:,idx) */
153:       DSGetMat(eps->ds,DS_MAT_Q,&Q);
154:       BVMultInPlace(eps->V,Q,eps->nconv,nv);
155:       ExtractSubmatrix(Q,eps->nconv,&Q1);
156:       BVMultInPlace(ctx->AV,Q1,0,nv-eps->nconv);
157:       MatDestroy(&Q);
158:       MatDestroy(&Q1);
159:       if (B) BVSetMatrix(eps->V,B,PETSC_FALSE);
160:     } else {
161:       /* No need to do Rayleigh-Ritz, just take diag(V'*A*V) */
162:       for (i=eps->nconv;i<nv;i++) {
163:         BVGetColumn(eps->V,i,&v);
164:         BVGetColumn(ctx->AV,i-eps->nconv,&av);
165:         MatMult(A,v,av);
166:         VecDot(av,v,eps->eigr+i);
167:         BVRestoreColumn(eps->V,i,&v);
168:         BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
169:       }
170:     }

172:     /* Compute gradient and check convergence */
173:     k = -1;
174:     for (i=eps->nconv;i<nv;i++) {
175:       BVGetColumn(eps->V,i,&v);
176:       BVGetColumn(ctx->AV,i-eps->nconv,&av);
177:       BVGetColumn(ctx->G,i-eps->nconv,&p);
178:       if (B) {
179:         BVGetColumn(ctx->W,i-eps->nconv,&bv);
180:         MatMult(B,v,bv);
181:         VecWAXPY(p,-eps->eigr[i],bv,av);
182:         BVRestoreColumn(ctx->W,i-eps->nconv,&bv);
183:       } else VecWAXPY(p,-eps->eigr[i],v,av);
184:       BVRestoreColumn(eps->V,i,&v);
185:       BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
186:       VecNorm(p,NORM_2,&resnorm);
187:       BVRestoreColumn(ctx->G,i-eps->nconv,&p);
188:       (*eps->converged)(eps,eps->eigr[i],0.0,resnorm,&eps->errest[i],eps->convergedctx);
189:       if (k==-1 && eps->errest[i] >= eps->tol) k = i;
190:     }
191:     if (k==-1) k = nv;
192:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);

194:     /* The next lines are necessary to avoid DS zeroing eigr */
195:     DSGetArray(eps->ds,DS_MAT_A,&C);
196:     for (i=eps->nconv;i<k;i++) C[i+i*ld] = eps->eigr[i];
197:     DSRestoreArray(eps->ds,DS_MAT_A,&C);

199:     if (eps->reason == EPS_CONVERGED_ITERATING) {

201:       /* Search direction */
202:       for (i=0;i<nv-eps->nconv;i++) {
203:         BVGetColumn(ctx->G,i,&v);
204:         STApply(eps->st,v,w);
205:         VecDot(w,v,&g);
206:         BVRestoreColumn(ctx->G,i,&v);
207:         beta = (!reset && eps->its>1)? g/gamma[i]: 0.0;
208:         gamma[i] = g;
209:         BVGetColumn(ctx->P,i,&v);
210:         VecAXPBY(v,1.0,beta,w);
211:         if (i+eps->nconv>0) {
212:           BVSetActiveColumns(eps->V,0,i+eps->nconv);
213:           BVOrthogonalizeVec(eps->V,v,NULL,NULL,NULL);
214:         }
215:         BVRestoreColumn(ctx->P,i,&v);
216:       }

218:       /* Minimization problem */
219:       for (i=eps->nconv;i<nv;i++) {
220:         BVGetColumn(eps->V,i,&v);
221:         BVGetColumn(ctx->AV,i-eps->nconv,&av);
222:         BVGetColumn(ctx->P,i-eps->nconv,&p);
223:         VecDot(av,v,&nu);
224:         VecDot(av,p,&pax);
225:         MatMult(A,p,w);
226:         VecDot(w,p,&pap);
227:         if (B) {
228:           BVGetColumn(ctx->W,i-eps->nconv,&bv);
229:           VecDot(bv,v,&mu);
230:           VecDot(bv,p,&pbx);
231:           BVRestoreColumn(ctx->W,i-eps->nconv,&bv);
232:           MatMult(B,p,w);
233:           VecDot(w,p,&pbp);
234:         } else {
235:           VecDot(v,v,&mu);
236:           VecDot(v,p,&pbx);
237:           VecDot(p,p,&pbp);
238:         }
239:         BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
240:         a = PetscRealPart(pap*pbx-pax*pbp);
241:         b = PetscRealPart(nu*pbp-mu*pap);
242:         c = PetscRealPart(mu*pax-nu*pbx);
243:         t = PetscMax(PetscMax(PetscAbsReal(a),PetscAbsReal(b)),PetscAbsReal(c));
244:         if (t!=0.0) { a /= t; b /= t; c /= t; }
245:         disc = b*b-4.0*a*c;
246:         d = PetscSqrtReal(PetscAbsReal(disc));
247:         if (b>=0.0 && a!=0.0) alpha = (b+d)/(2.0*a);
248:         else if (b!=d) alpha = 2.0*c/(b-d);
249:         else alpha = 0;
250:         /* Next iterate */
251:         if (alpha!=0.0) VecAXPY(v,alpha,p);
252:         BVRestoreColumn(eps->V,i,&v);
253:         BVRestoreColumn(ctx->P,i-eps->nconv,&p);
254:         BVOrthonormalizeColumn(eps->V,i,PETSC_TRUE,NULL,NULL);
255:       }
256:     }

258:     EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv);
259:     eps->nconv = k;
260:   }

262:   PetscFree(gamma);
263:   PetscFunctionReturn(0);
264: }

266: static PetscErrorCode EPSRQCGSetReset_RQCG(EPS eps,PetscInt nrest)
267: {
268:   EPS_RQCG *ctx = (EPS_RQCG*)eps->data;

270:   if (nrest==PETSC_DEFAULT) {
271:     ctx->nrest = 0;
272:     eps->state = EPS_STATE_INITIAL;
273:   } else {
275:     ctx->nrest = nrest;
276:   }
277:   PetscFunctionReturn(0);
278: }

280: /*@
281:    EPSRQCGSetReset - Sets the reset parameter of the RQCG iteration. Every
282:    nrest iterations, the solver performs a Rayleigh-Ritz projection step.

284:    Logically Collective on eps

286:    Input Parameters:
287: +  eps - the eigenproblem solver context
288: -  nrest - the number of iterations between resets

290:    Options Database Key:
291: .  -eps_rqcg_reset - Sets the reset parameter

293:    Level: advanced

295: .seealso: EPSRQCGGetReset()
296: @*/
297: PetscErrorCode EPSRQCGSetReset(EPS eps,PetscInt nrest)
298: {
301:   PetscTryMethod(eps,"EPSRQCGSetReset_C",(EPS,PetscInt),(eps,nrest));
302:   PetscFunctionReturn(0);
303: }

305: static PetscErrorCode EPSRQCGGetReset_RQCG(EPS eps,PetscInt *nrest)
306: {
307:   EPS_RQCG *ctx = (EPS_RQCG*)eps->data;

309:   *nrest = ctx->nrest;
310:   PetscFunctionReturn(0);
311: }

313: /*@
314:    EPSRQCGGetReset - Gets the reset parameter used in the RQCG method.

316:    Not Collective

318:    Input Parameter:
319: .  eps - the eigenproblem solver context

321:    Output Parameter:
322: .  nrest - the reset parameter

324:    Level: advanced

326: .seealso: EPSRQCGSetReset()
327: @*/
328: PetscErrorCode EPSRQCGGetReset(EPS eps,PetscInt *nrest)
329: {
332:   PetscUseMethod(eps,"EPSRQCGGetReset_C",(EPS,PetscInt*),(eps,nrest));
333:   PetscFunctionReturn(0);
334: }

336: PetscErrorCode EPSReset_RQCG(EPS eps)
337: {
338:   EPS_RQCG       *ctx = (EPS_RQCG*)eps->data;

340:   BVDestroy(&ctx->AV);
341:   BVDestroy(&ctx->W);
342:   BVDestroy(&ctx->P);
343:   BVDestroy(&ctx->G);
344:   ctx->allocsize = 0;
345:   PetscFunctionReturn(0);
346: }

348: PetscErrorCode EPSSetFromOptions_RQCG(PetscOptionItems *PetscOptionsObject,EPS eps)
349: {
350:   PetscBool      flg;
351:   PetscInt       nrest;

353:   PetscOptionsHead(PetscOptionsObject,"EPS RQCG Options");

355:     PetscOptionsInt("-eps_rqcg_reset","Reset parameter","EPSRQCGSetReset",20,&nrest,&flg);
356:     if (flg) EPSRQCGSetReset(eps,nrest);

358:   PetscOptionsTail();
359:   PetscFunctionReturn(0);
360: }

362: PetscErrorCode EPSDestroy_RQCG(EPS eps)
363: {
364:   PetscFree(eps->data);
365:   PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGSetReset_C",NULL);
366:   PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGGetReset_C",NULL);
367:   PetscFunctionReturn(0);
368: }

370: PetscErrorCode EPSView_RQCG(EPS eps,PetscViewer viewer)
371: {
372:   EPS_RQCG       *ctx = (EPS_RQCG*)eps->data;
373:   PetscBool      isascii;

375:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
376:   if (isascii) PetscViewerASCIIPrintf(viewer,"  reset every %" PetscInt_FMT " iterations\n",ctx->nrest);
377:   PetscFunctionReturn(0);
378: }

380: SLEPC_EXTERN PetscErrorCode EPSCreate_RQCG(EPS eps)
381: {
382:   EPS_RQCG       *rqcg;

384:   PetscNewLog(eps,&rqcg);
385:   eps->data = (void*)rqcg;

387:   eps->useds = PETSC_TRUE;
388:   eps->categ = EPS_CATEGORY_PRECOND;

390:   eps->ops->solve          = EPSSolve_RQCG;
391:   eps->ops->setup          = EPSSetUp_RQCG;
392:   eps->ops->setupsort      = EPSSetUpSort_Default;
393:   eps->ops->setfromoptions = EPSSetFromOptions_RQCG;
394:   eps->ops->destroy        = EPSDestroy_RQCG;
395:   eps->ops->reset          = EPSReset_RQCG;
396:   eps->ops->view           = EPSView_RQCG;
397:   eps->ops->backtransform  = EPSBackTransform_Default;
398:   eps->ops->setdefaultst   = EPSSetDefaultST_GMRES;

400:   PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGSetReset_C",EPSRQCGSetReset_RQCG);
401:   PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGGetReset_C",EPSRQCGGetReset_RQCG);
402:   PetscFunctionReturn(0);
403: }