Octopus
atomic_oct_m Module Reference

Data Types

type  valconf_t
 

Functions/Subroutines

subroutine, public valconf_copy (cout, cin)
 
subroutine, public write_valconf (c, s)
 
subroutine, public read_valconf (namespace, s, c)
 
subroutine, public valconf_unpolarized_to_polarized (conf)
 
subroutine, public atomhxc (namespace, functl, g, nspin, dens, v, extra)
 
subroutine, public atomxc (namespace, FUNCTL, AUTHOR, NR, MAXR, RMESH, NSPIN, DENS, EX, EC, DX, DC, V_XC)
 Finds total exchange-correlation energy and potential for a ! spherical electron density distribution. ! This version implements the Local (spin) Density Approximation and ! the Generalized-Gradient-Aproximation with the explicit mesh ! functional method of White & Bird, PRB 50, 4954 (1994). ! Gradients are defined by numerical derivatives, using 2*NN+1 mesh ! points, where NN is a parameter defined below ! Coded by L.C.Balbas and J.M.Soler. December 1996. Version 0.5. ! ! CHARACTER*(*) FUNCTL : Functional to be used: ! LDA or LSD => Local (spin) Density Approximation ! GGA => Generalized Gradient Corrections ! Uppercase is optional ! CHARACTER*(*) AUTHOR : Parametrization desired: ! CA or PZ => LSD Perdew & Zunger, PRB 23, 5075 (1981) ! PW92 => LSD Perdew & Wang, PRB, 45, 13244 (1992). This is ! the local density limit of the next: ! PBE => GGA Perdew, Burke & Ernzerhof, PRL 77, 3865 (1996) ! Uppercase is optional ! INTEGER NR : Number of radial mesh points ! INTEGER MAXR : Physical first dimension of RMESH, DENS and V_XC ! REAL*8 RMESH(MAXR) : Radial mesh points ! INTEGER NSPIN : NSPIN=1 => unpolarized; NSPIN=2 => polarized ! REAL*8 DENS(MAXR,NSPIN) : Total (NSPIN=1) or spin (NSPIN=2) electron ! density at mesh points ! ************************* OUTPUT ********************************** ! REAL*8 EX : Total exchange energy ! REAL*8 EC : Total correlation energy ! REAL*8 DX : IntegralOf( rho * (eps_x - v_x) ) ! REAL*8 DC : IntegralOf( rho * (eps_c - v_c) ) ! REAL*8 V_XC(MAXR,NSPIN): (Spin) exch-corr potential ! ************************ UNITS ************************************ ! Distances in atomic units (Bohr). ! Densities in atomic units (electrons/Bohr**3) ! Energy unit depending of parameter EUNIT below ! ********* ROUTINES CALLED ***************************************** ! GGAXC, LDAXC ! More...
 
subroutine, public vhrtre (rho, v, r, drdi, srdrdi, nr, a)
 VHRTRE CONSTRUCTS THE ELECTROSTATIC POTENTIAL DUE TO A SUPPLIED ! ELECTRON DENSITY. THE NUMEROV METHOD IS USED TO INTEGRATE ! POISSONS EQN. ! ! DESCRIPTION OF ARGUMENTS: ! RHO....4*PI*R**2 * THE ELECTRON DENSITY FOR WHICH WE CALCULATING ! THE ELECTROSTATIC POTENTIAL ! V......THE ELECTROSTATIC POTENTIAL DUE TO THE ELECTRON DENSITY ! RHO. THE CONSTANTS OF INTEGRATION ARE FIXED SO THAT THE ! POTENTIAL TENDS TO A CONSTANT AT THE ORIGIN AND TO ! 2*Q/R AT R=R(NR), WHERE Q IS THE INTEGRATED CHARGE ! CONTAINED IN RHO(R) ! R......THE RADIAL MESH R(I) = B*(EXP(A(I-1))-1) ! NR.....THE NUMBER OF RADIAL MESH POINTS ! DRDI...DR(I)/DI ! SRDRDI.SQRT(DR/DI) ! A......THE PARAMETER APPEARING IN R(I) = B*(EXP(A(I-1))-1) ! More...
 
subroutine, public egofv (namespace, h, s, n, e, g, y, l, z, a, b, rmax, nprin, nnode, dr, ierr)
 egofv determines the eigenenergy and wavefunction corresponding ! to a particular l, principal quantum number and boundary condition. ! ! two fundamental techniques are used to locate the solution: ! 1) node counting and bisection ! 2) variational estimate based on a slope discontinuity in psi ! the arguments are defined as follows: ! h,s: g = (h-e*s)*g ! nr: maximum allowed number of radial points ! e: e is the energy found ! ne: number of energies found ! l: the angular momentum ! ncor: the number of lower-energy state ! ! the individual energies are resolved by performing a fixed number ! of bisections after a given eigenvalue has been isolated ! More...
 
