Removed J1 from args

This commit is contained in:
Anthony Scemama 2023-08-03 02:12:15 +02:00
parent eec4af61e2
commit 803a59b586
3 changed files with 111 additions and 136 deletions

View File

@ -490,14 +490,11 @@ subroutine compute_r2_space_chol(nO,nV,t1,t2,tau,H_oo,H_vv,H_vo,r2,max_r2,gpu_da
!$omp end do
!$omp end parallel
double precision, allocatable :: J1(:,:,:,:)
allocate(J1(nO,nV,nV,nO))
double precision, allocatable :: K1(:,:,:,:)
allocate(K1(nO,nV,nO,nV))
call compute_r2_space_chol_gpu(nO,nV,cholesky_mo_num,gpu_data,t1,t2,tau, &
H_vv, g_occ, J1, K1, r2)
H_vv, g_occ, K1, r2)
!---
double precision, allocatable :: X_oovv(:,:,:,:)
@ -507,73 +504,8 @@ 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(:,:,:,:)
allocate(X_ovvo(nO,nV,nV,nO))
!$omp parallel &
!$omp private(u,v,gam,beta,i,a) &
!$omp default(shared)
do i = 1, nO
!$omp do
do a = 1, nV
do beta = 1, nV
do u = 1, nO
X_ovvo(u,beta,a,i) = (J1(u,a,beta,i) - 0.5d0 * K1(u,a,i,beta))
enddo
enddo
enddo
!$omp end do nowait
enddo
!$omp end parallel
deallocate(J1)
double precision, allocatable :: Y_voov(:,:,:,:)
allocate(Y_voov(nV,nO,nO,nV))
!$omp parallel &
!$omp private(u,v,gam,beta,i,a) &
!$omp default(shared)
!$omp do
do gam = 1, nV
do v = 1, nO
do i = 1, nO
do a = 1, nV
Y_voov(a,i,v,gam) = 2d0 * t2(i,v,a,gam) - t2(i,v,gam,a)
enddo
enddo
enddo
enddo
!$omp end do
!$omp end parallel
double precision, allocatable :: Z_ovov(:,:,:,:)
allocate(Z_ovov(nO,nV,nO,nV))
call dgemm('N','N', nO*nV,nO*nV,nV*nO, &
1d0, X_ovvo, size(X_ovvo,1) * size(X_ovvo,2), &
Y_voov, size(Y_voov,1) * size(Y_voov,2), &
0d0, Z_ovov, size(Z_ovov,1) * size(Z_ovov,2))
deallocate(X_ovvo,Y_voov)
!$omp parallel &
!$omp shared(nO,nV,r2,Z_ovov) &
!$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) + Z_ovov(u,beta,v,gam) + Z_ovov(v,gam,u,beta)
enddo
enddo
enddo
enddo
!$omp end do
!$omp end parallel
deallocate(Z_ovov)
double precision, allocatable :: Y_ovov(:,:,:,:), X_ovov(:,:,:,:)
allocate(X_ovov(nO,nV,nO,nV))
allocate(Y_ovov(nO,nV,nO,nV))

View File

