[1] | 1 | ! FORMAT_INPUT: Procedures to prepare data for input to the simulator |
---|
| 2 | ! Compiled/Modified: |
---|
| 3 | ! 08/28/2006 John Haynes (haynes@atmos.colostate.edu) |
---|
| 4 | ! |
---|
| 5 | ! irreg_to_grid (subroutine) |
---|
| 6 | ! order_data (subroutine) |
---|
| 7 | |
---|
| 8 | module format_input |
---|
| 9 | |
---|
| 10 | contains |
---|
| 11 | |
---|
| 12 | ! ---------------------------------------------------------------------------- |
---|
| 13 | ! SUBROUTINE IRREG_TO_GRID |
---|
| 14 | ! ---------------------------------------------------------------------------- |
---|
| 15 | subroutine irreg_to_grid(hgt_matrix,t_matrix,p_matrix,rh_matrix, & |
---|
| 16 | env_hgt_matrix,env_t_matrix,env_p_matrix,env_rh_matrix) |
---|
| 17 | use array_lib |
---|
| 18 | implicit none |
---|
| 19 | |
---|
| 20 | ! Purpose: |
---|
| 21 | ! Linearly interpolate sounding-level data to the hydrometeor-level |
---|
| 22 | ! resolution |
---|
| 23 | ! |
---|
| 24 | ! Inputs: |
---|
| 25 | ! [hgt_matrix] hydrometeor-level heights |
---|
| 26 | ! [env_hgt_matrix] sounding-level heights |
---|
| 27 | ! [env_t_matrix] sounding-level temperatures |
---|
| 28 | ! [env_p_matrix] sounding-level pressures |
---|
| 29 | ! [env_rh_matrix] sounding-level relative humidities |
---|
| 30 | ! |
---|
| 31 | ! Outputs: |
---|
| 32 | ! [t_matrix] hydrometeor-level temperatures |
---|
| 33 | ! [p_matrix] hydrometeor-level pressures |
---|
| 34 | ! [rh_matrix] hydrometeor-level relative humidities |
---|
| 35 | ! |
---|
| 36 | ! Created: |
---|
| 37 | ! 08/28/2006 John Haynes (haynes@atmos.colostate.edu) |
---|
| 38 | |
---|
| 39 | ! ----- INPUTS ----- |
---|
| 40 | real*8, dimension(:,:), intent(in) :: & |
---|
| 41 | hgt_matrix,env_hgt_matrix,env_t_matrix,env_p_matrix,env_rh_matrix |
---|
| 42 | |
---|
| 43 | ! ----- OUTPUTS ----- |
---|
| 44 | real*8, dimension(:,:), intent(out) :: & |
---|
| 45 | t_matrix,p_matrix,rh_matrix |
---|
| 46 | |
---|
| 47 | ! ----- INTERNAL ----- |
---|
| 48 | integer :: nprof, i |
---|
| 49 | integer,parameter :: KR8 = selected_real_kind(15,300) |
---|
| 50 | |
---|
| 51 | nprof = size(hgt_matrix,1) |
---|
| 52 | do i=1,nprof |
---|
| 53 | call lin_interpolate(env_t_matrix(i,:),env_hgt_matrix(i,:), & |
---|
| 54 | t_matrix(i,:),hgt_matrix(i,:),1000._KR8) |
---|
| 55 | call lin_interpolate(env_p_matrix(i,:),env_hgt_matrix(i,:), & |
---|
| 56 | p_matrix(i,:),hgt_matrix(i,:),1000._KR8) |
---|
| 57 | call lin_interpolate(env_rh_matrix(i,:),env_hgt_matrix(i,:), & |
---|
| 58 | rh_matrix(i,:),hgt_matrix(i,:),1000._KR8) |
---|
| 59 | enddo |
---|
| 60 | |
---|
| 61 | end subroutine irreg_to_grid |
---|
| 62 | |
---|
| 63 | ! ---------------------------------------------------------------------------- |
---|
| 64 | ! SUBROUTINE ORDER_DATA |
---|
| 65 | ! ---------------------------------------------------------------------------- |
---|
| 66 | subroutine order_data(hgt_matrix,hm_matrix,p_matrix,t_matrix, & |
---|
| 67 | rh_matrix,sfc_radar,hgt_reversed) |
---|
| 68 | implicit none |
---|
| 69 | |
---|
| 70 | ! Purpose: |
---|
| 71 | ! Ensure that input data is in top-down order/bottom-up order, |
---|
| 72 | ! for space-based/surface based radars, respectively |
---|
| 73 | ! |
---|
| 74 | ! Inputs: |
---|
| 75 | ! [hgt_matrix] heights |
---|
| 76 | ! [hm_matrix] mixing ratios |
---|
| 77 | ! [t_matrix] temperatures |
---|
| 78 | ! [p_matrix] pressures |
---|
| 79 | ! [rh_matrix] relative humidities |
---|
| 80 | ! [sfc_radar] 1=surface radar, 0=spaceborne |
---|
| 81 | ! |
---|
| 82 | ! Outputs: |
---|
| 83 | ! [hgt_matrix],[hm_matrix],[p_matrix,[t_matrix],[rh_matrix] in proper |
---|
| 84 | ! order for input to the radar simulator routine |
---|
| 85 | ! [hgt_reversed] T=heights were reordered,F=heights were not reordered |
---|
| 86 | ! |
---|
| 87 | ! Note: |
---|
| 88 | ! The order for all profiles is assumed to the same as the first profile. |
---|
| 89 | ! |
---|
| 90 | ! Created: |
---|
| 91 | ! 08/28/2006 John Haynes (haynes@atmos.colostate.edu) |
---|
| 92 | |
---|
| 93 | ! ----- INPUTS ----- |
---|
| 94 | integer, intent(in) :: sfc_radar |
---|
| 95 | |
---|
| 96 | ! ----- OUTPUTS ----- |
---|
| 97 | real*8, dimension(:,:), intent(inout) :: & |
---|
| 98 | hgt_matrix,p_matrix,t_matrix,rh_matrix |
---|
| 99 | real*8, dimension(:,:,:), intent(inout) :: & |
---|
| 100 | hm_matrix |
---|
| 101 | logical, intent(out) :: hgt_reversed |
---|
| 102 | |
---|
| 103 | ! ----- INTERNAL ----- |
---|
| 104 | integer :: ngate |
---|
| 105 | logical :: hgt_descending |
---|
| 106 | |
---|
| 107 | |
---|
| 108 | ngate = size(hgt_matrix,2) |
---|
| 109 | hgt_descending = hgt_matrix(1,1) > hgt_matrix(1,ngate) |
---|
| 110 | |
---|
| 111 | ! :: surface: heights must be ascending |
---|
| 112 | ! :: space-based: heights must be descending |
---|
| 113 | if ( & |
---|
| 114 | (sfc_radar == 1 .and. hgt_descending) .or. & |
---|
| 115 | (sfc_radar == 0 .and. (.not. hgt_descending)) & |
---|
| 116 | ) & |
---|
| 117 | then |
---|
| 118 | |
---|
| 119 | hgt_matrix(:,:) = hgt_matrix(:,ngate:1:-1) |
---|
| 120 | hm_matrix(:,:,:) = hm_matrix(:,:,ngate:1:-1) |
---|
| 121 | p_matrix(:,:) = p_matrix(:,ngate:1:-1) |
---|
| 122 | t_matrix(:,:) = t_matrix(:,ngate:1:-1) |
---|
| 123 | rh_matrix(:,:) = rh_matrix(:,ngate:1:-1) |
---|
| 124 | |
---|
| 125 | hgt_reversed = .true. |
---|
| 126 | else |
---|
| 127 | hgt_reversed = .false. |
---|
| 128 | endif |
---|
| 129 | |
---|
| 130 | end subroutine order_data |
---|
| 131 | |
---|
| 132 | end module format_input |
---|