TABLE OF CONTENTS


ABINIT/strconv [ Functions ]

[ Top ] [ Functions ]

NAME

 strconv

FUNCTION

 If original gprimd is input, convert from symmetric storage mode
 3x3 tensor in reduced coordinates "frac" to symmetric storage mode
 symmetric tensor in cartesian coordinates "cart".

INPUTS

  frac(6)=3x3 tensor in symmetric storage mode, reduced coordinates
  gprimd(3,3)=reciprocal space dimensional primitive translations (bohr^-1)

OUTPUT

  cart(6)=symmetric storage mode for symmetric 3x3 tensor in cartesian coords.

NOTES

 $cart(i,j)=G(i,a) G(j,b) frac(a,b)$
 "Symmetric" storage mode for 3x3 tensor is 6 element array with
 elements 11, 22, 33, 32, 31, and 21.
 "cart" may be same array as "frac".
 If rprimd transpose is input instead of gprimd, then convert tensor
 in cartesian coordinates to reduced coordinates

PARENTS

      ctocprj,d2frnl,mkcore,mkcore_paw,mkcore_wvl,nonlop_pl,nonlop_ylm
      stresssym

CHILDREN

SOURCE

165 subroutine strconv(frac,gprimd,cart)
166 
167  use defs_basis
168  use m_profiling_abi
169 
170 !This section has been created automatically by the script Abilint (TD).
171 !Do not modify the following lines by hand.
172 #undef ABI_FUNC
173 #define ABI_FUNC 'strconv'
174 !End of the abilint section
175 
176  implicit none
177 
178 !Arguments ------------------------------------
179 !arrays
180  real(dp),intent(in) :: frac(6),gprimd(3,3)
181  real(dp),intent(inout) :: cart(6) ! alias of frac   !vz_i
182 
183 !Local variables-------------------------------
184 !scalars
185  integer :: ii,jj
186 !arrays
187  real(dp) :: work1(3,3),work2(3,3)
188 
189 ! *************************************************************************
190 
191  work1(1,1)=frac(1)
192  work1(2,2)=frac(2)
193  work1(3,3)=frac(3)
194  work1(3,2)=frac(4) ; work1(2,3)=frac(4)
195  work1(3,1)=frac(5) ; work1(1,3)=frac(5)
196  work1(2,1)=frac(6) ; work1(1,2)=frac(6)
197 
198  do ii=1,3
199    work2(:,ii)=zero
200    do jj=1,3
201      work2(:,ii)=work2(:,ii)+gprimd(ii,jj)*work1(:,jj)
202    end do
203  end do
204 
205  do ii=1,3
206    work1(ii,:)=zero
207    do jj=1,3
208      work1(ii,:)=work1(ii,:)+gprimd(ii,jj)*work2(jj,:)
209    end do
210  end do
211 
212  cart(1)=work1(1,1)
213  cart(2)=work1(2,2)
214  cart(3)=work1(3,3)
215  cart(4)=work1(2,3)
216  cart(5)=work1(1,3)
217  cart(6)=work1(1,2)
218 
219 end subroutine strconv

ABINIT/stresssym [ Functions ]

[ Top ] [ Functions ]

NAME

 stresssym

FUNCTION

 For given order of point group, symmetrizes the stress tensor,
 in symmetrized storage mode and cartesian coordinates, using input
 3x3 symmetry operators in reduced coordinates.
 symmetrized tensor replaces input tensor.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR)
 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

 gprimd(3,3)=dimensional primitive translations for reciprocal space (bohr**-1)
 nsym=order of group.
 sym(3,3,nsym)=symmetry operators (usually symrec=expressed in terms
               of action on reciprocal lattice primitive translations);
               integers.

OUTPUT

 stress(6)=stress tensor, in cartesian coordinates, in symmetric storage mode

SIDE EFFECTS

PARENTS

      dfpt_nselt,dfpt_nstpaw,forstrnps,littlegroup_pert,pawgrnl,stress

CHILDREN

      matr3inv,strconv

SOURCE

 39 #if defined HAVE_CONFIG_H
 40 #include "config.h"
 41 #endif
 42 
 43 #include "abi_common.h"
 44 
 45 
 46 subroutine stresssym(gprimd,nsym,stress,sym)
 47 
 48  use defs_basis
 49  use m_profiling_abi
 50 
 51 !This section has been created automatically by the script Abilint (TD).
 52 !Do not modify the following lines by hand.
 53 #undef ABI_FUNC
 54 #define ABI_FUNC 'stresssym'
 55  use interfaces_32_util
 56  use interfaces_41_geometry, except_this_one => stresssym
 57 !End of the abilint section
 58 
 59  implicit none
 60 
 61 !Arguments ------------------------------------
 62 !scalars
 63  integer,intent(in) :: nsym
 64 !arrays
 65  integer,intent(in) :: sym(3,3,nsym)
 66  real(dp),intent(in) :: gprimd(3,3)
 67  real(dp),intent(inout) :: stress(6)
 68 
 69 !Local variables-------------------------------
 70 !scalars
 71  integer :: ii,isym,mu,nu
 72  real(dp) :: summ,tmp
 73 !arrays
 74  real(dp) :: rprimd(3,3),rprimdt(3,3),strfrac(6),tensor(3,3),tt(3,3)
 75 
 76 !*************************************************************************
 77 
 78 !Obtain matrix of real space dimensional primitive translations
 79 !(inverse tranpose of gprimd), and its transpose
 80  call matr3inv(gprimd,rprimd)
 81  rprimdt=transpose(rprimd)
 82 
 83 !Compute stress tensor in reduced coordinates
 84  call strconv(stress,rprimdt,strfrac)
 85 
 86 !Switch to full storage mode
 87  tensor(1,1)=strfrac(1)
 88  tensor(2,2)=strfrac(2)
 89  tensor(3,3)=strfrac(3)
 90  tensor(3,2)=strfrac(4)
 91  tensor(3,1)=strfrac(5)
 92  tensor(2,1)=strfrac(6)
 93  tensor(2,3)=tensor(3,2)
 94  tensor(1,3)=tensor(3,1)
 95  tensor(1,2)=tensor(2,1)
 96 
 97  do nu=1,3
 98    do mu=1,3
 99      tt(mu,nu)=tensor(mu,nu)/dble(nsym)
100      tensor(mu,nu)=0.0_dp
101    end do
102  end do
103 
104 !loop over all symmetry operations:
105  do isym=1,nsym
106    do mu=1,3
107      do nu=1,3
108        summ=0._dp
109        do ii=1,3
110          tmp=tt(ii,1)*sym(nu,1,isym)+tt(ii,2)*sym(nu,2,isym)+&
111 &         tt(ii,3)*sym(nu,3,isym)
112          summ=summ+sym(mu,ii,isym)*tmp
113        end do
114        tensor(mu,nu)=tensor(mu,nu)+summ
115      end do
116    end do
117  end do
118 
119 !Switch back to symmetric storage mode
120  strfrac(1)=tensor(1,1)
121  strfrac(2)=tensor(2,2)
122  strfrac(3)=tensor(3,3)
123  strfrac(4)=tensor(3,2)
124  strfrac(5)=tensor(3,1)
125  strfrac(6)=tensor(2,1)
126 
127 !Convert back stress tensor (symmetrized) in cartesian coordinates
128  call strconv(strfrac,gprimd,stress)
129 
130 end subroutine stresssym