subroutine yofe (e, de, dr, rmax, h, s, y, nmax, l, ncor, nnode, z, a, b)
 
subroutine nrmlzg (namespace, g, s, n)
 
subroutine bcorgn (e, h, s, l, zdr, y2)
 
subroutine bcrmax (e, dr, rmax, h, s, n, yn, a, b)
 
subroutine numin (e, h, s, y, n, nnode, yn, g, gsg, x, knk)
 
subroutine numout (e, h, s, y, ncor, knk, nnode, y2, g, gsg, x)
 

Variables

character(len=1), dimension(0:3), parameter spec_notation = (/ 's', 'p', 'd', 'f' /)
 Following stuff is for the "val_conf_t" data type, made to handle atomic configurations. More...
 
integer, parameter, public valconf_string_length = 80
 

Function/Subroutine Documentation

◆ valconf_copy()

subroutine, public atomic_oct_m::valconf_copy ( type(valconf_t), intent(out)  cout,
type(valconf_t), intent(in)  cin 
)

Definition at line 162 of file atomic.F90.

◆ write_valconf()

subroutine, public atomic_oct_m::write_valconf ( type(valconf_t), intent(in)  c,
character(len=valconf_string_length), intent(out)  s 
)

Definition at line 182 of file atomic.F90.

◆ read_valconf()

subroutine, public atomic_oct_m::read_valconf ( type(namespace_t), intent(in)  namespace,
character(len=valconf_string_length), intent(in)  s,
type(valconf_t), intent(out)  c 
)

Definition at line 198 of file atomic.F90.

◆ valconf_unpolarized_to_polarized()

subroutine, public atomic_oct_m::valconf_unpolarized_to_polarized ( type(valconf_t), intent(inout)  conf)

Definition at line 232 of file atomic.F90.

◆ atomhxc()

subroutine, public atomic_oct_m::atomhxc ( type(namespace_t), intent(in)  namespace,
character(len=*), intent(in)  functl,
type(logrid_t), intent(in)  g,
integer, intent(in)  nspin,
real(real64), dimension(g%nrval, nspin), intent(in)  dens,
real(real64), dimension(g%nrval, nspin), intent(out)  v,
real(real64), dimension(g%nrval), intent(in), optional  extra 
)

Definition at line 252 of file atomic.F90.

◆ atomxc()

subroutine, public atomic_oct_m::atomxc ( type(namespace_t), intent(in)  namespace,
character(len=*), intent(in)  FUNCTL,
character(len=*), intent(in)  AUTHOR,
integer, intent(in)  NR,
integer, intent(in)  MAXR,
real(real64), dimension(maxr), intent(in)  RMESH,
integer, intent(in)  NSPIN,
real(real64), dimension(maxr,nspin), intent(in)  DENS,
real(real64), intent(out)  EX,
real(real64), intent(out)  EC,
real(real64), intent(out)  DX,
real(real64), intent(out)  DC,
real(real64), dimension(maxr,nspin), intent(out)  V_XC 
)

Finds total exchange-correlation energy and potential for a ! spherical electron density distribution. ! This version implements the Local (spin) Density Approximation and ! the Generalized-Gradient-Aproximation with the explicit mesh ! functional method of White & Bird, PRB 50, 4954 (1994). ! Gradients are defined by numerical derivatives, using 2*NN+1 mesh ! points, where NN is a parameter defined below ! Coded by L.C.Balbas and J.M.Soler. December 1996. Version 0.5. ! ! CHARACTER*(*) FUNCTL : Functional to be used: ! LDA or LSD => Local (spin) Density Approximation ! GGA => Generalized Gradient Corrections ! Uppercase is optional ! CHARACTER*(*) AUTHOR : Parametrization desired: ! CA or PZ => LSD Perdew & Zunger, PRB 23, 5075 (1981) ! PW92 => LSD Perdew & Wang, PRB, 45, 13244 (1992). This is ! the local density limit of the next: ! PBE => GGA Perdew, Burke & Ernzerhof, PRL 77, 3865 (1996) ! Uppercase is optional ! INTEGER NR : Number of radial mesh points ! INTEGER MAXR : Physical first dimension of RMESH, DENS and V_XC ! REAL*8 RMESH(MAXR) : Radial mesh points ! INTEGER NSPIN : NSPIN=1 => unpolarized; NSPIN=2 => polarized ! REAL*8 DENS(MAXR,NSPIN) : Total (NSPIN=1) or spin (NSPIN=2) electron ! density at mesh points ! ************************* OUTPUT ********************************** ! REAL*8 EX : Total exchange energy ! REAL*8 EC : Total correlation energy ! REAL*8 DX : IntegralOf( rho * (eps_x - v_x) ) ! REAL*8 DC : IntegralOf( rho * (eps_c - v_c) ) ! REAL*8 V_XC(MAXR,NSPIN): (Spin) exch-corr potential ! ************************ UNITS ************************************ ! Distances in atomic units (Bohr). ! Densities in atomic units (electrons/Bohr**3) ! Energy unit depending of parameter EUNIT below ! ********* ROUTINES CALLED ***************************************** ! GGAXC, LDAXC !

