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