map_fragment_to_system_dipole_derivatives Subroutine

private subroutine map_fragment_to_system_dipole_derivatives(frag_dipole_derivs, monomers, sys_geom, sys_dipole_derivs)

Uses

  • proc~~map_fragment_to_system_dipole_derivatives~~UsesGraph proc~map_fragment_to_system_dipole_derivatives map_fragment_to_system_dipole_derivatives module~mqc_error mqc_error proc~map_fragment_to_system_dipole_derivatives->module~mqc_error module~mqc_physical_fragment mqc_physical_fragment proc~map_fragment_to_system_dipole_derivatives->module~mqc_physical_fragment 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 dipole derivatives to system coordinates with hydrogen cap redistribution

Dipole derivatives have shape (3, 3*N) where each column corresponds to the derivative of dipole w.r.t. one Cartesian coordinate of one atom. Bond connectivity is accessed via sys_geom%bonds

Arguments

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

(3, 3*natoms_frag)

integer, intent(in) :: monomers(:)
type(system_geometry_t), intent(in) :: sys_geom
real(kind=dp), intent(inout) :: sys_dipole_derivs(:,:)

(3, 3*total_atoms)


Calls

proc~~map_fragment_to_system_dipole_derivatives~~CallsGraph proc~map_fragment_to_system_dipole_derivatives map_fragment_to_system_dipole_derivatives proc~build_fragment_from_indices build_fragment_from_indices proc~map_fragment_to_system_dipole_derivatives->proc~build_fragment_from_indices proc~fragment_destroy physical_fragment_t%fragment_destroy proc~map_fragment_to_system_dipole_derivatives->proc~fragment_destroy proc~redistribute_cap_dipole_derivatives redistribute_cap_dipole_derivatives proc~map_fragment_to_system_dipole_derivatives->proc~redistribute_cap_dipole_derivatives 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~error_has_error error_t%error_has_error proc~build_fragment_from_indices->proc~error_has_error proc~fragment_compute_nelec physical_fragment_t%fragment_compute_nelec proc~build_fragment_from_indices->proc~fragment_compute_nelec 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 error error 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_dipole_derivatives~~CalledByGraph proc~map_fragment_to_system_dipole_derivatives map_fragment_to_system_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~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
type(error_t), private :: error
type(physical_fragment_t), private :: fragment

Source Code

   subroutine map_fragment_to_system_dipole_derivatives(frag_dipole_derivs, monomers, sys_geom, sys_dipole_derivs)
      !! Map fragment dipole derivatives to system coordinates with hydrogen cap redistribution
      !!
      !! Dipole derivatives have shape (3, 3*N) where each column corresponds to
      !! the derivative of dipole w.r.t. one Cartesian coordinate of one atom.
      !! Bond connectivity is accessed via sys_geom%bonds
      use mqc_physical_fragment, only: build_fragment_from_indices, redistribute_cap_dipole_derivatives
      use mqc_error, only: error_t
      real(dp), intent(in) :: frag_dipole_derivs(:, :)  !! (3, 3*natoms_frag)
      integer, intent(in) :: monomers(:)
      type(system_geometry_t), intent(in) :: sys_geom
      real(dp), intent(inout) :: sys_dipole_derivs(:, :)  !! (3, 3*total_atoms)

      type(physical_fragment_t) :: fragment
      type(error_t) :: error

      ! Zero out
      sys_dipole_derivs = 0.0_dp

      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)
         call redistribute_cap_dipole_derivatives(fragment, frag_dipole_derivs, sys_dipole_derivs)
         call fragment%destroy()
      else
         ! Code path for fragments without hydrogen caps
         ! Map fragment dipole derivative columns to system positions (fixed-size monomers only)
         block
            integer :: i_mon, i_atom
            integer :: frag_atom_idx, sys_atom_idx
            integer :: frag_col_start, sys_col_start
            integer :: n_monomers

            n_monomers = size(monomers)
            frag_atom_idx = 0

            ! Map each monomer's atoms
            do i_mon = 1, n_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
                  frag_col_start = (frag_atom_idx - 1)*3 + 1
                  sys_col_start = (sys_atom_idx - 1)*3 + 1

                  ! Copy the 3 columns (x, y, z derivatives) for this atom
                  sys_dipole_derivs(:, sys_col_start:sys_col_start + 2) = &
                     frag_dipole_derivs(:, frag_col_start:frag_col_start + 2)
               end do
            end do
         end block
      end if

   end subroutine map_fragment_to_system_dipole_derivatives