Actual source code: qarnoldi.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*

  3:    SLEPc quadratic eigensolver: "qarnoldi"

  5:    Method: Q-Arnoldi

  7:    Algorithm:

  9:        Quadratic Arnoldi with Krylov-Schur type restart.

 11:    References:

 13:        [1] K. Meerbergen, "The Quadratic Arnoldi method for the solution
 14:            of the quadratic eigenvalue problem", SIAM J. Matrix Anal.
 15:            Appl. 30(4):1462-1482, 2008.

 17:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 18:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 19:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

 21:    This file is part of SLEPc.

 23:    SLEPc is free software: you can redistribute it and/or modify it under  the
 24:    terms of version 3 of the GNU Lesser General Public License as published by
 25:    the Free Software Foundation.

 27:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 28:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 29:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 30:    more details.

 32:    You  should have received a copy of the GNU Lesser General  Public  License
 33:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 34:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 35: */

 37: #include <slepc-private/pepimpl.h>    /*I "slepcpep.h" I*/
 38: #include <petscblaslapack.h>

 40: typedef struct {
 41:   PetscReal keep;         /* restart parameter */
 42: } PEP_QARNOLDI;

 46: PetscErrorCode PEPSetUp_QArnoldi(PEP pep)
 47: {
 49:   PEP_QARNOLDI   *ctx = (PEP_QARNOLDI*)pep->data;
 50:   PetscBool      sinv,flg;

 53:   PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
 54:   if (!pep->max_it) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
 55:   if (!pep->which) {
 56:     PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
 57:     if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
 58:     else pep->which = PEP_LARGEST_MAGNITUDE;
 59:   }

 61:   if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver only available for quadratic problems");
 62:   if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
 63:   STGetTransform(pep->st,&flg);
 64:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag set, see STSetTransform()");

 66:   if (!ctx->keep) ctx->keep = 0.5;

 68:   PEPAllocateSolution(pep,0);
 69:   PEPSetWorkVecs(pep,4);

 71:   DSSetType(pep->ds,DSNHEP);
 72:   DSSetExtraRow(pep->ds,PETSC_TRUE);
 73:   DSAllocate(pep->ds,pep->ncv+1);
 74:   return(0);
 75: }

 79: /*
 80:   Compute a step of Classical Gram-Schmidt orthogonalization
 81: */
 82: static PetscErrorCode PEPQArnoldiCGS(PEP pep,PetscScalar *H,PetscBLASInt ldh,PetscScalar *h,PetscBLASInt j,BV V,Vec t,Vec v,Vec w,PetscReal *onorm,PetscReal *norm,PetscScalar *work)
 83: {
 85:   PetscBLASInt   ione = 1,j_1 = j+1;
 86:   PetscReal      x,y;
 87:   PetscScalar    dot,one = 1.0,zero = 0.0;

 90:   /* compute norm of v and w */
 91:   if (onorm) {
 92:     VecNorm(v,NORM_2,&x);
 93:     VecNorm(w,NORM_2,&y);
 94:     *onorm = PetscSqrtReal(x*x+y*y);
 95:   }

 97:   /* orthogonalize: compute h */
 98:   BVDotVec(V,v,h);
 99:   BVDotVec(V,w,work);
100:   if (j>0)
101:     PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&j_1,&j,&one,H,&ldh,work,&ione,&one,h,&ione));
102:   VecDot(w,t,&dot);
103:   h[j] += dot;

105:   /* orthogonalize: update v and w */
106:   BVMultVec(V,-1.0,1.0,v,h);
107:   if (j>0) {
108:     PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&j_1,&j,&one,H,&ldh,h,&ione,&zero,work,&ione));
109:     BVMultVec(V,-1.0,1.0,w,work);
110:   }
111:   VecAXPY(w,-h[j],t);

113:   /* compute norm of v and w */
114:   if (norm) {
115:     VecNorm(v,NORM_2,&x);
116:     VecNorm(w,NORM_2,&y);
117:     *norm = PetscSqrtReal(x*x+y*y);
118:   }
119:   return(0);
120: }

