Actual source code: lapack.c
slepc-3.16.0 2021-09-30
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 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 ierr,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;
29: eps->ncv = eps->n;
30: if (eps->mpd!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
31: if (eps->max_it==PETSC_DEFAULT) eps->max_it = 1;
32: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
33: if (eps->which==EPS_ALL && eps->inta!=eps->intb) SETERRQ(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 | EPS_FEATURE_STOPPING);
35: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION | EPS_FEATURE_CONVERGENCE);
36: EPSAllocateSolution(eps,0);
38: /* attempt to get dense representations of A and B separately */
39: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
40: if (isshift) {
41: STGetNumMatrices(eps->st,&nmat);
42: STGetMatrix(eps->st,0,&A);
43: MatHasOperation(A,MATOP_CREATE_SUBMATRICES,&flg);
44: if (flg) {
45: 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: PetscPopErrorHandler();
50: } else ierra = 1;
51: if (nmat>1) {
52: STGetMatrix(eps->st,1,&B);
53: MatHasOperation(B,MATOP_CREATE_SUBMATRICES,&flg);
54: if (flg) {
55: 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: 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) {
70: DSSetType(eps->ds,DSGHEP);
71: } else {
72: DSSetType(eps->ds,DSGNHEP); /* TODO: should be DSGHIEP */
73: }
74: } else {
75: DSSetType(eps->ds,DSGNHEP);
76: }
77: } else {
78: if (eps->ishermitian) {
79: DSSetType(eps->ds,DSHEP);
80: } else {
81: DSSetType(eps->ds,DSNHEP);
82: }
83: }
84: } else {
85: DSSetType(eps->ds,DSNHEP);
86: }
87: DSAllocate(eps->ds,eps->ncv);
88: DSGetLeadingDimension(eps->ds,&ld);
89: DSSetDimensions(eps->ds,eps->ncv,0,0);
91: if (denseok) {
92: STGetShift(eps->st,&shift);
93: if (shift != 0.0) {
94: if (nmat>1) {
95: MatAXPY(Adense,-shift,Bdense,SAME_NONZERO_PATTERN);
96: } else {
97: MatShift(Adense,-shift);
98: }
99: }
100: /* use dummy pc and ksp to avoid problems when B is not positive definite */
101: STGetKSP(eps->st,&ksp);
102: KSPSetType(ksp,KSPPREONLY);
103: KSPGetPC(ksp,&pc);
104: PCSetType(pc,PCNONE);
105: } else {
106: PetscInfo(eps,"Using slow explicit operator\n");
107: STGetOperator(eps->st,&shell);
108: MatComputeOperator(shell,MATDENSE,&OP);
109: STRestoreOperator(eps->st,&shell);
110: MatDestroy(&Adense);
111: MatCreateRedundantMatrix(OP,0,PETSC_COMM_SELF,MAT_INITIAL_MATRIX,&Adense);
112: MatDestroy(&OP);
113: }
115: /* fill DS matrices */
116: VecCreateSeqWithArray(PETSC_COMM_SELF,1,ld,NULL,&v);
117: DSGetArray(eps->ds,DS_MAT_A,&Ap);
118: for (i=0;i<ld;i++) {
119: VecPlaceArray(v,Ap+i*ld);
120: MatGetColumnVector(Adense,v,i);
121: VecResetArray(v);
122: }
123: DSRestoreArray(eps->ds,DS_MAT_A,&Ap);
124: if (denseok && eps->isgeneralized) {
125: DSGetArray(eps->ds,DS_MAT_B,&Bp);
126: for (i=0;i<ld;i++) {
127: VecPlaceArray(v,Bp+i*ld);
128: MatGetColumnVector(Bdense,v,i);
129: VecResetArray(v);
130: }
131: DSRestoreArray(eps->ds,DS_MAT_B,&Bp);
132: }
133: VecDestroy(&v);
134: DSSetState(eps->ds,DS_STATE_RAW);
135: MatDestroy(&Adense);
136: MatDestroy(&Bdense);
137: return(0);
138: }
140: PetscErrorCode EPSSolve_LAPACK(EPS eps)
141: {
143: PetscInt n=eps->n,i,low,high;
144: PetscScalar *array,*pX,*pY;
145: Vec v,w;
148: DSSolve(eps->ds,eps->eigr,eps->eigi);
149: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
150: DSSynchronize(eps->ds,eps->eigr,eps->eigi);
152: /* right eigenvectors */
153: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
154: DSGetArray(eps->ds,DS_MAT_X,&pX);
155: for (i=0;i<eps->ncv;i++) {
156: BVGetColumn(eps->V,i,&v);
157: VecGetOwnershipRange(v,&low,&high);
158: VecGetArray(v,&array);
159: PetscArraycpy(array,pX+i*n+low,high-low);
160: VecRestoreArray(v,&array);
161: BVRestoreColumn(eps->V,i,&v);
162: }
163: DSRestoreArray(eps->ds,DS_MAT_X,&pX);
165: /* left eigenvectors */
166: if (eps->twosided) {
167: DSVectors(eps->ds,DS_MAT_Y,NULL,NULL);
168: DSGetArray(eps->ds,DS_MAT_Y,&pY);
169: for (i=0;i<eps->ncv;i++) {
170: BVGetColumn(eps->W,i,&w);
171: VecGetOwnershipRange(w,&low,&high);
172: VecGetArray(w,&array);
173: PetscArraycpy(array,pY+i*n+low,high-low);
174: VecRestoreArray(w,&array);
175: BVRestoreColumn(eps->W,i,&w);
176: }
177: DSRestoreArray(eps->ds,DS_MAT_Y,&pY);
178: }
180: eps->nconv = eps->ncv;
181: eps->its = 1;
182: eps->reason = EPS_CONVERGED_TOL;
183: return(0);
184: }
186: SLEPC_EXTERN PetscErrorCode EPSCreate_LAPACK(EPS eps)
187: {
189: eps->useds = PETSC_TRUE;
190: eps->categ = EPS_CATEGORY_OTHER;
192: eps->ops->solve = EPSSolve_LAPACK;
193: eps->ops->setup = EPSSetUp_LAPACK;
194: eps->ops->setupsort = EPSSetUpSort_Default;
195: eps->ops->backtransform = EPSBackTransform_Default;
196: return(0);
197: }