!-------------------------------------- LICENCE BEGIN ------------------------------------
!Environment Canada - Atmospheric Science and Technology License/Disclaimer,
!                     version 3; Last Modified: May 7, 2008.
!This is free but copyrighted software; you can use/redistribute/modify it under the terms
!of the Environment Canada - Atmospheric Science and Technology License/Disclaimer
!version 3 or (at your option) any later version that should be found at:
!http://collaboration.cmc.ec.gc.ca/science/rpn.comm/license.html
!
!This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
!without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
!See the above mentioned License/Disclaimer for more details.
!You should have received a copy of the License/Disclaimer along with this software;
!if not, you can write to: EC-RPN COMM Group, 2121 TransCanada, suite 500, Dorval (Quebec),
!CANADA, H9P 1J3; or send e-mail to service.rpn@ec.gc.ca
!-------------------------------------- LICENCE END --------------------------------------
!
!*****************************************************************!! 
!  Author:  Bin He    *ARMA/MSC    Aug. 2005
! purposes:  
!  Create some utility functions to facilitate the MPI parallelization. 
!
!*****************************************************************!! 
!!!===================================================

 SUBROUTINE para_range1(n1,n2,mprocs,ibeg,iend)
 IMPLICIT NONE
 INTEGER,INTENT(IN) :: n1  
 INTEGER,INTENT(IN) :: n2  
 INTEGER,INTENT(IN) :: mprocs 
 INTEGER,DIMENSION(0:mprocs-1),INTENT(INOUT) :: ibeg 
 INTEGER,DIMENSION(0:mprocs-1),INTENT(INOUT) :: iend 
 INTEGER :: irank,size,iw1,iw2
 INTEGER,DIMENSION(mprocs) :: xpos
 size=n2-n1+1
 iw1 = size/mprocs
 iw2 = MOD(size,mprocs)
 xpos(1)=n1+iw1
 DO irank=2,mprocs-1
    xpos(irank)=iw1+xpos(irank-1)
    IF(irank<=iw2) xpos(irank)=xpos(irank)+1
 ENDDO
 ibeg(0)=n1
 iend(mprocs-1)=n2
 DO irank=1,mprocs-1
    ibeg(irank)=xpos(irank)
    iend(irank-1)=xpos(irank)
 ENDDO
 END SUBROUTINE para_range1
!
!=============================================================

SUBROUTINE tile_corner(myLDlat,myLDlon,myRUplat,myRUplon,npex,npey,myidx,myidy)  
 IMPLICIT NONE
 INTEGER,INTENT(IN) :: npex,npey,myidx,myidy
 INTEGER,INTENT(OUT) :: myLDlat,myLDlon,myRuplat,myRuplon

 !
 REAL*8 :: Dlon,Dlat

 IF(npex == 0 .or. npey == 0) stop "Error : npex or npey Zero" 
 Dlon=360.0/float(npex) 
 Dlat=180.0/float(npey) 
 myLDlon=INT(Dlon*myidx*100)
 if(myidx == 0) myLDlon=0
 myLDlat=INT(Dlat*myidy*100)
 if(myidy == 0) myLDlat=0

 myRUplon=INT(Dlon*(myidx+1)*100)
 if(myidx == (npex-1)) myRuplon=36000
 myRuplat=INT(DLAT*(myidy+1)*100)
 if(myidy == (npey-1)) myRUplat=18000
 print*,'myLDlon,myLDlat,myRupLon,myRuplat= ',myLDlon,myLDlat,myRupLon,myRuplat
END SUBROUTINE tile_corner
!=============================================================

!==================================================================

 SUBROUTINE para_disp(disp,mprocs,disp1,count) 1
    IMPLICIT NONE
    INTEGER,INTENT(IN) :: mprocs,disp1
    INTEGER,DIMENSION(mprocs),INTENT(IN) ::  count
    INTEGER,DIMENSION(mprocs),INTENT(OUT) :: disp
    INTEGER :: i
    disp(1)=disp1
    DO i=2,mprocs
       disp(i)=disp(i-1)+count(i-1)
    END DO
  END SUBROUTINE para_disp
