The Elk Code
 
Loading...
Searching...
No Matches
writevcl1223.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: writevcl1223
8! !INTERFACE:
9subroutine writevcl1223
10! !USES:
11use modmain
12use modmpi
13use modomp
14! !DESCRIPTION:
15! Generates Coulomb matrix elements of the type $V(1,2,2,3)$ and outputs them
16! to the file {\tt VCL1223.OUT}. Also writes the real diagonal of this matrix,
17! $V(1,2,2,1)$, to {\tt VCL1221.OUT}.
18!
19! !REVISION HISTORY:
20! Created 2008 (Sharma)
21!EOP
22!BOC
23implicit none
24! local variables
25integer ik,ist,recl,nthd
26! allocatable arrays
27real(8), allocatable :: vcl1221(:,:,:)
28complex(8), allocatable :: vcl1223(:,:,:,:)
29! determine record length for vcl1221 and open file
30allocate(vcl1221(nstsv,nstsv,nkpt))
31inquire(iolength=recl) vcl1221
32deallocate(vcl1221)
33open(260,file='VCL1221.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
34! determine record length for vcl1223 and open file
35allocate(vcl1223(nstsv,nstsv,nstsv,nkpt))
36inquire(iolength=recl) vcl1223
37deallocate(vcl1223)
38open(262,file='VCL1223.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
39call holdthd(nkptnr/np_mpi,nthd)
40!$OMP PARALLEL DEFAULT(SHARED) &
41!$OMP PRIVATE(vcl1221,vcl1223,ist) &
42!$OMP NUM_THREADS(nthd)
43allocate(vcl1221(nstsv,nstsv,nkpt),vcl1223(nstsv,nstsv,nstsv,nkpt))
44!$OMP DO SCHEDULE(DYNAMIC)
45do ik=1,nkptnr
46! distribute among MPI processes
47 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
48!$OMP CRITICAL(writevcl1223_)
49 write(*,'("Info(writevcl1223): ",I6," of ",I6," k-points")') ik,nkptnr
50!$OMP END CRITICAL(writevcl1223_)
51! calculate Coulomb matrix elements of the type V(1,2,2,3)
52 call genvcl1223(ik,vcl1223)
53! make a copy of the diagonal elements V(1,2,2,1)
54 do ist=1,nstsv
55 vcl1221(ist,:,:)=dble(vcl1223(ist,ist,:,:))
56 end do
57!$OMP CRITICAL(u260)
58 write(260,rec=ik) vcl1221
59!$OMP END CRITICAL(u260)
60!$OMP CRITICAL(u262)
61 write(262,rec=ik) vcl1223
62!$OMP END CRITICAL(u262)
63end do
64!$OMP END DO
65deallocate(vcl1221,vcl1223)
66!$OMP END PARALLEL
67call freethd(nthd)
68close(260)
69close(262)
70! synchronise MPI processes
71call mpi_barrier(mpicom,ierror)
72end subroutine
73!EOC
74
subroutine genvcl1223(ikp, vcl1223)
integer nkptnr
Definition modmain.f90:463
integer nkpt
Definition modmain.f90:461
integer nstsv
Definition modmain.f90:886
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 writevcl1223