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-2013, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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 <slepc-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: PetscBLASIntCast(PetscMax(7,eps->nev+PetscMin(eps->nev,6)),&tr->maxlan);
41: if (eps->ncv) {
42: if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
43: } else eps->ncv = tr->maxlan;
44: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
45: if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
47: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
49: if (eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is not available for generalized problems");
51: if (!eps->which) eps->which = EPS_LARGEST_REAL;
52: if (eps->which!=EPS_LARGEST_REAL && eps->which!=EPS_SMALLEST_REAL && eps->which!=EPS_TARGET_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
53: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
55: tr->restart = 0;
56: if (tr->maxlan+1-eps->ncv<=0) {
57: PetscBLASIntCast(tr->maxlan*(tr->maxlan+10),&tr->lwork);
58: } else {
59: PetscBLASIntCast(eps->nloc*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10),&tr->lwork);
60: }
61: PetscMalloc(tr->lwork*sizeof(PetscReal),&tr->work);
62: PetscLogObjectMemory(eps,tr->lwork*sizeof(PetscReal));
64: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
66: EPSAllocateSolution(eps);
68: /* dispatch solve method */
69: if (eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Left vectors not supported in this solver");
70: eps->ops->solve = EPSSolve_TRLAN;
71: return(0);
72: }
76: static PetscBLASInt MatMult_TRLAN(PetscBLASInt *n,PetscBLASInt *m,PetscReal *xin,PetscBLASInt *ldx,PetscReal *yout,PetscBLASInt *ldy)
77: {
79: Vec x,y;
80: PetscBLASInt i;
83: VecCreateMPIWithArray(PetscObjectComm((PetscObject)globaleps),1,*n,PETSC_DECIDE,NULL,&x);
84: VecCreateMPIWithArray(PetscObjectComm((PetscObject)globaleps),1,*n,PETSC_DECIDE,NULL,&y);
85: for (i=0;i<*m;i++) {
86: VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
87: VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
88: STApply(globaleps->st,x,y);
89: IPOrthogonalize(globaleps->ip,0,NULL,globaleps->nds,NULL,globaleps->defl,y,NULL,NULL,NULL);
90: VecResetArray(x);
91: VecResetArray(y);
92: }
93: VecDestroy(&x);
94: VecDestroy(&y);
95: return(0);
96: }
100: PetscErrorCode EPSSolve_TRLAN(EPS eps)
101: {
103: PetscInt i;
104: PetscBLASInt ipar[32],n,lohi,stat,ncv;
105: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
106: PetscScalar *pV;
109: PetscBLASIntCast(eps->ncv,&ncv);
110: PetscBLASIntCast(eps->nloc,&n);
112: if (eps->which==EPS_LARGEST_REAL || eps->which==EPS_TARGET_REAL) lohi = 1;
113: else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
114: else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
116: globaleps = eps;
118: ipar[0] = 0; /* stat: error flag */
119: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
120: PetscBLASIntCast(eps->nev,&ipar[2]); /* number of desired eigenpairs */
121: ipar[3] = 0; /* number of eigenpairs already converged */
122: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
123: ipar[5] = tr->restart; /* restarting scheme */
124: PetscBLASIntCast(eps->max_it,&ipar[6]); /* maximum number of MATVECs */
125: PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&ipar[7]);
126: ipar[8] = 0; /* verboseness */
127: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
128: ipar[10] = 1; /* use supplied starting vector */
129: ipar[11] = 0; /* checkpointing flag */
130: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
131: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
132: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
134: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
135: EPSGetStartVector(eps,0,eps->V[0],NULL);
136: VecGetArray(eps->V[0],&pV);
138: PetscStackCall("TRLan",TRLan_(MatMult_TRLAN,ipar,&n,&ncv,eps->eigr,pV,&n,tr->work,&tr->lwork));
140: VecRestoreArray(eps->V[0],&pV);
142: stat = ipar[0];
143: eps->nconv = ipar[3];
144: eps->its = ipar[25];
145: eps->reason = EPS_CONVERGED_TOL;
147: if (stat!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);
148: return(0);
149: }
153: PetscErrorCode EPSReset_TRLAN(EPS eps)
154: {
156: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
159: PetscFree(tr->work);
160: EPSFreeSolution(eps);
161: return(0);
162: }
166: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
167: {
171: PetscFree(eps->data);
172: return(0);
173: }
177: PETSC_EXTERN 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: }