mqc_cgto.f90 Source File

Data structures for cartesian contracted Gaussian type orbitals (CGTOs)


Files dependent on this one

sourcefile~~mqc_cgto.f90~~AfferentGraph sourcefile~mqc_cgto.f90 mqc_cgto.f90 sourcefile~mqc_basis_reader.f90 mqc_basis_reader.f90 sourcefile~mqc_basis_reader.f90->sourcefile~mqc_cgto.f90 sourcefile~mqc_physical_fragment.f90 mqc_physical_fragment.f90 sourcefile~mqc_physical_fragment.f90->sourcefile~mqc_cgto.f90 sourcefile~main.f90 main.f90 sourcefile~main.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_config_adapter.f90 mqc_config_adapter.f90 sourcefile~main.f90->sourcefile~mqc_config_adapter.f90 sourcefile~mqc_config_parser.f90 mqc_config_parser.F90 sourcefile~main.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_driver.f90 mqc_driver.f90 sourcefile~main.f90->sourcefile~mqc_driver.f90 sourcefile~mqc_calculation_interface.f90 mqc_calculation_interface.f90 sourcefile~mqc_calculation_interface.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_calculation_interface.f90->sourcefile~mqc_config_adapter.f90 sourcefile~mqc_calculation_interface.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_calculation_interface.f90->sourcefile~mqc_driver.f90 sourcefile~mqc_combinatorics.f90 mqc_combinatorics.f90 sourcefile~mqc_combinatorics.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_config_adapter.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_config_adapter.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_config_parser.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_config_parser_fragments.f90 mqc_config_parser_fragments.f90 sourcefile~mqc_config_parser_fragments.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_config_parser_fragments.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_driver.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_driver.f90->sourcefile~mqc_config_adapter.f90 sourcefile~mqc_driver.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_frag_utils.f90 mqc_frag_utils.f90 sourcefile~mqc_driver.f90->sourcefile~mqc_frag_utils.f90 sourcefile~mqc_many_body_expansion.f90 mqc_many_body_expansion.f90 sourcefile~mqc_driver.f90->sourcefile~mqc_many_body_expansion.f90 sourcefile~mqc_mbe.f90 mqc_mbe.f90 sourcefile~mqc_driver.f90->sourcefile~mqc_mbe.f90 sourcefile~mqc_mbe_fragment_distribution_scheme.f90 mqc_mbe_fragment_distribution_scheme.F90 sourcefile~mqc_driver.f90->sourcefile~mqc_mbe_fragment_distribution_scheme.f90 sourcefile~mqc_json_writer.f90 mqc_json_writer.f90 sourcefile~mqc_driver.f90->sourcefile~mqc_json_writer.f90 sourcefile~mqc_finite_differences.f90 mqc_finite_differences.f90 sourcefile~mqc_finite_differences.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_frag_utils.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_frag_utils.f90->sourcefile~mqc_combinatorics.f90 sourcefile~mqc_frag_utils.f90->sourcefile~mqc_config_adapter.f90 sourcefile~mqc_gmbe_utils.f90 mqc_gmbe_utils.f90 sourcefile~mqc_frag_utils.f90->sourcefile~mqc_gmbe_utils.f90 sourcefile~mqc_gmbe_fragment_distribution_scheme.f90 mqc_gmbe_fragment_distribution_scheme.f90 sourcefile~mqc_gmbe_fragment_distribution_scheme.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_gmbe_fragment_distribution_scheme.f90->sourcefile~mqc_mbe_fragment_distribution_scheme.f90 sourcefile~mqc_gmbe_utils.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_gmbe_utils.f90->sourcefile~mqc_combinatorics.f90 sourcefile~mqc_many_body_expansion.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_many_body_expansion.f90->sourcefile~mqc_config_adapter.f90 sourcefile~mqc_many_body_expansion.f90->sourcefile~mqc_gmbe_fragment_distribution_scheme.f90 sourcefile~mqc_many_body_expansion.f90->sourcefile~mqc_mbe_fragment_distribution_scheme.f90 sourcefile~mqc_mbe.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_mbe.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_mbe.f90->sourcefile~mqc_frag_utils.f90 sourcefile~mqc_mbe.f90->sourcefile~mqc_gmbe_utils.f90 sourcefile~mqc_mbe_io.f90 mqc_mbe_io.f90 sourcefile~mqc_mbe.f90->sourcefile~mqc_mbe_io.f90 sourcefile~mqc_mbe_fragment_distribution_scheme.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_mbe_fragment_distribution_scheme.f90->sourcefile~mqc_config_adapter.f90 sourcefile~mqc_mbe_fragment_distribution_scheme.f90->sourcefile~mqc_mbe.f90 sourcefile~mqc_mbe_fragment_distribution_scheme.f90->sourcefile~mqc_mbe_io.f90 sourcefile~mqc_method_base.f90 mqc_method_base.f90 sourcefile~mqc_mbe_fragment_distribution_scheme.f90->sourcefile~mqc_method_base.f90 sourcefile~mqc_method_factory.f90 mqc_method_factory.F90 sourcefile~mqc_mbe_fragment_distribution_scheme.f90->sourcefile~mqc_method_factory.f90 sourcefile~mqc_mbe_io.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_method_base.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_method_dft.f90 mqc_method_dft.f90 sourcefile~mqc_method_dft.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_method_dft.f90->sourcefile~mqc_method_base.f90 sourcefile~mqc_method_hf.f90 mqc_method_hf.f90 sourcefile~mqc_method_hf.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_method_hf.f90->sourcefile~mqc_method_base.f90 sourcefile~mqc_method_mcscf.f90 mqc_method_mcscf.f90 sourcefile~mqc_method_mcscf.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_method_mcscf.f90->sourcefile~mqc_method_base.f90 sourcefile~mqc_method_xtb.f90 mqc_method_xtb.f90 sourcefile~mqc_method_xtb.f90->sourcefile~mqc_physical_fragment.f90 sourcefile~mqc_method_xtb.f90->sourcefile~mqc_finite_differences.f90 sourcefile~mqc_method_xtb.f90->sourcefile~mqc_method_base.f90 sourcefile~mqc_config_parser_basic_sections.f90 mqc_config_parser_basic_sections.f90 sourcefile~mqc_config_parser_basic_sections.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_config_parser_calc_settings.f90 mqc_config_parser_calc_settings.f90 sourcefile~mqc_config_parser_calc_settings.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_config_parser_molecules.f90 mqc_config_parser_molecules.f90 sourcefile~mqc_config_parser_molecules.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_config_parser_structure.f90 mqc_config_parser_structure.f90 sourcefile~mqc_config_parser_structure.f90->sourcefile~mqc_config_parser.f90 sourcefile~mqc_json_writer.f90->sourcefile~mqc_mbe_io.f90 sourcefile~mqc_mbe_fragment_distribution_scheme_hessian.f90 mqc_mbe_fragment_distribution_scheme_hessian.F90 sourcefile~mqc_mbe_fragment_distribution_scheme_hessian.f90->sourcefile~mqc_finite_differences.f90 sourcefile~mqc_mbe_fragment_distribution_scheme_hessian.f90->sourcefile~mqc_mbe_fragment_distribution_scheme.f90 sourcefile~mqc_mbe_fragment_distribution_scheme_hessian.f90->sourcefile~mqc_method_base.f90 sourcefile~mqc_mbe_fragment_distribution_scheme_hessian.f90->sourcefile~mqc_method_factory.f90 sourcefile~mqc_mbe_mpi_fragment_distribution_scheme.f90 mqc_mbe_mpi_fragment_distribution_scheme.F90 sourcefile~mqc_mbe_mpi_fragment_distribution_scheme.f90->sourcefile~mqc_mbe_fragment_distribution_scheme.f90 sourcefile~mqc_method_factory.f90->sourcefile~mqc_method_base.f90 sourcefile~mqc_method_factory.f90->sourcefile~mqc_method_dft.f90 sourcefile~mqc_method_factory.f90->sourcefile~mqc_method_hf.f90 sourcefile~mqc_method_factory.f90->sourcefile~mqc_method_mcscf.f90 sourcefile~mqc_method_factory.f90->sourcefile~mqc_method_xtb.f90 sourcefile~mqc_serial_fragment_processor.f90 mqc_serial_fragment_processor.f90 sourcefile~mqc_serial_fragment_processor.f90->sourcefile~mqc_mbe_fragment_distribution_scheme.f90 sourcefile~mqc_unfragmented_workflow.f90 mqc_unfragmented_workflow.f90 sourcefile~mqc_unfragmented_workflow.f90->sourcefile~mqc_mbe_fragment_distribution_scheme.f90

