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 https://github.com/lobpcg/blopex     */
  4: /* @@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@ */
  5: /* This code was developed by Merico Argentati, Andrew Knyazev, Ilya Lashuk and Evgueni Ovtchinnikov */

  7: #include <slepcsys.h>
  8: #include <petscblaslapack.h>
  9: #include <interpreter.h>
 10: #include <temp_multivector.h>
 11: #include <fortran_matrix.h>
 12: #include "petsc-interface.h"

 14: #if !defined(PETSC_USE_COMPLEX)
 15: BlopexInt PETSC_dpotrf_interface (char *uplo,BlopexInt *n,double *a,BlopexInt * lda,BlopexInt *info)
 16: {
 17:   PetscBLASInt n_,lda_,info_;

 19:   /* type conversion */
 20:   n_ = *n;
 21:   lda_ = *lda;
 22:   info_ = *info;

 24:   LAPACKpotrf_(uplo,&n_,(PetscScalar*)a,&lda_,&info_);

 26:   *info = info_;
 27:   return 0;
 28: }

 30: 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)
 31: {
 32:   PetscBLASInt itype_,n_,lda_,ldb_,lwork_,info_;

 34:   itype_ = *itype;
 35:   n_ = *n;
 36:   lda_ = *lda;
 37:   ldb_ = *ldb;
 38:   lwork_ = *lwork;
 39:   info_ = *info;

 41:   LAPACKsygv_(&itype_,jobz,uplo,&n_,(PetscScalar*)a,&lda_,(PetscScalar*)b,&ldb_,(PetscScalar*)w,(PetscScalar*)work,&lwork_,&info_);

 43:   *info = info_;
 44:   return 0;
 45: }
 46: #else
 47: BlopexInt PETSC_zpotrf_interface (char *uplo,BlopexInt *n,komplex *a,BlopexInt* lda,BlopexInt *info)
 48: {
 49:   PetscBLASInt n_,lda_,info_;

 51:   /* type conversion */
 52:   n_ = *n;
 53:   lda_ = (PetscBLASInt)*lda;

 55:   LAPACKpotrf_(uplo,&n_,(PetscScalar*)a,&lda_,&info_);

 57:   *info = info_;
 58:   return 0;
 59: }

 61: 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)
 62: {
 63:   PetscBLASInt itype_,n_,lda_,ldb_,lwork_,info_;

 65:   itype_ = *itype;
 66:   n_ = *n;
 67:   lda_ = *lda;
 68:   ldb_ = *ldb;
 69:   lwork_ = *lwork;
 70:   info_ = *info;

 72:   LAPACKsygv_(&itype_,jobz,uplo,&n_,(PetscScalar*)a,&lda_,(PetscScalar*)b,&ldb_,(PetscReal*)w,(PetscScalar*)work,&lwork_,(PetscReal*)rwork,&info_);

 74:   *info = info_;
 75:   return 0;
 76: }
 77: #endif

 79: static void *PETSC_MimicVector(void *vvector)
 80: {
 81:   Vec temp;

 83:   PetscCallAbort(PETSC_COMM_SELF,VecDuplicate((Vec)vvector,&temp));
 84:   return (void*)temp;
 85: }

 87: static BlopexInt PETSC_DestroyVector(void *vvector)
 88: {
 89:   Vec v = (Vec)vvector;

 91:   PetscCall(VecDestroy(&v));
 92:   return 0;
 93: }

 95: static BlopexInt PETSC_InnerProd(void *x,void *y,void *result)
 96: {

 98:   PetscCall(VecDot((Vec)x,(Vec)y,(PetscScalar*)result));
 99:   return 0;
100: }

102: static BlopexInt PETSC_CopyVector(void *x,void *y)
103: {

105:   PetscCall(VecCopy((Vec)x,(Vec)y));
106:   return 0;
107: }

109: static BlopexInt PETSC_ClearVector(void *x)
110: {

112:   PetscCall(VecSet((Vec)x,0.0));
113:   return 0;
114: }

116: static BlopexInt PETSC_Axpy(void *alpha,void *x,void *y)
117: {

119:   PetscCall(VecAXPY((Vec)y,*(PetscScalar*)alpha,(Vec)x));
120:   return 0;
121: }

123: int PETSCSetupInterpreter(mv_InterfaceInterpreter *i)
124: {
125:   i->CreateVector = PETSC_MimicVector;
126:   i->DestroyVector = PETSC_DestroyVector;
127:   i->InnerProd = PETSC_InnerProd;
128:   i->CopyVector = PETSC_CopyVector;
129:   i->ClearVector = PETSC_ClearVector;
130:   i->Axpy = PETSC_Axpy;

132:   /* Multivector part */

134:   i->CreateMultiVector = mv_TempMultiVectorCreateFromSampleVector;
135:   i->CopyCreateMultiVector = mv_TempMultiVectorCreateCopy;
136:   i->DestroyMultiVector = mv_TempMultiVectorDestroy;

138:   i->Width = mv_TempMultiVectorWidth;
139:   i->Height = mv_TempMultiVectorHeight;
140:   i->SetMask = mv_TempMultiVectorSetMask;
141:   i->CopyMultiVector = mv_TempMultiVectorCopy;
142:   i->ClearMultiVector = mv_TempMultiVectorClear;
143:   i->SetRandomVectors = mv_TempMultiVectorSetRandom;
144:   i->Eval = mv_TempMultiVectorEval;

146: #if defined(PETSC_USE_COMPLEX)
147:   i->MultiInnerProd = mv_TempMultiVectorByMultiVector_complex;
148:   i->MultiInnerProdDiag = mv_TempMultiVectorByMultiVectorDiag_complex;
149:   i->MultiVecMat = mv_TempMultiVectorByMatrix_complex;
150:   i->MultiVecMatDiag = mv_TempMultiVectorByDiagonal_complex;
151:   i->MultiAxpy = mv_TempMultiVectorAxpy_complex;
152:   i->MultiXapy = mv_TempMultiVectorXapy_complex;
153: #else
154:   i->MultiInnerProd = mv_TempMultiVectorByMultiVector;
155:   i->MultiInnerProdDiag = mv_TempMultiVectorByMultiVectorDiag;
156:   i->MultiVecMat = mv_TempMultiVectorByMatrix;
157:   i->MultiVecMatDiag = mv_TempMultiVectorByDiagonal;
158:   i->MultiAxpy = mv_TempMultiVectorAxpy;
159:   i->MultiXapy = mv_TempMultiVectorXapy;
160: #endif

162:   return 0;
163: }