Actual source code: dsghep.c

slepc-3.12.2 2020-01-13
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: #include <slepc/private/dsimpl.h>
 12: #include <slepcblaslapack.h>

 14: PetscErrorCode DSAllocate_GHEP(DS ds,PetscInt ld)
 15: {

 19:   DSAllocateMat_Private(ds,DS_MAT_A);
 20:   DSAllocateMat_Private(ds,DS_MAT_B);
 21:   DSAllocateMat_Private(ds,DS_MAT_Q);
 22:   PetscFree(ds->perm);
 23:   PetscMalloc1(ld,&ds->perm);
 24:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
 25:   return(0);
 26: }

 28: PetscErrorCode DSView_GHEP(DS ds,PetscViewer viewer)
 29: {
 30:   PetscErrorCode    ierr;
 31:   PetscViewerFormat format;

 34:   PetscViewerGetFormat(viewer,&format);
 35:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
 36:   DSViewMat(ds,viewer,DS_MAT_A);
 37:   DSViewMat(ds,viewer,DS_MAT_B);
 38:   if (ds->state>DS_STATE_INTERMEDIATE) { DSViewMat(ds,viewer,DS_MAT_Q); }
 39:   if (ds->mat[DS_MAT_X]) { DSViewMat(ds,viewer,DS_MAT_X); }
 40:   return(0);
 41: }

 43: PetscErrorCode DSVectors_GHEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
 44: {
 45:   PetscScalar    *Q = ds->mat[DS_MAT_Q];
 46:   PetscInt       ld = ds->ld,i;

 50:   if (rnorm) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
 51:   switch (mat) {
 52:     case DS_MAT_X:
 53:     case DS_MAT_Y:
 54:       if (j) {
 55:         if (ds->state>=DS_STATE_CONDENSED) {
 56:           PetscArraycpy(ds->mat[mat]+(*j)*ld,Q+(*j)*ld,ld);
 57:         } else {
 58:           PetscArrayzero(ds->mat[mat]+(*j)*ld,ld);
 59:           *(ds->mat[mat]+(*j)+(*j)*ld) = 1.0;
 60:         }
 61:       } else {
 62:         if (ds->state>=DS_STATE_CONDENSED) {
 63:           PetscArraycpy(ds->mat[mat],Q,ld*ld);
 64:         } else {
 65:           PetscArrayzero(ds->mat[mat],ld*ld);
 66:           for (i=0;i<ds->n;i++) *(ds->mat[mat]+i+i*ld) = 1.0;
 67:         }
 68:       }
 69:       break;
 70:     case DS_MAT_U:
 71:     case DS_MAT_VT:
 72:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
 73:       break;
 74:     default:
 75:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
 76:   }
 77:   return(0);
 78: }

 80: PetscErrorCode DSSort_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
 81: {
 83:   PetscInt       n,l,i,*perm,ld=ds->ld;
 84:   PetscScalar    *A;

 87:   if (!ds->sc) return(0);
 88:   n = ds->n;
 89:   l = ds->l;
 90:   A  = ds->mat[DS_MAT_A];
 91:   perm = ds->perm;
 92:   for (i=l;i<n;i++) wr[i] = A[i+i*ld];
 93:   if (rr) {
 94:     DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
 95:   } else {
 96:     DSSortEigenvalues_Private(ds,wr,NULL,perm,PETSC_FALSE);
 97:   }
 98:   for (i=l;i<n;i++) A[i+i*ld] = wr[perm[i]];
 99:   for (i=l;i<n;i++) wr[i] = A[i+i*ld];
100:   DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
101:   return(0);
102: }

