Actual source code: trlan.c

  2: /*                       
  3:        This file implements a wrapper to the TRLAN package
  4: */
 5:  #include src/eps/impls/trlan/trlanp.h

  7: /* Nasty global variable to access EPS data from TRLan_ */
  8: static EPS globaleps;

 12: PetscErrorCode EPSSetUp_TRLAN(EPS eps)
 13: {
 15:   PetscInt       n;
 16:   EPS_TRLAN      *tr = (EPS_TRLAN *)eps->data;

 19:   VecGetSize(eps->vec_initial,&n);
 20:   if (eps->ncv) {
 21:     if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
 22:   }
 23:   else eps->ncv = eps->nev;
 24:   if (!eps->max_it) eps->max_it = PetscMax(100,n);
 25:   if (!eps->tol) eps->tol = 1.e-7;
 26: 
 27:   if (!eps->ishermitian)
 28:     SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");

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

 33:   tr->restart = 0;
 34:   VecGetLocalSize(eps->vec_initial,&n);
 35:   tr->maxlan = eps->nev+PetscMin(eps->nev,6);
 36:   if (tr->maxlan+1-eps->ncv<=0) tr->lwork = tr->maxlan*(tr->maxlan+10);
 37:   else tr->lwork = n*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10);
 38:   PetscMalloc(tr->lwork*sizeof(PetscReal),&tr->work);

 40:   EPSAllocateSolutionContiguous(eps);
 41:   return(0);
 42: }

 46: static int MatMult_TRLAN(int *n,int *m,PetscReal *xin,int *ldx,PetscReal *yout,int *ldy)
 47: {
 49:   Vec            x,y;
 50:   int            i;

 53:   VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,PETSC_NULL,&x);
 54:   VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,PETSC_NULL,&y);
 55:   for (i=0;i<*m;i++) {
 56:     VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
 57:     VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
 58:     STApply(globaleps->OP,x,y);
 59:     EPSOrthogonalize(globaleps,globaleps->nds,globaleps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 60:     VecResetArray(x);
 61:     VecResetArray(y);
 62:   }
 63:   return(0);
 64: }

 68: PetscErrorCode EPSSolve_TRLAN(EPS eps)
 69: {
 71:   PetscInt       nn;
 72:   int            ipar[32], i, n, lohi, stat;
 73:   EPS_TRLAN      *tr = (EPS_TRLAN *)eps->data;
 74:   PetscScalar    *pV;
 75: 

 78:   VecGetLocalSize(eps->vec_initial,&nn);
 79:   n = nn;
 80: 
 81:   if (eps->which==EPS_LARGEST_REAL) lohi = 1;
 82:   else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
 83:   else SETERRQ(1,"Wrong value of eps->which");

 85:   globaleps = eps;

 87:   ipar[0]  = 0;            /* stat: error flag */
 88:   ipar[1]  = lohi;         /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
 89:   ipar[2]  = eps->nev;     /* number of desired eigenpairs */
 90:   ipar[3]  = 0;            /* number of eigenpairs already converged */
 91:   ipar[4]  = tr->maxlan;   /* maximum Lanczos basis size */
 92:   ipar[5]  = tr->restart;  /* restarting scheme */
 93:   ipar[6]  = eps->max_it;  /* maximum number of MATVECs */
 94:   ipar[7]  = MPI_Comm_c2f(eps->comm);    /* communicator */
 95:   ipar[8]  = 0;            /* verboseness */
 96:   ipar[9]  = 99;           /* Fortran IO unit number used to write log messages */
 97:   ipar[10] = 1;            /* use supplied starting vector */
 98:   ipar[11] = 0;            /* checkpointing flag */
 99:   ipar[12] = 98;           /* Fortran IO unit number used to write checkpoint files */
100:   ipar[13] = 0;            /* number of flops per matvec per PE (not used) */
101:   tr->work[0] = eps->tol;  /* relative tolerance on residual norms */

103:   for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
104:   VecCopy(eps->vec_initial,eps->V[0]);
105:   VecGetArray(eps->V[0],&pV);

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

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

111:   stat        = ipar[0];
112:   eps->nconv  = ipar[3];
113:   eps->its    = ipar[25];
114:   eps->reason = EPS_CONVERGED_TOL;
115:   for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;
116: 
117:   if (stat!=0) { SETERRQ1(PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);}

119:   return(0);
120: }

124: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
125: {
127:   EPS_TRLAN      *tr = (EPS_TRLAN *)eps->data;

131:   PetscFree(tr->work);
132:   PetscFree(eps->data);
133:   EPSFreeSolutionContiguous(eps);
134:   return(0);
135: }

140: PetscErrorCode EPSCreate_TRLAN(EPS eps)
141: {
143:   EPS_TRLAN      *trlan;

146:   PetscNew(EPS_TRLAN,&trlan);
147:   PetscLogObjectMemory(eps,sizeof(EPS_TRLAN));
148:   eps->data                      = (void *) trlan;
149:   eps->ops->solve                = EPSSolve_TRLAN;
150:   eps->ops->setup                = EPSSetUp_TRLAN;
151:   eps->ops->destroy              = EPSDestroy_TRLAN;
152:   eps->ops->backtransform        = EPSBackTransform_Default;
153:   eps->ops->computevectors       = EPSComputeVectors_Default;
154:   eps->which = EPS_LARGEST_REAL;
155:   return(0);
156: }