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-2011, 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 <private/epsimpl.h>     /*I "slepceps.h" I*/
 26: #include <slepcblaslapack.h>

 28: typedef struct {
 29:   Mat OP,A,B;
 30: } EPS_LAPACK;

 34: PetscErrorCode EPSSetUp_LAPACK(EPS eps)
 35: {
 36:   PetscErrorCode ierr,ierra,ierrb;
 37:   EPS_LAPACK     *la = (EPS_LAPACK *)eps->data;
 38:   PetscBool      flg;
 39:   Mat            A,B;
 40:   PetscScalar    shift;
 41:   KSP            ksp;
 42:   PC             pc;
 43: 
 45:   eps->ncv = eps->n;
 46:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }

 48:   if (!eps->which) { EPSDefaultSetWhich(eps); }
 49:   if (eps->balance!=EPS_BALANCE_NONE)
 50:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Balancing not supported in Lapack solvers");

 52:   MatDestroy(&la->OP);
 53:   MatDestroy(&la->A);
 54:   MatDestroy(&la->B);

 56:   PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&flg);
 57:   STGetOperators(eps->OP,&A,&B);
 58: 
 59:   if (flg) {
 60:     la->OP = PETSC_NULL;
 61:     PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
 62:     ierra = SlepcMatConvertSeqDense(A,&la->A);
 63:     if (eps->isgeneralized) {
 64:       ierrb = SlepcMatConvertSeqDense(B,&la->B);
 65:     } else {
 66:       ierrb = 0;
 67:       la->B = PETSC_NULL;
 68:     }
 69:     PetscPopErrorHandler();
 70:     if (ierra == 0 && ierrb == 0) {
 71:       STGetShift(eps->OP,&shift);
 72:       if (shift != 0.0) {
 73:         MatShift(la->A,shift);
 74:       }
 75:       /* use dummy pc and ksp to avoid problems when B is not positive definite */
 76:       STGetKSP(eps->OP,&ksp);
 77:       KSPSetType(ksp,KSPPREONLY);
 78:       KSPGetPC(ksp,&pc);
 79:       PCSetType(pc,PCNONE);
 80:       EPSAllocateSolution(eps);
 81:       return(0);
 82:     }
 83:   }
 84:   PetscInfo(eps,"Using slow explicit operator\n");
 85:   la->A = PETSC_NULL;
 86:   la->B = PETSC_NULL;
 87:   STComputeExplicitOperator(eps->OP,&la->OP);
 88:   PetscTypeCompare((PetscObject)la->OP,MATSEQDENSE,&flg);
 89:   if (!flg) {
 90:     SlepcMatConvertSeqDense(la->OP,&la->OP);
 91:   }
 92:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
 93:   EPSAllocateSolution(eps);
 94:   return(0);
 95: }

 99: PetscErrorCode EPSSolve_LAPACK(EPS eps)
