Actual source code: petsc-interface.c
slepc-3.15.2 2021-09-20
1: /* @@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ */
2: /* @@@ BLOPEX (version 1.1) LGPL Version 2.1 or above.See www.gnu.org. */
3: /* @@@ Copyright 2010 BLOPEX team https://github.com/lobpcg/blopex */
4: /* @@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ */
5: /* This code was developed by Merico Argentati, Andrew Knyazev, Ilya Lashuk and Evgueni Ovtchinnikov */
7: #include <petscvec.h>
8: #include <petscblaslapack.h>
9: #include <interpreter.h>
10: #include <temp_multivector.h>
11: #include <fortran_matrix.h>
13: static PetscRandom LOBPCG_RandomContext = NULL;
15: #if !defined(PETSC_USE_COMPLEX)
16: BlopexInt PETSC_dpotrf_interface (char *uplo,BlopexInt *n,double *a,BlopexInt * lda,BlopexInt *info)
17: {
18: PetscBLASInt n_,lda_,info_;
20: /* type conversion */
21: n_ = *n;
22: lda_ = *lda;
23: info_ = *info;
25: LAPACKpotrf_(uplo,&n_,(PetscScalar*)a,&lda_,&info_);
27: *info = info_;
28: return 0;
29: }
31: BlopexInt PETSC_dsygv_interface (BlopexInt *itype,char *jobz,char *uplo,BlopexInt *n,double *a,BlopexInt *lda,double *b,BlopexInt *ldb,double *w,double *work,BlopexInt *lwork,BlopexInt *info)
32: {
33: PetscBLASInt itype_,n_,lda_,ldb_,lwork_,info_;
35: itype_ = *itype;
36: n_ = *n;
37: lda_ = *lda;
38: ldb_ = *ldb;
39: lwork_ = *lwork;
40: info_ = *info;
42: LAPACKsygv_(&itype_,jobz,uplo,&n_,(PetscScalar*)a,&lda_,(PetscScalar*)b,&ldb_,(PetscScalar*)w,(PetscScalar*)work,&lwork_,&info_);
44: *info = info_;
45: return 0;
46: }
47: #else
48: BlopexInt PETSC_zpotrf_interface (char *uplo,BlopexInt *n,komplex *a,BlopexInt* lda,BlopexInt *info)
49: {
50: PetscBLASInt n_,lda_,info_;
52: /* type conversion */
53: n_ = *n;
54: lda_ = (PetscBLASInt)*lda;
56: LAPACKpotrf_(uplo,&n_,(PetscScalar*)a,&lda_,&info_);
58: *info = info_;
59: return 0;
60: }
62: BlopexInt PETSC_zsygv_interface (BlopexInt *itype,char *jobz,char *uplo,BlopexInt *n,komplex *a,BlopexInt *lda,komplex *b,BlopexInt *ldb,double *w,komplex *work,BlopexInt *lwork,double *rwork,BlopexInt *info)
63: {
64: PetscBLASInt itype_,n_,lda_,ldb_,lwork_,info_;
66: itype_ = *itype;
67: n_ = *n;
68: lda_ = *lda;
69: ldb_ = *ldb;
70: lwork_ = *lwork;
71: info_ = *info;
73: LAPACKsygv_(&itype_,jobz,uplo,&n_,(PetscScalar*)a,&lda_,(PetscScalar*)b,&ldb_,(PetscReal*)w,(PetscScalar*)work,&lwork_,(PetscReal*)rwork,&info_);
75: *info = info_;
76: return 0;
77: }
78: #endif
80: void *PETSC_MimicVector(void *vvector)
81: {
82: PetscErrorCode ierr;
83: Vec temp;
85: VecDuplicate((Vec)vvector,&temp);CHKERRABORT(PETSC_COMM_SELF,ierr);
86: return (void*)temp;
87: }
89: BlopexInt PETSC_DestroyVector(void *vvector)
90: {
92: Vec v = (Vec)vvector;
94: VecDestroy(&v);
95: return 0;
96: }
98: BlopexInt PETSC_InnerProd(void *x,void *y,void *result)
99: {
102: VecDot((Vec)x,(Vec)y,(PetscScalar*)result);
103: return 0;
104: }
106: BlopexInt PETSC_CopyVector(void *x,void *y)
107: {
108: PetscErrorCode ierr;
110: VecCopy((Vec)x,(Vec)y);
111: return 0;
112: }
114: BlopexInt PETSC_ClearVector(void *x)
115: {
116: PetscErrorCode ierr;
118: VecSet((Vec)x,0.0);
119: return 0;
120: }
122: BlopexInt PETSC_SetRandomValues(void* v,BlopexInt seed)
123: {
126: /* note: without previous call to LOBPCG_InitRandomContext LOBPCG_RandomContext will be null,
127: and VecSetRandom will use internal petsc random context */
129: VecSetRandom((Vec)v,LOBPCG_RandomContext);
130: return 0;
131: }
133: BlopexInt PETSC_ScaleVector(double alpha,void *x)
134: {
137: VecScale((Vec)x,alpha);
138: return 0;
139: }
141: BlopexInt PETSC_Axpy(void *alpha,void *x,void *y)
142: {
145: VecAXPY((Vec)y,*(PetscScalar*)alpha,(Vec)x);
146: return 0;
147: }
149: BlopexInt PETSC_VectorSize(void *x)
150: {
151: PetscInt N;
152: VecGetSize((Vec)x,&N);
153: return N;
154: }
156: int LOBPCG_InitRandomContext(MPI_Comm comm,PetscRandom rand)
157: {
159: /* PetscScalar rnd_bound = 1.0; */
161: if (rand) {
162: PetscObjectReference((PetscObject)rand);
163: PetscRandomDestroy(&LOBPCG_RandomContext);
164: LOBPCG_RandomContext = rand;
165: } else {
166: PetscRandomCreate(comm,&LOBPCG_RandomContext);
167: }
168: return 0;
169: }
171: int LOBPCG_SetFromOptionsRandomContext(void)
172: {
174: PetscRandomSetFromOptions(LOBPCG_RandomContext);
176: #if defined(PETSC_USE_COMPLEX)
177: PetscRandomSetInterval(LOBPCG_RandomContext,(PetscScalar)PetscCMPLX(-1.0,-1.0),(PetscScalar)PetscCMPLX(1.0,1.0));
178: #else
179: PetscRandomSetInterval(LOBPCG_RandomContext,(PetscScalar)-1.0,(PetscScalar)1.0);
180: #endif
181: return 0;
182: }
184: int LOBPCG_DestroyRandomContext(void)
185: {
188: PetscRandomDestroy(&LOBPCG_RandomContext);
189: return 0;
190: }
192: int PETSCSetupInterpreter(mv_InterfaceInterpreter *i)
193: {
194: i->CreateVector = PETSC_MimicVector;
195: i->DestroyVector = PETSC_DestroyVector;
196: i->InnerProd = PETSC_InnerProd;
197: i->CopyVector = PETSC_CopyVector;
198: i->ClearVector = PETSC_ClearVector;
199: i->SetRandomValues = PETSC_SetRandomValues;
200: i->ScaleVector = PETSC_ScaleVector;
201: i->Axpy = PETSC_Axpy;
202: i->VectorSize = PETSC_VectorSize;
204: /* Multivector part */
206: i->CreateMultiVector = mv_TempMultiVectorCreateFromSampleVector;
207: i->CopyCreateMultiVector = mv_TempMultiVectorCreateCopy;
208: i->DestroyMultiVector = mv_TempMultiVectorDestroy;
210: i->Width = mv_TempMultiVectorWidth;
211: i->Height = mv_TempMultiVectorHeight;
212: i->SetMask = mv_TempMultiVectorSetMask;
213: i->CopyMultiVector = mv_TempMultiVectorCopy;
214: i->ClearMultiVector = mv_TempMultiVectorClear;
215: i->SetRandomVectors = mv_TempMultiVectorSetRandom;
216: i->Eval = mv_TempMultiVectorEval;
218: #if defined(PETSC_USE_COMPLEX)
219: i->MultiInnerProd = mv_TempMultiVectorByMultiVector_complex;
220: i->MultiInnerProdDiag = mv_TempMultiVectorByMultiVectorDiag_complex;
221: i->MultiVecMat = mv_TempMultiVectorByMatrix_complex;
222: i->MultiVecMatDiag = mv_TempMultiVectorByDiagonal_complex;
223: i->MultiAxpy = mv_TempMultiVectorAxpy_complex;
224: i->MultiXapy = mv_TempMultiVectorXapy_complex;
225: #else
226: i->MultiInnerProd = mv_TempMultiVectorByMultiVector;
227: i->MultiInnerProdDiag = mv_TempMultiVectorByMultiVectorDiag;
228: i->MultiVecMat = mv_TempMultiVectorByMatrix;
229: i->MultiVecMatDiag = mv_TempMultiVectorByDiagonal;
230: i->MultiAxpy = mv_TempMultiVectorAxpy;
231: i->MultiXapy = mv_TempMultiVectorXapy;
232: #endif
234: return 0;
235: }