Actual source code: scalapack.c
slepc-3.18.3 2023-03-24
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: 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: EPSCheckHermitianDefinite(eps);
30: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
32: eps->ncv = eps->n;
33: if (eps->mpd!=PETSC_DEFAULT) PetscInfo(eps,"Warning: parameter mpd ignored\n");
34: if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
35: if (!eps->which) EPSSetWhichEigenpairs_Default(eps);
37: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
38: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);
39: EPSAllocateSolution(eps,0);
41: /* convert matrices */
42: MatDestroy(&ctx->As);
43: MatDestroy(&ctx->Bs);
44: STGetNumMatrices(eps->st,&nmat);
45: STGetMatrix(eps->st,0,&A);
46: MatConvert(A,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->As);
47: if (nmat>1) {
48: STGetMatrix(eps->st,1,&B);
49: MatConvert(B,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->Bs);
50: }
51: STGetShift(eps->st,&shift);
52: if (shift != 0.0) {
53: if (nmat>1) MatAXPY(ctx->As,-shift,ctx->Bs,SAME_NONZERO_PATTERN);
54: else MatShift(ctx->As,-shift);
55: }
56: return 0;
57: }
59: PetscErrorCode EPSSolve_ScaLAPACK(EPS eps)
60: {
61: EPS_ScaLAPACK *ctx = (EPS_ScaLAPACK*)eps->data;
62: Mat A = ctx->As,B = ctx->Bs,Q,V;
63: Mat_ScaLAPACK *a = (Mat_ScaLAPACK*)A->data,*b,*q;
64: PetscReal rdummy=0.0,abstol=0.0,*gap=NULL,orfac=-1.0,*w = eps->errest; /* used to store real eigenvalues */
65: PetscScalar *work,minlwork[3];
66: PetscBLASInt i,m,info,idummy=0,lwork=-1,liwork=-1,minliwork,*iwork,*ifail=NULL,*iclustr=NULL,one=1;
67: #if defined(PETSC_USE_COMPLEX)
68: PetscReal *rwork,minlrwork[3];
69: PetscBLASInt lrwork=-1;
70: #endif
72: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
73: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
74: q = (Mat_ScaLAPACK*)Q->data;
76: if (B) {
78: b = (Mat_ScaLAPACK*)B->data;
79: PetscMalloc3(a->grid->nprow*a->grid->npcol,&gap,a->N,&ifail,2*a->grid->nprow*a->grid->npcol,&iclustr);
80: #if !defined(PETSC_USE_COMPLEX)
81: /* allocate workspace */
82: 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));
84: PetscBLASIntCast((PetscInt)minlwork[0],&lwork);
85: liwork = minliwork;
86: /* call computational routine */
87: PetscMalloc2(lwork,&work,liwork,&iwork);
88: 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));
90: PetscFree2(work,iwork);
91: #else
92: /* allocate workspace */
93: 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));
95: PetscBLASIntCast((PetscInt)PetscRealPart(minlwork[0]),&lwork);
96: PetscBLASIntCast((PetscInt)minlrwork[0],&lrwork);
97: lrwork += a->N*a->N;
98: liwork = minliwork;
99: /* call computational routine */
100: PetscMalloc3(lwork,&work,lrwork,&rwork,liwork,&iwork);
101: 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));
103: PetscFree3(work,rwork,iwork);
104: #endif
105: PetscFree3(gap,ifail,iclustr);
107: } else {
109: #if !defined(PETSC_USE_COMPLEX)
110: /* allocate workspace */
111: PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,minlwork,&lwork,&info));
113: PetscBLASIntCast((PetscInt)minlwork[0],&lwork);
114: PetscMalloc1(lwork,&work);
115: /* call computational routine */
116: PetscCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,work,&lwork,&info));
118: PetscFree(work);
119: #else
120: /* allocate workspace */
121: 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));
123: PetscBLASIntCast((PetscInt)PetscRealPart(minlwork[0]),&lwork);
124: lrwork = 4*a->N; /* PetscBLASIntCast((PetscInt)minlrwork[0],&lrwork); */
125: PetscMalloc2(lwork,&work,lrwork,&rwork);
126: /* call computational routine */
127: 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));
129: PetscFree2(work,rwork);
130: #endif
132: }
133: PetscFPTrapPop();
135: for (i=0;i<eps->ncv;i++) {
136: eps->eigr[i] = eps->errest[i];
137: eps->errest[i] = PETSC_MACHINE_EPSILON;
138: }
140: BVGetMat(eps->V,&V);
141: MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
142: BVRestoreMat(eps->V,&V);
143: MatDestroy(&Q);
145: eps->nconv = eps->ncv;
146: eps->its = 1;
147: eps->reason = EPS_CONVERGED_TOL;
148: return 0;
149: }
151: PetscErrorCode EPSDestroy_ScaLAPACK(EPS eps)
152: {
153: PetscFree(eps->data);
154: return 0;
155: }
157: PetscErrorCode EPSReset_ScaLAPACK(EPS eps)
158: {
159: EPS_ScaLAPACK *ctx = (EPS_ScaLAPACK*)eps->data;
161: MatDestroy(&ctx->As);
162: MatDestroy(&ctx->Bs);
163: return 0;
164: }
166: SLEPC_EXTERN PetscErrorCode EPSCreate_ScaLAPACK(EPS eps)
167: {
168: EPS_ScaLAPACK *ctx;
170: PetscNew(&ctx);
171: eps->data = (void*)ctx;
173: eps->categ = EPS_CATEGORY_OTHER;
175: eps->ops->solve = EPSSolve_ScaLAPACK;
176: eps->ops->setup = EPSSetUp_ScaLAPACK;
177: eps->ops->setupsort = EPSSetUpSort_Basic;
178: eps->ops->destroy = EPSDestroy_ScaLAPACK;
179: eps->ops->reset = EPSReset_ScaLAPACK;
180: eps->ops->backtransform = EPSBackTransform_Default;
181: eps->ops->setdefaultst = EPSSetDefaultST_NoFactor;
182: return 0;
183: }