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: {
14: PetscErrorCode ierr;
15: int 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 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: if (eps->isgeneralized)
34: SETERRQ(PETSC_ERR_SUP,"Requested method is not available for generalized problems");
36: tr->restart = 0;
37: VecGetLocalSize(eps->vec_initial,&n);
38: tr->maxlan = eps->nev+PetscMin(eps->nev,6);
39: if (tr->maxlan+1-eps->ncv<=0) tr->lwork = tr->maxlan*(tr->maxlan+10);
40: else tr->lwork = n*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10);
41: PetscMalloc(tr->lwork*sizeof(PetscReal),&tr->work);
43: EPSAllocateSolutionContiguous(eps);
44: return(0);
45: }
49: static int MatMult_TRLAN(int *n,int *m,PetscReal *xin,int *ldx,PetscReal *yout,int *ldy)
50: {
51: Vec x,y;
52: int i,ierr;
55: VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,PETSC_NULL,&x);
56: VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,PETSC_NULL,&y);
57: for (i=0;i<*m;i++) {
58: VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
59: VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
60: STApply(globaleps->OP,x,y);
61: EPSOrthogonalize(globaleps,globaleps->nds,globaleps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
62: }
63: return(0);
64: }
68: PetscErrorCode EPSSolve_TRLAN(EPS eps)
69: {
70: int ipar[32], i, n, lohi, stat, ierr;
71: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
72: PetscScalar *pV;
73:
76: VecGetLocalSize(eps->vec_initial,&n);
77:
78: if (eps->which==EPS_LARGEST_REAL) lohi = 1;
79: else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
80: else SETERRQ(1,"Wrong value of eps->which");
82: globaleps = eps;
84: ipar[0] = 0; /* stat: error flag */
85: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
86: ipar[2] = eps->nev; /* number of desired eigenpairs */
87: ipar[3] = 0; /* number of eigenpairs already converged */
88: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
89: ipar[5] = tr->restart; /* restarting scheme */
90: ipar[6] = eps->max_it; /* maximum number of MATVECs */
91: ipar[7] = MPI_Comm_c2f(eps->comm); /* communicator */
92: ipar[8] = 0; /* verboseness */
93: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
94: ipar[10] = 1; /* use supplied starting vector */
95: ipar[11] = 0; /* checkpointing flag */
96: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
97: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
98: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
100: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
101: VecCopy(eps->vec_initial,eps->V[0]);
102: VecGetArray(eps->V[0],&pV);
104: TRLan_ ( MatMult_TRLAN, ipar, &n, &eps->ncv, eps->eigr, pV, &n, tr->work, &tr->lwork );
106: VecRestoreArray( eps->V[0], &pV );
108: stat = ipar[0];
109: eps->nconv = ipar[3];
110: eps->its = ipar[25];
111: eps->reason = EPS_CONVERGED_TOL;
112: for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;
113:
114: if (stat!=0) { SETERRQ1(PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);}
116: return(0);
117: }
121: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
122: {
123: PetscErrorCode ierr;
124: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
128: if (tr->work) { PetscFree(tr->work); }
129: if (eps->data) {PetscFree(eps->data);}
130: EPSFreeSolutionContiguous(eps);
131: return(0);
132: }
134: EXTERN_C_BEGIN
137: PetscErrorCode EPSCreate_TRLAN(EPS eps)
138: {
139: PetscErrorCode ierr;
140: EPS_TRLAN *trlan;
143: PetscNew(EPS_TRLAN,&trlan);
144: PetscMemzero(trlan,sizeof(EPS_TRLAN));
145: PetscLogObjectMemory(eps,sizeof(EPS_TRLAN));
146: eps->data = (void *) trlan;
147: eps->ops->solve = EPSSolve_TRLAN;
148: eps->ops->setup = EPSSetUp_TRLAN;
149: eps->ops->destroy = EPSDestroy_TRLAN;
150: eps->ops->backtransform = EPSBackTransform_Default;
151: eps->ops->computevectors = EPSComputeVectors_Default;
152: eps->which = EPS_LARGEST_REAL;
153: return(0);
154: }
155: EXTERN_C_END