Actual source code: lapack.c
2: /*
3: This file implements a wrapper to the LAPACK eigenvalue subroutines.
4: Generalized problems are transformed to standard ones only if necessary.
5: */
6: #include src/eps/epsimpl.h
7: #include slepcblaslapack.h
9: typedef struct {
10: Mat OP,A,B;
11: } EPS_LAPACK;
15: PetscErrorCode EPSSetUp_LAPACK(EPS eps)
16: {
17: PetscErrorCode ierr,ierra,ierrb;
18: PetscInt N;
19: EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
20: PetscTruth flg;
21: Mat A,B;
22: PetscScalar shift;
23:
25: VecGetSize(eps->vec_initial,&N);
26: if (eps->nev<1 || eps->nev>N) SETERRQ(1,"Wrong value of nev");
27: eps->ncv = N;
29: if (la->OP) { MatDestroy(la->OP); }
30: if (la->A) { MatDestroy(la->A); }
31: if (la->B) { MatDestroy(la->B); }
33: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&flg);
34: STGetOperators(eps->OP,&A,&B);
35:
36: if (flg) {
37: la->OP = PETSC_NULL;
38: PetscPushErrorHandler(SlepcQuietErrorHandler,PETSC_NULL);
39: ierra = SlepcMatConvertSeqDense(A,&la->A);
40: if (eps->isgeneralized) {
41: ierrb = SlepcMatConvertSeqDense(B,&la->B);
42: } else {
43: ierrb = 0;
44: la->B = PETSC_NULL;
45: }
46: PetscPopErrorHandler();
47: if (ierra == 0 && ierrb == 0) {
48: STGetShift(eps->OP,&shift);
49: if (shift != 0.0) {
50: MatShift(la->A,shift);
51: }
52: EPSAllocateSolutionContiguous(eps);
53: return(0);
54: }
55: }
56: PetscInfo(eps,"Using slow explicit operator\n");
57: la->A = PETSC_NULL;
58: la->B = PETSC_NULL;
59: STComputeExplicitOperator(eps->OP,&la->OP);
60: PetscTypeCompare((PetscObject)la->OP,MATSEQDENSE,&flg);
61: if (!flg) {
62: SlepcMatConvertSeqDense(la->OP,&la->OP);
63: }
64: EPSAllocateSolutionContiguous(eps);
65: return(0);
66: }
70: PetscErrorCode EPSSolve_LAPACK(EPS eps)
71: {
73: PetscInt n,i,low,high;
74: PetscMPIInt size;
75: PetscScalar *array,*arrayb,*pV,*pW;
76: PetscReal *w;
77: EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
78: MPI_Comm comm = eps->comm;
79:
81: MPI_Comm_size(comm,&size);
82:
83: VecGetSize(eps->vec_initial,&n);
85: if (size == 1) {
86: VecGetArray(eps->V[0],&pV);
87: } else {
88: PetscMalloc(sizeof(PetscScalar)*n*n,&pV);
89: }
90: if (eps->solverclass == EPS_TWO_SIDE && (la->OP || !eps->ishermitian)) {
91: if (size == 1) {
92: VecGetArray(eps->W[0],&pW);
93: } else {
94: PetscMalloc(sizeof(PetscScalar)*n*n,&pW);
95: }
96: } else pW = PETSC_NULL;
97:
98:
99: if (la->OP) {
100: MatGetArray(la->OP,&array);
101: EPSDenseNHEP(n,array,eps->eigr,eps->eigi,pV,pW);
102: MatRestoreArray(la->OP,&array);
103: } else if (eps->ishermitian) {
104: #if defined(PETSC_USE_COMPLEX)
105: PetscMalloc(n*sizeof(PetscReal),&w);
106: #else
107: w = eps->eigr;
108: #endif
109: MatGetArray(la->A,&array);
110: if (!eps->isgeneralized) {
111: EPSDenseHEP(n,array,w,pV);
112: } else {
113: MatGetArray(la->B,&arrayb);
114: EPSDenseGHEP(n,array,arrayb,w,pV);
115: MatRestoreArray(la->B,&arrayb);
116: }
117: MatRestoreArray(la->A,&array);
118: #if defined(PETSC_USE_COMPLEX)
119: for (i=0;i<n;i++) {
120: eps->eigr[i] = w[i];
121: }
122: PetscFree(w);
123: #endif
124: for (i=0;i<n;i++) eps->eigi[i]=0.0;
125: } else {
126: MatGetArray(la->A,&array);
127: if (!eps->isgeneralized) {
128: EPSDenseNHEP(n,array,eps->eigr,eps->eigi,pV,pW);
129: } else {
130: MatGetArray(la->B,&arrayb);
131: EPSDenseGNHEP(n,array,arrayb,eps->eigr,eps->eigi,pV,pW);
132: MatRestoreArray(la->B,&arrayb);
133: }
134: MatRestoreArray(la->A,&array);
135: }
137: if (size == 1) {
138: VecRestoreArray(eps->V[0],&pV);
139: } else {
140: for (i=0; i<eps->ncv; i++) {
141: VecGetOwnershipRange(eps->V[i], &low, &high);
142: VecGetArray(eps->V[i], &array);
143: PetscMemcpy(array, pV+i*n+low, (high-low)*sizeof(PetscScalar));
144: VecRestoreArray(eps->V[i], &array);
145: }
146: PetscFree(pV);
147: }
148: if (pW) {
149: if (size == 1) {
150: VecRestoreArray(eps->W[0],&pW);
151: } else {
152: for (i=0; i<eps->ncv; i++) {
153: VecGetOwnershipRange(eps->W[i], &low, &high);
154: VecGetArray(eps->W[i], &array);
155: PetscMemcpy(array, pW+i*n+low, (high-low)*sizeof(PetscScalar));
156: VecRestoreArray(eps->W[i], &array);
157: }
158: PetscFree(pW);
159: }
160: }
162: eps->nconv = eps->ncv;
163: eps->its = 1;
164: eps->reason = EPS_CONVERGED_TOL;
165:
166: return(0);
167: }
171: PetscErrorCode EPSDestroy_LAPACK(EPS eps)
172: {
174: EPS_LAPACK *la = (EPS_LAPACK *)eps->data;
178: if (la->OP) { MatDestroy(la->OP); }
179: if (la->A) { MatDestroy(la->A); }
180: if (la->B) { MatDestroy(la->B); }
181: PetscFree(eps->data);
182: EPSFreeSolutionContiguous(eps);
183: return(0);
184: }
189: PetscErrorCode EPSCreate_LAPACK(EPS eps)
190: {
192: EPS_LAPACK *la;
195: PetscNew(EPS_LAPACK,&la);
196: PetscLogObjectMemory(eps,sizeof(EPS_LAPACK));
197: eps->data = (void *) la;
198: eps->ops->solve = EPSSolve_LAPACK;
199: eps->ops->solvets = EPSSolve_LAPACK;
200: eps->ops->setup = EPSSetUp_LAPACK;
201: eps->ops->destroy = EPSDestroy_LAPACK;
202: eps->ops->backtransform = EPSBackTransform_Default;
203: eps->ops->computevectors = EPSComputeVectors_Default;
204: return(0);
205: }