Actual source code: dsghep.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:       
  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>


 28: PetscErrorCode DSAllocate_GHEP(DS ds,PetscInt ld)
 29: {

 33:   DSAllocateMat_Private(ds,DS_MAT_A);
 34:   DSAllocateMat_Private(ds,DS_MAT_B);
 35:   DSAllocateMat_Private(ds,DS_MAT_Q);
 36:   PetscFree(ds->perm);
 37:   PetscMalloc(ld*sizeof(PetscInt),&ds->perm);
 38:   PetscLogObjectMemory(ds,ld*sizeof(PetscInt));
 39:   return(0);
 40: }

 44: PetscErrorCode DSView_GHEP(DS ds,PetscViewer viewer)
 45: {



 51:   DSViewMat_Private(ds,viewer,DS_MAT_A);
 52:   DSViewMat_Private(ds,viewer,DS_MAT_B);
 53:   if (ds->state>DS_STATE_INTERMEDIATE) {
 54:     DSViewMat_Private(ds,viewer,DS_MAT_Q);
 55:   }
 56:   return(0);
 57: }

 61: PetscErrorCode DSVectors_GHEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
 62: {
 63:   PetscScalar    *Q = ds->mat[DS_MAT_Q];
 64:   PetscInt       ld = ds->ld,i;

 68:   if (rnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
 69:   switch (mat) {
 70:     case DS_MAT_X:
 71:       if (j) {
 72:         if (ds->state>=DS_STATE_CONDENSED) {
 73:           PetscMemcpy(ds->mat[mat]+(*j)*ld,Q+(*j)*ld,ld*sizeof(PetscScalar));
 74:         } else {
 75:           PetscMemzero(ds->mat[mat]+(*j)*ld,ld*sizeof(PetscScalar));
 76:           *(ds->mat[mat]+(*j)+(*j)*ld) = 1.0;
 77:         }
 78:       } else {
 79:         if (ds->state>=DS_STATE_CONDENSED) {
 80:           PetscMemcpy(ds->mat[mat],Q,ld*ld*sizeof(PetscScalar));
 81:         } else {
 82:           PetscMemzero(ds->mat[mat],ld*ld*sizeof(PetscScalar));
 83:           for (i=0;i<ds->n;i++) *(ds->mat[mat]+i+i*ld) = 1.0;
 84:         }
 85:       }
 86:       break;
 87:     case DS_MAT_Y:
 88:     case DS_MAT_U:
 89:     case DS_MAT_VT:
 90:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
 91:       break;
 92:     default:
 93:       SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
 94:   }
 95:   return(0);
 96: }

100: PetscErrorCode DSNormalize_GHEP(DS ds,DSMatType mat,PetscInt col)
101: {
103:   switch (mat) {
104:     case DS_MAT_X:
105:     case DS_MAT_Q:
106:       /* All the matrices resulting from DSVectors and DSSolve are already B-normalized */
107:       break;
108:     case DS_MAT_Y:
109:     case DS_MAT_U:
110:     case DS_MAT_VT:
111:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
112:       break;
113:     default:
114:       SETERRQ(((PetscObject)ds)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
115:   }
116:   return(0);
117: }

121: PetscErrorCode DSSort_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
122: {
124:   PetscInt       n,l,i,*perm,ld=ds->ld;
125:   PetscScalar    *A;
126: 
128:   if (!ds->comp_fun) return(0);
129:   n = ds->n;
130:   l = ds->l;
131:   A  = ds->mat[DS_MAT_A];
132:   perm = ds->perm;
133:   for(i=l;i<n;i++) wr[i] = A[i+i*ld];
134:   if (rr) {
135:     DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
136:   } else {
137:     DSSortEigenvalues_Private(ds,wr,PETSC_NULL,perm,PETSC_FALSE);
138:   }
139:   for (i=l;i<n;i++) A[i+i*ld] = wr[perm[i]];
140:   for (i=l;i<n;i++) wr[i] = A[i+i*ld];
141:   DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
142:   return(0);
143: }

147: PetscErrorCode DSSolve_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
148: {
149: #if defined(SLEPC_MISSING_LAPACK_SYGVD)
151:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SYGVD - Lapack routine is unavailable");
152: #else
154:   PetscScalar    *work,*A,*B,*Q;
155:   PetscBLASInt   itype = 1,*iwork,info,n1,liwork,ld,lrwork=0,lwork;
156:   PetscInt       off,i;
157: #if defined(PETSC_USE_COMPLEX)
158:   PetscReal      *rwork,*rr;
159: #endif 

162:   n1 = PetscBLASIntCast(ds->n-ds->l);
163:   ld = PetscBLASIntCast(ds->ld);
164:   liwork = PetscBLASIntCast(5*ds->n+3);
165: #if defined(PETSC_USE_COMPLEX)
166:   lwork  = PetscBLASIntCast(ds->n*ds->n+2*ds->n);
167:   lrwork = PetscBLASIntCast(2*ds->n*ds->n+5*ds->n+1+n1);
168: #else
169:   lwork  = PetscBLASIntCast(2*ds->n*ds->n+6*ds->n+1);
170: #endif
171:   DSAllocateWork_Private(ds,lwork,lrwork,liwork);
172:   work = ds->work;
173:   iwork = ds->iwork;
174:   off = ds->l+ds->l*ld;
175:   A = ds->mat[DS_MAT_A];
176:   B = ds->mat[DS_MAT_B];
177:   Q = ds->mat[DS_MAT_Q];
178: #if defined(PETSC_USE_COMPLEX)
179:   rr = ds->rwork;
180:   rwork = ds->rwork + n1;
181:   lrwork = ds->lrwork - n1;
182:   LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,rr,work,&lwork,rwork,&lrwork,iwork,&liwork,&info);
183:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack ZHEGVD %d",info);
184:   for(i=0;i<n1;i++) wr[ds->l+i] = rr[i];
185: #else
186:   LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,wr+ds->l,work,&lwork,iwork,&liwork,&info);
187:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSYGVD %d",info);
188: #endif 
189:   PetscMemzero(Q+ds->l*ld,n1*ld*sizeof(PetscScalar));
190:   for(i=ds->l;i<ds->n;i++){
191:     PetscMemcpy(Q+ds->l+i*ld,A+ds->l+i*ld,n1*sizeof(PetscScalar));
192:   }
193:   PetscMemzero(B+ds->l*ld,n1*ld*sizeof(PetscScalar));
194:   PetscMemzero(A+ds->l*ld,n1*ld*sizeof(PetscScalar));
195:   for(i=ds->l;i<ds->n;i++){
196:     if (wi) wi[i] = 0.0;
197:     B[i+i*ld] = 1.0;
198:     A[i+i*ld] = wr[i];
199:   }
200:   return(0);
201: #endif 
202: }

204: EXTERN_C_BEGIN
207: PetscErrorCode DSCreate_GHEP(DS ds)
208: {
210:   ds->ops->allocate      = DSAllocate_GHEP;
211:   ds->ops->view          = DSView_GHEP;
212:   ds->ops->vectors       = DSVectors_GHEP;
213:   ds->ops->solve[0]      = DSSolve_GHEP;
214:   ds->ops->sort          = DSSort_GHEP;
215:   ds->ops->normalize     = DSNormalize_GHEP;
216:   return(0);
217: }
218: EXTERN_C_END