Actual source code: pepdefault.c

slepc-3.7.4 2017-05-17
Report Typos and Errors
  1: /*
  2:      This file contains some simple default routines for common PEP operations.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

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

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

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

 24: #include <slepc/private/pepimpl.h>     /*I "slepcpep.h" I*/

 28: /*@
 29:    PEPSetWorkVecs - Sets a number of work vectors into a PEP object.

 31:    Collective on PEP

 33:    Input Parameters:
 34: +  pep - polynomial eigensolver context
 35: -  nw  - number of work vectors to allocate

 37:    Developers Note:
 38:    This is PETSC_EXTERN because it may be required by user plugin PEP
 39:    implementations.

 41:    Level: developer
 42: @*/
 43: PetscErrorCode PEPSetWorkVecs(PEP pep,PetscInt nw)
 44: {
 46:   Vec            t;

 49:   if (pep->nwork < nw) {
 50:     VecDestroyVecs(pep->nwork,&pep->work);
 51:     pep->nwork = nw;
 52:     BVGetColumn(pep->V,0,&t);
 53:     VecDuplicateVecs(t,nw,&pep->work);
 54:     BVRestoreColumn(pep->V,0,&t);
 55:     PetscLogObjectParents(pep,nw,pep->work);
 56:   }
 57:   return(0);
 58: }

 62: /*
 63:   PEPConvergedRelative - Checks convergence relative to the eigenvalue.
 64: */
 65: PetscErrorCode PEPConvergedRelative(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 66: {
 67:   PetscReal w;

 70:   w = SlepcAbsEigenvalue(eigr,eigi);
 71:   *errest = res/w;
 72:   return(0);
 73: }

 77: /*
 78:   PEPConvergedNorm - Checks convergence relative to the matrix norms.
 79: */
 80: PetscErrorCode PEPConvergedNorm(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 81: {
 82:   PetscReal      w=0.0,t;
 83:   PetscInt       j;
 84:   PetscBool      flg;

 88:   /* initialization of matrix norms */
 89:   if (!pep->nrma[pep->nmat-1]) {
 90:     for (j=0;j<pep->nmat;j++) {
 91:       MatHasOperation(pep->A[j],MATOP_NORM,&flg);
 92:       if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"The convergence test related to the matrix norms requires a matrix norm operation");
 93:       MatNorm(pep->A[j],NORM_INFINITY,&pep->nrma[j]);
 94:     }
 95:   }
 96:   t = SlepcAbsEigenvalue(eigr,eigi);
 97:   for (j=pep->nmat-1;j>=0;j--) {
 98:     w = w*t+pep->nrma[j];
 99:   }
100:   *errest = res/w;
101:   return(0);
102: }

106: /*
107:   PEPConvergedAbsolute - Checks convergence absolutely.
108: */
109: PetscErrorCode PEPConvergedAbsolute(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
110: {
112:   *errest = res;
113:   return(0);
114: }

118: /*@C
119:    PEPStoppingBasic - Default routine to determine whether the outer eigensolver
120:    iteration must be stopped.

122:    Collective on PEP

124:    Input Parameters:
125: +  pep    - eigensolver context obtained from PEPCreate()
126: .  its    - current number of iterations
127: .  max_it - maximum number of iterations
128: .  nconv  - number of currently converged eigenpairs
129: .  nev    - number of requested eigenpairs
130: -  ctx    - context (not used here)

132:    Output Parameter:
133: .  reason - result of the stopping test

135:    Notes:
136:    A positive value of reason indicates that the iteration has finished successfully
137:    (converged), and a negative value indicates an error condition (diverged). If
138:    the iteration needs to be continued, reason must be set to PEP_CONVERGED_ITERATING
139:    (zero).

141:    PEPStoppingBasic() will stop if all requested eigenvalues are converged, or if
142:    the maximum number of iterations has been reached.

144:    Use PEPSetStoppingTest() to provide your own test instead of using this one.

146:    Level: advanced

148: .seealso: PEPSetStoppingTest(), PEPConvergedReason, PEPGetConvergedReason()
149: @*/
150: PetscErrorCode PEPStoppingBasic(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ctx)
151: {

155:   *reason = PEP_CONVERGED_ITERATING;
156:   if (nconv >= nev) {
157:     PetscInfo2(pep,"Polynomial eigensolver finished successfully: %D eigenpairs converged at iteration %D\n",nconv,its);
158:     *reason = PEP_CONVERGED_TOL;
159:   } else if (its >= max_it) {
160:     *reason = PEP_DIVERGED_ITS;
161:     PetscInfo1(pep,"Polynomial eigensolver iteration reached maximum number of iterations (%D)\n",its);
162:   }
163:   return(0);
164: }

168: PetscErrorCode PEPBackTransform_Default(PEP pep)
169: {

173:   STBackTransform(pep->st,pep->nconv,pep->eigr,pep->eigi);
174:   return(0);
175: }

179: PetscErrorCode PEPComputeVectors_Default(PEP pep)
180: {
182:   PetscInt       i;
183:   Vec            v;
184: #if !defined(PETSC_USE_COMPLEX)
185:   Vec            v1;
186: #endif

189:   PEPExtractVectors(pep);

191:   /* Fix eigenvectors if balancing was used */
192:   if ((pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) && pep->Dr && (pep->refine!=PEP_REFINE_MULTIPLE)) {
193:     for (i=0;i<pep->nconv;i++) {
194:       BVGetColumn(pep->V,i,&v);
195:       VecPointwiseMult(v,v,pep->Dr);
196:       BVRestoreColumn(pep->V,i,&v);
197:     }
198:   }

200:   /* normalization */
201:   for (i=0;i<pep->nconv;i++) {
202: #if !defined(PETSC_USE_COMPLEX)
203:     if (pep->eigi[i]!=0.0) {   /* first eigenvalue of a complex conjugate pair */
204:       BVGetColumn(pep->V,i,&v);
205:       BVGetColumn(pep->V,i+1,&v1);
206:       SlepcVecNormalize(v,v1,PETSC_TRUE,NULL);
207:       BVRestoreColumn(pep->V,i,&v);
208:       BVRestoreColumn(pep->V,i+1,&v1);
209:       i++;
210:     } else   /* real eigenvalue */
211: #endif
212:     {
213:       BVGetColumn(pep->V,i,&v);
214:       SlepcVecNormalize(v,NULL,PETSC_FALSE,NULL);
215:       BVRestoreColumn(pep->V,i,&v);
216:     }
217:   }
218:   return(0);
219: }

223: /*
224:    PEPKrylovConvergence - This is the analogue to EPSKrylovConvergence, but
225:    for polynomial Krylov methods.

227:    Differences:
228:    - Always non-symmetric
229:    - Does not check for STSHIFT
230:    - No correction factor
231:    - No support for true residual
232: */
233: PetscErrorCode PEPKrylovConvergence(PEP pep,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal beta,PetscInt *kout)
234: {
236:   PetscInt       k,newk,marker,inside;
237:   PetscScalar    re,im;
238:   PetscReal      resnorm;
239:   PetscBool      istrivial;

242:   RGIsTrivial(pep->rg,&istrivial);
243:   marker = -1;
244:   if (pep->trackall) getall = PETSC_TRUE;
245:   for (k=kini;k<kini+nits;k++) {
246:     /* eigenvalue */
247:     re = pep->eigr[k];
248:     im = pep->eigi[k];
249:     if (!istrivial) {
250:       STBackTransform(pep->st,1,&re,&im);
251:       RGCheckInside(pep->rg,1,&re,&im,&inside);
252:       if (marker==-1 && inside<0) marker = k;
253:       re = pep->eigr[k];
254:       im = pep->eigi[k];
255:     }
256:     newk = k;
257:     DSVectors(pep->ds,DS_MAT_X,&newk,&resnorm);
258:     resnorm *= beta;
259:     /* error estimate */
260:     (*pep->converged)(pep,re,im,resnorm,&pep->errest[k],pep->convergedctx);
261:     if (marker==-1 && pep->errest[k] >= pep->tol) marker = k;
262:     if (newk==k+1) {
263:       pep->errest[k+1] = pep->errest[k];
264:       k++;
265:     }
266:     if (marker!=-1 && !getall) break;
267:   }
268:   if (marker!=-1) k = marker;
269:   *kout = k;
270:   return(0);
271: }

275: /*
276:   PEPBuildDiagonalScaling - compute two diagonal matrices to be applied for balancing 
277:   in polynomial eigenproblems.
278: */
279: PetscErrorCode PEPBuildDiagonalScaling(PEP pep)
280: {
282:   PetscInt       it,i,j,k,nmat,nr,e,nz,lst,lend,nc=0,*cols,emax,emin,emaxl,eminl;
283:   const PetscInt *cidx,*ridx;
284:   Mat            M,*T,A;
285:   PetscMPIInt    n;
286:   PetscBool      cont=PETSC_TRUE,flg=PETSC_FALSE;
287:   PetscScalar    *array,*Dr,*Dl,t;
288:   PetscReal      l2,d,*rsum,*aux,*csum,w=1.0;
289:   MatStructure   str;
290:   MatInfo        info;

293:   l2 = 2*PetscLogReal(2.0);
294:   nmat = pep->nmat;
295:   PetscMPIIntCast(pep->n,&n);
296:   STGetMatStructure(pep->st,&str);
297:   PetscMalloc1(nmat,&T);
298:   for (k=0;k<nmat;k++) {
299:     STGetTOperators(pep->st,k,&T[k]);
300:   }
301:   /* Form local auxiliar matrix M */
302:   PetscObjectTypeCompareAny((PetscObject)T[0],&cont,MATMPIAIJ,MATSEQAIJ,"");
303:   if (!cont) SETERRQ(PetscObjectComm((PetscObject)T[0]),PETSC_ERR_SUP,"Only for MPIAIJ or SEQAIJ matrix types");
304:   PetscObjectTypeCompare((PetscObject)T[0],MATMPIAIJ,&cont);
305:   if (cont) {
306:     MatMPIAIJGetLocalMat(T[0],MAT_INITIAL_MATRIX,&M);
307:     flg = PETSC_TRUE; 
308:   } else {
309:     MatDuplicate(T[0],MAT_COPY_VALUES,&M);
310:   }
311:   MatGetInfo(M,MAT_LOCAL,&info);
312:   nz = (PetscInt)info.nz_used;
313:   MatSeqAIJGetArray(M,&array);
314:   for (i=0;i<nz;i++) {
315:     t = PetscAbsScalar(array[i]);
316:     array[i] = t*t;
317:   }
318:   MatSeqAIJRestoreArray(M,&array);
319:   for (k=1;k<nmat;k++) {
320:     if (flg) {
321:       MatMPIAIJGetLocalMat(T[k],MAT_INITIAL_MATRIX,&A);
322:     } else {
323:       if (str==SAME_NONZERO_PATTERN) {
324:         MatCopy(T[k],A,SAME_NONZERO_PATTERN);
325:       } else {
326:         MatDuplicate(T[k],MAT_COPY_VALUES,&A);
327:       }
328:     }
329:     MatGetInfo(A,MAT_LOCAL,&info);
330:     nz = (PetscInt)info.nz_used;
331:     MatSeqAIJGetArray(A,&array);
332:     for (i=0;i<nz;i++) {
333:       t = PetscAbsScalar(array[i]);
334:       array[i] = t*t;
335:     }
336:     MatSeqAIJRestoreArray(A,&array);
337:     w *= pep->slambda*pep->slambda*pep->sfactor;
338:     MatAXPY(M,w,A,str);
339:     if (flg || str!=SAME_NONZERO_PATTERN || k==nmat-2) {
340:       MatDestroy(&A);
341:     } 
342:   }
343:   MatGetRowIJ(M,0,PETSC_FALSE,PETSC_FALSE,&nr,&ridx,&cidx,&cont);
344:   if (!cont) SETERRQ(PetscObjectComm((PetscObject)T[0]), PETSC_ERR_SUP,"It is not possible to compute scaling diagonals for these PEP matrices");
345:   MatGetInfo(M,MAT_LOCAL,&info);
346:   nz = (PetscInt)info.nz_used;
347:   VecGetOwnershipRange(pep->Dl,&lst,&lend);
348:   PetscMalloc4(nr,&rsum,pep->n,&csum,pep->n,&aux,PetscMin(pep->n-lend+lst,nz),&cols);
349:   VecSet(pep->Dr,1.0);
350:   VecSet(pep->Dl,1.0);
351:   VecGetArray(pep->Dl,&Dl);
352:   VecGetArray(pep->Dr,&Dr);
353:   MatSeqAIJGetArray(M,&array);
354:   PetscMemzero(aux,pep->n*sizeof(PetscReal));
355:   for (j=0;j<nz;j++) {
356:     /* Search non-zero columns outsize lst-lend */
357:     if (aux[cidx[j]]==0 && (cidx[j]<lst || lend<=cidx[j])) cols[nc++] = cidx[j];
358:     /* Local column sums */
359:     aux[cidx[j]] += PetscAbsScalar(array[j]);
360:   }
361:   for (it=0;it<pep->sits && cont;it++) {
362:     emaxl = 0; eminl = 0;
363:     /* Column sum  */    
364:     if (it>0) { /* it=0 has been already done*/
365:       MatSeqAIJGetArray(M,&array);
366:       PetscMemzero(aux,pep->n*sizeof(PetscReal));
367:       for (j=0;j<nz;j++) aux[cidx[j]] += PetscAbsScalar(array[j]);
368:       MatSeqAIJRestoreArray(M,&array); 
369:     }
370:     MPI_Allreduce(aux,csum,n,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)pep->Dr));
371:     /* Update Dr */
372:     for (j=lst;j<lend;j++) {
373:       d = PetscLogReal(csum[j])/l2;
374:       e = -(PetscInt)((d < 0)?(d-0.5):(d+0.5));
375:       d = PetscPowReal(2.0,e);
376:       Dr[j-lst] *= d;
377:       aux[j] = d*d;
378:       emaxl = PetscMax(emaxl,e);
379:       eminl = PetscMin(eminl,e);
380:     }
381:     for (j=0;j<nc;j++) {
382:       d = PetscLogReal(csum[cols[j]])/l2;
383:       e = -(PetscInt)((d < 0)?(d-0.5):(d+0.5));
384:       d = PetscPowReal(2.0,e);
385:       aux[cols[j]] = d*d;
386:       emaxl = PetscMax(emaxl,e);
387:       eminl = PetscMin(eminl,e);
388:     }
389:     /* Scale M */
390:     MatSeqAIJGetArray(M,&array);
391:     for (j=0;j<nz;j++) {
392:       array[j] *= aux[cidx[j]];
393:     }
394:     MatSeqAIJRestoreArray(M,&array);
395:     /* Row sum */    
396:     PetscMemzero(rsum,nr*sizeof(PetscReal));
397:     MatSeqAIJGetArray(M,&array);
398:     for (i=0;i<nr;i++) {
399:       for (j=ridx[i];j<ridx[i+1];j++) rsum[i] += PetscAbsScalar(array[j]);
400:       /* Update Dl */
401:       d = PetscLogReal(rsum[i])/l2;
402:       e = -(PetscInt)((d < 0)?(d-0.5):(d+0.5));
403:       d = PetscPowReal(2.0,e);
404:       Dl[i] *= d;
405:       /* Scale M */
406:       for (j=ridx[i];j<ridx[i+1];j++) array[j] *= d*d;
407:       emaxl = PetscMax(emaxl,e);
408:       eminl = PetscMin(eminl,e);      
409:     }
410:     MatSeqAIJRestoreArray(M,&array);  
411:     /* Compute global max and min */
412:     MPI_Allreduce(&emaxl,&emax,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)pep->Dl));
413:     MPI_Allreduce(&eminl,&emin,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)pep->Dl));
414:     if (emax<=emin+2) cont = PETSC_FALSE;
415:   }
416:   VecRestoreArray(pep->Dr,&Dr);
417:   VecRestoreArray(pep->Dl,&Dl);
418:   /* Free memory*/
419:   MatDestroy(&M);
420:   PetscFree4(rsum,csum,aux,cols);
421:   PetscFree(T);
422:   return(0);
423: }

