compute_zpe Subroutine

public subroutine compute_zpe(frequencies, n_freqs, n_real, zpe_hartree, zpe_kcalmol)

Compute zero-point vibrational energy from frequencies.

ZPE = (1/2) * h * sum(nu_i) for all real frequencies. Imaginary frequencies (negative values) are skipped with a warning.

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: frequencies(:)

Vibrational frequencies in cm^-1

integer, intent(in) :: n_freqs

Total number of frequencies

integer, intent(out) :: n_real

Number of real (positive) frequencies used

real(kind=dp), intent(out) :: zpe_hartree

ZPE in Hartree

real(kind=dp), intent(out) :: zpe_kcalmol

ZPE in kcal/mol


Calls

proc~~compute_zpe~~CallsGraph proc~compute_zpe compute_zpe to_char to_char proc~compute_zpe->to_char warning warning proc~compute_zpe->warning

Called by

proc~~compute_zpe~~CalledByGraph proc~compute_zpe compute_zpe proc~compute_thermochemistry compute_thermochemistry proc~compute_thermochemistry->proc~compute_zpe proc~compute_mbe compute_mbe proc~compute_mbe->proc~compute_thermochemistry proc~print_vibrational_analysis print_vibrational_analysis proc~compute_mbe->proc~print_vibrational_analysis proc~gmbe_pie_coordinator gmbe_pie_coordinator proc~gmbe_pie_coordinator->proc~compute_thermochemistry proc~gmbe_pie_coordinator->proc~print_vibrational_analysis proc~hessian_coordinator hessian_coordinator proc~hessian_coordinator->proc~compute_thermochemistry proc~hessian_coordinator->proc~print_vibrational_analysis proc~print_vibrational_analysis->proc~compute_thermochemistry proc~serial_gmbe_pie_processor serial_gmbe_pie_processor proc~serial_gmbe_pie_processor->proc~compute_thermochemistry proc~serial_gmbe_pie_processor->proc~print_vibrational_analysis proc~unfragmented_calculation unfragmented_calculation proc~unfragmented_calculation->proc~compute_thermochemistry proc~unfragmented_calculation->proc~print_vibrational_analysis interface~hessian_coordinator hessian_coordinator interface~hessian_coordinator->proc~hessian_coordinator interface~unfragmented_calculation unfragmented_calculation interface~unfragmented_calculation->proc~unfragmented_calculation proc~compute_gmbe compute_gmbe proc~compute_gmbe->proc~print_vibrational_analysis proc~global_coordinator global_coordinator proc~global_coordinator->proc~compute_mbe proc~gmbe_run_distributed gmbe_context_t%gmbe_run_distributed proc~gmbe_run_distributed->proc~gmbe_pie_coordinator proc~gmbe_run_serial gmbe_context_t%gmbe_run_serial proc~gmbe_run_serial->proc~serial_gmbe_pie_processor proc~serial_fragment_processor serial_fragment_processor proc~serial_fragment_processor->proc~compute_mbe interface~global_coordinator global_coordinator interface~global_coordinator->proc~global_coordinator interface~serial_fragment_processor serial_fragment_processor interface~serial_fragment_processor->proc~serial_fragment_processor proc~distributed_unfragmented_hessian distributed_unfragmented_hessian proc~distributed_unfragmented_hessian->interface~hessian_coordinator proc~run_unfragmented_calculation run_unfragmented_calculation proc~run_unfragmented_calculation->interface~unfragmented_calculation interface~distributed_unfragmented_hessian distributed_unfragmented_hessian interface~distributed_unfragmented_hessian->proc~distributed_unfragmented_hessian proc~mbe_run_distributed mbe_context_t%mbe_run_distributed proc~mbe_run_distributed->interface~global_coordinator proc~mbe_run_serial mbe_context_t%mbe_run_serial proc~mbe_run_serial->interface~serial_fragment_processor proc~run_calculation run_calculation proc~run_calculation->proc~run_unfragmented_calculation

Variables

Type Visibility Attributes Name Initial
real(kind=dp), private :: freq_sum
integer, private :: i
integer, private :: n_imag

Source Code

   subroutine compute_zpe(frequencies, n_freqs, n_real, zpe_hartree, zpe_kcalmol)
      !! Compute zero-point vibrational energy from frequencies.
      !!
      !! ZPE = (1/2) * h * sum(nu_i) for all real frequencies.
      !! Imaginary frequencies (negative values) are skipped with a warning.
      real(dp), intent(in) :: frequencies(:)  !! Vibrational frequencies in cm^-1
      integer, intent(in) :: n_freqs          !! Total number of frequencies
      integer, intent(out) :: n_real          !! Number of real (positive) frequencies used
      real(dp), intent(out) :: zpe_hartree    !! ZPE in Hartree
      real(dp), intent(out) :: zpe_kcalmol    !! ZPE in kcal/mol

      integer :: i
      real(dp) :: freq_sum
      integer :: n_imag

      freq_sum = 0.0_dp
      n_real = 0
      n_imag = 0

      do i = 1, n_freqs
         if (frequencies(i) > IMAG_FREQ_THRESHOLD) then
            freq_sum = freq_sum + frequencies(i)
            n_real = n_real + 1
         else if (frequencies(i) < IMAG_FREQ_THRESHOLD) then
            n_imag = n_imag + 1
         end if
         ! Frequencies exactly at threshold (typically trans/rot modes ~0) are skipped
      end do

      if (n_imag > 0) then
         call logger%warning("Thermochemistry: "//trim(adjustl(to_char(n_imag)))// &
                             " imaginary frequency(ies) skipped")
      end if

      ! ZPE = 0.5 * sum(h * c * nu) where nu is in cm^-1
      ! In atomic units: ZPE = 0.5 * sum(nu * CM1_TO_KELVIN * KB_HARTREE)
      ! Actually: h*c*nu [cm^-1] = h*c*nu [J] = nu * CM1_TO_KELVIN * k_B [J]
      ! So ZPE [Hartree] = 0.5 * sum(nu) * CM1_TO_KELVIN * KB_HARTREE
      zpe_hartree = 0.5_dp*freq_sum*CM1_TO_KELVIN*KB_HARTREE
      zpe_kcalmol = zpe_hartree*HARTREE_TO_KCALMOL

   end subroutine compute_zpe