fm_reduce_an_head_send_broadcast Subroutine

private subroutine fm_reduce_an_head_send_broadcast(self, antot, antota, antotp)

Type Bound

fieldmat_type

Arguments

Type IntentOptional Attributes Name
class(fieldmat_type), intent(in) :: self
complex, intent(inout), dimension (-ntgrid:,:,:) :: antot
complex, intent(inout), dimension (-ntgrid:,:,:) :: antota
complex, intent(inout), dimension (-ntgrid:,:,:) :: antotp

Contents


Source Code

  subroutine fm_reduce_an_head_send_broadcast (self,antot, antota, antotp)
    use mp, only: broadcast, iproc, nbsend, nbrecv, initialise_requests, waitall, broadcast_sub
    use theta_grid, only: ntgrid
    use run_parameters, only: has_phi, has_apar, has_bpar
    use gs2_layouts, only: gf_lo, proc_id, idx, ik_idx, it_idx
    use kt_grids, only: kwork_filter
    implicit none
    class(fieldmat_type), intent(in) :: self
    complex, dimension (-ntgrid:,:,:), intent (in out) :: antot, antota, antotp

    integer :: igf, ik, it, ic, is, root, sendreqnum, recvreqnum, sender, receiver, tag, taglim
    integer, dimension (gf_lo%xypoints*3) :: send_handles
    !AJ This array is bigger than it needs to be but does mean we do not
    !AJ need to calculate how many supercells a give process is head for
    integer, dimension (gf_lo%naky*gf_lo%ntheta0*3) :: recv_handles
    integer :: tempit, tempis, tempismax
    complex, dimension(:,:,:,:,:), allocatable :: tempdata

    tempis = 0
    tempismax = 0
    do ik = 1,self%naky
       tempis = self%kyb(ik)%nsupercell
       if(tempis > tempismax) then
          tempismax = tempis
       end if
       
    end do

    allocate(tempdata(-ntgrid:ntgrid,gf_lo%ntheta0,self%nfield,gf_lo%naky,tempismax))

    call initialise_requests(send_handles)
    call initialise_requests(recv_handles)

    recvreqnum = 1

    taglim = self%naky*self%ntheta0+1

    !AJ Iterate through the supercells and each supercell that I am a head for 
    !AJ post a receive for data for a given element of the field array.
    do ik = 1,self%naky
       do is=1,self%kyb(ik)%nsupercell
          if(self%kyb(ik)%supercells(is)%is_head) then
             tempit = 1
             do ic=1,self%kyb(ik)%supercells(is)%ncell

                it = self%kyb(ik)%supercells(is)%cells(ic)%it_ind

                if(kwork_filter(it,ik)) cycle

                !AJ Calculate which process owns this it,ik point in gf_lo space
                !AJ i.e. which process will be sending the data to this supercell head
                sender = proc_id(gf_lo,idx(gf_lo,ik,it))

                if(sender /= iproc) then
                
                   tag = ik*self%ntheta0 + it
                
                   if(has_phi) then
                      call nbrecv(tempdata(:,tempit,1,ik,is),sender,tag,recv_handles(recvreqnum))
                      recvreqnum = recvreqnum + 1
                   end if

                   tag = tag + taglim
                   
                   if(has_apar) then
                      call nbrecv(tempdata(:,tempit,2,ik,is),sender,tag,recv_handles(recvreqnum))
                      recvreqnum = recvreqnum + 1
                   end if
                   
                   tag = tag + taglim
                   
                   if(has_bpar) then
                      call nbrecv(tempdata(:,tempit,3,ik,is),sender,tag,recv_handles(recvreqnum))
                      recvreqnum = recvreqnum + 1
                   end if

                else

                   if(has_phi) then
                      tempdata(:,tempit,1,ik,is) = antot(:,it,ik)
                   end if

                   if(has_apar) then
                      tempdata(:,tempit,2,ik,is) = antota(:,it,ik)
                   end if

                   if(has_bpar) then
                      tempdata(:,tempit,3,ik,is) = antotp(:,it,ik)
                   end if

                end if

                tempit = tempit + 1

             end do
          end if
       end do
    end do

    sendreqnum = 1

    !AJ Iterate through all the points I own in gf_lo and send that part of the field 
    !AJ arrays to the process that is the supercell head for the supercell that contains 
    !AJ that it, ik point.
    do igf = gf_lo%llim_proc, gf_lo%ulim_proc
       ik = ik_idx(gf_lo,igf)
       it = it_idx(gf_lo,igf)

       if(kwork_filter(it,ik)) cycle

       tag = ik*self%ntheta0 + it

       receiver = self%heads(it,ik)

       if(iproc /= receiver) then

          if(has_phi) then
             call nbsend(antot(:,it,ik),receiver,tag,send_handles(sendreqnum))
             sendreqnum = sendreqnum + 1
          end if
          
          tag = tag + taglim
          
          if(has_apar) then
             call nbsend(antota(:,it,ik),receiver,tag,send_handles(sendreqnum))
             sendreqnum = sendreqnum + 1
          end if
          
          tag = tag + taglim
          
          if(has_bpar) then
             call nbsend(antotp(:,it,ik),receiver,tag,send_handles(sendreqnum))
             sendreqnum = sendreqnum + 1
          end if
          
       end if
       
    end do

    call waitall(sendreqnum-1,send_handles)
    call waitall(recvreqnum-1,recv_handles)

    do ik = 1,self%naky
       if(.not.self%kyb(ik)%is_local) cycle
       do is=1,self%kyb(ik)%nsupercell
          if(.not. self%kyb(ik)%supercells(is)%is_local) cycle
          root = self%kyb(ik)%supercells(is)%head_iproc_pd
          
          if((self%kyb(ik)%supercells(is)%is_empty.or.self%kyb(ik)%supercells(is)%is_all_local)) cycle

          tempit = 1
          do ic=1,self%kyb(ik)%supercells(is)%ncell
             it = self%kyb(ik)%supercells(is)%cells(ic)%it_ind

             if(kwork_filter(it,ik)) cycle

             tempit = tempit + 1

          end do

          call broadcast_sub(tempdata(:,1:tempit,1:self%nfield,ik,is),root,self%kyb(ik)%supercells(is)%sc_sub_pd%id)

          if(has_phi) then
             tempit = 1
             do ic=1,self%kyb(ik)%supercells(is)%ncell
                it = self%kyb(ik)%supercells(is)%cells(ic)%it_ind

                if(kwork_filter(it,ik)) cycle

                antot(:,it,ik) = tempdata(:,tempit,1,ik,is)
                tempit = tempit + 1
             end do
          end if
             
          if(has_apar) then
             tempit = 1
             do ic=1,self%kyb(ik)%supercells(is)%ncell
                it = self%kyb(ik)%supercells(is)%cells(ic)%it_ind

                if(kwork_filter(it,ik)) cycle

                antota(:,it,ik) = tempdata(:,tempit,2,ik,is)
                tempit = tempit + 1
             end do
          end if
          
          if(has_bpar) then
             tempit = 1
             do ic=1,self%kyb(ik)%supercells(is)%ncell
                it = self%kyb(ik)%supercells(is)%cells(ic)%it_ind

                if(kwork_filter(it,ik)) cycle

                antotp(:,it,ik) = tempdata(:,tempit,3,ik,is)
                tempit = tempit + 1
             end do
          end if
             
       end do
    end do


    deallocate(tempdata)

  end subroutine fm_reduce_an_head_send_broadcast