427: /*
428:    PEPComputeScaleFactor - compute sfactor as described in [Betcke 2008].
429: */
430: PetscErrorCode PEPComputeScaleFactor(PEP pep)
431: {
433:   PetscBool      has0,has1,flg;
434:   PetscReal      norm0,norm1;
435:   Mat            T[2];
436:   PEPBasis       basis;
437:   PetscInt       i;

440:   if (pep->scale==PEP_SCALE_NONE || pep->scale==PEP_SCALE_DIAGONAL) {  /* no scalar scaling */
441:     pep->sfactor = 1.0;
442:     pep->dsfactor = 1.0;
443:     return(0);
444:   }
445:   if (pep->sfactor_set) return(0);  /* user provided value */
446:   pep->sfactor = 1.0;
447:   pep->dsfactor = 1.0;
448:   PEPGetBasis(pep,&basis);
449:   if (basis==PEP_BASIS_MONOMIAL) {
450:     STGetTransform(pep->st,&flg);
451:     if (flg) {
452:       STGetTOperators(pep->st,0,&T[0]);
453:       STGetTOperators(pep->st,pep->nmat-1,&T[1]);
454:     } else {
455:       T[0] = pep->A[0];
456:       T[1] = pep->A[pep->nmat-1];
457:     }
458:     if (pep->nmat>2) {
459:       MatHasOperation(T[0],MATOP_NORM,&has0);
460:       MatHasOperation(T[1],MATOP_NORM,&has1);
461:       if (has0 && has1) {
462:         MatNorm(T[0],NORM_INFINITY,&norm0);
463:         MatNorm(T[1],NORM_INFINITY,&norm1);
464:         pep->sfactor = PetscPowReal(norm0/norm1,1.0/(pep->nmat-1));
465:         pep->dsfactor = norm1;
466:         for (i=pep->nmat-2;i>0;i--) {
467:           STGetTOperators(pep->st,i,&T[1]);
468:           MatHasOperation(T[1],MATOP_NORM,&has1);
469:           if (has1) {
470:             MatNorm(T[1],NORM_INFINITY,&norm1);
471:             pep->dsfactor = pep->dsfactor*pep->sfactor+norm1;
472:           } else break;
473:         }
474:         if (has1) {
475:           pep->dsfactor = pep->dsfactor*pep->sfactor+norm0;
476:           pep->dsfactor = pep->nmat/pep->dsfactor;
477:         } else pep->dsfactor = 1.0;
478:       } 
479:     }
480:   } 
481:   return(0);
482: }

