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