Octopus
xc_fbe_oct_m Module Reference

Functions/Subroutines

subroutine, public x_fbe_calc (id, namespace, psolver, gr, st, space, ex, vxc)
 Interface to X(x_fbe_calc) Two possible run modes possible: adiabatic and Sturm-Liouville. In the first one, we assume no current and solve the local force-balance equation In the second case, we solve the Sturm-Liouville equation The energy is given by the virial relation. More...
 
subroutine solve_sturm_liouville (namespace, gr, st, space, fxc, vxc)
 Solve the Sturm-Liouville equation On entry, vxc is the adiabatic one, on exit, it is the solution of the Sturm-Liouville equation. More...
 
real(real64) function get_virial_energy (gr, nspin, fxc)
 Computes the energy from the force virial relation. More...
 
subroutine, public lda_c_fbe (st, n_blocks, l_dens, l_dedd, l_zk)
 Computes the local density correlation potential and energy obtained from the Colle-Salvetti approximation to the reduced density matrix, with a gradient expansion on the correlation force density. More...
 
subroutine, public fbe_c_lda_sl (namespace, psolver, gr, st, space, ec, vxc)
 Sturm-Liouville version of the FBE local-density correlation functional. More...
 
subroutine dx_fbe_calc (namespace, psolver, mesh, der, st, ex, vxc, fxc)
 
subroutine dfbe_x (namespace, mesh, der, psolver, st, isp, ex, vxc, fxc)
 This routine is adapted from the dslater routine. More...
 
subroutine zx_fbe_calc (namespace, psolver, mesh, der, st, ex, vxc, fxc)
 
subroutine zfbe_x (namespace, mesh, der, psolver, st, isp, ex, vxc, fxc)
 This routine is adapted from the zslater routine. More...
 

Variables

type(grid_t), pointer gr_aux => null()
 
real(real64), dimension(:), pointer rho_aux => null()
 
real(real64), dimension(:), allocatable diag_lapl
 diagonal of the laplacian More...
 

Function/Subroutine Documentation

◆ x_fbe_calc()

subroutine, public xc_fbe_oct_m::x_fbe_calc ( integer, intent(in)  id,
type(namespace_t), intent(in)  namespace,
type(poisson_t), intent(in)  psolver,
type(grid_t), intent(in)  gr,
type(states_elec_t), intent(inout)  st,
type(space_t), intent(in)  space,
real(real64), intent(inout)  ex,
real(real64), dimension(:,:), intent(inout), optional, contiguous  vxc 
)

Interface to X(x_fbe_calc) Two possible run modes possible: adiabatic and Sturm-Liouville. In the first one, we assume no current and solve the local force-balance equation In the second case, we solve the Sturm-Liouville equation The energy is given by the virial relation.

Parameters
[in,out]vxcvxc(grmeshnp, stdnspin)

Definition at line 164 of file xc_fbe.F90.

◆ solve_sturm_liouville()

subroutine xc_fbe_oct_m::solve_sturm_liouville ( type(namespace_t), intent(in)  namespace,
type(grid_t), intent(in), target  gr,
type(states_elec_t), intent(in), target  st,
type(space_t), intent(in)  space,
real(real64), dimension(:,:,:), intent(inout), contiguous  fxc,
real(real64), dimension(:,:), intent(inout), contiguous  vxc 
)
private

Solve the Sturm-Liouville equation On entry, vxc is the adiabatic one, on exit, it is the solution of the Sturm-Liouville equation.

Definition at line 219 of file xc_fbe.F90.

◆ get_virial_energy()

real(real64) function xc_fbe_oct_m::get_virial_energy ( type(grid_t), intent(in)  gr,
integer, intent(in)  nspin,
real(real64), dimension(:,:,:), intent(in)  fxc 
)
private

Computes the energy from the force virial relation.

Definition at line 327 of file xc_fbe.F90.

◆ lda_c_fbe()

subroutine, public xc_fbe_oct_m::lda_c_fbe ( type(states_elec_t), intent(in)  st,
integer, intent(in)  n_blocks,
real(real64), dimension(:,:), intent(in)  l_dens,
real(real64), dimension(:,:), intent(inout)  l_dedd,
real(real64), dimension(:), intent(inout), optional  l_zk 
)

