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-2011, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <private/epsimpl.h> /*I "slepceps.h" I*/
25: #include <../src/eps/impls/external/trlan/trlanp.h>
27: PetscErrorCode EPSSolve_TRLAN(EPS);
29: /* Nasty global variable to access EPS data from TRLan_ */
30: static EPS globaleps;
34: PetscErrorCode EPSSetUp_TRLAN(EPS eps)
35: {
37: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
40: tr->maxlan = PetscBLASIntCast(PetscMax(7,eps->nev+PetscMin(eps->nev,6)));
41: if (eps->ncv) {
42: if (eps->ncv<eps->nev) SETERRQ(((PetscObject)eps)->comm,1,"The value of ncv must be at least nev");
43: }
44: else eps->ncv = tr->maxlan;
45: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
46: if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
47:
48: if (!eps->ishermitian)
49: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
51: if (eps->isgeneralized)
52: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Requested method is not available for generalized problems");
54: if (!eps->which) eps->which = EPS_LARGEST_REAL;
55: if (eps->which!=EPS_LARGEST_REAL && eps->which!=EPS_SMALLEST_REAL && eps->which!=EPS_TARGET_REAL)
56: SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
58: tr->restart = 0;
59: if (tr->maxlan+1-eps->ncv<=0) { tr->lwork = PetscBLASIntCast(tr->maxlan*(tr->maxlan+10)); }
60: else { tr->lwork = PetscBLASIntCast(eps->nloc*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10)); }
61: PetscMalloc(tr->lwork*sizeof(PetscReal),&tr->work);
63: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
65: EPSAllocateSolution(eps);
67: /* dispatch solve method */
68: if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
69: eps->ops->solve = EPSSolve_TRLAN;
70: return(0);
71: }
75: static PetscBLASInt MatMult_TRLAN(PetscBLASInt *n,PetscBLASInt *m,PetscReal *xin,PetscBLASInt *ldx,PetscReal *yout,PetscBLASInt *ldy)
76: {
78: Vec x,y;
79: PetscBLASInt i;
82: VecCreateMPIWithArray(((PetscObject)globaleps)->comm,*n,PETSC_DECIDE,PETSC_NULL,&x);
83: VecCreateMPIWithArray(((PetscObject)globaleps)->comm,*n,PETSC_DECIDE,PETSC_NULL,&y);
84: for (i=0;i<*m;i++) {
85: VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
86: VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
87: STApply(globaleps->OP,x,y);
88: IPOrthogonalize(globaleps->ip,0,PETSC_NULL,globaleps->nds,PETSC_NULL,globaleps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
89: VecResetArray(x);
90: VecResetArray(y);
91: }
92: VecDestroy(&x);
93: VecDestroy(&y);
94: return(0);
95: }
99: PetscErrorCode EPSSolve_TRLAN(EPS eps)
100: {
102: PetscInt i;
103: PetscBLASInt ipar[32],n,lohi,stat,ncv;
104: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
105: PetscScalar *pV;
106:
108: ncv = PetscBLASIntCast(eps->ncv);
109: n = PetscBLASIntCast(eps->nloc);
110:
111: if (eps->which==EPS_LARGEST_REAL || eps->which==EPS_TARGET_REAL) lohi = 1;
112: else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
113: else SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
115: globaleps = eps;
117: ipar[0] = 0; /* stat: error flag */
118: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
119: ipar[2] = PetscBLASIntCast(eps->nev); /* number of desired eigenpairs */
120: ipar[3] = 0; /* number of eigenpairs already converged */
121: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
122: ipar[5] = tr->restart; /* restarting scheme */
123: ipar[6] = PetscBLASIntCast(eps->max_it); /* maximum number of MATVECs */
124: ipar[7] = PetscBLASIntCast(MPI_Comm_c2f(((PetscObject)eps)->comm)); /* communicator */
125: ipar[8] = 0; /* verboseness */
126: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
127: ipar[10] = 1; /* use supplied starting vector */
128: ipar[11] = 0; /* checkpointing flag */
129: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
130: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
131: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
133: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
134: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
135: VecGetArray(eps->V[0],&pV);
137: TRLan_(MatMult_TRLAN,ipar,&n,&ncv,eps->eigr,pV,&n,tr->work,&tr->lwork);
139: VecRestoreArray(eps->V[0],&pV);
141: stat = ipar[0];
142: eps->nconv = ipar[3];
143: eps->its = ipar[25];
144: eps->reason = EPS_CONVERGED_TOL;
145:
146: if (stat!=0) { SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);}
147: return(0);
148: }
152: PetscErrorCode EPSReset_TRLAN(EPS eps)
153: {
155: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
158: PetscFree(tr->work);
159: EPSFreeSolution(eps);
160: return(0);
161: }
165: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
166: {
170: PetscFree(eps->data);
171: return(0);
172: }
177: PetscErrorCode EPSCreate_TRLAN(EPS eps)
178: {
182: PetscNewLog(eps,EPS_TRLAN,&eps->data);
183: eps->ops->setup = EPSSetUp_TRLAN;
184: eps->ops->destroy = EPSDestroy_TRLAN;
185: eps->ops->reset = EPSReset_TRLAN;
186: eps->ops->backtransform = EPSBackTransform_Default;
187: eps->ops->computevectors = EPSComputeVectors_Default;
188: return(0);
189: }