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: static int EPSSetUp_LAPACK(EPS eps)
 13: {
 14:   int         ierr,size,rank,n;
 15:   EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
 16:   MPI_Comm    comm = eps->comm;
 17:   Mat         *T;
 18:   IS isrow, iscol;
 19:   PetscTruth flg;

 22:   EPSComputeExplicitOperator(eps,&la->BA);
 23:   MPI_Comm_size(comm,&size);
 24:   MPI_Comm_rank(comm,&rank);
 25: 
 26:   if (size > 1) {
 27:     /* assemble full matrix on every processor */
 28:     MatGetSize(la->BA,&n,&n);
 29:     ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &isrow);
 30:     ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &iscol);
 31:     MatGetSubMatrices(la->BA, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &T);
 32:     ISDestroy(isrow);
 33:     ISDestroy(iscol);
 34:     MatDestroy(la->BA);
 35:     la->BA = *T;
 36:   }
 37: 
 38:   PetscTypeCompare((PetscObject)la->BA, MATSEQDENSE, &flg);
 39:   if (!flg) {
 40:     /* convert matrix to MatSeqDense */
 41:     MatConvert(la->BA, MATSEQDENSE, &la->BA);
 42:   }

 44:   return(0);
 45: }

 49: static int EPSSetDefaults_LAPACK(EPS eps)
 50: {
 51:   int         ierr, N;

 54:   VecGetSize(eps->vec_initial,&N);
 55:   if (eps->nev<1 || eps->nev>N) SETERRQ(1,"Wrong value of nev");
 56:   if (eps->ncv) {
 57:     if (eps->ncv<=eps->nev) SETERRQ(1,"Wrong value of ncv");
 58: #ifndef PETSC_USE_COMPLEX
 59:   } else eps->ncv = PetscMin(eps->nev + 1, N);
 60: #else
 61:   } else eps->ncv = eps->nev;
 62: #endif

 64:   return(0);
 65: }

 69: static int  EPSSolve_LAPACK(EPS eps)
 70: {
 71:   int         ierr,n,size,rank,i,low,high;
 72:   PetscScalar *array,*pV;
 73:   EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
 74:   MPI_Comm    comm = eps->comm;
 75: 
 77:   MPI_Comm_size(comm,&size);
 78:   MPI_Comm_rank(comm,&rank);
 79: 
 80:   MatGetArray(la->BA,&array);
 81:   MatGetSize(la->BA,&n,&n);

 83:   if (!eps->dropvectors) {
 84:     if (size == 1) {
 85:       VecGetArray(eps->V[0],&pV);
 86:     } else {
 87:       PetscMalloc(sizeof(PetscScalar)*n,&pV);
 88:     }
 89:   } else pV = PETSC_NULL;
 90: 
 91:   EPSDenseNHEPSorted(n,array,eps->eigr,eps->eigi,pV,eps->ncv,eps->which);
 92: 
 93:   MatRestoreArray(la->BA,&array);

 95:   if (!eps->dropvectors) {
 96:     if (size == 1) {
 97:       VecRestoreArray(eps->V[0],&pV);
 98:     } else {
 99:       for (i=0; i<eps->ncv; i++) {
100:         VecGetOwnershipRange(eps->V[i], &low, &high);
101:         VecGetArray(eps->V[i], &array);
102:         PetscMemcpy(array, pV+i*n+low, (high-low)*sizeof(PetscScalar));
103:         VecRestoreArray(eps->V[i], &array);
104:       }
105:       PetscFree(pV);
106:     }
107:   }

109:   eps->nconv = eps->ncv;
110:   eps->its   = 1;
111:   eps->reason = EPS_CONVERGED_TOL;
112: 
113: #ifndef PETSC_USE_COMPLEX
114:   if (eps->eigi[eps->nconv - 1] >= 0.0) eps->nconv--;
115: #endif

117:   return(0);
118: }

122: int EPSDestroy_LAPACK(EPS eps)
123: {
124:   int         ierr,size;
125:   EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
126:   MPI_Comm    comm = eps->comm;

130:   MPI_Comm_size(comm,&size);
131:   if (la->BA) { MatDestroy(la->BA); }
132:   EPSDefaultDestroy(eps);
133:   return(0);
134: }

136: EXTERN_C_BEGIN
139: int EPSCreate_LAPACK(EPS eps)
140: {
141:   EPS_LAPACK *la;

145:   PetscNew(EPS_LAPACK,&la);
146:   PetscMemzero(la,sizeof(EPS_LAPACK));
147:   PetscLogObjectMemory(eps,sizeof(EPS_LAPACK));
148:   eps->data                      = (void *) la;
149:   eps->ops->setup                = EPSSetUp_LAPACK;
150:   eps->ops->setdefaults          = EPSSetDefaults_LAPACK;
151:   eps->ops->solve                = EPSSolve_LAPACK;
152:   eps->ops->destroy              = EPSDestroy_LAPACK;
153:   eps->ops->view                 = 0;
154:   eps->ops->backtransform        = EPSBackTransform_Default;
155:   return(0);
156: }
157: EXTERN_C_END