TABLE OF CONTENTS


ABINIT/mkeuler [ Functions ]

[ Top ] [ Functions ]

NAME

 mkeuler

FUNCTION

 For a given symmetry operation, determines the corresponding Euler angles

COPYRIGHT

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

  rot(3,3)= symmetry matrix

OUTPUT

  cosalp=  cos(alpha) with alpha=Euler angle 1
  cosbeta= cos(beta)  with beta =Euler angle 2
  cosgam=  cos(gamma) with gamma=Euler angle 3
  isn= error code (0 if the routine exit normally)
  sinalp= sin(alpha) with alpha=Euler angle 1
  singam= sin(gamma) with gamma=Euler angle 3

NOTES

  - This file comes from the file crystal_symmetry.f
    by N.A.W. Holzwarth and A. Tackett for the code pwpaw

PARENTS

      setsymrhoij

CHILDREN

SOURCE

 38 #if defined HAVE_CONFIG_H
 39 #include "config.h"
 40 #endif
 41 
 42 #include "abi_common.h"
 43 
 44 
 45 subroutine mkeuler(rot,cosbeta,cosalp,sinalp,cosgam,singam,isn)
 46 
 47  use defs_basis
 48  use m_profiling_abi
 49  use m_errors
 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 'mkeuler'
 55 !End of the abilint section
 56 
 57  implicit none
 58 
 59 !Arguments ---------------------------------------------
 60 !scalars
 61  integer,intent(out) :: isn
 62  real(dp),intent(out) :: cosalp,cosbeta,cosgam,sinalp,singam
 63 !arrays
 64  real(dp),intent(in) :: rot(3,3)
 65 
 66 !Local variables ---------------------------------------
 67 !scalars
 68  integer :: ier
 69  real(dp) :: check,sinbeta
 70  character(len=500) :: message
 71 
 72 ! *********************************************************************
 73 
 74  do isn= -1,1,2
 75    cosbeta=real(isn)*rot(3,3)
 76    if(abs(1._dp-cosbeta*cosbeta)<tol10) then
 77      sinbeta=zero
 78    else
 79      sinbeta=sqrt(1._dp-cosbeta*cosbeta)
 80    end if
 81    if (abs(sinbeta).gt.tol10)  then
 82      cosalp=isn*rot(3,1)/sinbeta
 83      sinalp=isn*rot(3,2)/sinbeta
 84      cosgam=-isn*rot(1,3)/sinbeta
 85      singam=isn*rot(2,3)/sinbeta
 86    else
 87      cosalp=isn*rot(1,1)/cosbeta
 88      sinalp=isn*rot(1,2)/cosbeta
 89      cosgam=one
 90      singam=zero
 91    end if
 92 
 93 !  Check matrix:
 94    ier=0
 95    check=cosalp*cosbeta*cosgam-sinalp*singam
 96    if (abs(check-isn*rot(1,1))>tol8) ier=ier+1
 97    check=sinalp*cosbeta*cosgam+cosalp*singam
 98    if (abs(check-isn*rot(1,2))>tol8) ier=ier+1
 99    check=-sinbeta*cosgam
100    if (abs(check-isn*rot(1,3))>tol8) ier=ier+1
101    check=-cosalp*cosbeta*singam-sinalp*cosgam
102    if (abs(check-isn*rot(2,1))>tol8) ier=ier+1
103    check=-sinalp*cosbeta*singam+cosalp*cosgam
104    if (abs(check-isn*rot(2,2))>tol8) ier=ier+1
105    check=sinbeta*singam
106    if (abs(check-isn*rot(2,3))>tol8) ier=ier+1
107    check=cosalp*sinbeta
108    if (abs(check-isn*rot(3,1))>tol8) ier=ier+1
109    check=sinalp*sinbeta
110    if (abs(check-isn*rot(3,2))>tol8) ier=ier+1
111    if (ier.eq.0) return
112  end do
113 
114  isn=0
115  write(message, '(7a)' )&
116 & 'Error during determination of symetries !',ch10,&
117 & 'Action: check your input file:',ch10,&
118 & '        unit cell vectors and/or atoms positions',ch10,&
119 & '        have to be given with a better precision...'
120  MSG_ERROR(message)
121 
122 end subroutine mkeuler