Actual source code: scalapack.c
slepc-3.15.1 2021-05-28
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 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: {
24: EPS_ScaLAPACK *ctx = (EPS_ScaLAPACK*)eps->data;
25: Mat A,B;
26: PetscInt nmat;
27: PetscBool isshift;
28: PetscScalar shift;
31: EPSCheckHermitianDefinite(eps);
32: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
33: if (!isshift) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support spectral transformations");
34: eps->ncv = eps->n;
35: if (eps->mpd!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
36: if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
37: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
38: if (eps->which==EPS_ALL && eps->inta!=eps->intb) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support interval computation");
39: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
40: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);
41: EPSAllocateSolution(eps,0);
43: /* convert matrices */
44: MatDestroy(&ctx->As);
45: MatDestroy(&ctx->Bs);
46: STGetNumMatrices(eps->st,&nmat);
47: STGetMatrix(eps->st,0,&A);
48: MatConvert(A,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->As);
49: if (nmat>1) {
50: STGetMatrix(eps->st,1,&B);
51: MatConvert(B,MATSCALAPACK,MAT_INITIAL_MATRIX,&ctx->Bs);
52: }
53: STGetShift(eps->st,&shift);
54: if (shift != 0.0) {
55: if (nmat>1) {
56: MatAXPY(ctx->As,-shift,ctx->Bs,SAME_NONZERO_PATTERN);
57: } else {
58: MatShift(ctx->As,-shift);
59: }
60: }
61: return(0);
62: }
64: PetscErrorCode EPSSolve_ScaLAPACK(EPS eps)
65: {
67: EPS_ScaLAPACK *ctx = (EPS_ScaLAPACK*)eps->data;
68: Mat A = ctx->As,B = ctx->Bs,Q,V;
69: Mat_ScaLAPACK *a = (Mat_ScaLAPACK*)A->data,*b,*q;
70: PetscReal rdummy=0.0,abstol=0.0,*gap=NULL,orfac=-1.0,*w = eps->errest; /* used to store real eigenvalues */
71: PetscScalar *work,minlwork[3];
72: PetscBLASInt i,m,info,idummy=0,lwork=-1,liwork=-1,minliwork,*iwork,*ifail=NULL,*iclustr=NULL,one=1;
73: #if defined(PETSC_USE_COMPLEX)
74: PetscReal *rwork,minlrwork[3];
75: PetscBLASInt lrwork=-1;
76: #endif
79: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
80: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
81: q = (Mat_ScaLAPACK*)Q->data;
83: if (B) {
85: b = (Mat_ScaLAPACK*)B->data;
86: PetscMalloc3(a->grid->nprow*a->grid->npcol,&gap,a->N,&ifail,2*a->grid->nprow*a->grid->npcol,&iclustr);
87: #if !defined(PETSC_USE_COMPLEX)
88: /* allocate workspace */
89: PetscStackCallBLAS("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));
91: PetscBLASIntCast((PetscInt)minlwork[0],&lwork);
92: liwork = minliwork;
93: /* call computational routine */
94: PetscMalloc2(lwork,&work,liwork,&iwork);
95: PetscStackCallBLAS("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));
97: PetscFree2(work,iwork);
98: #else
99: /* allocate workspace */
100: PetscStackCallBLAS("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));
102: PetscBLASIntCast((PetscInt)PetscRealPart(minlwork[0]),&lwork);
103: PetscBLASIntCast((PetscInt)minlrwork[0],&lrwork);
104: lrwork += a->N*a->N;
105: liwork = minliwork;
106: /* call computational routine */
107: PetscMalloc3(lwork,&work,lrwork,&rwork,liwork,&iwork);
108: PetscStackCallBLAS("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));
110: PetscFree3(work,rwork,iwork);
111: #endif
112: PetscFree3(gap,ifail,iclustr);
114: } else {
116: #if !defined(PETSC_USE_COMPLEX)
117: /* allocate workspace */
118: PetscStackCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,minlwork,&lwork,&info));
120: PetscBLASIntCast((PetscInt)minlwork[0],&lwork);
121: PetscMalloc1(lwork,&work);
122: /* call computational routine */
123: PetscStackCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,work,&lwork,&info));
125: PetscFree(work);
126: #else
127: /* allocate workspace */
128: PetscStackCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,minlwork,&lwork,minlrwork,&lrwork,&info));
130: PetscBLASIntCast((PetscInt)PetscRealPart(minlwork[0]),&lwork);
131: lrwork = 4*a->N; /* PetscBLASIntCast((PetscInt)minlrwork[0],&lrwork); */
132: PetscMalloc2(lwork,&work,lrwork,&rwork);
133: /* call computational routine */
134: PetscStackCallBLAS("SCALAPACKsyev",SCALAPACKsyev_("V","L",&a->N,a->loc,&one,&one,a->desc,w,q->loc,&one,&one,q->desc,work,&lwork,rwork,&lrwork,&info));
136: PetscFree2(work,rwork);
137: #endif
139: }
140: PetscFPTrapPop();
142: for (i=0;i<eps->ncv;i++) {
143: eps->eigr[i] = eps->errest[i];
144: eps->errest[i] = PETSC_MACHINE_EPSILON;
145: }
147: BVGetMat(eps->V,&V);
148: MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
149: BVRestoreMat(eps->V,&V);
150: MatDestroy(&Q);
152: eps->nconv = eps->ncv;
153: eps->its = 1;
154: eps->reason = EPS_CONVERGED_TOL;
155: return(0);
156: }
158: PetscErrorCode EPSDestroy_ScaLAPACK(EPS eps)
159: {
163: PetscFree(eps->data);
164: return(0);
165: }
167: PetscErrorCode EPSReset_ScaLAPACK(EPS eps)
168: {
170: EPS_ScaLAPACK *ctx = (EPS_ScaLAPACK*)eps->data;
173: MatDestroy(&ctx->As);
174: MatDestroy(&ctx->Bs);
175: return(0);
176: }
178: SLEPC_EXTERN PetscErrorCode EPSCreate_ScaLAPACK(EPS eps)
179: {
180: EPS_ScaLAPACK *ctx;
184: PetscNewLog(eps,&ctx);
185: eps->data = (void*)ctx;
187: eps->categ = EPS_CATEGORY_OTHER;
189: eps->ops->solve = EPSSolve_ScaLAPACK;
190: eps->ops->setup = EPSSetUp_ScaLAPACK;
191: eps->ops->setupsort = EPSSetUpSort_Basic;
192: eps->ops->destroy = EPSDestroy_ScaLAPACK;
193: eps->ops->reset = EPSReset_ScaLAPACK;
194: eps->ops->backtransform = EPSBackTransform_Default;
195: eps->ops->setdefaultst = EPSSetDefaultST_NoFactor;
196: return(0);
197: }