prismtrs_sort_add.F90

Go to the documentation of this file.
00001 !------------------------------------------------------------------------
00002 ! Copyright 2006-2010, CERFACS, Toulouse, France.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PRISMTrs_sort_add
00008 !
00009 ! !INTERFACE
00010       subroutine prismtrs_sort_add(add1, add2, weights, num_links, num_wts)
00011 
00012       USE PRISMDrv, dummy_INTERFACE => PRISMTrs_sort_add
00013       implicit none
00014 
00015 !-----------------------------------------------------------------------
00016 !
00017 !     Input and Output arrays
00018 !
00019 !-----------------------------------------------------------------------
00020 
00021       integer, intent(in) :: num_links, num_wts
00022       integer, intent(inout), dimension(num_links) :: 
00023               add1,      ! destination address array (num_links)
00024               add2        ! source      address array
00025 
00026       DOUBLE PRECISION, intent(inout), 
00027           dimension(num_wts, num_links) :: 
00028               weights     ! remapping weights (num_wts, num_links)
00029 !
00030 ! !DESCRIPTION
00031 !
00032 ! SUBROUTINE prismtrs_sort_add sorts address and weight arrays based 
00033 ! on the destination address with the source address as a secondary
00034 ! sorting criterion.  the method is a standard heap sort.This subroutine
00035 ! comes from the SCRIP 1.4 library available from Los Alamos National 
00036 ! Laboratory.
00037 !
00038 ! !REVISED HISTORY
00039 !   Date      Programmer   Description
00040 ! ----------  ----------   -----------
00041 ! 19/03/2007  S. Valcke    Duplicated from the SCRIP library
00042 !
00043 ! EOP
00044 !----------------------------------------------------------------------
00045 ! $Id: prismtrs_sort_add.F90 2082 2009-10-21 13:31:19Z hanke $
00046 ! $Author: hanke $
00047 !----------------------------------------------------------------------
00048 !
00049 !-----------------------------------------------------------------------
00050 !
00051 !     local variables
00052 !
00053   CHARACTER(LEN=len_cvs_string), SAVE  :: mycvs = 
00054      '$Id: prismtrs_sort_add.F90 2082 2009-10-21 13:31:19Z hanke $'
00055 !-----------------------------------------------------------------------
00056 
00057       integer :: 
00058                 add1_tmp, add2_tmp,  ! temp for addresses during swap
00059                 nwgt, 
00060                 lvl, final_lvl,    ! level indexes for heap sort levels
00061                 chk_lvl1, chk_lvl2, max_lvl
00062 
00063       double precision, dimension(SIZE(weights,DIM=1)) :: 
00064                 wgttmp              ! temp for holding wts during swap
00065 !
00066 #ifdef VERBOSE
00067   PRINT *, '| | | | | | | Enter PRISMTrs_sort_add'
00068   call psmile_flushstd
00069 #endif
00070 !-----------------------------------------------------------------------
00071 !
00072 !     start at the lowest level (N/2) of the tree and sift lower 
00073 !     values to the bottom of the tree, promoting the larger numbers
00074 !
00075 !-----------------------------------------------------------------------
00076 
00077       do lvl=num_links/2,1,-1
00078 
00079         final_lvl = lvl
00080         add1_tmp = add1(lvl)
00081         add2_tmp = add2(lvl)
00082         wgttmp(:) = weights(:,lvl)
00083 
00084         !***
00085         !*** loop until proper level is found for this link, or reach
00086         !*** bottom
00087         !***
00088 
00089         sift_loop1: do
00090 
00091           !***
00092           !*** find the largest of the two daughters
00093           !***
00094 
00095           chk_lvl1 = 2*final_lvl
00096           chk_lvl2 = 2*final_lvl+1
00097           if (chk_lvl1 .EQ. num_links) chk_lvl2 = chk_lvl1
00098 
00099           if ((add1(chk_lvl1) >  add1(chk_lvl2)) .OR. &
00100              ((add1(chk_lvl1) == add1(chk_lvl2)) .AND. &
00101               (add2(chk_lvl1) >  add2(chk_lvl2)))) then
00102             max_lvl = chk_lvl1
00103           else 
00104             max_lvl = chk_lvl2
00105           endif
00106 
00107           !***
00108           !*** if the parent is greater than both daughters,
00109           !*** the correct level has been found
00110           !***
00111 
00112           if ((add1_tmp .GT. add1(max_lvl)) .OR. &
00113              ((add1_tmp .EQ. add1(max_lvl)) .AND. &
00114               (add2_tmp .GT. add2(max_lvl)))) then
00115             add1(final_lvl) = add1_tmp
00116             add2(final_lvl) = add2_tmp
00117             weights(:,final_lvl) = wgttmp(:)
00118             exit sift_loop1
00119 
00120           !***
00121           !*** otherwise, promote the largest daughter and push
00122           !*** down one level in the tree.  if haven't reached
00123           !*** the end of the tree, repeat the process.  otherwise
00124           !*** store last values and exit the loop
00125           !***
00126 
00127           else 
00128             add1(final_lvl) = add1(max_lvl)
00129             add2(final_lvl) = add2(max_lvl)
00130             weights(:,final_lvl) = weights(:,max_lvl)
00131 
00132             final_lvl = max_lvl
00133             if (2*final_lvl > num_links) then
00134               add1(final_lvl) = add1_tmp
00135               add2(final_lvl) = add2_tmp
00136               weights(:,final_lvl) = wgttmp(:)
00137               exit sift_loop1
00138             endif
00139           endif
00140         end do sift_loop1
00141       end do
00142 
00143 !-----------------------------------------------------------------------
00144 !
00145 !     now that the heap has been sorted, strip off the top (largest)
00146 !     value and promote the values below
00147 !
00148 !-----------------------------------------------------------------------
00149 
00150       do lvl=num_links,3,-1
00151 
00152         !***
00153         !*** move the top value and insert it into the correct place
00154         !***
00155 
00156         add1_tmp = add1(lvl)
00157         add1(lvl) = add1(1)
00158 
00159         add2_tmp = add2(lvl)
00160         add2(lvl) = add2(1)
00161 
00162         wgttmp(:) = weights(:,lvl)
00163         weights(:,lvl) = weights(:,1)
00164 
00165         !***
00166         !*** as above this loop sifts the tmp values down until proper 
00167         !*** level is reached
00168         !***
00169 
00170         final_lvl = 1
00171 
00172         sift_loop2: do
00173 
00174           !***
00175           !*** find the largest of the two daughters
00176           !***
00177 
00178           chk_lvl1 = 2*final_lvl
00179           chk_lvl2 = 2*final_lvl+1
00180           if (chk_lvl2 >= lvl) chk_lvl2 = chk_lvl1
00181 
00182           if ((add1(chk_lvl1) >  add1(chk_lvl2)) .OR. &
00183              ((add1(chk_lvl1) == add1(chk_lvl2)) .AND. &
00184               (add2(chk_lvl1) >  add2(chk_lvl2)))) then
00185             max_lvl = chk_lvl1
00186           else 
00187             max_lvl = chk_lvl2
00188           endif
00189 
00190           !***
00191           !*** if the parent is greater than both daughters,
00192           !*** the correct level has been found
00193           !***
00194 
00195           if ((add1_tmp >  add1(max_lvl)) .OR. &
00196              ((add1_tmp == add1(max_lvl)) .AND. &
00197               (add2_tmp >  add2(max_lvl)))) then
00198             add1(final_lvl) = add1_tmp
00199             add2(final_lvl) = add2_tmp
00200             weights(:,final_lvl) = wgttmp(:)
00201             exit sift_loop2
00202 
00203           !***
00204           !*** otherwise, promote the largest daughter and push
00205           !*** down one level in the tree.  if haven't reached
00206           !*** the end of the tree, repeat the process.  otherwise
00207           !*** store last values and exit the loop
00208           !***
00209 
00210           else 
00211             add1(final_lvl) = add1(max_lvl)
00212             add2(final_lvl) = add2(max_lvl)
00213             weights(:,final_lvl) = weights(:,max_lvl)
00214 
00215             final_lvl = max_lvl
00216             if (2*final_lvl >= lvl) then
00217               add1(final_lvl) = add1_tmp
00218               add2(final_lvl) = add2_tmp
00219               weights(:,final_lvl) = wgttmp(:)
00220               exit sift_loop2
00221             endif
00222           endif
00223         end do sift_loop2
00224       end do
00225 
00226       !***
00227       !*** swap the last two entries
00228       !***
00229 
00230 
00231       add1_tmp = add1(2)
00232       add1(2)  = add1(1)
00233       add1(1)  = add1_tmp
00234 
00235       add2_tmp = add2(2)
00236       add2(2)  = add2(1)
00237       add2(1)  = add2_tmp
00238 
00239       wgttmp (:)   = weights(:,2)
00240       weights(:,2) = weights(:,1)
00241       weights(:,1) = wgttmp (:)
00242 !
00243 #ifdef VERBOSE
00244   PRINT *, '| | | | | | | Quit PRISMTrs_sort_add'
00245   call psmile_flushstd
00246 #endif
00247 !-----------------------------------------------------------------------
00248 
00249       end subroutine prismtrs_sort_add
00250 
00251 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00252 

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1