Actual source code: trlan.c
slepc-3.6.1 2015-09-03
1: /*
2: This file implements a wrapper to the TRLAN package
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, 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>
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: PetscBool istrivial;
38: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
41: PetscBLASIntCast(PetscMax(7,eps->nev+PetscMin(eps->nev,6)),&tr->maxlan);
42: if (eps->ncv) {
43: if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
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);
48: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
50: if (eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is not available for generalized problems");
52: if (!eps->which) eps->which = EPS_LARGEST_REAL;
53: 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");
54: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
56: tr->restart = 0;
57: if (tr->maxlan+1-eps->ncv<=0) {
58: PetscBLASIntCast(tr->maxlan*(tr->maxlan+10),&tr->lwork);
59: } else {
60: PetscBLASIntCast(eps->nloc*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10),&tr->lwork);
61: }
62: PetscMalloc1(tr->lwork,&tr->work);
63: PetscLogObjectMemory((PetscObject)eps,tr->lwork*sizeof(PetscReal));
65: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
66: RGIsTrivial(eps->rg,&istrivial);
67: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
69: EPSAllocateSolution(eps,0);
71: /* dispatch solve method */
72: eps->ops->solve = EPSSolve_TRLAN;
73: return(0);
74: }
78: static PetscBLASInt MatMult_TRLAN(PetscBLASInt *n,PetscBLASInt *m,PetscReal *xin,PetscBLASInt *ldx,PetscReal *yout,PetscBLASInt *ldy)
79: {
81: Vec x,y;
82: PetscBLASInt i;
85: VecCreateMPIWithArray(PetscObjectComm((PetscObject)globaleps),1,*n,PETSC_DECIDE,NULL,&x);
86: VecCreateMPIWithArray(PetscObjectComm((PetscObject)globaleps),1,*n,PETSC_DECIDE,NULL,&y);
87: for (i=0;i<*m;i++) {
88: VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
89: VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
90: STApply(globaleps->st,x,y);
91: BVOrthogonalizeVec(globaleps->V,y,NULL,NULL,NULL);
92: VecResetArray(x);
93: VecResetArray(y);
94: }
95: VecDestroy(&x);
96: VecDestroy(&y);
97: return(0);
98: }
102: PetscErrorCode EPSSolve_TRLAN(EPS eps)
103: {
105: PetscInt i;
106: PetscBLASInt ipar[32],n,lohi,stat,ncv;
107: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
108: PetscScalar *pV;
109: Vec v0;
112: PetscBLASIntCast(eps->ncv,&ncv);
113: PetscBLASIntCast(eps->nloc,&n);
115: if (eps->which==EPS_LARGEST_REAL || eps->which==EPS_TARGET_REAL) lohi = 1;
116: else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
117: else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
119: globaleps = eps;
121: ipar[0] = 0; /* stat: error flag */
122: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
123: PetscBLASIntCast(eps->nev,&ipar[2]); /* number of desired eigenpairs */
124: ipar[3] = 0; /* number of eigenpairs already converged */
125: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
126: ipar[5] = tr->restart; /* restarting scheme */
127: PetscBLASIntCast(eps->max_it,&ipar[6]); /* maximum number of MATVECs */
128: #if !defined(PETSC_HAVE_MPIUNI)
129: PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&ipar[7]);
130: #endif
131: ipar[8] = 0; /* verboseness */
132: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
133: ipar[10] = 1; /* use supplied starting vector */
134: ipar[11] = 0; /* checkpointing flag */
135: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
136: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
137: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
139: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
140: EPSGetStartVector(eps,0,NULL);
141: BVSetActiveColumns(eps->V,0,0); /* just for deflation space */
142: BVGetColumn(eps->V,0,&v0);
143: VecGetArray(v0,&pV);
145: PetscStackCall("TRLan",TRLan_(MatMult_TRLAN,ipar,&n,&ncv,eps->eigr,pV,&n,tr->work,&tr->lwork));
147: VecRestoreArray(v0,&pV);
148: BVRestoreColumn(eps->V,0,&v0);
150: stat = ipar[0];
151: eps->nconv = ipar[3];
152: eps->its = ipar[25];
153: eps->reason = EPS_CONVERGED_TOL;
155: if (stat!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);
156: return(0);
157: }
161: PetscErrorCode EPSReset_TRLAN(EPS eps)
162: {
164: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
167: PetscFree(tr->work);
168: return(0);
169: }
173: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
174: {
178: PetscFree(eps->data);
179: return(0);
180: }
184: PETSC_EXTERN PetscErrorCode EPSCreate_TRLAN(EPS eps)
185: {
186: EPS_TRLAN *ctx;
190: PetscNewLog(eps,&ctx);
191: eps->data = (void*)ctx;
193: eps->ops->setup = EPSSetUp_TRLAN;
194: eps->ops->destroy = EPSDestroy_TRLAN;
195: eps->ops->reset = EPSReset_TRLAN;
196: eps->ops->backtransform = EPSBackTransform_Default;
197: return(0);
198: }