Source Code

!! Data structures for cartesian contracted Gaussian type orbitals (CGTOs)
module mqc_cgto
   !! Defines data structures for cartesian contracted Gaussian type orbitals (CGTOs)
   use pic_types, only: dp
   implicit none
   private

   public :: cgto_type, atomic_basis_type, molecular_basis_type

   type :: cgto_type
      !! Contracted Gaussian type orbital (CGTO) data structure
      integer :: ang_mom
        !! Angular momentum quantum number (0=s, 1=p, 2=d, etc.)
      integer :: nfunc
        !! Number of primitive Gaussians in the contraction
      real(dp), allocatable :: exponents(:)
        !! Exponents (alpha values)
      real(dp), allocatable :: coefficients(:)
        !! Contraction coefficients
   contains
      procedure :: allocate_arrays => cgto_allocate_arrays
      procedure :: destroy => cgto_destroy
      procedure :: num_basis_functions => cgto_num_basis_functions
   end type cgto_type

   type :: atomic_basis_type
      !! Atomic basis set data structure
      character(len=:), allocatable :: element
      !! element symbol
      type(cgto_type), allocatable :: shells(:)
      !! array of contracted shells
      integer :: nshells
      !! number of shells in type
   contains
      procedure :: allocate_shells => allocate_basis_shells
      procedure :: destroy => atomic_basis_destroy
      procedure :: num_basis_functions => atomic_basis_num_basis_functions
   end type atomic_basis_type

   type :: molecular_basis_type
      !! Molecular basis set data structure (assembled basis)
      type(atomic_basis_type), allocatable :: elements(:)
      !! array of atomic basis types
      integer :: nelements
      !! total number of atoms/elements in a molecule
   contains
      procedure :: allocate_elements => basis_set_allocate_elements
      procedure :: destroy => basis_set_destroy
      procedure :: num_basis_functions => molecular_basis_num_basis_functions
   end type molecular_basis_type

