Perform complete vibrational analysis from Hessian matrix.
This is a convenience wrapper that computes: - Vibrational frequencies in cm⁻¹ - Reduced masses in amu - Force constants in Hartree/Bohr² (and optionally mdyne/Å) - Cartesian displacement vectors (normalized) - IR intensities in km/mol (if dipole_derivatives provided)
Optionally projects out translation/rotation modes.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=dp), | intent(in) | :: | hessian(:,:) |
Hessian matrix in Hartree/Bohr² (3N x 3N) |
||
| integer, | intent(in) | :: | element_numbers(:) |
Atomic numbers for each atom (N atoms) |
||
| real(kind=dp), | intent(out), | allocatable | :: | frequencies(:) |
Vibrational frequencies in cm⁻¹ |
|
| real(kind=dp), | intent(out), | allocatable | :: | reduced_masses(:) |
Reduced masses in amu |
|
| real(kind=dp), | intent(out), | allocatable | :: | force_constants(:) |
Force constants in Hartree/Bohr² |
|
| real(kind=dp), | intent(out), | allocatable | :: | cartesian_displacements(:,:) |
Cartesian displacement vectors (3N x 3N) |
|
| real(kind=dp), | intent(out), | optional, | allocatable | :: | eigenvalues_out(:) |
Raw eigenvalues from diagonalization |
| real(kind=dp), | intent(out), | optional, | allocatable | :: | eigenvectors_out(:,:) |
Mass-weighted eigenvectors |
| real(kind=dp), | intent(in), | optional | :: | coordinates(:,:) |
Atomic coordinates in Bohr (3, N) - required for projection |
|
| logical, | intent(in), | optional | :: | project_trans_rot |
If true, project out translation/rotation modes |
|
| real(kind=dp), | intent(out), | optional, | allocatable | :: | force_constants_mdyne(:) |
Force constants in mdyne/Å |
| real(kind=dp), | intent(in), | optional | :: | dipole_derivatives(:,:) |
Cartesian dipole derivatives (3, 3*N) in a.u. for IR intensities |
|
| real(kind=dp), | intent(out), | optional, | allocatable | :: | ir_intensities(:) |
IR intensities in km/mol |
| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| real(kind=dp), | private, | allocatable | :: | eigenvalues(:) | |||
| real(kind=dp), | private, | allocatable | :: | eigenvectors(:,:) |
subroutine compute_vibrational_analysis(hessian, element_numbers, frequencies, & reduced_masses, force_constants, & cartesian_displacements, & eigenvalues_out, eigenvectors_out, & coordinates, project_trans_rot, & force_constants_mdyne, & dipole_derivatives, ir_intensities) !! Perform complete vibrational analysis from Hessian matrix. !! !! This is a convenience wrapper that computes: !! - Vibrational frequencies in cm⁻¹ !! - Reduced masses in amu !! - Force constants in Hartree/Bohr² (and optionally mdyne/Å) !! - Cartesian displacement vectors (normalized) !! - IR intensities in km/mol (if dipole_derivatives provided) !! !! Optionally projects out translation/rotation modes. real(dp), intent(in) :: hessian(:, :) !! Hessian matrix in Hartree/Bohr² (3*N x 3*N) integer, intent(in) :: element_numbers(:) !! Atomic numbers for each atom (N atoms) real(dp), allocatable, intent(out) :: frequencies(:) !! Vibrational frequencies in cm⁻¹ real(dp), allocatable, intent(out) :: reduced_masses(:) !! Reduced masses in amu real(dp), allocatable, intent(out) :: force_constants(:) !! Force constants in Hartree/Bohr² real(dp), allocatable, intent(out) :: cartesian_displacements(:, :) !! Cartesian displacement vectors (3*N x 3*N) real(dp), allocatable, intent(out), optional :: eigenvalues_out(:) !! Raw eigenvalues from diagonalization real(dp), allocatable, intent(out), optional :: eigenvectors_out(:, :) !! Mass-weighted eigenvectors real(dp), intent(in), optional :: coordinates(:, :) !! Atomic coordinates in Bohr (3, N) - required for projection logical, intent(in), optional :: project_trans_rot !! If true, project out translation/rotation modes real(dp), allocatable, intent(out), optional :: force_constants_mdyne(:) !! Force constants in mdyne/Å real(dp), intent(in), optional :: dipole_derivatives(:, :) !! Cartesian dipole derivatives (3, 3*N) in a.u. for IR intensities real(dp), allocatable, intent(out), optional :: ir_intensities(:) !! IR intensities in km/mol real(dp), allocatable :: eigenvalues(:) real(dp), allocatable :: eigenvectors(:, :) ! First compute frequencies and eigenvectors call compute_vibrational_frequencies(hessian, element_numbers, frequencies, & eigenvalues_out=eigenvalues, & eigenvectors=eigenvectors, & coordinates=coordinates, & project_trans_rot=project_trans_rot) ! Compute reduced masses from eigenvectors call compute_reduced_masses(eigenvectors, element_numbers, reduced_masses) ! Compute force constants from eigenvalues and reduced masses call compute_force_constants(eigenvalues, reduced_masses, force_constants, & force_constants_mdyne) ! Compute Cartesian displacements from eigenvectors call compute_cartesian_displacements(eigenvectors, element_numbers, & cartesian_displacements) ! Compute IR intensities if dipole derivatives are provided if (present(dipole_derivatives) .and. present(ir_intensities)) then call compute_ir_intensities(dipole_derivatives, eigenvectors, element_numbers, & ir_intensities) end if ! Optionally return eigenvalues and eigenvectors if (present(eigenvalues_out)) then allocate (eigenvalues_out(size(eigenvalues))) eigenvalues_out = eigenvalues end if if (present(eigenvectors_out)) then allocate (eigenvectors_out(size(eigenvectors, 1), size(eigenvectors, 2))) eigenvectors_out = eigenvectors end if deallocate (eigenvalues, eigenvectors) end subroutine compute_vibrational_analysis