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