TABLE OF CONTENTS


ABINIT/matpointsym [ Functions ]

[ Top ] [ Functions ]

NAME

 matpointsym

FUNCTION

 For given order of point group, symmetrizes a 3x3 input matrix using the
 point symmetry of the input atom

COPYRIGHT

 Copyright (C) 2007-2018 ABINIT group (JWZ)
 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

 iatom=index of atom to symmetrize around
 natom=number of atoms in cell
 nsym=order of group
 rprimd(3,3)= real space primitive vectors
 symrel(3,3,nsym)=symmetry operators in terms of action on primitive translations
 tnons(3,nsym) = nonsymmorphic translations
 xred(3,natom)=locations of atoms in reduced coordinates

OUTPUT

SIDE EFFECTS

 mat3(3,3) = matrix to be symmetrized, in cartesian frame

PARENTS

      make_efg_el,make_efg_ion,make_efg_onsite

CHILDREN

      dgemm,dgemv,mati3inv,matrginv

SOURCE

 39 #if defined HAVE_CONFIG_H
 40 #include "config.h"
 41 #endif
 42 
 43 #include "abi_common.h"
 44 
 45 
 46 subroutine matpointsym(iatom,mat3,natom,nsym,rprimd,symrel,tnons,xred)
 47 
 48  use defs_basis
 49  use m_profiling_abi
 50  use m_linalg_interfaces
 51 
 52  use m_abilasi,         only : matrginv
 53 
 54 !This section has been created automatically by the script Abilint (TD).
 55 !Do not modify the following lines by hand.
 56 #undef ABI_FUNC
 57 #define ABI_FUNC 'matpointsym'
 58  use interfaces_32_util
 59 !End of the abilint section
 60 
 61  implicit none
 62 
 63 !Arguments ------------------------------------
 64 !scalars
 65  integer,intent(in) :: iatom,natom,nsym
 66 !arrays
 67  integer,intent(in) :: symrel(3,3,nsym)
 68  real(dp),intent(in) :: rprimd(3,3),tnons(3,nsym),xred(3,natom)
 69  real(dp),intent(inout) :: mat3(3,3)
 70 
 71 !Local variables-------------------------------
 72 !scalars
 73  integer :: cell_index,cell_indexp,ii,isym,nsym_point
 74  real(dp) :: xreddiff
 75 !arrays
 76  integer :: symrel_it(3,3)
 77  real(dp) :: mat3_tri(3,3),mat3_tri_sym(3,3),rprimd_inv(3,3),tmp_mat(3,3)
 78  real(dp) :: xredp(3)
 79 
 80 !**************************************************************************
 81 
 82 !DEBUG
 83 !write(std_out,*)' matpointsym : enter '
 84 !enddo
 85 !ENDDEBUG
 86 
 87 !copy rprimd input and construct inverse
 88  rprimd_inv = rprimd
 89  call matrginv(rprimd_inv,3,3)
 90 
 91 !transform input mat3 to triclinic frame with rprimd^{-1} * mat3 * rprimd
 92  call dgemm('N','N',3,3,3,one,rprimd_inv,3,mat3,3,zero,tmp_mat,3)
 93  call dgemm('N','N',3,3,3,one,tmp_mat,3,rprimd,3,zero,mat3_tri,3)
 94 
 95 !loop over symmetry elements to obtain symmetrized input matrix
 96  mat3_tri_sym = zero
 97  nsym_point = 0
 98  do isym = 1, nsym
 99 
100 ! skip any nonsymmorphic symmetry elements, want to consider point elements only
101    if(dot_product(tnons(:,isym),tnons(:,isym))>tol8) cycle
102 
103 ! for current symmetry element, find transformed reduced coordinates of target atom
104 ! via xredp = symrel * xred
105    call dgemv('N',3,3,one,dble(symrel(:,:,isym)),3,xred(:,iatom),1,zero,xredp,1)
106 
107 
108 ! shift xredp into the same unit cell as xred, for comparison
109 ! label cells as 0..1:0 1..2:1 2..3:2 and -1..0:-1 -2..-1:-2 and so forth
110    do ii = 1, 3
111 
112      cell_index = int(xred(ii,iatom))
113      if(xred(ii,iatom) < zero) cell_index = cell_index - 1
114      cell_indexp = int(xredp(ii))
115      if(xredp(ii) < zero) cell_indexp = cell_indexp - 1
116 
117      do while (cell_indexp < cell_index)
118        xredp(ii) = xredp(ii)+one
119        cell_indexp = cell_indexp + 1
120      end do
121      do while (cell_indexp > cell_index)
122        xredp(ii) = xredp(ii)-one
123        cell_indexp = cell_indexp - 1
124      end do
125 
126    end do
127 
128 ! now compare xredp to xred
129    xreddiff = dot_product(xredp-xred(:,iatom),xredp-xred(:,iatom))
130 
131    if (xreddiff < tol8) then
132 
133 !  accumulate symrel^{-1}*mat3_tri*symrel into mat3_tri_sym iff xredp = xred + L,
134 !  where is a lattice vector, so symrel leaves the target atom invariant
135 
136 !  mati3inv gives the inverse transpose of symrel
137      call mati3inv(symrel(:,:,isym),symrel_it)
138      call dgemm('N','N',3,3,3,one,mat3_tri,3,dble(symrel(:,:,isym)),3,zero,tmp_mat,3)
139      call dgemm('T','N',3,3,3,one,dble(symrel_it),3,tmp_mat,3,one,mat3_tri_sym,3)
140      nsym_point = nsym_point + 1
141    end if
142 
143  end do
144 
145 !normalize by number of point symmetry operations
146  mat3_tri_sym = mat3_tri_sym/dble(nsym_point)
147 
148 !transform mat3_tri_sym to cartesian frame with rprimd * mat3_tri_sym * rprimd^{-1}
149 
150  call dgemm('N','N',3,3,3,one,mat3_tri_sym,3,rprimd_inv,3,zero,tmp_mat,3)
151  call dgemm('N','N',3,3,3,one,rprimd,3,tmp_mat,3,zero,mat3,3)
152 
153 !DEBUG
154 !write(std_out,*)' matpointsym : exit '
155 !ENDDEBUG
156 
157 end subroutine matpointsym