result_isend Subroutine

public subroutine result_isend(result, comm, dest, tag, req)

Send calculation result over MPI (non-blocking) Sends SCF energy (non-blocking) and other components (blocking)

Arguments

Type IntentOptional Attributes Name
type(calculation_result_t), intent(in) :: result
type(comm_t), intent(in) :: comm
integer, intent(in) :: dest
integer, intent(in) :: tag
type(request_t), intent(out) :: req

Calls

proc~~result_isend~~CallsGraph proc~result_isend result_isend isend isend proc~result_isend->isend send send proc~result_isend->send

Called by

proc~~result_isend~~CalledByGraph proc~result_isend result_isend proc~node_coordinator node_coordinator proc~node_coordinator->proc~result_isend proc~node_worker node_worker proc~node_worker->proc~result_isend interface~node_coordinator node_coordinator interface~node_coordinator->proc~node_coordinator interface~node_worker node_worker interface~node_worker->proc~node_worker proc~gmbe_run_distributed gmbe_context_t%gmbe_run_distributed proc~gmbe_run_distributed->interface~node_coordinator proc~gmbe_run_distributed->interface~node_worker proc~mbe_run_distributed mbe_context_t%mbe_run_distributed proc~mbe_run_distributed->interface~node_coordinator proc~mbe_run_distributed->interface~node_worker

Source Code

   subroutine result_isend(result, comm, dest, tag, req)
      !! Send calculation result over MPI (non-blocking)
      !! Sends SCF energy (non-blocking) and other components (blocking)
      type(calculation_result_t), intent(in) :: result
      type(comm_t), intent(in) :: comm
      integer, intent(in) :: dest, tag
      type(request_t), intent(out) :: req

      ! Send SCF energy (non-blocking)
      call isend(comm, result%energy%scf, dest, tag, req)

      ! Send other energy components (blocking to avoid needing multiple request handles)
      call send(comm, result%energy%mp2%ss, dest, tag)
      call send(comm, result%energy%mp2%os, dest, tag)
      call send(comm, result%energy%cc%singles, dest, tag)
      call send(comm, result%energy%cc%doubles, dest, tag)
      call send(comm, result%energy%cc%triples, dest, tag)

      ! Send fragment metadata
      call send(comm, result%distance, dest, tag)

      ! Send gradient flag and data (blocking to avoid needing multiple request handles)
      call send(comm, result%has_gradient, dest, tag)
      if (result%has_gradient) then
         call send(comm, result%gradient, dest, tag)
      end if

      ! Send Hessian flag and data (blocking to avoid needing multiple request handles)
      call send(comm, result%has_hessian, dest, tag)
      if (result%has_hessian) then
         call send(comm, result%hessian, dest, tag)
      end if

      ! Send dipole flag and data (blocking to avoid needing multiple request handles)
      call send(comm, result%has_dipole, dest, tag)
      if (result%has_dipole) then
         call send(comm, result%dipole, dest, tag)
      end if

      ! Send dipole derivatives flag and data (blocking to avoid needing multiple request handles)
      call send(comm, result%has_dipole_derivatives, dest, tag)
      if (result%has_dipole_derivatives) then
         call send(comm, result%dipole_derivatives, dest, tag)
      end if
   end subroutine result_isend