Actual source code: lapack.c

  2: /*                       
  3:        This file implements a wrapper to the LAPACK eigenvalue subroutines.
  4:        Currently, only LAPACK routines for standard problems are used.
  5:        Generalized problems are transformed to standard ones.
  6: */
 7:  #include src/eps/impls/lapack/lapackp.h
 8:  #include slepcblaslapack.h

 12: PetscErrorCode EPSSetUp_LAPACK(EPS eps)
 13: {
 14:   PetscErrorCode ierr;
 15:   int            size,rank,n,N;
 16:   EPS_LAPACK     *la = (EPS_LAPACK *)eps->data;
 17:   MPI_Comm       comm = eps->comm;
 18:   Mat            *T;
 19:   IS             isrow, iscol;
 20:   PetscTruth     flg;

 23:   VecGetSize(eps->vec_initial,&N);
 24:   if (eps->nev<1 || eps->nev>N) SETERRQ(1,"Wrong value of nev");
 25:   if (eps->ncv) {
 26:     if (eps->ncv<=eps->nev) SETERRQ(1,"Wrong value of ncv");
 27: #ifndef PETSC_USE_COMPLEX
 28:   } else eps->ncv = PetscMin(eps->nev + 1, N);
 29: #else
 30:   } else eps->ncv = eps->nev;
 31: #endif

 33:   if (la->BA) { MatDestroy(la->BA); }
 34:   EPSComputeExplicitOperator(eps,&la->BA);
 35:   MPI_Comm_size(comm,&size);
 36:   MPI_Comm_rank(comm,&rank);
 37: 
 38:   if (size > 1) {
 39:     /* assemble full matrix on every processor */
 40:     MatGetSize(la->BA,&n,&n);
 41:     ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &isrow);
 42:     ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &iscol);
 43:     MatGetSubMatrices(la->BA, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &T);
 44:     ISDestroy(isrow);
 45:     ISDestroy(iscol);
 46:     MatDestroy(la->BA);
 47:     la->BA = *T;
 48:   }
 49: 
 50:   PetscTypeCompare((PetscObject)la->BA, MATSEQDENSE, &flg);
 51:   if (!flg) {
 52:     /* convert matrix to MatSeqDense */
 53:     MatConvert(la->BA, MATSEQDENSE, &la->BA);
 54:   }

 56:   EPSAllocateSolutionContiguous(eps);
 57:   return(0);
 58: }

 62: PetscErrorCode EPSSolve_LAPACK(EPS eps)
 63: {
 64:   PetscErrorCode ierr;
 65:   int            n,size,rank,i,low,high;
 66:   PetscScalar    *array,*pV;
 67:   EPS_LAPACK     *la = (EPS_LAPACK *)eps->data;
 68:   MPI_Comm       comm = eps->comm;
 69: 
 71:   MPI_Comm_size(comm,&size);
 72:   MPI_Comm_rank(comm,&rank);
 73: 
 74:   MatGetArray(la->BA,&array);
 75:   MatGetSize(la->BA,&n,&n);

 77:   if (size == 1) {
 78:     VecGetArray(eps->V[0],&pV);
 79:   } else {
 80:     PetscMalloc(sizeof(PetscScalar)*n,&pV);
 81:   }
 82: 
 83:   EPSDenseNHEPSorted(n,array,eps->eigr,eps->eigi,pV,eps->ncv,eps->which);
 84: 
 85:   MatRestoreArray(la->BA,&array);

 87:   if (size == 1) {
 88:     VecRestoreArray(eps->V[0],&pV);
 89:   } else {
 90:     for (i=0; i<eps->ncv; i++) {
 91:       VecGetOwnershipRange(eps->V[i], &low, &high);
 92:       VecGetArray(eps->V[i], &array);
 93:       PetscMemcpy(array, pV+i*n+low, (high-low)*sizeof(PetscScalar));
 94:       VecRestoreArray(eps->V[i], &array);
 95:     }
 96:     PetscFree(pV);
 97:   }

 99:   eps->nconv = eps->ncv;
100:   eps->its   = 1;
101:   eps->reason = EPS_CONVERGED_TOL;
102: 
103: #ifndef PETSC_USE_COMPLEX
104:   if (eps->eigi[eps->nconv - 1] >= 0.0 && eps->nconv != n) eps->nconv--;
105: #endif

107:   return(0);
108: }

112: PetscErrorCode EPSDestroy_LAPACK(EPS eps)
113: {
114:   PetscErrorCode ierr;
115:   int            size;
116:   EPS_LAPACK     *la = (EPS_LAPACK *)eps->data;
117:   MPI_Comm       comm = eps->comm;

121:   MPI_Comm_size(comm,&size);
122:   if (la->BA) { MatDestroy(la->BA); }
123:   if (eps->data) {PetscFree(eps->data);}
124:   EPSFreeSolutionContiguous(eps);
125:   return(0);
126: }

128: EXTERN_C_BEGIN
131: PetscErrorCode EPSCreate_LAPACK(EPS eps)
132: {
133:   PetscErrorCode ierr;
134:   EPS_LAPACK     *la;

137:   PetscNew(EPS_LAPACK,&la);
138:   PetscMemzero(la,sizeof(EPS_LAPACK));
139:   PetscLogObjectMemory(eps,sizeof(EPS_LAPACK));
140:   eps->data                      = (void *) la;
141:   eps->ops->solve                = EPSSolve_LAPACK;
142:   eps->ops->setup                = EPSSetUp_LAPACK;
143:   eps->ops->destroy              = EPSDestroy_LAPACK;
144:   eps->ops->backtransform        = EPSBackTransform_Default;
145:   eps->ops->computevectors       = EPSComputeVectors_Default;
146:   return(0);
147: }
148: EXTERN_C_END