The Elk Code
genhvec.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2010 Alexey I. Baranov.
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 genhvec
7 use modmain
8 use modpw
9 implicit none
10 ! local variables
11 logical lsym(48)
12 integer ih,jh,kh,lh,k
13 integer i1,i2,i3,iv(3)
14 integer nsym,isym,sym(3,3,48)
15 real(8) v1(3),v2(3),v3(3)
16 ! allocatable arrays
17 integer, allocatable :: idx(:),ivh0(:,:)
18 real(8), allocatable :: vhc0(:,:),hc0(:)
19 ! find the H-vector grid sizes
21 ! allocate global H-vector arrays
22 if (allocated(ivh)) deallocate(ivh)
23 allocate(ivh(3,nhtot))
24 if (allocated(mulh)) deallocate(mulh)
25 allocate(mulh(nhtot))
26 if (allocated(vhc)) deallocate(vhc)
27 allocate(vhc(3,nhtot))
28 if (allocated(hc)) deallocate(hc)
29 allocate(hc(nhtot))
30 ! allocate local arrays
31 allocate(idx(nhtot),ivh0(3,nhtot))
32 allocate(vhc0(3,nhtot),hc0(nhtot))
33 ih=0
34 do i1=inthv(1,1),inthv(2,1)
35  v1(:)=dble(i1)*bvec(:,1)
36  do i2=inthv(1,2),inthv(2,2)
37  v2(:)=v1(:)+dble(i2)*bvec(:,2)
38  do i3=inthv(1,3),inthv(2,3)
39  v3(:)=v2(:)+dble(i3)*bvec(:,3)
40  ih=ih+1
41 ! map from H-vector to (i1,i2,i3) index
42  ivh0(1,ih)=i1
43  ivh0(2,ih)=i2
44  ivh0(3,ih)=i3
45 ! H-vector in Cartesian coordinates
46  vhc0(:,ih)=v3(:)
47 ! length of each H-vector
48  hc0(ih)=sqrt(v3(1)**2+v3(2)**2+v3(3)**2)
49  end do
50  end do
51 end do
52 ! sort by vector length
53 call sortidx(nhtot,hc0,idx)
54 ! reorder arrays
55 do ih=1,nhtot
56  jh=idx(ih)
57  ivh(:,ih)=ivh0(:,jh)
58  hc(ih)=hc0(jh)
59  vhc(:,ih)=vhc0(:,jh)
60 end do
61 ! find the number of vectors with H < hmaxvr
62 nhvec=1
63 do ih=nhtot,1,-1
64  if (hc(ih) < hmaxvr) then
65  nhvec=ih
66  exit
67  end if
68 end do
69 ! find the subgroup of symmorphic, non-magnetic symmetries
70 lsym(:)=.false.
71 do isym=1,nsymcrys
72  if (tv0symc(isym).and.(lspnsymc(isym) == 1)) lsym(lsplsymc(isym))=.true.
73 end do
74 nsym=0
75 do isym=1,nsymlat
76  if (lsym(isym)) then
77  nsym=nsym+1
78  sym(:,:,nsym)=symlat(:,:,isym)
79  end if
80 end do
81 if (reduceh) then
82 ! find the subgroup of symmorphic, non-magnetic symmetries
83  lsym(:)=.false.
84  do isym=1,nsymcrys
85  if (tv0symc(isym).and.(lspnsymc(isym) == 1)) lsym(lsplsymc(isym))=.true.
86  end do
87  nsym=0
88  do isym=1,nsymlat
89  if (lsym(isym)) then
90  nsym=nsym+1
91  sym(:,:,nsym)=symlat(:,:,isym)
92  end if
93  end do
94 else
95 ! use only the identity element if no reduction is required
96  nsym=1
97 end if
98 ! reduce the H-vector set with the symmetries if required
99 if (nsym > 1) then
100  ivh0(:,1:nhvec)=ivh(:,1:nhvec)
101  hc0(1:nhvec)=hc(1:nhvec)
102  vhc0(:,1:nhvec)=vhc(:,1:nhvec)
103  kh=0
104  lh=nhvec
105  do ih=1,nhvec
106  do isym=1,nsym
107  call i3mtv(sym(:,:,isym),ivh0(:,ih),iv(:))
108  do jh=1,kh
109  k=abs(ivh(1,jh)-iv(1))+abs(ivh(2,jh)-iv(2))+abs(ivh(3,jh)-iv(3))
110  if (k == 0) then
111  ivh(:,lh)=ivh0(:,ih)
112  hc(lh)=hc0(ih)
113  vhc(:,lh)=vhc0(:,ih)
114  lh=lh-1
115  mulh(jh)=mulh(jh)+1
116  goto 10
117  end if
118  end do
119  end do
120  kh=kh+1
121  ivh(:,kh)=ivh0(:,ih)
122  hc(kh)=hc0(ih)
123  vhc(:,kh)=vhc0(:,ih)
124  mulh(kh)=1
125 10 continue
126  end do
127  nhvec=kh
128 else
129  mulh(:)=1
130 end if
131 deallocate(idx,ivh0,vhc0,hc0)
132 end subroutine
133 
integer, dimension(maxsymcrys) lspnsymc
Definition: modmain.f90:366
integer npfftg
Definition: modmain.f90:404
subroutine genhvec
Definition: genhvec.f90:7
integer nsymcrys
Definition: modmain.f90:358
integer, dimension(3, 3, 48) symlat
Definition: modmain.f90:344
pure subroutine i3mtv(a, x, y)
Definition: i3mtv.f90:10
logical, dimension(maxsymcrys) tv0symc
Definition: modmain.f90:362
real(8) hmaxvr
Definition: modpw.f90:14
subroutine gridsize(avec, gmaxvr, npfft, ngridg, ngtot, intgv)
Definition: gridsize.f90:10
real(8), dimension(:,:), allocatable vhc
Definition: modpw.f90:28
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
integer nsymlat
Definition: modmain.f90:342
logical reduceh
Definition: modpw.f90:12
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
integer nhtot
Definition: modpw.f90:18
pure subroutine sortidx(n, x, idx)
Definition: sortidx.f90:10
integer, dimension(:), allocatable mulh
Definition: modpw.f90:26
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
integer, dimension(2, 3) inthv
Definition: modpw.f90:20
real(8), dimension(:), allocatable hc
Definition: modpw.f90:30
Definition: modpw.f90:6
integer, dimension(:,:), allocatable ivh
Definition: modpw.f90:24
integer, dimension(3) ngridh
Definition: modpw.f90:16
integer nhvec
Definition: modpw.f90:22