Definition at line 361 of file atomic.F90.

◆ vhrtre()

subroutine, public atomic_oct_m::vhrtre ( real(real64), dimension(:), intent(in)  rho,
real(real64), dimension(:), intent(out)  v,
real(real64), dimension(:), intent(in)  r,
real(real64), dimension(:), intent(in)  drdi,
real(real64), dimension(:), intent(in)  srdrdi,
integer, intent(in)  nr,
real(real64), intent(in)  a 
)

VHRTRE CONSTRUCTS THE ELECTROSTATIC POTENTIAL DUE TO A SUPPLIED ! ELECTRON DENSITY. THE NUMEROV METHOD IS USED TO INTEGRATE ! POISSONS EQN. ! ! DESCRIPTION OF ARGUMENTS: ! RHO....4*PI*R**2 * THE ELECTRON DENSITY FOR WHICH WE CALCULATING ! THE ELECTROSTATIC POTENTIAL ! V......THE ELECTROSTATIC POTENTIAL DUE TO THE ELECTRON DENSITY ! RHO. THE CONSTANTS OF INTEGRATION ARE FIXED SO THAT THE ! POTENTIAL TENDS TO A CONSTANT AT THE ORIGIN AND TO ! 2*Q/R AT R=R(NR), WHERE Q IS THE INTEGRATED CHARGE ! CONTAINED IN RHO(R) ! R......THE RADIAL MESH R(I) = B*(EXP(A(I-1))-1) ! NR.....THE NUMBER OF RADIAL MESH POINTS ! DRDI...DR(I)/DI ! SRDRDI.SQRT(DR/DI) ! A......THE PARAMETER APPEARING IN R(I) = B*(EXP(A(I-1))-1) !

Definition at line 626 of file atomic.F90.

◆ egofv()

subroutine, public atomic_oct_m::egofv ( type(namespace_t), intent(in)  namespace,
real(real64), dimension(:), intent(in)  h,
real(real64), dimension(:), intent(in)  s,
integer, intent(in)  n,
real(real64), intent(inout)  e,
real(real64), dimension(:), intent(out)  g,
real(real64), dimension(:), intent(out)  y,
integer, intent(in)  l,
real(real64), intent(in)  z,
real(real64), intent(in)  a,
real(real64), intent(in)  b,
real(real64), intent(in)  rmax,
integer, intent(in)  nprin,
integer, intent(in)  nnode,
real(real64), intent(in)  dr,
integer, intent(out)  ierr 
)

egofv determines the eigenenergy and wavefunction corresponding ! to a particular l, principal quantum number and boundary condition. ! ! two fundamental techniques are used to locate the solution: ! 1) node counting and bisection ! 2) variational estimate based on a slope discontinuity in psi ! the arguments are defined as follows: ! h,s: g = (h-e*s)*g ! nr: maximum allowed number of radial points ! e: e is the energy found ! ne: number of energies found ! l: the angular momentum ! ncor: the number of lower-energy state ! ! the individual energies are resolved by performing a fixed number ! of bisections after a given eigenvalue has been isolated !

Parameters
[in]s(n)
[out]y(n)

Definition at line 764 of file atomic.F90.

◆ yofe()

