add_hydrogen_caps Subroutine

private subroutine add_hydrogen_caps(atoms_in_fragment, bonds, sys_geom, fragment, base_atom_count)

Add hydrogen caps to fragment for broken bonds Caps are placed at the position of the atom outside the fragment

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: atoms_in_fragment(:)

0-indexed atom indices in fragment

type(bond_t), intent(in) :: bonds(:)
type(system_geometry_t), intent(in) :: sys_geom
type(physical_fragment_t), intent(inout) :: fragment
integer, intent(in) :: base_atom_count

Number of non-cap atoms


Called by

proc~~add_hydrogen_caps~~CalledByGraph proc~add_hydrogen_caps add_hydrogen_caps proc~build_fragment_from_atom_list build_fragment_from_atom_list proc~build_fragment_from_atom_list->proc~add_hydrogen_caps proc~build_fragment_from_indices build_fragment_from_indices proc~build_fragment_from_indices->proc~add_hydrogen_caps proc~compute_gmbe compute_gmbe proc~compute_gmbe->proc~build_fragment_from_indices proc~process_intersection_derivatives process_intersection_derivatives proc~compute_gmbe->proc~process_intersection_derivatives proc~gmbe_pie_coordinator gmbe_pie_coordinator proc~gmbe_pie_coordinator->proc~build_fragment_from_atom_list proc~map_fragment_to_system_dipole_derivatives map_fragment_to_system_dipole_derivatives proc~map_fragment_to_system_dipole_derivatives->proc~build_fragment_from_indices proc~map_fragment_to_system_gradient map_fragment_to_system_gradient proc~map_fragment_to_system_gradient->proc~build_fragment_from_indices proc~map_fragment_to_system_hessian map_fragment_to_system_hessian proc~map_fragment_to_system_hessian->proc~build_fragment_from_indices proc~node_worker node_worker proc~node_worker->proc~build_fragment_from_atom_list proc~node_worker->proc~build_fragment_from_indices proc~process_intersection_derivatives->proc~build_fragment_from_atom_list proc~serial_fragment_processor serial_fragment_processor proc~serial_fragment_processor->proc~build_fragment_from_indices proc~compute_mbe compute_mbe proc~serial_fragment_processor->proc~compute_mbe proc~serial_gmbe_pie_processor serial_gmbe_pie_processor proc~serial_gmbe_pie_processor->proc~build_fragment_from_atom_list interface~node_worker node_worker interface~node_worker->proc~node_worker interface~serial_fragment_processor serial_fragment_processor interface~serial_fragment_processor->proc~serial_fragment_processor proc~compute_mbe->proc~map_fragment_to_system_dipole_derivatives proc~compute_mbe->proc~map_fragment_to_system_gradient proc~compute_mbe->proc~map_fragment_to_system_hessian proc~compute_mbe_dipole_derivatives compute_mbe_dipole_derivatives proc~compute_mbe->proc~compute_mbe_dipole_derivatives proc~compute_mbe_gradient compute_mbe_gradient proc~compute_mbe->proc~compute_mbe_gradient proc~compute_mbe_hessian compute_mbe_hessian proc~compute_mbe->proc~compute_mbe_hessian proc~compute_mbe_dipole_derivatives->proc~map_fragment_to_system_dipole_derivatives proc~compute_mbe_gradient->proc~map_fragment_to_system_gradient proc~compute_mbe_hessian->proc~map_fragment_to_system_hessian proc~gmbe_run_distributed gmbe_context_t%gmbe_run_distributed proc~gmbe_run_distributed->proc~gmbe_pie_coordinator proc~gmbe_run_distributed->interface~node_worker proc~gmbe_run_serial gmbe_context_t%gmbe_run_serial proc~gmbe_run_serial->proc~serial_gmbe_pie_processor proc~global_coordinator global_coordinator proc~global_coordinator->proc~compute_mbe proc~mbe_run_distributed mbe_context_t%mbe_run_distributed proc~mbe_run_distributed->interface~node_worker proc~mbe_run_serial mbe_context_t%mbe_run_serial proc~mbe_run_serial->interface~serial_fragment_processor interface~global_coordinator global_coordinator interface~global_coordinator->proc~global_coordinator

Variables

Type Visibility Attributes Name Initial
logical, private :: atom_i_in_frag
logical, private :: atom_j_in_frag
integer, private :: cap_idx
integer, private :: ibond

Source Code

   subroutine add_hydrogen_caps(atoms_in_fragment, bonds, sys_geom, fragment, base_atom_count)
      !! Add hydrogen caps to fragment for broken bonds
      !! Caps are placed at the position of the atom outside the fragment
      integer, intent(in) :: atoms_in_fragment(:)  !! 0-indexed atom indices in fragment
      type(bond_t), intent(in) :: bonds(:)
      type(system_geometry_t), intent(in) :: sys_geom
      type(physical_fragment_t), intent(inout) :: fragment
      integer, intent(in) :: base_atom_count  !! Number of non-cap atoms

      integer :: ibond, cap_idx
      logical :: atom_i_in_frag, atom_j_in_frag

      if (fragment%n_caps == 0) return

      cap_idx = 0
      do ibond = 1, size(bonds)
         if (.not. bonds(ibond)%is_broken) cycle

         atom_i_in_frag = any(atoms_in_fragment == bonds(ibond)%atom_i)
         atom_j_in_frag = any(atoms_in_fragment == bonds(ibond)%atom_j)

         if (atom_i_in_frag .and. .not. atom_j_in_frag) then
            ! atom_i is in fragment, atom_j is not → cap at position of atom_j
            cap_idx = cap_idx + 1
            fragment%element_numbers(base_atom_count + cap_idx) = 1  ! Hydrogen
            ! Place H at position of atom_j (1-indexed for coordinates array)
            fragment%coordinates(:, base_atom_count + cap_idx) = &
               sys_geom%coordinates(:, bonds(ibond)%atom_j + 1)
            fragment%cap_replaces_atom(cap_idx) = bonds(ibond)%atom_j

         else if (atom_j_in_frag .and. .not. atom_i_in_frag) then
            ! atom_j is in fragment, atom_i is not → cap at position of atom_i
            cap_idx = cap_idx + 1
            fragment%element_numbers(base_atom_count + cap_idx) = 1  ! Hydrogen
            ! Place H at position of atom_i (1-indexed for coordinates array)
            fragment%coordinates(:, base_atom_count + cap_idx) = &
               sys_geom%coordinates(:, bonds(ibond)%atom_i + 1)
            fragment%cap_replaces_atom(cap_idx) = bonds(ibond)%atom_i
         end if
      end do

   end subroutine add_hydrogen_caps