Actual source code: lapack.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 the LAPACK eigenvalue subroutines.
12: Generalized problems are transformed to standard ones only if necessary.
13: */
15: #include <slepc/private/epsimpl.h>
17: static PetscErrorCode EPSSetUp_LAPACK(EPS eps)
18: {
19: int ierra,ierrb;
20: PetscBool isshift,flg,denseok=PETSC_FALSE;
21: Mat A,B,OP,shell,Ar,Br,Adense=NULL,Bdense=NULL,Ads,Bds;
22: PetscScalar shift;
23: PetscInt nmat;
24: KSP ksp;
25: PC pc;
27: PetscFunctionBegin;
28: if (eps->nev==0) eps->nev = 1;
29: eps->ncv = eps->n;
30: if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
31: if (eps->max_it==PETSC_DETERMINE) eps->max_it = 1;
32: if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
33: PetscCheck(eps->which!=EPS_ALL || eps->inta==eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support interval computation");
34: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION);
35: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING);
36: PetscCall(EPSAllocateSolution(eps,0));
38: /* attempt to get dense representations of A and B separately */
39: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift));
40: if (isshift) {
41: PetscCall(STGetNumMatrices(eps->st,&nmat));
42: PetscCall(STGetMatrix(eps->st,0,&A));
43: PetscCall(MatHasOperation(A,MATOP_CREATE_SUBMATRICES,&flg));
44: if (flg) {
45: PetscCall(PetscPushErrorHandler(PetscReturnErrorHandler,NULL));
46: ierra = MatCreateRedundantMatrix(A,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Ar);
47: if (!ierra) ierra |= MatConvert(Ar,MATSEQDENSE,MAT_INITIAL_MATRIX,&Adense);
48: ierra |= MatDestroy(&Ar);
49: PetscCall(PetscPopErrorHandler());
50: } else ierra = 1;
51: if (nmat>1) {
52: PetscCall(STGetMatrix(eps->st,1,&B));
53: PetscCall(MatHasOperation(B,MATOP_CREATE_SUBMATRICES,&flg));
54: if (flg) {
55: PetscCall(PetscPushErrorHandler(PetscReturnErrorHandler,NULL));
56: ierrb = MatCreateRedundantMatrix(B,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Br);
57: if (!ierrb) ierrb |= MatConvert(Br,MATSEQDENSE,MAT_INITIAL_MATRIX,&Bdense);
58: ierrb |= MatDestroy(&Br);
59: PetscCall(PetscPopErrorHandler());
60: } else ierrb = 1;
61: } else ierrb = 0;
62: denseok = PetscNot(ierra || ierrb);
63: }
65: /* setup DS */
66: if (denseok) {
67: if (eps->isgeneralized) {
68: if (eps->ishermitian) {
69: if (eps->ispositive) PetscCall(DSSetType(eps->ds,DSGHEP));
70: else PetscCall(DSSetType(eps->ds,DSGNHEP)); /* TODO: should be DSGHIEP */
71: } else PetscCall(DSSetType(eps->ds,DSGNHEP));
72: } else {
73: if (eps->ishermitian) PetscCall(DSSetType(eps->ds,DSHEP));
74: else PetscCall(DSSetType(eps->ds,DSNHEP));
75: }
76: } else PetscCall(DSSetType(eps->ds,DSNHEP));
77: PetscCall(DSAllocate(eps->ds,eps->ncv));
78: PetscCall(DSSetDimensions(eps->ds,eps->ncv,0,0));
80: if (denseok) {
81: PetscCall(STGetShift(eps->st,&shift));
82: if (shift != 0.0) {
83: if (nmat>1) PetscCall(MatAXPY(Adense,-shift,Bdense,SAME_NONZERO_PATTERN));
84: else PetscCall(MatShift(Adense,-shift));
85: }
86: /* use dummy pc and ksp to avoid problems when B is not positive definite */
87: PetscCall(STGetKSP(eps->st,&ksp));
88: PetscCall(KSPSetType(ksp,KSPPREONLY));
89: PetscCall(KSPGetPC(ksp,&pc));
90: PetscCall(PCSetType(pc,PCNONE));
91: } else {
92: PetscCall(PetscInfo(eps,"Using slow explicit operator\n"));
93: PetscCall(STGetOperator(eps->st,&shell));
94: PetscCall(MatComputeOperator(shell,MATDENSE,&OP));
95: PetscCall(STRestoreOperator(eps->st,&shell));
96: PetscCall(MatDestroy(&Adense));
97: PetscCall(MatCreateRedundantMatrix(OP,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Adense));
98: PetscCall(MatDestroy(&OP));
99: }
101: /* fill DS matrices */
102: PetscCall(DSGetMat(eps->ds,DS_MAT_A,&Ads));
103: PetscCall(MatCopy(Adense,Ads,SAME_NONZERO_PATTERN));
104: PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&Ads));
105: if (denseok && eps->isgeneralized) {
106: PetscCall(DSGetMat(eps->ds,DS_MAT_B,&Bds));
107: PetscCall(MatCopy(Bdense,Bds,SAME_NONZERO_PATTERN));
108: PetscCall(DSRestoreMat(eps->ds,DS_MAT_B,&Bds));
109: }
110: PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
111: PetscCall(MatDestroy(&Adense));
112: PetscCall(MatDestroy(&Bdense));
113: PetscFunctionReturn(PETSC_SUCCESS);
114: }
116: static PetscErrorCode EPSSolve_LAPACK(EPS eps)
117: {
118: PetscInt n=eps->n,i,low,high;
119: PetscScalar *array,*pX,*pY;
120: Vec v,w;
122: PetscFunctionBegin;
123: PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
124: PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL));
125: PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));
127: /* right eigenvectors */
128: PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
129: PetscCall(DSGetArray(eps->ds,DS_MAT_X,&pX));
130: for (i=0;i<eps->ncv;i++) {
131: PetscCall(BVGetColumn(eps->V,i,&v));
132: PetscCall(VecGetOwnershipRange(v,&low,&high));
133: PetscCall(VecGetArray(v,&array));
134: PetscCall(PetscArraycpy(array,pX+i*n+low,high-low));
135: PetscCall(VecRestoreArray(v,&array));
136: PetscCall(BVRestoreColumn(eps->V,i,&v));
137: }
138: PetscCall(DSRestoreArray(eps->ds,DS_MAT_X,&pX));
140: /* left eigenvectors */
141: if (eps->twosided) {
142: PetscCall(DSVectors(eps->ds,DS_MAT_Y,NULL,NULL));
143: PetscCall(DSGetArray(eps->ds,DS_MAT_Y,&pY));
144: for (i=0;i<eps->ncv;i++) {
145: PetscCall(BVGetColumn(eps->W,i,&w));
146: PetscCall(VecGetOwnershipRange(w,&low,&high));
147: PetscCall(VecGetArray(w,&array));
148: PetscCall(PetscArraycpy(array,pY+i*n+low,high-low));
149: PetscCall(VecRestoreArray(w,&array));
150: PetscCall(BVRestoreColumn(eps->W,i,&w));
151: }
152: PetscCall(DSRestoreArray(eps->ds,DS_MAT_Y,&pY));
153: }
155: eps->nconv = eps->ncv;
156: eps->its = 1;
157: eps->reason = EPS_CONVERGED_TOL;
158: PetscFunctionReturn(PETSC_SUCCESS);
159: }
161: SLEPC_EXTERN PetscErrorCode EPSCreate_LAPACK(EPS eps)
162: {
163: PetscFunctionBegin;
164: eps->useds = PETSC_TRUE;
165: eps->categ = EPS_CATEGORY_OTHER;
167: eps->ops->solve = EPSSolve_LAPACK;
168: eps->ops->setup = EPSSetUp_LAPACK;
169: eps->ops->setupsort = EPSSetUpSort_Default;
170: eps->ops->backtransform = EPSBackTransform_Default;
171: PetscFunctionReturn(PETSC_SUCCESS);
172: }