subroutine atomic_oct_m::yofe ( real(real64), intent(inout)  e,
real(real64), intent(out)  de,
real(real64), intent(in)  dr,
real(real64), intent(in)  rmax,
real(real64), dimension(:), intent(in)  h,
real(real64), dimension(:), intent(in)  s,
real(real64), dimension(:), intent(out)  y,
integer, intent(in)  nmax,
integer, intent(in)  l,
integer, intent(in)  ncor,
integer, intent(out)  nnode,
real(real64), intent(in)  z,
real(real64), intent(in)  a,
real(real64), intent(in)  b 
)
private
Parameters
[in,out]eyofe integrates the radial Schroedinger eqn using the Numerov ! method. ! ! the arguments are defined as follows: ! e is the old energy (overwritten) by the new energy ! de is the e change predicted to elim the kink in psi ! dr is the log deriv (the boundary condition) ! gpp = (h-es)g (all diagonal in i (radius) ) ! y is the numerov independent variable y = g - gpp/12 ! n is the number of radial mesh points ! l is the angular momentum ! ncor is the number of states of lower energy ! nnode is 1 + the number of interior nodes in psi ! z is the atomic number ! a and b specify the radial mesh r(i)=(exp(a*(i-1))-1)*b !
[in]s(nmax)
[out]y(nmax)

Definition at line 904 of file atomic.F90.

◆ nrmlzg()

subroutine atomic_oct_m::nrmlzg ( type(namespace_t), intent(in)  namespace,
real(real64), dimension(:), intent(inout)  g,
real(real64), dimension(:), intent(in)  s,
integer, intent(in)  n 
)
private
Parameters
[in]namespacenrmlzg normalizes the supplied radial wavefunction ! ! the arguments are defined as follows: ! g is the radial wavefunction appropriate to the Numerov ! representation of the radial Schroedinger equation ! that is, the radial fcn r(r) = (drdi)**1/2 g(i) / r(i) ! gpp = (h-es)g (all diagonal in i (radius) ) ! n1 is the number of radial mesh points corresponding to ! the portion of the radial mesh on which the norm ! is defined ! n2 is the number of radial mesh points corresponding to ! the portion of the radial mesh on which the wavefunction ! is defined. For the intended use of this ! routine, n1 = nrval and n2 = nrcor ! a and b are the radial mesh parameters ! r(i) = ( exp(a*(i-1)) - 1 ) * b ! (dr/di = a*b at the origin) !

Definition at line 1016 of file atomic.F90.

◆ bcorgn()

subroutine atomic_oct_m::bcorgn ( real(real64), intent(in)  e,
real(real64), dimension(:), intent(in)  h,
real(real64), dimension(:), intent(in)  s,
integer, intent(in)  l,
real(real64), intent(in)  zdr,
real(real64), intent(out)  y2 
)
private

Definition at line 1083 of file atomic.F90.

◆ bcrmax()

subroutine atomic_oct_m::bcrmax ( real(real64), intent(in)  e,
real(real64), intent(in)  dr,
real(real64), intent(in)  rmax,
real(real64), dimension(:), intent(in)  h,
real(real64), dimension(:), intent(in)  s,
integer, intent(in)  n,
real(real64), intent(out)  yn,
real(real64), intent(in)  a,
real(real64), intent(in)  b 
)
private

Definition at line 1150 of file atomic.F90.

◆ numin()

subroutine atomic_oct_m::numin ( real(real64), intent(in)  e,
real(real64), dimension(:), intent(in)  h,
real(real64), dimension(:), intent(in)  s,
real(real64), dimension(:), intent(inout)  y,
integer, intent(in)  n,
integer, intent(out)  nnode,
real(real64), intent(in)  yn,
real(real64), intent(out)  g,
real(real64), intent(out)  gsg,
real(real64), intent(out)  x,
integer, intent(in)  knk 
)
private
Parameters
[in]s(n)
[in,out]y(n)

Definition at line 1178 of file atomic.F90.

◆ numout()

subroutine atomic_oct_m::numout ( real(real64), intent(in)  e,
real(real64), dimension(:), intent(in)  h,
real(real64), dimension(:), intent(in)  s,
real(real64), dimension(:), intent(inout)  y,
integer, intent(in)  ncor,
integer, intent(inout)  knk,
integer, intent(out)  nnode,
real(real64), intent(in)  y2,
real(real64), intent(out)  g,
real(real64), intent(out)  gsg,
real(real64), intent(out)  x 
)
private

Definition at line 1238 of file atomic.F90.

Variable Documentation

◆ spec_notation

character(len=1), dimension(0:3), parameter atomic_oct_m::spec_notation = (/ 's', 'p', 'd', 'f' /)
private

Following stuff is for the "val_conf_t" data type, made to handle atomic configurations.

Definition at line 141 of file atomic.F90.

◆ valconf_string_length

integer, parameter, public atomic_oct_m::valconf_string_length = 80

Definition at line 144 of file atomic.F90.