Actual source code: planso.c

  2: /*                       
  3:        This file implements a wrapper to the PLANSO package
  4: */
 5:  #include src/eps/impls/planso/plansop.h

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

 12: PetscErrorCode EPSSetUp_PLANSO(EPS eps)
 13: {
 15:   PetscInt       n;
 16:   EPS_PLANSO     *pl = (EPS_PLANSO *)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:     if (eps->ncv>n) SETERRQ(1,"The value of ncv cannot be larger than N");
 23:   } else eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),n);
 24:   eps->max_it = eps->ncv;
 25:   if (!eps->tol) eps->tol = 1.e-7;

 27:   if (!eps->ishermitian)
 28:     SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");

 30:   VecGetLocalSize(eps->vec_initial,&n);
 31:   pl->lwork = 5*n+1+4*eps->ncv+PetscMax(n,eps->ncv+1);
 32:   PetscFree(pl->work);
 33:   PetscMalloc(pl->lwork*sizeof(PetscReal),&pl->work);

 35:   EPSAllocateSolutionContiguous(eps);
 36:   return(0);
 37: }

 41: int PLANop_(int *n,PetscReal *s, PetscReal *q, PetscReal *p,MPI_Comm* c)
 42: {
 44:   Vec            x,y;

 47:   VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,(PetscScalar*)q,&x);
 48:   VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,(PetscScalar*)p,&y);
 49:   STApply(globaleps->OP,x,y);
 50:   EPSOrthogonalize(globaleps,globaleps->nds,globaleps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 51:   VecDestroy(x);
 52:   VecDestroy(y);
 53:   return(0);
 54: }

 58: int PLANopm_(int *n,PetscReal *q, PetscReal *s,MPI_Comm* c)
 59: {
 61:   Vec            x,y;

 64:   VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,(PetscScalar*)q,&x);
 65:   VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,(PetscScalar*)s,&y);
 66:   STApplyB(globaleps->OP,x,y);
 67:   VecDestroy(x);
 68:   VecDestroy(y);
 69:   return(0);
 70: }

 74: PetscErrorCode EPSSolve_PLANSO(EPS eps)
 75: {
 77:   PetscInt       nn;
 78:   int            i, n, msglvl, lohi,info;
 79:   PetscReal      condm;
 80:   EPS_PLANSO     *pl = (EPS_PLANSO *)eps->data;
 81:   MPI_Fint       fcomm;
 82: 

 85:   VecGetLocalSize(eps->vec_initial,&nn);
 86:   n = nn;
 87: 
 88:   if (eps->which==EPS_LARGEST_REAL) lohi = 1;
 89:   else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
 90:   else SETERRQ(1,"Wrong value of eps->which");

 92:   condm = 1.0;         /* estimated condition number: we have no information */
 93:   msglvl = 0;
 94:   globaleps = eps;
 95:   fcomm = MPI_Comm_c2f(eps->comm);

 97:   PLANdr2_( &n, &eps->ncv, &eps->nev, &lohi, &condm, &eps->tol, &eps->its, &eps->nconv,
 98:             eps->eigr, eps->eigi, pl->work, &pl->lwork, &info, &msglvl, &fcomm );

100:   for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;
101:   eps->reason = EPS_CONVERGED_TOL;
102: 
103:   if (info!=0) { SETERRQ1(PETSC_ERR_LIB,"Error in PLANSO (code=%d)",info);}

105:   return(0);
106: }

110: PetscErrorCode EPSDestroy_PLANSO(EPS eps)
111: {
113:   EPS_PLANSO     *pl = (EPS_PLANSO *)eps->data;

117:   PetscFree(pl->work);
118:   PetscFree(eps->data);
119:   EPSFreeSolutionContiguous(eps);
120:   return(0);
121: }

126: PetscErrorCode EPSCreate_PLANSO(EPS eps)
127: {
129:   EPS_PLANSO     *planso;

132:   PetscNew(EPS_PLANSO,&planso);
133:   PetscLogObjectMemory(eps,sizeof(EPS_PLANSO));
134:   eps->data                      = (void *) planso;
135:   eps->ops->solve                = EPSSolve_PLANSO;
136:   eps->ops->setup                = EPSSetUp_PLANSO;
137:   eps->ops->destroy              = EPSDestroy_PLANSO;
138:   eps->ops->backtransform        = EPSBackTransform_Default;
139:   eps->ops->computevectors       = EPSComputeVectors_Default;
140:   eps->which = EPS_LARGEST_REAL;
141:   return(0);
142: }