TABLE OF CONTENTS


ABINIT/symg [ Functions ]

[ Top ] [ Functions ]

NAME

 symg

FUNCTION

 Treat symmetries applied to the G vectors, in view of the application
 to symmetrization of the dielectric matrix.
 Generate a list of time-reversed G vectors, as well as a list
 of spatially-symmetric G vectors.

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (XG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

 kg_diel(3,npwdiel)=reduced planewave coordinates for the dielectric matrix.
 npwdiel=number of planewaves for the dielectric matrix
 nsym=number of symmetry
 symrel(3,3,nsym)=symmetry matrices in real space (integers)
 tnons(3,nsym)=reduced nonsymmorphic translations
 (symrel and tnons are in terms of real space primitive translations)

OUTPUT

 phdiel(2,npwdiel,nsym)=phase associated with point symmetries applied to G
 sym_g(npwdiel,nsym)=index list of symmetric G vectors
 (could save a bit of space by suppressing isym=1, since the
 corresponding symmetry is the identity)
 tmrev_g(npwdiel)=index list of inverted G vectors (time-reversed)

PARENTS

      suscep_stat

CHILDREN

SOURCE

 41 #if defined HAVE_CONFIG_H
 42 #include "config.h"
 43 #endif
 44 
 45 #include "abi_common.h"
 46 
 47 subroutine symg(kg_diel,npwdiel,nsym,phdiel,sym_g,symrel,tmrev_g,tnons)
 48 
 49  use defs_basis
 50  use m_profiling_abi
 51  use m_errors
 52 
 53 !This section has been created automatically by the script Abilint (TD).
 54 !Do not modify the following lines by hand.
 55 #undef ABI_FUNC
 56 #define ABI_FUNC 'symg'
 57 !End of the abilint section
 58 
 59  implicit none
 60 
 61 !Arguments ------------------------------------
 62 !scalars
 63  integer,intent(in) :: npwdiel,nsym
 64 !arrays
 65  integer,intent(in) :: kg_diel(3,npwdiel),symrel(3,3,nsym)
 66  integer,intent(out) :: sym_g(npwdiel,nsym),tmrev_g(npwdiel)
 67  real(dp),intent(in) :: tnons(3,nsym)
 68  real(dp),intent(out) :: phdiel(2,npwdiel,nsym)
 69 
 70 !Local variables-------------------------------
 71 !scalars
 72  integer :: g1,g2,g3,ipw,isym,j1,j2,j3,m1m,m1p,m2m,m2p,m3m,m3p,symmg,trevg
 73  real(dp) :: arg,tau1,tau2,tau3
 74  character(len=500) :: message
 75 !arrays
 76  integer,allocatable :: grid(:,:,:)
 77 
 78 ! *************************************************************************
 79 
 80 !Determines maximal bounds of the zone spanned by the planewaves
 81  m1m=0 ; m2m=0 ; m3m=0 ; m1p=0 ; m2p=0 ; m3p=0
 82  do ipw=1,npwdiel
 83    g1=kg_diel(1,ipw)
 84    g2=kg_diel(2,ipw)
 85    g3=kg_diel(3,ipw)
 86    if(g1<m1m)m1m=g1 ; if(g1>m1p)m1p=g1
 87    if(g2<m2m)m2m=g2 ; if(g2>m2p)m2p=g2
 88    if(g3<m3m)m3m=g3 ; if(g3>m3p)m3p=g3
 89  end do
 90 
 91 !Set up grid, that associate to each point the index of the
 92 !corresponding planewave, if there is one
 93  ABI_ALLOCATE(grid,(m1m:m1p,m2m:m2p,m3m:m3p))
 94  grid(:,:,:)=0
 95  do ipw=1,npwdiel
 96    g1=kg_diel(1,ipw)
 97    g2=kg_diel(2,ipw)
 98    g3=kg_diel(3,ipw)
 99    grid(g1,g2,g3)=ipw
100  end do
101 
102 !Set up tmrev_g and sym_g arrays
103  do ipw=1,npwdiel
104    g1=kg_diel(1,ipw)
105    g2=kg_diel(2,ipw)
106    g3=kg_diel(3,ipw)
107 
108 !  Treat first time-reversal symmetry
109    trevg=grid(-g1,-g2,-g3)
110    if(trevg==0)then
111      message = ' Do not find the time-reversed symmetric of a G-vector.'
112      MSG_BUG(message)
113    end if
114    tmrev_g(ipw)=trevg
115 
116 !  Treat now spatial symmetries
117    do isym=1,nsym
118 
119 !    Get rotated G vector Gj for each symmetry element
120 !    -- here we use the TRANSPOSE of symrel; assuming symrel expresses
121 !    the rotation in real space, the transpose is then appropriate
122 !    for G space symmetrization (according to Doug : see routine irrzg.f)
123      j1=symrel(1,1,isym)*g1+&
124 &     symrel(2,1,isym)*g2+symrel(3,1,isym)*g3
125      j2=symrel(1,2,isym)*g1+&
126 &     symrel(2,2,isym)*g2+symrel(3,2,isym)*g3
127      j3=symrel(1,3,isym)*g1+&
128 &     symrel(2,3,isym)*g2+symrel(3,3,isym)*g3
129      symmg=grid(j1,j2,j3)
130      if(symmg==0)then
131        message = ' Do not find the spatially symmetric of a G-vector.'
132        MSG_BUG(message)
133      end if
134      sym_g(ipw,isym)=symmg
135 
136 !    Get associated phase
137      tau1=tnons(1,isym)
138      tau2=tnons(2,isym)
139      tau3=tnons(3,isym)
140      if (abs(tau1)>tol12.or.abs(tau2)>tol12.or.abs(tau3)>tol12) then
141 !      compute exp(-2*Pi*I*G dot tau) using original G
142        arg=two_pi*(dble(g1)*tau1+dble(g2)*tau2+dble(g3)*tau3)
143        phdiel(1,ipw,isym)=cos(arg)
144        phdiel(2,ipw,isym)=-sin(arg)
145      else
146        phdiel(1,ipw,isym)=1._dp
147        phdiel(2,ipw,isym)=0._dp
148      end if
149 
150    end do
151 
152  end do
153 
154  ABI_DEALLOCATE(grid)
155 
156 end subroutine symg