The Elk Code
writevcl1221.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007-2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: writevcl1221
8 ! !INTERFACE:
9 subroutine writevcl1221
10 ! !USES:
11 use modmain
12 use modmpi
13 use modomp
14 ! !DESCRIPTION:
15 ! Generates Coulomb matrix elements of the type $V(1,2,2,1)$ and outputs them
16 ! to the file {\tt VCL1221.OUT}.
17 !
18 ! !REVISION HISTORY:
19 ! Created 2008 (Sharma)
20 !EOP
21 !BOC
22 implicit none
23 ! allocatable arrays
24 real(8), allocatable :: vcl1221(:,:,:)
25 integer recl,ik,nthd
26 ! determine record length for vcl1221 and open file
27 allocate(vcl1221(nstsv,nstsv,nkpt))
28 inquire(iolength=recl) vcl1221
29 deallocate(vcl1221)
30 open(260,file='VCL1221.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
31 call holdthd(nkptnr/np_mpi,nthd)
32 !$OMP PARALLEL DEFAULT(SHARED) &
33 !$OMP PRIVATE(vcl1221) &
34 !$OMP NUM_THREADS(nthd)
35 allocate(vcl1221(nstsv,nstsv,nkpt))
36 !$OMP DO SCHEDULE(DYNAMIC)
37 do ik=1,nkptnr
38 ! distribute among MPI processes
39  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
40 !$OMP CRITICAL(writevcl1221_)
41  write(*,'("Info(writevcl1221): ",I6," of ",I6," k-points")') ik,nkptnr
42 !$OMP END CRITICAL(writevcl1221_)
43 ! calculate Coulomb matrix elements of the type V(1,2,2,1)
44  call genvcl1221(ik,vcl1221)
45 !$OMP CRITICAL(u260)
46  write(260,rec=ik) vcl1221
47 !$OMP END CRITICAL(u260)
48 end do
49 !$OMP END DO
50 deallocate(vcl1221)
51 !$OMP END PARALLEL
52 call freethd(nthd)
53 close(260)
54 ! synchronise MPI processes
55 call mpi_barrier(mpicom,ierror)
56 end subroutine
57 !EOC
58 
integer nkpt
Definition: modmain.f90:461
Definition: modomp.f90:6
integer nkptnr
Definition: modmain.f90:463
integer np_mpi
Definition: modmpi.f90:13
integer nstsv
Definition: modmain.f90:889
subroutine writevcl1221
subroutine genvcl1221(ikp, vcl1221)
Definition: genvcl1221.f90:10
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