The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine oepvcl(vclcv,vclvv)
7use modmain
8use modmpi
9use modomp
10implicit none
11! arguments
12complex(8), intent(out) :: vclcv(ncrmax,natmtot,nstsv,nkpt)
13complex(8), intent(out) :: vclvv(nstsv,nstsv,nkpt)
14! local variables
15integer ik,ncv,nvv
16integer lp,nthd
17call holdthd(nkpt/np_mpi,nthd)
18!$OMP PARALLEL DO DEFAULT(SHARED) &
19!$OMP SCHEDULE(DYNAMIC) &
20!$OMP NUM_THREADS(nthd)
21do 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))
28end do
29!$OMP END PARALLEL DO
30call freethd(nthd)
31! broadcast matrix elements to all other processes
32ncv=ncrmax*natmtot*nstsv
33nvv=nstsv*nstsv
34do 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)
38end do
39end subroutine
40
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine oepvcl(vclcv, vclvv)
Definition oepvcl.f90:7
subroutine oepvclk(ikp, vclcv, vclvv)
Definition oepvclk.f90:7