57    integer :: nbatch_to_receive
 
   58    integer :: nbatch_to_send
 
   61    integer :: nblock_to_receive
 
   62    integer :: nblock_to_send
 
   77    class(states_elec_all_to_all_communications_t),  
intent(inout) :: this
 
   78    type(states_elec_t),   
intent(in) :: st
 
   79    integer,               
intent(in) :: task_from, task_to
 
   83    this%task_from = task_from
 
   84    this%task_to = task_to
 
   90    this%n_comms = max(this%nbatch_to_send, this%nbatch_to_receive)
 
   97    type(states_elec_t),   
intent(in)  :: st
 
   98    integer,               
intent(in)  :: task_from
 
   99    integer,               
intent(out) :: nblock_to_receive
 
  101    integer :: st_start, st_end, kpt_start, kpt_end, ib
 
  105    nbatch_to_receive = 0
 
  106    nblock_to_receive = 0
 
  109    if (task_from > -1) 
then 
  110      st_start   = st%st_kpt_task(task_from, 1)
 
  111      st_end     = st%st_kpt_task(task_from, 2)
 
  112      kpt_start  = st%st_kpt_task(task_from, 3)
 
  113      kpt_end    = st%st_kpt_task(task_from, 4)
 
  115      nblock_to_receive = 0
 
  116      do ib = 1, st%group%nblocks
 
  117        if (st%group%block_range(ib, 1) >= st_start .and. st%group%block_range(ib, 2) <= st_end) 
then 
  118          nblock_to_receive = nblock_to_receive + 1
 
  121      nbatch_to_receive = nblock_to_receive * (kpt_end-kpt_start+1)
 
  125      write(
message(1), 
'(a,i2,a,i2,a,i2)') 
'Task ', st%st_kpt_mpi_grp%rank, 
' will receive ', &
 
  126        nbatch_to_receive, 
' batches from task ', task_from
 
  135    type(states_elec_t),   
intent(in)  :: st
 
  136    integer,               
intent(in)  :: task_to
 
  137    integer,               
intent(out) :: nblock_to_send
 
  144    if (task_to > -1) 
then 
  145      nblock_to_send = (st%group%block_end-st%group%block_start+1)
 
  146      nbatch_to_send = nblock_to_send*(st%d%kpt%end-st%d%kpt%start+1)
 
  150      write(
message(1), 
'(a,i2,a,i2,a,i2)') 
'Task ', st%st_kpt_mpi_grp%rank, 
' will send ', nbatch_to_send, &
 
  151        ' batches to task ', task_to
 
  162    n_comms = this%n_comms
 
  169    nbatch_to_receive = this%nbatch_to_receive
 
  173  integer pure function states_elec_all_to_all_communications_get_nsend(this) result(nbatch_to_send)
 
  176    nbatch_to_send = this%nbatch_to_send
 
  182    type(states_elec_t),                             
intent(in)  :: st
 
  183    integer,                                         
intent(in)  :: icom
 
  184    integer,                                         
intent(in)  :: np
 
  185    class(wfs_elec_t), 
pointer,                      
intent(out) :: psib
 
  187    integer :: block_id, ib, ik
 
  192    block_id = mod(icom-1, this%nblock_to_receive)+1
 
  193    ik = int((icom-block_id)/this%nblock_to_receive) + st%st_kpt_task(this%task_from, 3)
 
  194    ib = block_id - 1 + st%group%iblock(st%st_kpt_task(this%task_from, 1))
 
  197      write(message(1), 
'(a,i2,a,i2,a,i2)') 
'Task ', st%st_kpt_mpi_grp%rank, 
' allocates memory for block ', &
 
  198        ib, 
' and k-point ', ik
 
  199      call messages_info(1, all_nodes=.
true.)
 
  202    call states_elec_parallel_allocate_batch(st, psib, np, ib, ik, packed=.
true.)
 
  210    type(states_elec_t),                             
intent(in)  :: st
 
  211    integer,                                         
intent(in)  :: icom
 
  212    integer,                                         
intent(out) :: ib
 
  213    integer,                                         
intent(out) :: ik
 
  218    ib = mod(icom-1, this%nblock_to_send) + 1
 
  219    ik = int((icom-ib)/this%nblock_to_send) + st%d%kpt%start
 
  220    ib = ib - 1 + st%group%iblock(st%st_start)
 
  223      write(message(1), 
'(a,i2,a,i2,a,i2)') 
'Task ', st%st_kpt_mpi_grp%rank, 
' will send the block ', &
 
  224        ib, 
' with k-point ', ik
 
  225      call messages_info(1, all_nodes=.
true.)
 
  234    type(states_elec_t),                             
intent(in)  :: st
 
  235    integer,                                         
intent(in)  :: icom
 
  236    integer,                                         
intent(out) :: ib
 
  237    integer,                                         
intent(out) :: ik
 
  242    ib = mod(icom-1, this%nblock_to_receive)+1
 
  243    ik = int((icom-ib)/this%nblock_to_receive) +  st%st_kpt_task(this%task_from, 3)
 
  244    ib = ib - 1 + st%group%iblock(st%st_kpt_task(this%task_from, 1))
 
  247      write(message(1), 
'(a,i2,a,i2,a,i2)') 
'Task ', st%st_kpt_mpi_grp%rank, 
' will receive the block ', &
 
  248        ib, 
' with k-point ', ik
 
  249      call messages_info(1, all_nodes=.
true.)
 
This module implements batches of mesh functions.
 
type(debug_t), save, public debug
 
This module defines the meshes, which are used in Octopus.
 
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
 
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
 
This module provides routines for communicating all batches in a ring-pattern scheme.
 
subroutine states_elec_all_to_all_communications_get_receive_indices(this, st, icom, ib, ik)
Given the icom step, returns the block and k-point indices to be received.
 
integer function states_elec_all_to_all_communications_eval_nsend(st, task_to, nblock_to_send)
How many batches we will send from task_send.
 
integer pure function states_elec_all_to_all_communications_get_nreceive(this)
Returns the number of receiv calls.
 
subroutine states_elec_all_to_all_communications_alloc_receive_batch(this, st, icom, np, psib)
Given the icom step, allocate the receiv buffer (wfs_elec_t)
 
integer pure function states_elec_all_to_all_communications_get_ncom(this)
Returns the number of communications.
 
subroutine states_elec_all_to_all_communications_get_send_indices(this, st, icom, ib, ik)
Given the icom step, returns the block and k-point indices to be sent.
 
integer function states_elec_all_to_all_communications_eval_nreceive(st, task_from, nblock_to_receive)
How many batches we will receive from task_from.
 
subroutine states_elec_all_to_all_communications_start(this, st, task_from, task_to)
Given a task to send to, and a task to receive from, initializes a states_elec_all_to_all_communicati...
 
integer pure function states_elec_all_to_all_communications_get_nsend(this)
Returns the number send calls.
 
This module provides routines for communicating states when using states parallelization.