count_hydrogen_caps Subroutine

private subroutine count_hydrogen_caps(atoms_in_fragment, bonds, n_caps)

Count how many hydrogen caps are needed for a fragment A cap is needed when exactly one atom of a broken bond is in the fragment

Arguments

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

0-indexed atom indices in fragment

type(bond_t), intent(in), optional :: bonds(:)
integer, intent(out) :: n_caps

Called by

proc~~count_hydrogen_caps~~CalledByGraph proc~count_hydrogen_caps count_hydrogen_caps proc~build_fragment_from_atom_list build_fragment_from_atom_list proc~build_fragment_from_atom_list->proc~count_hydrogen_caps proc~build_fragment_from_indices build_fragment_from_indices proc~build_fragment_from_indices->proc~count_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 :: ibond

Source Code

   subroutine count_hydrogen_caps(atoms_in_fragment, bonds, n_caps)
      !! Count how many hydrogen caps are needed for a fragment
      !! A cap is needed when exactly one atom of a broken bond is in the fragment
      integer, intent(in) :: atoms_in_fragment(:)  !! 0-indexed atom indices in fragment
      type(bond_t), intent(in), optional :: bonds(:)
      integer, intent(out) :: n_caps

      integer :: ibond
      logical :: atom_i_in_frag, atom_j_in_frag

      n_caps = 0
      if (.not. present(bonds)) return

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

         ! Check if exactly one atom of this bond is in the fragment
         atom_i_in_frag = any(atoms_in_fragment == bonds(ibond)%atom_i)
         atom_j_in_frag = any(atoms_in_fragment == bonds(ibond)%atom_j)

         ! Add cap only if one atom in fragment, other not (XOR condition)
         if ((atom_i_in_frag .and. .not. atom_j_in_frag) .or. &
             (.not. atom_i_in_frag .and. atom_j_in_frag)) then
            n_caps = n_caps + 1
         end if
      end do

   end subroutine count_hydrogen_caps