compute_rotational_constants Subroutine

public pure subroutine compute_rotational_constants(moments, is_linear, rot_const)

Convert moments of inertia to rotational constants in GHz.

B = h / (8 * pi^2 * I) where I is in SI units. For I in amu*Angstrom^2: B(GHz) = 505379.07 / I

Arguments

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

Moments in amu*Angstrom^2

logical, intent(in) :: is_linear

True if linear molecule

real(kind=dp), intent(out) :: rot_const(3)

Rotational constants in GHz


Called by

proc~~compute_rotational_constants~~CalledByGraph proc~compute_rotational_constants compute_rotational_constants proc~compute_thermochemistry compute_thermochemistry proc~compute_thermochemistry->proc~compute_rotational_constants 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
integer, private :: i

Source Code

   pure subroutine compute_rotational_constants(moments, is_linear, rot_const)
      !! Convert moments of inertia to rotational constants in GHz.
      !!
      !! B = h / (8 * pi^2 * I) where I is in SI units.
      !! For I in amu*Angstrom^2: B(GHz) = 505379.07 / I
      real(dp), intent(in) :: moments(3)      !! Moments in amu*Angstrom^2
      logical, intent(in) :: is_linear        !! True if linear molecule
      real(dp), intent(out) :: rot_const(3)   !! Rotational constants in GHz

      integer :: i

      rot_const = 0.0_dp

      if (is_linear) then
         ! For linear molecules, only one rotational constant matters
         ! Use the largest moment (moments are sorted ascending)
         if (moments(3) > LINEAR_THRESHOLD) then
            rot_const(1) = ROTCONST_AMUA2_TO_GHZ/moments(3)
         end if
      else
         ! For nonlinear molecules, compute all three
         do i = 1, 3
            if (moments(i) > LINEAR_THRESHOLD) then
               rot_const(i) = ROTCONST_AMUA2_TO_GHZ/moments(i)
            end if
         end do
      end if

   end subroutine compute_rotational_constants