!!!===================================================

 INTEGER FUNCTION get_loc_spsize(ntrbeg,ntrend,mxzw)
   IMPLICIT NONE
   INTEGER :: ntrbeg,ntrend
   INTEGER :: mxzw
   INTEGER :: jm,jm0,ins,jns,jn,count
   count=0
   DO jm0=ntrbeg,ntrend
      ins=1
      IF(jm0 == (mxzw-jm0)) ins=0
      DO jns=0,ins
         jm=(1-jns)*jm0+jns*(mxzw-jm0)
         DO jn=jm,mxzw
            count=count+1
         ENDDO
      ENDDO
   ENDDO
!!   count=count
   get_loc_spsize=count
 END FUNCTION get_loc_spsize
!!!===================================================

 SUBROUTINE get_loc_sppos(ntrbeg,ntrend,mxzw,nind)
   IMPLICIT NONE
   INTEGER,INTENT(IN) :: ntrbeg,ntrend
   INTEGER,INTENT(IN) :: mxzw
   INTEGER,DIMENSION(0:mxzw),INTENT(OUT) :: nind
   INTEGER :: pos,len,jm0,ntrbeg1,ntrend1
!
    nind=0
    nind(ntrbeg)=1
    jm0=ntrbeg
    DO jm0=ntrbeg+1,ntrend
       len=mxzw-jm0+2
       nind(jm0)=nind(jm0-1) + len
    ENDDO
    len=mxzw-ntrend+1
    ntrbeg1=mxzw-ntrend
    ntrend1=mxzw-ntrbeg
    IF(ntrbeg1 == ntrend) ntrbeg1=ntrbeg1+1
    nind(ntrbeg1) = nind(ntrend) + len
    DO jm0=ntrbeg1+1,ntrend1
       len=mxzw-jm0+2
       nind(jm0)=nind(jm0-1) + len
    ENDDO
 END SUBROUTINE get_loc_sppos
!=======================================================================

 SUBROUTINE get_disploc1(disploc,nlaloc,ntrbeg,ntrend,nind,nindloc,mxzw)
  IMPLICIT NONE
  INTEGER,INTENT(IN) :: nlaloc,ntrbeg,ntrend, mxzw
  INTEGER,DIMENSION(0:mxzw),INTENT(IN) :: nind 
  INTEGER,DIMENSION(0:mxzw),INTENT(IN) :: nindloc 
  INTEGER,DIMENSION(nlaloc),INTENT(OUT) :: disploc 
! define local variables  
   INTEGER :: jm,jm0,ins,jns,jn,indx,indx1 
   DO jm0=ntrbeg,ntrend
      ins=1
      IF(jm0 == (mxzw-jm0)) ins=0
      DO jns=0,ins
         jm=(1-jns)*jm0+jns*(mxzw-jm0)
         DO jn=jm,mxzw
            if(nind(jm) <  mxzw ) then 
              indx=nind(jm) + (jn -jm)  !! no imaginary part for the spectral coeff.  
            else
              indx=2*nind(jm) + 2*(jn -jm)-(mxzw+2) !! imaginary part of the spectra l coeff.  
             !! indx1=2*nindloc(jm)+2*(jn -jm)-(mxzw+2)!!imaginary part of the spectra l coeff.  
            endif  
             indx1=nindloc(jm) + (jn -jm)
             disploc(indx1)=indx   
         ENDDO
      ENDDO
   ENDDO
 END SUBROUTINE get_disploc1
!=======================================================================

 SUBROUTINE get_disploc0(disploc0,nlaloc,ntrbeg,ntrend,nind,nindloc,mxzw)
  IMPLICIT NONE
  INTEGER,INTENT(IN) :: nlaloc,ntrbeg,ntrend, mxzw
  INTEGER,DIMENSION(0:mxzw),INTENT(IN) :: nind,nindloc 
  INTEGER,DIMENSION(nlaloc),INTENT(OUT) :: disploc0 
