29  use, 
intrinsic :: iso_fortran_env
 
   51  subroutine x_slater_calc (namespace, gr, space, exxop, st, kpoints, ex, vxc)
 
   52    type(namespace_t),           
intent(in)    :: namespace
 
   53    type(grid_t),                
intent(in)    :: gr
 
   54    class(space_t),              
intent(in)    :: space
 
   55    type(exchange_operator_t),   
intent(in)    :: exxop
 
   56    type(states_elec_t),         
intent(inout) :: st
 
   57    type(kpoints_t),             
intent(in)    :: kpoints
 
   58    real(real64),                
intent(inout) :: ex
 
   59    real(real64), 
optional,             
intent(inout) :: vxc(:,:)
 
   64      call dslater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
 
   66      call zslater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
 
   86    real(real64),  
intent(in)  :: dens(:)
 
   87    real(real64),  
intent(out) :: alpha, betar, betai
 
   89    real(real64), 
parameter :: tiny = 1e-12_real64
 
   90    real(real64) :: mz, mm
 
   92    mz = dens(1) - dens(2)
 
   94    mm = 
sqrt(mz**2 + 
m_four*(dens(3)**2 + dens(4)**2))
 
  104    alpha = 
sqrt((mm + abs(mz))/(
m_two * mm))
 
  121  subroutine rotate_to_local(mat, alpha, betar, betai, alpha2, beta2, rot_mat)
 
  122    real(real64),  
intent(in)  :: mat(:)
 
  123    real(real64),  
intent(in)  :: alpha, betar, betai, alpha2, beta2
 
  124    real(real64),  
intent(out) :: rot_mat(:)
 
  126    complex(real64) :: cross
 
  128    rot_mat(1) = alpha2 * mat(1) + beta2 * mat(2) + 
m_two * alpha * (betar * mat(3) + betai * mat(4))
 
  129    rot_mat(2) = alpha2 * mat(2) + beta2 * mat(1) - 
m_two * alpha * (betar * mat(3) + betai * mat(4))
 
  130    cross = (cmplx(betar, betai, real64) )**2 * cmplx(mat(3), -mat(4), real64)
 
  131    rot_mat(3) = alpha2 * mat(3) + alpha * betar * (mat(2)-mat(1)) - real(cross)
 
  132    rot_mat(4) = alpha2 * mat(4) + alpha * betai * (mat(2)-mat(1)) - aimag(cross)
 
  142  subroutine rotate_to_global(mat, alpha, betar, betai, alpha2, beta2, rot_mat)
 
  143    real(real64),  
intent(in)  :: mat(:)
 
  144    real(real64),  
intent(in)  :: alpha, betar, betai, alpha2, beta2
 
  145    real(real64),  
intent(out) :: rot_mat(:)
 
  147    complex(real64) :: cross
 
  149    rot_mat(1) = alpha2 * mat(1) + beta2 * mat(2) - 
m_two * alpha * (betar * mat(3) + betai * mat(4))
 
  150    rot_mat(2) = alpha2 * mat(2) + beta2 * mat(1) + 
m_two * alpha * (betar * mat(3) + betai * mat(4))
 
  151    cross = (cmplx(betar, betai, real64) )**2 * cmplx(mat(3), -mat(4), real64)
 
  152    rot_mat(3) = alpha2 * mat(3) - alpha * betar * (mat(2)-mat(1)) - real(cross)
 
  153    rot_mat(4) = alpha2 * mat(4) - alpha * betai * (mat(2)-mat(1)) - aimag(cross)
 
  160#include "x_slater_inc.F90" 
  163#include "complex.F90" 
  164#include "x_slater_inc.F90" 
double sqrt(double __x) __attribute__((__nothrow__
 
real(real64), parameter, public m_two
 
real(real64), parameter, public m_zero
 
real(real64), parameter, public m_four
 
real(real64), parameter, public m_one
 
This module implements the underlying real-space grid.
 
This module defines various routines, operating on mesh functions.
 
pure logical function, public states_are_real(st)
 
This module handles spin dimensions of the states and the k-point distribution.
 
subroutine rotate_to_local(mat, alpha, betar, betai, alpha2, beta2, rot_mat)
Rotation to the local frame Given a matrix in spin space, this routine rotates according to the rotat...
 
subroutine rotate_to_global(mat, alpha, betar, betai, alpha2, beta2, rot_mat)
Rotation to the global (Cartesian) frame.
 
subroutine zslater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
 
subroutine, public x_slater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
Interface to X(slater_calc)
 
subroutine dslater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
 
subroutine get_rotation_matrix(dens, alpha, betar, betai)
This routine get the SU(2) matrix that diagonalize the spin-density matrix.