TABLE OF CONTENTS


ABINIT/create_mlms2jmj [ Functions ]

[ Top ] [ Functions ]

NAME

 create_mlms2jmj

FUNCTION

 For a given angular momentum lcor, give the rotation matrix msml2jmj

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

SIDE EFFECTS

  mlms2jmj= rotation matrix

NOTES

PARENTS

CHILDREN

SOURCE

 29 #if defined HAVE_CONFIG_H
 30 #include "config.h"
 31 #endif
 32 
 33 #include "abi_common.h"
 34 
 35 subroutine create_mlms2jmj(lcor,mlmstwojmj)
 36 
 37  use defs_basis
 38  use defs_datatypes
 39  use m_errors
 40  use m_profiling_abi
 41 
 42 !This section has been created automatically by the script Abilint (TD).
 43 !Do not modify the following lines by hand.
 44 #undef ABI_FUNC
 45 #define ABI_FUNC 'create_mlms2jmj'
 46 !End of the abilint section
 47 
 48  implicit none
 49 
 50 !Arguments ---------------------------------------------
 51 !scalars
 52  integer,intent(in) :: lcor
 53 !arrays
 54  complex(dpc),intent(out) :: mlmstwojmj(2*(2*lcor+1),2*(2*lcor+1))
 55 
 56 !Local variables ---------------------------------------
 57 !scalars
 58  integer :: jc1,jj,jm,ll,ml1,ms1
 59  real(dp) :: invsqrt2lp1,xj,xmj
 60  character(len=500) :: message
 61 !arrays
 62  integer, allocatable :: ind_msml(:,:)
 63  complex(dpc),allocatable :: mat_mlms2(:,:)
 64 
 65 !*********************************************************************
 66 
 67 !--------------- Built indices + allocations
 68  ll=lcor
 69  mlmstwojmj=czero
 70  ABI_ALLOCATE(ind_msml,(2,-ll:ll))
 71  ABI_ALLOCATE(mat_mlms2,(2*(2*lcor+1),2*(2*lcor+1)))
 72  mlmstwojmj=czero
 73  jc1=0
 74  do ms1=1,2
 75    do ml1=-ll,ll
 76      jc1=jc1+1
 77      ind_msml(ms1,ml1)=jc1
 78    end do
 79  end do
 80 !--------------- built mlmstwojmj
 81 !do jj=ll,ll+1    ! the physical value of j are ll-0.5,ll+0.5
 82 !xj(jj)=jj-0.5
 83  if(ll==0)then
 84    message=' ll should not be equal to zero !'
 85    MSG_BUG(message)
 86  end if
 87  jc1=0
 88  invsqrt2lp1=one/sqrt(float(2*lcor+1))
 89  do jj=ll,ll+1
 90    xj=float(jj)-half
 91    do jm=-jj,jj-1
 92      xmj=float(jm)+half
 93      jc1=jc1+1
 94      if(nint(xj+0.5)==ll+1) then
 95        if(nint(xmj+0.5)==ll+1)  then
 96          mlmstwojmj(ind_msml(2,ll),jc1)=1.0   !  J=L+0.5 and m_J=L+0.5
 97        else if(nint(xmj-0.5)==-ll-1) then
 98          mlmstwojmj(ind_msml(1,-ll),jc1)=1.0   !  J=L+0.5 and m_J=-L-0.5
 99        else
100          mlmstwojmj(ind_msml(2,nint(xmj-0.5)),jc1)=invsqrt2lp1*(sqrt(float(ll)+xmj+0.5))
101          mlmstwojmj(ind_msml(1,nint(xmj+0.5)),jc1)=invsqrt2lp1*(sqrt(float(ll)-xmj+0.5))
102        end if
103      end if
104      if(nint(xj+0.5)==ll) then
105        mlmstwojmj(ind_msml(1,nint(xmj+0.5)),jc1)=invsqrt2lp1*(sqrt(float(ll)+xmj+0.5))
106        mlmstwojmj(ind_msml(2,nint(xmj-0.5)),jc1)=-invsqrt2lp1*(sqrt(float(ll)-xmj+0.5))
107      end if
108    end do
109  end do
110 
111  ABI_DEALLOCATE(ind_msml)
112  ABI_DEALLOCATE(mat_mlms2)
113 
114  end subroutine create_mlms2jmj