xc_engine_t Derived Type

type, public :: xc_engine_t

@brief Main class which knows how to compute XC functional values, AO and MO values and gradients on a grid @details It is complemented with xc_consumer_t class to use calculation results


Inherits

type~~xc_engine_t~~InheritsGraph type~xc_engine_t xc_engine_t type~xc_libxc_t xc_libxc_t type~xc_engine_t->type~xc_libxc_t XCLib type~xc_lib_t xc_lib_t type~xc_libxc_t->type~xc_lib_t type~xc_pack_t xc_pack_t type~xc_lib_t->type~xc_pack_t ids

Components

Type Visibility Attributes Name Initial
real(kind=fp), public, allocatable :: xyzw(:,:)
real(kind=fp), public, allocatable :: aoMem_(:)
real(kind=fp), public, allocatable :: moMemA_(:)
real(kind=fp), public, allocatable :: moMemB_(:)
real(kind=fp), public, allocatable :: tmpWfAlpha(:)
real(kind=fp), public, allocatable :: tmpWfBeta(:)
integer, public, allocatable :: indices_p(:)
real(kind=fp), public, contiguous, pointer :: aoMem(:,:,:) => null()
real(kind=fp), public, contiguous, pointer :: moMemA(:,:,:) => null()
real(kind=fp), public, contiguous, pointer :: moMemB(:,:,:) => null()
real(kind=fp), public, contiguous, pointer :: aoV(:,:) => null()
real(kind=fp), public, contiguous, pointer :: moVA(:,:) => null()
real(kind=fp), public, contiguous, pointer :: moVB(:,:) => null()
real(kind=fp), public, contiguous, pointer :: aoG1(:,:,:) => null()
real(kind=fp), public, contiguous, pointer :: aoG2(:,:,:) => null()
real(kind=fp), public, contiguous, pointer :: moG1A(:,:,:) => null()
real(kind=fp), public, contiguous, pointer :: moG2A(:,:,:) => null()
real(kind=fp), public, contiguous, pointer :: moG1B(:,:,:) => null()
real(kind=fp), public, contiguous, pointer :: moG2B(:,:,:) => null()
real(kind=fp), public, contiguous, pointer :: wts(:) => null()
real(kind=fp), public, contiguous, pointer :: wfAlpha(:,:) => null()
real(kind=fp), public, contiguous, pointer :: wfBeta(:,:) => null()
real(kind=fp), public, contiguous, pointer :: wfAlpha_p(:,:) => null()
real(kind=fp), public, contiguous, pointer :: wfBeta_p(:,:) => null()
logical, public :: isGGA = .false.
logical, public :: needTau = .false.
logical, public :: hasBeta = .false.
logical, public :: isWFVecs = .true.
integer, public :: numAOs = 0
integer, public :: numAOs_p = 0
logical, public :: skip_p = .true.
integer, public :: numPts = 0
integer, public :: numAtoms = 0
integer, public :: maxPts = 0
integer, public :: maxAngMom = 0
integer, public :: nAODer = 0
integer, public :: nXCDer = 1
integer, public :: funTyp = 0
integer, public :: numAOVecs = 0
integer, public :: numTmpVec = 0
integer, public :: numOccAlpha = 0
integer, public :: numOccBeta = 0
real(kind=fp), public :: threshold = 1.0d-15
real(kind=fp), public :: ao_threshold = 1.0d-15
real(kind=fp), public :: ao_sparsity_ratio = 0.90d+0
type(xc_libxc_t), public, allocatable :: XCLib
integer, public :: dbgLevel = 0
real(kind=fp), public :: N_elec = 0.0
real(kind=fp), public :: E_kin = 0.0
real(kind=fp), public :: G_total(3) = 0.0
procedure(compute_density), public, pointer, pass :: compRho => null()
procedure(compute_density_grad), public, pointer, pass :: compDRho => null()
procedure(compute_density_tau), public, pointer, pass :: compTau => null()