! define local variables  
  INTEGER :: jm,jm0,ins,jns,jn,indx,indx1 
  disploc0=0
  DO jm0=ntrbeg,ntrend
     ins=1
     IF(jm0 == (mxzw-jm0)) ins=0
     DO jns=0,ins
        jm=(1-jns)*jm0+jns*(mxzw-jm0)
        DO jn=jm,mxzw
           indx1=nindloc(jm) +(jn-jm)
           indx=nind(jm) + (jn -jm)    
           disploc0(indx1)=indx   
        ENDDO
     ENDDO
  ENDDO
 END SUBROUTINE get_disploc0
!=======================================================================

 SUBROUTINE creat_pgdtype(imin,imax,ista,iend,jsta,jend,ksta,kend,&
                             maxzw,ioldtype,inewtype,lglb)
!=================================================================
!    3 Dimensional Case.
!=================================================================
 IMPLICIT NONE
 INTEGER,INTENT(IN) :: imin,imax
 INTEGER,INTENT(IN) :: ista,iend
 INTEGER,INTENT(IN) :: jsta,jend
 INTEGER,INTENT(IN) :: ksta,kend
 INTEGER,INTENT(IN) :: ioldtype
 INTEGER,INTENT(IN) :: maxzw
 INTEGER,INTENT(OUT) :: inewtype
 LOGICAL,INTENT(IN) :: lglb  
! Declair local variables. 
 INTEGER :: isize,ierr,itemp,itemp1,ilen,jlen,klen,ilenmx,rsd
 INTEGER :: k,j,count,ii0,izw,jm,ins,jns,MXZW,ilonr,iloni,jr,ji
 INTEGER :: jjsta,jjend
 INTEGER ,ALLOCATABLE,DIMENSION(:) :: disp,blklen
 ilen = iend -ista +1
 jlen = jend -jsta +1
 klen = kend -ksta +1
 ilenmx=imax - imin+1
 ALLOCATE(disp(4*ilenmx*jlen*klen+2))
 ALLOCATE(blklen(4*ilenmx*jlen*klen+2))
 MXZW=maxzw
 ii0=0
 rsd=1-imin
 IF(lglb) THEN
    jjsta=jsta
    jjend=jend
 ELSE
    jjsta=1
    jjend=jlen
 ENDIF
! construct a derived data type to the 3 D zoonal wave coefficients . 
  DO j=jjsta,jjend
     DO k=ksta,kend
        DO izw=ista,iend
           ins=1
           IF(izw == (MXZW - izw))ins=0
            DO jns=0,ins
               jm=(1-jns)*izw + jns*(MXZW-izw)
               ilonr=2*jm
               ilonr=2*jm +  rsd
               iloni=ilonr+1
               jr=2*ii0+1
               ji=jr+1
               disp(jr)=ilenmx*klen*(j-1)+ilenmx*(k-1)+ilonr
               disp(ji)=disp(jr)+1
               blklen(jr)=1
               blklen(ji)=1
               ii0=ii0+1
            ENDDO
         ENDDO
      ENDDO 
  ENDDO  
  count=ji
  CALL MPI_TYPE_INDEXED(count,blklen,disp,ioldtype,inewtype,ierr)
  CALL MPi_TYPE_COMMIT(inewtype,ierr)
!
  DEALLOCATE(disp,blklen)
  END SUBROUTINE creat_pgdtype
!=================================================================
!

 SUBROUTINE zwave_rtranspose(pgd,ibeg,iend,jbeglc,jendlc,kdim)
 IMPLICIT NONE
 INTEGER,INTENT(IN) :: ibeg,iend,kdim
 INTEGER,INTENT(IN) :: jbeglc,jendlc
 REAL*8,DIMENSION(ibeg:iend,kdim,jbeglc:jendlc),INTENT(INOUT) :: pgd