486: /*
487:    PEPBasisCoefficients - compute polynomial basis coefficients
488: */
489: PetscErrorCode PEPBasisCoefficients(PEP pep,PetscReal *pbc)
490: {
491:   PetscReal *ca,*cb,*cg;
492:   PetscInt  k,nmat=pep->nmat;
493:   
495:   ca = pbc;
496:   cb = pbc+nmat;
497:   cg = pbc+2*nmat;
498:   switch (pep->basis) {
499:   case PEP_BASIS_MONOMIAL:
500:     for (k=0;k<nmat;k++) {
501:       ca[k] = 1.0; cb[k] = 0.0; cg[k] = 0.0;
502:     }
503:     break;
504:   case PEP_BASIS_CHEBYSHEV1:
505:     ca[0] = 1.0; cb[0] = 0.0; cg[0] = 0.0;
506:     for (k=1;k<nmat;k++) {
507:       ca[k] = .5; cb[k] = 0.0; cg[k] = .5;
508:     }
509:     break;
510:   case PEP_BASIS_CHEBYSHEV2:
511:     ca[0] = .5; cb[0] = 0.0; cg[0] = 0.0;
512:     for (k=1;k<nmat;k++) {
513:       ca[k] = .5; cb[k] = 0.0; cg[k] = .5;
514:     }    
515:     break;
516:   case PEP_BASIS_LEGENDRE:
517:     ca[0] = 1.0; cb[0] = 0.0; cg[0] = 0.0;
518:     for (k=1;k<nmat;k++) {
519:       ca[k] = k+1; cb[k] = -2*k; cg[k] = k;
520:     }
521:     break;
522:   case PEP_BASIS_LAGUERRE:
523:     ca[0] = -1.0; cb[0] = 0.0; cg[0] = 0.0;
524:     for (k=1;k<nmat;k++) {
525:       ca[k] = -(k+1); cb[k] = 2*k+1; cg[k] = -k;
526:     }
527:     break;
528:   case PEP_BASIS_HERMITE:
529:     ca[0] = .5; cb[0] = 0.0; cg[0] = 0.0;
530:     for (k=1;k<nmat;k++) {
531:       ca[k] = .5; cb[k] = 0.0; cg[k] = -k;
532:     }
533:     break;
534:   default:
535:     SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'basis' value");
536:   }
537:   return(0);
538: }

