Actual source code: dsnep.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: #include <slepc-private/dsimpl.h> /*I "slepcds.h" I*/
23: #include <slepcblaslapack.h>
27: PetscErrorCode DSAllocate_NEP(DS ds,PetscInt ld)
28: {
32: if (!ds->nf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DSNEP requires passing some functions via DSSetFN()");
33: DSAllocateMat_Private(ds,DS_MAT_X);
34: PetscFree(ds->perm);
35: PetscMalloc(ld*sizeof(PetscInt),&ds->perm);
36: PetscLogObjectMemory(ds,ld*sizeof(PetscInt));
37: return(0);
38: }
42: PetscErrorCode DSView_NEP(DS ds,PetscViewer viewer)
43: {
44: PetscErrorCode ierr;
45: PetscViewerFormat format;
46: PetscInt i;
49: PetscViewerGetFormat(viewer,&format);
50: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
51: for (i=0;i<ds->nf;i++) {
52: FNView(ds->f[i],viewer);
53: DSViewMat_Private(ds,viewer,DSMatExtra[i]);
54: }
55: if (ds->state>DS_STATE_INTERMEDIATE) {
56: DSViewMat_Private(ds,viewer,DS_MAT_X);
57: }
58: return(0);
59: }
63: PetscErrorCode DSVectors_NEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
64: {
66: if (rnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
67: switch (mat) {
68: case DS_MAT_X:
69: break;
70: case DS_MAT_Y:
71: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
72: break;
73: default:
74: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
75: }
76: return(0);
77: }
81: PetscErrorCode DSNormalize_NEP(DS ds,DSMatType mat,PetscInt col)
82: {
84: PetscInt i,i0,i1;
85: PetscBLASInt ld,n,one = 1;
86: PetscScalar norm,*x;
89: switch (mat) {
90: case DS_MAT_X:
91: break;
92: case DS_MAT_Y:
93: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
94: break;
95: default:
96: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
97: }
98: PetscBLASIntCast(ds->n,&n);
99: PetscBLASIntCast(ds->ld,&ld);
100: DSGetArray(ds,mat,&x);
101: if (col < 0) {
102: i0 = 0; i1 = ds->n;
103: } else {
104: i0 = col; i1 = col+1;
105: }
106: for (i=i0;i<i1;i++) {
107: norm = BLASnrm2_(&n,&x[ld*i],&one);
108: norm = 1.0/norm;
109: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,&x[ld*i],&one));
110: }
111: return(0);
112: }
116: PetscErrorCode DSSort_NEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
117: {
119: PetscInt n,l,i,*perm,ld=ds->ld;
120: PetscScalar *A;
123: if (!ds->comparison) return(0);
124: n = ds->n;
125: l = ds->l;
126: A = ds->mat[DS_MAT_A];
127: perm = ds->perm;
128: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
129: if (rr) {
130: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
131: } else {
132: DSSortEigenvalues_Private(ds,wr,NULL,perm,PETSC_FALSE);
133: }
134: for (i=l;i<n;i++) A[i+i*ld] = wr[perm[i]];
135: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
136: DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
137: return(0);
138: }
142: PetscErrorCode DSSolve_NEP_SLP(DS ds,PetscScalar *wr,PetscScalar *wi)
143: {
144: #if defined(SLEPC_MISSING_LAPACK_GGEV)
146: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GGEV - Lapack routine is unavailable");
147: #else
149: PetscScalar *A,*B,*W,*X,*work,*alpha,*beta;
150: PetscScalar norm,sigma,lambda,mu,re,re2;
151: PetscBLASInt info,n,ld,lrwork=0,lwork,one=1;
152: PetscInt it,pos,j,maxit=100,result;
153: PetscReal tol;
154: #if defined(PETSC_USE_COMPLEX)
155: PetscReal *rwork;
156: #else
157: PetscReal *alphai,im,im2;
158: #endif
161: if (!ds->mat[DS_MAT_A]) {
162: DSAllocateMat_Private(ds,DS_MAT_A);
163: }
164: if (!ds->mat[DS_MAT_B]) {
165: DSAllocateMat_Private(ds,DS_MAT_B);
166: }
167: if (!ds->mat[DS_MAT_W]) {
168: DSAllocateMat_Private(ds,DS_MAT_W);
169: }
170: PetscBLASIntCast(ds->n,&n);
171: PetscBLASIntCast(ds->ld,&ld);
172: #if defined(PETSC_USE_COMPLEX)
173: PetscBLASIntCast(2*ds->n+2*ds->n,&lwork);
174: PetscBLASIntCast(8*ds->n,&lrwork);
175: #else
176: PetscBLASIntCast(3*ds->n+8*ds->n,&lwork);
177: #endif
178: DSAllocateWork_Private(ds,lwork,lrwork,0);
179: alpha = ds->work;
180: beta = ds->work + ds->n;
181: #if defined(PETSC_USE_COMPLEX)
182: work = ds->work + 2*ds->n;
183: lwork -= 2*ds->n;
184: #else
185: alphai = ds->work + 2*ds->n;
186: work = ds->work + 3*ds->n;
187: lwork -= 3*ds->n;
188: #endif
189: A = ds->mat[DS_MAT_A];
190: B = ds->mat[DS_MAT_B];
191: W = ds->mat[DS_MAT_W];
192: X = ds->mat[DS_MAT_X];
194: sigma = 0.0;
195: lambda = sigma;
196: tol = 1000*n*PETSC_MACHINE_EPSILON;
198: for (it=0;it<maxit;it++) {
200: /* evaluate T and T' */
201: DSComputeMatrix(ds,lambda,PETSC_FALSE,DS_MAT_A);
202: DSComputeMatrix(ds,lambda,PETSC_TRUE,DS_MAT_B);
204: /* % compute eigenvalue correction mu and eigenvector u */
205: #if defined(PETSC_USE_COMPLEX)
206: rwork = ds->rwork;
207: PetscStackCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,beta,NULL,&ld,W,&ld,work,&lwork,rwork,&info));
208: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack ZGGEV %d",info);
209: #else
210: PetscStackCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,alphai,beta,NULL,&ld,W,&ld,work,&lwork,&info));
211: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DGGEV %d",info);
212: #endif
214: /* find smallest eigenvalue */
215: j = 0;
216: if (beta[j]==0.0) re = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
217: else re = alpha[j]/beta[j];
218: #if !defined(PETSC_USE_COMPLEX)
219: if (beta[j]==0.0) im = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
220: else im = alphai[j]/beta[j];
221: #endif
222: pos = 0;
223: for (j=1;j<n;j++) {
224: if (beta[j]==0.0) re2 = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
225: else re2 = alpha[j]/beta[j];
226: #if !defined(PETSC_USE_COMPLEX)
227: if (beta[j]==0.0) im2 = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
228: else im2 = alphai[j]/beta[j];
229: SlepcCompareSmallestMagnitude(re,im,re2,im2,&result,NULL);
230: #else
231: SlepcCompareSmallestMagnitude(re,0.0,re2,0.0,&result,NULL);
232: #endif
233: if (result > 0) {
234: re = re2;
235: #if !defined(PETSC_USE_COMPLEX)
236: im = im2;
237: #endif
238: pos = j;
239: }
240: }
242: #if !defined(PETSC_USE_COMPLEX)
243: if (im!=0.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSNEP found a complex eigenvalue; try rerunning with complex scalars");
244: #endif
245: mu = alpha[pos];
246: PetscMemcpy(X,W+pos*ld,n*sizeof(PetscScalar));
247: norm = BLASnrm2_(&n,X,&one);
248: norm = 1.0/norm;
249: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X,&one));
251: /* correct eigenvalue approximation */
252: lambda = lambda - mu;
253: if (PetscAbsScalar(mu)<=tol) break;
254: }
256: wr[0] = lambda;
257: if (wi) wi[0] = 0.0;
259: if (it==maxit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"DSNEP did not converge");
260: return(0);
261: #endif
262: }
266: PETSC_EXTERN PetscErrorCode DSCreate_NEP(DS ds)
267: {
269: ds->ops->allocate = DSAllocate_NEP;
270: ds->ops->view = DSView_NEP;
271: ds->ops->vectors = DSVectors_NEP;
272: ds->ops->solve[0] = DSSolve_NEP_SLP;
273: ds->ops->sort = DSSort_NEP;
274: ds->ops->normalize = DSNormalize_NEP;
275: return(0);
276: }