!
 INTEGER :: ierr,size,lderive
 lderive=1
 SELECT CASE (lderive)
 CASE (1)
    size= (iend-ibeg+1)*kdim*(jendlc-jbeglc+1)
    CALL rpn_comm_allreduce(pgd,pgd,size,"mpi_double_precision","mpi_sum","EW",ierr)
 END SELECT  
 END SUBROUTINE zwave_rtranspose
!======================================================================

 SUBROUTINE get_glbgd(gdin,gdout,nibeg,niend,nkgdim,njbeg,njend,ibeglc,&,2
     iendlc,jbeglc,jendlc)
 USE procs_topo
 IMPLICIT NONE
!!
 INTEGER , INTENT(IN)  :: nibeg,niend
 INTEGER , INTENT(IN)  :: nkgdim
 INTEGER , INTENT(IN)  :: njbeg,njend
 INTEGER , INTENT(IN)  :: ibeglc,iendlc
 INTEGER , INTENT(IN)  :: jbeglc,jendlc
 REAL*8,DIMENSION(ibeglc:iendlc,nkgdim,jbeglc:jendlc),INTENT(IN) :: gdin
 REAL*8,DIMENSION(nibeg:niend,nkgdim,njbeg:njend),INTENT(OUT) :: gdout
 INTEGER,DIMENSION(npey) :: count,displs
!! --------------------------------------------------------------------
 INTEGER :: jendloce,size ,ier
!----------------------------------------------------------------------
 jendloce=jendlc
!!! Note:  since the last latitude may be overlapped with the starting latitute of the 
!!!       next process in the same communication group (comm_col), it is necessary to
!!    get rid of the last latitude for the process who is not the last in the colume. 
!! 
 if(myidy /= npey-1) jendloce=jendlc-1
 size=(niend-nibeg+1)*(nkgdim)*(jendloce-jbeglc+1)
 CALL rpn_comm_allgather(size,1,"mpi_integer",count,1,"mpi_integer","NS",ier)
 CALL para_disp(displs,npey,0,count)
 CALL rpn_comm_gatherv(gdin,size,"mpi_double_precision",gdout,count,displs,"mpi_double_precision",0,"NS",ier)
 size=(niend-nibeg+1)*(nkgdim)*(njend-njbeg+1)
 CALL rpn_comm_bcast(gdout,size,"mpi_double_precision",0,"NS",ier)
 END SUBROUTINE get_glbgd
!----------------------------------------------------------------------

 LOGICAL FUNCTION range_check(jm,ntrbeg,ntrend,ntrunc)
 IMPLICIT NONE
 INTEGER,INTENT(IN) :: jm
 INTEGER,INTENT(IN) :: ntrbeg,ntrend
 INTEGER,INTENT(IN) :: ntrunc
 IF( (jm <= ntrend .AND. jm >= ntrbeg) .or. (jm>= (Ntrunc-ntrend)&
     .and. jm <= (ntrunc-ntrbeg) )) THEN
    range_check = .true.
 ELSE
    range_check = .false.
 ENDIF
 END FUNCTION range_check
!----------------------------------------------------------------------

 SUBROUTINE  Abortmpi(infor,comm)
 CHARACTER (len = *),INTENT(IN) :: infor
 INTEGER,INTENT(IN) :: comm
 integer :: errorcode
 errorcode=1
 write(*,*) infor
 CALL MPI_ABORT(comm,errorcode)
 END SUBROUTINE Abortmpi
