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: {
 14:   PetscErrorCode ierr;
 15:   int            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:   }
 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;

 27: #if defined(PETSC_USE_COMPLEX)
 28:   SETERRQ(PETSC_ERR_SUP,"Requested method is not available for complex problems");
 29: #endif
 30:   if (!eps->ishermitian)
 31:     SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");

 33:   VecGetLocalSize(eps->vec_initial,&n);
 34:   pl->lwork = 5*n+1+4*eps->ncv+PetscMax(n,eps->ncv+1);
 35:   if (pl->work)  { PetscFree(pl->work); }
 36:   PetscMalloc(pl->lwork*sizeof(PetscReal),&pl->work);

 38:   EPSAllocateSolutionContiguous(eps);
 39:   return(0);
 40: }

 44: int PLANop_(int *n,PetscReal *s, PetscReal *q, PetscReal *p)
 45: {
 46:   PetscErrorCode ierr;
 47:   Vec            x,y;

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

 61: int PLANopm_(int *n,PetscReal *q, PetscReal *s)
 62: {
 63:   PetscErrorCode ierr;
 64:   Vec            x,y;

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

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

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

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

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

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

106:   return(0);
107: }

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

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

124: EXTERN_C_BEGIN
127: PetscErrorCode EPSCreate_PLANSO(EPS eps)
128: {
129:   PetscErrorCode ierr;
130:   EPS_PLANSO     *planso;

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