124: /*
125:   Compute a run of Q-Arnoldi iterations
126: */
127: static PetscErrorCode PEPQArnoldi(PEP pep,PetscScalar *H,PetscInt ldh,PetscInt k,PetscInt *M,Vec v,Vec w,PetscReal *beta,PetscBool *breakdown,PetscScalar *work)
128: {
129:   PetscErrorCode     ierr;
130:   PetscInt           i,j,l,m = *M;
131:   Vec                t = pep->work[2],u = pep->work[3];
132:   BVOrthogRefineType refinement;
133:   PetscReal          norm,onorm,eta;
134:   PetscScalar        *c = work + m;

137:   BVGetOrthogonalization(pep->V,NULL,&refinement,&eta);
138:   BVInsertVec(pep->V,k,v);
139:   for (j=k;j<m;j++) {
140:     /* apply operator */
141:     VecCopy(w,t);
142:     if (pep->Dr) {
143:       VecPointwiseMult(v,v,pep->Dr);
144:     }
145:     STMatMult(pep->st,0,v,u);
146:     VecCopy(t,v);
147:     if (pep->Dr) {
148:       VecPointwiseMult(t,t,pep->Dr);
149:     }
150:     STMatMult(pep->st,1,t,w);
151:     VecAXPY(u,pep->sfactor,w);
152:     STMatSolve(pep->st,u,w);
153:     VecScale(w,-1.0/(pep->sfactor*pep->sfactor));
154:     if (pep->Dr) {
155:       VecPointwiseDivide(w,w,pep->Dr);
156:     }
157:     VecCopy(v,t);
158:     BVSetActiveColumns(pep->V,0,j+1);

160:     /* orthogonalize */
161:     switch (refinement) {
162:       case BV_ORTHOG_REFINE_NEVER:
163:         PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,&norm,work);
164:         *breakdown = PETSC_FALSE;
165:         break;
166:       case BV_ORTHOG_REFINE_ALWAYS:
167:         PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,NULL,NULL,work);
168:         PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,&onorm,&norm,work);
169:         for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
170:         if (norm < eta * onorm) *breakdown = PETSC_TRUE;
171:         else *breakdown = PETSC_FALSE;
172:         break;
173:       case BV_ORTHOG_REFINE_IFNEEDED:
174:         PEPQArnoldiCGS(pep,H,ldh,H+ldh*j,j,pep->V,t,v,w,&onorm,&norm,work);
175:         /* ||q|| < eta ||h|| */
176:         l = 1;
177:         while (l<3 && norm < eta * onorm) {
178:           l++;
179:           onorm = norm;
180:           PEPQArnoldiCGS(pep,H,ldh,c,j,pep->V,t,v,w,NULL,&norm,work);
181:           for (i=0;i<=j;i++) H[ldh*j+i] += c[i];
182:         }
183:         if (norm < eta * onorm) *breakdown = PETSC_TRUE;
184:         else *breakdown = PETSC_FALSE;
185:         break;
186:       default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of ip->orth_ref");
187:     }
188:     VecScale(v,1.0/norm);
189:     VecScale(w,1.0/norm);

191:     H[j+1+ldh*j] = norm;
192:     if (j<m-1) {
193:       BVInsertVec(pep->V,j+1,v);
194:     }
195:   }
196:   *beta = norm;
197:   return(0);
198: }

