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

More on GPU

This commit is contained in:
Anthony Scemama 2023-07-17 14:51:50 +02:00
parent 79060ff6e1
commit 493f2cf7d8
3 changed files with 74 additions and 96 deletions

View File

@ -490,7 +490,7 @@ subroutine compute_r2_space_chol(nO,nV,t1,t2,tau,H_oo,H_vv,H_vo,r2,max_r2)
!$omp end parallel
call compute_r2_space_chol_gpu(nO,nV,cholesky_mo_num,t1,t2,tau, &
cc_space_v_ov_chol, cc_space_v_vo_chol, cc_space_v_vv_chol, &
cc_space_v_oo_chol, cc_space_v_ov_chol, cc_space_v_vo_chol, cc_space_v_vv_chol, &
cc_space_v_oooo, cc_space_v_vooo, cc_space_v_oovv, cc_space_v_vvoo, &
cc_space_f_vo, H_vv, g_occ, r2)
@ -498,92 +498,8 @@ subroutine compute_r2_space_chol(nO,nV,t1,t2,tau,H_oo,H_vv,H_vo,r2,max_r2)
!---
double precision, allocatable :: X_oovv(:,:,:,:)
double precision, allocatable :: X_vovv(:,:,:,:)
! allocate(X_vovv(nV,nO,nV,block_size))
! allocate(Y_oovv(nO,nO,nV,nV))
!
! do iblock = 1, nV, block_size
! do gam = iblock, min(nV, iblock+block_size-1)
! call dgemm('T','N',nV, nO*nV, cholesky_mo_num, 1.d0, &
! cc_space_v_vv_chol(1,1,gam), cholesky_mo_num, cc_space_v_ov_chol, &
! cholesky_mo_num, 0.d0, X_vovv(1,1,1,gam-iblock+1), nV)
!
! enddo
! call dgemm('N','N',nO,nO*nV*min(block_size, nV-iblock+1),nV, &
! 1d0, t1 , size(t1,1), &
! X_vovv, size(X_vovv,1), &
! 0d0, Y_oovv(1,1,1,iblock), size(Y_oovv,1))
!
! enddo
! deallocate(X_vovv)
!
! !$omp parallel &
! !$omp shared(nO,nV,r2,Y_oovv) &
! !$omp private(u,v,gam,beta) &
! !$omp default(none)
! !$omp do
! do gam = 1, nV
! do beta = 1, nV
! do v = 1, nO
! do u = 1, nO
! r2(u,v,beta,gam) = r2(u,v,beta,gam) + Y_oovv(v,u,beta,gam) + Y_oovv(u,v,gam,beta)
! enddo
! enddo
! enddo
! enddo
! !$omp end do
! !$omp end parallel
! deallocate(Y_oovv)
double precision, allocatable :: X_ovvo(:,:,:,:)
double precision, allocatable :: tcc(:,:,:), tcc2(:,:,:)
allocate(tcc2(cholesky_mo_num,nV,nO), X_ovvo(nO,nV,nV,nO))
allocate(tcc(cholesky_mo_num,nO,nV))
call dgemm('N','T', cholesky_mo_num*nV, nO, nV, 1.d0, &
cc_space_v_vv_chol, cholesky_mo_num*nV, t1, nO, &
0.d0, tcc2, cholesky_mo_num*nV)
call dgemm('N','N', cholesky_mo_num*nO, nV, nO, 1.d0, &
cc_space_v_oo_chol, cholesky_mo_num*nO, t1, nO, &
0.d0, tcc, cholesky_mo_num*nO)
call dgemm('T','N', nO*nV, nV*nO, cholesky_mo_num, 1.d0, &
tcc, cholesky_mo_num, tcc2, cholesky_mo_num, 0.d0, &
X_ovvo, nO*nV)
deallocate(tcc, tcc2)
!$omp parallel &
!$omp shared(nO,nV,r2,X_ovvo) &
!$omp private(u,v,gam,beta) &
!$omp default(none)
!$omp do
do gam = 1, nV
do beta = 1, nV
do v = 1, nO
do u = 1, nO
r2(u,v,beta,gam) = r2(u,v,beta,gam) - X_ovvo(u,beta,gam,v)
enddo
enddo
enddo
enddo
!$omp end do
!$omp do
do beta = 1, nV
do gam = 1, nV
do v = 1, nO
do u = 1, nO
r2(v,u,gam,beta) = r2(v,u,gam,beta) - X_ovvo(u,beta,gam,v)
enddo
enddo
enddo
enddo
!$omp end do
!$omp end parallel
deallocate(X_ovvo)
!-----
allocate(X_oovv(nO,nO,nV,nV))

