The Elk Code
oepvcl.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2006 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine oepvcl(vclcv,vclvv)
7 use modmain
8 use modmpi
9 use modomp
10 implicit none
11 ! arguments
12 complex(8), intent(out) :: vclcv(ncrmax,natmtot,nstsv,nkpt)
13 complex(8), intent(out) :: vclvv(nstsv,nstsv,nkpt)
14 ! local variables
15 integer ik,ncv,nvv
16 integer lp,nthd
17 call holdthd(nkpt/np_mpi,nthd)
18 !$OMP PARALLEL DO DEFAULT(SHARED) &
19 !$OMP SCHEDULE(DYNAMIC) &
20 !$OMP NUM_THREADS(nthd)
21 do ik=1,nkpt
22 ! distribute among MPI processes
23  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
24 !$OMP CRITICAL(oepvcl_)
25  write(*,'("Info(oepvcl): ",I6," of ",I6," k-points")') ik,nkpt
26 !$OMP END CRITICAL(oepvcl_)
27  call oepvclk(ik,vclcv(:,:,:,ik),vclvv(:,:,ik))
28 end do
29 !$OMP END PARALLEL DO
30 call freethd(nthd)
31 ! broadcast matrix elements to all other processes
32 ncv=ncrmax*natmtot*nstsv
33 nvv=nstsv*nstsv
34 do ik=1,nkpt
35  lp=mod(ik-1,np_mpi)
36  call mpi_bcast(vclcv(:,:,:,ik),ncv,mpi_double_complex,lp,mpicom,ierror)
37  call mpi_bcast(vclvv(:,:,ik),nvv,mpi_double_complex,lp,mpicom,ierror)
38 end do
39 end subroutine
40 
Definition: modomp.f90:6
subroutine oepvcl(vclcv, vclvv)
Definition: oepvcl.f90:7
integer np_mpi
Definition: modmpi.f90:13
subroutine oepvclk(ikp, vclcv, vclvv)
Definition: oepvclk.f90:7
Definition: modmpi.f90:6
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19