Actual source code: rqcg.c

slepc-3.16.2 2022-02-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: {
 40:   PetscInt       nmat;
 41:   EPS_RQCG       *ctx = (EPS_RQCG*)eps->data;

 44:   EPSCheckHermitianDefinite(eps);
 45:   EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
 46:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 47:   if (!eps->which) eps->which = EPS_SMALLEST_REAL;
 48:   if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(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:   EPSAllocateSolution(eps,0);
 55:   EPS_SetInnerProduct(eps);

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

 85: /*
 86:    ExtractSubmatrix - Returns B = A(k+1:end,k+1:end).
 87: */
 88: static PetscErrorCode ExtractSubmatrix(Mat A,PetscInt k,Mat *B)
 89: {
 90:   PetscErrorCode    ierr;
 91:   PetscInt          j,m,n;
 92:   PetscScalar       *pB;
 93:   const PetscScalar *pA;

 96:   MatGetSize(A,&m,&n);
 97:   MatCreateSeqDense(PETSC_COMM_SELF,m-k,n-k,NULL,B);
 98:   MatDenseGetArrayRead(A,&pA);
 99:   MatDenseGetArrayWrite(*B,&pB);
100:   for (j=k;j<n;j++) {
101:     PetscArraycpy(pB+(j-k)*(m-k),pA+j*m+k,m-k);
102:   }
103:   MatDenseRestoreArrayRead(A,&pA);
104:   MatDenseRestoreArrayWrite(*B,&pB);
105:   return(0);
106: }

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

120:   DSGetLeadingDimension(eps->ds,&ld);
121:   STGetNumMatrices(eps->st,&nmat);
122:   STGetMatrix(eps->st,0,&A);
123:   if (nmat>1) { STGetMatrix(eps->st,1,&B); }
124:   else B = NULL;
125:   PetscMalloc1(eps->mpd,&gamma);

127:   kini = eps->nini;
128:   while (eps->reason == EPS_CONVERGED_ITERATING) {
129:     eps->its++;
130:     nv = PetscMin(eps->nconv+eps->mpd,ncv);
131:     DSSetDimensions(eps->ds,nv,eps->nconv,0);
132:     for (;kini<nv;kini++) { /* Generate more initial vectors if necessary */
133:       BVSetRandomColumn(eps->V,kini);
134:       BVOrthonormalizeColumn(eps->V,kini,PETSC_TRUE,NULL,NULL);
135:     }
136:     reset = (eps->its>1 && (eps->its-1)%ctx->nrest==0)? PETSC_TRUE: PETSC_FALSE;

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

142:       /* Compute Rayleigh quotient */
143:       BVSetActiveColumns(eps->V,eps->nconv,nv);
144:       BVSetActiveColumns(ctx->AV,0,nv-eps->nconv);
145:       BVMatMult(eps->V,A,ctx->AV);
146:       DSGetArray(eps->ds,DS_MAT_A,&C);
147:       for (i=eps->nconv;i<nv;i++) {
148:         BVSetActiveColumns(eps->V,eps->nconv,i+1);
149:         BVGetColumn(ctx->AV,i-eps->nconv,&av);
150:         BVDotVec(eps->V,av,C+eps->nconv+i*ld);
151:         BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
152:         for (j=eps->nconv;j<i-1;j++) C[i+j*ld] = PetscConj(C[j+i*ld]);
153:       }
154:       DSRestoreArray(eps->ds,DS_MAT_A,&C);
155:       DSSetState(eps->ds,DS_STATE_RAW);

157:       /* Solve projected problem */
158:       DSSolve(eps->ds,eps->eigr,eps->eigi);
159:       DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
160:       DSSynchronize(eps->ds,eps->eigr,eps->eigi);

162:       /* Update vectors V(:,idx) = V * Y(:,idx) */
163:       DSGetMat(eps->ds,DS_MAT_Q,&Q);
164:       BVMultInPlace(eps->V,Q,eps->nconv,nv);
165:       ExtractSubmatrix(Q,eps->nconv,&Q1);
166:       BVMultInPlace(ctx->AV,Q1,0,nv-eps->nconv);
167:       MatDestroy(&Q);
168:       MatDestroy(&Q1);
169:       if (B) { BVSetMatrix(eps->V,B,PETSC_FALSE); }
170:     } else {
171:       /* No need to do Rayleigh-Ritz, just take diag(V'*A*V) */
172:       for (i=eps->nconv;i<nv;i++) {
173:         BVGetColumn(eps->V,i,&v);
174:         BVGetColumn(ctx->AV,i-eps->nconv,&av);
175:         MatMult(A,v,av);
176:         VecDot(av,v,eps->eigr+i);
177:         BVRestoreColumn(eps->V,i,&v);
178:         BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
179:       }
180:     }

182:     /* Compute gradient and check convergence */
183:     k = -1;
184:     for (i=eps->nconv;i<nv;i++) {
185:       BVGetColumn(eps->V,i,&v);
186:       BVGetColumn(ctx->AV,i-eps->nconv,&av);
187:       BVGetColumn(ctx->G,i-eps->nconv,&p);
188:       if (B) {
189:         BVGetColumn(ctx->W,i-eps->nconv,&bv);
190:         MatMult(B,v,bv);
191:         VecWAXPY(p,-eps->eigr[i],bv,av);
192:         BVRestoreColumn(ctx->W,i-eps->nconv,&bv);
193:       } else {
194:         VecWAXPY(p,-eps->eigr[i],v,av);
195:       }
196:       BVRestoreColumn(eps->V,i,&v);
197:       BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
198:       VecNorm(p,NORM_2,&resnorm);
199:       BVRestoreColumn(ctx->G,i-eps->nconv,&p);
200:       (*eps->converged)(eps,eps->eigr[i],0.0,resnorm,&eps->errest[i],eps->convergedctx);
201:       if (k==-1 && eps->errest[i] >= eps->tol) k = i;
202:     }
203:     if (k==-1) k = nv;
204:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);

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

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

213:       /* Search direction */
214:       for (i=0;i<nv-eps->nconv;i++) {
215:         BVGetColumn(ctx->G,i,&v);
216:         STApply(eps->st,v,w);
217:         VecDot(w,v,&g);
218:         BVRestoreColumn(ctx->G,i,&v);
219:         beta = (!reset && eps->its>1)? g/gamma[i]: 0.0;
220:         gamma[i] = g;
221:         BVGetColumn(ctx->P,i,&v);
222:         VecAXPBY(v,1.0,beta,w);
223:         if (i+eps->nconv>0) {
224:           BVSetActiveColumns(eps->V,0,i+eps->nconv);
225:           BVOrthogonalizeVec(eps->V,v,NULL,NULL,NULL);
226:         }
227:         BVRestoreColumn(ctx->P,i,&v);
228:       }

230:       /* Minimization problem */
231:       for (i=eps->nconv;i<nv;i++) {
232:         BVGetColumn(eps->V,i,&v);
233:         BVGetColumn(ctx->AV,i-eps->nconv,&av);
234:         BVGetColumn(ctx->P,i-eps->nconv,&p);
235:         VecDot(av,v,&nu);
236:         VecDot(av,p,&pax);
237:         MatMult(A,p,w);
238:         VecDot(w,p,&pap);
239:         if (B) {
240:           BVGetColumn(ctx->W,i-eps->nconv,&bv);
241:           VecDot(bv,v,&mu);
242:           VecDot(bv,p,&pbx);
243:           BVRestoreColumn(ctx->W,i-eps->nconv,&bv);
244:           MatMult(B,p,w);
245:           VecDot(w,p,&pbp);
246:         } else {
247:           VecDot(v,v,&mu);
248:           VecDot(v,p,&pbx);
249:           VecDot(p,p,&pbp);
250:         }
251:         BVRestoreColumn(ctx->AV,i-eps->nconv,&av);
252:         a = PetscRealPart(pap*pbx-pax*pbp);
253:         b = PetscRealPart(nu*pbp-mu*pap);
254:         c = PetscRealPart(mu*pax-nu*pbx);
255:         t = PetscMax(PetscMax(PetscAbsReal(a),PetscAbsReal(b)),PetscAbsReal(c));
256:         if (t!=0.0) { a /= t; b /= t; c /= t; }
257:         disc = b*b-4.0*a*c;
258:         d = PetscSqrtReal(PetscAbsReal(disc));
259:         if (b>=0.0 && a!=0.0) alpha = (b+d)/(2.0*a);
260:         else if (b!=d) alpha = 2.0*c/(b-d);
261:         else alpha = 0;
262:         /* Next iterate */
263:         if (alpha!=0.0) {
264:           VecAXPY(v,alpha,p);
265:         }
266:         BVRestoreColumn(eps->V,i,&v);
267:         BVRestoreColumn(ctx->P,i-eps->nconv,&p);
268:         BVOrthonormalizeColumn(eps->V,i,PETSC_TRUE,NULL,NULL);
269:       }
270:     }

