Actual source code: lapack.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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:    This file implements a wrapper to the LAPACK eigenvalue subroutines.
 12:    Generalized problems are transformed to standard ones only if necessary.
 13: */

 15: #include <slepc/private/epsimpl.h>

 17: PetscErrorCode EPSSetUp_LAPACK(EPS eps)
 18: {
 19:   PetscErrorCode ierr,ierra,ierrb;
 20:   PetscBool      isshift,flg,denseok=PETSC_FALSE;
 21:   Mat            A,B,OP,shell,Ar,Br,Adense=NULL,Bdense=NULL;
 22:   PetscScalar    shift,*Ap,*Bp;
 23:   PetscInt       i,ld,nmat;
 24:   KSP            ksp;
 25:   PC             pc;
 26:   Vec            v;

 29:   eps->ncv = eps->n;
 30:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 31:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 32:   if (eps->which==EPS_ALL && eps->inta!=eps->intb) SETERRQ(PetscObjectComm((PetscObject)eps),1,"This solver does not support interval computation");
 33:   if (eps->balance!=EPS_BALANCE_NONE) { PetscInfo(eps,"Warning: balancing ignored\n"); }
 34:   if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"User-defined stopping test not supported");
 35:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
 36:   EPSAllocateSolution(eps,0);

 38:   /* attempt to get dense representations of A and B separately */
 39:   PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
 40:   if (isshift) {
 41:     STGetNumMatrices(eps->st,&nmat);
 42:     STGetMatrix(eps->st,0,&A);
 43:     MatHasOperation(A,MATOP_CREATE_SUBMATRICES,&flg);
 44:     if (flg) {
 45:       PetscPushErrorHandler(PetscIgnoreErrorHandler,NULL);
 46:       ierra  = MatCreateRedundantMatrix(A,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Ar);
 47:       if (!ierra) { ierra |= MatConvert(Ar,MATSEQDENSE,MAT_INITIAL_MATRIX,&Adense); }
 48:       ierra |= MatDestroy(&Ar);
 49:       PetscPopErrorHandler();
 50:     } else ierra = 1;
 51:     if (nmat>1) {
 52:       STGetMatrix(eps->st,1,&B);
 53:       MatHasOperation(B,MATOP_CREATE_SUBMATRICES,&flg);
 54:       if (flg) {
 55:         PetscPushErrorHandler(PetscIgnoreErrorHandler,NULL);
 56:         ierrb  = MatCreateRedundantMatrix(B,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Br);
 57:         if (!ierrb) { ierrb |= MatConvert(Br,MATSEQDENSE,MAT_INITIAL_MATRIX,&Bdense); }
 58:         ierrb |= MatDestroy(&Br);
 59:         PetscPopErrorHandler();
 60:       } else ierrb = 1;
 61:     } else ierrb = 0;
 62:     denseok = PetscNot(ierra || ierrb);
 63:   }

 65:   /* setup DS */
 66:   if (denseok) {
 67:     if (eps->isgeneralized) {
 68:       if (eps->ishermitian) {
 69:         if (eps->ispositive) {
 70:           DSSetType(eps->ds,DSGHEP);
 71:         } else {
 72:           DSSetType(eps->ds,DSGNHEP); /* TODO: should be DSGHIEP */
 73:         }
 74:       } else {
 75:         DSSetType(eps->ds,DSGNHEP);
 76:       }
 77:     } else {
 78:       if (eps->ishermitian) {
 79:         DSSetType(eps->ds,DSHEP);
 80:       } else {
 81:         DSSetType(eps->ds,DSNHEP);
 82:       }
 83:     }
 84:   } else {
 85:     DSSetType(eps->ds,DSNHEP);
 86:   }
 87:   DSAllocate(eps->ds,eps->ncv);
 88:   DSGetLeadingDimension(eps->ds,&ld);
 89:   DSSetDimensions(eps->ds,eps->ncv,0,0,0);

 91:   if (denseok) {
 92:     STGetShift(eps->st,&shift);
 93:     if (shift != 0.0) {
 94:       MatShift(Adense,shift);
 95:     }
 96:     /* use dummy pc and ksp to avoid problems when B is not positive definite */
 97:     STGetKSP(eps->st,&ksp);
 98:     KSPSetType(ksp,KSPPREONLY);
 99:     KSPGetPC(ksp,&pc);
100:     PCSetType(pc,PCNONE);
101:   } else {
102:     PetscInfo(eps,"Using slow explicit operator\n");
103:     STGetOperator(eps->st,&shell);
104:     MatComputeExplicitOperator(shell,&OP);
105:     MatDestroy(&shell);
106:     MatDestroy(&Adense);
107:     MatCreateRedundantMatrix(OP,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Ar);
108:     MatDestroy(&OP);
109:     MatConvert(Ar,MATSEQDENSE,MAT_INITIAL_MATRIX,&Adense);
110:     MatDestroy(&Ar);
111:   }

