TABLE OF CONTENTS


ABINIT/prep_index_wavef_bandpp [ Functions ]

[ Top ] [ Functions ]

NAME

 prep_index_wavef_bandpp

FUNCTION

 this routine sorts the waves functions by bandpp and by processors 
 after the alltoall

COPYRIGHT

INPUTS

  nproc_band = number of processors below the band
  bandpp     = number of groups of couple of waves functions
  nspinor    = number of spin
  ndatarecv  = total number of values received by the processor and sended 
               by the other processors band
  recvcounts = number of values sended by each processor band and received 
               by the processor
  rdispls    = positions of the values received by the processor and  
               sended by each processor band

OUTPUT

  index_wavef_band = position of the sorted values

SIDE EFFECTS

PARENTS

      chebfi,prep_fourwf,prep_getghc,prep_nonlop

CHILDREN

SOURCE

 35 #if defined HAVE_CONFIG_H
 36 #include "config.h"
 37 #endif
 38 
 39 #include "abi_common.h"
 40 
 41 subroutine prep_index_wavef_bandpp(nproc_band,bandpp,&
 42                              nspinor,ndatarecv,&
 43                              recvcounts,rdispls,&
 44                              index_wavef_band)
 45 
 46  use defs_basis
 47  use m_profiling_abi
 48 
 49 !This section has been created automatically by the script Abilint (TD).
 50 !Do not modify the following lines by hand.
 51 #undef ABI_FUNC
 52 #define ABI_FUNC 'prep_index_wavef_bandpp'
 53 !End of the abilint section
 54 
 55  implicit none
 56 
 57 !Arguments ------------------------------------
 58 !scalars
 59  integer,intent(in) :: bandpp,ndatarecv,nproc_band,nspinor
 60 !arrays
 61  integer,intent(in) :: rdispls(nproc_band),recvcounts(nproc_band)
 62  integer,allocatable,intent(out) :: index_wavef_band(:)
 63 
 64 !Local variables-------------------------------
 65 !scalars
 66  integer :: delta,idebc,idebe,ifinc,ifine,iindex,iproc,kbandpp,nb
 67 
 68 ! *********************************************************************
 69 
 70 !DEBUG
 71 !write(std_out,*)' prep_index_wavef_banpp : enter '
 72 !write(std_out,*) 'ndatarecv = ', ndatarecv
 73 !write(std_out,*) 'rdispls(:) = ', rdispls(:) 
 74 !write(std_out,*) 'recvcounts(:) = ', recvcounts(:) 
 75 !ENDDEBUG
 76 
 77 
 78 !---------------------------------------------
 79 !Allocation
 80 !---------------------------------------------
 81  ABI_ALLOCATE(index_wavef_band ,(bandpp*nspinor*ndatarecv))
 82  index_wavef_band(:)   =0
 83  
 84 !---------------------------------------------
 85 !Calcul : loops on bandpp and processors band
 86 !---------------------------------------------
 87  nb = sum(recvcounts(1:nproc_band))
 88  do kbandpp=1,bandpp
 89    
 90    do iproc=1,nproc_band
 91      
 92      idebe = (rdispls(iproc) + 1)  + (kbandpp-1) * ndatarecv*nspinor
 93      ifine = idebe + recvcounts(iproc) -1
 94 
 95      if (iproc==1) then
 96        idebc =   (kbandpp-1)* recvcounts(iproc)*nspinor + 1  
 97      else
 98        idebc = (bandpp)  * sum(recvcounts(1:iproc-1))*nspinor &
 99        + (kbandpp-1)* recvcounts(iproc)*nspinor &
100        + 1 
101      end if
102      ifinc = idebc + recvcounts(iproc) -1
103      index_wavef_band(idebe:ifine) = (/( iindex,iindex=idebc,ifinc)/)
104      delta=ifine-idebe
105      if (nspinor==2) then
106        index_wavef_band(idebe+nb :idebe+nb +delta)=(/( iindex,iindex=ifinc+1,ifinc+1+delta)/)
107      end if
108    end do
109  end do
110 
111 end subroutine prep_index_wavef_bandpp