272:     EPSMonitor(eps,eps->its,k,eps->eigr,eps->eigi,eps->errest,nv);
273:     eps->nconv = k;
274:   }

276:   PetscFree(gamma);
277:   return(0);
278: }

280: static PetscErrorCode EPSRQCGSetReset_RQCG(EPS eps,PetscInt nrest)
281: {
282:   EPS_RQCG *ctx = (EPS_RQCG*)eps->data;

285:   if (nrest==PETSC_DEFAULT) {
286:     ctx->nrest = 0;
287:     eps->state = EPS_STATE_INITIAL;
288:   } else {
289:     if (nrest<=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reset parameter must be >0");
290:     ctx->nrest = nrest;
291:   }
292:   return(0);
293: }

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

299:    Logically Collective on eps

301:    Input Parameters:
302: +  eps - the eigenproblem solver context
303: -  nrest - the number of iterations between resets

305:    Options Database Key:
306: .  -eps_rqcg_reset - Sets the reset parameter

308:    Level: advanced

310: .seealso: EPSRQCGGetReset()
311: @*/
312: PetscErrorCode EPSRQCGSetReset(EPS eps,PetscInt nrest)
313: {

319:   PetscTryMethod(eps,"EPSRQCGSetReset_C",(EPS,PetscInt),(eps,nrest));
320:   return(0);
321: }

323: static PetscErrorCode EPSRQCGGetReset_RQCG(EPS eps,PetscInt *nrest)
324: {
325:   EPS_RQCG *ctx = (EPS_RQCG*)eps->data;

328:   *nrest = ctx->nrest;
329:   return(0);
330: }

332: /*@
333:    EPSRQCGGetReset - Gets the reset parameter used in the RQCG method.

335:    Not Collective

337:    Input Parameter:
338: .  eps - the eigenproblem solver context

340:    Output Parameter:
341: .  nrest - the reset parameter

343:    Level: advanced

345: .seealso: EPSRQCGSetReset()
346: @*/
347: PetscErrorCode EPSRQCGGetReset(EPS eps,PetscInt *nrest)
348: {

354:   PetscUseMethod(eps,"EPSRQCGGetReset_C",(EPS,PetscInt*),(eps,nrest));
355:   return(0);
356: }

358: PetscErrorCode EPSReset_RQCG(EPS eps)
359: {
361:   EPS_RQCG       *ctx = (EPS_RQCG*)eps->data;

364:   BVDestroy(&ctx->AV);
365:   BVDestroy(&ctx->W);
366:   BVDestroy(&ctx->P);
367:   BVDestroy(&ctx->G);
368:   ctx->allocsize = 0;
369:   return(0);
370: }

372: PetscErrorCode EPSSetFromOptions_RQCG(PetscOptionItems *PetscOptionsObject,EPS eps)
373: {
375:   PetscBool      flg;
376:   PetscInt       nrest;

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

381:     PetscOptionsInt("-eps_rqcg_reset","Reset parameter","EPSRQCGSetReset",20,&nrest,&flg);
382:     if (flg) { EPSRQCGSetReset(eps,nrest); }

384:   PetscOptionsTail();
385:   return(0);
386: }

388: PetscErrorCode EPSDestroy_RQCG(EPS eps)
389: {

393:   PetscFree(eps->data);
394:   PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGSetReset_C",NULL);
395:   PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGGetReset_C",NULL);
396:   return(0);
397: }

399: PetscErrorCode EPSView_RQCG(EPS eps,PetscViewer viewer)
400: {
402:   EPS_RQCG       *ctx = (EPS_RQCG*)eps->data;
403:   PetscBool      isascii;

406:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
407:   if (isascii) {
408:     PetscViewerASCIIPrintf(viewer,"  reset every %D iterations\n",ctx->nrest);
409:   }
410:   return(0);
411: }

413: SLEPC_EXTERN PetscErrorCode EPSCreate_RQCG(EPS eps)
414: {
415:   EPS_RQCG       *rqcg;

419:   PetscNewLog(eps,&rqcg);
420:   eps->data = (void*)rqcg;

422:   eps->useds = PETSC_TRUE;
423:   eps->categ = EPS_CATEGORY_PRECOND;

425:   eps->ops->solve          = EPSSolve_RQCG;
426:   eps->ops->setup          = EPSSetUp_RQCG;
427:   eps->ops->setupsort      = EPSSetUpSort_Default;
428:   eps->ops->setfromoptions = EPSSetFromOptions_RQCG;
429:   eps->ops->destroy        = EPSDestroy_RQCG;
430:   eps->ops->reset          = EPSReset_RQCG;
431:   eps->ops->view           = EPSView_RQCG;
432:   eps->ops->backtransform  = EPSBackTransform_Default;
433:   eps->ops->setdefaultst   = EPSSetDefaultST_GMRES;

435:   PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGSetReset_C",EPSRQCGSetReset_RQCG);
436:   PetscObjectComposeFunction((PetscObject)eps,"EPSRQCGGetReset_C",EPSRQCGGetReset_RQCG);
437:   return(0);
438: }