100: {
102:   PetscInt       n=eps->n,i,low,high;
103:   PetscMPIInt    size;
104:   PetscScalar    *array,*arrayb,*pV,*pW;
105:   PetscReal      *w;
106:   EPS_LAPACK     *la = (EPS_LAPACK *)eps->data;
107:   MPI_Comm       comm = ((PetscObject)eps)->comm;
108: 
110:   MPI_Comm_size(comm,&size);
111:   if (size == 1) {
112:     VecGetArray(eps->V[0],&pV);
113:   } else {
114:     PetscMalloc(sizeof(PetscScalar)*n*n,&pV);
115:   }
116:   if (eps->leftvecs) {
117:     if (size == 1) {
118:       VecGetArray(eps->W[0],&pW);
119:     } else {
120:       PetscMalloc(sizeof(PetscScalar)*n*n,&pW);
121:     }
122:   } else pW = PETSC_NULL;
123: 
124:   if (la->OP) {
125:     MatGetArray(la->OP,&array);
126:     EPSDenseNHEP(n,array,eps->eigr,eps->eigi,pV,pW);
127:     MatRestoreArray(la->OP,&array);
128:   } else if (eps->ishermitian) {
129: #if defined(PETSC_USE_COMPLEX)
130:     PetscMalloc(n*sizeof(PetscReal),&w);
131: #else
132:     w = eps->eigr;
133: #endif
134:     MatGetArray(la->A,&array);
135:     if (!eps->isgeneralized) {
136:       EPSDenseHEP(n,array,n,w,pV);
137:     } else {
138:       MatGetArray(la->B,&arrayb);
139:       EPSDenseGHEP(n,array,arrayb,w,pV);
140:       MatRestoreArray(la->B,&arrayb);
141:     }
142:     MatRestoreArray(la->A,&array);
143: #if defined(PETSC_USE_COMPLEX)
144:     for (i=0;i<n;i++) {
145:       eps->eigr[i] = w[i];
146:     }
147:     PetscFree(w);
148: #endif    
149:   } else {
150:     MatGetArray(la->A,&array);
151:     if (!eps->isgeneralized) {
152:       EPSDenseNHEP(n,array,eps->eigr,eps->eigi,pV,pW);
153:     } else {
154:       MatGetArray(la->B,&arrayb);
155:       EPSDenseGNHEP(n,array,arrayb,eps->eigr,eps->eigi,pV,pW);
156:       MatRestoreArray(la->B,&arrayb);
157:     }
158:     MatRestoreArray(la->A,&array);
159:   }

161:   if (size == 1) {
162:     VecRestoreArray(eps->V[0],&pV);
163:   } else {
164:     for (i=0; i<eps->ncv; i++) {
165:       VecGetOwnershipRange(eps->V[i],&low,&high);
166:       VecGetArray(eps->V[i],&array);
167:       PetscMemcpy(array,pV+i*n+low,(high-low)*sizeof(PetscScalar));
168:       VecRestoreArray(eps->V[i],&array);
169:     }
170:     PetscFree(pV);
171:   }
172:   if (pW) {
173:     if (size == 1) {
174:       VecRestoreArray(eps->W[0],&pW);
175:     } else {
176:       for (i=0; i<eps->ncv; i++) {
177:         VecGetOwnershipRange(eps->W[i],&low,&high);
178:         VecGetArray(eps->W[i],&array);
179:         PetscMemcpy(array,pW+i*n+low,(high-low)*sizeof(PetscScalar));
180:         VecRestoreArray(eps->W[i],&array);
181:       }
182:       PetscFree(pW);
183:     }
184:   }
185:   eps->nconv = eps->ncv;
186:   eps->its   = 1;
187:   eps->reason = EPS_CONVERGED_TOL;
188:   return(0);
189: }

193: PetscErrorCode EPSReset_LAPACK(EPS eps)
194: {
196:   EPS_LAPACK     *la = (EPS_LAPACK *)eps->data;

199:   MatDestroy(&la->OP);
200:   MatDestroy(&la->A);
201:   MatDestroy(&la->B);
202:   EPSFreeSolution(eps);
203:   return(0);
204: }

208: PetscErrorCode EPSDestroy_LAPACK(EPS eps)
209: {

213:   PetscFree(eps->data);
214:   return(0);
215: }

220: PetscErrorCode EPSCreate_LAPACK(EPS eps)
221: {

225:   PetscNewLog(eps,EPS_LAPACK,&eps->data);
226:   eps->ops->solve                = EPSSolve_LAPACK;
227:   eps->ops->setup                = EPSSetUp_LAPACK;
228:   eps->ops->destroy              = EPSDestroy_LAPACK;
229:   eps->ops->reset                = EPSReset_LAPACK;
230:   eps->ops->backtransform        = EPSBackTransform_Default;
231:   eps->ops->computevectors       = EPSComputeVectors_Default;
232:   return(0);
233: }