m_MergeSorts.F90 28 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195
  1. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  2. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  3. !-----------------------------------------------------------------------
  4. ! CVS m_MergeSorts.F90,v 1.3 2004-04-21 22:54:45 jacob Exp
  5. ! CVS MCT_2_8_0
  6. !BOP -------------------------------------------------------------------
  7. !
  8. ! !MODULE: m_MergeSorts - Tools for incremental indexed-sorting
  9. !
  10. ! !DESCRIPTION:
  11. !
  12. ! This tool module contains basic sorting procedures, that in
  13. ! addition to a couple of standard Fortran 90 statements in the
  14. ! array syntex, allow a full range sort or unsort operations.
  15. ! The main characteristics of the sorting algorithm used in this
  16. ! module are, a) stable, and b) index sorting.
  17. !
  18. ! !INTERFACE:
  19. module m_MergeSorts
  20. implicit none
  21. private ! except
  22. public :: IndexSet
  23. public :: IndexSort
  24. interface IndexSet
  25. module procedure setn_
  26. module procedure set_
  27. end interface
  28. interface IndexSort
  29. module procedure iSortn_
  30. module procedure rSortn_
  31. module procedure dSortn_
  32. module procedure cSortn_
  33. module procedure iSort_
  34. module procedure rSort_
  35. module procedure dSort_
  36. module procedure cSort_
  37. module procedure iSort1_
  38. module procedure rSort1_
  39. module procedure dSort1_
  40. module procedure cSort1_
  41. end interface
  42. ! !EXAMPLES:
  43. !
  44. ! ...
  45. ! integer, intent(in) :: No
  46. ! type(Observations), dimension(No), intent(inout) :: obs
  47. !
  48. ! integer, dimension(No) :: indx ! automatic array
  49. !
  50. ! call IndexSet(No,indx)
  51. ! call IndexSort(No,indx,obs(1:No)%lev,descend=.false.)
  52. ! call IndexSort(No,indx,obs(1:No)%lon,descend=.false.)
  53. ! call IndexSort(No,indx,obs(1:No)%lat,descend=.false.)
  54. ! call IndexSort(No,indx,obs(1:No)%kt,descend=.false.)
  55. ! call IndexSort(No,indx,obs(1:No)%ks,descend=.false.)
  56. ! call IndexSort(No,indx,obs(1:No)%kx,descend=.false.)
  57. ! call IndexSort(No,indx,obs(1:No)%kr,descend=.false.)
  58. !
  59. ! ! Sorting
  60. ! obs(1:No) = obs( (/ (indx(i),i=1,No) /) )
  61. ! ...
  62. ! ! Unsorting
  63. ! obs( (/ (indx(i),i=1,No) /) ) = obs(1:No)
  64. !
  65. ! !REVISION HISTORY:
  66. ! 15Mar00 - Jing Guo
  67. ! . Added interfaces without the explicit size
  68. ! . Added interfaces for two dimensional arrays
  69. ! 02Feb99 - Jing Guo <guo@thunder> - Added if(present(stat)) ...
  70. ! 04Jan99 - Jing Guo <guo@thunder> - revised
  71. ! 09Sep97 - Jing Guo <guo@thunder> - initial prototype/prolog/code
  72. !EOP ___________________________________________________________________
  73. character(len=*), parameter :: myname='MCT(MPEU)::m_MergeSorts'
  74. contains
  75. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  76. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  77. !BOP -------------------------------------------------------------------
  78. !
  79. ! !IROUTINE: setn_ - Initialize an array of data location indices
  80. !
  81. ! !DESCRIPTION:
  82. !
  83. ! !INTERFACE:
  84. subroutine setn_(n,indx)
  85. implicit none
  86. integer, intent(in) :: n ! size of indx(:)
  87. integer, dimension(n), intent(out) :: indx ! indices
  88. ! !REVISION HISTORY:
  89. ! 15Mar00 - Jing Guo
  90. ! . initial prototype/prolog/code
  91. ! . redefined for the original interface
  92. !EOP ___________________________________________________________________
  93. call set_(indx(1:n))
  94. end subroutine setn_
  95. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  96. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  97. !BOP -------------------------------------------------------------------
  98. !
  99. ! !IROUTINE: set_ - Initialize an array of data location indices
  100. !
  101. ! !DESCRIPTION:
  102. !
  103. ! !INTERFACE:
  104. subroutine set_(indx)
  105. implicit none
  106. integer, dimension(:), intent(out) :: indx ! indices
  107. ! !REVISION HISTORY:
  108. ! 15Mar00 - Jing Guo
  109. ! . Modified the interface, by removing the explicit size
  110. ! 09Sep97 - Jing Guo <guo@thunder> - initial prototype/prolog/code
  111. ! 04Jan99 - Jing Guo <guo@thunder> - revised prolog format
  112. !EOP ___________________________________________________________________
  113. integer :: i
  114. do i=1,size(indx)
  115. indx(i)=i
  116. end do
  117. end subroutine set_
  118. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  119. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  120. !BOP -------------------------------------------------------------------
  121. !
  122. ! !IROUTINE: iSortn_ - A stable merge index sorting of INTs.
  123. !
  124. ! !DESCRIPTION:
  125. !
  126. ! !INTERFACE:
  127. subroutine iSortn_(n,indx,keys,descend,stat)
  128. implicit none
  129. integer,intent(in) :: n
  130. integer, dimension(n), intent(inout) :: indx
  131. integer, dimension(n), intent(in) :: keys
  132. logical, optional, intent(in) :: descend
  133. integer, optional, intent(out) :: stat
  134. ! !REVISION HISTORY:
  135. ! 15Mar00 - Jing Guo
  136. ! . initial prototype/prolog/code
  137. ! . redefined for the original interface
  138. !EOP ___________________________________________________________________
  139. character(len=*),parameter :: myname_=myname//'::iSortn_'
  140. call iSort_(indx(1:n),keys(1:n),descend,stat)
  141. end subroutine iSortn_
  142. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  143. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  144. !BOP -------------------------------------------------------------------
  145. !
  146. ! !IROUTINE: rSortn_ - A stable merge index sorting REALs.
  147. !
  148. ! !DESCRIPTION:
  149. !
  150. ! !INTERFACE:
  151. subroutine rSortn_(n,indx,keys,descend,stat)
  152. use m_realkinds,only : SP
  153. implicit none
  154. integer,intent(in) :: n
  155. integer, dimension(n), intent(inout) :: indx
  156. real(SP),dimension(n), intent(in) :: keys
  157. logical, optional, intent(in) :: descend
  158. integer, optional, intent(out) :: stat
  159. ! !REVISION HISTORY:
  160. ! 15Mar00 - Jing Guo
  161. ! . initial prototype/prolog/code
  162. ! . redefined for the original interface
  163. !EOP ___________________________________________________________________
  164. character(len=*),parameter :: myname_=myname//'::rSortn_'
  165. call rSort_(indx(1:n),keys(1:n),descend,stat)
  166. end subroutine rSortn_
  167. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  168. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  169. !BOP -------------------------------------------------------------------
  170. !
  171. ! !IROUTINE: dSortn_ - A stable merge index sorting DOUBLEs.
  172. !
  173. ! !DESCRIPTION:
  174. !
  175. ! !INTERFACE:
  176. subroutine dSortn_(n,indx,keys,descend,stat)
  177. use m_realkinds,only : DP
  178. implicit none
  179. integer,intent(in) :: n
  180. integer, dimension(n), intent(inout) :: indx
  181. real(DP), dimension(n), intent(in) :: keys
  182. logical, optional, intent(in) :: descend
  183. integer, optional, intent(out) :: stat
  184. ! !REVISION HISTORY:
  185. ! 15Mar00 - Jing Guo
  186. ! . initial prototype/prolog/code
  187. ! . redefined for the original interface
  188. !EOP ___________________________________________________________________
  189. character(len=*),parameter :: myname_=myname//'::dSortn_'
  190. call dSort_(indx(1:n),keys(1:n),descend,stat)
  191. end subroutine dSortn_
  192. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  193. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  194. !BOP -------------------------------------------------------------------
  195. !
  196. ! !IROUTINE: cSortn_ - A stable merge index sorting of CHAR(*)s.
  197. !
  198. ! !DESCRIPTION:
  199. !
  200. ! !INTERFACE:
  201. subroutine cSortn_(n,indx,keys,descend,stat)
  202. implicit none
  203. integer,intent(in) :: n
  204. integer, dimension(n), intent(inout) :: indx
  205. character(len=*), dimension(n), intent(in) :: keys
  206. logical, optional, intent(in) :: descend
  207. integer, optional, intent(out) :: stat
  208. ! !REVISION HISTORY:
  209. ! 15Mar00 - Jing Guo
  210. ! . initial prototype/prolog/code
  211. ! . redefined for the original interface
  212. !EOP ___________________________________________________________________
  213. character(len=*),parameter :: myname_=myname//'::cSortn_'
  214. call cSort_(indx(1:n),keys(1:n),descend,stat)
  215. end subroutine cSortn_
  216. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  217. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  218. !BOP -------------------------------------------------------------------
  219. !
  220. ! !IROUTINE: iSort_ - A stable merge index sorting of INTs.
  221. !
  222. ! !DESCRIPTION:
  223. !
  224. ! !INTERFACE:
  225. subroutine iSort_(indx,keys,descend,stat)
  226. use m_stdio, only : stderr
  227. use m_die, only : die
  228. implicit none
  229. integer, dimension(:), intent(inout) :: indx
  230. integer, dimension(:), intent(in) :: keys
  231. logical, optional, intent(in) :: descend
  232. integer, optional, intent(out) :: stat
  233. ! !REVISION HISTORY:
  234. ! 15Mar00 - Jing Guo
  235. ! . Modified the interface, by removing the explicit size
  236. ! 02Feb99 - Jing Guo <guo@thunder> - Added if(present(stat)) ...
  237. ! 04Jan99 - Jing Guo <guo@thunder> - revised the prolog
  238. ! 09Sep97 - Jing Guo <guo@thunder> - initial prototype/prolog/code
  239. !EOP ___________________________________________________________________
  240. logical :: dsnd
  241. integer :: ierr
  242. integer, dimension(:),allocatable :: mtmp
  243. integer :: n
  244. character(len=*),parameter :: myname_=myname//'::iSort_'
  245. if(present(stat)) stat=0
  246. n=size(indx)
  247. allocate(mtmp(n),stat=ierr)
  248. if(ierr /= 0) then
  249. write(stderr,'(2a,i4)') myname_, &
  250. ': allocate(mtmp(:)) error, stat =',ierr
  251. if(.not.present(stat)) call die(myname_)
  252. stat=ierr
  253. return
  254. endif
  255. dsnd=.false.
  256. if(present(descend)) dsnd=descend
  257. call MergeSort_()
  258. deallocate(mtmp)
  259. contains
  260. subroutine MergeSort_()
  261. implicit none
  262. integer :: mstep,lstep
  263. integer :: lb,lm,le
  264. mstep=1
  265. do while(mstep < n)
  266. lstep=mstep*2
  267. lb=1
  268. do while(lb < n)
  269. lm=lb+mstep
  270. le=min(lm-1+mstep,n)
  271. call merge_(lb,lm,le)
  272. indx(lb:le)=mtmp(lb:le)
  273. lb=le+1
  274. end do
  275. mstep=lstep
  276. end do
  277. end subroutine MergeSort_
  278. subroutine merge_(lb,lm,le)
  279. integer,intent(in) :: lb,lm,le
  280. integer :: l1,l2,l
  281. l1=lb
  282. l2=lm
  283. do l=lb,le
  284. if(l2.gt.le) then
  285. mtmp(l)=indx(l1)
  286. l1=l1+1
  287. elseif(l1.ge.lm) then
  288. mtmp(l)=indx(l2)
  289. l2=l2+1
  290. else
  291. if(dsnd) then
  292. if(keys(indx(l1)) .ge. keys(indx(l2))) then
  293. mtmp(l)=indx(l1)
  294. l1=l1+1
  295. else
  296. mtmp(l)=indx(l2)
  297. l2=l2+1
  298. endif
  299. else
  300. if(keys(indx(l1)) .le. keys(indx(l2))) then
  301. mtmp(l)=indx(l1)
  302. l1=l1+1
  303. else
  304. mtmp(l)=indx(l2)
  305. l2=l2+1
  306. endif
  307. endif
  308. endif
  309. end do
  310. end subroutine merge_
  311. end subroutine iSort_
  312. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  313. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  314. !BOP -------------------------------------------------------------------
  315. !
  316. ! !IROUTINE: rSort_ - A stable merge index sorting REALs.
  317. !
  318. ! !DESCRIPTION:
  319. !
  320. ! !INTERFACE:
  321. subroutine rSort_(indx,keys,descend,stat)
  322. use m_stdio, only : stderr
  323. use m_die, only : die
  324. use m_realkinds,only : SP
  325. implicit none
  326. integer, dimension(:), intent(inout) :: indx
  327. real(SP),dimension(:), intent(in) :: keys
  328. logical, optional, intent(in) :: descend
  329. integer, optional, intent(out) :: stat
  330. ! !REVISION HISTORY:
  331. ! 15Mar00 - Jing Guo
  332. ! . Modified the interface, by removing the explicit size
  333. ! 02Feb99 - Jing Guo <guo@thunder> - Added if(present(stat)) ...
  334. ! 04Jan99 - Jing Guo <guo@thunder> - revised the prolog
  335. ! 09Sep97 - Jing Guo <guo@thunder> - initial prototype/prolog/code
  336. !EOP ___________________________________________________________________
  337. logical :: dsnd
  338. integer :: ierr
  339. integer, dimension(:),allocatable :: mtmp
  340. integer :: n
  341. character(len=*),parameter :: myname_=myname//'::rSort_'
  342. if(present(stat)) stat=0
  343. n=size(indx)
  344. allocate(mtmp(n),stat=ierr)
  345. if(ierr /= 0) then
  346. write(stderr,'(2a,i4)') myname_, &
  347. ': allocate(mtmp(:)) error, stat =',ierr
  348. if(.not.present(stat)) call die(myname_)
  349. stat=ierr
  350. return
  351. endif
  352. dsnd=.false.
  353. if(present(descend)) dsnd=descend
  354. call MergeSort_()
  355. deallocate(mtmp)
  356. contains
  357. subroutine MergeSort_()
  358. implicit none
  359. integer :: mstep,lstep
  360. integer :: lb,lm,le
  361. mstep=1
  362. do while(mstep < n)
  363. lstep=mstep*2
  364. lb=1
  365. do while(lb < n)
  366. lm=lb+mstep
  367. le=min(lm-1+mstep,n)
  368. call merge_(lb,lm,le)
  369. indx(lb:le)=mtmp(lb:le)
  370. lb=le+1
  371. end do
  372. mstep=lstep
  373. end do
  374. end subroutine MergeSort_
  375. subroutine merge_(lb,lm,le)
  376. integer,intent(in) :: lb,lm,le
  377. integer :: l1,l2,l
  378. l1=lb
  379. l2=lm
  380. do l=lb,le
  381. if(l2.gt.le) then
  382. mtmp(l)=indx(l1)
  383. l1=l1+1
  384. elseif(l1.ge.lm) then
  385. mtmp(l)=indx(l2)
  386. l2=l2+1
  387. else
  388. if(dsnd) then
  389. if(keys(indx(l1)) .ge. keys(indx(l2))) then
  390. mtmp(l)=indx(l1)
  391. l1=l1+1
  392. else
  393. mtmp(l)=indx(l2)
  394. l2=l2+1
  395. endif
  396. else
  397. if(keys(indx(l1)) .le. keys(indx(l2))) then
  398. mtmp(l)=indx(l1)
  399. l1=l1+1
  400. else
  401. mtmp(l)=indx(l2)
  402. l2=l2+1
  403. endif
  404. endif
  405. endif
  406. end do
  407. end subroutine merge_
  408. end subroutine rSort_
  409. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  410. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  411. !BOP -------------------------------------------------------------------
  412. !
  413. ! !IROUTINE: dSort_ - A stable merge index sorting DOUBLEs.
  414. !
  415. ! !DESCRIPTION:
  416. !
  417. ! !INTERFACE:
  418. subroutine dSort_(indx,keys,descend,stat)
  419. use m_stdio, only : stderr
  420. use m_die, only : die
  421. use m_realkinds,only : DP
  422. implicit none
  423. integer, dimension(:), intent(inout) :: indx
  424. real(DP), dimension(:), intent(in) :: keys
  425. logical, optional, intent(in) :: descend
  426. integer, optional, intent(out) :: stat
  427. ! !REVISION HISTORY:
  428. ! 15Mar00 - Jing Guo
  429. ! . Modified the interface, by removing the explicit size
  430. ! 02Feb99 - Jing Guo <guo@thunder> - Added if(present(stat)) ...
  431. ! 04Jan99 - Jing Guo <guo@thunder> - revised the prolog
  432. ! 09Sep97 - Jing Guo <guo@thunder> - initial prototype/prolog/code
  433. !EOP ___________________________________________________________________
  434. logical :: dsnd
  435. integer :: ierr
  436. integer, dimension(:),allocatable :: mtmp
  437. integer :: n
  438. character(len=*),parameter :: myname_=myname//'::dSort_'
  439. if(present(stat)) stat=0
  440. n=size(indx)
  441. allocate(mtmp(n),stat=ierr)
  442. if(ierr /= 0) then
  443. write(stderr,'(2a,i4)') myname_, &
  444. ': allocate(mtmp(:)) error, stat =',ierr
  445. if(.not.present(stat)) call die(myname_)
  446. stat=ierr
  447. return
  448. endif
  449. dsnd=.false.
  450. if(present(descend)) dsnd=descend
  451. call MergeSort_()
  452. deallocate(mtmp)
  453. contains
  454. subroutine MergeSort_()
  455. implicit none
  456. integer :: mstep,lstep
  457. integer :: lb,lm,le
  458. mstep=1
  459. do while(mstep < n)
  460. lstep=mstep*2
  461. lb=1
  462. do while(lb < n)
  463. lm=lb+mstep
  464. le=min(lm-1+mstep,n)
  465. call merge_(lb,lm,le)
  466. indx(lb:le)=mtmp(lb:le)
  467. lb=le+1
  468. end do
  469. mstep=lstep
  470. end do
  471. end subroutine MergeSort_
  472. subroutine merge_(lb,lm,le)
  473. integer,intent(in) :: lb,lm,le
  474. integer :: l1,l2,l
  475. l1=lb
  476. l2=lm
  477. do l=lb,le
  478. if(l2.gt.le) then
  479. mtmp(l)=indx(l1)
  480. l1=l1+1
  481. elseif(l1.ge.lm) then
  482. mtmp(l)=indx(l2)
  483. l2=l2+1
  484. else
  485. if(dsnd) then
  486. if(keys(indx(l1)) .ge. keys(indx(l2))) then
  487. mtmp(l)=indx(l1)
  488. l1=l1+1
  489. else
  490. mtmp(l)=indx(l2)
  491. l2=l2+1
  492. endif
  493. else
  494. if(keys(indx(l1)) .le. keys(indx(l2))) then
  495. mtmp(l)=indx(l1)
  496. l1=l1+1
  497. else
  498. mtmp(l)=indx(l2)
  499. l2=l2+1
  500. endif
  501. endif
  502. endif
  503. end do
  504. end subroutine merge_
  505. end subroutine dSort_
  506. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  507. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  508. !BOP -------------------------------------------------------------------
  509. !
  510. ! !IROUTINE: cSort_ - A stable merge index sorting of CHAR(*)s.
  511. !
  512. ! !DESCRIPTION:
  513. !
  514. ! !INTERFACE:
  515. subroutine cSort_(indx,keys,descend,stat)
  516. use m_stdio, only : stderr
  517. use m_die, only : die
  518. implicit none
  519. integer, dimension(:), intent(inout) :: indx
  520. character(len=*), dimension(:), intent(in) :: keys
  521. logical, optional, intent(in) :: descend
  522. integer, optional, intent(out) :: stat
  523. ! !REVISION HISTORY:
  524. ! 15Mar00 - Jing Guo
  525. ! . Modified the interface, by removing the explicit size
  526. ! 02Feb99 - Jing Guo <guo@thunder> - Added if(present(stat)) ...
  527. ! 04Jan99 - Jing Guo <guo@thunder> - revised the prolog
  528. ! 09Sep97 - Jing Guo <guo@thunder> - initial prototype/prolog/code
  529. !EOP ___________________________________________________________________
  530. logical :: dsnd
  531. integer :: ierr
  532. integer, dimension(:),allocatable :: mtmp
  533. integer :: n
  534. character(len=*),parameter :: myname_=myname//'::cSort_'
  535. if(present(stat)) stat=0
  536. n=size(indx)
  537. allocate(mtmp(n),stat=ierr)
  538. if(ierr /= 0) then
  539. write(stderr,'(2a,i4)') myname_, &
  540. ': allocate(mtmp(:)) error, stat =',ierr
  541. if(.not.present(stat)) call die(myname_)
  542. stat=ierr
  543. return
  544. endif
  545. dsnd=.false.
  546. if(present(descend)) dsnd=descend
  547. call MergeSort_()
  548. deallocate(mtmp)
  549. contains
  550. subroutine MergeSort_()
  551. implicit none
  552. integer :: mstep,lstep
  553. integer :: lb,lm,le
  554. mstep=1
  555. do while(mstep < n)
  556. lstep=mstep*2
  557. lb=1
  558. do while(lb < n)
  559. lm=lb+mstep
  560. le=min(lm-1+mstep,n)
  561. call merge_(lb,lm,le)
  562. indx(lb:le)=mtmp(lb:le)
  563. lb=le+1
  564. end do
  565. mstep=lstep
  566. end do
  567. end subroutine MergeSort_
  568. subroutine merge_(lb,lm,le)
  569. integer,intent(in) :: lb,lm,le
  570. integer :: l1,l2,l
  571. l1=lb
  572. l2=lm
  573. do l=lb,le
  574. if(l2.gt.le) then
  575. mtmp(l)=indx(l1)
  576. l1=l1+1
  577. elseif(l1.ge.lm) then
  578. mtmp(l)=indx(l2)
  579. l2=l2+1
  580. else
  581. if(dsnd) then
  582. if(keys(indx(l1)) .ge. keys(indx(l2))) then
  583. mtmp(l)=indx(l1)
  584. l1=l1+1
  585. else
  586. mtmp(l)=indx(l2)
  587. l2=l2+1
  588. endif
  589. else
  590. if(keys(indx(l1)) .le. keys(indx(l2))) then
  591. mtmp(l)=indx(l1)
  592. l1=l1+1
  593. else
  594. mtmp(l)=indx(l2)
  595. l2=l2+1
  596. endif
  597. endif
  598. endif
  599. end do
  600. end subroutine merge_
  601. end subroutine cSort_
  602. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  603. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  604. !BOP -------------------------------------------------------------------
  605. !
  606. ! !IROUTINE: iSort1_ - A stable merge index sorting of INTs.
  607. !
  608. ! !DESCRIPTION:
  609. !
  610. ! !INTERFACE:
  611. subroutine iSort1_(indx,keys,ikey,descend,stat)
  612. use m_stdio, only : stderr
  613. use m_die, only : die
  614. implicit none
  615. integer, dimension(:), intent(inout) :: indx
  616. integer, dimension(:,:), intent(in) :: keys
  617. integer,intent(in) :: ikey
  618. logical, optional, intent(in) :: descend
  619. integer, optional, intent(out) :: stat
  620. ! !REVISION HISTORY:
  621. ! 15Mar00 - Jing Guo
  622. ! . initial prototype/prolog/code
  623. ! . Copied code from iSort_
  624. ! . Extended the interface and the algorithm to handle
  625. ! 2-d arrays with an index.
  626. !EOP ___________________________________________________________________
  627. logical :: dsnd
  628. integer :: ierr
  629. integer, dimension(:),allocatable :: mtmp
  630. integer :: n
  631. character(len=*),parameter :: myname_=myname//'::iSort1_'
  632. if(present(stat)) stat=0
  633. n=size(indx)
  634. allocate(mtmp(n),stat=ierr)
  635. if(ierr /= 0) then
  636. write(stderr,'(2a,i4)') myname_, &
  637. ': allocate(mtmp(:)) error, stat =',ierr
  638. if(.not.present(stat)) call die(myname_)
  639. stat=ierr
  640. return
  641. endif
  642. dsnd=.false.
  643. if(present(descend)) dsnd=descend
  644. call MergeSort_()
  645. deallocate(mtmp)
  646. contains
  647. subroutine MergeSort_()
  648. implicit none
  649. integer :: mstep,lstep
  650. integer :: lb,lm,le
  651. mstep=1
  652. do while(mstep < n)
  653. lstep=mstep*2
  654. lb=1
  655. do while(lb < n)
  656. lm=lb+mstep
  657. le=min(lm-1+mstep,n)
  658. call merge_(lb,lm,le)
  659. indx(lb:le)=mtmp(lb:le)
  660. lb=le+1
  661. end do
  662. mstep=lstep
  663. end do
  664. end subroutine MergeSort_
  665. subroutine merge_(lb,lm,le)
  666. integer,intent(in) :: lb,lm,le
  667. integer :: l1,l2,l
  668. l1=lb
  669. l2=lm
  670. do l=lb,le
  671. if(l2.gt.le) then
  672. mtmp(l)=indx(l1)
  673. l1=l1+1
  674. elseif(l1.ge.lm) then
  675. mtmp(l)=indx(l2)
  676. l2=l2+1
  677. else
  678. if(dsnd) then
  679. if(keys(ikey,indx(l1)) .ge. keys(ikey,indx(l2))) then
  680. mtmp(l)=indx(l1)
  681. l1=l1+1
  682. else
  683. mtmp(l)=indx(l2)
  684. l2=l2+1
  685. endif
  686. else
  687. if(keys(ikey,indx(l1)) .le. keys(ikey,indx(l2))) then
  688. mtmp(l)=indx(l1)
  689. l1=l1+1
  690. else
  691. mtmp(l)=indx(l2)
  692. l2=l2+1
  693. endif
  694. endif
  695. endif
  696. end do
  697. end subroutine merge_
  698. end subroutine iSort1_
  699. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  700. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  701. !BOP -------------------------------------------------------------------
  702. !
  703. ! !IROUTINE: rSort1_ - A stable merge index sorting REALs.
  704. !
  705. ! !DESCRIPTION:
  706. !
  707. ! !INTERFACE:
  708. subroutine rSort1_(indx,keys,ikey,descend,stat)
  709. use m_stdio, only : stderr
  710. use m_die, only : die
  711. use m_realkinds,only : SP
  712. implicit none
  713. integer, dimension(:), intent(inout) :: indx
  714. real(SP),dimension(:,:), intent(in) :: keys
  715. integer,intent(in) :: ikey
  716. logical, optional, intent(in) :: descend
  717. integer, optional, intent(out) :: stat
  718. ! !REVISION HISTORY:
  719. ! 15Mar00 - Jing Guo
  720. ! . initial prototype/prolog/code
  721. ! . Copied code from rSort_
  722. ! . Extended the interface and the algorithm to handle
  723. ! 2-d arrays with an index.
  724. !EOP ___________________________________________________________________
  725. logical :: dsnd
  726. integer :: ierr
  727. integer, dimension(:),allocatable :: mtmp
  728. integer :: n
  729. character(len=*),parameter :: myname_=myname//'::rSort1_'
  730. if(present(stat)) stat=0
  731. n=size(indx)
  732. allocate(mtmp(n),stat=ierr)
  733. if(ierr /= 0) then
  734. write(stderr,'(2a,i4)') myname_, &
  735. ': allocate(mtmp(:)) error, stat =',ierr
  736. if(.not.present(stat)) call die(myname_)
  737. stat=ierr
  738. return
  739. endif
  740. dsnd=.false.
  741. if(present(descend)) dsnd=descend
  742. call MergeSort_()
  743. deallocate(mtmp)
  744. contains
  745. subroutine MergeSort_()
  746. implicit none
  747. integer :: mstep,lstep
  748. integer :: lb,lm,le
  749. mstep=1
  750. do while(mstep < n)
  751. lstep=mstep*2
  752. lb=1
  753. do while(lb < n)
  754. lm=lb+mstep
  755. le=min(lm-1+mstep,n)
  756. call merge_(lb,lm,le)
  757. indx(lb:le)=mtmp(lb:le)
  758. lb=le+1
  759. end do
  760. mstep=lstep
  761. end do
  762. end subroutine MergeSort_
  763. subroutine merge_(lb,lm,le)
  764. integer,intent(in) :: lb,lm,le
  765. integer :: l1,l2,l
  766. l1=lb
  767. l2=lm
  768. do l=lb,le
  769. if(l2.gt.le) then
  770. mtmp(l)=indx(l1)
  771. l1=l1+1
  772. elseif(l1.ge.lm) then
  773. mtmp(l)=indx(l2)
  774. l2=l2+1
  775. else
  776. if(dsnd) then
  777. if(keys(ikey,indx(l1)) .ge. keys(ikey,indx(l2))) then
  778. mtmp(l)=indx(l1)
  779. l1=l1+1
  780. else
  781. mtmp(l)=indx(l2)
  782. l2=l2+1
  783. endif
  784. else
  785. if(keys(ikey,indx(l1)) .le. keys(ikey,indx(l2))) then
  786. mtmp(l)=indx(l1)
  787. l1=l1+1
  788. else
  789. mtmp(l)=indx(l2)
  790. l2=l2+1
  791. endif
  792. endif
  793. endif
  794. end do
  795. end subroutine merge_
  796. end subroutine rSort1_
  797. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  798. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  799. !BOP -------------------------------------------------------------------
  800. !
  801. ! !IROUTINE: dSort1_ - A stable merge index sorting DOUBLEs.
  802. !
  803. ! !DESCRIPTION:
  804. !
  805. ! !INTERFACE:
  806. subroutine dSort1_(indx,keys,ikey,descend,stat)
  807. use m_stdio, only : stderr
  808. use m_die, only : die
  809. use m_realkinds,only : DP
  810. implicit none
  811. integer, dimension(:), intent(inout) :: indx
  812. real(DP), dimension(:,:), intent(in) :: keys
  813. integer,intent(in) :: ikey
  814. logical, optional, intent(in) :: descend
  815. integer, optional, intent(out) :: stat
  816. ! !REVISION HISTORY:
  817. ! 15Mar00 - Jing Guo
  818. ! . initial prototype/prolog/code
  819. ! . Copied code from dSort_
  820. ! . Extended the interface and the algorithm to handle
  821. ! 2-d arrays with an index.
  822. !EOP ___________________________________________________________________
  823. logical :: dsnd
  824. integer :: ierr
  825. integer, dimension(:),allocatable :: mtmp
  826. integer :: n
  827. character(len=*),parameter :: myname_=myname//'::dSort1_'
  828. if(present(stat)) stat=0
  829. n=size(indx)
  830. allocate(mtmp(n),stat=ierr)
  831. if(ierr /= 0) then
  832. write(stderr,'(2a,i4)') myname_, &
  833. ': allocate(mtmp(:)) error, stat =',ierr
  834. if(.not.present(stat)) call die(myname_)
  835. stat=ierr
  836. return
  837. endif
  838. dsnd=.false.
  839. if(present(descend)) dsnd=descend
  840. call MergeSort_()
  841. deallocate(mtmp)
  842. contains
  843. subroutine MergeSort_()
  844. implicit none
  845. integer :: mstep,lstep
  846. integer :: lb,lm,le
  847. mstep=1
  848. do while(mstep < n)
  849. lstep=mstep*2
  850. lb=1
  851. do while(lb < n)
  852. lm=lb+mstep
  853. le=min(lm-1+mstep,n)
  854. call merge_(lb,lm,le)
  855. indx(lb:le)=mtmp(lb:le)
  856. lb=le+1
  857. end do
  858. mstep=lstep
  859. end do
  860. end subroutine MergeSort_
  861. subroutine merge_(lb,lm,le)
  862. integer,intent(in) :: lb,lm,le
  863. integer :: l1,l2,l
  864. l1=lb
  865. l2=lm
  866. do l=lb,le
  867. if(l2.gt.le) then
  868. mtmp(l)=indx(l1)
  869. l1=l1+1
  870. elseif(l1.ge.lm) then
  871. mtmp(l)=indx(l2)
  872. l2=l2+1
  873. else
  874. if(dsnd) then
  875. if(keys(ikey,indx(l1)) .ge. keys(ikey,indx(l2))) then
  876. mtmp(l)=indx(l1)
  877. l1=l1+1
  878. else
  879. mtmp(l)=indx(l2)
  880. l2=l2+1
  881. endif
  882. else
  883. if(keys(ikey,indx(l1)) .le. keys(ikey,indx(l2))) then
  884. mtmp(l)=indx(l1)
  885. l1=l1+1
  886. else
  887. mtmp(l)=indx(l2)
  888. l2=l2+1
  889. endif
  890. endif
  891. endif
  892. end do
  893. end subroutine merge_
  894. end subroutine dSort1_
  895. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  896. ! NASA/GSFC, Data Assimilation Office, Code 910.3, GEOS/DAS !
  897. !BOP -------------------------------------------------------------------
  898. !
  899. ! !IROUTINE: cSort1_ - A stable merge index sorting of CHAR(*)s.
  900. !
  901. ! !DESCRIPTION:
  902. !
  903. ! !INTERFACE:
  904. subroutine cSort1_(indx,keys,ikey,descend,stat)
  905. use m_stdio, only : stderr
  906. use m_die, only : die
  907. implicit none
  908. integer, dimension(:), intent(inout) :: indx
  909. character(len=*), dimension(:,:), intent(in) :: keys
  910. integer,intent(in) :: ikey
  911. logical, optional, intent(in) :: descend
  912. integer, optional, intent(out) :: stat
  913. ! !REVISION HISTORY:
  914. ! 15Mar00 - Jing Guo
  915. ! . initial prototype/prolog/code
  916. ! . Copied code from cSort_
  917. ! . Extended the interface and the algorithm to handle
  918. ! 2-d arrays with an index.
  919. !EOP ___________________________________________________________________
  920. logical :: dsnd
  921. integer :: ierr
  922. integer, dimension(:),allocatable :: mtmp
  923. integer :: n
  924. character(len=*),parameter :: myname_=myname//'::cSort1_'
  925. if(present(stat)) stat=0
  926. n=size(indx)
  927. allocate(mtmp(n),stat=ierr)
  928. if(ierr /= 0) then
  929. write(stderr,'(2a,i4)') myname_, &
  930. ': allocate(mtmp(:)) error, stat =',ierr
  931. if(.not.present(stat)) call die(myname_)
  932. stat=ierr
  933. return
  934. endif
  935. dsnd=.false.
  936. if(present(descend)) dsnd=descend
  937. call MergeSort_()
  938. deallocate(mtmp)
  939. contains
  940. subroutine MergeSort_()
  941. implicit none
  942. integer :: mstep,lstep
  943. integer :: lb,lm,le
  944. mstep=1
  945. do while(mstep < n)
  946. lstep=mstep*2
  947. lb=1
  948. do while(lb < n)
  949. lm=lb+mstep
  950. le=min(lm-1+mstep,n)
  951. call merge_(lb,lm,le)
  952. indx(lb:le)=mtmp(lb:le)
  953. lb=le+1
  954. end do
  955. mstep=lstep
  956. end do
  957. end subroutine MergeSort_
  958. subroutine merge_(lb,lm,le)
  959. integer,intent(in) :: lb,lm,le
  960. integer :: l1,l2,l
  961. l1=lb
  962. l2=lm
  963. do l=lb,le
  964. if(l2.gt.le) then
  965. mtmp(l)=indx(l1)
  966. l1=l1+1
  967. elseif(l1.ge.lm) then
  968. mtmp(l)=indx(l2)
  969. l2=l2+1
  970. else
  971. if(dsnd) then
  972. if(keys(ikey,indx(l1)) .ge. keys(ikey,indx(l2))) then
  973. mtmp(l)=indx(l1)
  974. l1=l1+1
  975. else
  976. mtmp(l)=indx(l2)
  977. l2=l2+1
  978. endif
  979. else
  980. if(keys(ikey,indx(l1)) .le. keys(ikey,indx(l2))) then
  981. mtmp(l)=indx(l1)
  982. l1=l1+1
  983. else
  984. mtmp(l)=indx(l2)
  985. l2=l2+1
  986. endif
  987. endif
  988. endif
  989. end do
  990. end subroutine merge_
  991. end subroutine cSort1_
  992. !-----------------------------------------------------------------------
  993. end module m_MergeSorts
  994. !.