113:   /* fill DS matrices */
114:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,ld,NULL,&v);
115:   DSGetArray(eps->ds,DS_MAT_A,&Ap);
116:   for (i=0;i<ld;i++) {
117:     VecPlaceArray(v,Ap+i*ld);
118:     MatGetColumnVector(Adense,v,i);
119:     VecResetArray(v);
120:   }
121:   DSRestoreArray(eps->ds,DS_MAT_A,&Ap);
122:   if (denseok && eps->isgeneralized) {
123:     DSGetArray(eps->ds,DS_MAT_B,&Bp);
124:     for (i=0;i<ld;i++) {
125:       VecPlaceArray(v,Bp+i*ld);
126:       MatGetColumnVector(Bdense,v,i);
127:       VecResetArray(v);
128:     }
129:     DSRestoreArray(eps->ds,DS_MAT_B,&Bp);
130:   }
131:   VecDestroy(&v);
132:   DSSetState(eps->ds,DS_STATE_RAW);
133:   MatDestroy(&Adense);
134:   MatDestroy(&Bdense);
135:   return(0);
136: }

138: PetscErrorCode EPSSolve_LAPACK(EPS eps)
139: {
141:   PetscInt       n=eps->n,i,low,high;
142:   PetscScalar    *array,*pX,*pY;
143:   Vec            v,w;

146:   DSSolve(eps->ds,eps->eigr,eps->eigi);
147:   DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
148:   DSSynchronize(eps->ds,eps->eigr,eps->eigi);

150:   /* right eigenvectors */
151:   DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
152:   DSGetArray(eps->ds,DS_MAT_X,&pX);
153:   for (i=0;i<eps->ncv;i++) {
154:     BVGetColumn(eps->V,i,&v);
155:     VecGetOwnershipRange(v,&low,&high);
156:     VecGetArray(v,&array);
157:     PetscMemcpy(array,pX+i*n+low,(high-low)*sizeof(PetscScalar));
158:     VecRestoreArray(v,&array);
159:     BVRestoreColumn(eps->V,i,&v);
160:   }
161:   DSRestoreArray(eps->ds,DS_MAT_X,&pX);

163:   /* left eigenvectors */
164:   if (eps->twosided) {
165:     DSVectors(eps->ds,DS_MAT_Y,NULL,NULL);
166:     DSGetArray(eps->ds,DS_MAT_Y,&pY);
167:     for (i=0;i<eps->ncv;i++) {
168:       BVGetColumn(eps->W,i,&w);
169:       VecGetOwnershipRange(w,&low,&high);
170:       VecGetArray(w,&array);
171:       PetscMemcpy(array,pY+i*n+low,(high-low)*sizeof(PetscScalar));
172:       VecRestoreArray(w,&array);
173:       BVRestoreColumn(eps->W,i,&w);
174:     }
175:     DSRestoreArray(eps->ds,DS_MAT_Y,&pY);
176:   }

178:   eps->nconv  = eps->ncv;
179:   eps->its    = 1;
180:   eps->reason = EPS_CONVERGED_TOL;
181:   return(0);
182: }

184: SLEPC_EXTERN PetscErrorCode EPSCreate_LAPACK(EPS eps)
185: {
187:   eps->useds = PETSC_TRUE;
188:   eps->hasts = PETSC_TRUE;
189:   eps->categ = EPS_CATEGORY_OTHER;

191:   eps->ops->solve          = EPSSolve_LAPACK;
192:   eps->ops->setup          = EPSSetUp_LAPACK;
193:   eps->ops->backtransform  = EPSBackTransform_Default;
194:   return(0);
195: }