Actual source code: lapack.c
slepc-3.18.0 2022-10-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: PetscErrorCode EPSSetUp_LAPACK(EPS eps)
18: {
19: PetscErrorCode 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: eps->ncv = eps->n;
28: if (eps->mpd!=PETSC_DEFAULT) PetscInfo(eps,"Warning: parameter mpd ignored\n");
29: if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
30: if (!eps->which) EPSSetWhichEigenpairs_Default(eps);
32: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_STOPPING);
33: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);
34: EPSAllocateSolution(eps,0);
36: /* attempt to get dense representations of A and B separately */
37: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
38: if (isshift) {
39: STGetNumMatrices(eps->st,&nmat);
40: STGetMatrix(eps->st,0,&A);
41: MatHasOperation(A,MATOP_CREATE_SUBMATRICES,&flg);
42: if (flg) {
43: PetscPushErrorHandler(PetscReturnErrorHandler,NULL);
44: ierra = MatCreateRedundantMatrix(A,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Ar);
45: if (!ierra) ierra |= MatConvert(Ar,MATSEQDENSE,MAT_INITIAL_MATRIX,&Adense);
46: ierra |= MatDestroy(&Ar);
47: PetscPopErrorHandler();
48: } else ierra = 1;
49: if (nmat>1) {
50: STGetMatrix(eps->st,1,&B);
51: MatHasOperation(B,MATOP_CREATE_SUBMATRICES,&flg);
52: if (flg) {
53: PetscPushErrorHandler(PetscReturnErrorHandler,NULL);
54: ierrb = MatCreateRedundantMatrix(B,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Br);
55: if (!ierrb) ierrb |= MatConvert(Br,MATSEQDENSE,MAT_INITIAL_MATRIX,&Bdense);
56: ierrb |= MatDestroy(&Br);
57: PetscPopErrorHandler();
58: } else ierrb = 1;
59: } else ierrb = 0;
60: denseok = PetscNot(ierra || ierrb);
61: }
63: /* setup DS */
64: if (denseok) {
65: if (eps->isgeneralized) {
66: if (eps->ishermitian) {
67: if (eps->ispositive) DSSetType(eps->ds,DSGHEP);
68: else DSSetType(eps->ds,DSGNHEP); /* TODO: should be DSGHIEP */
69: } else DSSetType(eps->ds,DSGNHEP);
70: } else {
71: if (eps->ishermitian) DSSetType(eps->ds,DSHEP);
72: else DSSetType(eps->ds,DSNHEP);
73: }
74: } else DSSetType(eps->ds,DSNHEP);
75: DSAllocate(eps->ds,eps->ncv);
76: DSSetDimensions(eps->ds,eps->ncv,0,0);
78: if (denseok) {
79: STGetShift(eps->st,&shift);
80: if (shift != 0.0) {
81: if (nmat>1) MatAXPY(Adense,-shift,Bdense,SAME_NONZERO_PATTERN);
82: else MatShift(Adense,-shift);
83: }
84: /* use dummy pc and ksp to avoid problems when B is not positive definite */
85: STGetKSP(eps->st,&ksp);
86: KSPSetType(ksp,KSPPREONLY);
87: KSPGetPC(ksp,&pc);
88: PCSetType(pc,PCNONE);
89: } else {
90: PetscInfo(eps,"Using slow explicit operator\n");
91: STGetOperator(eps->st,&shell);
92: MatComputeOperator(shell,MATDENSE,&OP);
93: STRestoreOperator(eps->st,&shell);
94: MatDestroy(&Adense);
95: MatCreateRedundantMatrix(OP,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Adense);
96: MatDestroy(&OP);
97: }
99: /* fill DS matrices */
100: DSGetMat(eps->ds,DS_MAT_A,&Ads);
101: MatCopy(Adense,Ads,SAME_NONZERO_PATTERN);
102: DSRestoreMat(eps->ds,DS_MAT_A,&Ads);
103: if (denseok && eps->isgeneralized) {
104: DSGetMat(eps->ds,DS_MAT_B,&Bds);
105: MatCopy(Bdense,Bds,SAME_NONZERO_PATTERN);
106: DSRestoreMat(eps->ds,DS_MAT_B,&Bds);
107: }
108: DSSetState(eps->ds,DS_STATE_RAW);
109: MatDestroy(&Adense);
110: MatDestroy(&Bdense);
111: return 0;
112: }
114: PetscErrorCode EPSSolve_LAPACK(EPS eps)
115: {
116: PetscInt n=eps->n,i,low,high;
117: PetscScalar *array,*pX,*pY;
118: Vec v,w;
120: DSSolve(eps->ds,eps->eigr,eps->eigi);
121: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
122: DSSynchronize(eps->ds,eps->eigr,eps->eigi);
124: /* right eigenvectors */
125: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
126: DSGetArray(eps->ds,DS_MAT_X,&pX);
127: for (i=0;i<eps->ncv;i++) {
128: BVGetColumn(eps->V,i,&v);
129: VecGetOwnershipRange(v,&low,&high);
130: VecGetArray(v,&array);
131: PetscArraycpy(array,pX+i*n+low,high-low);
132: VecRestoreArray(v,&array);
133: BVRestoreColumn(eps->V,i,&v);
134: }
135: DSRestoreArray(eps->ds,DS_MAT_X,&pX);
137: /* left eigenvectors */
138: if (eps->twosided) {
139: DSVectors(eps->ds,DS_MAT_Y,NULL,NULL);
140: DSGetArray(eps->ds,DS_MAT_Y,&pY);
141: for (i=0;i<eps->ncv;i++) {
142: BVGetColumn(eps->W,i,&w);
143: VecGetOwnershipRange(w,&low,&high);
144: VecGetArray(w,&array);
145: PetscArraycpy(array,pY+i*n+low,high-low);
146: VecRestoreArray(w,&array);
147: BVRestoreColumn(eps->W,i,&w);
148: }
149: DSRestoreArray(eps->ds,DS_MAT_Y,&pY);
150: }
152: eps->nconv = eps->ncv;
153: eps->its = 1;
154: eps->reason = EPS_CONVERGED_TOL;
155: return 0;
156: }
158: SLEPC_EXTERN PetscErrorCode EPSCreate_LAPACK(EPS eps)
159: {
160: eps->useds = PETSC_TRUE;
161: eps->categ = EPS_CATEGORY_OTHER;
163: eps->ops->solve = EPSSolve_LAPACK;
164: eps->ops->setup = EPSSetUp_LAPACK;
165: eps->ops->setupsort = EPSSetUpSort_Default;
166: eps->ops->backtransform = EPSBackTransform_Default;
167: return 0;
168: }