Computes the local density correlation potential and energy obtained from the Colle-Salvetti approximation to the reduced density matrix, with a gradient expansion on the correlation force density.

The spin-polarized case does not have an energy and uses the approximated potential that neglects \nabla \zeta. The energy is therefore wrong for the spin-polarized case

Definition at line 362 of file xc_fbe.F90.

◆ fbe_c_lda_sl()

subroutine, public xc_fbe_oct_m::fbe_c_lda_sl ( type(namespace_t), intent(in)  namespace,
type(poisson_t), intent(in)  psolver,
type(grid_t), intent(in)  gr,
type(states_elec_t), intent(inout)  st,
type(space_t), intent(in)  space,
real(real64), intent(inout)  ec,
real(real64), dimension(:,:), intent(inout), optional, contiguous  vxc 
)

Sturm-Liouville version of the FBE local-density correlation functional.

Parameters
[in,out]vxcvxc(grmeshnp, stdnspin)

Definition at line 430 of file xc_fbe.F90.

◆ dx_fbe_calc()

subroutine xc_fbe_oct_m::dx_fbe_calc ( type(namespace_t), intent(in)  namespace,
type(poisson_t), intent(in)  psolver,
class(mesh_t), intent(in)  mesh,
type(derivatives_t), intent(in)  der,
type(states_elec_t), intent(inout)  st,
real(real64), intent(inout)  ex,
real(real64), dimension(:,:), intent(inout), optional, contiguous  vxc,
real(real64), dimension(:,:,:), intent(out), optional, contiguous  fxc 
)
private
Parameters
[in,out]vxcvxc(grmeshnp, stdnspin)
[out]fxcthe exchange force density

Definition at line 579 of file xc_fbe.F90.

◆ dfbe_x()

subroutine xc_fbe_oct_m::dfbe_x ( type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
type(derivatives_t), intent(in)  der,
type(poisson_t), intent(in)  psolver,
type(states_elec_t), intent(inout), target  st,
integer, intent(in)  isp,
real(real64), intent(inout)  ex,
real(real64), dimension(:), intent(inout), optional, contiguous  vxc,
real(real64), dimension(:,:), intent(inout), optional, contiguous  fxc 
)
private

This routine is adapted from the dslater routine.

Parameters
[in,out]vxcvxc(grmeshnp)
[in,out]fxcthe exchange force density

Definition at line 636 of file xc_fbe.F90.

◆ zx_fbe_calc()

subroutine xc_fbe_oct_m::zx_fbe_calc ( type(namespace_t), intent(in)  namespace,
type(poisson_t), intent(in)  psolver,
class(mesh_t), intent(in)  mesh,
type(derivatives_t), intent(in)  der,
type(states_elec_t), intent(inout)  st,
real(real64), intent(inout)  ex,
real(real64), dimension(:,:), intent(inout), optional, contiguous  vxc,
real(real64), dimension(:,:,:), intent(out), optional, contiguous  fxc 
)
private
Parameters
[in,out]vxcvxc(grmeshnp, stdnspin)
[out]fxcthe exchange force density

Definition at line 941 of file xc_fbe.F90.

◆ zfbe_x()

subroutine xc_fbe_oct_m::zfbe_x ( type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
type(derivatives_t), intent(in)  der,
type(poisson_t), intent(in)  psolver,
type(states_elec_t), intent(inout), target  st,
integer, intent(in)  isp,
real(real64), intent(inout)  ex,
real(real64), dimension(:), intent(inout), optional, contiguous  vxc,
real(real64), dimension(:,:), intent(inout), optional, contiguous  fxc 
)
private

This routine is adapted from the zslater routine.

Parameters
[in,out]vxcvxc(grmeshnp)
[in,out]fxcthe exchange force density

Definition at line 998 of file xc_fbe.F90.

Variable Documentation

◆ gr_aux

type(grid_t), pointer xc_fbe_oct_m::gr_aux => null()
private

Definition at line 151 of file xc_fbe.F90.

◆ rho_aux

real(real64), dimension(:), pointer xc_fbe_oct_m::rho_aux => null()
private

Definition at line 152 of file xc_fbe.F90.

◆ diag_lapl

real(real64), dimension(:), allocatable xc_fbe_oct_m::diag_lapl
private

diagonal of the laplacian

Definition at line 153 of file xc_fbe.F90.