dispersion_model Derived Type

type, public, abstract :: dispersion_model

Abstract base dispersion model to evaluate C6 coefficients


Components

Type Visibility Attributes Name Initial
real(kind=wp), public, allocatable :: aiw(:,:,:)

Reference dynamic polarizibilities

real(kind=wp), public, allocatable :: c6(:,:,:,:)

Reference C6 coefficients

real(kind=wp), public, allocatable :: cn(:,:)

Reference coordination numbers

real(kind=wp), public, allocatable :: en(:)

Electronegativity

real(kind=wp), public, allocatable :: eta(:)

Chemical hardness

real(kind=wp), public :: ga

Charge scaling height

real(kind=wp), public :: gc

Charge scaling steepness

integer, public :: ncoup

Number of atoms coupled to by pairwise parameters

integer, public, allocatable :: ngw(:,:)

Number of Gaussian weights for each reference

real(kind=wp), public, allocatable :: q(:,:)

Reference partial charges

real(kind=wp), public, allocatable :: r4r2(:)

Expectation values for C8 extrapolation

real(kind=wp), public, allocatable :: rcov(:)

Covalent radii for coordination number

integer, public, allocatable :: ref(:)

Number of reference systems

real(kind=wp), public, allocatable :: zeff(:)

Effective nuclear charges


Type-Bound Procedures

procedure(get_atomic_c6), public, deferred :: get_atomic_c6

Evaluate C6 coefficient

  • subroutine get_atomic_c6(self, mol, gwvec, gwdcn, gwdq, c6, dc6dcn, dc6dq) Prototype

    Calculate atomic dispersion coefficients and their derivatives w.r.t. the coordination numbers and atomic partial charges.

    Arguments

    Type IntentOptional Attributes Name
    class(dispersion_model), intent(in) :: self

    Instance of the dispersion model

    class(structure_type), intent(in) :: mol

    Molecular structure data

    real(kind=wp), intent(in) :: gwvec(:,:,:)

    Weighting function for the atomic reference systems

    real(kind=wp), intent(in), optional :: gwdcn(:,:,:)

    Derivative of the weighting function w.r.t. the coordination number

    real(kind=wp), intent(in), optional :: gwdq(:,:,:)

    Derivative of the weighting function w.r.t. the partial charge

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

    C6 coefficients for all atom pairs.

    real(kind=wp), intent(out), optional :: dc6dcn(:,:)

    Derivative of the C6 w.r.t. the coordination number

    real(kind=wp), intent(out), optional :: dc6dq(:,:)

    Derivative of the C6 w.r.t. the partial charge

procedure(get_polarizibilities), public, deferred :: get_polarizibilities

Evaluate atomic polarizibilities

  • subroutine get_polarizibilities(self, mol, gwvec, gwdcn, gwdq, alpha, dadcn, dadq) Prototype

    Calculate atomic polarizibilities and their derivatives w.r.t. the coordination numbers and atomic partial charges.

    Arguments

    Type IntentOptional Attributes Name
    class(dispersion_model), intent(in) :: self

    Instance of the dispersion model

    class(structure_type), intent(in) :: mol

    Molecular structure data

    real(kind=wp), intent(in) :: gwvec(:,:,:)

    Weighting function for the atomic reference systems

    real(kind=wp), intent(in), optional :: gwdcn(:,:,:)

    Derivative of the weighting function w.r.t. the coordination number

    real(kind=wp), intent(in), optional :: gwdq(:,:,:)

    Derivative of the weighting function w.r.t. the partial charge

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

    Static polarizibilities for all atoms.

    real(kind=wp), intent(out), optional :: dadcn(:)

    Derivative of the polarizibility w.r.t. the coordination number

    real(kind=wp), intent(out), optional :: dadq(:)

    Derivative of the polarizibility w.r.t. the partial charge

procedure(weight_references), public, deferred :: weight_references

Generate weights for all reference systems

  • subroutine weight_references(self, mol, cn, q, gwvec, gwdcn, gwdq) Prototype

    Calculate the weights of the reference system and the derivatives w.r.t. coordination number for later use.

    Arguments

    Type IntentOptional Attributes Name
    class(dispersion_model), intent(in) :: self

    Instance of the dispersion model

    class(structure_type), intent(in) :: mol

    Molecular structure data

    real(kind=wp), intent(in) :: cn(:)

    Coordination number of every atom: [nat]

    real(kind=wp), intent(in) :: q(:)

    Partial charge of every atom: [nat]

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

    weighting for the atomic reference systems: [nref, nat, ncoup]

    real(kind=wp), intent(out), optional :: gwdcn(:,:,:)

    derivative of the weighting function w.r.t. the coordination number: [nref, nat, ncoup]

    real(kind=wp), intent(out), optional :: gwdq(:,:,:)

    derivative of the weighting function w.r.t. the charge scaling: [nref, nat, ncoup]

Source Code

   type, abstract :: dispersion_model

      !> Number of atoms coupled to by pairwise parameters
      integer :: ncoup

      !> Charge scaling height
      real(wp) :: ga

      !> Charge scaling steepness
      real(wp) :: gc

      !> Effective nuclear charges
      real(wp), allocatable :: zeff(:)

      !> Chemical hardness
      real(wp), allocatable :: eta(:)

      !> Electronegativity
      real(wp), allocatable :: en(:)

      !> Covalent radii for coordination number
      real(wp), allocatable :: rcov(:)

      !> Expectation values for C8 extrapolation
      real(wp), allocatable :: r4r2(:)

      !> Number of reference systems
      integer, allocatable :: ref(:)

      !> Number of Gaussian weights for each reference
      integer, allocatable :: ngw(:, :)

      !> Reference coordination numbers
      real(wp), allocatable :: cn(:, :)

      !> Reference partial charges
      real(wp), allocatable :: q(:, :)

      !> Reference dynamic polarizibilities
      real(wp), allocatable :: aiw(:, :, :)

      !> Reference C6 coefficients
      real(wp), allocatable :: c6(:, :, :, :)

   contains

      !> Generate weights for all reference systems
      procedure(weight_references), deferred :: weight_references

      !> Evaluate C6 coefficient
      procedure(get_atomic_c6), deferred :: get_atomic_c6

      !> Evaluate atomic polarizibilities
      procedure(get_polarizibilities), deferred :: get_polarizibilities

   end type dispersion_model