Actual source code: lapack.c
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-2009, Universidad Politecnica de Valencia, Spain
9: This file is part of SLEPc.
10:
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 private/epsimpl.h
26: #include slepcblaslapack.h
28: typedef struct {
29: Mat OP,A,B;
30: } EPS_LAPACK;
34: PetscErrorCode EPSSetUp_LAPACK(EPS eps)
35: {
36: PetscErrorCode ierr,ierra,ierrb;
37: PetscInt N;
38: EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
39: PetscTruth flg;
40: Mat A,B;
41: PetscScalar shift;
42: KSP ksp;
43: PC pc;
44:
46: VecGetSize(eps->vec_initial,&N);
47: eps->ncv = N;
48: if (eps->mpd) PetscInfo(eps,"Warning: parameter mpd ignored\n");
50: if (la->OP) { MatDestroy(la->OP); }
51: if (la->A) { MatDestroy(la->A); }
52: if (la->B) { MatDestroy(la->B); }
54: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&flg);
55: STGetOperators(eps->OP,&A,&B);
56:
57: if (flg) {
58: la->OP = PETSC_NULL;
59: PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
60: ierra = SlepcMatConvertSeqDense(A,&la->A);
61: if (eps->isgeneralized) {
62: ierrb = SlepcMatConvertSeqDense(B,&la->B);
63: } else {
64: ierrb = 0;
65: la->B = PETSC_NULL;
66: }
67: PetscPopErrorHandler();
68: if (ierra == 0 && ierrb == 0) {
69: STGetShift(eps->OP,&shift);
70: if (shift != 0.0) {
71: MatShift(la->A,shift);
72: }
73: /* use dummy pc and ksp to avoid problems when B is not positive definite */
74: STGetKSP(eps->OP,&ksp);
75: KSPSetType(ksp,KSPPREONLY);
76: KSPGetPC(ksp,&pc);
77: PCSetType(pc,PCNONE);
78: EPSAllocateSolution(eps);
79: return(0);
80: }
81: }
82: PetscInfo(eps,"Using slow explicit operator\n");
83: la->A = PETSC_NULL;
84: la->B = PETSC_NULL;
85: STComputeExplicitOperator(eps->OP,&la->OP);
86: PetscTypeCompare((PetscObject)la->OP,MATSEQDENSE,&flg);
87: if (!flg) {
88: SlepcMatConvertSeqDense(la->OP,&la->OP);
89: }
90: if (eps->extraction) {
91: PetscInfo(eps,"Warning: extraction type ignored\n");
92: }
93: EPSAllocateSolution(eps);
94: return(0);
95: }
99: PetscErrorCode EPSSolve_LAPACK(EPS eps)
100: {
102: PetscInt n,i,low,high;
103: PetscMPIInt size;
104: PetscScalar *array,*arrayb,*pV,*pW;
105: PetscReal *w;
106: EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
107: MPI_Comm comm = ((PetscObject)eps)->comm;
108:
110: MPI_Comm_size(comm,&size);
111:
112: VecGetSize(eps->vec_initial,&n);
114: if (size == 1) {
115: VecGetArray(eps->V[0],&pV);
116: } else {
117: PetscMalloc(sizeof(PetscScalar)*n*n,&pV);
118: }
119: if (eps->solverclass == EPS_TWO_SIDE && (la->OP || !eps->ishermitian)) {
120: if (size == 1) {
121: VecGetArray(eps->W[0],&pW);
122: } else {
123: PetscMalloc(sizeof(PetscScalar)*n*n,&pW);
124: }
125: } else pW = PETSC_NULL;
126:
127:
128: if (la->OP) {
129: MatGetArray(la->OP,&array);
130: EPSDenseNHEP(n,array,eps->eigr,eps->eigi,pV,pW);
131: MatRestoreArray(la->OP,&array);
132: } else if (eps->ishermitian) {
133: #if defined(PETSC_USE_COMPLEX)
134: PetscMalloc(n*sizeof(PetscReal),&w);
135: #else
136: w = eps->eigr;
137: #endif
138: MatGetArray(la->A,&array);
139: if (!eps->isgeneralized) {
140: EPSDenseHEP(n,array,n,w,pV);
141: } else {
142: MatGetArray(la->B,&arrayb);
143: EPSDenseGHEP(n,array,arrayb,w,pV);
144: MatRestoreArray(la->B,&arrayb);
145: }
146: MatRestoreArray(la->A,&array);
147: #if defined(PETSC_USE_COMPLEX)
148: for (i=0;i<n;i++) {
149: eps->eigr[i] = w[i];
150: }
151: PetscFree(w);
152: #endif
153: } else {
154: MatGetArray(la->A,&array);
155: if (!eps->isgeneralized) {
156: EPSDenseNHEP(n,array,eps->eigr,eps->eigi,pV,pW);
157: } else {
158: MatGetArray(la->B,&arrayb);
159: EPSDenseGNHEP(n,array,arrayb,eps->eigr,eps->eigi,pV,pW);
160: MatRestoreArray(la->B,&arrayb);
161: }
162: MatRestoreArray(la->A,&array);
163: }
165: if (size == 1) {
166: VecRestoreArray(eps->V[0],&pV);
167: } else {
168: for (i=0; i<eps->ncv; i++) {
169: VecGetOwnershipRange(eps->V[i], &low, &high);
170: VecGetArray(eps->V[i], &array);
171: PetscMemcpy(array, pV+i*n+low, (high-low)*sizeof(PetscScalar));
172: VecRestoreArray(eps->V[i], &array);
173: }
174: PetscFree(pV);
175: }
176: if (pW) {
177: if (size == 1) {
178: VecRestoreArray(eps->W[0],&pW);
179: } else {
180: for (i=0; i<eps->ncv; i++) {
181: VecGetOwnershipRange(eps->W[i], &low, &high);
182: VecGetArray(eps->W[i], &array);
183: PetscMemcpy(array, pW+i*n+low, (high-low)*sizeof(PetscScalar));
184: VecRestoreArray(eps->W[i], &array);
185: }
186: PetscFree(pW);
187: }
188: }
190: eps->nconv = eps->ncv;
191: eps->its = 1;
192: eps->reason = EPS_CONVERGED_TOL;
193:
194: return(0);
195: }
199: PetscErrorCode EPSDestroy_LAPACK(EPS eps)
200: {
202: EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
206: if (la->OP) { MatDestroy(la->OP); }
207: if (la->A) { MatDestroy(la->A); }
208: if (la->B) { MatDestroy(la->B); }
209: PetscFree(eps->data);
210: EPSFreeSolution(eps);
211: return(0);
212: }
217: PetscErrorCode EPSCreate_LAPACK(EPS eps)
218: {
220: EPS_LAPACK *la;
223: PetscNew(EPS_LAPACK,&la);
224: PetscLogObjectMemory(eps,sizeof(EPS_LAPACK));
225: eps->data = (void *) la;
226: eps->ops->solve = EPSSolve_LAPACK;
227: eps->ops->solvets = EPSSolve_LAPACK;
228: eps->ops->setup = EPSSetUp_LAPACK;
229: eps->ops->destroy = EPSDestroy_LAPACK;
230: eps->ops->backtransform = EPSBackTransform_Default;
231: eps->ops->computevectors = EPSComputeVectors_Default;
232: return(0);
233: }