Actual source code: petsc-interface.c
1: /* @@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ */
2: /* @@@ BLOPEX (version 1.1) LGPL Version 2.1 or above.See www.gnu.org. */
3: /* @@@ Copyright 2010 BLOPEX team http://code.google.com/p/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 "blopex_interpreter.h"
10: #include "blopex_temp_multivector.h"
12: static PetscRandom LOBPCG_RandomContext = NULL;
14: typedef struct {
15: double real,imag;
16: } komplex;
18: BlopexInt PETSC_dpotrf_interface (char *uplo,BlopexInt *n,double *a,BlopexInt * lda,BlopexInt *info)
19: {
20: PetscBLASInt n_,lda_,info_;
22: /* type conversion */
23: n_ = *n;
24: lda_ = *lda;
25: info_ = *info;
27: LAPACKpotrf_(uplo,&n_,(PetscScalar*)a,&lda_,&info_);
29: *info = info_;
30: return 0;
31: }
33: BlopexInt PETSC_zpotrf_interface (char *uplo,BlopexInt *n,komplex *a,BlopexInt* lda,BlopexInt *info)
34: {
35: PetscBLASInt n_,lda_,info_;
37: /* type conversion */
38: n_ = *n;
39: lda_ = (PetscBLASInt)*lda;
41: LAPACKpotrf_(uplo,&n_,(PetscScalar*)a,&lda_,&info_);
43: *info = info_;
44: return 0;
45: }
47: 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)
48: {
49: #if !defined(PETSC_USE_COMPLEX)
50: PetscBLASInt itype_,n_,lda_,ldb_,lwork_,info_;
52: itype_ = *itype;
53: n_ = *n;
54: lda_ = *lda;
55: ldb_ = *ldb;
56: lwork_ = *lwork;
57: info_ = *info;
59: LAPACKsygv_(&itype_,jobz,uplo,&n_,(PetscScalar*)a,&lda_,(PetscScalar*)b,&ldb_,(PetscScalar*)w,(PetscScalar*)work,&lwork_,&info_);
61: *info = info_;
62: #endif
63: return 0;
64: }
66: 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)
67: {
68: #if defined(PETSC_USE_COMPLEX)
69: PetscBLASInt itype_,n_,lda_,ldb_,lwork_,info_;
71: itype_ = *itype;
72: n_ = *n;
73: lda_ = *lda;
74: ldb_ = *ldb;
75: lwork_ = *lwork;
76: info_ = *info;
78: LAPACKsygv_(&itype_,jobz,uplo,&n_,(PetscScalar*)a,&lda_,(PetscScalar*)b,&ldb_,(PetscReal*)w,(PetscScalar*)work,&lwork_,(PetscReal*)rwork,&info_);
80: *info = info_;
81: #endif
82: return 0;
83: }
85: void *PETSC_MimicVector(void *vvector)
86: {
87: PetscErrorCode ierr;
88: Vec temp;
90: VecDuplicate((Vec)vvector,&temp);CHKERRABORT(PETSC_COMM_SELF,ierr);
91: return (void*)temp;
92: }
94: BlopexInt PETSC_DestroyVector(void *vvector)
95: {
97: Vec v = (Vec)vvector;
99: VecDestroy(&v);
100: return 0;
101: }
103: BlopexInt PETSC_InnerProd(void *x,void *y,void *result)
104: {
107: VecDot((Vec)x,(Vec)y,(PetscScalar*)result);
108: return 0;
109: }
111: BlopexInt PETSC_CopyVector(void *x,void *y)
112: {
113: PetscErrorCode ierr;
115: VecCopy((Vec)x,(Vec)y);
116: return 0;
117: }
119: BlopexInt PETSC_ClearVector(void *x)
120: {
121: PetscErrorCode ierr;
123: VecSet((Vec)x,0.0);
124: return 0;
125: }
127: BlopexInt PETSC_SetRandomValues(void* v,BlopexInt seed)
128: {
131: /* note: without previous call to LOBPCG_InitRandomContext LOBPCG_RandomContext will be null,
132: and VecSetRandom will use internal petsc random context */
134: VecSetRandom((Vec)v,LOBPCG_RandomContext);
135: return 0;
136: }
138: BlopexInt PETSC_ScaleVector(double alpha,void *x)
139: {
142: VecScale((Vec)x,alpha);
143: return 0;
144: }
146: BlopexInt PETSC_Axpy(void *alpha,void *x,void *y)
147: {
150: VecAXPY((Vec)y,*(PetscScalar*)alpha,(Vec)x);
151: return 0;
152: }
154: BlopexInt PETSC_VectorSize(void *x)
155: {
156: PetscInt N;
157: VecGetSize((Vec)x,&N);
158: return N;
159: }
161: int LOBPCG_InitRandomContext(MPI_Comm comm,PetscRandom rand)
162: {
164: /* PetscScalar rnd_bound = 1.0; */
166: if (rand) {
167: PetscObjectReference((PetscObject)rand);
168: PetscRandomDestroy(&LOBPCG_RandomContext);
169: LOBPCG_RandomContext = rand;
170: } else {
171: PetscRandomCreate(comm,&LOBPCG_RandomContext);
172: }
173: return 0;
174: }
176: int LOBPCG_SetFromOptionsRandomContext(void)
177: {
179: PetscRandomSetFromOptions(LOBPCG_RandomContext);
181: #if defined(PETSC_USE_COMPLEX)
182: PetscRandomSetInterval(LOBPCG_RandomContext,(PetscScalar)-1.0-1.0*PETSC_i,(PetscScalar)1.0+1.0*PETSC_i);
183: #else
184: PetscRandomSetInterval(LOBPCG_RandomContext,(PetscScalar)-1.0,(PetscScalar)1.0);
185: #endif
186: return 0;
187: }
189: int LOBPCG_DestroyRandomContext(void)
190: {
193: PetscRandomDestroy(&LOBPCG_RandomContext);
194: return 0;
195: }
197: int PETSCSetupInterpreter(mv_InterfaceInterpreter *i)
198: {
199: i->CreateVector = PETSC_MimicVector;
200: i->DestroyVector = PETSC_DestroyVector;
201: i->InnerProd = PETSC_InnerProd;
202: i->CopyVector = PETSC_CopyVector;
203: i->ClearVector = PETSC_ClearVector;
204: i->SetRandomValues = PETSC_SetRandomValues;
205: i->ScaleVector = PETSC_ScaleVector;
206: i->Axpy = PETSC_Axpy;
207: i->VectorSize = PETSC_VectorSize;
209: /* Multivector part */
211: i->CreateMultiVector = mv_TempMultiVectorCreateFromSampleVector;
212: i->CopyCreateMultiVector = mv_TempMultiVectorCreateCopy;
213: i->DestroyMultiVector = mv_TempMultiVectorDestroy;
215: i->Width = mv_TempMultiVectorWidth;
216: i->Height = mv_TempMultiVectorHeight;
217: i->SetMask = mv_TempMultiVectorSetMask;
218: i->CopyMultiVector = mv_TempMultiVectorCopy;
219: i->ClearMultiVector = mv_TempMultiVectorClear;
220: i->SetRandomVectors = mv_TempMultiVectorSetRandom;
221: i->Eval = mv_TempMultiVectorEval;
223: #if defined(PETSC_USE_COMPLEX)
224: i->MultiInnerProd = mv_TempMultiVectorByMultiVector_complex;
225: i->MultiInnerProdDiag = mv_TempMultiVectorByMultiVectorDiag_complex;
226: i->MultiVecMat = mv_TempMultiVectorByMatrix_complex;
227: i->MultiVecMatDiag = mv_TempMultiVectorByDiagonal_complex;
228: i->MultiAxpy = mv_TempMultiVectorAxpy_complex;
229: i->MultiXapy = mv_TempMultiVectorXapy_complex;
230: #else
231: i->MultiInnerProd = mv_TempMultiVectorByMultiVector;
232: i->MultiInnerProdDiag = mv_TempMultiVectorByMultiVectorDiag;
233: i->MultiVecMat = mv_TempMultiVectorByMatrix;
234: i->MultiVecMatDiag = mv_TempMultiVectorByDiagonal;
235: i->MultiAxpy = mv_TempMultiVectorAxpy;
236: i->MultiXapy = mv_TempMultiVectorXapy;
237: #endif
239: return 0;
240: }