redistribute_cap_gradients Subroutine

public subroutine redistribute_cap_gradients(fragment, fragment_gradient, system_gradient)

Redistribute hydrogen cap gradients to original atoms

This subroutine handles gradient redistribution for fragments with hydrogen caps. Hydrogen caps are virtual atoms added at broken bonds - their gradients represent forces at the bond breakpoint and must be transferred to the original atoms they replace.

Algorithm: 1. For real atoms (indices 1 to n_atoms - n_caps): Accumulate gradient to system using local_to_global mapping 2. For hydrogen caps (indices n_atoms - n_caps + 1 to n_atoms): Add cap gradient to the original atom it replaces (from cap_replaces_atom)

Example: Fragment: [C, C, H_cap] where H_cap replaces atom 5 in system Fragment gradient: [(3,1), (3,2), (3,3)] - Atoms 1,2: accumulate to system using local_to_global - Atom 3 (cap): add gradient to system atom 5 (cap_replaces_atom(1) + 1)

Arguments

Type IntentOptional Attributes Name
type(physical_fragment_t), intent(in) :: fragment
real(kind=dp), intent(in) :: fragment_gradient(:,:)

(3, n_atoms_fragment)

real(kind=dp), intent(inout) :: system_gradient(:,:)

(3, n_atoms_system)


Called by

proc~~redistribute_cap_gradients~~CalledByGraph proc~redistribute_cap_gradients redistribute_cap_gradients proc~compute_gmbe compute_gmbe proc~compute_gmbe->proc~redistribute_cap_gradients 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~redistribute_cap_gradients proc~map_fragment_to_system_gradient map_fragment_to_system_gradient proc~map_fragment_to_system_gradient->proc~redistribute_cap_gradients proc~process_intersection_derivatives->proc~redistribute_cap_gradients proc~serial_gmbe_pie_processor serial_gmbe_pie_processor proc~serial_gmbe_pie_processor->proc~redistribute_cap_gradients proc~compute_mbe compute_mbe proc~compute_mbe->proc~map_fragment_to_system_gradient proc~compute_mbe_gradient compute_mbe_gradient proc~compute_mbe->proc~compute_mbe_gradient proc~compute_mbe_gradient->proc~map_fragment_to_system_gradient 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~global_coordinator global_coordinator proc~global_coordinator->proc~compute_mbe 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~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

Variables

Type Visibility Attributes Name Initial
integer, private :: global_idx
integer, private :: global_original_idx
integer, private :: i
integer, private :: i_cap
integer, private :: local_cap_idx
integer, private :: local_idx
integer, private :: n_real_atoms

Source Code

   subroutine redistribute_cap_gradients(fragment, fragment_gradient, system_gradient)
      !! Redistribute hydrogen cap gradients to original atoms
      !!
      !! This subroutine handles gradient redistribution for fragments with hydrogen caps.
      !! Hydrogen caps are virtual atoms added at broken bonds - their gradients represent
      !! forces at the bond breakpoint and must be transferred to the original atoms they replace.
      !!
      !! Algorithm:
      !!   1. For real atoms (indices 1 to n_atoms - n_caps):
      !!      Accumulate gradient to system using local_to_global mapping
      !!   2. For hydrogen caps (indices n_atoms - n_caps + 1 to n_atoms):
      !!      Add cap gradient to the original atom it replaces (from cap_replaces_atom)
      !!
      !! Example:
      !!   Fragment: [C, C, H_cap] where H_cap replaces atom 5 in system
      !!   Fragment gradient: [(3,1), (3,2), (3,3)]
      !!   - Atoms 1,2: accumulate to system using local_to_global
      !!   - Atom 3 (cap): add gradient to system atom 5 (cap_replaces_atom(1) + 1)
      type(physical_fragment_t), intent(in) :: fragment
      real(dp), intent(in) :: fragment_gradient(:, :)   !! (3, n_atoms_fragment)
      real(dp), intent(inout) :: system_gradient(:, :)  !! (3, n_atoms_system)

      integer :: i, local_idx, global_idx
      integer :: i_cap, local_cap_idx, global_original_idx
      integer :: n_real_atoms

      n_real_atoms = fragment%n_atoms - fragment%n_caps

      ! Accumulate gradients for real atoms using local→global mapping
      do i = 1, n_real_atoms
         global_idx = fragment%local_to_global(i)
         system_gradient(:, global_idx) = system_gradient(:, global_idx) + fragment_gradient(:, i)
      end do

      ! Redistribute cap gradients to original atoms they replace
      if (fragment%n_caps > 0) then
         do i_cap = 1, fragment%n_caps
            local_cap_idx = n_real_atoms + i_cap
            ! cap_replaces_atom is 0-indexed, add 1 for Fortran arrays
            global_original_idx = fragment%cap_replaces_atom(i_cap) + 1

            ! Add cap gradient to the atom it replaces
            system_gradient(:, global_original_idx) = system_gradient(:, global_original_idx) + &
                                                      fragment_gradient(:, local_cap_idx)
         end do
      end if

   end subroutine redistribute_cap_gradients