Type-Bound Procedures

procedure, public :: init

  • private subroutine init(self, xco)

    @brief Initialize xc_engine_t instance @param[in] numAOs number of atomic orbitals in a basis @param[in] nAt number of atoms in a system @param[in] maxAngMom maximum angular momentum of basis functions and their derivatives @param[in] maxPts maximum known number of non-zero points in a slice @param[in] limPts maximum possible number of points (i.e. max(nRad*nAng)) in a slice @param[in] nDer degree of energy derivative needed @param[in] hasBeta .TRUE. if open-shell calculation @param[in] isGGA .TRUE. if GGA/metaGGA functional @param[in] needTau .TRUE. if metaGGA functional @param[in] vec_or_dens .TRUE./.FALSE. - wavefunction is MO vectors/density @param[in] nOccAlpha number of occupied orbitals, alpha spin @param[in] nOccBeta number of occupied orbitals, beta spin @param[in] wfAlpha wavefunction, alpha spin @param[in] wfBeta wavefunction, beta spin @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t), intent(inout), target :: self
    type(xc_options_t), intent(in), target :: xco

procedure, public :: echo => echoVars

  • private subroutine echoVars(self)

    @brief Print parameters of the xc_engine_t instance @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: self

procedure, public :: getStats

  • private subroutine getStats(self, E_xc, E_exch, E_corr, N_elec, E_kin, G_total)

    @brief Get debug statistics @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: self
    real(kind=fp), intent(out), optional :: E_xc
    real(kind=fp), intent(out), optional :: E_exch
    real(kind=fp), intent(out), optional :: E_corr
    real(kind=fp), intent(out), optional :: N_elec
    real(kind=fp), intent(out), optional :: E_kin
    real(kind=fp), intent(out), optional :: G_total(3)

procedure, public :: resetPointers

  • private subroutine resetPointers(self, numPts)

    @brief Adjust internal memory storage for a given number of grid points @param[in] numPts number of grid points @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: self
    integer, intent(in) :: numPts

procedure, public :: resetOrbPointers

  • private subroutine resetOrbPointers(self)

    @brief Adjust internal AO/MO memory storage for a given number of grid points @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t), target :: self

procedure, public :: resetXCPointers

  • private subroutine resetXCPointers(self)

    @brief Adjust XC memory storage for a given number of grid points @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t), target :: self

procedure, public :: compAOs

  • private subroutine compAOs(self, basis, nDer, xyz)

    @brief Compute atomic orbital values/gradient/hessian in a grid point @param[in] iPtIn index of the point in self%xyzw array @param[in] iPtOut index of the point in AO/MO arrays @param[out] nnz number of non-zero AOs in the point @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: self
    type(basis_set), intent(in) :: basis
    integer, intent(in) :: nDer
    real(kind=fp), intent(in) :: xyz(:,:)

procedure, public :: pruneAOs

  • private subroutine pruneAOs(self, skip)

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t), target :: self
    logical :: skip

procedure, public :: resetPrunedPointers

  • private subroutine resetPrunedPointers(self)

    @brief Adjust XC memory storage for a given number of pruned grid points @author Konstantin Komarov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t), target :: self

procedure, public :: compMOs

  • private subroutine compMOs(self)

    @brief Transform AOs to "MOs" @details Multiply AO vector to the MO coefficient matrix or density matrix. True MOs are only obtained in the former case. @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: self

procedure, public :: compRMOs

  • private subroutine compRMOs(xce, da, mo)

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: xce
    real(kind=fp) :: da(:,:,:)
    real(kind=fp) :: mo(:,:,:)

procedure, public :: compRMOGs

  • private subroutine compRMOGs(xce, da, moG1)

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: xce
    real(kind=fp) :: da(:,:,:)
    real(kind=fp) :: moG1(:,:,:,:)

