@brief Generate transformation to spherical harmonics basis
@details This is used by the MINI basis during the Huckel guess, to convert any d or f shells to spherical harmonics.
@author MWS: 4/2013 rewrite returns pure s,p,d,f spherical
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp) | :: | w(l1co,l1co) | ||||
integer | :: | l1co | ||||
integer | :: | l0co | ||||
logical | :: | notsp | ||||
type(basis_set), | intent(in) | :: | basis |
subroutine mksphar(w,l1co,l0co,notsp,basis) use messages, only: show_message, WITH_ABORT use basis_tools, only: basis_set ! implicit none ! ! type(basis_set), intent(in) :: basis logical :: notsp real(kind=dp) :: w(l1co,l1co) integer :: l1co, l0co ! real(kind=dp), parameter :: RT0304 = SQRT( 3.0D+00/ 4.0D+00) real(kind=dp), parameter :: RT0920 = SQRT( 9.0D+00/20.0D+00) integer :: sph, ncont, n, cart, ndim, l ! ! Transform to spherical harmonic basis w = 0 sph = 1 ncont = 0 notsp = .false. do n = 1, basis%nshell notsp = notsp .or. (basis%am(n) > 1) ! ! S,P,L shell is already spherical harmonics. if (basis%am(n) <= 1) then cart = basis%ao_offset(n) ndim = basis%naos(n) do l = 1, ndim w(cart+l-1, sph+l-1) = 1.0d0 enddo sph = sph + ndim end if ! if (basis%am(n) == 2) then cart = basis%ao_offset(n) ! true D orbitals, from eg irrep of Oh w(cart ,sph ) = rt0304 w(cart+1,sph ) = -rt0304 w(cart ,sph+1) = -0.5d0 w(cart+1,sph+1) = -0.5d0 w(cart+2,sph+1) = 1.0d0 ! true d orbitals, from t2g irrep of oh w(cart+3,sph+2) = 1.0d0 w(cart+4,sph+3) = 1.0d0 w(cart+5,sph+4) = 1.0d0 ncont = ncont+1 sph = sph+5 end if ! if (basis%am(n) == 3) then cart = basis%ao_offset(n) ! T1U (IN OH) COMBINATIONS w(cart ,sph ) = 1.0d0 w(cart+5,sph ) = -rt0920 w(cart+7,sph ) = -rt0920 w(cart+1,sph+1) = 1.0d0 w(cart+3,sph+1) = -rt0920 w(cart+8,sph+1) = -rt0920 w(cart+2,sph+2) = 1.0d0 w(cart+4,sph+2) = -rt0920 w(cart+6,sph+2) = -rt0920 ! t2u (in oh) combinations w(cart+5,sph+3) = rt0304 w(cart+7,sph+3) = -rt0304 w(cart+3,sph+4) = rt0304 w(cart+8,sph+4) = -rt0304 w(cart+4,sph+5) = rt0304 w(cart+6,sph+5) = -rt0304 ! a2u (in oh) combination w(cart+9,sph+6) = 1.0d0 ncont = ncont+3 sph = sph+7 end if ! if(basis%am(n) > 3) then CALL show_message('MKSPHAR: CALLED FOR BASIS WITH G/H/I AOS', WITH_ABORT) end if end do l0co = l1co - ncont end subroutine mksphar