redistribute_cap_dipole_derivatives Subroutine

public subroutine redistribute_cap_dipole_derivatives(fragment, fragment_dipole_derivs, system_dipole_derivs)

Redistribute hydrogen cap dipole derivatives to original atoms

Dipole derivatives have shape (3, 3*N_atoms) where each column corresponds to the derivative of the 3 dipole components w.r.t. one Cartesian coordinate. The mapping is similar to the column dimension of the Hessian.

Algorithm: 1. For real atoms: accumulate to system using local_to_global mapping 2. For hydrogen caps: add to the original atom the cap replaces

Arguments

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

(3, 3*n_atoms_fragment)

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

(3, 3*n_atoms_system)


Called by

proc~~redistribute_cap_dipole_derivatives~~CalledByGraph proc~redistribute_cap_dipole_derivatives redistribute_cap_dipole_derivatives proc~gmbe_pie_coordinator gmbe_pie_coordinator proc~gmbe_pie_coordinator->proc~redistribute_cap_dipole_derivatives proc~map_fragment_to_system_dipole_derivatives map_fragment_to_system_dipole_derivatives proc~map_fragment_to_system_dipole_derivatives->proc~redistribute_cap_dipole_derivatives proc~serial_gmbe_pie_processor serial_gmbe_pie_processor proc~serial_gmbe_pie_processor->proc~redistribute_cap_dipole_derivatives proc~compute_mbe compute_mbe proc~compute_mbe->proc~map_fragment_to_system_dipole_derivatives proc~compute_mbe_dipole_derivatives compute_mbe_dipole_derivatives proc~compute_mbe->proc~compute_mbe_dipole_derivatives proc~compute_mbe_dipole_derivatives->proc~map_fragment_to_system_dipole_derivatives 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_i
integer, private :: global_original_idx
integer, private :: i
integer, private :: i_cap
integer, private :: icart
integer, private :: local_cap_idx
integer, private :: local_i
integer, private :: n_real_atoms

Source Code

   subroutine redistribute_cap_dipole_derivatives(fragment, fragment_dipole_derivs, system_dipole_derivs)
      !! Redistribute hydrogen cap dipole derivatives to original atoms
      !!
      !! Dipole derivatives have shape (3, 3*N_atoms) where each column corresponds to
      !! the derivative of the 3 dipole components w.r.t. one Cartesian coordinate.
      !! The mapping is similar to the column dimension of the Hessian.
      !!
      !! Algorithm:
      !!   1. For real atoms: accumulate to system using local_to_global mapping
      !!   2. For hydrogen caps: add to the original atom the cap replaces
      type(physical_fragment_t), intent(in) :: fragment
      real(dp), intent(in) :: fragment_dipole_derivs(:, :)   !! (3, 3*n_atoms_fragment)
      real(dp), intent(inout) :: system_dipole_derivs(:, :)  !! (3, 3*n_atoms_system)

      integer :: i, local_i, global_i, icart
      integer :: i_cap, local_cap_idx, global_original_idx
      integer :: n_real_atoms

      n_real_atoms = fragment%n_atoms - fragment%n_caps

      ! Accumulate dipole derivative columns for real atoms
      do i = 1, n_real_atoms
         local_i = i
         global_i = fragment%local_to_global(local_i)
         if (global_i <= 0) cycle

         do icart = 1, 3
            system_dipole_derivs(:, (global_i - 1)*3 + icart) = &
               system_dipole_derivs(:, (global_i - 1)*3 + icart) + &
               fragment_dipole_derivs(:, (local_i - 1)*3 + icart)
         end do
      end do

      ! Redistribute cap contributions to their original atoms
      if (fragment%n_caps > 0 .and. allocated(fragment%cap_replaces_atom)) 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

            do icart = 1, 3
               system_dipole_derivs(:, (global_original_idx - 1)*3 + icart) = &
                  system_dipole_derivs(:, (global_original_idx - 1)*3 + icart) + &
                  fragment_dipole_derivs(:, (local_cap_idx - 1)*3 + icart)
            end do
         end do
      end if

   end subroutine redistribute_cap_dipole_derivatives