123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620 |
- !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
- ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
- !-----------------------------------------------------------------------
- ! CVS m_rankMerge.F90,v 1.3 2004-04-21 22:54:48 jacob Exp
- ! CVS MCT_2_8_0
- !BOP -------------------------------------------------------------------
- !
- ! !MODULE: m_rankMerge - A merging tool through ranking
- !
- ! !DESCRIPTION:
- !
- ! !INTERFACE:
- module m_rankMerge
- implicit none
- private ! except
- public :: rankSet ! set inital ranks
- public :: rankMerge ! merge two ranks
- public :: IndexedRankMerge ! index-merge two array segments
- interface rankSet; module procedure set_; end interface
- interface rankMerge; module procedure &
- imerge_, & ! rank-merging two integer arrays
- rmerge_, & ! rank-merging two real arrays
- dmerge_, & ! rank-merging two dble arrays
- uniq_ ! merging to rank arrays
- end interface
- interface IndexedRankMerge; module procedure &
- iindexmerge_, & ! merging two index arrays of integers
- rindexmerge_, & ! merging two index arrays of reals
- dindexmerge_ ! merging two index arrays of dbles
- end interface
- ! !REVISION HISTORY:
- ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
- ! - initial prototype/prolog/code
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname='MCT(MPEU)::m_rankMerge'
- contains
- !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
- ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: set_ - set initial ranking
- !
- ! !DESCRIPTION:
- !
- ! !INTERFACE:
- subroutine set_(rank)
- implicit none
- integer,dimension(:),intent(out) :: rank
- ! !REVISION HISTORY:
- ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
- ! - initial prototype/prolog/code
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::set_'
- integer :: i
- do i=1,size(rank)
- rank(i)=0
- end do
- end subroutine set_
- !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
- ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: imerge_ - merge two sorted integer arrays by ranking
- !
- ! !DESCRIPTION:
- !
- ! !INTERFACE:
- subroutine imerge_(value_i,value_j,krank_i,krank_j,descend)
- implicit none
- integer,dimension(:),intent(in) :: value_j ! value of j-vec
- integer,dimension(:),intent(in) :: value_i ! value of i-vec
- integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
- integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
- logical,optional,intent(in) :: descend
- ! !REVISION HISTORY:
- ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
- ! - initial prototype/prolog/code
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::imerge_'
- integer :: ni,nj
- logical :: descend_
- logical :: geti
- integer :: value_sv,value
- integer :: krank
- integer :: i,j
- descend_=.false.
- if(present(descend)) descend_=descend
-
- ni=size(krank_i)
- nj=size(krank_j)
-
- i=1
- j=1
- krank=0 ! a preset rank value
- value_sv=0
- do
- geti=j>nj
- if(geti) then ! .eqv. j>nj
- if(i>ni) exit ! i>ni
- value = value_i(i)
- else ! .eqv. j<=nj
- geti = i<=ni
- if(geti) then ! .eqv. i<=ni
- value = value_i(i)
- geti = krank_i(i) <= krank_j(j)
- if(krank_i(i)==krank_j(j)) then
- geti = value_i(i)<=value_j(j)
- if(descend_) geti = value_i(i)>=value_j(j)
- endif
- endif
- if(.not.geti) value = value_j(j)
- endif
- if(krank==0 .or. value /= value_sv) then
- krank=krank+1 ! the next rank value
- value_sv=value
- endif
-
- if(geti) then
- krank_i(i)=krank
- i=i+1
- else
- krank_j(j)=krank
- j=j+1
- endif
- end do
- end subroutine imerge_
- !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
- ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: rmerge_ - merge two sorted real arrays by ranking
- !
- ! !DESCRIPTION:
- !
- ! !INTERFACE:
- subroutine rmerge_(value_i,value_j,krank_i,krank_j,descend)
- use m_realkinds, only : SP
- implicit none
- real(SP),dimension(:),intent(in) :: value_i ! value of i-vec
- real(SP),dimension(:),intent(in) :: value_j ! value of j-vec
- integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
- integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
- logical,optional,intent(in) :: descend
- ! !REVISION HISTORY:
- ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
- ! - initial prototype/prolog/code
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::rmerge_'
- integer :: ni,nj
- logical :: descend_
- logical :: geti
- real(SP) :: value_sv,value
- integer :: krank
- integer :: i,j
-
- descend_=.false.
- if(present(descend)) descend_=descend
- ni=size(krank_i)
- nj=size(krank_j)
-
- i=1
- j=1
- krank=0 ! a preset rank value
- value_sv=0
- do
- geti=j>nj
- if(geti) then ! .eqv. j>nj
- if(i>ni) exit ! i>ni
- value = value_i(i)
- else ! .eqv. j<=nj
- geti = i<=ni
- if(geti) then ! .eqv. i<=ni
- value = value_i(i)
- geti = krank_i(i) <= krank_j(j)
- if(krank_i(i)==krank_j(j)) then
- geti = value_i(i)<=value_j(j)
- if(descend_) geti = value_i(i)>=value_j(j)
- endif
- endif
- if(.not.geti) value = value_j(j)
- endif
-
- if(krank==0 .or. value /= value_sv) then
- krank=krank+1 ! the next rank value
- value_sv=value
- endif
-
- if(geti) then
- krank_i(i)=krank
- i=i+1
- else
- krank_j(j)=krank
- j=j+1
- endif
- end do
- end subroutine rmerge_
- !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
- ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: dmerge_ - merge two sorted real arrays by ranking
- !
- ! !DESCRIPTION:
- !
- ! !INTERFACE:
- subroutine dmerge_(value_i,value_j,krank_i,krank_j,descend)
- use m_realkinds, only : DP
- implicit none
- real(DP),dimension(:),intent(in) :: value_i ! value of i-vec
- real(DP),dimension(:),intent(in) :: value_j ! value of j-vec
- integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
- integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
- logical,optional,intent(in) :: descend
- ! !REVISION HISTORY:
- ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
- ! - initial prototype/prolog/code
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::dmerge_'
- integer :: ni,nj
- logical :: descend_
- logical :: geti
- real(DP):: value_sv,value
- integer :: krank
- integer :: i,j
-
- descend_=.false.
- if(present(descend)) descend_=descend
- ni=size(krank_i)
- nj=size(krank_j)
-
- i=1
- j=1
- krank=0 ! a preset rank value
- value_sv=0
- do
- geti=j>nj
- if(geti) then ! .eqv. j>nj
- if(i>ni) exit ! i>ni
- value = value_i(i)
- else ! .eqv. j<=nj
- geti = i<=ni
- if(geti) then ! .eqv. i<=ni
- value = value_i(i)
- geti = krank_i(i) <= krank_j(j)
- if(krank_i(i)==krank_j(j)) then
- geti = value_i(i)<=value_j(j)
- if(descend_) geti = value_i(i)>=value_j(j)
- endif
- endif
- if(.not.geti) value = value_j(j)
- endif
-
- if(krank==0 .or. value /= value_sv) then
- krank=krank+1 ! the next rank value
- value_sv=value
- endif
-
- if(geti) then
- krank_i(i)=krank
- i=i+1
- else
- krank_j(j)=krank
- j=j+1
- endif
- end do
- end subroutine dmerge_
- !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
- ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: iindexmerge_ - merge two sorted integer arrays by ranking
- !
- ! !DESCRIPTION:
- !
- ! !INTERFACE:
- subroutine iindexmerge_(indx_i,indx_j,value,krank_i,krank_j,descend)
- implicit none
- integer,dimension(:),intent(in) :: indx_i ! of the i-vec
- integer,dimension(:),intent(in) :: indx_j ! of the j-vec
- integer,dimension(:),intent(in) :: value ! of the full
- integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
- integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
- logical,optional,intent(in) :: descend
- ! !REVISION HISTORY:
- ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
- ! - initial prototype/prolog/code
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::iindexmerge_'
- integer :: ni,nj
- logical :: descend_
- logical :: geti
- integer :: value_sv,value_
- integer :: krank
- integer :: i,j,li,lj
- descend_=.false.
- if(present(descend)) descend_=descend
-
- ni=size(krank_i)
- nj=size(krank_j)
-
- i=1
- j=1
- krank=0 ! a preset rank value
- value_sv=0
- do
- geti=j>nj
- if(geti) then ! .eqv. j>nj
- if(i>ni) exit ! i>ni
- li=indx_i(i)
- value_ = value(li)
- else ! .eqv. j<=nj
- lj=indx_j(j)
- geti = i<=ni
- if(geti) then ! .eqv. i<=ni
- li=indx_i(i)
- value_ = value(li)
- geti = krank_i(i) <= krank_j(j)
- if(krank_i(i)==krank_j(j)) then
- geti = value(li)<=value(lj)
- if(descend_) geti = value(li)>=value(lj)
- endif
- endif
- if(.not.geti) value_ = value(lj)
- endif
- if(krank==0 .or. value_ /= value_sv) then
- krank=krank+1 ! the next rank value
- value_sv=value_
- endif
-
- if(geti) then
- krank_i(i)=krank
- i=i+1
- else
- krank_j(j)=krank
- j=j+1
- endif
- end do
- end subroutine iindexmerge_
- !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
- ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: rindexmerge_ - merge two sorted real arrays by ranking
- !
- ! !DESCRIPTION:
- !
- ! !INTERFACE:
- subroutine rindexmerge_(indx_i,indx_j,value,krank_i,krank_j,descend)
- use m_realkinds,only : SP
- implicit none
- integer,dimension(:),intent(in) :: indx_i ! of the i-vec
- integer,dimension(:),intent(in) :: indx_j ! of the j-vec
- real(SP),dimension(:),intent(in) :: value ! of the full
- integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
- integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
- logical,optional,intent(in) :: descend
- ! !REVISION HISTORY:
- ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
- ! - initial prototype/prolog/code
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::rindexmerge_'
- integer :: ni,nj
- logical :: descend_
- logical :: geti
- real(SP):: value_sv,value_
- integer :: krank
- integer :: i,j,li,lj
-
- descend_=.false.
- if(present(descend)) descend_=descend
- ni=size(krank_i)
- nj=size(krank_j)
-
- i=1
- j=1
- krank=0 ! a preset rank value
- value_sv=0
- do
- geti=j>nj
- if(geti) then ! .eqv. j>nj
- if(i>ni) exit ! i>ni
- li=indx_i(i)
- value_ = value(li)
- else ! .eqv. j<=nj
- lj=indx_j(j)
- geti = i<=ni
- if(geti) then ! .eqv. i<=ni
- li=indx_i(i)
- value_ = value(li)
- geti = krank_i(i) <= krank_j(j)
- if(krank_i(i)==krank_j(j)) then
- geti = value(li)<=value(lj)
- if(descend_) geti = value(li)>=value(lj)
- endif
- endif
- if(.not.geti) value_ = value(lj)
- endif
-
- if(krank==0 .or. value_ /= value_sv) then
- krank=krank+1 ! the next rank value
- value_sv=value_
- endif
-
- if(geti) then
- krank_i(i)=krank
- i=i+1
- else
- krank_j(j)=krank
- j=j+1
- endif
- end do
- end subroutine rindexmerge_
- !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
- ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: dindexmerge_ - merge two sorted real arrays by ranking
- !
- ! !DESCRIPTION:
- !
- ! !INTERFACE:
- subroutine dindexmerge_(indx_i,indx_j,value,krank_i,krank_j,descend)
- use m_realkinds,only : DP
- implicit none
- integer,dimension(:),intent(in) :: indx_i ! of the i-vec
- integer,dimension(:),intent(in) :: indx_j ! of the j-vec
- real(DP),dimension(:),intent(in) :: value ! of the full
- integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
- integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
- logical,optional,intent(in) :: descend
- ! !REVISION HISTORY:
- ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
- ! - initial prototype/prolog/code
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::dindexmerge_'
- integer :: ni,nj
- logical :: descend_
- logical :: geti
- real(DP):: value_sv,value_
- integer :: krank
- integer :: i,j,li,lj
-
- descend_=.false.
- if(present(descend)) descend_=descend
- ni=size(krank_i)
- nj=size(krank_j)
-
- i=1
- j=1
- krank=0 ! a preset rank value
- value_sv=0
- do
- geti=j>nj
- if(geti) then ! .eqv. j>nj
- if(i>ni) exit ! i>ni
- li=indx_i(i)
- value_ = value(li)
- else ! .eqv. j<=nj
- lj=indx_j(j)
- geti = i<=ni
- if(geti) then ! .eqv. i<=ni
- li=indx_i(i)
- value_ = value(li)
- geti = krank_i(i) <= krank_j(j)
- if(krank_i(i)==krank_j(j)) then
- geti = value(li)<=value(lj)
- if(descend_) geti = value(li)>=value(lj)
- endif
- endif
- if(.not.geti) value_ = value(lj)
- endif
-
- if(krank==0 .or. value_ /= value_sv) then
- krank=krank+1 ! the next rank value
- value_sv=value_
- endif
-
- if(geti) then
- krank_i(i)=krank
- i=i+1
- else
- krank_j(j)=krank
- j=j+1
- endif
- end do
- end subroutine dindexmerge_
- !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
- ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
- !BOP -------------------------------------------------------------------
- !
- ! !IROUTINE: uniq_ - merge two rank arrays with unique rank values
- !
- ! !DESCRIPTION:
- !
- ! !INTERFACE:
- subroutine uniq_(krank_i,krank_j)
- implicit none
- integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
- integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
- ! !REVISION HISTORY:
- ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
- ! - initial prototype/prolog/code
- !EOP ___________________________________________________________________
- character(len=*),parameter :: myname_=myname//'::uniq_'
- integer :: ni,nj
- integer :: i,j
- integer :: krank
- logical :: geti
- ni=size(krank_i)
- nj=size(krank_j)
- i=1
- j=1
- krank=0
- do
- geti=j>nj
- if(geti) then ! .eqv. j>nj
- if(i>ni) exit ! i>ni
- else ! .eqv. j<=nj
- geti = i<=ni
- if(geti) geti = krank_i(i) <= krank_j(j) ! if(i<=ni) ..
- endif
-
- krank=krank+1 ! the next rank value
- if(geti) then
- krank_i(i)=krank
- i=i+1
- else
- krank_j(j)=krank
- j=j+1
- endif
- end do
- end subroutine uniq_
- end module m_rankMerge
|