@ -16,7 +16,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
double* tau,
double* H_vv,
double* g_occ,
double* J1,
double* K1,
double* r2)
{
@ -24,6 +23,8 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
int ngpus = 1;
cudaGetDeviceCount(&ngpus);
double* J1 = malloc(nO*nV*nV*nO*sizeof(double));
#pragma omp parallel num_threads(ngpus)
{
int m,n,k, lda, ldb, ldc;
@ -233,36 +234,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
C = d_r2; ldc = nO*nO;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_N, nO*nO, nV*nV, &alpha, A, lda, &beta, B, ldb, C, ldc);
/*
double * Y_oovv = malloc(nO*nO*nV*nV*sizeof(double));
lda=nO*nO;
cublasGetMatrix(nO*nO, nV*nV, sizeof(double), d_Y_oovv, lda, Y_oovv, lda);
cudaFree(d_Y_oovv);
double * r2_tmp = malloc(nO*nO*nV*nV*sizeof(double));
lda=nO*nO;
cublasGetMatrix(nO*nO, nV*nV, sizeof(double), d_r2, lda, r2_tmp, lda);
for (int j=0 ; j<nV ; ++j) {
for (int i=0 ; i<nV ; ++i) {
double * xx = &(r2_tmp[nO*nO*(i+nV*j)]);
const double * yy = &(Y_oovv[nO*nO*(j+nV*i)]);
for (int k=0 ; k<nO ; ++k) {
for (int l=0 ; l<nO ; ++l) {
xx[l + k*nO] += yy[k + l*nO];
}
}
}
}
free(Y_oovv);
lda=nO*nO;
cublasSetMatrix(nO*nO, nV*nV, sizeof(double), r2_tmp, lda, d_r2, lda);
free(r2_tmp);
*/
//--
for (int i=0 ; i<nV ; ++i) {
cudaStreamCreate(&(stream[i]));
}
@ -283,7 +254,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
}
cublasSetStream(handle, NULL);
cudaFree(d_Y_oovv);
//--
}
// g_occ
@ -321,35 +291,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
C = d_r2; ldc = nO*nO;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_N, nO*nO, nV*nV, &alpha, A, lda, &beta, B, ldb, C, ldc);
/*
double * X_oovv = malloc(nO*nO*nV*nV*sizeof(double));
lda=nO*nO;
cublasGetMatrix(nO*nO, nV*nV, sizeof(double), d_X_oovv, lda, X_oovv, lda);
cudaFree(d_X_oovv);
double * r2_tmp = malloc(nO*nO*nV*nV*sizeof(double));
lda=nO*nO;
cublasGetMatrix(nO*nO, nV*nV, sizeof(double), d_r2, lda, r2_tmp, lda);
for (int j=0 ; j<nV ; ++j) {
for (int i=0 ; i<nV ; ++i) {
double * xx = &(r2_tmp[nO*nO*(i+nV*j)]);
const double * yy = &(X_oovv[nO*nO*(j+nV*i)]);
for (int k=0 ; k<nO ; ++k) {
for (int l=0 ; l<nO ; ++l) {
xx[l + k*nO] -= yy[k + l*nO];
}
}
}
}
free(X_oovv);
lda=nO*nO;
cublasSetMatrix(nO*nO, nV*nV, sizeof(double), r2_tmp, lda, d_r2, lda);
free(r2_tmp);
*/
//--
for (int i=0 ; i<nV ; ++i) {
cudaStreamCreate(&(stream[i]));
}
@ -369,8 +310,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
}
cublasSetStream(handle, NULL);
cudaFree(d_X_oovv);
//--
}
#pragma omp section
@ -724,7 +663,6 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
cublasSetStream(handle, NULL);
cudaFree(d_J1_tmp);
//---
double* d_X_voov;
cudaMalloc((void **)&d_X_voov, nV*nO*nO*nV * sizeof(double));
@ -983,7 +921,7 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
}
cublasSetStream(handle, NULL);
cudaFree(d_K1tmp);
// cudaFree(d_Z);
cudaFree(d_Z);
lda = nO*nV;
cublasGetMatrix(nO*nV, nO*nV, sizeof(double), d_K1, lda, K1, lda);
@ -993,6 +931,112 @@ void compute_r2_space_chol_gpu(const int nO, const int nV, const int cholesky_mo
} // end sections
#pragma omp barrier
#pragma omp sections
{
#pragma omp section
{
double* d_K1;
lda = nO*nV;
cudaMalloc((void **)&d_K1, nO*nV*nO*nV * sizeof(double));
cublasSetMatrix(lda, nO*nV, sizeof(double), K1, lda, d_K1, lda);
double* d_J1;
lda = nO*nV;
cudaMalloc((void **)&d_J1, nO*nV*nV*nO * sizeof(double));
cublasSetMatrix(lda, nV*nO, sizeof(double), J1, lda, d_J1, lda);
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.5;
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_K1[nO*nV*(i+nO*b)]); ldb = nO;
C = &(d_X_ovvo[nO*(b+nV*nV*i)]); ldc = nO*nV;
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_J1);
double* d_Y_voov;
cudaMalloc((void **)&d_Y_voov, nV*nO*nO*nV * sizeof(double));
for (int i=0 ; i<nV ; ++i) {
cudaStreamCreate(&(stream[i]));
}
alpha = 2.0;
beta = -1.0;
for (int v=0 ; v<nO ; ++v) {
for (int g=0 ; g<nV ; ++g) {
cublasSetStream(handle, stream[g]);
A = &(d_t2[nO*(v+nO*nV*g)]); lda = nO*nO;
B = &(d_t2[nO*(v+nO*g)]); ldb = nO*nO*nV;
C = &(d_Y_voov[nV*nO*(v+nO*g)]); ldc = nV;
cublasDgeam(handle, CUBLAS_OP_T, CUBLAS_OP_T, nV, nO, &alpha, A, lda, &beta, B, ldb, C, ldc);
}
}
for (int i=0 ; i<nV ; ++i) {
cudaStreamDestroy(stream[i]);
}
cublasSetStream(handle, NULL);
double* d_Z_ovov;
cudaMalloc((void **)&d_Z_ovov, nO*nV*nO*nV * sizeof(double));
alpha = 1.0;
beta = 0.0;
m=nO*nV; n=nO*nV; k=nV*nO;
A=d_X_ovvo; lda=nO*nV;
B=d_Y_voov; ldb=nV*nO;
C=d_Z_ovov; ldc=nO*nV;
cublasDgemm(handle, CUBLAS_OP_N, CUBLAS_OP_N, m, n, k, &alpha, A, lda, B, ldb, &beta, C, ldc);
cudaFree(d_X_ovvo);
cudaFree(d_Y_voov);
for (int i=0 ; i<nV ; ++i) {
cudaStreamCreate(&(stream[i]));
}
alpha = 1.0;
beta = 1.0;
for (int b=0 ; b<nV ; ++b) {
for (int g=0 ; g<nV ; ++g) {
cublasSetStream(handle, stream[g]);
A = &(d_r2[nO*nO*(b+nV*g)]); lda = nO;
B = &(d_Z_ovov[nO*(b+nV*nO*g)]); ldb = nO*nV;
C = &(d_r2[nO*nO*(b+nV*g)]); ldc = nO;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_N, nO, nO, &alpha, A, lda, &beta, B, ldb, C, ldc);
}
for (int g=0 ; g<nV ; ++g) {
cublasSetStream(handle, stream[g]);
A = &(d_r2[nO*nO*(b+nV*g)]); lda = nO;
B = &(d_Z_ovov[nO*(g+nV*nO*b)]); ldb = nO*nV;
C = &(d_r2[nO*nO*(b+nV*g)]); ldc = nO;
cublasDgeam(handle, CUBLAS_OP_N, CUBLAS_OP_T, nO, nO, &alpha, A, lda, &beta, B, ldb, C, ldc);
}
}
for (int i=0 ; i<nV ; ++i) {
cudaStreamDestroy(stream[i]);
}
cublasSetStream(handle, NULL);
cudaFree(d_Z_ovov);
}
} // end sections
lda = cholesky_mo_num * nV;

View File

@ -27,7 +27,7 @@ module gpu_module
end function
subroutine compute_r2_space_chol_gpu(nO,nV,cholesky_mo_num, gpu_data, t1, t2, tau,&
H_vv, g_occ, J1, K1, r2) bind(C)
H_vv, g_occ, K1, r2) bind(C)
import c_int, c_double, c_ptr
integer(c_int), intent(in), value :: nO, nV, cholesky_mo_num
type(c_ptr), value :: gpu_data
@ -36,7 +36,6 @@ module gpu_module
real(c_double), intent(in) :: tau(nO,nO,nV,nV)
real(c_double), intent(in) :: H_vv(nV,nV)
real(c_double), intent(in) :: g_occ(nO,nO)
real(c_double) :: J1(nO,nV,nV,nO)
real(c_double) :: K1(nO,nV,nO,nV)
real(c_double), intent(out) :: r2(nO,nO,nV,nV)
end subroutine