Actual source code: default.c

  1: /*
  2:      This file contains some simple default routines for common operations.  
  3: */
 4:  #include src/eps/epsimpl.h
 5:  #include slepcblaslapack.h

  9: PetscErrorCode EPSDestroy_Default(EPS eps)
 10: {
 11:   PetscErrorCode ierr;

 15:   if (eps->data) {PetscFree(eps->data);}

 17:   /* free work vectors */
 18:   EPSDefaultFreeWork(eps);
 19:   EPSFreeSolution(eps);
 20:   return(0);
 21: }

 25: PetscErrorCode EPSBackTransform_Default(EPS eps)
 26: {
 27:   PetscErrorCode ierr;
 28:   int            i;

 32:   for (i=0;i<eps->nconv;i++) {
 33:     STBackTransform(eps->OP,&eps->eigr[i],&eps->eigi[i]);
 34:   }
 35:   return(0);
 36: }

 40: /*
 41:   EPSComputeVectors_Default - Compute eigenvectors from the vectors
 42:   provided by the eigensolver. This version just copies the vectors
 43:   and is intended for solvers such as power that provide the eigenvector.
 44:  */
 45: PetscErrorCode EPSComputeVectors_Default(EPS eps)
 46: {
 47:   PetscErrorCode ierr;
 48:   int            i;

 51:   for (i=0;i<eps->nconv;i++) {
 52:     VecCopy(eps->V[i],eps->AV[i]);
 53:   }
 54:   eps->evecsavailable = PETSC_TRUE;
 55:   return(0);
 56: }

 60: /*
 61:   EPSComputeVectors_Schur - Compute eigenvectors from the vectors
 62:   provided by the eigensolver. This version is intended for solvers 
 63:   that provide Schur vectors. Given the partial Schur decomposition
 64:   OP*V=V*T, the following steps are performed:
 65:       1) compute eigenvectors of T: T*Y=Y*D
 66:       2) compute eigenvectors of OP: X=V*Y
 67:  */
 68: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
 69: {
 70:   PetscErrorCode ierr;
 71:   int            i,mout,info,ncv=eps->ncv;
 72:   PetscScalar    *Y,*work;
 73: #if defined(PETSC_USE_COMPLEX)
 74:   PetscReal      *rwork;
 75: #endif
 76: 
 78:   if (eps->ishermitian) {
 79:     EPSComputeVectors_Default(eps);
 80:     return(0);
 81:   }

 83: #if defined(PETSC_BLASLAPACK_ESSL_ONLY)
 84:   SETERRQ(PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
 85: #endif 

 87:   PetscMalloc(ncv*ncv*sizeof(PetscScalar),&Y);
 88:   PetscMalloc(3*ncv*sizeof(PetscScalar),&work);
 89: #if defined(PETSC_USE_COMPLEX)
 90:   PetscMalloc(ncv*sizeof(PetscReal),&rwork);
 91: #endif

 93: #if !defined(PETSC_USE_COMPLEX)
 94:   LAtrevc_("R","A",PETSC_NULL,&ncv,eps->T,&ncv,PETSC_NULL,&ncv,Y,&ncv,&ncv,&mout,work,&info,1,1);
 95: #else
 96:   LAtrevc_("R","A",PETSC_NULL,&ncv,eps->T,&ncv,PETSC_NULL,&ncv,Y,&ncv,&ncv,&mout,work,rwork,&info,1,1);
 97: #endif
 98:   if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREVC %i",info);

100:   for (i=0;i<eps->nconv;i++) {
101:     VecCopy(eps->V[i],eps->AV[i]);
102:   }
103:   EPSReverseProjection(eps,eps->AV,Y,0,ncv,eps->work);
104: 
105:   PetscFree(Y);
106:   PetscFree(work);
107: #if defined(PETSC_USE_COMPLEX)
108:   PetscFree(rwork);
109: #endif
110:   eps->evecsavailable = PETSC_TRUE;
111:   return(0);
112: }

116: /*
117:   EPSDefaultGetWork - Gets a number of work vectors.

119:   Input Parameters:
120: + eps  - eigensolver context
121: - nw   - number of work vectors to allocate

123:   Notes:
124:   Call this only if no work vectors have been allocated.

126:  */
127: PetscErrorCode EPSDefaultGetWork(EPS eps, int nw)
128: {
129:   PetscErrorCode ierr;


133:   if (eps->nwork != nw) {
134:     if (eps->nwork > 0) {
135:       VecDestroyVecs(eps->work,eps->nwork);
136:     }
137:     eps->nwork = nw;
138:     VecDuplicateVecs(eps->vec_initial,nw,&eps->work);
139:     PetscLogObjectParents(eps,nw,eps->work);
140:   }
141: 
142:   return(0);
143: }

147: /*
148:   EPSDefaultFreeWork - Free work vectors.

150:   Input Parameters:
151: . eps  - eigensolver context

153:  */
154: PetscErrorCode EPSDefaultFreeWork(EPS eps)
155: {
156:   PetscErrorCode ierr;

160:   if (eps->work)  {
161:     VecDestroyVecs(eps->work,eps->nwork);
162:   }
163:   return(0);
164: }