Actual source code: trlan.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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 "trlan.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: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
29: EPSCheckHermitian(eps);
30: EPSCheckStandard(eps);
31: PetscBLASIntCast(PetscMax(7,eps->nev+PetscMin(eps->nev,6)),&tr->maxlan);
32: if (eps->ncv!=PETSC_DEFAULT) {
33: if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev");
34: } else eps->ncv = tr->maxlan;
35: if (eps->mpd!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
36: if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(1000,eps->n);
38: if (!eps->which) eps->which = EPS_LARGEST_REAL;
39: if (eps->which!=EPS_SMALLEST_REAL && eps->which!=EPS_LARGEST_REAL && eps->which!=EPS_TARGET_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only smallest, largest or target real eigenvalues");
40: EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING);
41: EPSCheckIgnored(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_EXTRACTION);
43: tr->restart = 0;
44: if (tr->maxlan+1-eps->ncv<=0) {
45: PetscBLASIntCast(tr->maxlan*(tr->maxlan+10),&tr->lwork);
46: } else {
47: PetscBLASIntCast(eps->nloc*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10),&tr->lwork);
48: }
49: if (tr->work) { PetscFree(tr->work); }
50: PetscMalloc1(tr->lwork,&tr->work);
51: PetscLogObjectMemory((PetscObject)eps,tr->lwork*sizeof(PetscReal));
53: EPSAllocateSolution(eps,0);
54: return(0);
55: }
57: static PetscBLASInt MatMult_TRLAN(PetscBLASInt *n,PetscBLASInt *m,PetscReal *xin,PetscBLASInt *ldx,PetscReal *yout,PetscBLASInt *ldy)
58: {
60: Vec x=globaldata.x,y=globaldata.y;
61: EPS eps=globaldata.eps;
62: PetscBLASInt i;
65: for (i=0;i<*m;i++) {
66: VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx));
67: VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy));
68: STApply(eps->st,x,y);
69: BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL);
70: VecResetArray(x);
71: VecResetArray(y);
72: }
73: return(0);
74: }
76: PetscErrorCode EPSSolve_TRLAN(EPS eps)
77: {
79: PetscInt i;
80: PetscBLASInt ipar[32],n,lohi,stat,ncv;
81: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
82: PetscScalar *pV;
83: Vec v0;
84: Mat A;
85: #if !defined(PETSC_HAVE_MPIUNI)
86: MPI_Fint fcomm;
87: #endif
90: PetscBLASIntCast(eps->ncv,&ncv);
91: PetscBLASIntCast(eps->nloc,&n);
93: if (eps->which==EPS_LARGEST_REAL || eps->which==EPS_TARGET_REAL) lohi = 1;
94: else if (eps->which==EPS_SMALLEST_REAL) lohi = -1;
95: else SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"Wrong value of eps->which");
97: globaldata.eps = eps;
98: STGetMatrix(eps->st,0,&A);
99: MatCreateVecsEmpty(A,&globaldata.x,&globaldata.y);
101: ipar[0] = 0; /* stat: error flag */
102: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
103: PetscBLASIntCast(eps->nev,&ipar[2]); /* number of desired eigenpairs */
104: ipar[3] = 0; /* number of eigenpairs already converged */
105: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
106: ipar[5] = tr->restart; /* restarting scheme */
107: PetscBLASIntCast(eps->max_it,&ipar[6]); /* maximum number of MATVECs */
108: #if !defined(PETSC_HAVE_MPIUNI)
109: fcomm = MPI_Comm_c2f(PetscObjectComm((PetscObject)eps));
110: ipar[7] = fcomm;
111: #endif
112: ipar[8] = 0; /* verboseness */
113: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
114: ipar[10] = 1; /* use supplied starting vector */
115: ipar[11] = 0; /* checkpointing flag */
116: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
117: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
118: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
120: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
121: EPSGetStartVector(eps,0,NULL);
122: BVSetActiveColumns(eps->V,0,0); /* just for deflation space */
123: BVGetColumn(eps->V,0,&v0);
124: VecGetArray(v0,&pV);
126: PetscStackCall("TRLan",TRLan_(MatMult_TRLAN,ipar,&n,&ncv,eps->eigr,pV,&n,tr->work,&tr->lwork));
128: VecRestoreArray(v0,&pV);
129: BVRestoreColumn(eps->V,0,&v0);
131: stat = ipar[0];
132: eps->nconv = ipar[3];
133: eps->its = ipar[25];
134: eps->reason = EPS_CONVERGED_TOL;
136: VecDestroy(&globaldata.x);
137: VecDestroy(&globaldata.y);
138: if (stat!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in TRLAN (code=%d)",stat);
139: return(0);
140: }
142: PetscErrorCode EPSReset_TRLAN(EPS eps)
143: {
145: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
148: PetscFree(tr->work);
149: return(0);
150: }
152: PetscErrorCode EPSDestroy_TRLAN(EPS eps)
153: {
157: PetscFree(eps->data);
158: return(0);
159: }
161: SLEPC_EXTERN PetscErrorCode EPSCreate_TRLAN(EPS eps)
162: {
163: EPS_TRLAN *ctx;
167: PetscNewLog(eps,&ctx);
168: eps->data = (void*)ctx;
170: eps->ops->solve = EPSSolve_TRLAN;
171: eps->ops->setup = EPSSetUp_TRLAN;
172: eps->ops->setupsort = EPSSetUpSort_Basic;
173: eps->ops->destroy = EPSDestroy_TRLAN;
174: eps->ops->reset = EPSReset_TRLAN;
175: eps->ops->backtransform = EPSBackTransform_Default;
176: return(0);
177: }