Actual source code: dsghep.c
slepc-3.14.2 2021-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: default:
74: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
75: }
76: return(0);
77: }
79: PetscErrorCode DSSort_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
80: {
82: PetscInt n,l,i,*perm,ld=ds->ld;
83: PetscScalar *A;
86: if (!ds->sc) return(0);
87: n = ds->n;
88: l = ds->l;
89: A = ds->mat[DS_MAT_A];
90: perm = ds->perm;
91: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
92: if (rr) {
93: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
94: } else {
95: DSSortEigenvalues_Private(ds,wr,NULL,perm,PETSC_FALSE);
96: }
97: for (i=l;i<n;i++) A[i+i*ld] = wr[perm[i]];
98: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
99: DSPermuteColumns_Private(ds,l,n,DS_MAT_Q,perm);
100: return(0);
101: }
103: PetscErrorCode DSSolve_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
104: {
106: PetscScalar *work,*A,*B,*Q;
107: PetscBLASInt itype = 1,*iwork,info,n1,liwork,ld,lrwork=0,lwork;
108: PetscInt off,i;
109: #if defined(PETSC_USE_COMPLEX)
110: PetscReal *rwork,*rr;
111: #endif
114: PetscBLASIntCast(ds->n-ds->l,&n1);
115: PetscBLASIntCast(ds->ld,&ld);
116: PetscBLASIntCast(5*ds->n+3,&liwork);
117: #if defined(PETSC_USE_COMPLEX)
118: PetscBLASIntCast(ds->n*ds->n+2*ds->n,&lwork);
119: PetscBLASIntCast(2*ds->n*ds->n+5*ds->n+1+n1,&lrwork);
120: #else
121: PetscBLASIntCast(2*ds->n*ds->n+6*ds->n+1,&lwork);
122: #endif
123: DSAllocateWork_Private(ds,lwork,lrwork,liwork);
124: work = ds->work;
125: iwork = ds->iwork;
126: off = ds->l+ds->l*ld;
127: A = ds->mat[DS_MAT_A];
128: B = ds->mat[DS_MAT_B];
129: Q = ds->mat[DS_MAT_Q];
130: #if defined(PETSC_USE_COMPLEX)
131: rr = ds->rwork;
132: rwork = ds->rwork + n1;
133: lrwork = ds->lrwork - n1;
134: PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,rr,work,&lwork,rwork,&lrwork,iwork,&liwork,&info));
135: for (i=0;i<n1;i++) wr[ds->l+i] = rr[i];
136: #else
137: PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,wr+ds->l,work,&lwork,iwork,&liwork,&info));
138: #endif
139: SlepcCheckLapackInfo("sygvd",info);
140: PetscArrayzero(Q+ds->l*ld,n1*ld);
141: for (i=ds->l;i<ds->n;i++) {
142: PetscArraycpy(Q+ds->l+i*ld,A+ds->l+i*ld,n1);
143: }
144: PetscArrayzero(B+ds->l*ld,n1*ld);
145: PetscArrayzero(A+ds->l*ld,n1*ld);
146: for (i=ds->l;i<ds->n;i++) {
147: if (wi) wi[i] = 0.0;
148: B[i+i*ld] = 1.0;
149: A[i+i*ld] = wr[i];
150: }
151: return(0);
152: }
154: PetscErrorCode DSSynchronize_GHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
155: {
157: PetscInt ld=ds->ld,l=ds->l,k;
158: PetscMPIInt n,rank,off=0,size,ldn;
161: k = 2*(ds->n-l)*ld;
162: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
163: if (eigr) k += (ds->n-l);
164: DSAllocateWork_Private(ds,k,0,0);
165: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
166: PetscMPIIntCast(ds->n-l,&n);
167: PetscMPIIntCast(ld*(ds->n-l),&ldn);
168: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
169: if (!rank) {
170: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
171: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
172: if (ds->state>DS_STATE_RAW) {
173: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
174: }
175: if (eigr) {
176: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
177: }
178: }
179: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
180: if (rank) {
181: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
182: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
183: if (ds->state>DS_STATE_RAW) {
184: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
185: }
186: if (eigr) {
187: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
188: }
189: }
190: return(0);
191: }
193: PetscErrorCode DSHermitian_GHEP(DS ds,DSMatType m,PetscBool *flg)
194: {
196: if (m==DS_MAT_A || m==DS_MAT_B) *flg = PETSC_TRUE;
197: else *flg = PETSC_FALSE;
198: return(0);
199: }
201: SLEPC_EXTERN PetscErrorCode DSCreate_GHEP(DS ds)
202: {
204: ds->ops->allocate = DSAllocate_GHEP;
205: ds->ops->view = DSView_GHEP;
206: ds->ops->vectors = DSVectors_GHEP;
207: ds->ops->solve[0] = DSSolve_GHEP;
208: ds->ops->sort = DSSort_GHEP;
209: ds->ops->synchronize = DSSynchronize_GHEP;
210: ds->ops->hermitian = DSHermitian_GHEP;
211: return(0);
212: }