!!-------------------------------------------------------------------

 SUBROUTINE mpiabort(infor) 1
 CHARACTER (len = *),INTENT(IN) :: infor
 write(*,*) infor
 stop 
 END SUBROUTINE MPIABORT 


 SUBROUTINE data_distr2d1_rpncomm(gdg,nibeg,niend,njbeg,njend,nkbeg,nkend,gdloc,&
                 ibeg,iend,jbeg,jend,kbeg,kend)
  IMPLICIT NONE
  INTEGER,INTENT(IN) :: nibeg,niend  ! Global index
  INTEGER,INTENT(IN) :: ibeg,iend    ! local index
  INTEGER,INTENT(IN) :: njbeg,njend
  INTEGER,INTENT(IN) :: jbeg,jend
  INTEGER,INTENT(IN) :: nkbeg,nkend
  INTEGER,INTENT(IN) :: kbeg,kend
  REAL*8,DIMENSION(nibeg:niend,nkbeg:nkend,njbeg:njend),INTENT(IN) :: gdg
  REAL*8,DIMENSION(ibeg:iend,kbeg:kend,jbeg:jend),INTENT(OUT) :: gdloc
!
!   Declair the local variables.
  INTEGER :: gibeg,giend,gjbeg,gjend,nig,njg,libeg,liend,ljbeg,ljend 
  INTEGER :: ierr
!
!
  gibeg=nibeg
  giend=(niend-nibeg+1)*(nkend-nkbeg+1)+nibeg  
  gjbeg=njbeg
  gjend=njend 

  libeg=gibeg
  liend=giend 

  ljbeg=1
  ljend=jend-jbeg+1 

  nig=(niend-nibeg+1)*(nkend-nkbeg+1) 
  njg=(njend-njbeg+1) 

  CALL rpn_comm_dist(gdg,gibeg,giend,gjbeg,gjend,nig,njg,1,0,0,2,gdloc,&
                    libeg,liend,ljbeg, ljend,0,1,.false.,.false.,ierr)    

 END SUBROUTINE data_distr2d1_rpncomm
!!==========================================================

  SUBROUTINE data_distr2d1(gdg,nibeg,niend,njbeg,njend,nkbeg,nkend,gdloc,&,1
             ibeg,iend,jbeg,jend,kbeg,kend)
  USE procs_topo 
  IMPLICIT NONE
  INTEGER,INTENT(IN) :: nibeg,niend  ! Global index
  INTEGER,INTENT(IN) :: ibeg,iend    ! local index
  INTEGER,INTENT(IN) :: njbeg,njend
  INTEGER,INTENT(IN) :: jbeg,jend
  INTEGER,INTENT(IN) :: nkbeg,nkend
  INTEGER,INTENT(IN) :: kbeg,kend
  REAL*8,DIMENSION(nibeg:niend,nkbeg:nkend,njbeg:njend),INTENT(IN) :: gdg
  REAL*8,DIMENSION(ibeg:iend,kbeg:kend,jbeg:jend),INTENT(OUT) :: gdloc
!   Declair the local variables.
  INTEGER :: ierr
  INTEGER,DIMENSION(npey) :: count,disp

  INTEGER :: newtype,ierr,sizeloc ,root

!
  INTEGER :: i,j,jj,backsize
!
  root=0
  sizeloc=(iend-ibeg+1)*(jend-jbeg+1)*(kend-kbeg+1)
 
  CALL rpn_comm_ALLGATHER(sizeloc,1,"mpi_integer",count,1,"mpi_integer","NS",ierr)
  backsize=(iend-ibeg+1)*(kend-kbeg+1)
  disp(1)=0
  DO j=2,npey
     disp(j)=disp(j-1) + count(j-1) -backsize
  ENDDO
  IF(myidx == 0 ) &
      CALL rpn_comm_SCATTERV(gdg,count,disp,"mpi_double_precision",gdloc,sizeloc,&
           "mpi_double_precision",root,"NS",ierr)
  CALL rpn_comm_BARRIER("grid",ierr) 
  IF(npex >1 ) &
     CALL rpn_comm_BCAST(gdloc,sizeloc,"mpi_double_precision",0,"EW",ierr)

  print*,'======END of data_distr2d1========='
  END SUBROUTINE data_distr2d1
!!==========================================================