Actual source code: dsghep.c
slepc-3.7.3 2016-09-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, 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>
23: #include <slepcblaslapack.h>
27: PetscErrorCode DSAllocate_GHEP(DS ds,PetscInt ld)
28: {
32: DSAllocateMat_Private(ds,DS_MAT_A);
33: DSAllocateMat_Private(ds,DS_MAT_B);
34: DSAllocateMat_Private(ds,DS_MAT_Q);
35: PetscFree(ds->perm);
36: PetscMalloc1(ld,&ds->perm);
37: PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
38: return(0);
39: }
43: PetscErrorCode DSView_GHEP(DS ds,PetscViewer viewer)
44: {
48: DSViewMat(ds,viewer,DS_MAT_A);
49: DSViewMat(ds,viewer,DS_MAT_B);
50: if (ds->state>DS_STATE_INTERMEDIATE) {
51: DSViewMat(ds,viewer,DS_MAT_Q);
52: }
53: return(0);
54: }
58: PetscErrorCode DSVectors_GHEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
59: {
60: PetscScalar *Q = ds->mat[DS_MAT_Q];
61: PetscInt ld = ds->ld,i;
65: if (rnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
66: switch (mat) {
67: case DS_MAT_X:
68: case DS_MAT_Y:
69: if (j) {
70: if (ds->state>=DS_STATE_CONDENSED) {
71: PetscMemcpy(ds->mat[mat]+(*j)*ld,Q+(*j)*ld,ld*sizeof(PetscScalar));
72: } else {
73: PetscMemzero(ds->mat[mat]+(*j)*ld,ld*sizeof(PetscScalar));
74: *(ds->mat[mat]+(*j)+(*j)*ld) = 1.0;
75: }
76: } else {
77: if (ds->state>=DS_STATE_CONDENSED) {
78: PetscMemcpy(ds->mat[mat],Q,ld*ld*sizeof(PetscScalar));
79: } else {
80: PetscMemzero(ds->mat[mat],ld*ld*sizeof(PetscScalar));
81: for (i=0;i<ds->n;i++) *(ds->mat[mat]+i+i*ld) = 1.0;
82: }
83: }
84: break;
85: case DS_MAT_U:
86: case DS_MAT_VT:
87: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
88: break;
89: default:
90: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
91: }
92: return(0);
93: }
97: PetscErrorCode DSNormalize_GHEP(DS ds,DSMatType mat,PetscInt col)
98: {
100: PetscInt i,i0,i1;
101: PetscBLASInt ld,n,one = 1;
102: PetscScalar norm,*x;
105: switch (mat) {
106: case DS_MAT_X:
107: case DS_MAT_Y:
108: case DS_MAT_Q:
109: break;
110: case DS_MAT_U:
111: case DS_MAT_VT:
112: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
113: break;
114: default:
115: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
116: }
117: /* All the matrices resulting from DSVectors and DSSolve are B-normalized,
118: but function returns 2-normalized vectors. */
119: PetscBLASIntCast(ds->n,&n);
120: PetscBLASIntCast(ds->ld,&ld);
121: DSGetArray(ds,mat,&x);
122: if (col < 0) {
123: i0 = 0; i1 = ds->n;
124: } else {
125: i0 = col; i1 = col+1;
126: }
127: for (i=i0;i<i1;i++) {
128: norm = BLASnrm2_(&n,&x[ld*i],&one);
129: norm = 1.0/norm;
130: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,&x[ld*i],&one));
131: }
132: return(0);
133: }
137: PetscErrorCode DSSort_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
138: {
140: PetscInt n,l,i,*perm,ld=ds->ld;
141: PetscScalar *A;
144: if (!ds->sc) return(0);
145: n = ds->n;
146: l = ds->l;
147: A = ds->mat[DS_MAT_A];
148: perm = ds->perm;
149: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
150: if (rr) {
151: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
152: } else {
153: DSSortEigenvalues_Private(ds,wr,NULL,perm,PETSC_FALSE);
154: }
155: for (i=l;i<n;i++) A[i+i*ld] = wr[perm[i]];
156: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
157: DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
158: return(0);
159: }
163: PetscErrorCode DSSolve_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
164: {
165: #if defined(SLEPC_MISSING_LAPACK_SYGVD)
167: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SYGVD - Lapack routine is unavailable");
168: #else
170: PetscScalar *work,*A,*B,*Q;
171: PetscBLASInt itype = 1,*iwork,info,n1,liwork,ld,lrwork=0,lwork;
172: PetscInt off,i;
173: #if defined(PETSC_USE_COMPLEX)
174: PetscReal *rwork,*rr;
175: #endif
178: PetscBLASIntCast(ds->n-ds->l,&n1);
179: PetscBLASIntCast(ds->ld,&ld);
180: PetscBLASIntCast(5*ds->n+3,&liwork);
181: #if defined(PETSC_USE_COMPLEX)
182: PetscBLASIntCast(ds->n*ds->n+2*ds->n,&lwork);
183: PetscBLASIntCast(2*ds->n*ds->n+5*ds->n+1+n1,&lrwork);
184: #else
185: PetscBLASIntCast(2*ds->n*ds->n+6*ds->n+1,&lwork);
186: #endif
187: DSAllocateWork_Private(ds,lwork,lrwork,liwork);
188: work = ds->work;
189: iwork = ds->iwork;
190: off = ds->l+ds->l*ld;
191: A = ds->mat[DS_MAT_A];
192: B = ds->mat[DS_MAT_B];
193: Q = ds->mat[DS_MAT_Q];
194: #if defined(PETSC_USE_COMPLEX)
195: rr = ds->rwork;
196: rwork = ds->rwork + n1;
197: lrwork = ds->lrwork - n1;
198: PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,rr,work,&lwork,rwork,&lrwork,iwork,&liwork,&info));
199: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack ZHEGVD %d",info);
200: for (i=0;i<n1;i++) wr[ds->l+i] = rr[i];
201: #else
202: PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,wr+ds->l,work,&lwork,iwork,&liwork,&info));
203: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSYGVD %d",info);
204: #endif
205: PetscMemzero(Q+ds->l*ld,n1*ld*sizeof(PetscScalar));
206: for (i=ds->l;i<ds->n;i++) {
207: PetscMemcpy(Q+ds->l+i*ld,A+ds->l+i*ld,n1*sizeof(PetscScalar));
208: }
209: PetscMemzero(B+ds->l*ld,n1*ld*sizeof(PetscScalar));
210: PetscMemzero(A+ds->l*ld,n1*ld*sizeof(PetscScalar));
211: for (i=ds->l;i<ds->n;i++) {
212: if (wi) wi[i] = 0.0;
213: B[i+i*ld] = 1.0;
214: A[i+i*ld] = wr[i];
215: }
216: return(0);
217: #endif
218: }
222: PETSC_EXTERN PetscErrorCode DSCreate_GHEP(DS ds)
223: {
225: ds->ops->allocate = DSAllocate_GHEP;
226: ds->ops->view = DSView_GHEP;
227: ds->ops->vectors = DSVectors_GHEP;
228: ds->ops->solve[0] = DSSolve_GHEP;
229: ds->ops->sort = DSSort_GHEP;
230: ds->ops->normalize = DSNormalize_GHEP;
231: return(0);
232: }