map_fragment_to_system_gradient Subroutine

private subroutine map_fragment_to_system_gradient(frag_grad, monomers, sys_geom, sys_grad, world_comm)

Uses

  • proc~~map_fragment_to_system_gradient~~UsesGraph proc~map_fragment_to_system_gradient map_fragment_to_system_gradient module~mqc_error mqc_error proc~map_fragment_to_system_gradient->module~mqc_error module~mqc_physical_fragment mqc_physical_fragment proc~map_fragment_to_system_gradient->module~mqc_physical_fragment pic_logger pic_logger proc~map_fragment_to_system_gradient->pic_logger module~mqc_physical_fragment->module~mqc_error module~mqc_cgto mqc_cgto module~mqc_physical_fragment->module~mqc_cgto module~mqc_elements mqc_elements module~mqc_physical_fragment->module~mqc_elements module~mqc_geometry mqc_geometry module~mqc_physical_fragment->module~mqc_geometry module~mqc_physical_constants mqc_physical_constants module~mqc_physical_fragment->module~mqc_physical_constants module~mqc_xyz_reader mqc_xyz_reader module~mqc_physical_fragment->module~mqc_xyz_reader pic_types pic_types module~mqc_physical_fragment->pic_types module~mqc_cgto->pic_types module~mqc_elements->pic_types pic_ascii pic_ascii module~mqc_elements->pic_ascii module~mqc_geometry->pic_types module~mqc_physical_constants->pic_types module~mqc_xyz_reader->module~mqc_error module~mqc_xyz_reader->module~mqc_geometry module~mqc_xyz_reader->pic_types

Map fragment gradient to system gradient coordinates with hydrogen cap redistribution

This function rebuilds the fragment to get local→global mappings and cap information, then redistributes gradients including hydrogen caps to their original atoms. Bond connectivity is accessed via sys_geom%bonds.

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: frag_grad(:,:)

(3, natoms_frag)

integer, intent(in) :: monomers(:)

Monomer indices in fragment

type(system_geometry_t), intent(in) :: sys_geom
real(kind=dp), intent(inout) :: sys_grad(:,:)

(3, total_atoms)

type(comm_t), intent(in), optional :: world_comm

MPI communicator for abort


Calls

proc~~map_fragment_to_system_gradient~~CallsGraph proc~map_fragment_to_system_gradient map_fragment_to_system_gradient abort_comm abort_comm proc~map_fragment_to_system_gradient->abort_comm configuration configuration proc~map_fragment_to_system_gradient->configuration debug debug proc~map_fragment_to_system_gradient->debug error error proc~map_fragment_to_system_gradient->error proc~build_fragment_from_indices build_fragment_from_indices proc~map_fragment_to_system_gradient->proc~build_fragment_from_indices proc~error_get_full_trace error_t%error_get_full_trace proc~map_fragment_to_system_gradient->proc~error_get_full_trace proc~error_has_error error_t%error_has_error proc~map_fragment_to_system_gradient->proc~error_has_error proc~fragment_destroy physical_fragment_t%fragment_destroy proc~map_fragment_to_system_gradient->proc~fragment_destroy proc~redistribute_cap_gradients redistribute_cap_gradients proc~map_fragment_to_system_gradient->proc~redistribute_cap_gradients proc~build_fragment_from_indices->proc~error_has_error proc~add_hydrogen_caps add_hydrogen_caps proc~build_fragment_from_indices->proc~add_hydrogen_caps proc~calculate_monomer_distance calculate_monomer_distance proc~build_fragment_from_indices->proc~calculate_monomer_distance proc~check_duplicate_atoms check_duplicate_atoms proc~build_fragment_from_indices->proc~check_duplicate_atoms proc~count_hydrogen_caps count_hydrogen_caps proc~build_fragment_from_indices->proc~count_hydrogen_caps proc~error_add_context error_t%error_add_context proc~build_fragment_from_indices->proc~error_add_context proc~fragment_compute_nelec physical_fragment_t%fragment_compute_nelec proc~build_fragment_from_indices->proc~fragment_compute_nelec proc~error_get_full_trace->proc~error_has_error proc~basis_set_destroy molecular_basis_type%basis_set_destroy proc~fragment_destroy->proc~basis_set_destroy proc~atomic_basis_destroy atomic_basis_type%atomic_basis_destroy proc~basis_set_destroy->proc~atomic_basis_destroy proc~to_angstrom to_angstrom proc~calculate_monomer_distance->proc~to_angstrom proc~check_duplicate_atoms->error proc~element_number_to_symbol element_number_to_symbol proc~check_duplicate_atoms->proc~element_number_to_symbol proc~error_set error_t%error_set proc~check_duplicate_atoms->proc~error_set to_char to_char proc~check_duplicate_atoms->to_char proc~cgto_destroy cgto_type%cgto_destroy proc~atomic_basis_destroy->proc~cgto_destroy

