Actual source code: peputils.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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:    Common subroutines for several PEP solvers
 12: */

 14: #include <slepc/private/pepimpl.h>
 15: #include <slepcblaslapack.h>

 17: /*
 18:   Computes T_j=phy_idx(T) for a matrix T.
 19:     Tp (if j>0) and Tpp (if j>1) are the evaluation
 20:     of phy_(j-1) and phy(j-2)respectively.
 21: */
 22: PetscErrorCode PEPEvaluateBasisMat(PEP pep,PetscInt k,PetscScalar *T,PetscInt ldt,PetscInt idx,PetscScalar *Tpp,PetscInt ldtpp,PetscScalar *Tp,PetscInt ldtp,PetscScalar *Tj,PetscInt ldtj)
 23: {
 24:   PetscInt       i;
 25:   PetscReal      *ca,*cb,*cg;
 26:   PetscScalar    g,a;
 27:   PetscBLASInt   k_,ldt_,ldtj_,ldtp_;

 29:   if (idx==0) {
 30:     for (i=0;i<k;i++) {
 31:       PetscArrayzero(Tj+i*ldtj,k);
 32:       Tj[i+i*ldtj] = 1.0;
 33:     }
 34:   } else {
 35:     PetscBLASIntCast(ldt,&ldt_);
 36:     PetscBLASIntCast(ldtj,&ldtj_);
 37:     PetscBLASIntCast(ldtp,&ldtp_);
 38:     PetscBLASIntCast(k,&k_);
 39:     ca = pep->pbc; cb = pep->pbc+pep->nmat; cg = pep->pbc+2*pep->nmat;
 40:     for (i=0;i<k;i++) T[i*ldt+i] -= cb[idx-1];
 41:     if (idx>1) {
 42:       for (i=0;i<k;i++) PetscArraycpy(Tj+i*ldtj,Tpp+i*ldtpp,k);
 43:     }
 44:     a = 1/ca[idx-1];
 45:     g = (idx==1)?0.0:-cg[idx-1]/ca[idx-1];
 46:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&k_,&k_,&k_,&a,T,&ldt_,Tp,&ldtp_,&g,Tj,&ldtj_));
 47:     for (i=0;i<k;i++) T[i*ldt+i] += cb[idx-1];
 48:   }
 49:   PetscFunctionReturn(0);
 50: }

 52: /*
 53:    PEPEvaluateBasis - evaluate the polynomial basis on a given parameter sigma
 54: */
 55: PetscErrorCode PEPEvaluateBasis(PEP pep,PetscScalar sigma,PetscScalar isigma,PetscScalar *vals,PetscScalar *ivals)
 56: {
 57:   PetscInt   nmat=pep->nmat,k;
 58:   PetscReal  *a=pep->pbc,*b=pep->pbc+nmat,*g=pep->pbc+2*nmat;

 60:   if (ivals) for (k=0;k<nmat;k++) ivals[k] = 0.0;
 61:   vals[0] = 1.0;
 62:   vals[1] = (sigma-b[0])/a[0];
 63: #if !defined(PETSC_USE_COMPLEX)
 64:   if (ivals) ivals[1] = isigma/a[0];
 65: #endif
 66:   for (k=2;k<nmat;k++) {
 67:     vals[k] = ((sigma-b[k-1])*vals[k-1]-g[k-1]*vals[k-2])/a[k-1];
 68:     if (ivals) vals[k] -= isigma*ivals[k-1]/a[k-1];
 69: #if !defined(PETSC_USE_COMPLEX)
 70:     if (ivals) ivals[k] = ((sigma-b[k-1])*ivals[k-1]+isigma*vals[k-1]-g[k-1]*ivals[k-2])/a[k-1];
 71: #endif
 72:   }
 73:   PetscFunctionReturn(0);
 74: }

 76: /*
 77:    PEPEvaluateBasisDerivative - evaluate the derivative of the polynomial basis on a given parameter sigma
 78: */
 79: PetscErrorCode PEPEvaluateBasisDerivative(PEP pep,PetscScalar sigma,PetscScalar isigma,PetscScalar *vals,PetscScalar *ivals)
 80: {
 81:   PetscInt       nmat=pep->nmat,k;
 82:   PetscReal      *a=pep->pbc,*b=pep->pbc+nmat,*g=pep->pbc+2*nmat;

 84:   PEPEvaluateBasis(pep,sigma,isigma,vals,ivals);
 85:   for (k=nmat-1;k>0;k--) {
 86:     vals[k] = vals[k-1];
 87:     if (ivals) ivals[k] = ivals[k-1];
 88:   }
 89:   vals[0] = 0.0;
 90:   vals[1] = vals[1]/a[0];
 91: #if !defined(PETSC_USE_COMPLEX)
 92:   if (ivals) ivals[1] = ivals[1]/a[0];
 93: #endif
 94:   for (k=2;k<nmat;k++) {
 95:     vals[k] += (sigma-b[k-1])*vals[k-1]-g[k-1]*vals[k-2];
 96:     if (ivals) vals[k] -= isigma*ivals[k-1];
 97:     vals[k] /= a[k-1];
 98: #if !defined(PETSC_USE_COMPLEX)
 99:     if (ivals) {
100:       ivals[k] += (sigma-b[k-1])*ivals[k-1]+isigma*vals[k-1]-g[k-1]*ivals[k-2];
101:       ivals[k] /= a[k-1];
102:     }
103: #endif
104:   }
105:   PetscFunctionReturn(0);
106: }