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-2010, Universidad 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
 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:   PetscTruth     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) eps->which = EPS_LARGEST_MAGNITUDE;
 49:   if (eps->balance!=EPS_BALANCE_NONE)
 50:     SETERRQ(PETSC_ERR_SUP,"Balancing not supported in Lapack solvers");

 52:   if (la->OP) { MatDestroy(la->OP); }
 53:   if (la->A) { MatDestroy(la->A); }
 54:   if (la->B) { 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) {
 93:      PetscInfo(eps,"Warning: extraction type ignored\n");
 94:   }
 95:   EPSAllocateSolution(eps);
 96:   return(0);
 97: }

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

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

188:   eps->nconv = eps->ncv;
189:   eps->its   = 1;
190:   eps->reason = EPS_CONVERGED_TOL;
191: 
192:   return(0);
193: }

197: PetscErrorCode EPSDestroy_LAPACK(EPS eps)
198: {
200:   EPS_LAPACK     *la = (EPS_LAPACK *)eps->data;

204:   if (la->OP) { MatDestroy(la->OP); }
205:   if (la->A) { MatDestroy(la->A); }
206:   if (la->B) { MatDestroy(la->B); }
207:   PetscFree(eps->data);
208:   EPSFreeSolution(eps);
209:   return(0);
210: }

215: PetscErrorCode EPSCreate_LAPACK(EPS eps)
216: {
218:   EPS_LAPACK     *la;

221:   PetscNew(EPS_LAPACK,&la);
222:   PetscLogObjectMemory(eps,sizeof(EPS_LAPACK));
223:   eps->data                      = (void *) la;
224:   eps->ops->solve                = EPSSolve_LAPACK;
225:   eps->ops->setup                = EPSSetUp_LAPACK;
226:   eps->ops->destroy              = EPSDestroy_LAPACK;
227:   eps->ops->backtransform        = EPSBackTransform_Default;
228:   eps->ops->computevectors       = EPSComputeVectors_Default;
229:   return(0);
230: }