Called by

proc~~map_fragment_to_system_gradient~~CalledByGraph proc~map_fragment_to_system_gradient map_fragment_to_system_gradient 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~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 :: current_log_level
type(error_t), private :: error
integer, private :: frag_atom_idx
type(physical_fragment_t), private :: fragment
integer, private :: i_atom
integer, private :: i_mon
integer, private :: sys_atom_idx

Source Code

   subroutine map_fragment_to_system_gradient(frag_grad, monomers, sys_geom, sys_grad, world_comm)
      !! Map fragment gradient to system gradient coordinates with hydrogen cap redistribution
      !!
      !! This function rebuilds the fragment to get local→global mappings and cap information,
      !! then redistributes gradients including hydrogen caps to their original atoms.
      !! Bond connectivity is accessed via sys_geom%bonds.
      use mqc_physical_fragment, only: build_fragment_from_indices, redistribute_cap_gradients
      use mqc_error, only: error_t
      use pic_logger, only: verbose_level
      real(dp), intent(in) :: frag_grad(:, :)  !! (3, natoms_frag)
      integer, intent(in) :: monomers(:)  !! Monomer indices in fragment
      type(system_geometry_t), intent(in) :: sys_geom
      real(dp), intent(inout) :: sys_grad(:, :)  !! (3, total_atoms)
      type(comm_t), intent(in), optional :: world_comm  !! MPI communicator for abort

      type(physical_fragment_t) :: fragment
      type(error_t) :: error
      integer :: i_mon, i_atom, frag_atom_idx, sys_atom_idx
      integer :: current_log_level

      ! Explicitly zero out the entire sys_grad array
      sys_grad = 0.0_dp

      ! Debug output
      call logger%configuration(level=current_log_level)
      if (current_log_level >= debug_level) then
         block
            character(len=256) :: debug_msg
            write (debug_msg, '(a,i0,a,*(i0,1x))') "  Mapping fragment with ", size(monomers), " monomers: ", monomers
            call logger%debug(trim(debug_msg))
            write (debug_msg, '(a,i0,a)') "  Fragment has ", size(frag_grad, 2), " atoms"
            call logger%debug(trim(debug_msg))
         end block
      end if

      if (allocated(sys_geom%bonds)) then
         ! Rebuild fragment to get local→global mapping and cap information
         call build_fragment_from_indices(sys_geom, monomers, fragment, error, sys_geom%bonds)
         if (error%has_error()) then
            call logger%error(error%get_full_trace())
            if (present(world_comm)) then
               call abort_comm(world_comm, 1)
            else
               error stop "Failed to build fragment in gradient mapping"
            end if
         end if

         ! Use new gradient redistribution with cap handling
         call redistribute_cap_gradients(fragment, frag_grad, sys_grad)

         ! Clean up
         call fragment%destroy()
      else
         ! Old code path for fragments without hydrogen caps
         ! Map fragment gradient to system positions (fixed-size monomers only)
         frag_atom_idx = 0
         do i_mon = 1, size(monomers)
            do i_atom = 1, sys_geom%atoms_per_monomer
               frag_atom_idx = frag_atom_idx + 1
               sys_atom_idx = (monomers(i_mon) - 1)*sys_geom%atoms_per_monomer + i_atom

               if (current_log_level >= debug_level .and. i_atom == 1) then
                  block
                     character(len=256) :: debug_msg
                     write (debug_msg, '(a,i0,a,i0,a,i0)') &
                        "    Monomer ", monomers(i_mon), ": frag atoms ", &
                        frag_atom_idx, " -> sys atom ", sys_atom_idx
                     call logger%debug(trim(debug_msg))
                  end block
               end if

               sys_grad(:, sys_atom_idx) = frag_grad(:, frag_atom_idx)
            end do
         end do
      end if

   end subroutine map_fragment_to_system_gradient