Actual source code: pepdefault.c

slepc-3.5.2 2014-10-10
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-2014, 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:   PEPConvergedEigRelative - Checks convergence relative to the eigenvalue.
 64: */
 65: PetscErrorCode PEPConvergedEigRelative(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:   PEPConvergedNormRelative - Checks convergence relative to the matrix norms.
 79: */
 80: PetscErrorCode PEPConvergedNormRelative(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 81: {
 83:   PetscInt       i;
 84:   PetscReal      w=0.0;
 85:   PetscScalar    t[20],*vals=t,*ivals=NULL;
 86: #if !defined(PETSC_USE_COMPLEX)
 87:   PetscScalar    it[20];
 88: #endif

 91: #if !defined(PETSC_USE_COMPLEX)
 92:   ivals = it;
 93: #endif
 94:   if (pep->nmat>20) {
 95: #if !defined(PETSC_USE_COMPLEX)
 96:     PetscMalloc2(pep->nmat,&vals,pep->nmat,&ivals);
 97: #else
 98:     PetscMalloc1(pep->nmat,&vals);
 99: #endif
100:   }
101:   PEPEvaluateBasis(pep,eigr,eigi,vals,ivals);
102:   for (i=0;i<pep->nmat;i++) w += SlepcAbsEigenvalue(vals[i],ivals[i])*pep->nrma[i];
103:   *errest = res/w;
104:   if (pep->nmat>20) {
105: #if !defined(PETSC_USE_COMPLEX)
106:     PetscFree2(vals,ivals);
107: #else
108:     PetscFree(vals);
109: #endif
110:   }
111:   return(0);
112: }

116: /*
117:   PEPConvergedAbsolute - Checks convergence absolutely.
118: */
119: PetscErrorCode PEPConvergedAbsolute(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
120: {
122:   *errest = res;
123:   return(0);
124: }

128: PetscErrorCode PEPComputeVectors_Schur(PEP pep)
129: {
131:   PetscInt       n,i;
132:   Mat            Z;
133:   Vec            v;
134: #if !defined(PETSC_USE_COMPLEX)
135:   Vec            v1;
136:   PetscScalar    tmp;
137:   PetscReal      norm,normi;
138: #endif

141:   DSGetDimensions(pep->ds,&n,NULL,NULL,NULL,NULL);
142:   DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
143:   DSGetMat(pep->ds,DS_MAT_X,&Z);
144:   BVSetActiveColumns(pep->V,0,n);
145:   BVMultInPlace(pep->V,Z,0,n);
146:   MatDestroy(&Z);

148:   /* Fix eigenvectors if balancing was used */
149:   if ((pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) && pep->Dr && (pep->refine!=PEP_REFINE_MULTIPLE)) {
150:     for (i=0;i<n;i++) {
151:       BVGetColumn(pep->V,i,&v);
152:       VecPointwiseMult(v,v,pep->Dr);
153:       BVRestoreColumn(pep->V,i,&v);
154:     }
155:   }

157:   /* normalization */
158:   for (i=0;i<n;i++) {
159: #if !defined(PETSC_USE_COMPLEX)
160:     if (pep->eigi[i] != 0.0) {
161:       BVGetColumn(pep->V,i,&v);
162:       BVGetColumn(pep->V,i+1,&v1);
163:       VecNorm(v,NORM_2,&norm);
164:       VecNorm(v1,NORM_2,&normi);
165:       tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
166:       VecScale(v,tmp);
167:       VecScale(v1,tmp);
168:       BVRestoreColumn(pep->V,i,&v);
169:       BVRestoreColumn(pep->V,i+1,&v1);
170:       i++;
171:     } else
172: #endif
173:     {
174:       BVGetColumn(pep->V,i,&v);
175:       VecNormalize(v,NULL);
176:       BVRestoreColumn(pep->V,i,&v);
177:     }
178:   }
179:   return(0);
180: }

184: PetscErrorCode PEPComputeVectors_Indefinite(PEP pep)
185: {
187:   PetscInt       n,i;
188:   Mat            Z;
189:   Vec            v;
190: #if !defined(PETSC_USE_COMPLEX)
191:   Vec            v1;
192:   PetscScalar    tmp;
193:   PetscReal      norm,normi;
194: #endif

197:   DSGetDimensions(pep->ds,&n,NULL,NULL,NULL,NULL);
198:   DSVectors(pep->ds,DS_MAT_X,NULL,NULL);
199:   DSGetMat(pep->ds,DS_MAT_X,&Z);
200:   BVSetActiveColumns(pep->V,0,n);
201:   BVMultInPlace(pep->V,Z,0,n);
202:   MatDestroy(&Z);

204:   /* normalization */
205:   for (i=0;i<n;i++) {
206: #if !defined(PETSC_USE_COMPLEX)
207:     if (pep->eigi[i] != 0.0) {
208:       BVGetColumn(pep->V,i,&v);
209:       BVGetColumn(pep->V,i+1,&v1);
210:       VecNorm(v,NORM_2,&norm);
211:       VecNorm(v1,NORM_2,&normi);
212:       tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
213:       VecScale(v,tmp);
214:       VecScale(v1,tmp);
215:       BVRestoreColumn(pep->V,i,&v);
216:       BVRestoreColumn(pep->V,i+1,&v1);
217:       i++;
218:     } else
219: #endif
220:     {
221:       BVGetColumn(pep->V,i,&v);
222:       VecNormalize(v,NULL);
223:       BVRestoreColumn(pep->V,i,&v);
224:     }
225:   }
226:   return(0);
227: }

231: /*
232:    PEPKrylovConvergence - This is the analogue to EPSKrylovConvergence, but
233:    for polynomial Krylov methods.

235:    Differences:
236:    - Always non-symmetric
237:    - Does not check for STSHIFT
238:    - No correction factor
239:    - No support for true residual
240: */
241: PetscErrorCode PEPKrylovConvergence(PEP pep,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal beta,PetscInt *kout)
242: {
244:   PetscInt       k,newk,marker,ld,inside;
245:   PetscScalar    re,im;
246:   PetscReal      resnorm;
247:   PetscBool      istrivial;

250:   RGIsTrivial(pep->rg,&istrivial);
251:   DSGetLeadingDimension(pep->ds,&ld);
252:   marker = -1;
253:   if (pep->trackall) getall = PETSC_TRUE;
254:   for (k=kini;k<kini+nits;k++) {
255:     /* eigenvalue */
256:     re = pep->eigr[k];
257:     im = pep->eigi[k];
258:     if (!istrivial || pep->conv==PEP_CONV_NORM) {
259:       STBackTransform(pep->st,1,&re,&im);
260:     }
261:     if (!istrivial) {
262:       RGCheckInside(pep->rg,1,&re,&im,&inside);
263:       if (marker==-1 && inside<=0) marker = k;
264:       if (!pep->conv==PEP_CONV_NORM) {  /* make sure pep->converged below uses the right value */
265:         re = pep->eigr[k];
266:         im = pep->eigi[k];
267:       }
268:     }
269:     newk = k;
270:     DSVectors(pep->ds,DS_MAT_X,&newk,&resnorm);
271:     resnorm *= beta;
272:     /* error estimate */
273:     (*pep->converged)(pep,re,im,resnorm,&pep->errest[k],pep->convergedctx);
274:     if (marker==-1 && pep->errest[k] >= pep->tol) marker = k;
275:     if (newk==k+1) {
276:       pep->errest[k+1] = pep->errest[k];
277:       k++;
278:     }
279:     if (marker!=-1 && !getall) break;
280:   }
281:   if (marker!=-1) k = marker;
282:   *kout = k;
283:   return(0);
284: }

288: /*
289:   PEPBuildDiagonalScaling - compute two diagonal matrices to be applied for balancing 
290:   in polynomial eigenproblems.
291: */
292: PetscErrorCode PEPBuildDiagonalScaling(PEP pep)
293: {
295:   PetscInt       it,i,j,k,nmat,nr,e,nz,lst,lend,nc=0,*cols,emax,emin,emaxl,eminl;
296:   const PetscInt *cidx,*ridx;
297:   Mat            M,*T,A;
298:   PetscMPIInt    n;
299:   PetscBool      cont=PETSC_TRUE,flg=PETSC_FALSE;
300:   PetscScalar    *array,*Dr,*Dl,t;
301:   PetscReal      l2,d,*rsum,*aux,*csum,w=1.0;
302:   MatStructure   str;
303:   MatInfo        info;

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

440: /*
441:    PEPComputeScaleFactor - compute sfactor as described in [Betcke 2008].
442: */
443: PetscErrorCode PEPComputeScaleFactor(PEP pep)
444: {
446:   PetscBool      has0,has1,flg;
447:   PetscReal      norm0,norm1;
448:   Mat            T[2];
449:   PEPBasis       basis;

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

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

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