compute_mbe_delta Function

private function compute_mbe_delta(fragment_idx, fragment, lookup, energies, delta_energies, n, world_comm) result(delta_E)

Uses

    • pic_io
  • proc~~compute_mbe_delta~~UsesGraph proc~compute_mbe_delta compute_mbe_delta pic_io pic_io proc~compute_mbe_delta->pic_io

Bottom-up computation of n-body correction (non-recursive, uses pre-computed subset deltas) deltaE(i1,i2,…,in) = E(i1,i2,…,in) - sum of all subset deltaE values All subsets must have been computed already (guaranteed by processing fragments in order)

Arguments

Type IntentOptional Attributes Name
integer(kind=int64), intent(in) :: fragment_idx

Index of this fragment (already known)

integer, intent(in) :: fragment(:)
type(fragment_lookup_t), intent(in) :: lookup

Pre-built hash table for lookups

real(kind=dp), intent(in) :: energies(:)

Pre-computed delta values

real(kind=dp), intent(in) :: delta_energies(:)

Pre-computed delta values

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

MPI communicator for abort

Return Value real(kind=dp)


Calls

proc~~compute_mbe_delta~~CallsGraph proc~compute_mbe_delta compute_mbe_delta abort_comm abort_comm proc~compute_mbe_delta->abort_comm error error proc~compute_mbe_delta->error proc~fragment_lookup_find fragment_lookup_t%fragment_lookup_find proc~compute_mbe_delta->proc~fragment_lookup_find proc~get_next_combination get_next_combination proc~compute_mbe_delta->proc~get_next_combination fnv_1a_hash fnv_1a_hash proc~fragment_lookup_find->fnv_1a_hash proc~arrays_equal_internal arrays_equal_internal proc~fragment_lookup_find->proc~arrays_equal_internal sort sort proc~fragment_lookup_find->sort

Called by

proc~~compute_mbe_delta~~CalledByGraph proc~compute_mbe_delta compute_mbe_delta proc~compute_mbe compute_mbe proc~compute_mbe->proc~compute_mbe_delta 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
logical, private :: has_next
integer, private :: i
integer, private :: indices(MAX_MBE_LEVEL)
integer, private :: subset(MAX_MBE_LEVEL)
integer(kind=int64), private :: subset_idx
integer, private :: subset_size

Source Code

   function compute_mbe_delta(fragment_idx, fragment, lookup, energies, delta_energies, n, world_comm) result(delta_E)
      !! Bottom-up computation of n-body correction (non-recursive, uses pre-computed subset deltas)
      !! deltaE(i1,i2,...,in) = E(i1,i2,...,in) - sum of all subset deltaE values
      !! All subsets must have been computed already (guaranteed by processing fragments in order)
      integer(int64), intent(in) :: fragment_idx  !! Index of this fragment (already known)
      integer, intent(in) :: fragment(:), n
      type(fragment_lookup_t), intent(in) :: lookup  !! Pre-built hash table for lookups
      real(dp), intent(in) :: energies(:), delta_energies(:)  !! Pre-computed delta values
      type(comm_t), intent(in), optional :: world_comm  !! MPI communicator for abort
      real(dp) :: delta_E

      integer :: subset_size, i
      integer :: indices(MAX_MBE_LEVEL), subset(MAX_MBE_LEVEL)  ! Stack arrays to avoid heap contention
      integer(int64) :: subset_idx
      logical :: has_next

      ! Start with the full n-mer energy
      delta_E = energies(fragment_idx)

      ! Subtract all proper subsets (size 1 to n-1)
      do subset_size = 1, n - 1
         ! Initialize first combination
         do i = 1, subset_size
            indices(i) = i
         end do

         ! Loop through all combinations
         do
            ! Build current subset
            do i = 1, subset_size
               subset(i) = fragment(indices(i))
            end do

            ! Look up subset index
            subset_idx = lookup%find(subset(1:subset_size), subset_size)
            if (subset_idx < 0) then
               block
                  use pic_io, only: to_char
                  character(len=512) :: error_msg
                  integer :: j
                  write (error_msg, '(a,i0,a,*(i0,1x))') "Subset not found! Fragment idx=", fragment_idx, &
                     " seeking subset: ", (subset(j), j=1, subset_size)
                  call logger%error(trim(error_msg))
                  write (error_msg, '(a,*(i0,1x))') "  Full fragment: ", (fragment(j), j=1, n)
                  call logger%error(trim(error_msg))
                  if (present(world_comm)) then
                     call abort_comm(world_comm, 1)
                  else
                     error stop "Subset not found in bottom-up MBE!"
                  end if
               end block
            end if

            ! Subtract pre-computed delta energy
            delta_E = delta_E - delta_energies(subset_idx)

            ! Get next combination
            call get_next_combination(indices, subset_size, n, has_next)
            if (.not. has_next) exit
         end do
      end do

   end function compute_mbe_delta