Actual source code: dsghep.c
slepc-3.22.1 2024-10-28
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: static PetscErrorCode DSAllocate_GHEP(DS ds,PetscInt ld)
15: {
16: PetscFunctionBegin;
17: PetscCall(DSAllocateMat_Private(ds,DS_MAT_A));
18: PetscCall(DSAllocateMat_Private(ds,DS_MAT_B));
19: PetscCall(DSAllocateMat_Private(ds,DS_MAT_Q));
20: PetscCall(PetscFree(ds->perm));
21: PetscCall(PetscMalloc1(ld,&ds->perm));
22: PetscFunctionReturn(PETSC_SUCCESS);
23: }
25: static PetscErrorCode DSView_GHEP(DS ds,PetscViewer viewer)
26: {
27: PetscViewerFormat format;
29: PetscFunctionBegin;
30: PetscCall(PetscViewerGetFormat(viewer,&format));
31: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(PETSC_SUCCESS);
32: PetscCall(DSViewMat(ds,viewer,DS_MAT_A));
33: PetscCall(DSViewMat(ds,viewer,DS_MAT_B));
34: if (ds->state>DS_STATE_INTERMEDIATE) PetscCall(DSViewMat(ds,viewer,DS_MAT_Q));
35: if (ds->omat[DS_MAT_X]) PetscCall(DSViewMat(ds,viewer,DS_MAT_X));
36: PetscFunctionReturn(PETSC_SUCCESS);
37: }
39: static PetscErrorCode DSVectors_GHEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
40: {
41: PetscScalar *Z;
42: const PetscScalar *Q;
43: PetscInt ld = ds->ld;
45: PetscFunctionBegin;
46: PetscCheck(!rnorm,PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
47: switch (mat) {
48: case DS_MAT_X:
49: case DS_MAT_Y:
50: if (j) {
51: PetscCall(MatDenseGetArray(ds->omat[mat],&Z));
52: if (ds->state>=DS_STATE_CONDENSED) {
53: PetscCall(MatDenseGetArrayRead(ds->omat[DS_MAT_Q],&Q));
54: PetscCall(PetscArraycpy(Z+(*j)*ld,Q+(*j)*ld,ld));
55: PetscCall(MatDenseRestoreArrayRead(ds->omat[DS_MAT_Q],&Q));
56: } else {
57: PetscCall(PetscArrayzero(Z+(*j)*ld,ld));
58: Z[(*j)+(*j)*ld] = 1.0;
59: }
60: PetscCall(MatDenseRestoreArray(ds->omat[mat],&Z));
61: } else {
62: if (ds->state>=DS_STATE_CONDENSED) PetscCall(MatCopy(ds->omat[DS_MAT_Q],ds->omat[mat],SAME_NONZERO_PATTERN));
63: else PetscCall(DSSetIdentity(ds,mat));
64: }
65: break;
66: case DS_MAT_U:
67: case DS_MAT_V:
68: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
69: default:
70: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
71: }
72: PetscFunctionReturn(PETSC_SUCCESS);
73: }
75: static PetscErrorCode DSSort_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *k)
76: {
77: PetscInt n,l,i,*perm,ld=ds->ld;
78: PetscScalar *A;
80: PetscFunctionBegin;
81: if (!ds->sc) PetscFunctionReturn(PETSC_SUCCESS);
82: n = ds->n;
83: l = ds->l;
84: PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A));
85: perm = ds->perm;
86: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
87: if (rr) PetscCall(DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE));
88: else PetscCall(DSSortEigenvalues_Private(ds,wr,NULL,perm,PETSC_FALSE));
89: for (i=l;i<n;i++) A[i+i*ld] = wr[perm[i]];
90: for (i=l;i<n;i++) wr[i] = A[i+i*ld];
91: PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A));
92: PetscCall(DSPermuteColumns_Private(ds,l,n,n,DS_MAT_Q,perm));
93: PetscFunctionReturn(PETSC_SUCCESS);
94: }
96: static PetscErrorCode DSSolve_GHEP(DS ds,PetscScalar *wr,PetscScalar *wi)
97: {
98: PetscScalar *work,*A,*B,*Q;
99: PetscBLASInt itype = 1,*iwork,info,n1,liwork,ld,lrwork=0,lwork;
100: PetscInt off,i;
101: #if defined(PETSC_USE_COMPLEX)
102: PetscReal *rwork,*rr;
103: #endif
105: PetscFunctionBegin;
106: PetscCall(PetscBLASIntCast(ds->n-ds->l,&n1));
107: PetscCall(PetscBLASIntCast(ds->ld,&ld));
108: PetscCall(PetscBLASIntCast(5*ds->n+3,&liwork));
109: #if defined(PETSC_USE_COMPLEX)
110: PetscCall(PetscBLASIntCast(ds->n*ds->n+2*ds->n,&lwork));
111: PetscCall(PetscBLASIntCast(2*ds->n*ds->n+5*ds->n+1+n1,&lrwork));
112: #else
113: PetscCall(PetscBLASIntCast(2*ds->n*ds->n+6*ds->n+1,&lwork));
114: #endif
115: PetscCall(DSAllocateWork_Private(ds,lwork,lrwork,liwork));
116: work = ds->work;
117: iwork = ds->iwork;
118: off = ds->l+ds->l*ld;
119: PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A));
120: PetscCall(MatDenseGetArray(ds->omat[DS_MAT_B],&B));
121: PetscCall(MatDenseGetArray(ds->omat[DS_MAT_Q],&Q));
122: #if defined(PETSC_USE_COMPLEX)
123: rr = ds->rwork;
124: rwork = ds->rwork + n1;
125: lrwork = ds->lrwork - n1;
126: PetscCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,rr,work,&lwork,rwork,&lrwork,iwork,&liwork,&info));
127: for (i=0;i<n1;i++) wr[ds->l+i] = rr[i];
128: #else
129: PetscCallBLAS("LAPACKsygvd",LAPACKsygvd_(&itype,"V","U",&n1,A+off,&ld,B+off,&ld,wr+ds->l,work,&lwork,iwork,&liwork,&info));
130: #endif
131: SlepcCheckLapackInfo("sygvd",info);
132: PetscCall(PetscArrayzero(Q+ds->l*ld,n1*ld));
133: for (i=ds->l;i<ds->n;i++) PetscCall(PetscArraycpy(Q+ds->l+i*ld,A+ds->l+i*ld,n1));
134: PetscCall(PetscArrayzero(B+ds->l*ld,n1*ld));
135: PetscCall(PetscArrayzero(A+ds->l*ld,n1*ld));
136: for (i=ds->l;i<ds->n;i++) {
137: if (wi) wi[i] = 0.0;
138: B[i+i*ld] = 1.0;
139: A[i+i*ld] = wr[i];
140: }
141: PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A));
142: PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_B],&B));
143: PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q));
144: PetscFunctionReturn(PETSC_SUCCESS);
145: }
147: #if !defined(PETSC_HAVE_MPIUNI)
148: static PetscErrorCode DSSynchronize_GHEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
149: {
150: PetscScalar *A,*B,*Q;
151: PetscInt ld=ds->ld,l=ds->l,k;
152: PetscMPIInt n,rank,off=0,size,ldn;
154: PetscFunctionBegin;
155: k = 2*(ds->n-l)*ld;
156: if (ds->state>DS_STATE_RAW) k += (ds->n-l)*ld;
157: if (eigr) k += (ds->n-l);
158: PetscCall(DSAllocateWork_Private(ds,k,0,0));
159: PetscCall(PetscMPIIntCast(k*sizeof(PetscScalar),&size));
160: PetscCall(PetscMPIIntCast(ds->n-l,&n));
161: PetscCall(PetscMPIIntCast(ld*(ds->n-l),&ldn));
162: PetscCall(MatDenseGetArray(ds->omat[DS_MAT_A],&A));
163: PetscCall(MatDenseGetArray(ds->omat[DS_MAT_B],&B));
164: if (ds->state>DS_STATE_RAW) PetscCall(MatDenseGetArray(ds->omat[DS_MAT_Q],&Q));
165: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank));
166: if (!rank) {
167: PetscCallMPI(MPI_Pack(A+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds)));
168: PetscCallMPI(MPI_Pack(B+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds)));
169: if (ds->state>DS_STATE_RAW) PetscCallMPI(MPI_Pack(Q+l*ld,ldn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds)));
170: if (eigr) PetscCallMPI(MPI_Pack(eigr+l,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds)));
171: }
172: PetscCallMPI(MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds)));
173: if (rank) {
174: PetscCallMPI(MPI_Unpack(ds->work,size,&off,A+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds)));
175: PetscCallMPI(MPI_Unpack(ds->work,size,&off,B+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds)));
176: if (ds->state>DS_STATE_RAW) PetscCallMPI(MPI_Unpack(ds->work,size,&off,Q+l*ld,ldn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds)));
177: if (eigr) PetscCallMPI(MPI_Unpack(ds->work,size,&off,eigr+l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds)));
178: }
179: PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_A],&A));
180: PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_B],&B));
181: if (ds->state>DS_STATE_RAW) PetscCall(MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q));
182: PetscFunctionReturn(PETSC_SUCCESS);
183: }
184: #endif
186: static PetscErrorCode DSHermitian_GHEP(DS ds,DSMatType m,PetscBool *flg)
187: {
188: PetscFunctionBegin;
189: if (m==DS_MAT_A || m==DS_MAT_B) *flg = PETSC_TRUE;
190: else *flg = PETSC_FALSE;
191: PetscFunctionReturn(PETSC_SUCCESS);
192: }
194: /*MC
195: DSGHEP - Dense Generalized Hermitian Eigenvalue Problem.
197: Level: beginner
199: Notes:
200: The problem is expressed as A*X = B*X*Lambda, where both A and B are
201: real symmetric (or complex Hermitian) and B is positive-definite. Lambda
202: is a diagonal matrix whose diagonal elements are the arguments of DSSolve().
203: After solve, A is overwritten with Lambda, and B is overwritten with I.
205: No intermediate state is implemented, nor compact storage.
207: Used DS matrices:
208: + DS_MAT_A - first problem matrix
209: . DS_MAT_B - second problem matrix
210: - DS_MAT_Q - matrix of B-orthogonal eigenvectors, which is equal to X
212: Implemented methods:
213: . 0 - Divide and Conquer (_sygvd)
215: .seealso: DSCreate(), DSSetType(), DSType
216: M*/
217: SLEPC_EXTERN PetscErrorCode DSCreate_GHEP(DS ds)
218: {
219: PetscFunctionBegin;
220: ds->ops->allocate = DSAllocate_GHEP;
221: ds->ops->view = DSView_GHEP;
222: ds->ops->vectors = DSVectors_GHEP;
223: ds->ops->solve[0] = DSSolve_GHEP;
224: ds->ops->sort = DSSort_GHEP;
225: #if !defined(PETSC_HAVE_MPIUNI)
226: ds->ops->synchronize = DSSynchronize_GHEP;
227: #endif
228: ds->ops->hermitian = DSHermitian_GHEP;
229: PetscFunctionReturn(PETSC_SUCCESS);
230: }