Line data Source code
1 : ! %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2 : ! Copyright (c) 2015, Regents of the University of Colorado
3 : ! All rights reserved.
4 : !
5 : ! Redistribution and use in source and binary forms, with or without modification, are
6 : ! permitted provided that the following conditions are met:
7 : !
8 : ! 1. Redistributions of source code must retain the above copyright notice, this list of
9 : ! conditions and the following disclaimer.
10 : !
11 : ! 2. Redistributions in binary form must reproduce the above copyright notice, this list
12 : ! of conditions and the following disclaimer in the documentation and/or other
13 : ! materials provided with the distribution.
14 : !
15 : ! 3. Neither the name of the copyright holder nor the names of its contributors may be
16 : ! used to endorse or promote products derived from this software without specific prior
17 : ! written permission.
18 : !
19 : ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY
20 : ! EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF
21 : ! MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
22 : ! THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
23 : ! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT
24 : ! OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
25 : ! INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
26 : ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
27 : ! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
28 : !
29 : ! History
30 : ! 11/2005: John Haynes - Created
31 : ! 09/2006 placed into subroutine form (Roger Marchand,JMH)
32 : ! 08/2007 added equivalent volume spheres, Z and N scalling most distrubtion types (Roger Marchand)
33 : ! 01/2008 'Do while' to determine if hydrometeor(s) present in volume
34 : ! changed for vectorization purposes (A. Bodas-Salcedo)
35 : !
36 : ! 07/2010 V3.0 ... Modified to load or save scale factors to disk as a Look-Up Table (LUT)
37 : ! ... All hydrometeor and radar simulator properties now included in hp structure
38 : ! ... hp structure should be initialized by call to radar_simulator_init prior
39 : ! ... to calling this subroutine.
40 : ! Also ... Support of Morrison 2-moment style microphyscis (Np_matrix) added
41 : ! ... Changes implement by Roj Marchand following work by Laura Fowler
42 : !
43 : ! 10/2011 Modified ngate loop to go in either direction depending on flag
44 : ! hp%radar_at_layer_one. This affects the direction in which attenuation is summed.
45 : !
46 : ! Also removed called to AVINT for gas and hydrometeor attenuation and replaced with simple
47 : ! summation. (Roger Marchand)
48 : ! May 2015 - D. Swales - Modified for COSPv2.0
49 : ! %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
50 : module quickbeam
51 : USE COSP_KINDS, ONLY: wp
52 : USE MOD_COSP_CONFIG, ONLY: R_UNDEF,cloudsat_histRef,use_vgrid,vgrid_zl,vgrid_zu,&
53 : pClass_noPrecip, pClass_Rain1, pClass_Rain2, pClass_Rain3,&
54 : pClass_Snow1, pClass_Snow2, pClass_Mixed1, pClass_Mixed2, &
55 : pClass_Rain4, pClass_default, Zenonbinval, Zbinvallnd, &
56 : N_HYDRO,nCloudsatPrecipClass,cloudsat_preclvl
57 :
58 : USE MOD_COSP_STATS, ONLY: COSP_LIDAR_ONLY_CLOUD,hist1D,COSP_CHANGE_VERTICAL_GRID
59 : implicit none
60 :
61 : integer,parameter :: &
62 : maxhclass = 20, & ! Qucikbeam maximum number of hydrometeor classes.
63 : nRe_types = 550, & ! Quickbeam maximum number or Re size bins allowed in N and Z_scaled look up table.
64 : nd = 85, & ! Qucikbeam number of discrete particles used in construction DSDs.
65 : mt_ntt = 39, & ! Quickbeam number of temperatures in mie LUT.
66 : Re_BIN_LENGTH = 10, & ! Quickbeam minimum Re interval in scale LUTs
67 : Re_MAX_BIN = 250 ! Quickbeam maximum Re interval in scale LUTs
68 : real(wp),parameter :: &
69 : dmin = 0.1, & ! Quickbeam minimum size of discrete particle
70 : dmax = 10000. ! Quickbeam maximum size of discrete particle
71 :
72 : !djs logical :: radar_at_layer_one ! If true radar is assume to be at the edge
73 : ! of the first layer, if the first layer is the
74 : ! surface than a ground-based radar. If the
75 : ! first layer is the top-of-atmosphere, then
76 : ! a space borne radar.
77 :
78 : ! ##############################################################################################
79 : type radar_cfg
80 : ! Radar properties
81 : real(wp) :: freq,k2
82 : integer :: nhclass ! Number of hydrometeor classes in use
83 : integer :: use_gas_abs, do_ray
84 : logical :: radar_at_layer_one ! If true radar is assume to be at the edge
85 : ! of the first layer, if the first layer is the
86 : ! surface than a ground-based radar. If the
87 : ! first layer is the top-of-atmosphere, then
88 : ! a space borne radar.
89 :
90 : ! Variables used to store Z scale factors
91 : character(len=240) :: scale_LUT_file_name
92 : logical :: load_scale_LUTs, update_scale_LUTs
93 : logical, dimension(maxhclass,nRe_types) :: N_scale_flag
94 : logical, dimension(maxhclass,mt_ntt,nRe_types) :: Z_scale_flag,Z_scale_added_flag
95 : real(wp),dimension(maxhclass,mt_ntt,nRe_types) :: Ze_scaled,Zr_scaled,kr_scaled
96 : real(wp),dimension(maxhclass,nd,nRe_types) :: fc, rho_eff
97 : real(wp),dimension(Re_MAX_BIN) :: base_list,step_list
98 :
99 : end type radar_cfg
100 :
101 : contains
102 : ! ######################################################################################
103 : ! SUBROUTINE quickbeam_subcolumn
104 : ! ######################################################################################
105 92880 : subroutine quickbeam_subcolumn(rcfg,nprof,ngate,hgt_matrix,z_vol,kr_vol,g_vol,dBZe,Ze_non)
106 :
107 : ! INPUTS
108 : type(radar_cfg),intent(inout) :: &
109 : rcfg ! Derived type for radar simulator setup
110 : integer,intent(in) :: &
111 : nprof, & ! Number of hydrometeor profiles
112 : ngate ! Number of vertical layers
113 : real(wp),intent(in),dimension(nprof,ngate) :: &
114 : hgt_matrix, & ! Height of hydrometeors (km)
115 : z_vol, & ! Effective reflectivity factor (mm^6/m^3)
116 : kr_vol, & ! Attenuation coefficient hydro (dB/km)
117 : g_vol ! Attenuation coefficient gases (dB/km)
118 :
119 : ! OUTPUTS
120 : real(wp), intent(out),dimension(nprof,ngate) :: &
121 : Ze_non, & ! Radar reflectivity without attenuation (dBZ)
122 : dBZe ! Effective radar reflectivity factor (dBZ)
123 :
124 : ! LOCAL VARIABLES
125 : integer :: k,pr,start_gate,end_gate,d_gate
126 : real(wp),dimension(nprof,ngate) :: &
127 185760 : Ze_ray, & ! Rayleigh reflectivity (dBZ)
128 185760 : g_to_vol, & ! Gaseous atteunation, radar to vol (dB)
129 185760 : a_to_vol, & ! Hydromets attenuation, radar to vol (dB)
130 185760 : z_ray ! Reflectivity factor, Rayleigh only (mm^6/m^3)
131 :
132 : ! Load scaling matricies from disk -- but only the first time this subroutine is called
133 92880 : if(rcfg%load_scale_LUTs) then
134 0 : call load_scale_LUTs(rcfg)
135 0 : rcfg%load_scale_LUTs=.false.
136 0 : rcfg%Z_scale_added_flag = .false. ! will be set true if scaling Look Up Tables are modified during run
137 : endif
138 :
139 : ! Initialization
140 130366800 : g_to_vol = 0._wp
141 130366800 : a_to_vol = 0._wp
142 :
143 : ! Loop over each range gate (ngate) ... starting with layer closest to the radar !
144 92880 : if(rcfg%radar_at_layer_one) then
145 : start_gate = 1
146 : end_gate = ngate
147 : d_gate = 1
148 : else
149 0 : start_gate = ngate
150 0 : end_gate = 1
151 0 : d_gate = -1
152 : endif
153 7894800 : do k=start_gate,end_gate,d_gate
154 : ! Loop over each profile (nprof)
155 130366800 : do pr=1,nprof
156 : ! Attenuation due to hydrometeors between radar and volume
157 :
158 : ! NOTE old scheme integrates attenuation only for the layers ABOVE
159 : ! the current layer ... i.e. 1 to k-1 rather than 1 to k ...
160 : ! which may be a problem. ROJ
161 : ! in the new scheme I assign half the attenuation to the current layer
162 122472000 : if(d_gate==1) then
163 : ! dheight calcuations assumes hgt_matrix points are the cell mid-points.
164 122472000 : if (k>2) then
165 : ! add to previous value to half of above layer + half of current layer
166 239112000 : a_to_vol(pr,k)= a_to_vol(pr,k-1) + &
167 358668000 : (kr_vol(pr,k-1)+kr_vol(pr,k))*(hgt_matrix(pr,k-1)-hgt_matrix(pr,k))
168 : else
169 2916000 : a_to_vol(pr,k)= kr_vol(pr,k)*(hgt_matrix(pr,k)-hgt_matrix(pr,k+1))
170 : endif
171 : else ! d_gate==-1
172 0 : if(k<ngate) then
173 : ! Add to previous value half of above layer + half of current layer
174 0 : a_to_vol(pr,k) = a_to_vol(pr,k+1) + &
175 0 : (kr_vol(pr,k+1)+kr_vol(pr,k))*(hgt_matrix(pr,k+1)-hgt_matrix(pr,k))
176 : else
177 0 : a_to_vol(pr,k)= kr_vol(pr,k)*(hgt_matrix(pr,k)-hgt_matrix(pr,k-1))
178 : endif
179 : endif
180 :
181 : ! Attenuation due to gaseous absorption between radar and volume
182 130273920 : if ((rcfg%use_gas_abs == 1) .or. (rcfg%use_gas_abs == 2 .and. pr .eq. 1)) then
183 122472000 : if (d_gate==1) then
184 122472000 : if (k>1) then
185 : ! Add to previous value to half of above layer + half of current layer
186 242028000 : g_to_vol(pr,k) = g_to_vol(pr,k-1) + &
187 363042000 : 0.5*(g_vol(pr,k-1)+g_vol(pr,k))*(hgt_matrix(pr,k-1)-hgt_matrix(pr,k))
188 : else
189 1458000 : g_to_vol(pr,k)= 0.5_wp*g_vol(pr,k)*(hgt_matrix(pr,k)-hgt_matrix(pr,k+1))
190 : endif
191 : else ! d_gate==-1
192 0 : if (k<ngate) then
193 : ! Add to previous value to half of above layer + half of current layer
194 0 : g_to_vol(pr,k) = g_to_vol(pr,k+1) + &
195 0 : 0.5_wp*(g_vol(pr,k+1)+g_vol(pr,k))*(hgt_matrix(pr,k+1)-hgt_matrix(pr,k))
196 : else
197 0 : g_to_vol(pr,k)= 0.5_wp*g_vol(pr,k)*(hgt_matrix(pr,k)-hgt_matrix(pr,k-1))
198 : endif
199 : endif
200 0 : elseif(rcfg%use_gas_abs == 2) then
201 : ! Using value calculated for the first column
202 0 : g_to_vol(pr,k) = g_to_vol(1,k)
203 0 : elseif (rcfg%use_gas_abs == 0) then
204 0 : g_to_vol(pr,k) = 0._wp
205 : endif
206 : enddo ! End loop over pr (profile)
207 : enddo ! End loop of k (range gate)
208 :
209 : ! Compute Rayleigh reflectivity, and full, attenuated reflectivity
210 92880 : if(rcfg%do_ray == 1) then
211 0 : where(z_ray(1:nprof,1:ngate) > 0._wp)
212 : Ze_ray(1:nprof,1:ngate) = 10._wp*log10(z_ray(1:nprof,1:ngate))
213 : elsewhere
214 : Ze_Ray(1:nprof,1:ngate) = 0._wp
215 : endwhere
216 : else
217 130366800 : Ze_ray(1:nprof,1:ngate) = R_UNDEF
218 : end if
219 :
220 651462480 : where(z_vol(1:nprof,1:ngate) > 0._wp)
221 : Ze_non(1:nprof,1:ngate) = 10._wp*log10(z_vol(1:nprof,1:ngate))
222 : dBZe(1:nprof,1:ngate) = Ze_non(1:nprof,1:ngate)-a_to_vol(1:nprof,1:ngate)-g_to_vol(1:nprof,1:ngate)
223 : elsewhere
224 : dBZe(1:nprof,1:ngate) = R_UNDEF
225 : Ze_non(1:nprof,1:ngate) = R_UNDEF
226 : end where
227 :
228 : ! Save any updates made
229 92880 : if (rcfg%update_scale_LUTs) call save_scale_LUTs(rcfg)
230 :
231 92880 : end subroutine quickbeam_subcolumn
232 : ! ######################################################################################
233 : ! SUBROUTINE quickbeam_column
234 : ! ######################################################################################
235 9288 : subroutine quickbeam_column(npoints, ncolumns, nlevels, llm, DBZE_BINS, platform, &
236 9288 : Ze_tot, Ze_tot_non, land, surfelev, t2m, fracPrecipIce, zlev, zlev_half, cfad_ze, &
237 9288 : cloudsat_precip_cover, cloudsat_pia)
238 : ! Inputs
239 : integer,intent(in) :: &
240 : npoints, & ! Number of horizontal grid points
241 : ncolumns, & ! Number of subcolumns
242 : nlevels, & ! Number of vertical layers in OLD grid
243 : llm, & ! Number of vertical layers in NEW grid
244 : DBZE_BINS ! Number of bins for cfad.
245 : character(len=*),intent(in) :: &
246 : platform ! Name of platform (e.g. cloudsat)
247 : real(wp),dimension(Npoints),intent(in) :: &
248 : land, & ! Land/Sea mask. (1/0)
249 : surfelev, & ! Surface Elevation (m)
250 : t2m ! Near-surface temperature
251 : real(wp),dimension(Npoints,Ncolumns),intent(in) :: &
252 : fracPrecipIce ! Fraction of precipitation which is frozen. (1)
253 : real(wp),intent(in),dimension(npoints,ncolumns,Nlevels) :: &
254 : Ze_tot, & ! Effective reflectivity factor (dBZ)
255 : Ze_tot_non ! Effective reflectivity factor w/o attenuation (dBZ)
256 : real(wp),intent(in),dimension(npoints,Nlevels) :: &
257 : zlev ! Model full levels
258 : real(wp),intent(in),dimension(npoints,Nlevels+1) :: &
259 : zlev_half ! Model half levels
260 :
261 : ! Outputs
262 : real(wp),intent(inout),dimension(npoints,DBZE_BINS,llm) :: &
263 : cfad_ze !
264 : real(wp),dimension(Npoints,nCloudsatPrecipClass),intent(out) :: &
265 : cloudsat_precip_cover ! Model precip rate in by CloudSat precip flag
266 : real(wp),dimension(Npoints),intent(out) :: &
267 : cloudsat_pia ! Cloudsat path integrated attenuation
268 :
269 : ! Local variables
270 : integer :: i,j
271 : real(wp) :: zstep
272 18576 : real(wp),dimension(npoints,ncolumns,llm) :: ze_toti,ze_noni
273 : logical :: lcloudsat = .false.
274 :
275 : ! Which platforms to create diagnostics for?
276 9288 : if (platform .eq. 'cloudsat') lcloudsat=.true.
277 :
278 : ! Create Cloudsat diagnostics.
279 9288 : if (lcloudsat) then
280 9288 : if (use_vgrid) then
281 : ! Regrid in the vertical (*NOTE* This routine requires SFC-2-TOA ordering, so flip
282 : ! inputs and outputs to maintain TOA-2-SFC ordering convention in COSP2.)
283 : call cosp_change_vertical_grid(Npoints,Ncolumns,Nlevels,zlev(:,nlevels:1:-1),&
284 9288 : zlev_half(:,nlevels:1:-1),Ze_tot(:,:,nlevels:1:-1),llm,vgrid_zl(llm:1:-1),&
285 220277232 : vgrid_zu(llm:1:-1),Ze_toti(:,:,llm:1:-1),log_units=.true.)
286 :
287 : ! Effective reflectivity histogram
288 155088 : do i=1,Npoints
289 5987088 : do j=1,llm
290 151777800 : cfad_ze(i,:,j) = hist1D(Ncolumns,Ze_toti(i,:,j),DBZE_BINS,cloudsat_histRef)
291 : enddo
292 : enddo
293 93433608 : where(cfad_ze .ne. R_UNDEF) cfad_ze = cfad_ze/Ncolumns
294 :
295 : ! Compute cloudsat near-surface precipitation diagnostics
296 : ! First, regrid in the vertical Ze_tot_non.
297 : call cosp_change_vertical_grid(Npoints,Ncolumns,Nlevels,zlev(:,nlevels:1:-1),&
298 0 : zlev_half(:,nlevels:1:-1),Ze_tot_non(:,:,nlevels:1:-1),llm,vgrid_zl(llm:1:-1),&
299 220267944 : vgrid_zu(llm:1:-1),Ze_noni(:,:,llm:1:-1),log_units=.true.)
300 : ! Compute the zstep distance between two atmopsheric layers
301 9288 : zstep = vgrid_zl(1)-vgrid_zl(2)
302 : ! Now call routine to generate diagnostics.
303 : call cloudsat_precipOccurence(Npoints, Ncolumns, llm, N_HYDRO, Ze_toti, Ze_noni, &
304 9288 : land, surfelev, t2m, fracPrecipIce, cloudsat_precip_cover, cloudsat_pia, zstep)
305 : else
306 : ! Effective reflectivity histogram
307 0 : do i=1,Npoints
308 0 : do j=1,llm
309 0 : cfad_ze(i,:,j) = hist1D(Ncolumns,Ze_tot(i,:,j),DBZE_BINS,cloudsat_histRef)
310 : enddo
311 : enddo
312 0 : where(cfad_ze .ne. R_UNDEF) cfad_ze = cfad_ze/Ncolumns
313 : endif
314 : endif
315 :
316 9288 : end subroutine quickbeam_column
317 : ! ##############################################################################################
318 : ! ##############################################################################################
319 : ! ######################################################################################
320 : ! SUBROUTINE cloudsat_precipOccurence
321 : !
322 : ! Notes from July 2016: Add precip flag also looped over subcolumns
323 : ! Modified by Tristan L'Ecuyer (TSL) to add precipitation flagging
324 : ! based on CloudSat's 2C-PRECIP-COLUMN algorithm (Haynes et al, JGR, 2009).
325 : ! To mimic the satellite algorithm, this code applies thresholds to non-attenuated
326 : ! reflectivities, Ze_non, consistent with those outlined in Haynes et al, JGR (2009).
327 : !
328 : ! Procedures/Notes:
329 : !
330 : ! (1) If the 2-way attenuation exceeds 40 dB, the pixel will be flagged as 'heavy rain'
331 : ! consistent with the multiple-scattering analysis of Battaglia et al, JGR (2008).
332 : ! (2) Rain, snow, and mixed precipitation scenes are separated according to the fraction
333 : ! of the total precipitation hydrometeor mixing ratio that exists as ice.
334 : ! (3) The entire analysis is applied to the range gate from 480-960 m to be consistent with
335 : ! CloudSat's 720 m ground-clutter.
336 : ! (4) Only a single flag is assigned to each model grid since there is no variation in
337 : ! hydrometeor contents across a single model level. Unlike CFADs, whose variation enters
338 : ! due to differening attenuation corrections from hydrometeors aloft, the non-attenuated
339 : ! reflectivities used in the computation of this flag cannot vary across sub-columns.
340 : !
341 : ! radar_prec_flag = 1-Rain possible 2-Rain probable 3-Rain certain
342 : ! 4-Snow possible 5-Snow certain
343 : ! 6-Mixed possible 7-Mixed certain
344 : ! 8-Heavy Rain
345 : ! 9- default value
346 : !
347 : ! Modified by Dustin Swales (University of Colorado) for use with COSP2.
348 : ! *NOTE* All inputs (Ze_out, Ze_non_out, fracPrecipIce) are at a single level from the
349 : ! statistical output grid used by Cloudsat. This level is controlled by the
350 : ! parameter cloudsat_preclvl, defined in src/cosp_config.F90
351 : ! ######################################################################################
352 9288 : subroutine cloudsat_precipOccurence(Npoints, Ncolumns, llm, Nhydro, Ze_out, Ze_non_out, &
353 9288 : land, surfelev, t2m, fracPrecipIce, cloudsat_precip_cover, cloudsat_pia, zstep)
354 :
355 : ! Inputs
356 : integer,intent(in) :: &
357 : Npoints, & ! Number of columns
358 : Ncolumns, & ! Numner of subcolumns
359 : Nhydro, & ! Number of hydrometeor types
360 : llm ! Number of levels
361 : real(wp),dimension(Npoints),intent(in) :: &
362 : land, & ! Land/Sea mask. (1/0)
363 : surfelev, & ! Surface Elevation (m)
364 : t2m ! Near-surface temperature
365 : real(wp),dimension(Npoints,Ncolumns,llm),intent(in) :: &
366 : Ze_out, & ! Effective reflectivity factor (dBZ)
367 : Ze_non_out ! Effective reflectivity factor, w/o attenuation (dBZ)
368 : real(wp),dimension(Npoints,Ncolumns),intent(in) :: &
369 : fracPrecipIce ! Fraction of precipitation which is frozen. (1)
370 : real(wp),intent(in) :: &
371 : zstep ! Distance between two atmopsheric layers (m)
372 :
373 : ! Outputs
374 : real(wp),dimension(Npoints,nCloudsatPrecipClass),intent(out) :: &
375 : cloudsat_precip_cover ! Model precip rate in by CloudSat precip flag
376 : real(wp),dimension(Npoints),intent(out) :: &
377 : cloudsat_pia ! Cloudsat path integrated attenuation
378 :
379 : ! Local variables
380 : integer,dimension(Npoints,Ncolumns) :: &
381 18576 : cloudsat_pflag, & ! Subcolumn precipitation flag
382 18576 : cloudsat_precip_pia ! Subcolumn path integrated attenutation.
383 : integer,dimension(Npoints) :: &
384 9288 : cloudsat_preclvl_index ! Altitude index for precip flags calculation
385 : ! in 40-level grid (one layer above surfelev)
386 : integer :: pr,i,k,m,j
387 : real(wp) :: Zmax
388 :
389 : ! Initialize
390 1560168 : cloudsat_pflag(:,:) = pClass_default
391 1560168 : cloudsat_precip_pia(:,:) = 0._wp
392 1560168 : cloudsat_precip_cover(:,:) = 0._wp
393 155088 : cloudsat_pia(:) = 0._wp
394 155088 : cloudsat_preclvl_index(:) = 0._wp
395 :
396 : ! Computing altitude index for precip flags calculation
397 155088 : cloudsat_preclvl_index(:) = cloudsat_preclvl - floor( surfelev(:)/zstep )
398 :
399 : ! ######################################################################################
400 : ! SUBCOLUMN processing
401 : ! ######################################################################################
402 155088 : do i=1, Npoints
403 1613088 : do pr=1,Ncolumns
404 : ! Compute precipitation flag
405 : ! ################################################################################
406 : ! 1) Oceanic points.
407 : ! ################################################################################
408 1458000 : if (land(i) .eq. 0) then
409 :
410 : ! 1a) Compute the PIA in all profiles containing hydrometeors
411 983910 : if ( (Ze_non_out(i,pr,cloudsat_preclvl_index(i)).gt.-100) .and. (Ze_out(i,pr,cloudsat_preclvl_index(i)).gt.-100) ) then
412 439746 : if ( (Ze_non_out(i,pr,cloudsat_preclvl_index(i)).lt.100) .and. (Ze_out(i,pr,cloudsat_preclvl_index(i)).lt.100) ) then
413 439746 : cloudsat_precip_pia(i,pr) = Ze_non_out(i,pr,cloudsat_preclvl_index(i)) - Ze_out(i,pr,cloudsat_preclvl_index(i))
414 : endif
415 : endif
416 :
417 : ! Snow
418 983910 : if(fracPrecipIce(i,pr).gt.0.9) then
419 111533 : if(Ze_non_out(i,pr,cloudsat_preclvl_index(i)).gt.Zenonbinval(2)) then
420 4862 : cloudsat_pflag(i,pr) = pClass_Snow2 ! TSL: Snow certain
421 : endif
422 111533 : if(Ze_non_out(i,pr,cloudsat_preclvl_index(i)).gt.Zenonbinval(4).and. &
423 : Ze_non_out(i,pr,cloudsat_preclvl_index(i)).le.Zenonbinval(2)) then
424 15267 : cloudsat_pflag(i,pr) = pClass_Snow1 ! TSL: Snow possible
425 : endif
426 : endif
427 :
428 : ! Mixed
429 983910 : if(fracPrecipIce(i,pr).gt.0.1.and.fracPrecipIce(i,pr).le.0.9) then
430 78191 : if(Ze_non_out(i,pr,cloudsat_preclvl_index(i)).gt.Zenonbinval(2)) then
431 7013 : cloudsat_pflag(i,pr) = pClass_Mixed2 ! TSL: Mixed certain
432 : endif
433 78191 : if(Ze_non_out(i,pr,cloudsat_preclvl_index(i)).gt.Zenonbinval(4).and. &
434 : Ze_non_out(i,pr,cloudsat_preclvl_index(i)).le.Zenonbinval(2)) then
435 17738 : cloudsat_pflag(i,pr) = pClass_Mixed1 ! TSL: Mixed possible
436 : endif
437 : endif
438 :
439 : ! Rain
440 983910 : if(fracPrecipIce(i,pr).le.0.1) then
441 794186 : if(Ze_non_out(i,pr,cloudsat_preclvl_index(i)).gt.Zenonbinval(1)) then
442 29377 : cloudsat_pflag(i,pr) = pClass_Rain3 ! TSL: Rain certain
443 : endif
444 794186 : if(Ze_non_out(i,pr,cloudsat_preclvl_index(i)).gt.Zenonbinval(3).and. &
445 : Ze_non_out(i,pr,cloudsat_preclvl_index(i)).le.Zenonbinval(1)) then
446 4676 : cloudsat_pflag(i,pr) = pClass_Rain2 ! TSL: Rain probable
447 : endif
448 794186 : if(Ze_non_out(i,pr,cloudsat_preclvl_index(i)).gt.Zenonbinval(4).and. &
449 : Ze_non_out(i,pr,cloudsat_preclvl_index(i)).le.Zenonbinval(3)) then
450 17522 : cloudsat_pflag(i,pr) = pClass_Rain1 ! TSL: Rain possible
451 : endif
452 794186 : if(cloudsat_precip_pia(i,pr).gt.40) then
453 351 : cloudsat_pflag(i,pr) = pClass_Rain4 ! TSL: Heavy Rain
454 : endif
455 : endif
456 :
457 : ! No precipitation
458 983910 : if(Ze_non_out(i,pr,cloudsat_preclvl_index(i)).le.-15) then
459 887455 : cloudsat_pflag(i,pr) = pClass_noPrecip ! TSL: Not Raining
460 : endif
461 : endif ! Ocean points
462 :
463 : ! ################################################################################
464 : ! 2) Land points.
465 : ! *NOTE* For land points we go up a layer higher, so cloudsat_preclvl_index(i)-1
466 : !
467 : ! ################################################################################
468 1603800 : if (land(i) .eq. 1) then
469 : ! 2a) Compute the PIA in all profiles containing hydrometeors
470 474090 : if ( (Ze_non_out(i,pr,cloudsat_preclvl_index(i)-1).gt.-100) .and. (Ze_out(i,pr,cloudsat_preclvl_index(i)-1).gt.-100) ) then
471 192650 : if ( (Ze_non_out(i,pr,cloudsat_preclvl_index(i)-1).lt.100) .and. (Ze_out(i,pr,cloudsat_preclvl_index(i)-1).lt.100) ) then
472 192650 : cloudsat_precip_pia(i,pr) = Ze_non_out(i,pr,cloudsat_preclvl_index(i)-1) - Ze_out(i,pr,cloudsat_preclvl_index(i)-1)
473 : endif
474 : endif
475 :
476 : ! Find Zmax, the maximum reflectivity value in the attenuated profile (Ze_out);
477 19911780 : Zmax=maxval(Ze_out(i,pr,:))
478 :
479 : ! Snow (T<273)
480 474090 : if(t2m(i) .lt. 273._wp) then
481 194080 : if(Ze_out(i,pr,cloudsat_preclvl_index(i)-1) .gt. Zbinvallnd(5)) then
482 129 : cloudsat_pflag(i,pr) = pClass_Snow2 ! JEK: Snow certain
483 : endif
484 194080 : if(Ze_out(i,pr,cloudsat_preclvl_index(i)-1) .gt. Zbinvallnd(6) .and. &
485 : Ze_out(i,pr,cloudsat_preclvl_index(i)-1).le.Zbinvallnd(5)) then
486 7407 : cloudsat_pflag(i,pr) = pClass_Snow1 ! JEK: Snow possible
487 : endif
488 : endif
489 :
490 : ! Mized phase (273<T<275)
491 474090 : if(t2m(i) .ge. 273._wp .and. t2m(i) .le. 275._wp) then
492 12120 : if ((Zmax .gt. Zbinvallnd(1) .and. cloudsat_precip_pia(i,pr).gt.30) .or. &
493 12120 : (Ze_out(i,pr,cloudsat_preclvl_index(i)-1) .gt. Zbinvallnd(4))) then
494 81 : cloudsat_pflag(i,pr) = pClass_Mixed2 ! JEK: Mixed certain
495 : endif
496 12120 : if ((Ze_out(i,pr,cloudsat_preclvl_index(i)-1) .gt. Zbinvallnd(6) .and. &
497 24240 : Ze_out(i,pr,cloudsat_preclvl_index(i)-1) .le. Zbinvallnd(4)) .and. &
498 : (Zmax .gt. Zbinvallnd(5)) ) then
499 57 : cloudsat_pflag(i,pr) = pClass_Mixed1 ! JEK: Mixed possible
500 : endif
501 : endif
502 :
503 : ! Rain (T>275)
504 474090 : if(t2m(i) .gt. 275) then
505 267890 : if ((Zmax .gt. Zbinvallnd(1) .and. cloudsat_precip_pia(i,pr).gt.30) .or. &
506 267890 : (Ze_out(i,pr,cloudsat_preclvl_index(i)-1) .gt. Zbinvallnd(2))) then
507 5383 : cloudsat_pflag(i,pr) = pClass_Rain3 ! JEK: Rain certain
508 : endif
509 267890 : if((Ze_out(i,pr,cloudsat_preclvl_index(i)-1) .gt. Zbinvallnd(6)) .and. &
510 : (Zmax .gt. Zbinvallnd(3))) then
511 8466 : cloudsat_pflag(i,pr) = pClass_Rain2 ! JEK: Rain probable
512 : endif
513 267890 : if((Ze_out(i,pr,cloudsat_preclvl_index(i)-1) .gt. Zbinvallnd(6)) .and. &
514 : (Zmax.lt.Zbinvallnd(3))) then
515 4292 : cloudsat_pflag(i,pr) = pClass_Rain1 ! JEK: Rain possible
516 : endif
517 267890 : if(cloudsat_precip_pia(i,pr).gt.40) then
518 257 : cloudsat_pflag(i,pr) = pClass_Rain4 ! JEK: Heavy Rain
519 : endif
520 : endif
521 :
522 : ! No precipitation
523 474090 : if(Ze_out(i,pr,cloudsat_preclvl_index(i)-1) .le. -15) then
524 452589 : cloudsat_pflag(i,pr) = pClass_noPrecip ! JEK: Not Precipitating
525 : endif
526 : endif ! Land points
527 : enddo ! Sub-columns
528 : enddo ! Gridpoints
529 :
530 : ! ######################################################################################
531 : ! COLUMN processing
532 : ! ######################################################################################
533 :
534 : ! Aggregate subcolumns
535 155088 : do i=1,Npoints
536 : ! Gridmean precipitation fraction for each precipitation type
537 1603800 : do k=1,nCloudsatPrecipClass
538 14536083 : if (any(cloudsat_pflag(i,:) .eq. k-1)) then
539 2043504 : cloudsat_precip_cover(i,k) = count(cloudsat_pflag(i,:) .eq. k-1)
540 : endif
541 : enddo
542 :
543 : ! Gridmean path integrated attenuation (pia)
544 1613088 : cloudsat_pia(i)=sum(cloudsat_precip_pia(i,:))
545 : enddo
546 :
547 : ! Normalize by number of subcolumns
548 : where ((cloudsat_precip_cover /= R_UNDEF).and.(cloudsat_precip_cover /= 0.0)) &
549 1560168 : cloudsat_precip_cover = cloudsat_precip_cover / Ncolumns
550 : where ((cloudsat_pia/= R_UNDEF).and.(cloudsat_pia/= 0.0)) &
551 155088 : cloudsat_pia = cloudsat_pia / Ncolumns
552 :
553 9288 : end subroutine cloudsat_precipOccurence
554 :
555 : ! ##############################################################################################
556 : ! ##############################################################################################
557 0 : subroutine load_scale_LUTs(rcfg)
558 :
559 : type(radar_cfg), intent(inout) :: rcfg
560 : logical :: LUT_file_exists
561 : integer :: i,j,k,ind
562 :
563 : ! Load scale LUT from file
564 : inquire(file=trim(rcfg%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat', &
565 0 : exist=LUT_file_exists)
566 :
567 0 : if(.not.LUT_file_exists) then
568 0 : write(*,*) '*************************************************'
569 0 : write(*,*) 'Warning: Could NOT FIND radar LUT file: ', &
570 0 : trim(rcfg%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat'
571 0 : write(*,*) 'Will calculated LUT values as needed'
572 0 : write(*,*) '*************************************************'
573 0 : return
574 : else
575 : OPEN(unit=12,file=trim(rcfg%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat',&
576 : form='unformatted', &
577 : err= 89, &
578 : access='DIRECT',&
579 0 : recl=28)
580 0 : write(*,*) 'Loading radar LUT file: ', &
581 0 : trim(rcfg%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat'
582 :
583 0 : do i=1,maxhclass
584 0 : do j=1,mt_ntt
585 0 : do k=1,nRe_types
586 0 : ind = i+(j-1)*maxhclass+(k-1)*(nRe_types*mt_ntt)
587 0 : read(12,rec=ind) rcfg%Z_scale_flag(i,j,k), &
588 0 : rcfg%Ze_scaled(i,j,k), &
589 0 : rcfg%Zr_scaled(i,j,k), &
590 0 : rcfg%kr_scaled(i,j,k)
591 : enddo
592 : enddo
593 : enddo
594 0 : close(unit=12)
595 0 : return
596 : endif
597 :
598 0 : 89 write(*,*) 'Error: Found but could NOT READ radar LUT file: ', &
599 0 : trim(rcfg%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat'
600 :
601 : end subroutine load_scale_LUTs
602 :
603 : ! ##############################################################################################
604 : ! ##############################################################################################
605 0 : subroutine save_scale_LUTs(rcfg)
606 : type(radar_cfg), intent(inout) :: rcfg
607 : logical :: LUT_file_exists
608 : integer :: i,j,k,ind
609 :
610 : inquire(file=trim(rcfg%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat', &
611 0 : exist=LUT_file_exists)
612 :
613 : OPEN(unit=12,file=trim(rcfg%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat',&
614 0 : form='unformatted',err= 99,access='DIRECT',recl=28)
615 :
616 0 : write(*,*) 'Creating or Updating radar LUT file: ', &
617 0 : trim(rcfg%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat'
618 :
619 0 : do i=1,maxhclass
620 0 : do j=1,mt_ntt
621 0 : do k=1,nRe_types
622 0 : ind = i+(j-1)*maxhclass+(k-1)*(nRe_types*mt_ntt)
623 0 : if(.not.LUT_file_exists .or. rcfg%Z_scale_added_flag(i,j,k)) then
624 0 : rcfg%Z_scale_added_flag(i,j,k)=.false.
625 0 : write(12,rec=ind) rcfg%Z_scale_flag(i,j,k), &
626 0 : rcfg%Ze_scaled(i,j,k), &
627 0 : rcfg%Zr_scaled(i,j,k), &
628 0 : rcfg%kr_scaled(i,j,k)
629 : endif
630 : enddo
631 : enddo
632 : enddo
633 0 : close(unit=12)
634 0 : return
635 :
636 0 : 99 write(*,*) 'Error: Unable to create/update radar LUT file: ', &
637 0 : trim(rcfg%scale_LUT_file_name) // '_radar_Z_scale_LUT.dat'
638 0 : return
639 :
640 : end subroutine save_scale_LUTs
641 : ! ##############################################################################################
642 : ! ##############################################################################################
643 0 : subroutine quickbeam_init()
644 :
645 :
646 0 : end subroutine quickBeam_init
647 : ! ##############################################################################################
648 : ! ##############################################################################################
649 :
650 :
651 0 : end module quickbeam
652 :
653 :
|