generic, public :: compRRho => compRRho_ab, compRRho_a

  • private subroutine compRRho_ab(xce, mo, rho)

    @brief Compute electronic density in a grid point, density-driven calculation @param[in] xce XC engine, parameters @param[in] mo "Molecular orbitals" @param[out] rho electronic density @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: xce
    real(kind=fp), intent(in) :: mo(:,:,:,:)
    real(kind=fp), intent(out) :: rho(:,:,:)
  • private subroutine compRRho_a(xce, mo, rho)

    @brief Compute electronic density in a grid point, density-driven calculation @param[in] xce XC engine, parameters @param[in] mo "Molecular orbitals" @param[out] rho electronic density @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: xce
    real(kind=fp), intent(in) :: mo(:,:,:)
    real(kind=fp), intent(out) :: rho(:,:)

generic, public :: compRDRho => compRDRho_ab, compRDRho_a

  • private subroutine compRDRho_ab(xce, mo, drrho)

    @brief Compute electronic density gradient in a grid point, density-driven calculation @param[in] xce XC engine, parameters @param[in] mo "Molecular orbitals" @param[out] drho density directional derivative (along X, Y, and Z axes) @param[out] drrho dRho/d[x,y,z] vector and its dot product with drho @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: xce
    real(kind=fp), intent(in) :: mo(:,:,:,:)
    real(kind=fp), intent(out) :: drrho(:,:,:,:)
  • private subroutine compRDRho_a(xce, mo, drrho)

    @brief Compute electronic density gradient in a grid point, density-driven calculation @param[in] xce XC engine, parameters @param[in] mo "Molecular orbitals" @param[out] drho density directional derivative (along X, Y, and Z axes) @param[out] drrho dRho/d[x,y,z] vector and its dot product with drho @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: xce
    real(kind=fp), intent(in) :: mo(:,:,:)
    real(kind=fp), intent(out) :: drrho(:,:,:)

generic, public :: compRTau => compRTau_ab, compRTau_a

  • private subroutine compRTau_ab(xce, moG1, rtau)

    @brief Compute tau: (MO)' times (AO)' @param[in] xce XC engine, parameters @param[in] mog1 MO directional derivatives @param[out] rtau kinetic energy density @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: xce
    real(kind=fp), intent(in) :: moG1(:,:,:,:,:)
    real(kind=fp), intent(out) :: rtau(:,:,:)
  • private subroutine compRTau_a(xce, moG1, tau)

    @brief Compute tau: (MO)' times (AO)' @param[in] xce XC engine, parameters @param[in] mog1 MO directional derivatives @param[out] tau kinetic energy density @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: xce
    real(kind=fp), intent(in) :: moG1(:,:,:,:)
    real(kind=fp), intent(out) :: tau(:,:)

procedure, public :: compRhoAll

  • private subroutine compRhoAll(self, skip)

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: self
    logical, intent(out) :: skip

procedure, public :: compXC

  • private subroutine compXC(self, functional, skip)

    @brief Compute XC contribution to the energy @param[inout] bfGrad array of gradient contributinos per basis function @param[inout] exec XC energy integral @param[inout] ecorl correlation energy integral @param[inout] totele density integral == number of electrons @param[inout] totkin kinetic energy integral @param[inout] togradxyz density gradient integral @author Vladimir Mironov

    Arguments

    Type IntentOptional Attributes Name
    class(xc_engine_t) :: self
    type(functional_t) :: functional
    logical :: skip

Source Code

  type :: xc_engine_t
