Actual source code: lapack.c
slepc-3.7.1 2016-05-27
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-2016, 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->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"User-defined stopping test not supported");
46: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
47: EPSAllocateSolution(eps,0);
49: /* attempt to get dense representations of A and B separately */
50: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isshift);
51: if (isshift) {
52: STGetNumMatrices(eps->st,&nmat);
53: STGetOperators(eps->st,0,&A);
54: if (nmat>1) { STGetOperators(eps->st,1,&B); }
55: PetscPushErrorHandler(PetscIgnoreErrorHandler,NULL);
56: ierra = SlepcMatConvertSeqDense(A,&Adense);
57: if (nmat>1) {
58: ierrb = SlepcMatConvertSeqDense(B,&Bdense);
59: } else {
60: ierrb = 0;
61: }
62: PetscPopErrorHandler();
63: denseok = PetscNot(ierra || ierrb);
64: }
66: /* setup DS */
67: if (denseok) {
68: if (eps->isgeneralized) {
69: if (eps->ishermitian) {
70: if (eps->ispositive) {
71: DSSetType(eps->ds,DSGHEP);
72: } else {
73: DSSetType(eps->ds,DSGNHEP); /* TODO: should be DSGHIEP */
74: }
75: } else {
76: DSSetType(eps->ds,DSGNHEP);
77: }
78: } else {
79: if (eps->ishermitian) {
80: DSSetType(eps->ds,DSHEP);
81: } else {
82: DSSetType(eps->ds,DSNHEP);
83: }
84: }
85: } else {
86: DSSetType(eps->ds,DSNHEP);
87: }
88: DSAllocate(eps->ds,eps->ncv);
89: DSGetLeadingDimension(eps->ds,&ld);
90: DSSetDimensions(eps->ds,eps->ncv,0,0,0);
92: if (denseok) {
93: STGetShift(eps->st,&shift);
94: if (shift != 0.0) {
95: MatShift(Adense,shift);
96: }
97: /* use dummy pc and ksp to avoid problems when B is not positive definite */
98: STGetKSP(eps->st,&ksp);
99: KSPSetType(ksp,KSPPREONLY);
100: KSPGetPC(ksp,&pc);
101: PCSetType(pc,PCNONE);
102: } else {
103: PetscInfo(eps,"Using slow explicit operator\n");
104: STComputeExplicitOperator(eps->st,&OP);
105: MatDestroy(&Adense);
106: SlepcMatConvertSeqDense(OP,&Adense);
107: }
109: /* fill DS matrices */
110: VecCreateSeqWithArray(PETSC_COMM_SELF,1,ld,NULL,&v);
111: DSGetArray(eps->ds,DS_MAT_A,&Ap);
112: for (i=0;i<ld;i++) {
113: VecPlaceArray(v,Ap+i*ld);
114: MatGetColumnVector(Adense,v,i);
115: VecResetArray(v);
116: }
117: DSRestoreArray(eps->ds,DS_MAT_A,&Ap);
118: if (denseok && eps->isgeneralized) {
119: DSGetArray(eps->ds,DS_MAT_B,&Bp);
120: for (i=0;i<ld;i++) {
121: VecPlaceArray(v,Bp+i*ld);
122: MatGetColumnVector(Bdense,v,i);
123: VecResetArray(v);
124: }
125: DSRestoreArray(eps->ds,DS_MAT_B,&Bp);
126: }
127: VecDestroy(&v);
128: MatDestroy(&Adense);
129: if (!denseok) { MatDestroy(&OP); }
130: if (denseok && eps->isgeneralized) { MatDestroy(&Bdense); }
131: return(0);
132: }
136: PetscErrorCode EPSSolve_LAPACK(EPS eps)
137: {
139: PetscInt n=eps->n,i,low,high;
140: PetscScalar *array,*pX;
141: Vec v;
144: DSSolve(eps->ds,eps->eigr,eps->eigi);
145: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
147: /* right eigenvectors */
148: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
149: DSGetArray(eps->ds,DS_MAT_X,&pX);
150: for (i=0;i<eps->ncv;i++) {
151: BVGetColumn(eps->V,i,&v);
152: VecGetOwnershipRange(v,&low,&high);
153: VecGetArray(v,&array);
154: PetscMemcpy(array,pX+i*n+low,(high-low)*sizeof(PetscScalar));
155: VecRestoreArray(v,&array);
156: BVRestoreColumn(eps->V,i,&v);
157: }
158: DSRestoreArray(eps->ds,DS_MAT_X,&pX);
160: eps->nconv = eps->ncv;
161: eps->its = 1;
162: eps->reason = EPS_CONVERGED_TOL;
163: return(0);
164: }
168: PETSC_EXTERN PetscErrorCode EPSCreate_LAPACK(EPS eps)
169: {
171: eps->ops->solve = EPSSolve_LAPACK;
172: eps->ops->setup = EPSSetUp_LAPACK;
173: eps->ops->backtransform = EPSBackTransform_Default;
174: return(0);
175: }