The Elk Code
symrfmt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 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 symrfmt(nrmt_,nrmti_,npmt_,ld,rfmt)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: nrmt_(nspecies),nrmti_(nspecies),npmt_(nspecies),ld
11 real(8), intent(inout) :: rfmt(ld,natmtot)
12 ! local variables
13 integer is,ia,ja,ias,jas
14 integer nr,nri,np,isym,lspl
15 real(8) t0
16 ! allocatable arrays
17 real(8), allocatable :: rfmt1(:,:),rfmt2(:)
18 allocate(rfmt1(ld,natmmax),rfmt2(ld))
19 t0=1.d0/dble(nsymcrys)
20 do is=1,nspecies
21  nr=nrmt_(is)
22  nri=nrmti_(is)
23  np=npmt_(is)
24 ! make a copy of the input function
25  do ia=1,natoms(is)
26  ias=idxas(ia,is)
27  rfmt1(1:np,ia)=rfmt(1:np,ias)
28  end do
29 ! loop over atoms
30  do ia=1,natoms(is)
31 ! only symmetrise first equivalent atom (rotate into others)
32  if (.not.tfeqat(ia,is)) cycle
33  ias=idxas(ia,is)
34  rfmt(1:np,ias)=0.d0
35 ! loop over crystal symmetries
36  do isym=1,nsymcrys
37 ! index to spatial rotation lattice symmetry
38  lspl=lsplsymc(isym)
39 ! equivalent atom index (symmetry rotates atom ja into atom ia)
40  ja=ieqatom(ia,is,isym)
41 ! apply the rotation to the muffin-tin function
42  call rotrfmt(symlatc(:,:,lspl),nr,nri,rfmt1(:,ja),rfmt2)
43 ! accumulate in original function array
44  rfmt(1:np,ias)=rfmt(1:np,ias)+rfmt2(1:np)
45  end do
46 ! normalise
47  rfmt(1:np,ias)=t0*rfmt(1:np,ias)
48 ! rotate into equivalent atoms
49  do ja=1,natoms(is)
50  if (eqatoms(ia,ja,is).and.(ia /= ja)) then
51  isym=findloc(ieqatom(ia,is,1:nsymcrys),ja,1)
52  jas=idxas(ja,is)
53 ! inverse symmetry (which rotates atom ia into atom ja)
54  lspl=isymlat(lsplsymc(isym))
55 ! rotate symmetrised function into equivalent muffin-tin
56  call rotrfmt(symlatc(:,:,lspl),nr,nri,rfmt(:,ias),rfmt(:,jas))
57  end if
58  end do
59  end do
60 end do
61 deallocate(rfmt1,rfmt2)
62 end subroutine
63 
subroutine symrfmt(nrmt_, nrmti_, npmt_, ld, rfmt)
Definition: symrfmt.f90:7
integer natmmax
Definition: modmain.f90:38
logical, dimension(:,:), allocatable tfeqat
Definition: modmain.f90:372
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
integer nsymcrys
Definition: modmain.f90:358
integer, dimension(48) isymlat
Definition: modmain.f90:348
logical, dimension(:,:,:), allocatable eqatoms
Definition: modmain.f90:370
integer, dimension(:,:,:), allocatable ieqatom
Definition: modmain.f90:368
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
real(8), dimension(3, 3, 48) symlatc
Definition: modmain.f90:350
subroutine rotrfmt(rot, nr, nri, rfmt1, rfmt2)
Definition: rotrfmt.f90:7
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36