!    private
    real(KIND=fp), allocatable :: xyzw(:, :)
    real(KIND=fp), allocatable :: aoMem_(:) !< AO memory
    real(KIND=fp), allocatable :: moMemA_(:) !< MO memory (alpha)
    real(KIND=fp), allocatable :: moMemB_(:) !< MO memory (beta)
    real(KIND=fp), allocatable :: tmpWfAlpha(:) !< tmp alpha spin wavefunction
    real(KIND=fp), allocatable :: tmpWfBeta(:) !< tmp beta spin wavefunction
    integer, allocatable :: indices_p(:) !< AO significant indices

    real(KIND=fp), contiguous, pointer :: &
        aoMem(:, :, :) => null() & !< AO memory
      , moMemA(:, :, :) => null() & !< MO memory (alpha)
      , moMemB(:, :, :) => null() & !< MO memory (beta)
      , aoV(:, :) => null() & !< AO values
      , moVA(:, :) => null() & !< MO values (alpha)
      , moVB(:, :) => null() & !< MO values (beta)
      , aoG1(:, :, :) => null() & !< AO gradient
      , aoG2(:, :, :) => null() & !< AO 2nd der.
      , moG1A(:, :, :) => null() & !< MO gradient (alpha)
      , moG2A(:, :, :) => null() & !< MO 2nd der. (alpha)
      , moG1B(:, :, :) => null() & !< MO gradient (beta)
      , moG2B(:, :, :) => null() & !< MO 2nd der. (beta)
      , wts(:) => null() & !< weights
      , wfAlpha(:, :) => null() & !< alpha spin wavefunction
      , wfBeta(:, :) => null() & !< beta spin wavefunction
      , wfAlpha_p(:, :) => null() & !< pruned alpha spin wavefunction
      , wfBeta_p(:, :) => null()   !< pruned beta spin wavefunction

    logical :: isGGA = .false.
    logical :: needTau = .false.
    logical :: hasBeta = .false.
    logical :: isWFVecs = .true.   !< .TRUE.  - wfA and wfB are MO vectors
                                   !< .FALSE. - wfA and wfB are densities
    integer :: numAOs = 0 !< number of AOs
    integer :: numAOs_p = 0 !< number of pruned AOs
    logical :: skip_p = .true. !< skip if no pruned numAOs
    integer :: numPts = 0
    integer :: numAtoms = 0
    integer :: maxPts = 0
    integer :: maxAngMom = 0
    integer :: nAODer = 0
    integer :: nXCDer = 1
    integer :: funTyp = 0   !< 0 - LDA, 1 - GGA, 2 - MGGA

    integer :: numAOVecs = 0
    integer :: numTmpVec = 0

    integer :: numOccAlpha = 0
    integer :: numOccBeta = 0

    real(kind=fp) :: threshold  = 1.0d-15
    real(kind=fp) :: ao_threshold  = 1.0d-15
    real(kind=fp) :: ao_sparsity_ratio = 0.90d+0 !< Cut off if more than 90% pruned AOs
                                                 !< (skip_p becomes False).

    type(xc_libxc_t), allocatable :: XCLib

    integer       :: dbgLevel   = 0
    real(kind=fp) :: N_elec     = 0.0
    real(kind=fp) :: E_kin      = 0.0
    real(kind=fp) :: G_total(3) = 0.0

    procedure(compute_density), pointer, pass :: compRho => null()
    procedure(compute_density_grad), pointer, pass :: compDRho => null()
    procedure(compute_density_tau), pointer, pass :: compTau => null()

  contains
    procedure :: init
    procedure :: echo => echoVars
    procedure :: getStats
    procedure :: resetPointers
    procedure :: resetOrbPointers
    procedure :: resetXCPointers
    procedure :: compAOs
    procedure :: pruneAOs
    procedure :: resetPrunedPointers
    procedure :: compMOs
    procedure :: compRMOs
    procedure :: compRMOGs

    generic :: compRRho  => compRRho_ab, compRRho_a
    generic :: compRDRho => compRDRho_ab, compRDRho_a
    generic :: compRTau => compRTau_ab, compRTau_a

    procedure :: compRhoAll

    procedure :: compXC

    procedure, private :: compRRho_ab
    procedure, private :: compRDRho_ab
    procedure, private :: compRTau_ab
    procedure, private :: compRRho_a
    procedure, private :: compRDRho_a
    procedure, private :: compRTau_a

  end type xc_engine_t