| Line | Branch | Exec | Source |
|---|---|---|---|
| 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 */ | ||
| 6 | |||
| 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" | ||
| 13 | |||
| 14 | #if !defined(PETSC_USE_COMPLEX) | ||
| 15 | 1618 | BlopexInt PETSC_dpotrf_interface (char *uplo,BlopexInt *n,double *a,BlopexInt * lda,BlopexInt *info) | |
| 16 | { | ||
| 17 | 1618 | PetscBLASInt n_,lda_,info_; | |
| 18 | |||
| 19 | /* type conversion */ | ||
| 20 | 1618 | n_ = *n; | |
| 21 | 1618 | lda_ = *lda; | |
| 22 | 1618 | info_ = *info; | |
| 23 | |||
| 24 | 1618 | LAPACKpotrf_(uplo,&n_,(PetscScalar*)a,&lda_,&info_); | |
| 25 | |||
| 26 | 1618 | *info = info_; | |
| 27 | 1618 | return 0; | |
| 28 | } | ||
| 29 | |||
| 30 | 833 | 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 | 833 | PetscBLASInt itype_,n_,lda_,ldb_,lwork_,info_; | |
| 33 | |||
| 34 | 833 | itype_ = *itype; | |
| 35 | 833 | n_ = *n; | |
| 36 | 833 | lda_ = *lda; | |
| 37 | 833 | ldb_ = *ldb; | |
| 38 | 833 | lwork_ = *lwork; | |
| 39 | 833 | info_ = *info; | |
| 40 | |||
| 41 | 833 | LAPACKsygv_(&itype_,jobz,uplo,&n_,(PetscScalar*)a,&lda_,(PetscScalar*)b,&ldb_,(PetscScalar*)w,(PetscScalar*)work,&lwork_,&info_); | |
| 42 | |||
| 43 | 833 | *info = info_; | |
| 44 | 833 | return 0; | |
| 45 | } | ||
| 46 | #else | ||
| 47 | 1554 | BlopexInt PETSC_zpotrf_interface (char *uplo,BlopexInt *n,komplex *a,BlopexInt* lda,BlopexInt *info) | |
| 48 | { | ||
| 49 | 1554 | PetscBLASInt n_,lda_,info_; | |
| 50 | |||
| 51 | /* type conversion */ | ||
| 52 | 1554 | n_ = *n; | |
| 53 | 1554 | lda_ = (PetscBLASInt)*lda; | |
| 54 | |||
| 55 | 1554 | LAPACKpotrf_(uplo,&n_,(PetscScalar*)a,&lda_,&info_); | |
| 56 | |||
| 57 | 1554 | *info = info_; | |
| 58 | 1554 | return 0; | |
| 59 | } | ||
| 60 | |||
| 61 | 801 | 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 | 801 | PetscBLASInt itype_,n_,lda_,ldb_,lwork_,info_; | |
| 64 | |||
| 65 | 801 | itype_ = *itype; | |
| 66 | 801 | n_ = *n; | |
| 67 | 801 | lda_ = *lda; | |
| 68 | 801 | ldb_ = *ldb; | |
| 69 | 801 | lwork_ = *lwork; | |
| 70 | 801 | info_ = *info; | |
| 71 | |||
| 72 | 801 | LAPACKsygv_(&itype_,jobz,uplo,&n_,(PetscScalar*)a,&lda_,(PetscScalar*)b,&ldb_,(PetscReal*)w,(PetscScalar*)work,&lwork_,(PetscReal*)rwork,&info_); | |
| 73 | |||
| 74 | 801 | *info = info_; | |
| 75 | 801 | return 0; | |
| 76 | } | ||
| 77 | #endif | ||
| 78 | |||
| 79 | 1482 | static void *PETSC_MimicVector(void *vvector) | |
| 80 | { | ||
| 81 | 1482 | Vec temp; | |
| 82 | |||
| 83 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1482 | PetscCallAbort(PETSC_COMM_SELF,VecDuplicate((Vec)vvector,&temp)); |
| 84 | 1482 | return (void*)temp; | |
| 85 | } | ||
| 86 | |||
| 87 | 1740 | static BlopexInt PETSC_DestroyVector(void *vvector) | |
| 88 | { | ||
| 89 | 1740 | Vec v = (Vec)vvector; | |
| 90 | |||
| 91 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1740 | PetscCall(VecDestroy(&v)); |
| 92 | return 0; | ||
| 93 | } | ||
| 94 | |||
| 95 | 182906 | static BlopexInt PETSC_InnerProd(void *x,void *y,void *result) | |
| 96 | { | ||
| 97 | |||
| 98 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
182906 | PetscCall(VecDot((Vec)x,(Vec)y,(PetscScalar*)result)); |
| 99 | return 0; | ||
| 100 | } | ||
| 101 | |||
| 102 | 56324 | static BlopexInt PETSC_CopyVector(void *x,void *y) | |
| 103 | { | ||
| 104 | |||
| 105 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
56324 | PetscCall(VecCopy((Vec)x,(Vec)y)); |
| 106 | return 0; | ||
| 107 | } | ||
| 108 | |||
| 109 | 74522 | static BlopexInt PETSC_ClearVector(void *x) | |
| 110 | { | ||
| 111 | |||
| 112 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
74522 | PetscCall(VecSet((Vec)x,0.0)); |
| 113 | return 0; | ||
| 114 | } | ||
| 115 | |||
| 116 | 277012 | static BlopexInt PETSC_Axpy(void *alpha,void *x,void *y) | |
| 117 | { | ||
| 118 | |||
| 119 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
277012 | PetscCall(VecAXPY((Vec)y,*(PetscScalar*)alpha,(Vec)x)); |
| 120 | return 0; | ||
| 121 | } | ||
| 122 | |||
| 123 | 42 | int PETSCSetupInterpreter(mv_InterfaceInterpreter *i) | |
| 124 | { | ||
| 125 | 42 | i->CreateVector = PETSC_MimicVector; | |
| 126 | 42 | i->DestroyVector = PETSC_DestroyVector; | |
| 127 | 42 | i->InnerProd = PETSC_InnerProd; | |
| 128 | 42 | i->CopyVector = PETSC_CopyVector; | |
| 129 | 42 | i->ClearVector = PETSC_ClearVector; | |
| 130 | 42 | i->Axpy = PETSC_Axpy; | |
| 131 | |||
| 132 | /* Multivector part */ | ||
| 133 | |||
| 134 | 42 | i->CreateMultiVector = mv_TempMultiVectorCreateFromSampleVector; | |
| 135 | 42 | i->CopyCreateMultiVector = mv_TempMultiVectorCreateCopy; | |
| 136 | 42 | i->DestroyMultiVector = mv_TempMultiVectorDestroy; | |
| 137 | |||
| 138 | 42 | i->Width = mv_TempMultiVectorWidth; | |
| 139 | 42 | i->Height = mv_TempMultiVectorHeight; | |
| 140 | 42 | i->SetMask = mv_TempMultiVectorSetMask; | |
| 141 | 42 | i->CopyMultiVector = mv_TempMultiVectorCopy; | |
| 142 | 42 | i->ClearMultiVector = mv_TempMultiVectorClear; | |
| 143 | 42 | i->SetRandomVectors = mv_TempMultiVectorSetRandom; | |
| 144 | 42 | i->Eval = mv_TempMultiVectorEval; | |
| 145 | |||
| 146 | #if defined(PETSC_USE_COMPLEX) | ||
| 147 | 21 | i->MultiInnerProd = mv_TempMultiVectorByMultiVector_complex; | |
| 148 | 21 | i->MultiInnerProdDiag = mv_TempMultiVectorByMultiVectorDiag_complex; | |
| 149 | 21 | i->MultiVecMat = mv_TempMultiVectorByMatrix_complex; | |
| 150 | 21 | i->MultiVecMatDiag = mv_TempMultiVectorByDiagonal_complex; | |
| 151 | 21 | i->MultiAxpy = mv_TempMultiVectorAxpy_complex; | |
| 152 | 21 | i->MultiXapy = mv_TempMultiVectorXapy_complex; | |
| 153 | #else | ||
| 154 | 21 | i->MultiInnerProd = mv_TempMultiVectorByMultiVector; | |
| 155 | 21 | i->MultiInnerProdDiag = mv_TempMultiVectorByMultiVectorDiag; | |
| 156 | 21 | i->MultiVecMat = mv_TempMultiVectorByMatrix; | |
| 157 | 21 | i->MultiVecMatDiag = mv_TempMultiVectorByDiagonal; | |
| 158 | 21 | i->MultiAxpy = mv_TempMultiVectorAxpy; | |
| 159 | 21 | i->MultiXapy = mv_TempMultiVectorXapy; | |
| 160 | #endif | ||
| 161 | |||
| 162 | 42 | return 0; | |
| 163 | } | ||
| 164 |