104: PetscErrorCode DSSolve_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
105: {
106: #if defined(SLEPC_MISSING_LAPACK_SYGVD)
108:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SYGVD - Lapack routine is unavailable");
109: #else
111:   PetscScalar    *work,*A,*B,*Q;
112:   PetscBLASInt   itype = 1,*iwork,info,n1,liwork,ld,lrwork=0,lwork;
113:   PetscInt       off,i;
114: #if defined(PETSC_USE_COMPLEX)
115:   PetscReal      *rwork,*rr;
116: #endif

119:   PetscBLASIntCast(ds->n-ds->l,&n1);
120:   PetscBLASIntCast(ds->ld,&ld);
121:   PetscBLASIntCast(5*ds->n+3,&liwork);
122: #if defined(PETSC_USE_COMPLEX)
123:   PetscBLASIntCast(ds->n*ds->n+2*ds->n,&lwork);
124:   PetscBLASIntCast(2*ds->n*ds->n+5*ds->n+1+n1,&lrwork);
125: #else
126:   PetscBLASIntCast(2*ds->n*ds->n+6*ds->n+1,&lwork);
127: #endif
128:   DSAllocateWork_Private(ds,lwork,lrwork,liwork);
129:   work = ds->work;
130:   iwork = ds->iwork;
131:   off = ds->l+ds->l*ld;
132:   A = ds->mat[DS_MAT_A];
133:   B = ds->mat[DS_MAT_B];
134:   Q = ds->mat[DS_MAT_Q];
135: #if defined(PETSC_USE_COMPLEX)
136:   rr = ds->rwork;
137:   rwork = ds->rwork + n1;
138:   lrwork = ds->lrwork - n1;
139:   PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,rr,work,&lwork,rwork,&lrwork,iwork,&liwork,&info));
140:   for (i=0;i<n1;i++) wr[ds->l+i] = rr[i];
141: #else
142:   PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,wr+ds->l,work,&lwork,iwork,&liwork,&info));
143: #endif
144:   SlepcCheckLapackInfo("sygvd",info);
145:   PetscArrayzero(Q+ds->l*ld,n1*ld);
146:   for (i=ds->l;i<ds->n;i++) {
147:     PetscArraycpy(Q+ds->l+i*ld,A+ds->l+i*ld,n1);
148:   }
149:   PetscArrayzero(B+ds->l*ld,n1*ld);
150:   PetscArrayzero(A+ds->l*ld,n1*ld);
151:   for (i=ds->l;i<ds->n;i++) {
152:     if (wi) wi[i] = 0.0;
153:     B[i+i*ld] = 1.0;
154:     A[i+i*ld] = wr[i];
155:   }
156:   return(0);
157: #endif
158: }

160: PetscErrorCode DSSynchronize_GHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
161: {
163:   PetscInt       ld=ds->ld,l=ds->l,k;
164:   PetscMPIInt    n,rank,off=0,size,ldn;

167:   k = 2*(ds->n-l)*ld;
168:   if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
169:   if (eigr) k += (ds->n-l);
170:   DSAllocateWork_Private(ds,k,0,0);
171:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
172:   PetscMPIIntCast(ds->n-l,&n);
173:   PetscMPIIntCast(ld*(ds->n-l),&ldn);
174:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
175:   if (!rank) {
176:     MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
177:     MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
178:     if (ds->state>DS_STATE_RAW) {
179:       MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
180:     }
181:     if (eigr) {
182:       MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
183:     }
184:   }
185:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
186:   if (rank) {
187:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
188:     MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
189:     if (ds->state>DS_STATE_RAW) {
190:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
191:     }
192:     if (eigr) {
193:       MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
194:     }
195:   }
196:   return(0);
197: }

199: PetscErrorCode DSHermitian_GHEP(DS ds,DSMatType m,PetscBool *flg)
200: {
202:   if (m==DS_MAT_A || m==DS_MAT_B) *flg = PETSC_TRUE;
203:   else *flg = PETSC_FALSE;
204:   return(0);
205: }

207: SLEPC_EXTERN 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->synchronize   = DSSynchronize_GHEP;
216:   ds->ops->hermitian     = DSHermitian_GHEP;
217:   return(0);
218: }