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-2009, Universidad 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 src/eps/impls/trlan/trlanp.h
26: /* Nasty global variable to access EPS data from TRLan_ */
27: static EPS globaleps;
31: PetscErrorCode EPSSetUp_TRLAN(EPS eps)
32: {
34: PetscInt n;
35: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
38: VecGetSize(eps->vec_initial,&n);
39: if (eps->ncv) {
40: if (eps->ncv<eps->nev) SETERRQ(1,"The value of ncv must be at least nev");
41: }
42: else eps->ncv = eps->nev;
43: if (eps->mpd) PetscInfo(eps,"Warning: parameter mpd ignored\n");
44: if (!eps->max_it) eps->max_it = PetscMax(1000,n);
45:
46: if (!eps->ishermitian)
47: SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
49: if (eps->isgeneralized)
50: SETERRQ(PETSC_ERR_SUP,"Requested method is not available for generalized problems");
52: tr->restart = 0;
53: VecGetLocalSize(eps->vec_initial,&n);
54: tr->maxlan = PetscBLASIntCast(eps->nev+PetscMin(eps->nev,6));
55: if (tr->maxlan+1-eps->ncv<=0) { tr->lwork = PetscBLASIntCast(tr->maxlan*(tr->maxlan+10)); }
56: else { tr->lwork = PetscBLASIntCast(n*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10)); }
57: PetscMalloc(tr->lwork*sizeof(PetscReal),&tr->work);
59: if (eps->extraction) {
60: PetscInfo(eps,"Warning: extraction type ignored\n");
61: }
63: EPSAllocateSolution(eps);
64: EPSDefaultGetWork(eps,1);
65: return(0);
66: }
70: static PetscBLASInt MatMult_TRLAN(PetscBLASInt *n,PetscBLASInt *m,PetscReal *xin,PetscBLASInt *ldx,PetscReal *yout,PetscBLASInt *ldy)
71: {
73: Vec x,y;
74: PetscBLASInt i;
77: VecCreateMPIWithArray(((PetscObject)globaleps)->comm,*n,PETSC_DECIDE,PETSC_NULL,&x);
78: VecCreateMPIWithArray(((PetscObject)globaleps)->comm,*n,PETSC_DECIDE,PETSC_NULL,&y);
79: for (i=0;i<*m;i++) {
80: VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
81: VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
82: STApply(globaleps->OP,x,y);
83: IPOrthogonalize(globaleps->ip,globaleps->nds,PETSC_NULL,globaleps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL,globaleps->work[0],PETSC_NULL);
84: VecResetArray(x);
85: VecResetArray(y);
86: }
87: VecDestroy(x);
88: VecDestroy(y);
89: return(0);
90: }
94: PetscErrorCode EPSSolve_TRLAN(EPS eps)
95: {
97: PetscInt i,nn;
98: PetscBLASInt ipar[32], n, lohi, stat, ncv;
99: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
100: PetscScalar *pV;
101:
104: ncv = PetscBLASIntCast(eps->ncv);
105: VecGetLocalSize(eps->vec_initial,&nn);
106: n = PetscBLASIntCast(nn);
107:
108: if (eps->which==EPS_LARGEST_REAL) lohi = 1;
109: else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
110: else SETERRQ(1,"Wrong value of eps->which");
112: globaleps = eps;
114: ipar[0] = 0; /* stat: error flag */
115: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
116: ipar[2] = PetscBLASIntCast(eps->nev); /* number of desired eigenpairs */
117: ipar[3] = 0; /* number of eigenpairs already converged */
118: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
119: ipar[5] = tr->restart; /* restarting scheme */
120: ipar[6] = PetscBLASIntCast(eps->max_it); /* maximum number of MATVECs */
121: ipar[7] = PetscBLASIntCast(MPI_Comm_c2f(((PetscObject)eps)->comm)); /* communicator */
122: ipar[8] = 0; /* verboseness */
123: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
124: ipar[10] = 1; /* use supplied starting vector */
125: ipar[11] = 0; /* checkpointing flag */
126: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
127: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
128: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
130: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
131: EPSGetStartVector(eps,0,eps->V[0],PETSC_NULL);
132: VecGetArray(eps->V[0],&pV);
134: TRLan_ ( MatMult_TRLAN, ipar, &n, &ncv, eps->eigr, pV, &n, tr->work, &tr->lwork );
136: VecRestoreArray( eps->V[0], &pV );
138: stat = ipar[0];
139: eps->nconv = ipar[3];
140: eps->its = ipar[25];
141: eps->reason = EPS_CONVERGED_TOL;
142:
143: if (stat!=0) { SETERRQ1(PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);}
145: return(0);
146: }
150: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
151: {
153: EPS_TRLAN *tr = (EPS_TRLAN *)eps->data;
157: PetscFree(tr->work);
158: PetscFree(eps->data);
159: EPSFreeSolution(eps);
160: return(0);
161: }
166: PetscErrorCode EPSCreate_TRLAN(EPS eps)
167: {
169: EPS_TRLAN *trlan;
172: PetscNew(EPS_TRLAN,&trlan);
173: PetscLogObjectMemory(eps,sizeof(EPS_TRLAN));
174: eps->data = (void *) trlan;
175: eps->ops->solve = EPSSolve_TRLAN;
176: eps->ops->setup = EPSSetUp_TRLAN;
177: eps->ops->destroy = EPSDestroy_TRLAN;
178: eps->ops->backtransform = EPSBackTransform_Default;
179: eps->ops->computevectors = EPSComputeVectors_Default;
180: eps->which = EPS_LARGEST_REAL;
181: return(0);
182: }