Actual source code: lapack.c
slepc-3.22.2 2024-12-02
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: eps->ncv = eps->n;
29: if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
30: if (eps->max_it==PETSC_DETERMINE) eps->max_it = 1;
31: if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
32: PetscCheck(eps->which!=EPS_ALL || eps->inta==eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support interval computation");
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: PetscCall(EPSAllocateSolution(eps,0));
37: /* attempt to get dense representations of A and B separately */
38: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift));
39: if (isshift) {
40: PetscCall(STGetNumMatrices(eps->st,&nmat));
41: PetscCall(STGetMatrix(eps->st,0,&A));
42: PetscCall(MatHasOperation(A,MATOP_CREATE_SUBMATRICES,&flg));
43: if (flg) {
44: PetscCall(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: PetscCall(PetscPopErrorHandler());
49: } else ierra = 1;
50: if (nmat>1) {
51: PetscCall(STGetMatrix(eps->st,1,&B));
52: PetscCall(MatHasOperation(B,MATOP_CREATE_SUBMATRICES,&flg));
53: if (flg) {
54: PetscCall(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: PetscCall(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) PetscCall(DSSetType(eps->ds,DSGHEP));
69: else PetscCall(DSSetType(eps->ds,DSGNHEP)); /* TODO: should be DSGHIEP */
70: } else PetscCall(DSSetType(eps->ds,DSGNHEP));
71: } else {
72: if (eps->ishermitian) PetscCall(DSSetType(eps->ds,DSHEP));
73: else PetscCall(DSSetType(eps->ds,DSNHEP));
74: }
75: } else PetscCall(DSSetType(eps->ds,DSNHEP));
76: PetscCall(DSAllocate(eps->ds,eps->ncv));
77: PetscCall(DSSetDimensions(eps->ds,eps->ncv,0,0));
79: if (denseok) {
80: PetscCall(STGetShift(eps->st,&shift));
81: if (shift != 0.0) {
82: if (nmat>1) PetscCall(MatAXPY(Adense,-shift,Bdense,SAME_NONZERO_PATTERN));
83: else PetscCall(MatShift(Adense,-shift));
84: }
85: /* use dummy pc and ksp to avoid problems when B is not positive definite */
86: PetscCall(STGetKSP(eps->st,&ksp));
87: PetscCall(KSPSetType(ksp,KSPPREONLY));
88: PetscCall(KSPGetPC(ksp,&pc));
89: PetscCall(PCSetType(pc,PCNONE));
90: } else {
91: PetscCall(PetscInfo(eps,"Using slow explicit operator\n"));
92: PetscCall(STGetOperator(eps->st,&shell));
93: PetscCall(MatComputeOperator(shell,MATDENSE,&OP));
94: PetscCall(STRestoreOperator(eps->st,&shell));
95: PetscCall(MatDestroy(&Adense));
96: PetscCall(MatCreateRedundantMatrix(OP,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Adense));
97: PetscCall(MatDestroy(&OP));
98: }
100: /* fill DS matrices */
101: PetscCall(DSGetMat(eps->ds,DS_MAT_A,&Ads));
102: PetscCall(MatCopy(Adense,Ads,SAME_NONZERO_PATTERN));
103: PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&Ads));
104: if (denseok && eps->isgeneralized) {
105: PetscCall(DSGetMat(eps->ds,DS_MAT_B,&Bds));
106: PetscCall(MatCopy(Bdense,Bds,SAME_NONZERO_PATTERN));
107: PetscCall(DSRestoreMat(eps->ds,DS_MAT_B,&Bds));
108: }
109: PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
110: PetscCall(MatDestroy(&Adense));
111: PetscCall(MatDestroy(&Bdense));
112: PetscFunctionReturn(PETSC_SUCCESS);
113: }
115: static PetscErrorCode EPSSolve_LAPACK(EPS eps)
116: {
117: PetscInt n=eps->n,i,low,high;
118: PetscScalar *array,*pX,*pY;
119: Vec v,w;
121: PetscFunctionBegin;
122: PetscCall(DSSolve(eps->ds,eps->eigr,eps->eigi));
123: PetscCall(DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL));
124: PetscCall(DSSynchronize(eps->ds,eps->eigr,eps->eigi));
126: /* right eigenvectors */
127: PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
128: PetscCall(DSGetArray(eps->ds,DS_MAT_X,&pX));
129: for (i=0;i<eps->ncv;i++) {
130: PetscCall(BVGetColumn(eps->V,i,&v));
131: PetscCall(VecGetOwnershipRange(v,&low,&high));
132: PetscCall(VecGetArray(v,&array));
133: PetscCall(PetscArraycpy(array,pX+i*n+low,high-low));
134: PetscCall(VecRestoreArray(v,&array));
135: PetscCall(BVRestoreColumn(eps->V,i,&v));
136: }
137: PetscCall(DSRestoreArray(eps->ds,DS_MAT_X,&pX));
139: /* left eigenvectors */
140: if (eps->twosided) {
141: PetscCall(DSVectors(eps->ds,DS_MAT_Y,NULL,NULL));
142: PetscCall(DSGetArray(eps->ds,DS_MAT_Y,&pY));
143: for (i=0;i<eps->ncv;i++) {
144: PetscCall(BVGetColumn(eps->W,i,&w));
145: PetscCall(VecGetOwnershipRange(w,&low,&high));
146: PetscCall(VecGetArray(w,&array));
147: PetscCall(PetscArraycpy(array,pY+i*n+low,high-low));
148: PetscCall(VecRestoreArray(w,&array));
149: PetscCall(BVRestoreColumn(eps->W,i,&w));
150: }
151: PetscCall(DSRestoreArray(eps->ds,DS_MAT_Y,&pY));
152: }
154: eps->nconv = eps->ncv;
155: eps->its = 1;
156: eps->reason = EPS_CONVERGED_TOL;
157: PetscFunctionReturn(PETSC_SUCCESS);
158: }
160: SLEPC_EXTERN PetscErrorCode EPSCreate_LAPACK(EPS eps)
161: {
162: PetscFunctionBegin;
163: eps->useds = PETSC_TRUE;
164: eps->categ = EPS_CATEGORY_OTHER;
166: eps->ops->solve = EPSSolve_LAPACK;
167: eps->ops->setup = EPSSetUp_LAPACK;
168: eps->ops->setupsort = EPSSetUpSort_Default;
169: eps->ops->backtransform = EPSBackTransform_Default;
170: PetscFunctionReturn(PETSC_SUCCESS);
171: }