View File

@ -63,6 +63,7 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
double* t1,
double* t2,
double* tau,
double* cc_space_v_oo_chol,
double* cc_space_v_ov_chol,
double* cc_space_v_vo_chol,
double* cc_space_v_vv_chol,
@ -96,9 +97,10 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
double* d_tau;
double* d_r2;
double* d_cc_space_v_vv_chol;
double* d_cc_space_v_vo_chol;
double* d_cc_space_v_oo_chol;
double* d_cc_space_v_ov_chol;
double* d_cc_space_v_vo_chol;
double* d_cc_space_v_vv_chol;
double* d_t1;
double* d_t2;
double* d_tmp_cc;
@ -116,6 +118,10 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cudaMalloc((void **)&d_cc_space_v_vv_chol, lda * nV * sizeof(double));
cublasSetMatrix(cholesky_mo_num*nV, nV, sizeof(double), cc_space_v_vv_chol, lda, d_cc_space_v_vv_chol, lda);
lda = cholesky_mo_num * nO;
cudaMalloc((void **)&d_cc_space_v_oo_chol, lda * nO * sizeof(double));
cublasSetMatrix(cholesky_mo_num*nO, nO, sizeof(double), cc_space_v_oo_chol, lda, d_cc_space_v_oo_chol, lda);
lda = cholesky_mo_num * nO;
cudaMalloc((void **)&d_cc_space_v_ov_chol, lda * nV * sizeof(double));
cublasSetMatrix(cholesky_mo_num*nO, nV, sizeof(double), cc_space_v_ov_chol, lda, d_cc_space_v_ov_chol, lda);
@ -360,9 +366,9 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
lda = nO*nO;
cudaMalloc((void **)&d_Y_oovv, nO*nO*nV*nV * sizeof(double));
for (size_t iblock=0 ; iblock<nV ; iblock += BLOCK_SIZE) {
size_t mbs = nV < iblock+BLOCK_SIZE ? nV : iblock+BLOCK_SIZE;
for (size_t gam=iblock ; gam<mbs ; ++gam) {
for (int iblock=0 ; iblock<nV ; iblock += BLOCK_SIZE) {
int mbs = nV < iblock+BLOCK_SIZE ? nV : iblock+BLOCK_SIZE;
for (int gam=iblock ; gam<mbs ; ++gam) {
alpha = 1.0;
beta = 0.0;
m=nV; n=nO*nV; k=cholesky_mo_num;
@ -396,7 +402,62 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
}
}
cudaFree(d_Y_oovv);
}
#pragma omp section
{
double* d_tcc2;
lda = cholesky_mo_num;
cudaMalloc((void **)&d_tcc2, cholesky_mo_num*nV*nO * sizeof(double));
alpha = 1.0;
beta = 0.0;
m=cholesky_mo_num*nV; n=nO; k=nV;
A=d_cc_space_v_vv_chol; lda=cholesky_mo_num*nV;
B=d_t1; ldb=nO;
C=d_tcc2; ldc=cholesky_mo_num*nV;
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;
beta = 0.0;
m=cholesky_mo_num*nO; n=nV; k=nO;
A=d_cc_space_v_oo_chol; lda=cholesky_mo_num*nO;
B=d_t1; ldb=nO;
C=d_tcc; ldc=cholesky_mo_num*nO;
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;
beta = 0.0;
m=nO*nV; n=nV*nO; k=cholesky_mo_num;
A=d_tcc; lda=cholesky_mo_num;
B=d_tcc2; ldb=cholesky_mo_num;
C=d_X_ovvo; ldc=nO*nV;
cublasDgemm(handle, CUBLAS_OP_T, CUBLAS_OP_N, m, n, k, &alpha, A, lda, B, ldb, &beta, C, ldc);
cudaFree(d_tcc);
cudaFree(d_tcc2);
for(int gam = 0; gam < nV; gam++){
for(int bet = 0; bet < nV; bet++){
alpha = 1.0;
beta = -1.0;
A = &(d_r2[nO*nO*(bet+nV*gam)]); lda = nO;
B = &(d_X_ovvo[nO*(bet+nV*gam)]); ldb = nO*nV*nV;
C = &(d_r2[nO*nO*(bet+nV*gam)]); ldc = nO;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_N, nO, nO, &alpha, A, lda, &beta, B, ldb, C, ldc);
B = &(d_X_ovvo[nO*(gam+nV*bet)]); ldb = nO*nV*nV;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_T, nO, nO, &alpha, A, lda, &beta, B, ldb, C, ldc);
}
}
cudaFree(d_X_ovvo);
}
}
@ -420,7 +481,7 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cudaMalloc((void**)&d_tmpB1, nV*BLOCK_SIZE*nV*sizeof(double));
#pragma omp for
for (size_t gam=0 ; gam<nV ; ++gam)
for (int gam=0 ; gam<nV ; ++gam)
{
double* d_tmp_cc_ = &(d_tmp_cc[gam*nV*cholesky_mo_num]);
double* d_cc_space_v_vv_chol_ = &(d_cc_space_v_vv_chol[gam*nV*cholesky_mo_num]);
@ -432,9 +493,9 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
C = d_tmp_cc2 ; ldc = cholesky_mo_num;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_N, cholesky_mo_num, nV, &alpha, A, lda, &beta, B, ldb, C, ldc);
for (size_t iblock=0 ; iblock<nV ; iblock += BLOCK_SIZE)
for (int iblock=0 ; iblock<nV ; iblock += BLOCK_SIZE)
{
const size_t mbs = BLOCK_SIZE < nV-iblock ? BLOCK_SIZE : nV-iblock;
const int mbs = BLOCK_SIZE < nV-iblock ? BLOCK_SIZE : nV-iblock;
alpha=-1.0; beta=0.0;
m=nV*mbs; n=nV; k=cholesky_mo_num;
@ -486,7 +547,7 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasGetMatrix(nO*nO, nV*nV, sizeof(double), d_r2, lda, r2_tmp, lda);
#pragma omp critical
{
for (size_t i=0 ; i<nO*nO*nV*nV ; ++i) {
for (size_t i=0 ; i<(size_t) nO*nO*nV*nV ; ++i) {
r2[i] += r2_tmp[i];
}
}
@ -496,7 +557,7 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasDestroy(handle);
}
for (size_t i=0 ; i<nO*nO*nV*nV ; ++i)
for (size_t i=0 ; i<(size_t) nO*nO*nV*nV ; ++i)
{
r2[i] += cc_space_v_oovv[i];
}

View File

@ -5,7 +5,7 @@ module gpu_module
interface
subroutine compute_r2_space_chol_gpu(nO,nV,cholesky_mo_num, t1, t2, tau,&
cc_space_v_ov_chol, cc_space_v_vo_chol, cc_space_v_vv_chol, &
cc_space_v_oo_chol, cc_space_v_ov_chol, cc_space_v_vo_chol, cc_space_v_vv_chol, &
cc_space_v_oooo, cc_space_v_vooo, cc_space_v_oovv, cc_space_v_vvoo, &
cc_space_f_vo, H_vv, g_occ, r2) bind(C)
import c_int, c_double
@ -13,6 +13,7 @@ module gpu_module
real(c_double), intent(in) :: t1(nO,nV)
real(c_double), intent(in) :: t2(nO,nO,nV,nV)
real(c_double), intent(in) :: tau(nO,nO,nV,nV)
real(c_double), intent(in) :: cc_space_v_oo_chol(cholesky_mo_num,nO,nO)
real(c_double), intent(in) :: cc_space_v_ov_chol(cholesky_mo_num,nO,nV)
real(c_double), intent(in) :: cc_space_v_vo_chol(cholesky_mo_num,nV,nO)
real(c_double), intent(in) :: cc_space_v_vv_chol(cholesky_mo_num,nV,nV)