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: static int EPSSetUp_TRLAN(EPS eps)
13: {
14: int ierr, n;
15: EPS_TRLAN *tr = (EPS_TRLAN *)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: if (eps->isgeneralized)
25: SETERRQ(PETSC_ERR_SUP,"Requested method is not available for generalized problems");
27: tr->restart = 0;
28: VecGetLocalSize(eps->vec_initial,&n);
29: tr->maxlan = eps->nev+PetscMin(eps->nev,6);
30: if (tr->maxlan+1-eps->ncv<=0) tr->lwork = tr->maxlan*(tr->maxlan+10);
31: else tr->lwork = n*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10);
32: PetscMalloc(tr->lwork*sizeof(PetscReal),&tr->work);
34: return(0);
35: }
39: static int EPSSetDefaults_TRLAN(EPS eps)
40: {
41: int ierr, N;
44: VecGetSize(eps->vec_initial,&N);
45: if (eps->ncv) {
46: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
47: }
48: else eps->ncv = eps->nev;
49: if (!eps->max_it) eps->max_it = PetscMax(100,N);
50: if (!eps->tol) eps->tol = 1.e-7;
51: return(0);
52: }
56: static int MatMult_TRLAN(int *n,int *m,PetscReal *xin,int *ldx,PetscReal *yout,int *ldy)
57: {
58: Vec x,y;
59: int i,ierr;
62: VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,PETSC_NULL,&x);
63: VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,PETSC_NULL,&y);
64: for (i=0;i<*m;i++) {
65: VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
66: VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
67: STApply(globaleps->OP,x,y);
68: }
69: return(0);
70: }
74: static int EPSSolve_TRLAN(EPS eps)
75: {
76: int ipar[32], i, n, lohi, stat, ierr;
77: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
78: PetscScalar *pV;
79:
82: VecGetLocalSize(eps->vec_initial,&n);
83:
84: if (eps->which==EPS_LARGEST_REAL) lohi = 1;
85: else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
86: else SETERRQ(1,"Wrong value of eps->which");
88: globaleps = eps;
90: ipar[0] = 0; /* stat: error flag */
91: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
92: ipar[2] = eps->nev; /* number of desired eigenpairs */
93: ipar[3] = 0; /* number of eigenpairs already converged */
94: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
95: ipar[5] = tr->restart; /* restarting scheme */
96: ipar[6] = eps->max_it; /* maximum number of MATVECs */
97: ipar[7] = MPI_Comm_c2f(eps->comm); /* communicator */
98: ipar[8] = 0; /* verboseness */
99: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
100: ipar[10] = 1; /* use supplied starting vector */
101: ipar[11] = 0; /* checkpointing flag */
102: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
103: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
104: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
106: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
107: VecCopy(eps->vec_initial,eps->V[0]);
108: VecGetArray(eps->V[0],&pV);
110: TRLan_ ( MatMult_TRLAN, ipar, &n, &eps->ncv, eps->eigr, pV, &n, tr->work, &tr->lwork );
112: VecRestoreArray( eps->V[0], &pV );
114: stat = ipar[0];
115: eps->nconv = ipar[3];
116: eps->its = ipar[25];
117: eps->reason = EPS_CONVERGED_TOL;
118: for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;
119:
120: if (stat!=0) { SETERRQ1(PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);}
122: return(0);
123: }
127: int EPSDestroy_TRLAN(EPS eps)
128: {
129: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
130: int ierr;
134: if (tr->work) { PetscFree(tr->work); }
135: if (eps->data) { PetscFree(eps->data); }
136: return(0);
137: }
139: EXTERN_C_BEGIN
142: int EPSCreate_TRLAN(EPS eps)
143: {
144: EPS_TRLAN *trlan;
145: int ierr;
148: PetscNew(EPS_TRLAN,&trlan);
149: PetscMemzero(trlan,sizeof(EPS_TRLAN));
150: PetscLogObjectMemory(eps,sizeof(EPS_TRLAN));
151: eps->data = (void *) trlan;
152: eps->ops->setup = EPSSetUp_TRLAN;
153: eps->ops->setdefaults = EPSSetDefaults_TRLAN;
154: eps->ops->solve = EPSSolve_TRLAN;
155: eps->ops->destroy = EPSDestroy_TRLAN;
156: eps->ops->view = 0;
157: eps->ops->backtransform = EPSBackTransform_Default;
158: eps->which = EPS_LARGEST_REAL;
159: return(0);
160: }
161: EXTERN_C_END