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".

COPYRIGHT

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

  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

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

CHILDREN

SOURCE

 40 #if defined HAVE_CONFIG_H
 41 #include "config.h"
 42 #endif
 43 
 44 #include "abi_common.h"
 45 
 46 
 47 subroutine strconv(frac,gprimd,cart)
 48 
 49  use defs_basis
 50  use m_profiling_abi
 51 
 52 !This section has been created automatically by the script Abilint (TD).
 53 !Do not modify the following lines by hand.
 54 #undef ABI_FUNC
 55 #define ABI_FUNC 'strconv'
 56 !End of the abilint section
 57 
 58  implicit none
 59 
 60 !Arguments ------------------------------------
 61 !arrays
 62  real(dp),intent(in) :: frac(6),gprimd(3,3)
 63  real(dp),intent(inout) :: cart(6) ! alias of frac   !vz_i
 64 
 65 !Local variables-------------------------------
 66 !scalars
 67  integer :: ii,jj
 68 !arrays
 69  real(dp) :: work1(3,3),work2(3,3)
 70 
 71 ! *************************************************************************
 72 
 73  work1(1,1)=frac(1)
 74  work1(2,2)=frac(2)
 75  work1(3,3)=frac(3)
 76  work1(3,2)=frac(4) ; work1(2,3)=frac(4)
 77  work1(3,1)=frac(5) ; work1(1,3)=frac(5)
 78  work1(2,1)=frac(6) ; work1(1,2)=frac(6)
 79 
 80  do ii=1,3
 81    work2(:,ii)=zero
 82    do jj=1,3
 83      work2(:,ii)=work2(:,ii)+gprimd(ii,jj)*work1(:,jj)
 84    end do
 85  end do
 86 
 87  do ii=1,3
 88    work1(ii,:)=zero
 89    do jj=1,3
 90      work1(ii,:)=work1(ii,:)+gprimd(ii,jj)*work2(jj,:)
 91    end do
 92  end do
 93 
 94  cart(1)=work1(1,1)
 95  cart(2)=work1(2,2)
 96  cart(3)=work1(3,3)
 97  cart(4)=work1(2,3)
 98  cart(5)=work1(1,3)
 99  cart(6)=work1(1,2)
100 
101 end subroutine strconv