1
0
mirror of https://gitlab.com/scemama/qp_plugins_scemama.git synced 2024-11-07 06:33:40 +01:00

Removed compute_J1

This commit is contained in:
Anthony Scemama 2023-08-02 23:50:28 +02:00
parent 9d5df395aa
commit 2763954c19
2 changed files with 77 additions and 108 deletions

View File

@ -506,9 +506,6 @@ subroutine compute_r2_space_chol(nO,nV,t1,t2,tau,H_oo,H_vv,H_vo,r2,max_r2,gpu_da
double precision, allocatable :: X_vovo(:,:,:,:), Y_oovo(:,:,:,:)
call compute_J1_chol(nO,nV,t1,t2,cc_space_v_ovvo,cc_space_v_ovoo, &
cc_space_v_vvoo,J1)
double precision, allocatable :: K1(:,:,:,:)
allocate(K1(nO,nV,nO,nV))
call compute_K1_chol(nO,nV,t1,t2,cc_space_v_ovoo,cc_space_v_vvoo, &
@ -724,92 +721,6 @@ subroutine compute_r2_space_chol(nO,nV,t1,t2,tau,H_oo,H_vv,H_vo,r2,max_r2,gpu_da
end
subroutine compute_J1_chol(nO,nV,t1,t2,v_ovvo,v_ovoo,v_vvoo,J1)
implicit none
integer, intent(in) :: nO,nV
double precision, intent(in) :: t1(nO, nV)
double precision, intent(in) :: t2(nO, nO, nV, nV)
double precision, intent(in) :: v_ovvo(nO,nV,nV,nO), v_ovoo(nO,nV,nO,nO)
double precision, intent(in) :: v_vvoo(nV,nV,nO,nO)
double precision, intent(inout) :: J1(nO, nV, nV, nO)
integer :: a,tmp_a,b,k,l,c,d,tmp_c,tmp_d,i,j,u,v, beta, gam
double precision, allocatable :: X_ovoo(:,:,:,:), Y_ovov(:,:,:,:)
allocate(X_ovoo(nO,nV,nO,nO),Y_ovov(nO,nV,nO,nV))
deallocate(X_ovoo)
double precision, allocatable :: tmp_cc(:,:,:), J1_tmp(:,:,:,:)
!TODO: I am here
!- cc_space_v_vvoo(a,b,i,j) * (0.5d0 * t2(u,j,b,beta) + t1(u,b) * t1(j,beta)) &
double precision, allocatable :: X_voov(:,:,:,:), Z_ovvo(:,:,:,:)
allocate(X_voov(nV,nO,nO,nV), Z_ovvo(nO,nV,nV,nO))
double precision, allocatable :: X_ovvo(:,:,:,:), Y_vovo(:,:,:,:)
allocate(X_ovvo(nO,nV,nV,nO))
allocate(Y_vovo(nV,nO,nV,nO))
!$omp parallel &
!$omp shared(nO,nV,J1,Z_ovvo,t2,Y_vovo,v_vvoo,X_ovvo) &
!$omp private(i,beta,a,u,j,b) &
!$omp default(none)
!+ 0.5d0 * (2d0 * cc_space_v_vvoo(a,b,i,j) - cc_space_v_vvoo(b,a,i,j)) * t2(u,j,beta,b)
do j = 1, nO
!$omp do
do b = 1, nV
do i = 1, nO
do a = 1, nV
Y_vovo(a,i,b,j) = 0.5d0 * (2d0 * v_vvoo(a,b,i,j) - v_vvoo(b,a,i,j))
enddo
enddo
enddo
!$omp end do nowait
enddo
do j = 1, nO
!$omp do
do b = 1, nV
do beta = 1, nV
do u = 1, nO
X_ovvo(u,beta,b,j) = t2(u,j,beta,b)
enddo
enddo
enddo
!$omp end do nowait
enddo
!$omp end parallel
call dgemm('N','T',nO*nV,nV*nO,nV*nO, &
1d0, X_ovvo, size(X_ovvo,1) * size(X_ovvo,2), &
Y_vovo, size(Y_vovo,1) * size(Y_vovo,2), &
0d0, Z_ovvo, size(Z_ovvo,1) * size(Z_ovvo,2))
!$omp parallel &
!$omp shared(nO,nV,J1,Z_ovvo) &
!$omp private(i,beta,a,u) &
!$omp default(none)
do i = 1, nO
!$omp do
do beta = 1, nV
do a = 1, nV
do u = 1, nO
J1(u,a,beta,i) = J1(u,a,beta,i) + Z_ovvo(u,beta,a,i)
enddo
enddo
enddo
!$omp end do nowait
enddo
!$omp end parallel
deallocate(X_ovvo,Z_ovvo)
end
! K1
subroutine compute_K1_chol(nO,nV,t1,t2,v_ovoo,v_vvoo,v_ovov,K1)

View File

@ -302,7 +302,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasDgemm(handle, CUBLAS_OP_N, CUBLAS_OP_N, m, n, k, &alpha, A, lda, B, ldb, &beta, C, ldc);
double* d_X_oovv;
lda = nO*nO;
cudaMalloc((void **)&d_X_oovv, nO*nO*nV*nV * sizeof(double));
alpha = 1.0;
@ -376,11 +375,9 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
#pragma omp section
{
double* d_X_vovv;
lda = nV*nO;
cudaMalloc((void **)&d_X_vovv, nV*nO*nV*BLOCK_SIZE * sizeof(double));
double* d_Y_oovv;
lda = nO*nO;
cudaMalloc((void **)&d_Y_oovv, nO*nO*nV*nV * sizeof(double));
for (int iblock=0 ; iblock<nV ; iblock += BLOCK_SIZE) {
@ -443,7 +440,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
#pragma omp section
{
double* d_tcc2;
lda = cholesky_mo_num;
cudaMalloc((void **)&d_tcc2, cholesky_mo_num*nV*nO * sizeof(double));
alpha = 1.0;
@ -455,7 +451,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasDgemm(handle, CUBLAS_OP_N, CUBLAS_OP_T, m, n, k, &alpha, A, lda, B, ldb, &beta, C, ldc);
double* d_tcc;
lda = cholesky_mo_num;
cudaMalloc((void **)&d_tcc, cholesky_mo_num*nO*nV * sizeof(double));
alpha = 1.0;
@ -467,7 +462,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasDgemm(handle, CUBLAS_OP_N, CUBLAS_OP_N, m, n, k, &alpha, A, lda, B, ldb, &beta, C, ldc);
double* d_X_ovvo;
lda = nO*nV;
cudaMalloc((void **)&d_X_ovvo, nO*nV*nV*nO * sizeof(double));
alpha = 1.0;
@ -511,7 +505,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
#pragma omp section
{
double* d_X_oovv;
lda = nO*nO;
cudaMalloc((void **)&d_X_oovv, nO*nO*nV*nV * sizeof(double));
alpha = 1.0;
@ -546,7 +539,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasSetStream(handle, NULL);
double* d_X_vovo;
lda = nV*nO;
cudaMalloc((void **)&d_X_vovo, nV*nO*nV*nO * sizeof(double));
alpha = 0.0;
@ -569,7 +561,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasSetStream(handle, NULL);
double* d_Y_oovo;
lda = nO*nO;
cudaMalloc((void **)&d_Y_oovo, nO*nO*nV*nO * sizeof(double));
alpha = 1.0;
@ -635,7 +626,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
{
double* d_J1;
lda = nO*nV;
cudaMalloc((void **)&d_J1, nO*nV*nV*nO * sizeof(double));
alpha = 1.0;
@ -647,7 +637,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
double* d_X_ovoo;
lda = nO*nV;
cudaMalloc((void **)&d_X_ovoo, nO*nV*nO*nO * sizeof(double));
alpha = 0.0;
beta = 1.0;
@ -670,7 +659,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
double* d_Y_ovov;
lda = nO*nV;
cudaMalloc((void **)&d_Y_ovov, nO*nV*nO*nV * sizeof(double));
alpha = 1.0;
@ -702,7 +690,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasSetStream(handle, NULL);
double* d_tmp_cc;
lda = cholesky_mo_num;
cudaMalloc((void **)&d_tmp_cc, cholesky_mo_num*nV*nO * sizeof(double));
alpha = 1.0;
@ -714,7 +701,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasDgemm(handle, CUBLAS_OP_N, CUBLAS_OP_T, m, n, k, &alpha, A, lda, B, ldb, &beta, C, ldc);
double* d_J1_tmp;
lda = nV*nO;
cudaMalloc((void **)&d_J1_tmp, nV*nO*nV*nO * sizeof(double));
alpha = 1.0;
@ -747,7 +733,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
//---
double* d_X_voov;
lda = nV*nO;
cudaMalloc((void **)&d_X_voov, nV*nO*nO*nV * sizeof(double));
for (int i=0 ; i<nV ; ++i) {
@ -781,7 +766,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasSetStream(handle, NULL);
double* d_Z_ovvo;
lda = nO*nV;
cudaMalloc((void **)&d_Z_ovvo, nO*nV*nV*nO * sizeof(double));
alpha = -1.0;
@ -812,6 +796,80 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cudaStreamDestroy(stream[i]);
}
cublasSetStream(handle, NULL);
double* d_Y_vovo;
cudaMalloc((void **)&d_Y_vovo, nV*nO*nV*nO * sizeof(double));
for (int i=0 ; i<nO ; ++i) {
cudaStreamCreate(&(stream[i]));
}
alpha = 1.0;
beta = -0.5;
for (int j=0 ; j<nO ; ++j) {
for (int i=0 ; i<nO ; ++i) {
cublasSetStream(handle, stream[i]);
A = &(d_cc_space_v_vvoo[nV*nV*(i+nO*j)]); lda = nV;
B = &(d_cc_space_v_vvoo[nV*nV*(i+nO*j)]); ldb = nV;
C = &(d_Y_vovo[nV*(i+nO*nV*j)]); ldc = nV*nO;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_T, nV, nV, &alpha, A, lda, &beta, B, ldb, C, ldc);
}
}
for (int i=0 ; i<nO ; ++i) {
cudaStreamDestroy(stream[i]);
}
cublasSetStream(handle, NULL);
double* d_X_ovvo;
cudaMalloc((void **)&d_X_ovvo, nO*nV*nV*nO * sizeof(double));
for (int i=0 ; i<nV ; ++i) {
cudaStreamCreate(&(stream[i]));
}
alpha = 1.0;
beta = 0.0;
for (int j=0 ; j<nO ; ++j) {
for (int b=0 ; b<nV ; ++b) {
cublasSetStream(handle, stream[b]);
A = &(d_t2[nO*(j+nO*nV*b)]); lda = nO*nO;
B = &(d_t2[nO*(j+nO*nV*b)]); ldb = nO*nO;
C = &(d_X_ovvo[nO*nV*(b+nV*j)]); ldc = nO;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_N, nO, nV, &alpha, A, lda, &beta, B, ldb, C, ldc);
}
}
for (int i=0 ; i<nV ; ++i) {
cudaStreamDestroy(stream[i]);
}
cublasSetStream(handle, NULL);
alpha = 1.0;
beta = 0.0;
m=nO*nV; n=nV*nO; k=nV*nO;
A=d_X_ovvo; lda=nO*nV;
B=d_Y_vovo; ldb=nV*nO;
C=d_Z_ovvo; ldc=nO*nV;
cublasDgemm(handle, CUBLAS_OP_N, CUBLAS_OP_T, m, n, k, &alpha, A, lda, B, ldb, &beta, C, ldc);
cudaFree(d_X_ovvo);
cudaFree(d_Y_vovo);
for (int i=0 ; i<nV ; ++i) {
cudaStreamCreate(&(stream[i]));
}
alpha = 1.0;
beta = 1.0;
for (int i=0 ; i<nO ; ++i) {
for (int b=0 ; b<nV ; ++b) {
cublasSetStream(handle, stream[b]);
A = &(d_J1[nO*nV*(b+nV*i)]); lda = nO;
B = &(d_Z_ovvo[nO*(b+nV*nV*i)]); ldb = nO*nV;
C = &(d_J1[nO*nV*(b+nV*i)]); ldc = nO;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_N, nO, nV, &alpha, A, lda, &beta, B, ldb, C, ldc);
}
}
for (int i=0 ; i<nV ; ++i) {
cudaStreamDestroy(stream[i]);
}
cublasSetStream(handle, NULL);
cudaFree(d_Z_ovvo);
lda = nO*nV;
@ -820,9 +878,9 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
}
#pragma omp section
{
}
// #pragma omp section
// {
// }
} // end sections