File: writevclijjk.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (76 lines) | stat: -rw-r--r-- 2,172 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76

! Copyright (C) 2007-2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: writevclijjk
! !INTERFACE:
subroutine writevclijjk
! !USES:
use modmain
use modmpi
use modomp
! !DESCRIPTION:
!   Generates Coulomb matrix elements of the type $(i-jj-k)$ and outputs them to
!   the file {\tt VCLIJJK.OUT}. Also writes the real diagonal of this matrix,
!   $(i-jj-i)$, to {\tt VCLIJJI.OUT}.
!
! !REVISION HISTORY:
!   Created 2008 (Sharma)
!EOP
!BOC
implicit none
! local variables
integer ik,ist,recl,nthd
! allocatable arrays
real(8), allocatable :: vclijji(:,:,:)
complex(8), allocatable :: vclijjk(:,:,:,:)
! determine record length for vclijji and open file
allocate(vclijji(nstsv,nstsv,nkpt))
inquire(iolength=recl) vclijji
deallocate(vclijji)
open(260,file='VCLIJJI.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
! determine record length for vclijjk and open file
allocate(vclijjk(nstsv,nstsv,nstsv,nkpt))
inquire(iolength=recl) vclijjk
deallocate(vclijjk)
open(262,file='VCLIJJK.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
call omp_hold(nkptnr/np_mpi,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(vclijji,vclijjk,ist) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ik=1,nkptnr
! distribute among MPI processes
  if (mod(ik-1,np_mpi).ne.lp_mpi) cycle
  allocate(vclijji(nstsv,nstsv,nkpt))
  allocate(vclijjk(nstsv,nstsv,nstsv,nkpt))
!$OMP CRITICAL(writevclijjk_)
  write(*,'("Info(writevclijjk): ",I6," of ",I6," k-points")') ik,nkptnr
!$OMP END CRITICAL(writevclijjk_)
! calculate Coulomb matrix elements of the type (i-jj-k)
  call genvclijjk(ik,vclijjk)
! make a copy of the diagonal elements (i-jj-i)
  do ist=1,nstsv
    vclijji(ist,:,:)=dble(vclijjk(ist,ist,:,:))
  end do
!$OMP CRITICAL(u260)
  write(260,rec=ik) vclijji
!$OMP END CRITICAL(u260)
!$OMP CRITICAL(u262)
  write(262,rec=ik) vclijjk
!$OMP END CRITICAL(u262)
  deallocate(vclijji,vclijjk)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
close(260)
close(262)
! synchronise MPI processes
call mpi_barrier(mpicom,ierror)
return
end subroutine
!EOC