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: static int EPSSetUp_PLANSO(EPS eps)
 13: {
 14:   int        ierr, n;
 15:   EPS_PLANSO *pl = (EPS_PLANSO *)eps->data;

 18: #if defined(PETSC_USE_COMPLEX)
 19:   SETERRQ(PETSC_ERR_SUP,"Requested method is not available for complex problems");
 20: #endif
 21:   if (!eps->ishermitian)
 22:     SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");

 24:   VecGetLocalSize(eps->vec_initial,&n);
 25:   pl->lwork = 5*n+1+4*eps->ncv+PetscMax(n,eps->ncv+1);
 26:   PetscMalloc(pl->lwork*sizeof(PetscReal),&pl->work);

 28:   return(0);
 29: }

 33: static int EPSSetDefaults_PLANSO(EPS eps)
 34: {
 35:   int         ierr, N;

 38:   VecGetSize(eps->vec_initial,&N);
 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->max_it) eps->max_it = PetscMax(100,N);
 44:   if (!eps->tol) eps->tol = 1.e-7;
 45:   return(0);
 46: }

 50: int  PLANop_(int *n,PetscReal *s, PetscReal *q, PetscReal *p)
 51: {
 52:   Vec    x,y;
 53:   int    ierr;

 56:   VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,(PetscScalar*)q,&x);
 57:   VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,(PetscScalar*)p,&y);
 58:   STApply(globaleps->OP,x,y);
 59:   return(0);
 60: }

 64: int  PLANopm_(int *n,PetscReal *q, PetscReal *s)
 65: {
 66:   Vec    x,y;
 67:   int    ierr;

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

 78: static int  EPSSolve_PLANSO(EPS eps)
 79: {
 80:   int        i, n, msglvl, lohi, ierr;
 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, &ierr, &msglvl, &fcomm );

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

106:   return(0);
107: }

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

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

123: EXTERN_C_BEGIN
126: int EPSCreate_PLANSO(EPS eps)
127: {
128:   EPS_PLANSO *planso;
129:   int        ierr;

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