202: PetscErrorCode PEPSolve_QArnoldi(PEP pep)
203: {
205:   PEP_QARNOLDI   *ctx = (PEP_QARNOLDI*)pep->data;
206:   PetscInt       j,k,l,lwork,nv,ld,newn;
207:   Vec            v=pep->work[0],w=pep->work[1];
208:   Mat            Q;
209:   PetscScalar    *S,*work;
210:   PetscReal      beta=0.0,norm,x,y;
211:   PetscBool      breakdown=PETSC_FALSE;

214:   DSGetLeadingDimension(pep->ds,&ld);
215:   lwork = 7*pep->ncv;
216:   PetscMalloc1(lwork,&work);

218:   /* Get the starting Arnoldi vector */
219:   if (pep->nini==0) {
220:     BVSetRandomColumn(pep->V,0,pep->rand);
221:   }
222:   /* w is always a random vector */
223:   BVSetRandomColumn(pep->V,1,pep->rand);
224:   BVCopyVec(pep->V,0,v);
225:   BVCopyVec(pep->V,1,w);
226:   VecNorm(v,NORM_2,&x);
227:   VecNorm(w,NORM_2,&y);
228:   norm = PetscSqrtReal(x*x+y*y);
229:   VecScale(v,1.0/norm);
230:   VecScale(w,1.0/norm);

232:    /* Restart loop */
233:   l = 0;
234:   while (pep->reason == PEP_CONVERGED_ITERATING) {
235:     pep->its++;

237:     /* Compute an nv-step Arnoldi factorization */
238:     nv = PetscMin(pep->nconv+pep->mpd,pep->ncv);
239:     DSGetArray(pep->ds,DS_MAT_A,&S);
240:     PEPQArnoldi(pep,S,ld,pep->nconv+l,&nv,v,w,&beta,&breakdown,work);
241:     DSRestoreArray(pep->ds,DS_MAT_A,&S);
242:     DSSetDimensions(pep->ds,nv,0,pep->nconv,pep->nconv+l);
243:     if (l==0) {
244:       DSSetState(pep->ds,DS_STATE_INTERMEDIATE);
245:     } else {
246:       DSSetState(pep->ds,DS_STATE_RAW);
247:     }
248:     BVSetActiveColumns(pep->V,pep->nconv,nv);

250:     /* Solve projected problem */
251:     DSSolve(pep->ds,pep->eigr,pep->eigi);
252:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
253:     DSUpdateExtraRow(pep->ds);

255:     /* Check convergence */
256:     PEPKrylovConvergence(pep,PETSC_FALSE,pep->nconv,nv-pep->nconv,beta,&k);
257:     if (pep->its >= pep->max_it) pep->reason = PEP_DIVERGED_ITS;
258:     if (k >= pep->nev) pep->reason = PEP_CONVERGED_TOL;

260:     /* Update l */
261:     if (pep->reason != PEP_CONVERGED_ITERATING || breakdown) l = 0;
262:     else l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));

264:     if (pep->reason == PEP_CONVERGED_ITERATING) {
265:       if (breakdown) {
266:         /* Stop if breakdown */
267:         PetscInfo2(pep,"Breakdown Quadratic Arnoldi method (it=%D norm=%g)\n",pep->its,(double)beta);
268:         pep->reason = PEP_DIVERGED_BREAKDOWN;
269:       } else {
270:         /* Prepare the Rayleigh quotient for restart */
271:         DSTruncate(pep->ds,k+l);
272:         DSGetDimensions(pep->ds,&newn,NULL,NULL,NULL,NULL);
273:         l = newn-k;
274:       }
275:     }
276:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
277:     DSGetMat(pep->ds,DS_MAT_Q,&Q);
278:     BVMultInPlace(pep->V,Q,pep->nconv,k+l);
279:     MatDestroy(&Q);

281:     pep->nconv = k;
282:     PEPMonitor(pep,pep->its,pep->nconv,pep->eigr,pep->eigi,pep->errest,nv);
283:   }

285:   for (j=0;j<pep->nconv;j++) {
286:     pep->eigr[j] *= pep->sfactor;
287:     pep->eigi[j] *= pep->sfactor;
288:   }

290:   /* truncate Schur decomposition and change the state to raw so that
291:      DSVectors() computes eigenvectors from scratch */
292:   DSSetDimensions(pep->ds,pep->nconv,0,0,0);
293:   DSSetState(pep->ds,DS_STATE_RAW);
294:   PetscFree(work);
295:   return(0);
296: }

300: static PetscErrorCode PEPQArnoldiSetRestart_QArnoldi(PEP pep,PetscReal keep)
301: {
302:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

305:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
306:   else {
307:     if (keep<0.1 || keep>0.9) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
308:     ctx->keep = keep;
309:   }
310:   return(0);
311: }

