get_pairwise_dispersion Subroutine

public subroutine get_pairwise_dispersion(mol, disp, param, cutoff, energy2, energy3)

Wrapper to handle the evaluation of pairwise representation of the dispersion energy

Arguments

Type IntentOptional Attributes Name
class(structure_type), intent(in) :: mol

Molecular structure data

class(d4_model), intent(in) :: disp

Dispersion model

class(damping_param), intent(in) :: param

Damping parameters

type(realspace_cutoff), intent(in) :: cutoff

Realspace cutoffs

real(kind=wp), intent(out) :: energy2(:,:)

Pairwise representation of additive dispersion energy

real(kind=wp), intent(out) :: energy3(:,:)

Pairwise representation of non-additive dispersion energy


Source Code

subroutine get_pairwise_dispersion(mol, disp, param, cutoff, energy2, energy3)
   !DEC$ ATTRIBUTES DLLEXPORT :: get_pairwise_dispersion

   !> Molecular structure data
   class(structure_type), intent(in) :: mol

   !> Dispersion model
   class(d4_model), intent(in) :: disp

   !> Damping parameters
   class(damping_param), intent(in) :: param

   !> Realspace cutoffs
   type(realspace_cutoff), intent(in) :: cutoff

   !> Pairwise representation of additive dispersion energy
   real(wp), intent(out) :: energy2(:, :)

   !> Pairwise representation of non-additive dispersion energy
   real(wp), intent(out) :: energy3(:, :)

   integer :: mref
   real(wp), allocatable :: cn(:), q(:), gwvec(:, :), c6(:, :), lattr(:, :)

   mref = maxval(disp%ref)

   allocate(cn(mol%nat))
   call get_lattice_points(mol%periodic, mol%lattice, cutoff%cn, lattr)
   call get_coordination_number(mol, lattr, cutoff%cn, disp%rcov, disp%en, cn)

   allocate(q(mol%nat))
   call get_charges(mol, q)

   allocate(gwvec(mref, mol%nat))
   call disp%weight_references(mol, cn, q, gwvec)

   allocate(c6(mol%nat, mol%nat))
   call disp%get_atomic_c6(mol, gwvec, c6=c6)

   energy2(:, :) = 0.0_wp
   energy3(:, :) = 0.0_wp
   call get_lattice_points(mol%periodic, mol%lattice, cutoff%disp2, lattr)
   call param%get_pairwise_dispersion2(mol, lattr, cutoff%disp2, disp%r4r2, &
      & c6, energy2)

   q(:) = 0.0_wp
   call disp%weight_references(mol, cn, q, gwvec)
   call disp%get_atomic_c6(mol, gwvec, c6=c6)

   call get_lattice_points(mol%periodic, mol%lattice, cutoff%disp3, lattr)
   call param%get_pairwise_dispersion3(mol, lattr, cutoff%disp3, disp%r4r2, &
      & c6, energy3)

end subroutine get_pairwise_dispersion