Actual source code: dsghep.c
slepc-3.16.2 2022-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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;
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: DSSetIdentity(ds,mat);
66: }
67: }
68: break;
69: case DS_MAT_U:
70: case DS_MAT_V:
71: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
72: default:
73: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
74: }
75: return(0);
76: }
78: PetscErrorCode DSSort_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
79: {
81: PetscInt n,l,i,*perm,ld=ds->ld;
82: PetscScalar *A;
85: if (!ds->sc) return(0);
86: n = ds->n;
87: l = ds->l;
88: A = ds->mat[DS_MAT_A];
89: perm = ds->perm;
90: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
91: if (rr) {
92: DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
93: } else {
94: DSSortEigenvalues_Private(ds,wr,NULL,perm,PETSC_FALSE);
95: }
96: for (i=l;i<n;i++) A[i+i*ld] = wr[perm[i]];
97: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
98: DSPermuteColumns_Private(ds,l,n,n,DS_MAT_Q,perm);
99: return(0);
100: }
102: PetscErrorCode DSSolve_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
103: {
105: PetscScalar *work,*A,*B,*Q;
106: PetscBLASInt itype = 1,*iwork,info,n1,liwork,ld,lrwork=0,lwork;
107: PetscInt off,i;
108: #if defined(PETSC_USE_COMPLEX)
109: PetscReal *rwork,*rr;
110: #endif
113: PetscBLASIntCast(ds->n-ds->l,&n1);
114: PetscBLASIntCast(ds->ld,&ld);
115: PetscBLASIntCast(5*ds->n+3,&liwork);
116: #if defined(PETSC_USE_COMPLEX)
117: PetscBLASIntCast(ds->n*ds->n+2*ds->n,&lwork);
118: PetscBLASIntCast(2*ds->n*ds->n+5*ds->n+1+n1,&lrwork);
119: #else
120: PetscBLASIntCast(2*ds->n*ds->n+6*ds->n+1,&lwork);
121: #endif
122: DSAllocateWork_Private(ds,lwork,lrwork,liwork);
123: work = ds->work;
124: iwork = ds->iwork;
125: off = ds->l+ds->l*ld;
126: A = ds->mat[DS_MAT_A];
127: B = ds->mat[DS_MAT_B];
128: Q = ds->mat[DS_MAT_Q];
129: #if defined(PETSC_USE_COMPLEX)
130: rr = ds->rwork;
131: rwork = ds->rwork + n1;
132: lrwork = ds->lrwork - n1;
133: PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,rr,work,&lwork,rwork,&lrwork,iwork,&liwork,&info));
134: for (i=0;i<n1;i++) wr[ds->l+i] = rr[i];
135: #else
136: PetscStackCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,wr+ds->l,work,&lwork,iwork,&liwork,&info));
137: #endif
138: SlepcCheckLapackInfo("sygvd",info);
139: PetscArrayzero(Q+ds->l*ld,n1*ld);
140: for (i=ds->l;i<ds->n;i++) {
141: PetscArraycpy(Q+ds->l+i*ld,A+ds->l+i*ld,n1);
142: }
143: PetscArrayzero(B+ds->l*ld,n1*ld);
144: PetscArrayzero(A+ds->l*ld,n1*ld);
145: for (i=ds->l;i<ds->n;i++) {
146: if (wi) wi[i] = 0.0;
147: B[i+i*ld] = 1.0;
148: A[i+i*ld] = wr[i];
149: }
150: return(0);
151: }
153: PetscErrorCode DSSynchronize_GHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
154: {
156: PetscInt ld=ds->ld,l=ds->l,k;
157: PetscMPIInt n,rank,off=0,size,ldn;
160: k = 2*(ds->n-l)*ld;
161: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
162: if (eigr) k += (ds->n-l);
163: DSAllocateWork_Private(ds,k,0,0);
164: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
165: PetscMPIIntCast(ds->n-l,&n);
166: PetscMPIIntCast(ld*(ds->n-l),&ldn);
167: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
168: if (!rank) {
169: MPI_Pack(ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
170: MPI_Pack(ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
171: if (ds->state>DS_STATE_RAW) {
172: MPI_Pack(ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
173: }
174: if (eigr) {
175: MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
176: }
177: }
178: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
179: if (rank) {
180: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_A]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
181: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_B]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
182: if (ds->state>DS_STATE_RAW) {
183: MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Q]+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
184: }
185: if (eigr) {
186: MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
187: }
188: }
189: return(0);
190: }
192: PetscErrorCode DSHermitian_GHEP(DS ds,DSMatType m,PetscBool *flg)
193: {
195: if (m==DS_MAT_A || m==DS_MAT_B) *flg = PETSC_TRUE;
196: else *flg = PETSC_FALSE;
197: return(0);
198: }
200: /*MC
201: DSGHEP - Dense Generalized Hermitian Eigenvalue Problem.
203: Level: beginner
205: Notes:
206: The problem is expressed as A*X = B*X*Lambda, where both A and B are
207: real symmetric (or complex Hermitian) and B is positive-definite. Lambda
208: is a diagonal matrix whose diagonal elements are the arguments of DSSolve().
209: After solve, A is overwritten with Lambda, and B is overwritten with I.
211: No intermediate state is implemented, nor compact storage.
213: Used DS matrices:
214: + DS_MAT_A - first problem matrix
215: . DS_MAT_B - second problem matrix
216: - DS_MAT_Q - matrix of B-orthogonal eigenvectors, which is equal to X
218: Implemented methods:
219: . 0 - Divide and Conquer (_sygvd)
221: .seealso: DSCreate(), DSSetType(), DSType
222: M*/
223: SLEPC_EXTERN PetscErrorCode DSCreate_GHEP(DS ds)
224: {
226: ds->ops->allocate = DSAllocate_GHEP;
227: ds->ops->view = DSView_GHEP;
228: ds->ops->vectors = DSVectors_GHEP;
229: ds->ops->solve[0] = DSSolve_GHEP;
230: ds->ops->sort = DSSort_GHEP;
231: ds->ops->synchronize = DSSynchronize_GHEP;
232: ds->ops->hermitian = DSHermitian_GHEP;
233: return(0);
234: }