9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-12-06 11:43:30 +01:00

Removed omp_set_nested

This commit is contained in:
Anthony Scemama 2021-05-24 21:55:14 +02:00
parent 1685e41258
commit e7dd374ab9
5 changed files with 80 additions and 58 deletions

View File

@ -6,6 +6,23 @@ BEGIN_PROVIDER [ integer, ao_prim_num_max ]
ao_prim_num_max = maxval(ao_prim_num)
END_PROVIDER
BEGIN_PROVIDER [ integer, ao_shell, (ao_num) ]
implicit none
BEGIN_DOC
! Index of the shell to which the AO corresponds
END_DOC
integer :: i, j, k, n
k=0
do i=1,shell_num
n = shell_ang_mom(i)+1
do j=1,(n*(n+1))/2
k = k+1
ao_shell(k) = i
enddo
enddo
END_PROVIDER
BEGIN_PROVIDER [ double precision, ao_coef_normalized, (ao_num,ao_prim_num_max) ]
&BEGIN_PROVIDER [ double precision, ao_coef_normalization_factor, (ao_num) ]
implicit none
@ -30,7 +47,8 @@ END_PROVIDER
! Normalization of the primitives
if (primitives_normalized) then
do j=1,ao_prim_num(i)
call overlap_gaussian_xyz(C_A,C_A,ao_expo(i,j),ao_expo(i,j),powA,powA,overlap_x,overlap_y,overlap_z,norm,nz)
call overlap_gaussian_xyz(C_A,C_A,ao_expo(i,j),ao_expo(i,j), &
powA,powA,overlap_x,overlap_y,overlap_z,norm,nz)
ao_coef_normalized(i,j) = ao_coef(i,j)/sqrt(norm)
enddo
else
@ -198,38 +216,11 @@ END_PROVIDER
END_DOC
do i = 1, nucl_num
Nucl_num_shell_Aos(i) = 0
do j = 1, Nucl_N_Aos(i)
if(ao_l(Nucl_Aos(i,j))==0)then
! S type function
Nucl_num_shell_Aos(i)+=1
Nucl_list_shell_Aos(i,Nucl_num_shell_Aos(i))=Nucl_Aos(i,j)
elseif(ao_l(Nucl_Aos(i,j))==1)then
! P type function
if(ao_power(Nucl_Aos(i,j),1)==1)then
if (ao_power(Nucl_Aos(i,j),1) == ao_l(Nucl_Aos(i,j))) then
Nucl_num_shell_Aos(i)+=1
Nucl_list_shell_Aos(i,Nucl_num_shell_Aos(i))=Nucl_Aos(i,j)
endif
elseif(ao_l(Nucl_Aos(i,j))==2)then
! D type function
if(ao_power(Nucl_Aos(i,j),1)==2)then
Nucl_num_shell_Aos(i)+=1
Nucl_list_shell_Aos(i,Nucl_num_shell_Aos(i))=Nucl_Aos(i,j)
endif
elseif(ao_l(Nucl_Aos(i,j))==3)then
! F type function
if(ao_power(Nucl_Aos(i,j),1)==3)then
Nucl_num_shell_Aos(i)+=1
Nucl_list_shell_Aos(i,Nucl_num_shell_Aos(i))=Nucl_Aos(i,j)
endif
elseif(ao_l(Nucl_Aos(i,j))==4)then
! G type function
if(ao_power(Nucl_Aos(i,j),1)==4)then
Nucl_num_shell_Aos(i)+=1
Nucl_list_shell_Aos(i,Nucl_num_shell_Aos(i))=Nucl_Aos(i,j)
endif
endif
enddo
enddo

View File

@ -17,7 +17,7 @@ interface: ezfio, provider
type: double precision
doc: Number of primitives per |AO|
size: (basis.shell_num)
interface: ezfio, provider
interface: ezfio
[prim_num]
type: integer

View File

