check_duplicate_atoms Subroutine

public subroutine check_duplicate_atoms(fragment, error)

Uses

    • pic_logger
    • pic_io
  • proc~~check_duplicate_atoms~~UsesGraph proc~check_duplicate_atoms check_duplicate_atoms pic_io pic_io proc~check_duplicate_atoms->pic_io pic_logger pic_logger proc~check_duplicate_atoms->pic_logger

Validate that fragment has no spatially overlapping atoms Checks if any two atoms are too close together (< 0.01 Bohr) This catches bugs in geometry construction or fragment building

Arguments

Type IntentOptional Attributes Name
type(physical_fragment_t), intent(in) :: fragment
type(error_t), intent(out) :: error

Calls

proc~~check_duplicate_atoms~~CallsGraph proc~check_duplicate_atoms check_duplicate_atoms 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

Called by

proc~~check_duplicate_atoms~~CalledByGraph proc~check_duplicate_atoms check_duplicate_atoms proc~build_fragment_from_atom_list build_fragment_from_atom_list proc~build_fragment_from_atom_list->proc~check_duplicate_atoms proc~build_fragment_from_indices build_fragment_from_indices proc~build_fragment_from_indices->proc~check_duplicate_atoms proc~unfragmented_calculation unfragmented_calculation proc~unfragmented_calculation->proc~check_duplicate_atoms interface~unfragmented_calculation unfragmented_calculation interface~unfragmented_calculation->proc~unfragmented_calculation 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~run_unfragmented_calculation run_unfragmented_calculation proc~run_unfragmented_calculation->interface~unfragmented_calculation 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 proc~run_calculation run_calculation proc~run_calculation->proc~run_unfragmented_calculation interface~global_coordinator global_coordinator interface~global_coordinator->proc~global_coordinator proc~compute_energy_and_forces compute_energy_and_forces proc~compute_energy_and_forces->proc~run_calculation proc~run_multi_molecule_calculations run_multi_molecule_calculations proc~run_multi_molecule_calculations->proc~run_calculation program~main main program~main->proc~run_calculation

Variables

Type Visibility Attributes Name Initial
real(kind=dp), private, parameter :: MIN_ATOM_DISTANCE = 0.01_dp

Bohr - atoms closer than this are overlapping

real(kind=dp), private :: distance
real(kind=dp), private :: dx
real(kind=dp), private :: dy
real(kind=dp), private :: dz
integer, private :: i
integer, private :: j
integer, private :: n_atoms

Source Code

   subroutine check_duplicate_atoms(fragment, error)
      !! Validate that fragment has no spatially overlapping atoms
      !! Checks if any two atoms are too close together (< 0.01 Bohr)
      !! This catches bugs in geometry construction or fragment building
      use pic_logger, only: logger => global_logger
      use pic_io, only: to_char

      type(physical_fragment_t), intent(in) :: fragment
      type(error_t), intent(out) :: error

      integer :: i, j, n_atoms
      real(dp) :: distance, dx, dy, dz
      real(dp), parameter :: MIN_ATOM_DISTANCE = 0.01_dp  !! Bohr - atoms closer than this are overlapping

      ! Only check non-cap atoms (caps can be close to replaced atoms)
      n_atoms = fragment%n_atoms - fragment%n_caps

      if (n_atoms < 2) return

      do i = 1, n_atoms - 1
         do j = i + 1, n_atoms
            dx = fragment%coordinates(1, i) - fragment%coordinates(1, j)
            dy = fragment%coordinates(2, i) - fragment%coordinates(2, j)
            dz = fragment%coordinates(3, i) - fragment%coordinates(3, j)
            distance = sqrt(dx*dx + dy*dy + dz*dz)

            if (distance < MIN_ATOM_DISTANCE) then
               ! Build detailed error message
               call error%set(ERROR_VALIDATION, &
                              "Fragment contains overlapping atoms "//to_char(i)//" and "//to_char(j)// &
                              " (distance: "//to_char(distance)//" Bohr). "// &
                              "This indicates bad input geometry or a bug in fragment construction.")

               ! Log detailed information for debugging
               call logger%error("ERROR: Fragment contains overlapping atoms!")
               call logger%error("  Atoms "//to_char(i)//" and "//to_char(j)//" are too close together")
               call logger%error("  Distance: "//to_char(distance)//" Bohr ("// &
                                 to_char(distance*0.529177_dp)//" Angstrom)")
               call logger%error("  Atom "//to_char(i)//": "// &
                                 element_number_to_symbol(fragment%element_numbers(i))// &
                                 " at ("//to_char(fragment%coordinates(1, i))//", "// &
                                 to_char(fragment%coordinates(2, i))//", "// &
                                 to_char(fragment%coordinates(3, i))//") Bohr")
               call logger%error("  Atom "//to_char(j)//": "// &
                                 element_number_to_symbol(fragment%element_numbers(j))// &
                                 " at ("//to_char(fragment%coordinates(1, j))//", "// &
                                 to_char(fragment%coordinates(2, j))//", "// &
                                 to_char(fragment%coordinates(3, j))//") Bohr")
               return
            end if
         end do
      end do
   end subroutine check_duplicate_atoms