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-2010, Universidad 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 src/eps/impls/external/trlan/trlanp.h

 26: PetscErrorCode EPSSolve_TRLAN(EPS);

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

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

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

 49:   if (eps->isgeneralized)
 50:     SETERRQ(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)
 54:     SETERRQ(1,"Wrong value of eps->which");

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

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

 66:   EPSAllocateSolution(eps);

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

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

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

100: PetscErrorCode EPSSolve_TRLAN(EPS eps)
101: {
103:   PetscInt       i;
104:   PetscBLASInt   ipar[32], n, lohi, stat, ncv;
105:   EPS_TRLAN      *tr = (EPS_TRLAN *)eps->data;
106:   PetscScalar    *pV;
107: 

110:   ncv = PetscBLASIntCast(eps->ncv);
111:   n = PetscBLASIntCast(eps->nloc);
112: 
113:   if (eps->which==EPS_LARGEST_REAL) lohi = 1;
114:   else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
115:   else SETERRQ(1,"Wrong value of eps->which");

117:   globaleps = eps;

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

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

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

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

143:   stat        = ipar[0];
144:   eps->nconv  = ipar[3];
145:   eps->its    = ipar[25];
146:   eps->reason = EPS_CONVERGED_TOL;
147: 
148:   if (stat!=0) { SETERRQ1(PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);}

150:   return(0);
151: }

155: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
156: {
158:   EPS_TRLAN      *tr = (EPS_TRLAN *)eps->data;

162:   PetscFree(tr->work);
163:   PetscFree(eps->data);
164:   EPSFreeSolution(eps);
165:   return(0);
166: }

171: PetscErrorCode EPSCreate_TRLAN(EPS eps)
172: {
174:   EPS_TRLAN      *trlan;

177:   PetscNew(EPS_TRLAN,&trlan);
178:   PetscLogObjectMemory(eps,sizeof(EPS_TRLAN));
179:   eps->data                      = (void *) trlan;
180:   eps->ops->setup                = EPSSetUp_TRLAN;
181:   eps->ops->destroy              = EPSDestroy_TRLAN;
182:   eps->ops->backtransform        = EPSBackTransform_Default;
183:   eps->ops->computevectors       = EPSComputeVectors_Default;
184:   return(0);
185: }