Actual source code: dsghep.c
slepc-3.11.2 2019-07-30
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) {
39: DSViewMat(ds,viewer,DS_MAT_Q);
40: }
41: if (ds->mat[DS_MAT_X]) {
42: DSViewMat(ds,viewer,DS_MAT_X);
43: }
44: return(0);
45: }
47: PetscErrorCode DSVectors_GHEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
48: {
49: PetscScalar *Q = ds->mat[DS_MAT_Q];
50: PetscInt ld = ds->ld,i;
54: if (rnorm) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
55: switch (mat) {
56: case DS_MAT_X:
57: case DS_MAT_Y:
58: if (j) {
59: if (ds->state>=DS_STATE_CONDENSED) {
60: PetscMemcpy(ds->mat[mat]+(*j)*ld,Q+(*j)*ld,ld*sizeof(PetscScalar));
61: } else {
62: PetscMemzero(ds->mat[mat]+(*j)*ld,ld*sizeof(PetscScalar));
63: *(ds->mat[mat]+(*j)+(*j)*ld) = 1.0;
64: }
65: } else {
66: if (ds->state>=DS_STATE_CONDENSED) {
67: PetscMemcpy(ds->mat[mat],Q,ld*ld*sizeof(PetscScalar));
68: } else {
69: PetscMemzero(ds->mat[mat],ld*ld*sizeof(PetscScalar));
70: for (i=0;i<ds->n;i++) *(ds->mat[mat]+i+i*ld) = 1.0;
71: }
72: }
73: break;
74: case DS_MAT_U:
75: case DS_MAT_VT:
76: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
77: break;
78: default:
79: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
80: }
81: return(0);
82: }
84: PetscErrorCode DSSort_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
85: {
87: PetscInt n,l,i,*perm,ld=ds->ld;
88: PetscScalar *A;
91: if (!ds->sc) return(0);
92: n = ds->n;
93: l = ds->l;
94: A = ds->mat[DS_MAT_A];
95: perm = ds->perm;
96: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
97: if (rr) {
98: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
99: } else {
100: DSSortEigenvalues_Private(ds,wr,NULL,perm,PETSC_FALSE);
101: }
102: for (i=l;i<n;i++) A[i+i*ld] = wr[perm[i]];
103: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
104: DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
105: return(0);
106: }
108: PetscErrorCode DSSolve_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
109: {
110: #if defined(SLEPC_MISSING_LAPACK_SYGVD)
112: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SYGVD - Lapack routine is unavailable");
113: #else
115: PetscScalar *work,*A,*B,*Q;
116: PetscBLASInt itype = 1,*iwork,info,n1,liwork,ld,lrwork=0,lwork;
117: PetscInt off,i;
118: #if defined(PETSC_USE_COMPLEX)
119: PetscReal *rwork,*rr;
120: #endif
123: PetscBLASIntCast(ds->n-ds->l,&n1);
124: PetscBLASIntCast(ds->ld,&ld);
125: PetscBLASIntCast(5*ds->n+3,&liwork);
126: #if defined(PETSC_USE_COMPLEX)
127: PetscBLASIntCast(ds->n*ds->n+2*ds->n,&lwork);
128: PetscBLASIntCast(2*ds->n*ds->n+5*ds->n+1+n1,&lrwork);
129: #else
130: PetscBLASIntCast(2*ds->n*ds->n+6*ds->n+1,&lwork);
131: #endif
132: DSAllocateWork_Private(ds,lwork,lrwork,liwork);
133: work = ds->work;
134: iwork = ds->iwork;
135: off = ds->l+ds->l*ld;
136: A = ds->mat[DS_MAT_A];
137: B = ds->mat[DS_MAT_B];
138: Q = ds->mat[DS_MAT_Q];
139: #if defined(PETSC_USE_COMPLEX)
140: rr = ds->rwork;
141: rwork = ds->rwork + n1;
142: lrwork = ds->lrwork - n1;
143: PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,rr,work,&lwork,rwork,&lrwork,iwork,&liwork,&info));
144: for (i=0;i<n1;i++) wr[ds->l+i] = rr[i];
145: #else
146: PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,wr+ds->l,work,&lwork,iwork,&liwork,&info));
147: #endif
148: SlepcCheckLapackInfo("sygvd",info);
149: PetscMemzero(Q+ds->l*ld,n1*ld*sizeof(PetscScalar));
150: for (i=ds->l;i<ds->n;i++) {
151: PetscMemcpy(Q+ds->l+i*ld,A+ds->l+i*ld,n1*sizeof(PetscScalar));
152: }
153: PetscMemzero(B+ds->l*ld,n1*ld*sizeof(PetscScalar));
154: PetscMemzero(A+ds->l*ld,n1*ld*sizeof(PetscScalar));
155: for (i=ds->l;i<ds->n;i++) {
156: if (wi) wi[i] = 0.0;
157: B[i+i*ld] = 1.0;
158: A[i+i*ld] = wr[i];
159: }
160: return(0);
161: #endif
162: }
164: PetscErrorCode DSSynchronize_GHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
165: {
167: PetscInt ld=ds->ld,l=ds->l,k;
168: PetscMPIInt n,rank,off=0,size,ldn;
171: k = 2*(ds->n-l)*ld;
172: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
173: if (eigr) k += (ds->n-l);
174: if (eigi) k += (ds->n-l);
175: DSAllocateWork_Private(ds,k,0,0);
176: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
177: PetscMPIIntCast(ds->n-l,&n);
178: PetscMPIIntCast(ld*(ds->n-l),&ldn);
179: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
180: if (!rank) {
181: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
182: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
183: if (ds->state>DS_STATE_RAW) {
184: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
185: }
186: if (eigr) {
187: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
188: }
189: if (eigi) {
190: MPI_Pack(eigi+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
191: }
192: }
193: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
194: if (rank) {
195: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
196: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
197: if (ds->state>DS_STATE_RAW) {
198: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
199: }
200: if (eigr) {
201: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
202: }
203: if (eigi) {
204: MPI_Unpack(ds->work,size,&off,eigi+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
205: }
206: }
207: return(0);
208: }
210: PetscErrorCode DSHermitian_GHEP(DS ds,DSMatType m,PetscBool *flg)
211: {
213: if (m==DS_MAT_A || m==DS_MAT_B) *flg = PETSC_TRUE;
214: else *flg = PETSC_FALSE;
215: return(0);
216: }
218: SLEPC_EXTERN PetscErrorCode DSCreate_GHEP(DS ds)
219: {
221: ds->ops->allocate = DSAllocate_GHEP;
222: ds->ops->view = DSView_GHEP;
223: ds->ops->vectors = DSVectors_GHEP;
224: ds->ops->solve[0] = DSSolve_GHEP;
225: ds->ops->sort = DSSort_GHEP;
226: ds->ops->synchronize = DSSynchronize_GHEP;
227: ds->ops->hermitian = DSHermitian_GHEP;
228: return(0);
229: }