315: /*@
316:    PEPQArnoldiSetRestart - Sets the restart parameter for the Q-Arnoldi
317:    method, in particular the proportion of basis vectors that must be kept
318:    after restart.

320:    Logically Collective on PEP

322:    Input Parameters:
323: +  pep  - the eigenproblem solver context
324: -  keep - the number of vectors to be kept at restart

326:    Options Database Key:
327: .  -pep_qarnoldi_restart - Sets the restart parameter

329:    Notes:
330:    Allowed values are in the range [0.1,0.9]. The default is 0.5.

332:    Level: advanced

334: .seealso: PEPQArnoldiGetRestart()
335: @*/
336: PetscErrorCode PEPQArnoldiSetRestart(PEP pep,PetscReal keep)
337: {

343:   PetscTryMethod(pep,"PEPQArnoldiSetRestart_C",(PEP,PetscReal),(pep,keep));
344:   return(0);
345: }

349: static PetscErrorCode PEPQArnoldiGetRestart_QArnoldi(PEP pep,PetscReal *keep)
350: {
351:   PEP_QARNOLDI *ctx = (PEP_QARNOLDI*)pep->data;

354:   *keep = ctx->keep;
355:   return(0);
356: }

360: /*@
361:    PEPQArnoldiGetRestart - Gets the restart parameter used in the Q-Arnoldi method.

363:    Not Collective

365:    Input Parameter:
366: .  pep - the eigenproblem solver context

368:    Output Parameter:
369: .  keep - the restart parameter

371:    Level: advanced

373: .seealso: PEPQArnoldiSetRestart()
374: @*/
375: PetscErrorCode PEPQArnoldiGetRestart(PEP pep,PetscReal *keep)
376: {

382:   PetscTryMethod(pep,"PEPQArnoldiGetRestart_C",(PEP,PetscReal*),(pep,keep));
383:   return(0);
384: }

388: PetscErrorCode PEPSetFromOptions_QArnoldi(PEP pep)
389: {
391:   PetscBool      flg;
392:   PetscReal      keep;

395:   PetscOptionsHead("PEP Q-Arnoldi Options");
396:   PetscOptionsReal("-pep_qarnoldi_restart","Proportion of vectors kept after restart","PEPQArnoldiSetRestart",0.5,&keep,&flg);
397:   if (flg) {
398:     PEPQArnoldiSetRestart(pep,keep);
399:   }
400:   PetscOptionsTail();
401:   return(0);
402: }

406: PetscErrorCode PEPView_QArnoldi(PEP pep,PetscViewer viewer)
407: {
409:   PEP_QARNOLDI   *ctx = (PEP_QARNOLDI*)pep->data;
410:   PetscBool      isascii;

413:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
414:   if (isascii) {
415:     PetscViewerASCIIPrintf(viewer,"  Q-Arnoldi: %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
416:   }
417:   return(0);
418: }

422: PetscErrorCode PEPDestroy_QArnoldi(PEP pep)
423: {

427:   PetscFree(pep->data);
428:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",NULL);
429:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",NULL);
430:   return(0);
431: }

435: PETSC_EXTERN PetscErrorCode PEPCreate_QArnoldi(PEP pep)
436: {
437:   PEP_QARNOLDI   *ctx;

441:   PetscNewLog(pep,&ctx);
442:   pep->data = (void*)ctx;

444:   pep->ops->solve          = PEPSolve_QArnoldi;
445:   pep->ops->setup          = PEPSetUp_QArnoldi;
446:   pep->ops->setfromoptions = PEPSetFromOptions_QArnoldi;
447:   pep->ops->destroy        = PEPDestroy_QArnoldi;
448:   pep->ops->view           = PEPView_QArnoldi;
449:   pep->ops->computevectors = PEPComputeVectors_Schur;
450:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiSetRestart_C",PEPQArnoldiSetRestart_QArnoldi);
451:   PetscObjectComposeFunction((PetscObject)pep,"PEPQArnoldiGetRestart_C",PEPQArnoldiGetRestart_QArnoldi);
452:   return(0);
453: }