Actual source code: trlan.c
1: /*
2: This file implements a wrapper to the TRLAN package
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include src/eps/impls/trlan/trlanp.h
15: /* Nasty global variable to access EPS data from TRLan_ */
16: static EPS globaleps;
20: PetscErrorCode EPSSetUp_TRLAN(EPS eps)
21: {
23: PetscInt n;
24: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
27: VecGetSize(eps->vec_initial,&n);
28: if (eps->ncv) {
29: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
30: }
31: else eps->ncv = eps->nev;
32: if (!eps->max_it) eps->max_it = PetscMax(1000,n);
33:
34: if (!eps->ishermitian)
35: SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
37: if (eps->isgeneralized)
38: SETERRQ(PETSC_ERR_SUP,"Requested method is not available for generalized problems");
40: tr->restart = 0;
41: VecGetLocalSize(eps->vec_initial,&n);
42: tr->maxlan = eps->nev+PetscMin(eps->nev,6);
43: if (tr->maxlan+1-eps->ncv<=0) tr->lwork = tr->maxlan*(tr->maxlan+10);
44: else tr->lwork = n*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10);
45: PetscMalloc(tr->lwork*sizeof(PetscReal),&tr->work);
47: EPSAllocateSolutionContiguous(eps);
48: EPSDefaultGetWork(eps,1);
49: return(0);
50: }
54: static int MatMult_TRLAN(int *n,int *m,PetscReal *xin,int *ldx,PetscReal *yout,int *ldy)
55: {
57: Vec x,y;
58: int i;
61: VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,PETSC_NULL,&x);
62: VecCreateMPIWithArray(globaleps->comm,*n,PETSC_DECIDE,PETSC_NULL,&y);
63: for (i=0;i<*m;i++) {
64: VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
65: VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
66: STApply(globaleps->OP,x,y);
67: IPOrthogonalize(globaleps->ip,globaleps->nds,PETSC_NULL,globaleps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL,globaleps->work[0]);
68: VecResetArray(x);
69: VecResetArray(y);
70: }
71: VecDestroy(x);
72: VecDestroy(y);
73: return(0);
74: }
78: PetscErrorCode EPSSolve_TRLAN(EPS eps)
79: {
81: PetscInt nn;
82: int ipar[32], i, n, lohi, stat;
83: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
84: PetscScalar *pV;
85:
88: VecGetLocalSize(eps->vec_initial,&nn);
89: n = nn;
90:
91: if (eps->which==EPS_LARGEST_REAL) lohi = 1;
92: else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
93: else SETERRQ(1,"Wrong value of eps->which");
95: globaleps = eps;
97: ipar[0] = 0; /* stat: error flag */
98: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
99: ipar[2] = eps->nev; /* number of desired eigenpairs */
100: ipar[3] = 0; /* number of eigenpairs already converged */
101: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
102: ipar[5] = tr->restart; /* restarting scheme */
103: ipar[6] = eps->max_it; /* maximum number of MATVECs */
104: ipar[7] = MPI_Comm_c2f(eps->comm); /* communicator */
105: ipar[8] = 0; /* verboseness */
106: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
107: ipar[10] = 1; /* use supplied starting vector */
108: ipar[11] = 0; /* checkpointing flag */
109: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
110: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
111: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
113: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
114: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
115: VecGetArray(eps->V[0],&pV);
117: TRLan_ ( MatMult_TRLAN, ipar, &n, &eps->ncv, eps->eigr, pV, &n, tr->work, &tr->lwork );
119: VecRestoreArray( eps->V[0], &pV );
121: stat = ipar[0];
122: eps->nconv = ipar[3];
123: eps->its = ipar[25];
124: eps->reason = EPS_CONVERGED_TOL;
125:
126: if (stat!=0) { SETERRQ1(PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);}
128: return(0);
129: }
133: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
134: {
136: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
140: PetscFree(tr->work);
141: PetscFree(eps->data);
142: EPSFreeSolutionContiguous(eps);
143: return(0);
144: }
149: PetscErrorCode EPSCreate_TRLAN(EPS eps)
150: {
152: EPS_TRLAN *trlan;
155: PetscNew(EPS_TRLAN,&trlan);
156: PetscLogObjectMemory(eps,sizeof(EPS_TRLAN));
157: eps->data = (void *) trlan;
158: eps->ops->solve = EPSSolve_TRLAN;
159: eps->ops->setup = EPSSetUp_TRLAN;
160: eps->ops->destroy = EPSDestroy_TRLAN;
161: eps->ops->backtransform = EPSBackTransform_Default;
162: eps->ops->computevectors = EPSComputeVectors_Default;
163: eps->which = EPS_LARGEST_REAL;
164: return(0);
165: }