[3491] | 1 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 2 | ! |
---|
| 3 | ! Module pour remplir la variable "cospIN", calculs des proprietes optiques subcolumn |
---|
| 4 | ! pour les differents simulateurs |
---|
| 5 | ! |
---|
| 6 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 7 | |
---|
| 8 | MODULE LMDZ_COSP_SUBSAMPLE_AND_OPTICS_MOD |
---|
| 9 | |
---|
| 10 | use cosp_kinds, only: wp |
---|
| 11 | USE MOD_COSP_CONFIG |
---|
| 12 | USE mod_quickbeam_optics, only: size_distribution,quickbeam_optics,gases |
---|
| 13 | use quickbeam, only: radar_cfg |
---|
| 14 | use MOD_COSP, only: cosp_optical_inputs,cosp_column_inputs |
---|
| 15 | USE mod_rng, ONLY: rng_state, init_rng |
---|
| 16 | USE mod_scops, ONLY: scops |
---|
| 17 | USE mod_prec_scops, ONLY: prec_scops |
---|
| 18 | USE MOD_COSP_UTILS, ONLY: cosp_precip_mxratio |
---|
| 19 | use cosp_optics, ONLY: cosp_simulator_optics,lidar_optics,modis_optics, & |
---|
| 20 | modis_optics_partition |
---|
| 21 | use mod_cosp_stats, ONLY: COSP_CHANGE_VERTICAL_GRID |
---|
| 22 | use lmdz_cosp_read_outputkeys, ONLY: cosp_config |
---|
| 23 | |
---|
| 24 | |
---|
| 25 | implicit none |
---|
| 26 | |
---|
| 27 | ! Indices to address arrays of LS and CONV hydrometeors |
---|
| 28 | integer,parameter :: & |
---|
| 29 | I_LSCLIQ = 1, & ! Large-scale (stratiform) liquid |
---|
| 30 | I_LSCICE = 2, & ! Large-scale (stratiform) ice |
---|
| 31 | I_LSRAIN = 3, & ! Large-scale (stratiform) rain |
---|
| 32 | I_LSSNOW = 4, & ! Large-scale (stratiform) snow |
---|
| 33 | I_CVCLIQ = 5, & ! Convective liquid |
---|
| 34 | I_CVCICE = 6, & ! Convective ice |
---|
| 35 | I_CVRAIN = 7, & ! Convective rain |
---|
| 36 | I_CVSNOW = 8, & ! Convective snow |
---|
| 37 | I_LSGRPL = 9 ! Large-scale (stratiform) groupel |
---|
| 38 | |
---|
| 39 | ! Stratiform and convective clouds in frac_out. |
---|
| 40 | integer, parameter :: & |
---|
| 41 | I_LSC = 1, & ! Large-scale clouds |
---|
| 42 | I_CVC = 2 ! Convective clouds |
---|
| 43 | |
---|
| 44 | ! Microphysical settings for the precipitation flux to mixing ratio conversion |
---|
| 45 | real(wp),parameter,dimension(N_HYDRO) :: & |
---|
| 46 | ! LSL LSI LSR LSS CVL CVI CVR CVS LSG |
---|
| 47 | N_ax = (/-1., -1., 8.e6, 3.e6, -1., -1., 8.e6, 3.e6, 4.e6/),& |
---|
| 48 | N_bx = (/-1., -1., 0.0, 0.0, -1., -1., 0.0, 0.0, 0.0/),& |
---|
| 49 | alpha_x = (/-1., -1., 0.0, 0.0, -1., -1., 0.0, 0.0, 0.0/),& |
---|
| 50 | c_x = (/-1., -1., 842.0, 4.84, -1., -1., 842.0, 4.84, 94.5/),& |
---|
| 51 | d_x = (/-1., -1., 0.8, 0.25, -1., -1., 0.8, 0.25, 0.5/),& |
---|
| 52 | g_x = (/-1., -1., 0.5, 0.5, -1., -1., 0.5, 0.5, 0.5/),& |
---|
| 53 | a_x = (/-1., -1., 524.0, 52.36, -1., -1., 524.0, 52.36, 209.44/),& |
---|
| 54 | b_x = (/-1., -1., 3.0, 3.0, -1., -1., 3.0, 3.0, 3.0/),& |
---|
| 55 | gamma_1 = (/-1., -1., 17.83725, 8.284701, -1., -1., 17.83725, 8.284701, 11.63230/),& |
---|
| 56 | gamma_2 = (/-1., -1., 6.0, 6.0, -1., -1., 6.0, 6.0, 6.0/),& |
---|
| 57 | gamma_3 = (/-1., -1., 2.0, 2.0, -1., -1., 2.0, 2.0, 2.0/),& |
---|
| 58 | gamma_4 = (/-1., -1., 6.0, 6.0, -1., -1., 6.0, 6.0, 6.0/) |
---|
| 59 | |
---|
| 60 | |
---|
| 61 | contains |
---|
| 62 | |
---|
| 63 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 64 | ! SUBROUTINE subsample_and_optics |
---|
| 65 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 66 | subroutine subsample_and_optics(cfg, nPoints, nLevels, nColumns, nHydro, overlap, & |
---|
| 67 | use_precipitation_fluxes, lidar_ice_type, sd, tca, cca, fl_lsrainIN, fl_lssnowIN, & |
---|
| 68 | fl_lsgrplIN, fl_ccrainIN, fl_ccsnowIN, mr_lsliq, mr_lsice, mr_ccliq, mr_ccice, & |
---|
| 69 | reffIN, dtau_c, dtau_s, dem_c, dem_s, cospstateIN, cospIN) |
---|
| 70 | ! Inputs |
---|
| 71 | integer,intent(in) :: nPoints, nLevels, nColumns, nHydro, overlap, lidar_ice_type |
---|
| 72 | real(wp),intent(in),dimension(nPoints,nLevels) :: tca,cca,mr_lsliq,mr_lsice,mr_ccliq, & |
---|
| 73 | mr_ccice,dtau_c,dtau_s,dem_c,dem_s,fl_lsrainIN,fl_lssnowIN,fl_lsgrplIN,fl_ccrainIN,& |
---|
| 74 | fl_ccsnowIN |
---|
| 75 | real(wp),intent(in),dimension(nPoints,nLevels,nHydro) :: reffIN |
---|
| 76 | logical,intent(in) :: use_precipitation_fluxes |
---|
| 77 | type(size_distribution),intent(inout) :: sd |
---|
| 78 | type(cosp_config),intent(in) :: cfg ! Configuration options |
---|
| 79 | |
---|
| 80 | ! Outputs |
---|
| 81 | type(cosp_optical_inputs),intent(inout) :: cospIN |
---|
| 82 | type(cosp_column_inputs),intent(inout) :: cospstateIN |
---|
| 83 | |
---|
| 84 | type(rng_state),allocatable,dimension(:) :: rngs ! Seeds for random number generator |
---|
| 85 | integer,dimension(:),allocatable :: seed |
---|
| 86 | integer,dimension(:),allocatable :: cloudsat_preclvl_index !PREC_BUG |
---|
| 87 | integer :: i,j,k |
---|
| 88 | real(wp),dimension(:,:), allocatable :: & |
---|
| 89 | ls_p_rate, cv_p_rate, frac_ls, frac_cv, prec_ls, prec_cv,g_vol |
---|
| 90 | real(wp),dimension(:,:,:), allocatable :: & |
---|
| 91 | frac_prec, MODIS_cloudWater, MODIS_cloudIce, fracPrecipIce, fracPrecipIce_statGrid,& |
---|
| 92 | MODIS_watersize,MODIS_iceSize, MODIS_opticalThicknessLiq,MODIS_opticalThicknessIce |
---|
| 93 | real(wp),dimension(:,:,:,:),allocatable :: & |
---|
| 94 | mr_hydro, Reff, Np |
---|
| 95 | real(wp),dimension(nPoints,nLevels) :: & |
---|
| 96 | column_frac_out, column_prec_out, fl_lsrain, fl_lssnow, fl_lsgrpl, fl_ccrain, fl_ccsnow |
---|
| 97 | ! real(wp),dimension(nPoints,nColumns,Nlvgrid_local) :: tempOut |
---|
| 98 | logical :: cmpGases=.true. |
---|
| 99 | |
---|
| 100 | if (Ncolumns .gt. 1) then |
---|
| 101 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 102 | ! Generate subcolumns for clouds (SCOPS) and precipitation type (PREC_SCOPS) |
---|
| 103 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 104 | ! RNG used for subcolumn generation |
---|
| 105 | allocate(rngs(nPoints),seed(nPoints)) |
---|
| 106 | seed(:)=0 |
---|
| 107 | seed = int(cospstateIN%phalf(:,Nlevels+1)) ! In case of NPoints=1 |
---|
| 108 | ! *NOTE* Chunking will change the seed |
---|
| 109 | if (NPoints .gt. 1) seed=int((cospstateIN%phalf(:,Nlevels+1)-minval(cospstateIN%phalf(:,Nlevels+1)))/ & |
---|
| 110 | (maxval(cospstateIN%phalf(:,Nlevels+1))-minval(cospstateIN%phalf(:,Nlevels+1)))*100000) + 1 |
---|
| 111 | call init_rng(rngs, seed) |
---|
| 112 | |
---|
| 113 | ! Call scops |
---|
| 114 | call scops(NPoints,Nlevels,Ncolumns,rngs,tca,cca,overlap,cospIN%frac_out,0) |
---|
| 115 | deallocate(seed,rngs) |
---|
| 116 | |
---|
| 117 | ! Sum up precipitation rates |
---|
| 118 | allocate(ls_p_rate(nPoints,nLevels),cv_p_rate(nPoints,Nlevels)) |
---|
| 119 | if(use_precipitation_fluxes) then |
---|
| 120 | ls_p_rate(:,1:nLevels) = fl_lsrainIN + fl_lssnowIN + fl_lsgrplIN |
---|
| 121 | cv_p_rate(:,1:nLevels) = fl_ccrainIN + fl_ccsnowIN |
---|
| 122 | else |
---|
| 123 | ls_p_rate(:,1:nLevels) = 0 ! mixing_ratio(rain) + mixing_ratio(snow) + mixing_ratio (groupel) |
---|
| 124 | cv_p_rate(:,1:nLevels) = 0 ! mixing_ratio(rain) + mixing_ratio(snow) |
---|
| 125 | endif |
---|
| 126 | |
---|
| 127 | ! Call PREC_SCOPS |
---|
| 128 | allocate(frac_prec(nPoints,nColumns,nLevels)) |
---|
| 129 | call prec_scops(nPoints,nLevels,nColumns,ls_p_rate,cv_p_rate,cospIN%frac_out,frac_prec) |
---|
| 130 | deallocate(ls_p_rate,cv_p_rate) |
---|
| 131 | |
---|
| 132 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 133 | ! Compute fraction in each gridbox for precipitation and cloud type. |
---|
| 134 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 135 | ! Allocate |
---|
| 136 | allocate(frac_ls(nPoints,nLevels),prec_ls(nPoints,nLevels), & |
---|
| 137 | frac_cv(nPoints,nLevels),prec_cv(nPoints,nLevels)) |
---|
| 138 | |
---|
| 139 | ! Initialize |
---|
| 140 | frac_ls(1:nPoints,1:nLevels) = 0._wp |
---|
| 141 | prec_ls(1:nPoints,1:nLevels) = 0._wp |
---|
| 142 | frac_cv(1:nPoints,1:nLevels) = 0._wp |
---|
| 143 | prec_cv(1:nPoints,1:nLevels) = 0._wp |
---|
| 144 | do j=1,nPoints |
---|
| 145 | do k=1,nLevels |
---|
| 146 | do i=1,nColumns |
---|
| 147 | if (cospIN%frac_out(j,i,k) .eq. 1) frac_ls(j,k) = frac_ls(j,k)+1._wp |
---|
| 148 | if (cospIN%frac_out(j,i,k) .eq. 2) frac_cv(j,k) = frac_cv(j,k)+1._wp |
---|
| 149 | if (frac_prec(j,i,k) .eq. 1) prec_ls(j,k) = prec_ls(j,k)+1._wp |
---|
| 150 | if (frac_prec(j,i,k) .eq. 2) prec_cv(j,k) = prec_cv(j,k)+1._wp |
---|
| 151 | if (frac_prec(j,i,k) .eq. 3) prec_cv(j,k) = prec_cv(j,k)+1._wp |
---|
| 152 | if (frac_prec(j,i,k) .eq. 3) prec_ls(j,k) = prec_ls(j,k)+1._wp |
---|
| 153 | enddo |
---|
| 154 | frac_ls(j,k)=frac_ls(j,k)/nColumns |
---|
| 155 | frac_cv(j,k)=frac_cv(j,k)/nColumns |
---|
| 156 | prec_ls(j,k)=prec_ls(j,k)/nColumns |
---|
| 157 | prec_cv(j,k)=prec_cv(j,k)/nColumns |
---|
| 158 | enddo |
---|
| 159 | enddo |
---|
| 160 | |
---|
| 161 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 162 | ! Assign gridmean mixing-ratios (mr_XXXXX), effective radius (ReffIN) and number |
---|
| 163 | ! concentration (not defined) to appropriate sub-column. Here we are using scops. |
---|
| 164 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 165 | allocate(mr_hydro(nPoints,nColumns,nLevels,nHydro), & |
---|
| 166 | Reff(nPoints,nColumns,nLevels,nHydro), & |
---|
| 167 | Np(nPoints,nColumns,nLevels,nHydro)) |
---|
| 168 | |
---|
| 169 | ! Initialize |
---|
| 170 | mr_hydro(:,:,:,:) = 0._wp |
---|
| 171 | Reff(:,:,:,:) = 0._wp |
---|
| 172 | Np(:,:,:,:) = 0._wp |
---|
| 173 | do k=1,nColumns |
---|
| 174 | ! Subcolumn cloud fraction |
---|
| 175 | column_frac_out = cospIN%frac_out(:,k,:) |
---|
| 176 | |
---|
| 177 | ! LS clouds |
---|
| 178 | where (column_frac_out == I_LSC) |
---|
| 179 | mr_hydro(:,k,:,I_LSCLIQ) = mr_lsliq |
---|
| 180 | mr_hydro(:,k,:,I_LSCICE) = mr_lsice |
---|
| 181 | Reff(:,k,:,I_LSCLIQ) = ReffIN(:,:,I_LSCLIQ) |
---|
| 182 | Reff(:,k,:,I_LSCICE) = ReffIN(:,:,I_LSCICE) |
---|
| 183 | ! CONV clouds |
---|
| 184 | elsewhere (column_frac_out == I_CVC) |
---|
| 185 | mr_hydro(:,k,:,I_CVCLIQ) = mr_ccliq |
---|
| 186 | mr_hydro(:,k,:,I_CVCICE) = mr_ccice |
---|
| 187 | Reff(:,k,:,I_CVCLIQ) = ReffIN(:,:,I_CVCLIQ) |
---|
| 188 | Reff(:,k,:,I_CVCICE) = ReffIN(:,:,I_CVCICE) |
---|
| 189 | end where |
---|
| 190 | |
---|
| 191 | ! Subcolumn precipitation |
---|
| 192 | column_prec_out = frac_prec(:,k,:) |
---|
| 193 | |
---|
| 194 | ! LS Precipitation |
---|
| 195 | where ((column_prec_out == 1) .or. (column_prec_out == 3) ) |
---|
| 196 | Reff(:,k,:,I_LSRAIN) = ReffIN(:,:,I_LSRAIN) |
---|
| 197 | Reff(:,k,:,I_LSSNOW) = ReffIN(:,:,I_LSSNOW) |
---|
| 198 | Reff(:,k,:,I_LSGRPL) = ReffIN(:,:,I_LSGRPL) |
---|
| 199 | ! CONV precipitation |
---|
| 200 | elsewhere ((column_prec_out == 2) .or. (column_prec_out == 3)) |
---|
| 201 | Reff(:,k,:,I_CVRAIN) = ReffIN(:,:,I_CVRAIN) |
---|
| 202 | Reff(:,k,:,I_CVSNOW) = ReffIN(:,:,I_CVSNOW) |
---|
| 203 | end where |
---|
| 204 | enddo |
---|
| 205 | |
---|
| 206 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 207 | ! Convert the subcolumn mixing ratio and precipitation fluxes from gridbox mean |
---|
| 208 | ! values to fraction-based values. |
---|
| 209 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 210 | ! Initialize |
---|
| 211 | fl_lsrain(:,:) = 0._wp |
---|
| 212 | fl_lssnow(:,:) = 0._wp |
---|
| 213 | fl_lsgrpl(:,:) = 0._wp |
---|
| 214 | fl_ccrain(:,:) = 0._wp |
---|
| 215 | fl_ccsnow(:,:) = 0._wp |
---|
| 216 | do k=1,nLevels |
---|
| 217 | do j=1,nPoints |
---|
| 218 | ! In-cloud mixing ratios. |
---|
| 219 | if (frac_ls(j,k) .ne. 0.) then |
---|
| 220 | mr_hydro(j,:,k,I_LSCLIQ) = mr_hydro(j,:,k,I_LSCLIQ)/frac_ls(j,k) |
---|
| 221 | mr_hydro(j,:,k,I_LSCICE) = mr_hydro(j,:,k,I_LSCICE)/frac_ls(j,k) |
---|
| 222 | endif |
---|
| 223 | if (frac_cv(j,k) .ne. 0.) then |
---|
| 224 | mr_hydro(j,:,k,I_CVCLIQ) = mr_hydro(j,:,k,I_CVCLIQ)/frac_cv(j,k) |
---|
| 225 | mr_hydro(j,:,k,I_CVCICE) = mr_hydro(j,:,k,I_CVCICE)/frac_cv(j,k) |
---|
| 226 | endif |
---|
| 227 | ! Precipitation |
---|
| 228 | if (use_precipitation_fluxes) then |
---|
| 229 | if (prec_ls(j,k) .ne. 0.) then |
---|
| 230 | fl_lsrain(j,k) = fl_lsrainIN(j,k)/prec_ls(j,k) |
---|
| 231 | fl_lssnow(j,k) = fl_lssnowIN(j,k)/prec_ls(j,k) |
---|
| 232 | fl_lsgrpl(j,k) = fl_lsgrplIN(j,k)/prec_ls(j,k) |
---|
| 233 | endif |
---|
| 234 | if (prec_cv(j,k) .ne. 0.) then |
---|
| 235 | fl_ccrain(j,k) = fl_ccrainIN(j,k)/prec_cv(j,k) |
---|
| 236 | fl_ccsnow(j,k) = fl_ccsnowIN(j,k)/prec_cv(j,k) |
---|
| 237 | endif |
---|
| 238 | else |
---|
| 239 | if (prec_ls(j,k) .ne. 0.) then |
---|
| 240 | mr_hydro(j,:,k,I_LSRAIN) = mr_hydro(j,:,k,I_LSRAIN)/prec_ls(j,k) |
---|
| 241 | mr_hydro(j,:,k,I_LSSNOW) = mr_hydro(j,:,k,I_LSSNOW)/prec_ls(j,k) |
---|
| 242 | mr_hydro(j,:,k,I_LSGRPL) = mr_hydro(j,:,k,I_LSGRPL)/prec_ls(j,k) |
---|
| 243 | endif |
---|
| 244 | if (prec_cv(j,k) .ne. 0.) then |
---|
| 245 | mr_hydro(j,:,k,I_CVRAIN) = mr_hydro(j,:,k,I_CVRAIN)/prec_cv(j,k) |
---|
| 246 | mr_hydro(j,:,k,I_CVSNOW) = mr_hydro(j,:,k,I_CVSNOW)/prec_cv(j,k) |
---|
| 247 | endif |
---|
| 248 | endif |
---|
| 249 | enddo |
---|
| 250 | enddo |
---|
| 251 | |
---|
| 252 | deallocate(frac_ls,prec_ls,frac_cv,prec_cv) |
---|
| 253 | |
---|
| 254 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 255 | ! Convert precipitation fluxes to mixing ratios |
---|
| 256 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 257 | if (use_precipitation_fluxes) then |
---|
| 258 | ! LS rain |
---|
| 259 | call cosp_precip_mxratio(nPoints, nLevels, nColumns, cospstateIN%pfull, & |
---|
| 260 | cospstateIN%at, frac_prec, 1._wp, n_ax(I_LSRAIN), n_bx(I_LSRAIN), & |
---|
| 261 | alpha_x(I_LSRAIN), c_x(I_LSRAIN), d_x(I_LSRAIN), g_x(I_LSRAIN), & |
---|
| 262 | a_x(I_LSRAIN), b_x(I_LSRAIN), gamma_1(I_LSRAIN), gamma_2(I_LSRAIN), & |
---|
| 263 | gamma_3(I_LSRAIN), gamma_4(I_LSRAIN), fl_lsrain, & |
---|
| 264 | mr_hydro(:,:,:,I_LSRAIN), Reff(:,:,:,I_LSRAIN)) |
---|
| 265 | ! LS snow |
---|
| 266 | call cosp_precip_mxratio(nPoints, nLevels, nColumns, cospstateIN%pfull, & |
---|
| 267 | cospstateIN%at, frac_prec, 1._wp, n_ax(I_LSSNOW), n_bx(I_LSSNOW), & |
---|
| 268 | alpha_x(I_LSSNOW), c_x(I_LSSNOW), d_x(I_LSSNOW), g_x(I_LSSNOW), & |
---|
| 269 | a_x(I_LSSNOW), b_x(I_LSSNOW), gamma_1(I_LSSNOW), gamma_2(I_LSSNOW), & |
---|
| 270 | gamma_3(I_LSSNOW), gamma_4(I_LSSNOW), fl_lssnow, & |
---|
| 271 | mr_hydro(:,:,:,I_LSSNOW), Reff(:,:,:,I_LSSNOW)) |
---|
| 272 | ! CV rain |
---|
| 273 | call cosp_precip_mxratio(nPoints, nLevels, nColumns, cospstateIN%pfull, & |
---|
| 274 | cospstateIN%at, frac_prec, 2._wp, n_ax(I_CVRAIN), n_bx(I_CVRAIN), & |
---|
| 275 | alpha_x(I_CVRAIN), c_x(I_CVRAIN), d_x(I_CVRAIN), g_x(I_CVRAIN), & |
---|
| 276 | a_x(I_CVRAIN), b_x(I_CVRAIN), gamma_1(I_CVRAIN), gamma_2(I_CVRAIN), & |
---|
| 277 | gamma_3(I_CVRAIN), gamma_4(I_CVRAIN), fl_ccrain, & |
---|
| 278 | mr_hydro(:,:,:,I_CVRAIN), Reff(:,:,:,I_CVRAIN)) |
---|
| 279 | ! CV snow |
---|
| 280 | call cosp_precip_mxratio(nPoints, nLevels, nColumns, cospstateIN%pfull, & |
---|
| 281 | cospstateIN%at, frac_prec, 2._wp, n_ax(I_CVSNOW), n_bx(I_CVSNOW), & |
---|
| 282 | alpha_x(I_CVSNOW), c_x(I_CVSNOW), d_x(I_CVSNOW), g_x(I_CVSNOW), & |
---|
| 283 | a_x(I_CVSNOW), b_x(I_CVSNOW), gamma_1(I_CVSNOW), gamma_2(I_CVSNOW), & |
---|
| 284 | gamma_3(I_CVSNOW), gamma_4(I_CVSNOW), fl_ccsnow, & |
---|
| 285 | mr_hydro(:,:,:,I_CVSNOW), Reff(:,:,:,I_CVSNOW)) |
---|
| 286 | ! LS groupel. |
---|
| 287 | call cosp_precip_mxratio(nPoints, nLevels, nColumns, cospstateIN%pfull, & |
---|
| 288 | cospstateIN%at, frac_prec, 1._wp, n_ax(I_LSGRPL), n_bx(I_LSGRPL), & |
---|
| 289 | alpha_x(I_LSGRPL), c_x(I_LSGRPL), d_x(I_LSGRPL), g_x(I_LSGRPL), & |
---|
| 290 | a_x(I_LSGRPL), b_x(I_LSGRPL), gamma_1(I_LSGRPL), gamma_2(I_LSGRPL), & |
---|
| 291 | gamma_3(I_LSGRPL), gamma_4(I_LSGRPL), fl_lsgrpl, & |
---|
| 292 | mr_hydro(:,:,:,I_LSGRPL), Reff(:,:,:,I_LSGRPL)) |
---|
| 293 | deallocate(frac_prec) |
---|
| 294 | endif |
---|
| 295 | |
---|
| 296 | else |
---|
| 297 | cospIN%frac_out(:,:,:) = 1 |
---|
| 298 | allocate(mr_hydro(nPoints,1,nLevels,nHydro),Reff(nPoints,1,nLevels,nHydro), & |
---|
| 299 | Np(nPoints,1,nLevels,nHydro)) |
---|
| 300 | mr_hydro(:,1,:,I_LSCLIQ) = mr_lsliq |
---|
| 301 | mr_hydro(:,1,:,I_LSCICE) = mr_lsice |
---|
| 302 | mr_hydro(:,1,:,I_CVCLIQ) = mr_ccliq |
---|
| 303 | mr_hydro(:,1,:,I_CVCICE) = mr_ccice |
---|
| 304 | Reff(:,1,:,:) = ReffIN |
---|
| 305 | endif |
---|
| 306 | |
---|
| 307 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 308 | ! 11 micron emissivity |
---|
| 309 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 310 | if (cfg%Lisccp) then |
---|
| 311 | call cosp_simulator_optics(nPoints,nColumns,nLevels,cospIN%frac_out,dem_c,dem_s, & |
---|
| 312 | cospIN%emiss_11) |
---|
| 313 | endif |
---|
| 314 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 315 | ! 0.67 micron optical depth |
---|
| 316 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 317 | if (cfg%Lisccp .or. cfg%Lmisr .or. cfg%Lmodis) then |
---|
| 318 | call cosp_simulator_optics(nPoints,nColumns,nLevels,cospIN%frac_out,dtau_c,dtau_s, & |
---|
| 319 | cospIN%tau_067) |
---|
| 320 | endif |
---|
| 321 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 322 | ! LIDAR Polarized optics |
---|
| 323 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 324 | if (cfg%Lcalipso) then |
---|
| 325 | call lidar_optics(nPoints, nColumns, nLevels, 4, lidar_ice_type, 532, .false., & |
---|
| 326 | mr_hydro(:,:,:,I_LSCLIQ), mr_hydro(:,:,:,I_LSCICE), mr_hydro(:,:,:,I_CVCLIQ), & |
---|
| 327 | mr_hydro(:,:,:,I_CVCICE), ReffIN(:,:,I_LSCLIQ), ReffIN(:,:,I_LSCICE), & |
---|
| 328 | ReffIN(:,:,I_CVCLIQ), ReffIN(:,:,I_CVCICE), cospstateIN%pfull, & |
---|
| 329 | cospstateIN%phalf, cospstateIN%at, cospIN%beta_mol_calipso, & |
---|
| 330 | cospIN%betatot_calipso, cospIN%tau_mol_calipso, cospIN%tautot_calipso, & |
---|
| 331 | cospIN%tautot_S_liq, cospIN%tautot_S_ice, cospIN%betatot_ice_calipso, & |
---|
| 332 | cospIN%betatot_liq_calipso, cospIN%tautot_ice_calipso, cospIN%tautot_liq_calipso) |
---|
| 333 | endif |
---|
| 334 | |
---|
| 335 | if (cfg%LgrLidar532) then |
---|
| 336 | call lidar_optics(nPoints, nColumns, nLevels, 4, lidar_ice_type, 532, .true., & |
---|
| 337 | mr_hydro(:,:,:,I_LSCLIQ), mr_hydro(:,:,:,I_LSCICE), mr_hydro(:,:,:,I_CVCLIQ), & |
---|
| 338 | mr_hydro(:,:,:,I_CVCICE), ReffIN(:,:,I_LSCLIQ), ReffIN(:,:,I_LSCICE), & |
---|
| 339 | ReffIN(:,:,I_CVCLIQ), ReffIN(:,:,I_CVCICE), cospstateIN%pfull, & |
---|
| 340 | cospstateIN%phalf, cospstateIN%at, cospIN%beta_mol_grLidar532, & |
---|
| 341 | cospIN%betatot_grLidar532, cospIN%tau_mol_grLidar532, cospIN%tautot_grLidar532) |
---|
| 342 | endif |
---|
| 343 | |
---|
| 344 | if (cfg%Latlid) then |
---|
| 345 | call lidar_optics(nPoints, nColumns, nLevels, 4, lidar_ice_type, 355, .false., & |
---|
| 346 | mr_hydro(:,:,:,I_LSCLIQ), mr_hydro(:,:,:,I_LSCICE), mr_hydro(:,:,:,I_CVCLIQ), & |
---|
| 347 | mr_hydro(:,:,:,I_CVCICE), ReffIN(:,:,I_LSCLIQ), ReffIN(:,:,I_LSCICE), & |
---|
| 348 | ReffIN(:,:,I_CVCLIQ), ReffIN(:,:,I_CVCICE), cospstateIN%pfull, & |
---|
| 349 | cospstateIN%phalf, cospstateIN%at, cospIN%beta_mol_atlid, cospIN%betatot_atlid,& |
---|
| 350 | cospIN%tau_mol_atlid, cospIN%tautot_atlid) |
---|
| 351 | endif |
---|
| 352 | |
---|
| 353 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 354 | ! CLOUDSAT RADAR OPTICS |
---|
| 355 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 356 | if (cfg%Lcloudsat) then |
---|
| 357 | |
---|
| 358 | ! Compute gaseous absorption (assume identical for each subcolun) |
---|
| 359 | allocate(g_vol(nPoints,nLevels)) |
---|
| 360 | g_vol(:,:)=0._wp |
---|
| 361 | do i=1,nPoints |
---|
| 362 | do j=1,nLevels |
---|
| 363 | if (cospIN%rcfg_cloudsat%use_gas_abs == 1 .or. (cospIN%rcfg_cloudsat%use_gas_abs == 2 .and. j .eq. 1)) then |
---|
| 364 | g_vol(i,j) = gases(cospstateIN%pfull(i,j), cospstateIN%at(i,j),cospstateIN%qv(i,j),cospIN%rcfg_cloudsat%freq) |
---|
| 365 | endif |
---|
| 366 | cospIN%g_vol_cloudsat(i,:,j)=g_vol(i,j) |
---|
| 367 | end do |
---|
| 368 | end do |
---|
| 369 | |
---|
| 370 | ! Loop over all subcolumns |
---|
| 371 | allocate(fracPrecipIce(nPoints,nColumns,nLevels)) |
---|
| 372 | fracPrecipIce(:,:,:) = 0._wp |
---|
| 373 | do k=1,nColumns |
---|
| 374 | call quickbeam_optics(sd, cospIN%rcfg_cloudsat, nPoints, nLevels, R_UNDEF, & |
---|
| 375 | mr_hydro(:,k,:,1:nHydro)*1000._wp, Reff(:,k,:,1:nHydro)*1.e6_wp,& |
---|
| 376 | Np(:,k,:,1:nHydro), cospstateIN%pfull, cospstateIN%at, & |
---|
| 377 | cospstateIN%qv, cospIN%z_vol_cloudsat(1:nPoints,k,:), & |
---|
| 378 | cospIN%kr_vol_cloudsat(1:nPoints,k,:)) |
---|
| 379 | |
---|
| 380 | ! At each model level, what fraction of the precipitation is frozen? |
---|
| 381 | where(mr_hydro(:,k,:,I_LSRAIN) .gt. 0 .or. mr_hydro(:,k,:,I_LSSNOW) .gt. 0 .or. & |
---|
| 382 | mr_hydro(:,k,:,I_CVRAIN) .gt. 0 .or. mr_hydro(:,k,:,I_CVSNOW) .gt. 0 .or. & |
---|
| 383 | mr_hydro(:,k,:,I_LSGRPL) .gt. 0) |
---|
| 384 | fracPrecipIce(:,k,:) = (mr_hydro(:,k,:,I_LSSNOW) + mr_hydro(:,k,:,I_CVSNOW) + & |
---|
| 385 | mr_hydro(:,k,:,I_LSGRPL)) / & |
---|
| 386 | (mr_hydro(:,k,:,I_LSSNOW) + mr_hydro(:,k,:,I_CVSNOW) + mr_hydro(:,k,:,I_LSGRPL) + & |
---|
| 387 | mr_hydro(:,k,:,I_LSRAIN) + mr_hydro(:,k,:,I_CVRAIN)) |
---|
| 388 | elsewhere |
---|
| 389 | fracPrecipIce(:,k,:) = 0._wp |
---|
| 390 | endwhere |
---|
| 391 | enddo |
---|
| 392 | |
---|
| 393 | ! Regrid frozen fraction to Cloudsat/Calipso statistical grid |
---|
| 394 | allocate(fracPrecipIce_statGrid(nPoints,nColumns,Nlvgrid)) |
---|
| 395 | fracPrecipIce_statGrid(:,:,:) = 0._wp |
---|
| 396 | call cosp_change_vertical_grid(Npoints, Ncolumns, Nlevels, cospstateIN%hgt_matrix(:,Nlevels:1:-1), & |
---|
| 397 | cospstateIN%hgt_matrix_half(:,Nlevels:1:-1), fracPrecipIce(:,:,Nlevels:1:-1), Nlvgrid, & |
---|
| 398 | ! vgrid_zl(Nlvgrid:1:-1), vgrid_zu(Nlvgrid:1:-1), fracPrecipIce_statGrid) !!! ORIGINAL |
---|
| 399 | vgrid_zl(Nlvgrid:1:-1), vgrid_zu(Nlvgrid:1:-1), fracPrecipIce_statGrid(:,:,Nlvgrid:1:-1)) !DEBUG fracPrecipIce_statGrid |
---|
| 400 | |
---|
| 401 | ! Find proper layer to compute precip flags in Cloudsat/Calipso statistical grid !PREC_BUG |
---|
| 402 | allocate(cloudsat_preclvl_index(nPoints)) !PREC_BUG |
---|
| 403 | cloudsat_preclvl_index(:) = 0._wp !PREC_BUG |
---|
| 404 | ! Computing altitude index for precip flags calculation (2nd layer above surfelev) !PREC_BUG |
---|
| 405 | cloudsat_preclvl_index(:) = 39 - floor( cospstateIN%surfelev(:)/480. ) !PREC_BUG |
---|
| 406 | |
---|
| 407 | ! For near-surface diagnostics, we only need the frozen fraction at one layer. |
---|
| 408 | do i=1,nPoints !PREC_BUG |
---|
| 409 | cospIN%fracPrecipIce(i,:) = fracPrecipIce_statGrid(i,:,cloudsat_preclvl_index(i)) !PREC_BUG |
---|
| 410 | enddo !PREC_BUG |
---|
| 411 | ! cospIN%fracPrecipIce(:,:) = fracPrecipIce_statGrid(:,:,cloudsat_preclvl) !!! ORIGINAL |
---|
| 412 | ! cospIN%fracPrecipIce(:,:) = fracPrecipIce_statGrid(:,:,2) !!! TEST ARTEM |
---|
| 413 | |
---|
| 414 | endif |
---|
| 415 | |
---|
| 416 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 417 | ! MODIS optics |
---|
| 418 | !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 419 | if (cfg%Lmodis) then |
---|
| 420 | allocate(MODIS_cloudWater(nPoints,nColumns,nLevels), & |
---|
| 421 | MODIS_cloudIce(nPoints,nColumns,nLevels), & |
---|
| 422 | MODIS_waterSize(nPoints,nColumns,nLevels), & |
---|
| 423 | MODIS_iceSize(nPoints,nColumns,nLevels), & |
---|
| 424 | MODIS_opticalThicknessLiq(nPoints,nColumns,nLevels), & |
---|
| 425 | MODIS_opticalThicknessIce(nPoints,nColumns,nLevels)) |
---|
| 426 | ! Cloud water |
---|
| 427 | call cosp_simulator_optics(nPoints,nColumns,nLevels,cospIN%frac_out, & |
---|
| 428 | mr_hydro(:,:,:,I_CVCLIQ),mr_hydro(:,:,:,I_LSCLIQ),MODIS_cloudWater) |
---|
| 429 | ! Cloud ice |
---|
| 430 | call cosp_simulator_optics(nPoints,nColumns,nLevels,cospIN%frac_out, & |
---|
| 431 | mr_hydro(:,:,:,I_CVCICE),mr_hydro(:,:,:,I_LSCICE),MODIS_cloudIce) |
---|
| 432 | ! Water droplet size |
---|
| 433 | call cosp_simulator_optics(nPoints,nColumns,nLevels,cospIN%frac_out, & |
---|
| 434 | Reff(:,:,:,I_CVCLIQ),Reff(:,:,:,I_LSCLIQ),MODIS_waterSize) |
---|
| 435 | ! Ice crystal size |
---|
| 436 | call cosp_simulator_optics(nPoints,nColumns,nLevels,cospIN%frac_out, & |
---|
| 437 | Reff(:,:,:,I_CVCICE),Reff(:,:,:,I_LSCICE),MODIS_iceSize) |
---|
| 438 | |
---|
| 439 | ! Partition optical thickness into liquid and ice parts |
---|
| 440 | call modis_optics_partition(nPoints, nLevels, nColumns, MODIS_cloudWater, & |
---|
| 441 | MODIS_cloudIce, MODIS_waterSize, MODIS_iceSize, cospIN%tau_067, & |
---|
| 442 | MODIS_opticalThicknessLiq, MODIS_opticalThicknessIce) |
---|
| 443 | |
---|
| 444 | ! Compute assymetry parameter and single scattering albedo |
---|
| 445 | call modis_optics(nPoints, nLevels, nColumns, MODIS_opticalThicknessLiq, & |
---|
| 446 | MODIS_waterSize*1.0e6_wp, MODIS_opticalThicknessIce, & |
---|
| 447 | MODIS_iceSize*1.0e6_wp, cospIN%fracLiq, cospIN%asym, cospIN%ss_alb) |
---|
| 448 | |
---|
| 449 | ! Deallocate memory |
---|
| 450 | deallocate(MODIS_cloudWater,MODIS_cloudIce,MODIS_WaterSize,MODIS_iceSize, & |
---|
| 451 | MODIS_opticalThicknessLiq,MODIS_opticalThicknessIce,mr_hydro, & |
---|
| 452 | Np,Reff) |
---|
| 453 | endif |
---|
| 454 | end subroutine subsample_and_optics |
---|
| 455 | |
---|
| 456 | END MODULE LMDZ_COSP_SUBSAMPLE_AND_OPTICS_MOD |
---|