Actual source code: scalapack.c
slepc-3.23.1 2025-05-01
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 eigensolvers in ScaLAPACK.
12: */
14: #include <slepc/private/epsimpl.h>
15: #include <slepc/private/slepcscalapack.h>
17: typedef struct {
18: Mat As,Bs; /* converted matrices */
19: } EPS_ScaLAPACK;
21: static PetscErrorCode EPSSetUp_ScaLAPACK(EPS eps)
22: {
23: EPS_ScaLAPACK *ctx = (EPS_ScaLAPACK*)eps->data;
24: Mat A,B;
25: PetscInt nmat;
26: PetscBool isshift;
27: PetscScalar shift;
29: PetscFunctionBegin;
30: EPSCheckHermitianDefinite(eps);
31: EPSCheckNotStructured(eps);
32: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift));
33: PetscCheck(isshift,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations");
34: if (eps->nev==0) eps->nev = 1;
35: eps->ncv = eps->n;
36: if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
37: if (eps->max_it==PETSC_DETERMINE) eps->max_it = 1;
38: if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
39: PetscCheck(eps->which!=EPS_ALL || eps->inta==eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support interval computation");
40: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION);
41: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING);
42: PetscCall(EPSAllocateSolution(eps,0));
44: /* convert matrices */
45: PetscCall(MatDestroy(&ctx->As));
46: PetscCall(MatDestroy(&ctx->Bs));
47: PetscCall(STGetNumMatrices(eps->st,&nmat));
48: PetscCall(STGetMatrix(eps->st,0,&A));
49: PetscCall(MatConvert(A,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->As));
50: if (nmat>1) {
51: PetscCall(STGetMatrix(eps->st,1,&B));
52: PetscCall(MatConvert(B,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->Bs));
53: }
54: PetscCall(STGetShift(eps->st,&shift));
55: if (shift != 0.0) {
56: if (nmat>1) PetscCall(MatAXPY(ctx->As,-shift,ctx->Bs,SAME_NONZERO_PATTERN));
57: else PetscCall(MatShift(ctx->As,-shift));
58: }
59: PetscFunctionReturn(PETSC_SUCCESS);
60: }
62: static PetscErrorCode EPSSolve_ScaLAPACK(EPS eps)
63: {
64: EPS_ScaLAPACK *ctx = (EPS_ScaLAPACK*)eps->data;
65: Mat A = ctx->As,B = ctx->Bs,Q,V;
66: Mat_ScaLAPACK *a = (Mat_ScaLAPACK*)A->data,*b,*q;
67: PetscReal rdummy=0.0,abstol=0.0,*gap=NULL,orfac=-1.0,*w = eps->errest; /* used to store real eigenvalues */
68: PetscScalar *work,minlwork[3];
69: PetscBLASInt i,m,info,idummy=0,lwork=-1,liwork=-1,minliwork,*iwork,*ifail=NULL,*iclustr=NULL,one=1;
70: #if defined(PETSC_USE_COMPLEX)
71: PetscReal *rwork,minlrwork[3];
72: PetscBLASInt lrwork=-1;
73: #endif
75: PetscFunctionBegin;
76: PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q));
77: PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
78: q = (Mat_ScaLAPACK*)Q->data;
80: if (B) {
82: b = (Mat_ScaLAPACK*)B->data;
83: PetscCall(PetscMalloc3(a->grid->nprow*a->grid->npcol,&gap,a->N,&ifail,2*a->grid->nprow*a->grid->npcol,&iclustr));
84: #if !defined(PETSC_USE_COMPLEX)
85: /* allocate workspace */
86: PetscCallBLAS("SCALAPACKsygvx",SCALAPACKsygvx_(&one,"V","A","L",&a->N,a->loc,&one,&one,a->desc,b->loc,&one,&one,b->desc,&rdummy,&rdummy,&idummy,&idummy,&abstol,&m,&idummy,w,&orfac,q->loc,&one,&one,q->desc,minlwork,&lwork,&minliwork,&liwork,ifail,iclustr,gap,&info));
87: PetscCheckScaLapackInfo("sygvx",info);
88: PetscCall(PetscBLASIntCast((PetscInt)minlwork[0],&lwork));
89: liwork = minliwork;
90: /* call computational routine */
91: PetscCall(PetscMalloc2(lwork,&work,liwork,&iwork));
92: PetscCallBLAS("SCALAPACKsygvx",SCALAPACKsygvx_(&one,"V","A","L",&a->N,a->loc,&one,&one,a->desc,b->loc,&one,&one,b->desc,&rdummy,&rdummy,&idummy,&idummy,&abstol,&m,&idummy,w,&orfac,q->loc,&one,&one,q->desc,work,&lwork,iwork,&liwork,ifail,iclustr,gap,&info));
93: PetscCheckScaLapackInfo("sygvx",info);
94: PetscCall(PetscFree2(work,iwork));
95: #else
96: /* allocate workspace */
97: PetscCallBLAS("SCALAPACKsygvx",SCALAPACKsygvx_(&one,"V","A","L",&a->N,a->loc,&one,&one,a->desc,b->loc,&one,&one,b->desc,&rdummy,&rdummy,&idummy,&idummy,&abstol,&m,&idummy,w,&orfac,q->loc,&one,&one,q->desc,minlwork,&lwork,minlrwork,&lrwork,&minliwork,&liwork,ifail,iclustr,gap,&info));
98: PetscCheckScaLapackInfo("sygvx",info);
99: PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(minlwork[0]),&lwork));
100: PetscCall(PetscBLASIntCast((PetscInt)minlrwork[0],&lrwork));
101: lrwork += a->N*a->N;
102: liwork = minliwork;
103: /* call computational routine */
104: PetscCall(PetscMalloc3(lwork,&work,lrwork,&rwork,liwork,&iwork));
105: PetscCallBLAS("SCALAPACKsygvx",SCALAPACKsygvx_(&one,"V","A","L",&a->N,a->loc,&one,&one,a->desc,b->loc,&one,&one,b->desc,&rdummy,&rdummy,&idummy,&idummy,&abstol,&m,&idummy,w,&orfac,q->loc,&one,&one,q->desc,work,&lwork,rwork,&lrwork,iwork,&liwork,ifail,iclustr,gap,&info));
106: PetscCheckScaLapackInfo("sygvx",info);
107: PetscCall(PetscFree3(work,rwork,iwork));
108: #endif
109: PetscCall(PetscFree3(gap,ifail,iclustr));
111: } else {
113: #if !defined(PETSC_USE_COMPLEX)
114: /* allocate workspace */
115: PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,minlwork,&lwork,&info));
116: PetscCheckScaLapackInfo("syev",info);
117: PetscCall(PetscBLASIntCast((PetscInt)minlwork[0],&lwork));
118: PetscCall(PetscMalloc1(lwork,&work));
119: /* call computational routine */
120: PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,work,&lwork,&info));
121: PetscCheckScaLapackInfo("syev",info);
122: PetscCall(PetscFree(work));
123: #else
124: /* allocate workspace */
125: PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,minlwork,&lwork,minlrwork,&lrwork,&info));
126: PetscCheckScaLapackInfo("syev",info);
127: PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(minlwork[0]),&lwork));
128: lrwork = 4*a->N; /* PetscCall(PetscBLASIntCast((PetscInt)minlrwork[0],&lrwork)); */
129: PetscCall(PetscMalloc2(lwork,&work,lrwork,&rwork));
130: /* call computational routine */
131: PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,work,&lwork,rwork,&lrwork,&info));
132: PetscCheckScaLapackInfo("syev",info);
133: PetscCall(PetscFree2(work,rwork));
134: #endif
136: }
137: PetscCall(PetscFPTrapPop());
139: for (i=0;i<eps->ncv;i++) {
140: eps->eigr[i] = eps->errest[i];
141: eps->errest[i] = PETSC_MACHINE_EPSILON;
142: }
144: PetscCall(BVGetMat(eps->V,&V));
145: PetscCall(MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V));
146: PetscCall(BVRestoreMat(eps->V,&V));
147: PetscCall(MatDestroy(&Q));
149: eps->nconv = eps->ncv;
150: eps->its = 1;
151: eps->reason = EPS_CONVERGED_TOL;
152: PetscFunctionReturn(PETSC_SUCCESS);
153: }
155: static PetscErrorCode EPSDestroy_ScaLAPACK(EPS eps)
156: {
157: PetscFunctionBegin;
158: PetscCall(PetscFree(eps->data));
159: PetscFunctionReturn(PETSC_SUCCESS);
160: }
162: static PetscErrorCode EPSReset_ScaLAPACK(EPS eps)
163: {
164: EPS_ScaLAPACK *ctx = (EPS_ScaLAPACK*)eps->data;
166: PetscFunctionBegin;
167: PetscCall(MatDestroy(&ctx->As));
168: PetscCall(MatDestroy(&ctx->Bs));
169: PetscFunctionReturn(PETSC_SUCCESS);
170: }
172: SLEPC_EXTERN PetscErrorCode EPSCreate_ScaLAPACK(EPS eps)
173: {
174: EPS_ScaLAPACK *ctx;
176: PetscFunctionBegin;
177: PetscCall(PetscNew(&ctx));
178: eps->data = (void*)ctx;
180: eps->categ = EPS_CATEGORY_OTHER;
182: eps->ops->solve = EPSSolve_ScaLAPACK;
183: eps->ops->setup = EPSSetUp_ScaLAPACK;
184: eps->ops->setupsort = EPSSetUpSort_Basic;
185: eps->ops->destroy = EPSDestroy_ScaLAPACK;
186: eps->ops->reset = EPSReset_ScaLAPACK;
187: eps->ops->backtransform = EPSBackTransform_Default;
188: eps->ops->setdefaultst = EPSSetDefaultST_NoFactor;
189: PetscFunctionReturn(PETSC_SUCCESS);
190: }