TABLE OF CONTENTS


ABINIT/create_slm2ylm [ Functions ]

[ Top ] [ Functions ]

NAME

 create_slm2ylm

FUNCTION

 For a given angular momentum lcor, compute slm2ylm

COPYRIGHT

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

  lcor= angular momentum, size of the matrix is 2(2*lcor+1)

OUTPUT

  slm2ylm(2lcor+1,2lcor+1) = rotation matrix.

NOTES

  usefull only in ndij==4

PARENTS

CHILDREN

SOURCE

30 #if defined HAVE_CONFIG_H
31 #include "config.h"
32 #endif
33 
34 #include "abi_common.h"
35 
36 subroutine create_slm2ylm(lcor,slmtwoylm)
37 
38  use defs_basis
39  use defs_datatypes
40  use m_errors
41  use m_profiling_abi
42 
43 !This section has been created automatically by the script Abilint (TD).
44 !Do not modify the following lines by hand.
45 #undef ABI_FUNC
46 #define ABI_FUNC 'create_slm2ylm'
47 !End of the abilint section
48 
49  implicit none
50 
51 !Arguments ---------------------------------------------
52 !scalars
53  integer,intent(in) :: lcor
54 !arrays
55  complex(dpc),intent(out) :: slmtwoylm(2*lcor+1,2*lcor+1)
56 
57 !Local variables ---------------------------------------
58 !scalars
59  integer :: jm,ll,mm,im
60  real(dp),parameter :: invsqrt2=one/sqrt2
61  real(dp) :: onem
62 !arrays
63 
64 ! *********************************************************************
65 
66  ll=lcor
67  slmtwoylm=czero
68  do im=1,2*ll+1
69    mm=im-ll-1;jm=-mm+ll+1
70    onem=dble((-1)**mm)
71    if (mm> 0) then
72      slmtwoylm(im,im)= cmplx(onem*invsqrt2,zero,kind=dp)
73      slmtwoylm(jm,im)= cmplx(invsqrt2,     zero,kind=dp)
74    end if
75    if (mm==0) then
76      slmtwoylm(im,im)=cone
77    end if
78    if (mm< 0) then
79      slmtwoylm(im,im)= cmplx(zero,     invsqrt2,kind=dp)
80      slmtwoylm(jm,im)=-cmplx(zero,onem*invsqrt2,kind=dp)
81    end if
82  end do
83 
84 end subroutine create_slm2ylm