Actual source code: rqcg.c

slepc-3.14.0 2020-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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: {
 91:   PetscInt       j,m,n;
 92:   PetscScalar    *pA,*pB;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

298:    Logically Collective on eps

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

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

307:    Level: advanced

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

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

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

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

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

334:    Not Collective

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

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

342:    Level: advanced

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

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

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

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

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

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

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

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

387: PetscErrorCode EPSDestroy_RQCG(EPS eps)
388: {

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

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

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

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

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

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

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

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