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