Actual source code: trlan.c

  1: /*
  2:        This file implements a wrapper to the TRLAN package

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.
  9:       
 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc-private/epsimpl.h>          /*I "slepceps.h" I*/
 25: #include <../src/eps/impls/external/trlan/trlanp.h>

 27: PetscErrorCode EPSSolve_TRLAN(EPS);

 29: /* Nasty global variable to access EPS data from TRLan_ */
 30: static EPS globaleps;

 34: PetscErrorCode EPSSetUp_TRLAN(EPS eps)
 35: {
 37:   EPS_TRLAN      *tr = (EPS_TRLAN *)eps->data;

 40:   tr->maxlan = PetscBLASIntCast(PetscMax(7,eps->nev+PetscMin(eps->nev,6)));
 41:   if (eps->ncv) {
 42:     if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
 43:   }
 44:   else eps->ncv = tr->maxlan;
 45:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 46:   if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
 47: 
 48:   if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");

 50:   if (eps->isgeneralized) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Requested method is not available for generalized problems");

 52:   if (!eps->which) eps->which = EPS_LARGEST_REAL;
 53:   if (eps->which!=EPS_LARGEST_REAL && eps->which!=EPS_SMALLEST_REAL && eps->which!=EPS_TARGET_REAL) SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
 54:   if (eps->arbit_func) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");

 56:   tr->restart = 0;
 57:   if (tr->maxlan+1-eps->ncv<=0) { tr->lwork = PetscBLASIntCast(tr->maxlan*(tr->maxlan+10)); }
 58:   else { tr->lwork = PetscBLASIntCast(eps->nloc*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10)); }
 59:   PetscMalloc(tr->lwork*sizeof(PetscReal),&tr->work);

 61:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }

 63:   EPSAllocateSolution(eps);

 65:   /* dispatch solve method */
 66:   if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
 67:   eps->ops->solve = EPSSolve_TRLAN;
 68:   return(0);
 69: }

 73: static PetscBLASInt MatMult_TRLAN(PetscBLASInt *n,PetscBLASInt *m,PetscReal *xin,PetscBLASInt *ldx,PetscReal *yout,PetscBLASInt *ldy)
 74: {
 76:   Vec            x,y;
 77:   PetscBLASInt   i;

 80:   VecCreateMPIWithArray(((PetscObject)globaleps)->comm,1,*n,PETSC_DECIDE,PETSC_NULL,&x);
 81:   VecCreateMPIWithArray(((PetscObject)globaleps)->comm,1,*n,PETSC_DECIDE,PETSC_NULL,&y);
 82:   for (i=0;i<*m;i++) {
 83:     VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
 84:     VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
 85:     STApply(globaleps->OP,x,y);
 86:     IPOrthogonalize(globaleps->ip,0,PETSC_NULL,globaleps->nds,PETSC_NULL,globaleps->defl,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 87:     VecResetArray(x);
 88:     VecResetArray(y);
 89:   }
 90:   VecDestroy(&x);
 91:   VecDestroy(&y);
 92:   return(0);
 93: }

 97: PetscErrorCode EPSSolve_TRLAN(EPS eps)
 98: {
100:   PetscInt       i;
101:   PetscBLASInt   ipar[32],n,lohi,stat,ncv;
102:   EPS_TRLAN      *tr = (EPS_TRLAN *)eps->data;
103:   PetscScalar    *pV;
104: 
106:   ncv = PetscBLASIntCast(eps->ncv);
107:   n = PetscBLASIntCast(eps->nloc);
108: 
109:   if (eps->which==EPS_LARGEST_REAL || eps->which==EPS_TARGET_REAL) lohi = 1;
110:   else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
111:   else SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");

113:   globaleps = eps;

115:   ipar[0]  = 0;            /* stat: error flag */
116:   ipar[1]  = lohi;         /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
117:   ipar[2]  = PetscBLASIntCast(eps->nev); /* number of desired eigenpairs */
118:   ipar[3]  = 0;            /* number of eigenpairs already converged */
119:   ipar[4]  = tr->maxlan;   /* maximum Lanczos basis size */
120:   ipar[5]  = tr->restart;  /* restarting scheme */
121:   ipar[6]  = PetscBLASIntCast(eps->max_it); /* maximum number of MATVECs */
122:   ipar[7]  = PetscBLASIntCast(MPI_Comm_c2f(((PetscObject)eps)->comm)); /* communicator */
123:   ipar[8]  = 0;            /* verboseness */
124:   ipar[9]  = 99;           /* Fortran IO unit number used to write log messages */
125:   ipar[10] = 1;            /* use supplied starting vector */
126:   ipar[11] = 0;            /* checkpointing flag */
127:   ipar[12] = 98;           /* Fortran IO unit number used to write checkpoint files */
128:   ipar[13] = 0;            /* number of flops per matvec per PE (not used) */
129:   tr->work[0] = eps->tol;  /* relative tolerance on residual norms */

131:   for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
132:   EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
133:   VecGetArray(eps->V[0],&pV);

135:   TRLan_(MatMult_TRLAN,ipar,&n,&ncv,eps->eigr,pV,&n,tr->work,&tr->lwork);

137:   VecRestoreArray(eps->V[0],&pV);

139:   stat        = ipar[0];
140:   eps->nconv  = ipar[3];
141:   eps->its    = ipar[25];
142:   eps->reason = EPS_CONVERGED_TOL;
143: 
144:   if (stat!=0) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);
145:   return(0);
146: }

150: PetscErrorCode EPSReset_TRLAN(EPS eps)
151: {
153:   EPS_TRLAN      *tr = (EPS_TRLAN *)eps->data;

156:   PetscFree(tr->work);
157:   EPSFreeSolution(eps);
158:   return(0);
159: }

163: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
164: {

168:   PetscFree(eps->data);
169:   return(0);
170: }

172: EXTERN_C_BEGIN
175: PetscErrorCode EPSCreate_TRLAN(EPS eps)
176: {

180:   PetscNewLog(eps,EPS_TRLAN,&eps->data);
181:   eps->ops->setup                = EPSSetUp_TRLAN;
182:   eps->ops->destroy              = EPSDestroy_TRLAN;
183:   eps->ops->reset                = EPSReset_TRLAN;
184:   eps->ops->backtransform        = EPSBackTransform_Default;
185:   eps->ops->computevectors       = EPSComputeVectors_Default;
186:   return(0);
187: }
188: EXTERN_C_END