10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-11-19 04:22:32 +01:00

Removed dependency to cblas

This commit is contained in:
Anthony Scemama 2021-04-28 00:19:46 +02:00
parent 83f2f996d1
commit 8da63eeb7d
4 changed files with 26 additions and 19 deletions

View File

@ -254,11 +254,11 @@ subroutine run_pt2_slave_large(thread,iproc,energy)
call omp_unset_lock(global_selection_buffer_lock) call omp_unset_lock(global_selection_buffer_lock)
if ( iproc == 1 ) then if ( iproc == 1 ) then
call omp_set_lock(global_selection_buffer_lock) call omp_set_lock(global_selection_buffer_lock)
call push_pt2_results_async_send(zmq_socket_push, i_generator, pt2_data, global_selection_buffer, task_id, 1,sending) call push_pt2_results_async_send(zmq_socket_push, (/i_generator/), (/pt2_data/), global_selection_buffer, (/task_id/), 1,sending)
global_selection_buffer%cur = 0 global_selection_buffer%cur = 0
call omp_unset_lock(global_selection_buffer_lock) call omp_unset_lock(global_selection_buffer_lock)
else else
call push_pt2_results_async_send(zmq_socket_push, i_generator, pt2_data, b, task_id, 1,sending) call push_pt2_results_async_send(zmq_socket_push, (/i_generator/), (/pt2_data/), b, (/task_id/), 1,sending)
endif endif
call pt2_dealloc(pt2_data) call pt2_dealloc(pt2_data)

View File

@ -46,4 +46,18 @@ module cfunctions
real (kind=C_DOUBLE ),intent(out) :: csftodetmatrix(rowsmax,colsmax) real (kind=C_DOUBLE ),intent(out) :: csftodetmatrix(rowsmax,colsmax)
end subroutine getCSFtoDETTransformationMatrix end subroutine getCSFtoDETTransformationMatrix
end interface end interface
end module cfunctions end module cfunctions
subroutine f_dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC) &
bind(C, name='f_dgemm')
use iso_c_binding
implicit none
character, intent(in), value :: TRANSA, TRANSB
integer, intent(in), value :: M,N,K,LDA,LDB,LDC
double precision, intent(in), value :: ALPHA, BETA
double precision, intent(in) :: A(LDA,*), B(LDB,*)
double precision, intent(out) :: C(LDC,*)
call dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
end subroutine

View File

@ -267,6 +267,7 @@ void callBlasMatxMat(double *A, int rowA, int colA, double *B, int rowB, int col
double alpha = 1.0; double alpha = 1.0;
double beta = 0.0; double beta = 0.0;
int val = 0; int val = 0;
if (transA) val |= 0x1; if (transA) val |= 0x1;
if (transB) val |= 0x2; if (transB) val |= 0x2;
@ -275,15 +276,13 @@ void callBlasMatxMat(double *A, int rowA, int colA, double *B, int rowB, int col
m = rowA; m = rowA;
n = colB; n = colB;
k = colA; k = colA;
cblas_dgemm(CblasRowMajor, CblasNoTrans, CblasNoTrans, f_dgemm('N', 'N', n, m, k, alpha, B, n, A, k, beta, C, n);
m, n, k, alpha, A, k, B, n, beta, C, n);
break; break;
case 1: // transA, notransB case 1: // transA, notransB
m = colA; m = colA;
n = colB; n = colB;
k = rowA; k = rowA;
cblas_dgemm(CblasRowMajor, CblasTrans, CblasNoTrans, f_dgemm('N', 'T', n, m, k, alpha, B, n, A, colA, beta, C, n);
m, n, k, alpha, A, colA, B, n, beta, C, n);
break; break;
case 2: // notransA, transB case 2: // notransA, transB
//m = rowA; //m = rowA;
@ -292,15 +291,13 @@ void callBlasMatxMat(double *A, int rowA, int colA, double *B, int rowB, int col
m = rowA; m = rowA;
n = rowB; n = rowB;
k = colA; k = colA;
cblas_dgemm(CblasRowMajor, CblasNoTrans, CblasTrans, f_dgemm('T', 'N', n, m, k, alpha, B, colB, A, k, beta, C, n);
m, n, k, alpha, A, k, B, colB, beta, C, n);
break; break;
case 3: // transA, transB case 3: // transA, transB
m = colA; m = colA;
n = rowB; n = rowB;
k = rowA; k = rowA;
cblas_dgemm(CblasRowMajor, CblasTrans, CblasTrans, f_dgemm('T', 'T', n, m, k, alpha, B, colB, A, colA, beta, C, n);
m, n, k, alpha, A, colA, B, colB, beta, C, n);
break; break;
default: default:
printf("Impossible !!!!\n"); printf("Impossible !!!!\n");

View File

@ -22,14 +22,10 @@ struct bin_tree {
int NBF; int NBF;
}; };
typedef enum {CblasRowMajor=101, CblasColMajor=102} CBLAS_LAYOUT; void f_dgemm(const char transb, const char transa, const int n, const int m, const int k,
typedef enum {CblasNoTrans=111, CblasTrans=112, CblasConjTrans=113} CBLAS_TRANSPOSE; const double alpha, const double* B, const int ldb, const double* A,
typedef CBLAS_LAYOUT CBLAS_ORDER; const int lda, const double beta, double* C, const int ldc);
void cblas_dgemm(CBLAS_LAYOUT layout, CBLAS_TRANSPOSE TransA,
CBLAS_TRANSPOSE TransB, const int M, const int N,
const int K, const double alpha, const double *A,
const int lda, const double *B, const int ldb,
const double beta, double *C, const int ldc);
#define MAX_SOMO 32 #define MAX_SOMO 32