m_rankMerge.F90 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620
  1. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  2. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  3. !-----------------------------------------------------------------------
  4. ! CVS m_rankMerge.F90,v 1.3 2004-04-21 22:54:48 jacob Exp
  5. ! CVS MCT_2_8_0
  6. !BOP -------------------------------------------------------------------
  7. !
  8. ! !MODULE: m_rankMerge - A merging tool through ranking
  9. !
  10. ! !DESCRIPTION:
  11. !
  12. ! !INTERFACE:
  13. module m_rankMerge
  14. implicit none
  15. private ! except
  16. public :: rankSet ! set inital ranks
  17. public :: rankMerge ! merge two ranks
  18. public :: IndexedRankMerge ! index-merge two array segments
  19. interface rankSet; module procedure set_; end interface
  20. interface rankMerge; module procedure &
  21. imerge_, & ! rank-merging two integer arrays
  22. rmerge_, & ! rank-merging two real arrays
  23. dmerge_, & ! rank-merging two dble arrays
  24. uniq_ ! merging to rank arrays
  25. end interface
  26. interface IndexedRankMerge; module procedure &
  27. iindexmerge_, & ! merging two index arrays of integers
  28. rindexmerge_, & ! merging two index arrays of reals
  29. dindexmerge_ ! merging two index arrays of dbles
  30. end interface
  31. ! !REVISION HISTORY:
  32. ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
  33. ! - initial prototype/prolog/code
  34. !EOP ___________________________________________________________________
  35. character(len=*),parameter :: myname='MCT(MPEU)::m_rankMerge'
  36. contains
  37. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  38. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  39. !BOP -------------------------------------------------------------------
  40. !
  41. ! !IROUTINE: set_ - set initial ranking
  42. !
  43. ! !DESCRIPTION:
  44. !
  45. ! !INTERFACE:
  46. subroutine set_(rank)
  47. implicit none
  48. integer,dimension(:),intent(out) :: rank
  49. ! !REVISION HISTORY:
  50. ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
  51. ! - initial prototype/prolog/code
  52. !EOP ___________________________________________________________________
  53. character(len=*),parameter :: myname_=myname//'::set_'
  54. integer :: i
  55. do i=1,size(rank)
  56. rank(i)=0
  57. end do
  58. end subroutine set_
  59. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  60. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  61. !BOP -------------------------------------------------------------------
  62. !
  63. ! !IROUTINE: imerge_ - merge two sorted integer arrays by ranking
  64. !
  65. ! !DESCRIPTION:
  66. !
  67. ! !INTERFACE:
  68. subroutine imerge_(value_i,value_j,krank_i,krank_j,descend)
  69. implicit none
  70. integer,dimension(:),intent(in) :: value_j ! value of j-vec
  71. integer,dimension(:),intent(in) :: value_i ! value of i-vec
  72. integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
  73. integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
  74. logical,optional,intent(in) :: descend
  75. ! !REVISION HISTORY:
  76. ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
  77. ! - initial prototype/prolog/code
  78. !EOP ___________________________________________________________________
  79. character(len=*),parameter :: myname_=myname//'::imerge_'
  80. integer :: ni,nj
  81. logical :: descend_
  82. logical :: geti
  83. integer :: value_sv,value
  84. integer :: krank
  85. integer :: i,j
  86. descend_=.false.
  87. if(present(descend)) descend_=descend
  88. ni=size(krank_i)
  89. nj=size(krank_j)
  90. i=1
  91. j=1
  92. krank=0 ! a preset rank value
  93. value_sv=0
  94. do
  95. geti=j>nj
  96. if(geti) then ! .eqv. j>nj
  97. if(i>ni) exit ! i>ni
  98. value = value_i(i)
  99. else ! .eqv. j<=nj
  100. geti = i<=ni
  101. if(geti) then ! .eqv. i<=ni
  102. value = value_i(i)
  103. geti = krank_i(i) <= krank_j(j)
  104. if(krank_i(i)==krank_j(j)) then
  105. geti = value_i(i)<=value_j(j)
  106. if(descend_) geti = value_i(i)>=value_j(j)
  107. endif
  108. endif
  109. if(.not.geti) value = value_j(j)
  110. endif
  111. if(krank==0 .or. value /= value_sv) then
  112. krank=krank+1 ! the next rank value
  113. value_sv=value
  114. endif
  115. if(geti) then
  116. krank_i(i)=krank
  117. i=i+1
  118. else
  119. krank_j(j)=krank
  120. j=j+1
  121. endif
  122. end do
  123. end subroutine imerge_
  124. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  125. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  126. !BOP -------------------------------------------------------------------
  127. !
  128. ! !IROUTINE: rmerge_ - merge two sorted real arrays by ranking
  129. !
  130. ! !DESCRIPTION:
  131. !
  132. ! !INTERFACE:
  133. subroutine rmerge_(value_i,value_j,krank_i,krank_j,descend)
  134. use m_realkinds, only : SP
  135. implicit none
  136. real(SP),dimension(:),intent(in) :: value_i ! value of i-vec
  137. real(SP),dimension(:),intent(in) :: value_j ! value of j-vec
  138. integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
  139. integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
  140. logical,optional,intent(in) :: descend
  141. ! !REVISION HISTORY:
  142. ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
  143. ! - initial prototype/prolog/code
  144. !EOP ___________________________________________________________________
  145. character(len=*),parameter :: myname_=myname//'::rmerge_'
  146. integer :: ni,nj
  147. logical :: descend_
  148. logical :: geti
  149. real(SP) :: value_sv,value
  150. integer :: krank
  151. integer :: i,j
  152. descend_=.false.
  153. if(present(descend)) descend_=descend
  154. ni=size(krank_i)
  155. nj=size(krank_j)
  156. i=1
  157. j=1
  158. krank=0 ! a preset rank value
  159. value_sv=0
  160. do
  161. geti=j>nj
  162. if(geti) then ! .eqv. j>nj
  163. if(i>ni) exit ! i>ni
  164. value = value_i(i)
  165. else ! .eqv. j<=nj
  166. geti = i<=ni
  167. if(geti) then ! .eqv. i<=ni
  168. value = value_i(i)
  169. geti = krank_i(i) <= krank_j(j)
  170. if(krank_i(i)==krank_j(j)) then
  171. geti = value_i(i)<=value_j(j)
  172. if(descend_) geti = value_i(i)>=value_j(j)
  173. endif
  174. endif
  175. if(.not.geti) value = value_j(j)
  176. endif
  177. if(krank==0 .or. value /= value_sv) then
  178. krank=krank+1 ! the next rank value
  179. value_sv=value
  180. endif
  181. if(geti) then
  182. krank_i(i)=krank
  183. i=i+1
  184. else
  185. krank_j(j)=krank
  186. j=j+1
  187. endif
  188. end do
  189. end subroutine rmerge_
  190. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  191. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  192. !BOP -------------------------------------------------------------------
  193. !
  194. ! !IROUTINE: dmerge_ - merge two sorted real arrays by ranking
  195. !
  196. ! !DESCRIPTION:
  197. !
  198. ! !INTERFACE:
  199. subroutine dmerge_(value_i,value_j,krank_i,krank_j,descend)
  200. use m_realkinds, only : DP
  201. implicit none
  202. real(DP),dimension(:),intent(in) :: value_i ! value of i-vec
  203. real(DP),dimension(:),intent(in) :: value_j ! value of j-vec
  204. integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
  205. integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
  206. logical,optional,intent(in) :: descend
  207. ! !REVISION HISTORY:
  208. ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
  209. ! - initial prototype/prolog/code
  210. !EOP ___________________________________________________________________
  211. character(len=*),parameter :: myname_=myname//'::dmerge_'
  212. integer :: ni,nj
  213. logical :: descend_
  214. logical :: geti
  215. real(DP):: value_sv,value
  216. integer :: krank
  217. integer :: i,j
  218. descend_=.false.
  219. if(present(descend)) descend_=descend
  220. ni=size(krank_i)
  221. nj=size(krank_j)
  222. i=1
  223. j=1
  224. krank=0 ! a preset rank value
  225. value_sv=0
  226. do
  227. geti=j>nj
  228. if(geti) then ! .eqv. j>nj
  229. if(i>ni) exit ! i>ni
  230. value = value_i(i)
  231. else ! .eqv. j<=nj
  232. geti = i<=ni
  233. if(geti) then ! .eqv. i<=ni
  234. value = value_i(i)
  235. geti = krank_i(i) <= krank_j(j)
  236. if(krank_i(i)==krank_j(j)) then
  237. geti = value_i(i)<=value_j(j)
  238. if(descend_) geti = value_i(i)>=value_j(j)
  239. endif
  240. endif
  241. if(.not.geti) value = value_j(j)
  242. endif
  243. if(krank==0 .or. value /= value_sv) then
  244. krank=krank+1 ! the next rank value
  245. value_sv=value
  246. endif
  247. if(geti) then
  248. krank_i(i)=krank
  249. i=i+1
  250. else
  251. krank_j(j)=krank
  252. j=j+1
  253. endif
  254. end do
  255. end subroutine dmerge_
  256. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  257. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  258. !BOP -------------------------------------------------------------------
  259. !
  260. ! !IROUTINE: iindexmerge_ - merge two sorted integer arrays by ranking
  261. !
  262. ! !DESCRIPTION:
  263. !
  264. ! !INTERFACE:
  265. subroutine iindexmerge_(indx_i,indx_j,value,krank_i,krank_j,descend)
  266. implicit none
  267. integer,dimension(:),intent(in) :: indx_i ! of the i-vec
  268. integer,dimension(:),intent(in) :: indx_j ! of the j-vec
  269. integer,dimension(:),intent(in) :: value ! of the full
  270. integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
  271. integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
  272. logical,optional,intent(in) :: descend
  273. ! !REVISION HISTORY:
  274. ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
  275. ! - initial prototype/prolog/code
  276. !EOP ___________________________________________________________________
  277. character(len=*),parameter :: myname_=myname//'::iindexmerge_'
  278. integer :: ni,nj
  279. logical :: descend_
  280. logical :: geti
  281. integer :: value_sv,value_
  282. integer :: krank
  283. integer :: i,j,li,lj
  284. descend_=.false.
  285. if(present(descend)) descend_=descend
  286. ni=size(krank_i)
  287. nj=size(krank_j)
  288. i=1
  289. j=1
  290. krank=0 ! a preset rank value
  291. value_sv=0
  292. do
  293. geti=j>nj
  294. if(geti) then ! .eqv. j>nj
  295. if(i>ni) exit ! i>ni
  296. li=indx_i(i)
  297. value_ = value(li)
  298. else ! .eqv. j<=nj
  299. lj=indx_j(j)
  300. geti = i<=ni
  301. if(geti) then ! .eqv. i<=ni
  302. li=indx_i(i)
  303. value_ = value(li)
  304. geti = krank_i(i) <= krank_j(j)
  305. if(krank_i(i)==krank_j(j)) then
  306. geti = value(li)<=value(lj)
  307. if(descend_) geti = value(li)>=value(lj)
  308. endif
  309. endif
  310. if(.not.geti) value_ = value(lj)
  311. endif
  312. if(krank==0 .or. value_ /= value_sv) then
  313. krank=krank+1 ! the next rank value
  314. value_sv=value_
  315. endif
  316. if(geti) then
  317. krank_i(i)=krank
  318. i=i+1
  319. else
  320. krank_j(j)=krank
  321. j=j+1
  322. endif
  323. end do
  324. end subroutine iindexmerge_
  325. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  326. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  327. !BOP -------------------------------------------------------------------
  328. !
  329. ! !IROUTINE: rindexmerge_ - merge two sorted real arrays by ranking
  330. !
  331. ! !DESCRIPTION:
  332. !
  333. ! !INTERFACE:
  334. subroutine rindexmerge_(indx_i,indx_j,value,krank_i,krank_j,descend)
  335. use m_realkinds,only : SP
  336. implicit none
  337. integer,dimension(:),intent(in) :: indx_i ! of the i-vec
  338. integer,dimension(:),intent(in) :: indx_j ! of the j-vec
  339. real(SP),dimension(:),intent(in) :: value ! of the full
  340. integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
  341. integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
  342. logical,optional,intent(in) :: descend
  343. ! !REVISION HISTORY:
  344. ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
  345. ! - initial prototype/prolog/code
  346. !EOP ___________________________________________________________________
  347. character(len=*),parameter :: myname_=myname//'::rindexmerge_'
  348. integer :: ni,nj
  349. logical :: descend_
  350. logical :: geti
  351. real(SP):: value_sv,value_
  352. integer :: krank
  353. integer :: i,j,li,lj
  354. descend_=.false.
  355. if(present(descend)) descend_=descend
  356. ni=size(krank_i)
  357. nj=size(krank_j)
  358. i=1
  359. j=1
  360. krank=0 ! a preset rank value
  361. value_sv=0
  362. do
  363. geti=j>nj
  364. if(geti) then ! .eqv. j>nj
  365. if(i>ni) exit ! i>ni
  366. li=indx_i(i)
  367. value_ = value(li)
  368. else ! .eqv. j<=nj
  369. lj=indx_j(j)
  370. geti = i<=ni
  371. if(geti) then ! .eqv. i<=ni
  372. li=indx_i(i)
  373. value_ = value(li)
  374. geti = krank_i(i) <= krank_j(j)
  375. if(krank_i(i)==krank_j(j)) then
  376. geti = value(li)<=value(lj)
  377. if(descend_) geti = value(li)>=value(lj)
  378. endif
  379. endif
  380. if(.not.geti) value_ = value(lj)
  381. endif
  382. if(krank==0 .or. value_ /= value_sv) then
  383. krank=krank+1 ! the next rank value
  384. value_sv=value_
  385. endif
  386. if(geti) then
  387. krank_i(i)=krank
  388. i=i+1
  389. else
  390. krank_j(j)=krank
  391. j=j+1
  392. endif
  393. end do
  394. end subroutine rindexmerge_
  395. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  396. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  397. !BOP -------------------------------------------------------------------
  398. !
  399. ! !IROUTINE: dindexmerge_ - merge two sorted real arrays by ranking
  400. !
  401. ! !DESCRIPTION:
  402. !
  403. ! !INTERFACE:
  404. subroutine dindexmerge_(indx_i,indx_j,value,krank_i,krank_j,descend)
  405. use m_realkinds,only : DP
  406. implicit none
  407. integer,dimension(:),intent(in) :: indx_i ! of the i-vec
  408. integer,dimension(:),intent(in) :: indx_j ! of the j-vec
  409. real(DP),dimension(:),intent(in) :: value ! of the full
  410. integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
  411. integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
  412. logical,optional,intent(in) :: descend
  413. ! !REVISION HISTORY:
  414. ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
  415. ! - initial prototype/prolog/code
  416. !EOP ___________________________________________________________________
  417. character(len=*),parameter :: myname_=myname//'::dindexmerge_'
  418. integer :: ni,nj
  419. logical :: descend_
  420. logical :: geti
  421. real(DP):: value_sv,value_
  422. integer :: krank
  423. integer :: i,j,li,lj
  424. descend_=.false.
  425. if(present(descend)) descend_=descend
  426. ni=size(krank_i)
  427. nj=size(krank_j)
  428. i=1
  429. j=1
  430. krank=0 ! a preset rank value
  431. value_sv=0
  432. do
  433. geti=j>nj
  434. if(geti) then ! .eqv. j>nj
  435. if(i>ni) exit ! i>ni
  436. li=indx_i(i)
  437. value_ = value(li)
  438. else ! .eqv. j<=nj
  439. lj=indx_j(j)
  440. geti = i<=ni
  441. if(geti) then ! .eqv. i<=ni
  442. li=indx_i(i)
  443. value_ = value(li)
  444. geti = krank_i(i) <= krank_j(j)
  445. if(krank_i(i)==krank_j(j)) then
  446. geti = value(li)<=value(lj)
  447. if(descend_) geti = value(li)>=value(lj)
  448. endif
  449. endif
  450. if(.not.geti) value_ = value(lj)
  451. endif
  452. if(krank==0 .or. value_ /= value_sv) then
  453. krank=krank+1 ! the next rank value
  454. value_sv=value_
  455. endif
  456. if(geti) then
  457. krank_i(i)=krank
  458. i=i+1
  459. else
  460. krank_j(j)=krank
  461. j=j+1
  462. endif
  463. end do
  464. end subroutine dindexmerge_
  465. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  466. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  467. !BOP -------------------------------------------------------------------
  468. !
  469. ! !IROUTINE: uniq_ - merge two rank arrays with unique rank values
  470. !
  471. ! !DESCRIPTION:
  472. !
  473. ! !INTERFACE:
  474. subroutine uniq_(krank_i,krank_j)
  475. implicit none
  476. integer,dimension(:),intent(inout) :: krank_i ! rank of i-vec
  477. integer,dimension(:),intent(inout) :: krank_j ! rank of j-vec
  478. ! !REVISION HISTORY:
  479. ! 13Mar00 - Jing Guo <guo@dao.gsfc.nasa.gov>
  480. ! - initial prototype/prolog/code
  481. !EOP ___________________________________________________________________
  482. character(len=*),parameter :: myname_=myname//'::uniq_'
  483. integer :: ni,nj
  484. integer :: i,j
  485. integer :: krank
  486. logical :: geti
  487. ni=size(krank_i)
  488. nj=size(krank_j)
  489. i=1
  490. j=1
  491. krank=0
  492. do
  493. geti=j>nj
  494. if(geti) then ! .eqv. j>nj
  495. if(i>ni) exit ! i>ni
  496. else ! .eqv. j<=nj
  497. geti = i<=ni
  498. if(geti) geti = krank_i(i) <= krank_j(j) ! if(i<=ni) ..
  499. endif
  500. krank=krank+1 ! the next rank value
  501. if(geti) then
  502. krank_i(i)=krank
  503. i=i+1
  504. else
  505. krank_j(j)=krank
  506. j=j+1
  507. endif
  508. end do
  509. end subroutine uniq_
  510. end module m_rankMerge