Actual source code: lapack.c
slepc-3.6.1 2015-09-03
1: /*
2: This file implements a wrapper to the LAPACK eigenvalue subroutines.
3: Generalized problems are transformed to standard ones only if necessary.
5: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6: SLEPc - Scalable Library for Eigenvalue Problem Computations
7: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
9: This file is part of SLEPc.
11: SLEPc is free software: you can redistribute it and/or modify it under the
12: terms of version 3 of the GNU Lesser General Public License as published by
13: the Free Software Foundation.
15: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
16: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
18: more details.
20: You should have received a copy of the GNU Lesser General Public License
21: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: */
25: #include <slepc/private/epsimpl.h>
29: PetscErrorCode EPSSetUp_LAPACK(EPS eps)
30: {
31: PetscErrorCode ierr,ierra,ierrb;
32: PetscBool isshift,denseok=PETSC_FALSE;
33: Mat A,B,OP,Adense=NULL,Bdense=NULL;
34: PetscScalar shift,*Ap,*Bp;
35: PetscInt i,ld,nmat;
36: KSP ksp;
37: PC pc;
38: Vec v;
41: eps->ncv = eps->n;
42: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
43: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
44: if (eps->balance!=EPS_BALANCE_NONE) { PetscInfo(eps,"Warning: balancing ignored\n"); }
45: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
46: EPSAllocateSolution(eps,0);
48: /* attempt to get dense representations of A and B separately */
49: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
50: if (isshift) {
51: STGetNumMatrices(eps->st,&nmat);
52: STGetOperators(eps->st,0,&A);
53: if (nmat>1) { STGetOperators(eps->st,1,&B); }
54: PetscPushErrorHandler(PetscIgnoreErrorHandler,NULL);
55: ierra = SlepcMatConvertSeqDense(A,&Adense);
56: if (eps->isgeneralized) {
57: ierrb = SlepcMatConvertSeqDense(B,&Bdense);
58: } else {
59: ierrb = 0;
60: }
61: PetscPopErrorHandler();
62: denseok = (ierra == 0 && ierrb == 0)? PETSC_TRUE: PETSC_FALSE;
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,0);
91: if (denseok) {
92: STGetShift(eps->st,&shift);
93: if (shift != 0.0) {
94: MatShift(Adense,shift);
95: }
96: /* use dummy pc and ksp to avoid problems when B is not positive definite */
97: STGetKSP(eps->st,&ksp);
98: KSPSetType(ksp,KSPPREONLY);
99: KSPGetPC(ksp,&pc);
100: PCSetType(pc,PCNONE);
101: } else {
102: PetscInfo(eps,"Using slow explicit operator\n");
103: STComputeExplicitOperator(eps->st,&OP);
104: MatDestroy(&Adense);
105: SlepcMatConvertSeqDense(OP,&Adense);
106: }
108: /* fill DS matrices */
109: VecCreateSeqWithArray(PETSC_COMM_SELF,1,ld,NULL,&v);
110: DSGetArray(eps->ds,DS_MAT_A,&Ap);
111: for (i=0;i<ld;i++) {
112: VecPlaceArray(v,Ap+i*ld);
113: MatGetColumnVector(Adense,v,i);
114: VecResetArray(v);
115: }
116: DSRestoreArray(eps->ds,DS_MAT_A,&Ap);
117: if (denseok && eps->isgeneralized) {
118: DSGetArray(eps->ds,DS_MAT_B,&Bp);
119: for (i=0;i<ld;i++) {
120: VecPlaceArray(v,Bp+i*ld);
121: MatGetColumnVector(Bdense,v,i);
122: VecResetArray(v);
123: }
124: DSRestoreArray(eps->ds,DS_MAT_B,&Bp);
125: }
126: VecDestroy(&v);
127: MatDestroy(&Adense);
128: if (!denseok) { MatDestroy(&OP); }
129: if (denseok && eps->isgeneralized) { MatDestroy(&Bdense); }
130: return(0);
131: }
135: PetscErrorCode EPSSolve_LAPACK(EPS eps)
136: {
138: PetscInt n=eps->n,i,low,high;
139: PetscScalar *array,*pX;
140: Vec v;
143: DSSolve(eps->ds,eps->eigr,eps->eigi);
144: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
146: /* right eigenvectors */
147: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
148: DSGetArray(eps->ds,DS_MAT_X,&pX);
149: for (i=0;i<eps->ncv;i++) {
150: BVGetColumn(eps->V,i,&v);
151: VecGetOwnershipRange(v,&low,&high);
152: VecGetArray(v,&array);
153: PetscMemcpy(array,pX+i*n+low,(high-low)*sizeof(PetscScalar));
154: VecRestoreArray(v,&array);
155: BVRestoreColumn(eps->V,i,&v);
156: }
157: DSRestoreArray(eps->ds,DS_MAT_X,&pX);
159: eps->nconv = eps->ncv;
160: eps->its = 1;
161: eps->reason = EPS_CONVERGED_TOL;
162: return(0);
163: }
167: PETSC_EXTERN PetscErrorCode EPSCreate_LAPACK(EPS eps)
168: {
170: eps->ops->solve = EPSSolve_LAPACK;
171: eps->ops->setup = EPSSetUp_LAPACK;
172: eps->ops->backtransform = EPSBackTransform_Default;
173: return(0);
174: }