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