@ -1,32 +1,62 @@
BEGIN_PROVIDER [ double precision, basis_normalization_factor, (shell_num) ]
BEGIN_PROVIDER [ double precision, shell_normalization_factor , (shell_num) ]
implicit none
BEGIN_DOC
! Normalization factors of the shells
! Number of primitives per |AO|
END_DOC
double precision :: norm,overlap_x,overlap_y,overlap_z,C_A(3), c
integer :: l, powA(3), nz
integer :: i,j,k
nz=100
C_A(1) = 0.d0
C_A(2) = 0.d0
C_A(3) = 0.d0
do i=1,shell_num
powA(1) = shell_ang_mom(i)
powA(2) = 0
powA(3) = 0
logical :: has
PROVIDE ezfio_filename
if (mpi_master) then
if (size(shell_normalization_factor) == 0) return
! Normalization of the contracted basis functions
norm = 0.d0
do j=shell_prim_index(i), shell_prim_index(i)+shell_prim_num(i)-1
do k=shell_prim_index(i), shell_prim_index(i)+shell_prim_num(i)-1
call overlap_gaussian_xyz(C_A,C_A,shell_prim_expo(j),shell_prim_expo(k), &
powA,powA,overlap_x,overlap_y,overlap_z,c,nz)
norm = norm+c*shell_prim_coef(j)*shell_prim_coef(k)
enddo
enddo
basis_normalization_factor(i) = basis_normalization_factor(i) * sqrt(norm)
call ezfio_has_basis_shell_normalization_factor(has)
if (has) then
write(6,'(A)') '.. >>>>> [ IO READ: shell_normalization_factor ] <<<<< ..'
call ezfio_get_basis_shell_normalization_factor(shell_normalization_factor)
else
enddo
double precision :: norm,overlap_x,overlap_y,overlap_z,C_A(3), c
integer :: l, powA(3), nz
integer :: i,j,k
nz=100
C_A(1) = 0.d0
C_A(2) = 0.d0
C_A(3) = 0.d0
do i=1,shell_num
powA(1) = shell_ang_mom(i)
powA(2) = 0
powA(3) = 0
! Normalization of the contracted basis functions
norm = 0.d0
do j=shell_prim_index(i), shell_prim_index(i)+shell_prim_num(i)-1
do k=shell_prim_index(i), shell_prim_index(i)+shell_prim_num(i)-1
call overlap_gaussian_xyz(C_A,C_A,shell_prim_expo(j),shell_prim_expo(k),&
powA,powA,overlap_x,overlap_y,overlap_z,c,nz)
norm = norm+c*shell_prim_coef(j)*shell_prim_coef(k)
enddo
enddo
shell_normalization_factor(i) = dsqrt(norm)
enddo
endif
endif
IRP_IF MPI_DEBUG
print *, irp_here, mpi_rank
call MPI_BARRIER(MPI_COMM_WORLD, ierr)
IRP_ENDIF
IRP_IF MPI
include 'mpif.h'
integer :: ierr
call MPI_BCAST( shell_normalization_factor, (shell_num), MPI_DOUBLE_PRECISION, 0, MPI_COMM_WORLD, ierr)
if (ierr /= MPI_SUCCESS) then
stop 'Unable to read shell_normalization_factor with MPI'
endif
IRP_ENDIF
call write_time(6)
END_PROVIDER

View File

@ -286,7 +286,7 @@ subroutine ZMQ_pt2(E, pt2_data, pt2_data_err, relative_error, N_in)
call write_int(6,nproc_target,'Number of threads for PT2')
call write_double(6,mem,'Memory (Gb)')
call omp_set_nested(.false.)
call omp_set_max_active_levels(1)
print '(A)', '========== ======================= ===================== ===================== ==========='
@ -313,6 +313,7 @@ subroutine ZMQ_pt2(E, pt2_data, pt2_data_err, relative_error, N_in)
endif
!$OMP END PARALLEL
call end_parallel_job(zmq_to_qp_run_socket, zmq_socket_pull, 'pt2')
call omp_set_max_active_levels(8)
print '(A)', '========== ======================= ===================== ===================== ==========='

View File

@ -72,7 +72,7 @@ subroutine run_dress_slave(thread,iproce,energy)
provide psi_energy
ending = dress_N_cp+1
ntask_tbd = 0
call omp_set_nested(.true.)
call omp_set_max_active_levels(8)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(interesting, breve_delta_m, task_id) &
@ -84,7 +84,7 @@ subroutine run_dress_slave(thread,iproce,energy)
zmq_socket_push = new_zmq_push_socket(thread)
integer, external :: connect_to_taskserver
!$OMP CRITICAL
call omp_set_nested(.false.)
call omp_set_max_active_levels(1)
if (connect_to_taskserver(zmq_to_qp_run_socket,worker_id,thread) == -1) then
print *, irp_here, ': Unable to connect to task server'
stop -1
@ -296,7 +296,7 @@ subroutine run_dress_slave(thread,iproce,energy)
!$OMP END CRITICAL
!$OMP END PARALLEL
call omp_set_nested(.false.)
call omp_set_max_active_levels(1)
! do i=0,dress_N_cp+1
! call omp_destroy_lock(lck_sto(i))
! end do