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