Actual source code: trlan.c
slepc-3.22.2 2024-12-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: static PetscErrorCode EPSSetUp_TRLAN(EPS eps)
24: {
25: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
27: PetscFunctionBegin;
28: EPSCheckHermitian(eps);
29: EPSCheckStandard(eps);
30: EPSCheckNotStructured(eps);
31: PetscCall(PetscBLASIntCast(PetscMax(7,eps->nev+PetscMin(eps->nev,6)),&tr->maxlan));
32: if (eps->ncv!=PETSC_DETERMINE) {
33: PetscCheck(eps->ncv>=eps->nev,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_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
36: if (eps->max_it==PETSC_DETERMINE) eps->max_it = PetscMax(1000,eps->n);
38: if (!eps->which) eps->which = EPS_LARGEST_REAL;
39: PetscCheck(eps->which==EPS_SMALLEST_REAL || eps->which==EPS_LARGEST_REAL || eps->which==EPS_TARGET_REAL,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) PetscCall(PetscBLASIntCast(tr->maxlan*(tr->maxlan+10),&tr->lwork));
45: else PetscCall(PetscBLASIntCast(eps->nloc*(tr->maxlan+1-eps->ncv) + tr->maxlan*(tr->maxlan+10),&tr->lwork));
46: if (tr->work) PetscCall(PetscFree(tr->work));
47: PetscCall(PetscMalloc1(tr->lwork,&tr->work));
49: PetscCall(EPSAllocateSolution(eps,0));
50: PetscFunctionReturn(PETSC_SUCCESS);
51: }
53: static PetscBLASInt MatMult_TRLAN(PetscBLASInt *n,PetscBLASInt *m,PetscReal *xin,PetscBLASInt *ldx,PetscReal *yout,PetscBLASInt *ldy)
54: {
55: Vec x=globaldata.x,y=globaldata.y;
56: EPS eps=globaldata.eps;
57: PetscBLASInt i;
59: PetscFunctionBegin;
60: for (i=0;i<*m;i++) {
61: PetscCall(VecPlaceArray(x,(PetscScalar*)xin+i*(*ldx)));
62: PetscCall(VecPlaceArray(y,(PetscScalar*)yout+i*(*ldy)));
63: PetscCall(STApply(eps->st,x,y));
64: PetscCall(BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL));
65: PetscCall(VecResetArray(x));
66: PetscCall(VecResetArray(y));
67: }
68: PetscFunctionReturn(PETSC_SUCCESS);
69: }
71: static PetscErrorCode EPSSolve_TRLAN(EPS eps)
72: {
73: PetscInt i;
74: PetscBLASInt ipar[32],n,lohi,stat,ncv;
75: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
76: PetscScalar *pV;
77: Vec v0;
78: Mat A;
79: #if !defined(PETSC_HAVE_MPIUNI)
80: MPI_Fint fcomm;
81: #endif
83: PetscFunctionBegin;
84: PetscCall(PetscBLASIntCast(eps->ncv,&ncv));
85: PetscCall(PetscBLASIntCast(eps->nloc,&n));
87: PetscCheck(eps->which==EPS_LARGEST_REAL || eps->which==EPS_TARGET_REAL || eps->which==EPS_SMALLEST_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"Wrong value of eps->which");
88: lohi = (eps->which==EPS_SMALLEST_REAL)? -1: 1;
90: globaldata.eps = eps;
91: PetscCall(STGetMatrix(eps->st,0,&A));
92: PetscCall(MatCreateVecsEmpty(A,&globaldata.x,&globaldata.y));
94: ipar[0] = 0; /* stat: error flag */
95: ipar[1] = lohi; /* smallest (lohi<0) or largest eigenvalues (lohi>0) */
96: PetscCall(PetscBLASIntCast(eps->nev,&ipar[2])); /* number of desired eigenpairs */
97: ipar[3] = 0; /* number of eigenpairs already converged */
98: ipar[4] = tr->maxlan; /* maximum Lanczos basis size */
99: ipar[5] = tr->restart; /* restarting scheme */
100: PetscCall(PetscBLASIntCast(eps->max_it,&ipar[6])); /* maximum number of MATVECs */
101: #if !defined(PETSC_HAVE_MPIUNI)
102: fcomm = MPI_Comm_c2f(PetscObjectComm((PetscObject)eps));
103: ipar[7] = fcomm;
104: #endif
105: ipar[8] = 0; /* verboseness */
106: ipar[9] = 99; /* Fortran IO unit number used to write log messages */
107: ipar[10] = 1; /* use supplied starting vector */
108: ipar[11] = 0; /* checkpointing flag */
109: ipar[12] = 98; /* Fortran IO unit number used to write checkpoint files */
110: ipar[13] = 0; /* number of flops per matvec per PE (not used) */
111: tr->work[0] = eps->tol; /* relative tolerance on residual norms */
113: for (i=0;i<eps->ncv;i++) eps->eigr[i]=0.0;
114: PetscCall(EPSGetStartVector(eps,0,NULL));
115: PetscCall(BVSetActiveColumns(eps->V,0,0)); /* just for deflation space */
116: PetscCall(BVGetColumn(eps->V,0,&v0));
117: PetscCall(VecGetArray(v0,&pV));
119: PetscStackCallExternalVoid("TRLan",TRLan_(MatMult_TRLAN,ipar,&n,&ncv,eps->eigr,pV,&n,tr->work,&tr->lwork));
121: PetscCall(VecRestoreArray(v0,&pV));
122: PetscCall(BVRestoreColumn(eps->V,0,&v0));
124: stat = ipar[0];
125: eps->nconv = ipar[3];
126: eps->its = ipar[25];
127: eps->reason = EPS_CONVERGED_TOL;
129: PetscCall(VecDestroy(&globaldata.x));
130: PetscCall(VecDestroy(&globaldata.y));
131: PetscCheck(stat==0,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in TRLAN (code=%" PetscBLASInt_FMT ")",stat);
132: PetscFunctionReturn(PETSC_SUCCESS);
133: }
135: static PetscErrorCode EPSReset_TRLAN(EPS eps)
136: {
137: EPS_TRLAN *tr = (EPS_TRLAN*)eps->data;
139: PetscFunctionBegin;
140: PetscCall(PetscFree(tr->work));
141: PetscFunctionReturn(PETSC_SUCCESS);
142: }
144: static PetscErrorCode EPSDestroy_TRLAN(EPS eps)
145: {
146: PetscFunctionBegin;
147: PetscCall(PetscFree(eps->data));
148: PetscFunctionReturn(PETSC_SUCCESS);
149: }
151: SLEPC_EXTERN PetscErrorCode EPSCreate_TRLAN(EPS eps)
152: {
153: EPS_TRLAN *ctx;
155: PetscFunctionBegin;
156: PetscCall(PetscNew(&ctx));
157: eps->data = (void*)ctx;
159: eps->ops->solve = EPSSolve_TRLAN;
160: eps->ops->setup = EPSSetUp_TRLAN;
161: eps->ops->setupsort = EPSSetUpSort_Basic;
162: eps->ops->destroy = EPSDestroy_TRLAN;
163: eps->ops->reset = EPSReset_TRLAN;
164: eps->ops->backtransform = EPSBackTransform_Default;
165: PetscFunctionReturn(PETSC_SUCCESS);
166: }