542: /*
543:    PEPEvaluateBasis - evaluate the polynomial basis on a given parameter sigma
544: */
545: PetscErrorCode PEPEvaluateBasis(PEP pep,PetscScalar sigma,PetscScalar isigma,PetscScalar *vals,PetscScalar *ivals)
546: {
547:   PetscInt   nmat=pep->nmat,k;
548:   PetscReal  *a=pep->pbc,*b=pep->pbc+nmat,*g=pep->pbc+2*nmat;
549:   
551:   if (ivals) for (k=0;k<nmat;k++) ivals[k] = 0.0;
552:   vals[0] = 1.0;
553:   vals[1] = (sigma-b[0])/a[0];
554: #if !defined(PETSC_USE_COMPLEX)
555:   if (ivals) ivals[1] = isigma/a[0];
556: #endif
557:   for (k=2;k<nmat;k++) {
558:     vals[k] = ((sigma-b[k-1])*vals[k-1]-g[k-1]*vals[k-2])/a[k-1];
559:     if (ivals) vals[k] -= isigma*ivals[k-1]/a[k-1];
560: #if !defined(PETSC_USE_COMPLEX)
561:     if (ivals) ivals[k] = ((sigma-b[k-1])*ivals[k-1]+isigma*vals[k-1]-g[k-1]*ivals[k-2])/a[k-1];
562: #endif
563:   }
564:   return(0);
565: }