Actual source code: lapack.c

  1: /*
  2:        This file implements a wrapper to the LAPACK eigenvalue subroutines.
  3:        Generalized problems are transformed to standard ones only if necessary.

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

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

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

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

 25: #include <slepc-private/epsimpl.h>     /*I "slepceps.h" I*/
 26: #include <slepcblaslapack.h>

 30: PetscErrorCode EPSSetUp_LAPACK(EPS eps)
 31: {
 32:   PetscErrorCode ierr,ierra,ierrb;
 33:   PetscBool      isshift,denseok=PETSC_FALSE;
 34:   Mat            A,B,OP,Adense,Bdense;
 35:   PetscScalar    shift,*Ap,*Bp;
 36:   PetscInt       i,ld;
 37:   KSP            ksp;
 38:   PC             pc;
 39:   Vec            v;
 40: 
 42:   eps->ncv = eps->n;
 43:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 44:   if (!eps->which) { EPSDefaultSetWhich(eps); }
 45:   if (eps->balance!=EPS_BALANCE_NONE) { PetscInfo(eps,"Warning: balancing ignored\n"); }
 46:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
 47:   EPSAllocateSolution(eps);

 49:   /* attempt to get dense representations of A and B separately */
 50:   PetscObjectTypeCompare((PetscObject)eps->OP,STSHIFT,&isshift);
 51:   if (isshift) {
 52:     STGetOperators(eps->OP,&A,&B);
 53:     PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
 54:     ierra = SlepcMatConvertSeqDense(A,&Adense);
 55:     if (eps->isgeneralized) {
 56:       ierrb = SlepcMatConvertSeqDense(B,&Bdense);
 57:     } else {
 58:       ierrb = 0;
 59:     }
 60:     PetscPopErrorHandler();
 61:     denseok = (ierra == 0 && ierrb == 0)? PETSC_TRUE: PETSC_FALSE;
 62:   } else Adense = PETSC_NULL;

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

 90:   if (denseok) {
 91:     STGetShift(eps->OP,&shift);
 92:     if (shift != 0.0) {
 93:       MatShift(Adense,shift);
 94:     }
 95:     /* use dummy pc and ksp to avoid problems when B is not positive definite */
 96:     STGetKSP(eps->OP,&ksp);
 97:     KSPSetType(ksp,KSPPREONLY);
 98:     KSPGetPC(ksp,&pc);
 99:     PCSetType(pc,PCNONE);
100:   } else {
101:     PetscInfo(eps,"Using slow explicit operator\n");
102:     STComputeExplicitOperator(eps->OP,&OP);
103:     MatDestroy(&Adense);
104:     SlepcMatConvertSeqDense(OP,&Adense);
105:   }

107:   /* fill DS matrices */
108:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,ld,PETSC_NULL,&v);
109:   DSGetArray(eps->ds,DS_MAT_A,&Ap);
110:   for (i=0;i<ld;i++) {
111:     VecPlaceArray(v,Ap+i*ld);
112:     MatGetColumnVector(Adense,v,i);
113:     VecResetArray(v);
114:   }
115:   DSRestoreArray(eps->ds,DS_MAT_A,&Ap);
116:   if (denseok && eps->isgeneralized) {
117:     DSGetArray(eps->ds,DS_MAT_B,&Bp);
118:     for (i=0;i<ld;i++) {
119:       VecPlaceArray(v,Bp+i*ld);
120:       MatGetColumnVector(Bdense,v,i);
121:       VecResetArray(v);
122:     }
123:     DSRestoreArray(eps->ds,DS_MAT_B,&Bp);
124:   }
125:   VecDestroy(&v);
126:   MatDestroy(&Adense);
127:   if (!denseok) { MatDestroy(&OP); }
128:   if (denseok && eps->isgeneralized) { MatDestroy(&Bdense); }
129:   return(0);
130: }

134: PetscErrorCode EPSSolve_LAPACK(EPS eps)
135: {
137:   PetscInt       n=eps->n,i,low,high;
138:   PetscScalar    *array,*pX,*pY;
139: 
141:   DSSolve(eps->ds,eps->eigr,eps->eigi);
142:   DSSort(eps->ds,eps->eigr,eps->eigi,PETSC_NULL,PETSC_NULL,PETSC_NULL);

144:   /* right eigenvectors */
145:   DSVectors(eps->ds,DS_MAT_X,PETSC_NULL,PETSC_NULL);
146:   DSGetArray(eps->ds,DS_MAT_X,&pX);
147:   for (i=0;i<eps->ncv;i++) {
148:     VecGetOwnershipRange(eps->V[i],&low,&high);
149:     VecGetArray(eps->V[i],&array);
150:     PetscMemcpy(array,pX+i*n+low,(high-low)*sizeof(PetscScalar));
151:     VecRestoreArray(eps->V[i],&array);
152:   }
153:   DSRestoreArray(eps->ds,DS_MAT_X,&pX);

155:   /* left eigenvectors */
156:   if (eps->leftvecs) {
157:     DSVectors(eps->ds,DS_MAT_Y,PETSC_NULL,PETSC_NULL);
158:     DSGetArray(eps->ds,DS_MAT_Y,&pY);
159:     for (i=0;i<eps->ncv;i++) {
160:       VecGetOwnershipRange(eps->W[i],&low,&high);
161:       VecGetArray(eps->W[i],&array);
162:       PetscMemcpy(array,pY+i*n+low,(high-low)*sizeof(PetscScalar));
163:       VecRestoreArray(eps->W[i],&array);
164:     }
165:     DSRestoreArray(eps->ds,DS_MAT_Y,&pY);
166:   }
167:   eps->nconv  = eps->ncv;
168:   eps->its    = 1;
169:   eps->reason = EPS_CONVERGED_TOL;
170:   return(0);
171: }

175: PetscErrorCode EPSReset_LAPACK(EPS eps)
176: {

180:   EPSFreeSolution(eps);
181:   return(0);
182: }

184: EXTERN_C_BEGIN
187: PetscErrorCode EPSCreate_LAPACK(EPS eps)
188: {
190:   eps->ops->solve                = EPSSolve_LAPACK;
191:   eps->ops->setup                = EPSSetUp_LAPACK;
192:   eps->ops->reset                = EPSReset_LAPACK;
193:   eps->ops->backtransform        = EPSBackTransform_Default;
194:   eps->ops->computevectors       = EPSComputeVectors_Default;
195:   return(0);
196: }
197: EXTERN_C_END