contains

   pure subroutine cgto_allocate_arrays(self, nfunc)
      !! Allocate arrays for exponents and coefficients in a CGTO
      class(cgto_type), intent(inout) :: self
      integer, intent(in) :: nfunc

      self%nfunc = nfunc
      allocate (self%exponents(nfunc))
      allocate (self%coefficients(nfunc))
   end subroutine cgto_allocate_arrays

   pure subroutine cgto_destroy(self)
      !! Clean up allocated memory in a CGTO
      class(cgto_type), intent(inout) :: self

      if (allocated(self%exponents)) deallocate (self%exponents)
      if (allocated(self%coefficients)) deallocate (self%coefficients)
      self%nfunc = 0
      self%ang_mom = 0
   end subroutine cgto_destroy

   pure subroutine allocate_basis_shells(self, nshells)
      !! Allocate array of shells in an atomic basis
      class(atomic_basis_type), intent(inout) :: self
      integer, intent(in) :: nshells

      self%nshells = nshells
      allocate (self%shells(nshells))
   end subroutine allocate_basis_shells

   pure subroutine atomic_basis_destroy(self)
      !! Clean up allocated memory in an atomic basis
      class(atomic_basis_type), intent(inout) :: self
      integer :: i

      if (allocated(self%shells)) then
         do i = 1, self%nshells
            call self%shells(i)%destroy()
         end do
         deallocate (self%shells)
      end if
      if (allocated(self%element)) deallocate (self%element)
      self%nshells = 0
   end subroutine atomic_basis_destroy

   pure subroutine basis_set_allocate_elements(self, nelements)
      !! Allocate array of atomic basis elements in a molecular basis set
      class(molecular_basis_type), intent(inout) :: self
      integer, intent(in) :: nelements

      self%nelements = nelements
      allocate (self%elements(nelements))

   end subroutine basis_set_allocate_elements

   pure subroutine basis_set_destroy(self)
      !! Clean up allocated memory in a molecular basis set
      class(molecular_basis_type), intent(inout) :: self
      integer :: i

      if (allocated(self%elements)) then
         do i = 1, self%nelements
            call self%elements(i)%destroy()
         end do
         deallocate (self%elements)
      end if

      self%nelements = 0
   end subroutine basis_set_destroy

   pure function cgto_num_basis_functions(self) result(nbf)
      !! Get number of basis functions in a shell (Cartesian)
      class(cgto_type), intent(in) :: self
      integer :: nbf

      ! Cartesian: (ang_mom+1)*(ang_mom+2)/2
      nbf = (self%ang_mom + 1)*(self%ang_mom + 2)/2
   end function cgto_num_basis_functions

   pure function atomic_basis_num_basis_functions(self) result(nbf)
      !! Get total number of basis functions for an atom
      class(atomic_basis_type), intent(in) :: self
      integer :: nbf
      integer :: ishell

      nbf = 0
      do ishell = 1, self%nshells
         nbf = nbf + self%shells(ishell)%num_basis_functions()
      end do
   end function atomic_basis_num_basis_functions

   pure function molecular_basis_num_basis_functions(self) result(nbf)
      !! Get total number of basis functions for the molecule
      class(molecular_basis_type), intent(in) :: self
      integer :: nbf
      integer :: iatom

      nbf = 0
      do iatom = 1, self%nelements
         nbf = nbf + self%elements(iatom)%num_basis_functions()
      end do
   end function molecular_basis_num_basis_functions

end module mqc_cgto