source: LMDZ5/branches/LMDZ5V1.0-dev/libf/cosp/cosp_radar.F90 @ 4249

Last change on this file since 4249 was 1279, checked in by Laurent Fairhead, 15 years ago

Merged LMDZ4-dev branch changes r1241:1278 into the trunk
Running trunk and LMDZ4-dev in LMDZOR configuration on local
machine (sequential) and SX8 (4-proc) yields identical results
(restart and restartphy are identical binarily)
Log history from r1241 to r1278 is available by switching to
source:LMDZ4/branches/LMDZ4-dev-20091210

File size: 8.8 KB
Line 
1! (c) British Crown Copyright 2008, the Met Office.
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without modification, are permitted
5! provided that the following conditions are met:
6!
7!     * Redistributions of source code must retain the above copyright notice, this list
8!       of conditions and the following disclaimer.
9!     * Redistributions in binary form must reproduce the above copyright notice, this list
10!       of conditions and the following disclaimer in the documentation and/or other materials
11!       provided with the distribution.
12!     * Neither the name of the Met Office nor the names of its contributors may be used
13!       to endorse or promote products derived from this software without specific prior written
14!       permission.
15!
16! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR
17! IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND
18! FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR
19! CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
20! DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
21! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER
22! IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
23! OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
24
25MODULE MOD_COSP_RADAR
26  USE MOD_COSP_CONSTANTS
27  USE MOD_COSP_TYPES
28  use radar_simulator_types
29  use array_lib
30  use atmos_lib
31  use format_input
32  IMPLICIT NONE
33 
34  INTERFACE
35    subroutine radar_simulator(freq,k2,do_ray,use_gas_abs,use_mie_table,mt, &
36        nhclass,hp,nprof,ngate,nsizes,D,hgt_matrix,hm_matrix,re_matrix,p_matrix,t_matrix, &
37        rh_matrix,Ze_non,Ze_ray,h_atten_to_vol,g_atten_to_vol,dBZe, &
38        g_to_vol_in,g_to_vol_out)
39 
40        use m_mrgrnk
41        use array_lib
42        use math_lib
43        use optics_lib
44        use radar_simulator_types
45        implicit none
46        ! ----- INPUTS ----- 
47        type(mie), intent(in) :: mt
48        type(class_param) :: hp
49        real*8, intent(in) :: freq,k2
50        integer, intent(in) ::  do_ray,use_gas_abs,use_mie_table, &
51            nhclass,nprof,ngate,nsizes
52        real*8, dimension(nsizes), intent(in) :: D
53        real*8, dimension(nprof,ngate), intent(in) :: hgt_matrix, p_matrix, &
54            t_matrix,rh_matrix
55        real*8, dimension(nhclass,nprof,ngate), intent(in) :: hm_matrix
56        real*8, dimension(nhclass,nprof,ngate), intent(inout) :: re_matrix
57        ! ----- OUTPUTS -----
58        real*8, dimension(nprof,ngate), intent(out) :: Ze_non,Ze_ray, &
59            g_atten_to_vol,dBZe,h_atten_to_vol   
60        ! ----- OPTIONAL -----
61        real*8, optional, dimension(ngate,nprof) :: &
62            g_to_vol_in,g_to_vol_out
63     end subroutine radar_simulator
64  END INTERFACE
65
66CONTAINS
67
68!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
69!------------------- SUBROUTINE COSP_RADAR ------------------------
70!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
71SUBROUTINE COSP_RADAR(gbx,sgx,sghydro,z)
72  IMPLICIT NONE
73
74  ! Arguments
75  type(cosp_gridbox),intent(in) :: gbx  ! Gridbox info
76  type(cosp_subgrid),intent(in) :: sgx  ! Subgrid info
77  type(cosp_sghydro),intent(in) :: sghydro  ! Subgrid info for hydrometeors
78  type(cosp_sgradar),intent(inout) :: z ! Output from simulator, subgrid
79
80  ! Local variables
81  integer :: &
82  nsizes                        ! num of discrete drop sizes
83
84  real*8 :: &
85  freq, &                       ! radar frequency (GHz)
86  k2                            ! |K|^2, -1=use frequency dependent default
87 
88  real*8, dimension(:,:), allocatable :: &
89  g_to_vol ! integrated atten due to gases, r>v (dB)
90 
91  real*8, dimension(:,:), allocatable :: &
92  Ze_non, &                     ! radar reflectivity withOUT attenuation (dBZ)
93  Ze_ray, &                     ! Rayleigh reflectivity (dBZ)
94  h_atten_to_vol, &             ! attenuation by hydromets, radar to vol (dB)
95  g_atten_to_vol, &             ! gaseous atteunation, radar to vol (dB)
96  dBZe, &                       ! effective radar reflectivity factor (dBZ)
97  hgt_matrix, &                 ! height of hydrometeors (km)
98  t_matrix, &                   !temperature (k)
99  p_matrix, &                   !pressure (hPa)
100  rh_matrix                     !relative humidity (%)
101 
102  real*8, dimension(:,:,:), allocatable :: &
103  hm_matrix, &                  ! hydrometeor mixing ratio (g kg^-1)
104  re_matrix
105
106  integer, parameter :: one = 1
107  logical :: hgt_reversed
108  integer :: pr,i,j,k,unt
109
110! ----- main program settings ------
111
112  freq = gbx%radar_freq
113  k2 = gbx%k2
114 
115  !
116  ! note:  intitialization section that was here has been relocated to SUBROUTINE CONSTRUCT_COSP_GRIDBOX by roj, Feb 2008
117  !
118  mt_ttl=gbx%mt_ttl  ! these variables really should be moved into the mt structure rather than kept as global arrays.
119  mt_tti=gbx%mt_tti
120
121  ! Inputs to Quickbeam
122  allocate(hgt_matrix(gbx%Npoints,gbx%Nlevels),p_matrix(gbx%Npoints,gbx%Nlevels), &
123           t_matrix(gbx%Npoints,gbx%Nlevels),rh_matrix(gbx%Npoints,gbx%Nlevels))
124  allocate(hm_matrix(gbx%Nhydro,gbx%Npoints,gbx%Nlevels))
125  allocate(re_matrix(gbx%Nhydro,gbx%Npoints,gbx%Nlevels))
126
127  ! Outputs from Quickbeam
128  allocate(Ze_non(gbx%Npoints,gbx%Nlevels))
129  allocate(Ze_ray(gbx%Npoints,gbx%Nlevels))
130  allocate(h_atten_to_vol(gbx%Npoints,gbx%Nlevels))
131  allocate(g_atten_to_vol(gbx%Npoints,gbx%Nlevels))
132  allocate(dBZe(gbx%Npoints,gbx%Nlevels))
133 
134  ! Optional argument. It is computed and returned in the first call to
135  ! radar_simulator, and passed as input in the rest
136  allocate(g_to_vol(gbx%Nlevels,gbx%Npoints))
137 
138  p_matrix   = gbx%p/100.0     ! From Pa to hPa
139  hgt_matrix = gbx%zlev/1000.0 ! From m to km
140  t_matrix   = gbx%T-273.15    ! From K to C
141  rh_matrix  = gbx%q
142  re_matrix  = 0.0
143 
144  ! Quickbeam assumes the first row is closest to the radar
145  call order_data(hgt_matrix,hm_matrix,p_matrix,t_matrix, &
146      rh_matrix,gbx%surface_radar,hgt_reversed)
147 
148  ! ----- loop over subcolumns -----
149  do pr=1,sgx%Ncolumns
150      !  atmospheric profiles are the same within the same gridbox
151      !  only hydrometeor profiles will be different
152      if (hgt_reversed) then 
153         do i=1,gbx%Nhydro 
154            hm_matrix(i,:,:) = sghydro%mr_hydro(:,pr,gbx%Nlevels:1:-1,i)*1000.0 ! Units from kg/kg to g/kg
155            if (gbx%use_reff) then
156              re_matrix(i,:,:) = sghydro%Reff(:,pr,gbx%Nlevels:1:-1,i)*1.e6     ! Units from m to micron
157            endif
158         enddo 
159      else 
160         do i=1,gbx%Nhydro
161            hm_matrix(i,:,:) = sghydro%mr_hydro(:,pr,:,i)*1000.0 ! Units from kg/kg to g/kg
162            if (gbx%use_reff) then
163              re_matrix(i,:,:) = sghydro%Reff(:,pr,:,i)*1.e6       ! Units from m to micron
164            endif
165         enddo
166      endif 
167
168      !   ----- call radar simulator -----
169      if (pr == 1) then ! Compute gaseous attenuation for all profiles
170         j=0
171         if (gbx%Npoints == 53) then
172           unt=10
173           j=1
174         endif
175         if (gbx%Npoints == 153) then
176           unt=11
177           j=101
178         endif
179         call radar_simulator(freq,k2,gbx%do_ray,gbx%use_gas_abs,gbx%use_mie_tables,gbx%mt, &    !  v0.2: mt changed to gbx%mt, roj
180           gbx%Nhydro,gbx%hp,gbx%Npoints,gbx%Nlevels,gbx%nsizes,gbx%D, &                         !  v0.2: hp->gbx%hp, D->gbx%d, nsizes->gbx%nsizes, roj
181           hgt_matrix,hm_matrix,re_matrix,p_matrix,t_matrix,rh_matrix, &
182           Ze_non,Ze_ray,h_atten_to_vol,g_atten_to_vol,dBZe,g_to_vol_out=g_to_vol)
183      else ! Use gaseous atteunuation for pr = 1
184         call radar_simulator(freq,k2,gbx%do_ray,gbx%use_gas_abs,gbx%use_mie_tables,gbx%mt, &
185           gbx%Nhydro,gbx%hp,gbx%Npoints,gbx%Nlevels,gbx%nsizes,gbx%D, &
186           hgt_matrix,hm_matrix,re_matrix,p_matrix,t_matrix,rh_matrix, &
187           Ze_non,Ze_ray,h_atten_to_vol,g_atten_to_vol,dBZe,g_to_vol_in=g_to_vol)
188      endif
189      ! ----- BEGIN output section -----
190      ! spaceborne radar : from TOA to SURFACE
191      if (gbx%surface_radar == 1) then
192        z%Ze_tot(:,pr,:)=dBZe(:,:)
193      else if (gbx%surface_radar == 0) then ! Spaceborne
194        z%Ze_tot(:,pr,:)=dBZe(:,gbx%Nlevels:1:-1)
195      endif
196
197  enddo !pr
198 
199  ! Change undefined value to one defined in COSP
200  where (z%Ze_tot == -999.0) z%Ze_tot = R_UNDEF
201
202  deallocate(hgt_matrix,p_matrix,t_matrix,rh_matrix)
203  deallocate(hm_matrix,re_matrix, &
204      Ze_non,Ze_ray,h_atten_to_vol,g_atten_to_vol,dBZe)
205  deallocate(g_to_vol)
206 
207  ! deallocate(mt_ttl,mt_tti)   !v0.2: roj feb 2008 can not be done here,
208                                !these variables now part of gbx structure and dealocated later
209
210END SUBROUTINE COSP_RADAR
211
212END MODULE MOD_COSP_RADAR
Note: See TracBrowser for help on using the repository browser.