Actual source code: trlan.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: This file implements a wrapper to the TRLAN package
12: */
14: #include <slepc/private/epsimpl.h>
15: #include <../src/eps/impls/external/trlan/trlanp.h>
17: /* Nasty global variable to access EPS data from TRLan_ */
18: static struct {
19: EPS eps;
20: Vec x,y;
21: } globaldata;
23: PetscErrorCode EPSSetUp_TRLAN(EPS eps)
24: {
26: PetscBool istrivial;
27: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
30: PetscBLASIntCast(PetscMax(7,eps->nev+PetscMin(eps->nev,6)),&tr->maxlan);
31: if (eps->ncv) {
32: if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
33: } else eps->ncv = tr->maxlan;
34: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
35: if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
37: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
39: if (eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is not available for generalized problems");
41: if (!eps->which) eps->which = EPS_LARGEST_REAL;
42: 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");
43: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
44: if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");
46: tr->restart = 0;
47: if (tr->maxlan+1-eps->ncv<=0) {
48: PetscBLASIntCast(tr->maxlan*(tr->maxlan+10),&tr->lwork);
49: } else {
50: PetscBLASIntCast(eps->nloc*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10),&tr->lwork);
51: }
52: if (tr->work) { PetscFree(tr->work); }
53: PetscMalloc1(tr->lwork,&tr->work);
54: PetscLogObjectMemory((PetscObject)eps,tr->lwork*sizeof(PetscReal));
56: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
57: RGIsTrivial(eps->rg,&istrivial);
58: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
60: EPSAllocateSolution(eps,0);
61: return(0);
62: }
64: static PetscBLASInt MatMult_TRLAN(PetscBLASInt *n,PetscBLASInt *m,PetscReal *xin,PetscBLASInt *ldx,PetscReal *yout,PetscBLASInt *ldy)
65: {
67: Vec x=globaldata.x,y=globaldata.y;
68: EPS eps=globaldata.eps;
69: PetscBLASInt i;
72: for (i=0;i<*m;i++) {
73: VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
74: VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
75: STApply(eps->st,x,y);
76: BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL);
77: VecResetArray(x);
78: VecResetArray(y);
79: }
80: return(0);
81: }
83: PetscErrorCode EPSSolve_TRLAN(EPS eps)
84: {
86: PetscInt i;
87: PetscBLASInt ipar[32],n,lohi,stat,ncv;
88: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
89: PetscScalar *pV;
90: Vec v0;
91: Mat A;
94: PetscBLASIntCast(eps->ncv,&ncv);
95: PetscBLASIntCast(eps->nloc,&n);
97: if (eps->which==EPS_LARGEST_REAL || eps->which==EPS_TARGET_REAL) lohi = 1;
98: else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
99: else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
101: globaldata.eps = eps;
102: STGetMatrix(eps->st,0,&A);
103: MatCreateVecsEmpty(A,&globaldata.x,&globaldata.y);
105: ipar[0] = 0; /* stat: error flag */
106: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
107: PetscBLASIntCast(eps->nev,&ipar[2]); /* number of desired eigenpairs */
108: ipar[3] = 0; /* number of eigenpairs already converged */
109: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
110: ipar[5] = tr->restart; /* restarting scheme */
111: PetscBLASIntCast(eps->max_it,&ipar[6]); /* maximum number of MATVECs */
112: #if !defined(PETSC_HAVE_MPIUNI)
113: PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&ipar[7]);
114: #endif
115: ipar[8] = 0; /* verboseness */
116: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
117: ipar[10] = 1; /* use supplied starting vector */
118: ipar[11] = 0; /* checkpointing flag */
119: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
120: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
121: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
123: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
124: EPSGetStartVector(eps,0,NULL);
125: BVSetActiveColumns(eps->V,0,0); /* just for deflation space */
126: BVGetColumn(eps->V,0,&v0);
127: VecGetArray(v0,&pV);
129: PetscStackCall("TRLan",TRLan_(MatMult_TRLAN,ipar,&n,&ncv,eps->eigr,pV,&n,tr->work,&tr->lwork));
131: VecRestoreArray(v0,&pV);
132: BVRestoreColumn(eps->V,0,&v0);
134: stat = ipar[0];
135: eps->nconv = ipar[3];
136: eps->its = ipar[25];
137: eps->reason = EPS_CONVERGED_TOL;
139: VecDestroy(&globaldata.x);
140: VecDestroy(&globaldata.y);
141: if (stat!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);
142: return(0);
143: }
145: PetscErrorCode EPSReset_TRLAN(EPS eps)
146: {
148: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
151: PetscFree(tr->work);
152: return(0);
153: }
155: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
156: {
160: PetscFree(eps->data);
161: return(0);
162: }
164: SLEPC_EXTERN PetscErrorCode EPSCreate_TRLAN(EPS eps)
165: {
166: EPS_TRLAN *ctx;
170: PetscNewLog(eps,&ctx);
171: eps->data = (void*)ctx;
173: eps->ops->solve = EPSSolve_TRLAN;
174: eps->ops->setup = EPSSetUp_TRLAN;
175: eps->ops->destroy = EPSDestroy_TRLAN;
176: eps->ops->reset = EPSReset_TRLAN;
177: eps->ops->backtransform = EPSBackTransform_Default;
178: return(0);
179: }