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-2007, Universidad Politecnica de Valencia, Spain
9: This file is part of SLEPc. See the README file for conditions of use
10: and additional information.
11: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
12: */
14: #include src/eps/epsimpl.h
15: #include slepcblaslapack.h
17: typedef struct {
18: Mat OP,A,B;
19: } EPS_LAPACK;
23: PetscErrorCode EPSSetUp_LAPACK(EPS eps)
24: {
25: PetscErrorCode ierr,ierra,ierrb;
26: PetscInt N;
27: EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
28: PetscTruth flg;
29: Mat A,B;
30: PetscScalar shift;
31: KSP ksp;
32: PC pc;
33:
35: VecGetSize(eps->vec_initial,&N);
36: eps->ncv = N;
38: if (la->OP) { MatDestroy(la->OP); }
39: if (la->A) { MatDestroy(la->A); }
40: if (la->B) { MatDestroy(la->B); }
42: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&flg);
43: STGetOperators(eps->OP,&A,&B);
44:
45: if (flg) {
46: la->OP = PETSC_NULL;
47: PetscPushErrorHandler(PetscIgnoreErrorHandler,PETSC_NULL);
48: ierra = SlepcMatConvertSeqDense(A,&la->A);
49: if (eps->isgeneralized) {
50: ierrb = SlepcMatConvertSeqDense(B,&la->B);
51: } else {
52: ierrb = 0;
53: la->B = PETSC_NULL;
54: }
55: PetscPopErrorHandler();
56: if (ierra == 0 && ierrb == 0) {
57: STGetShift(eps->OP,&shift);
58: if (shift != 0.0) {
59: MatShift(la->A,shift);
60: }
61: /* use dummy pc and ksp to avoid problems when B is not positive definite */
62: STGetKSP(eps->OP,&ksp);
63: KSPSetType(ksp,KSPPREONLY);
64: KSPGetPC(ksp,&pc);
65: PCSetType(pc,PCNONE);
66: EPSAllocateSolutionContiguous(eps);
67: return(0);
68: }
69: }
70: PetscInfo(eps,"Using slow explicit operator\n");
71: la->A = PETSC_NULL;
72: la->B = PETSC_NULL;
73: STComputeExplicitOperator(eps->OP,&la->OP);
74: PetscTypeCompare((PetscObject)la->OP,MATSEQDENSE,&flg);
75: if (!flg) {
76: SlepcMatConvertSeqDense(la->OP,&la->OP);
77: }
78: EPSAllocateSolutionContiguous(eps);
79: return(0);
80: }
84: PetscErrorCode EPSSolve_LAPACK(EPS eps)
85: {
87: PetscInt n,i,low,high;
88: PetscMPIInt size;
89: PetscScalar *array,*arrayb,*pV,*pW;
90: PetscReal *w;
91: EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
92: MPI_Comm comm = eps->comm;
93:
95: MPI_Comm_size(comm,&size);
96:
97: VecGetSize(eps->vec_initial,&n);
99: if (size == 1) {
100: VecGetArray(eps->V[0],&pV);
101: } else {
102: PetscMalloc(sizeof(PetscScalar)*n*n,&pV);
103: }
104: if (eps->solverclass == EPS_TWO_SIDE && (la->OP || !eps->ishermitian)) {
105: if (size == 1) {
106: VecGetArray(eps->W[0],&pW);
107: } else {
108: PetscMalloc(sizeof(PetscScalar)*n*n,&pW);
109: }
110: } else pW = PETSC_NULL;
111:
112:
113: if (la->OP) {
114: MatGetArray(la->OP,&array);
115: EPSDenseNHEP(n,array,eps->eigr,eps->eigi,pV,pW);
116: MatRestoreArray(la->OP,&array);
117: } else if (eps->ishermitian) {
118: #if defined(PETSC_USE_COMPLEX)
119: PetscMalloc(n*sizeof(PetscReal),&w);
120: #else
121: w = eps->eigr;
122: #endif
123: MatGetArray(la->A,&array);
124: if (!eps->isgeneralized) {
125: EPSDenseHEP(n,array,n,w,pV);
126: } else {
127: MatGetArray(la->B,&arrayb);
128: EPSDenseGHEP(n,array,arrayb,w,pV);
129: MatRestoreArray(la->B,&arrayb);
130: }
131: MatRestoreArray(la->A,&array);
132: #if defined(PETSC_USE_COMPLEX)
133: for (i=0;i<n;i++) {
134: eps->eigr[i] = w[i];
135: }
136: PetscFree(w);
137: #endif
138: } else {
139: MatGetArray(la->A,&array);
140: if (!eps->isgeneralized) {
141: EPSDenseNHEP(n,array,eps->eigr,eps->eigi,pV,pW);
142: } else {
143: MatGetArray(la->B,&arrayb);
144: EPSDenseGNHEP(n,array,arrayb,eps->eigr,eps->eigi,pV,pW);
145: MatRestoreArray(la->B,&arrayb);
146: }
147: MatRestoreArray(la->A,&array);
148: }
150: if (size == 1) {
151: VecRestoreArray(eps->V[0],&pV);
152: } else {
153: for (i=0; i<eps->ncv; i++) {
154: VecGetOwnershipRange(eps->V[i], &low, &high);
155: VecGetArray(eps->V[i], &array);
156: PetscMemcpy(array, pV+i*n+low, (high-low)*sizeof(PetscScalar));
157: VecRestoreArray(eps->V[i], &array);
158: }
159: PetscFree(pV);
160: }
161: if (pW) {
162: if (size == 1) {
163: VecRestoreArray(eps->W[0],&pW);
164: } else {
165: for (i=0; i<eps->ncv; i++) {
166: VecGetOwnershipRange(eps->W[i], &low, &high);
167: VecGetArray(eps->W[i], &array);
168: PetscMemcpy(array, pW+i*n+low, (high-low)*sizeof(PetscScalar));
169: VecRestoreArray(eps->W[i], &array);
170: }
171: PetscFree(pW);
172: }
173: }
175: eps->nconv = eps->ncv;
176: eps->its = 1;
177: eps->reason = EPS_CONVERGED_TOL;
178:
179: return(0);
180: }
184: PetscErrorCode EPSDestroy_LAPACK(EPS eps)
185: {
187: EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
191: if (la->OP) { MatDestroy(la->OP); }
192: if (la->A) { MatDestroy(la->A); }
193: if (la->B) { MatDestroy(la->B); }
194: PetscFree(eps->data);
195: EPSFreeSolutionContiguous(eps);
196: return(0);
197: }
202: PetscErrorCode EPSCreate_LAPACK(EPS eps)
203: {
205: EPS_LAPACK *la;
208: PetscNew(EPS_LAPACK,&la);
209: PetscLogObjectMemory(eps,sizeof(EPS_LAPACK));
210: eps->data = (void *) la;
211: eps->ops->solve = EPSSolve_LAPACK;
212: eps->ops->solvets = EPSSolve_LAPACK;
213: eps->ops->setup = EPSSetUp_LAPACK;
214: eps->ops->destroy = EPSDestroy_LAPACK;
215: eps->ops->backtransform = EPSBackTransform_Default;
216: eps->ops->computevectors = EPSComputeVectors_Default;
217: return(0);
218: }