source: LMDZ6/trunk/libf/phylmd/cosp/mod_cosp_radar.F90 @ 3625

Last change on this file since 3625 was 3233, checked in by idelkadi, 7 years ago

Nettoyage du code :

  • changement de nomes des routines pour avoir les memes nome des modules
  • corrections
  • dos2unix appliquee aux fichiers
  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 7.5 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 MOD_COSP_UTILS
29  use radar_simulator_types
30  use array_lib
31  use atmos_lib
32  use format_input
33  IMPLICIT NONE
34
35  INTERFACE
36    subroutine radar_simulator(hp,nprof,ngate,undef, &
37        hgt_matrix,hm_matrix,re_matrix,Np_matrix, &
38        p_matrix,t_matrix,rh_matrix, &
39        Ze_non,Ze_ray,g_to_vol,a_to_vol,dBZe, &
40        g_to_vol_in,g_to_vol_out)
41
42        use m_mrgrnk
43        use array_lib
44        use math_lib
45        use optics_lib
46        use radar_simulator_types
47        implicit none
48
49        ! ----- INPUTS ----- 
50        type(class_param) :: hp
51
52        integer, intent(in) :: nprof,ngate
53
54        real undef
55        real*8, dimension(nprof,ngate), intent(in) :: hgt_matrix, p_matrix, &
56            t_matrix,rh_matrix
57        real*8, dimension(hp%nhclass,nprof,ngate), intent(in) :: hm_matrix
58        real*8, dimension(hp%nhclass,nprof,ngate), intent(inout) :: re_matrix
59        real*8, dimension(hp%nhclass,nprof,ngate), intent(inout) :: Np_matrix
60
61        ! ----- OUTPUTS -----
62        real*8, dimension(nprof,ngate), intent(out) :: Ze_non,Ze_ray, &
63            g_to_vol,dBZe,a_to_vol
64        ! ----- OPTIONAL -----
65        real*8, optional, dimension(nprof,ngate) :: &
66            g_to_vol_in,g_to_vol_out
67     end subroutine radar_simulator
68  END INTERFACE
69
70CONTAINS
71
72!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
73!------------------- SUBROUTINE COSP_RADAR ------------------------
74!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
75SUBROUTINE COSP_RADAR(gbx,sgx,sghydro,z)
76  IMPLICIT NONE
77
78  ! Arguments
79  type(cosp_gridbox),intent(inout) :: gbx  ! Gridbox info
80  type(cosp_subgrid),intent(in) :: sgx  ! Subgrid info
81  type(cosp_sghydro),intent(in) :: sghydro  ! Subgrid info for hydrometeors
82  type(cosp_sgradar),intent(inout) :: z ! Output from simulator, subgrid
83
84  ! Local variables
85  integer :: &
86  nsizes            ! num of discrete drop sizes
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, &          ! effective radius (microns).   Optional. 0 ==> use Np_matrix or defaults
105  Np_matrix         ! total number concentration (kg^-1).   Optional 0==> use defaults
106
107  integer, parameter :: one = 1
108  ! logical :: hgt_reversed
109  logical :: hgt_descending
110  integer :: pr,i,j,k,unt,ngate
111
112! ----- main program settings ------
113
114  ! Inputs to Quickbeam
115  allocate(hgt_matrix(gbx%Npoints,gbx%Nlevels),p_matrix(gbx%Npoints,gbx%Nlevels), &
116           t_matrix(gbx%Npoints,gbx%Nlevels),rh_matrix(gbx%Npoints,gbx%Nlevels))
117  allocate(hm_matrix(gbx%Nhydro,gbx%Npoints,gbx%Nlevels))
118  allocate(re_matrix(gbx%Nhydro,gbx%Npoints,gbx%Nlevels))
119  allocate(Np_matrix(gbx%Nhydro,gbx%Npoints,gbx%Nlevels))
120
121  ! Outputs from Quickbeam
122  allocate(Ze_non(gbx%Npoints,gbx%Nlevels))
123  allocate(Ze_ray(gbx%Npoints,gbx%Nlevels))
124  allocate(h_atten_to_vol(gbx%Npoints,gbx%Nlevels))
125  allocate(g_atten_to_vol(gbx%Npoints,gbx%Nlevels))
126  allocate(dBZe(gbx%Npoints,gbx%Nlevels))
127
128  ! Optional argument. It is computed and returned in the first call to
129  ! radar_simulator, and passed as input in the rest
130  allocate(g_to_vol(gbx%Npoints,gbx%Nlevels))
131
132  ! Even if there is no unit conversion, they are needed for type conversion
133  p_matrix   = gbx%p/100.0     ! From Pa to hPa
134  hgt_matrix = gbx%zlev/1000.0 ! From m to km
135  t_matrix   = gbx%T
136  rh_matrix  = gbx%q
137  re_matrix  = 0.0
138
139
140  ! set flag denoting position of radar relative to hgt_matrix orientation
141          ngate = size(hgt_matrix,2)
142
143          hgt_descending = hgt_matrix(1,1) > hgt_matrix(1,ngate)
144
145          if ( &
146             (gbx%surface_radar == 1 .and. hgt_descending) .or.  &
147             (gbx%surface_radar == 0 .and. (.not. hgt_descending)) &
148             ) &
149          then
150            gbx%hp%radar_at_layer_one = .false.
151          else
152            gbx%hp%radar_at_layer_one = .true.
153          endif
154
155  ! ----- loop over subcolumns -----
156  do pr=1,sgx%Ncolumns
157
158      !  NOTE:
159      !  atmospheric profiles are the same within the same gridbox
160      !  only hydrometeor profiles will be different for each subgridbox
161
162         do i=1,gbx%Nhydro
163            hm_matrix(i,:,:) = sghydro%mr_hydro(:,pr,:,i)*1000.0 ! Units from kg/kg to g/kg
164            if (gbx%use_reff) then
165              re_matrix(i,:,:) = sghydro%Reff(:,pr,:,i)*1.e6       ! Units from m to micron
166              Np_matrix(i,:,:) = sghydro%Np(:,pr,:,i)              ! Units [#/kg]
167            endif
168         enddo
169
170      !   ----- call radar simulator -----
171      if (pr == 1) then ! Compute gaseous attenuation for all profiles
172         call radar_simulator(gbx%hp,gbx%Npoints,gbx%Nlevels,R_UNDEF, &
173           hgt_matrix,hm_matrix,re_matrix,Np_matrix, &
174           p_matrix,t_matrix,rh_matrix, &
175           Ze_non,Ze_ray,h_atten_to_vol,g_atten_to_vol,dBZe,g_to_vol_out=g_to_vol)
176      else ! Use gaseous atteunuation for pr = 1
177         call radar_simulator(gbx%hp,gbx%Npoints,gbx%Nlevels,R_UNDEF, &
178           hgt_matrix,hm_matrix,re_matrix,Np_matrix, &
179           p_matrix,t_matrix,rh_matrix, &
180           Ze_non,Ze_ray,h_atten_to_vol,g_atten_to_vol,dBZe,g_to_vol_in=g_to_vol)
181      endif
182
183      ! store caluculated dBZe values for later output/processing
184      z%Ze_tot(:,pr,:)=dBZe(:,:)
185  enddo !pr
186
187  deallocate(hgt_matrix,p_matrix,t_matrix,rh_matrix)
188  deallocate(hm_matrix,re_matrix, &
189      Ze_non,Ze_ray,h_atten_to_vol,g_atten_to_vol,dBZe)
190  deallocate(g_to_vol)
191END SUBROUTINE COSP_RADAR
192
193END MODULE MOD_COSP_RADAR
Note: See TracBrowser for help on using the repository browser.