The Elk Code
getcfgq.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2011 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 subroutine getcfgq(fname,vpl,ng,m,cf)
7 use modmain
8 use modramdisk
9 implicit none
10 ! arguments
11 character(*), intent(in) :: fname
12 real(8), intent(in) :: vpl(3)
13 integer, intent(in) :: ng,m
14 complex(8), intent(out) :: cf(ng,ng,m)
15 ! local variables
16 logical tgs
17 integer isym,iq,lspl,ilspl
18 integer igq,jgq,igp,jgp,i,j
19 integer recl,ng_,m_
20 real(8) vql_(3),si(3,3)
21 real(8) vgql(3),v(3),t1
22 complex(8) z1
23 ! automatic arrays
24 logical done(ng)
25 integer map(ng)
26 ! allocatable arrays
27 real(8), allocatable :: vgpl(:,:)
28 complex(8), allocatable :: cf_(:,:,:),x(:)
29 ! find the equivalent reduced q-point and symmetry which rotates vql to vpl
30 call findqpt(vpl,isym,iq)
31 !$OMP CRITICAL(u245)
32 ! read from RAM disk if required
33 if (ramdisk) then
34  call getrd(trim(fname),iq,tgs,v1=vql_,n1=ng_,n2=m_,nzv=ng*ng*m,zva=cf)
35  if (tgs) goto 10
36 end if
37 ! find the record length
38 inquire(iolength=recl) vql_,ng_,m_,cf
39 open(245,file=trim(fname),form='UNFORMATTED',access='DIRECT',recl=recl)
40 read(245,rec=iq) vql_,ng_,m_,cf
41 close(245)
42 10 continue
43 !$OMP END CRITICAL(u245)
44 t1=abs(vql(1,iq)-vql_(1))+abs(vql(2,iq)-vql_(2))+abs(vql(3,iq)-vql_(3))
45 if (t1 > epslat) then
46  write(*,*)
47  write(*,'("Error(getcfgq): differing vectors for q-point ",I8)') iq
48  write(*,'(" current : ",3G18.10)') vql(:,iq)
49  write(*,'(" file : ",3G18.10)') vql_
50  write(*,'(" in file ",A)') trim(fname)
51  write(*,*)
52  stop
53 end if
54 if (ng /= ng_) then
55  write(*,*)
56  write(*,'("Error(getcfgq): differing ng for q-point ",I8)') iq
57  write(*,'(" current : ",I8)') ng
58  write(*,'(" file : ",I8)') ng_
59  write(*,'(" in file ",A)') trim(fname)
60  write(*,*)
61  stop
62 end if
63 if (m /= m_) then
64  write(*,*)
65  write(*,'("Error(getcfgq): differing m for q-point ",I8)') iq
66  write(*,'(" current : ",I8)') m
67  write(*,'(" file : ",I8)') m_
68  write(*,'(" in file ",A)') trim(fname)
69  write(*,*)
70  stop
71 end if
72 ! if p = q then return
73 t1=abs(vpl(1)-vql(1,iq))+abs(vpl(2)-vql(2,iq))+abs(vpl(3)-vql(3,iq))
74 if (t1 < epslat) return
75 ! allocate local arrays
76 allocate(vgpl(3,ng),cf_(ng,ng,m),x(ng))
77 ! perform translation operation and store in temporary array
78 if (tv0symc(isym)) then
79 ! translation vector is zero
80  cf_(1:ng,1:ng,1:m)=cf(1:ng,1:ng,1:m)
81 else
82 ! non-zero translation vector gives a phase factor
83  v(1:3)=vtcsymc(1:3,isym)
84  do igq=1,ng
85  t1=vgc(1,igq)*v(1)+vgc(2,igq)*v(2)+vgc(3,igq)*v(3)
86  x(igq)=cmplx(cos(t1),-sin(t1),8)
87  end do
88  do jgq=1,ng
89  z1=conjg(x(jgq))
90  do igq=1,ng
91  cf_(igq,jgq,1:m)=z1*x(igq)*cf(igq,jgq,1:m)
92  end do
93  end do
94 end if
95 ! index to spatial rotation in lattice point group
96 lspl=lsplsymc(isym)
97 ! the inverse of the spatial symmetry
98 ilspl=isymlat(lspl)
99 si(1:3,1:3)=dble(symlat(1:3,1:3,ilspl))
100 ! find the map from {G+q} to {G+p}
101 map(1:ng)=0
102 do igp=1,ng
103  vgpl(1:3,igp)=dble(ivg(1:3,igp))+vpl(1:3)
104 end do
105 done(1:ng)=.false.
106 i=1
107 do igq=1,ng
108  vgql(1:3)=dble(ivg(1:3,igq))+vql(1:3,iq)
109  call r3mtv(si,vgql,v)
110  j=0
111  do igp=i,ng
112  if (done(igp)) cycle
113  if (j == 0) j=igp
114  t1=abs(v(1)-vgpl(1,igp))+abs(v(2)-vgpl(2,igp))+abs(v(3)-vgpl(3,igp))
115  if (t1 < epslat) then
116  map(igp)=igq
117  done(igp)=.true.
118  if (igp == j) j=j+1
119  exit
120  end if
121  end do
122  if (j > 0) i=j
123 end do
124 ! rotate correlation function (passive transformation)
125 do jgp=1,ng
126  jgq=map(jgp)
127  do igp=1,ng
128  igq=map(igp)
129  if ((igq == 0).or.(jgq == 0)) then
130  cf(igp,jgp,1:m)=0.d0
131  else
132  cf(igp,jgp,1:m)=cf_(igq,jgq,1:m)
133  end if
134  end do
135 end do
136 deallocate(vgpl,cf_,x)
137 end subroutine
138 
pure subroutine r3mtv(a, x, y)
Definition: r3mtv.f90:10
logical ramdisk
Definition: modramdisk.f90:9
integer, dimension(48) isymlat
Definition: modmain.f90:348
type(file_t), dimension(:), allocatable, private file
Definition: modramdisk.f90:29
integer, dimension(3, 3, 48) symlat
Definition: modmain.f90:344
subroutine getcfgq(fname, vpl, ng, m, cf)
Definition: getcfgq.f90:7
logical, dimension(maxsymcrys) tv0symc
Definition: modmain.f90:362
subroutine getrd(fname, irec, tgs, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
Definition: modramdisk.f90:214
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
subroutine findqpt(vpl, isym, iq)
Definition: findqpt.f90:7
real(8), dimension(:,:), allocatable vql
Definition: modmain.f90:545
integer, dimension(:,:), allocatable ivg
Definition: modmain.f90:400
real(8) epslat
Definition: modmain.f90:24
real(8), dimension(3, maxsymcrys) vtcsymc
Definition: modmain.f90:360