TABLE OF CONTENTS


ABINIT/matr3eigval [ Functions ]

[ Top ] [ Functions ]

NAME

 matr3eigval

FUNCTION

 Find the eigenvalues of a real symmetric 3x3 matrix,
 entered in full storage mode.

COPYRIGHT

 Copyright (C) 2003-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 .

INPUTS

  matr(3,3)=real symmetric 3x3 matrix

OUTPUT

  eigval(3)=three eigenvalues

PARENTS

      chkdilatmx

CHILDREN

      zhpev

SOURCE

30 #if defined HAVE_CONFIG_H
31 #include "config.h"
32 #endif
33 
34 #include "abi_common.h"
35 
36 
37 subroutine matr3eigval(eigval,matr)
38 
39  use defs_basis
40 
41 !This section has been created automatically by the script Abilint (TD).
42 !Do not modify the following lines by hand.
43 #undef ABI_FUNC
44 #define ABI_FUNC 'matr3eigval'
45 !End of the abilint section
46 
47  implicit none
48 
49 !Arguments ------------------------------------
50 !arrays
51  real(dp),intent(in) :: matr(3,3)
52  real(dp),intent(out) :: eigval(3)
53 
54 !Local variables-------------------------------
55 !scalars
56  integer :: ier
57 !arrays
58  real(dp) :: eigvec(2,3,3),matrx(2,6),zhpev1(2,2*3-1),zhpev2(3*3-2)
59 
60 ! *************************************************************************
61 
62  matrx(1,1)=matr(1,1)
63  matrx(1,2)=matr(1,2)
64  matrx(1,3)=matr(2,2)
65  matrx(1,4)=matr(1,3)
66  matrx(1,5)=matr(2,3)
67  matrx(1,6)=matr(3,3)
68  matrx(2,:)=zero
69 
70  call ZHPEV ('V','U',3,matrx,eigval,eigvec,3,zhpev1,zhpev2,ier)
71 !write(std_out,*)' eigval=',eigval
72 
73 end subroutine matr3eigval