[253] | 1 | Subroutine aeropacity(ngrid,nlayer,nq,pplay,pplev,pq, & |
---|
| 2 | aerosol,reffrad,QREFvis3d,QREFir3d,tau_col,cloudfrac,totcloudfrac,clearsky) |
---|
[135] | 3 | |
---|
[726] | 4 | use radinc_h, only : L_TAUMAX,naerkind |
---|
| 5 | use aerosol_mod |
---|
[787] | 6 | USE comgeomfi_h |
---|
[858] | 7 | USE tracer_h, only: noms,rho_co2,rho_ice |
---|
[135] | 8 | |
---|
| 9 | implicit none |
---|
| 10 | |
---|
| 11 | !================================================================== |
---|
| 12 | ! |
---|
| 13 | ! Purpose |
---|
| 14 | ! ------- |
---|
| 15 | ! Compute aerosol optical depth in each gridbox. |
---|
| 16 | ! |
---|
| 17 | ! Authors |
---|
| 18 | ! ------- |
---|
| 19 | ! F. Forget |
---|
| 20 | ! F. Montmessin (water ice scheme) |
---|
| 21 | ! update J.-B. Madeleine (2008) |
---|
| 22 | ! dust removal, simplification by Robin Wordsworth (2009) |
---|
| 23 | ! |
---|
| 24 | ! Input |
---|
| 25 | ! ----- |
---|
| 26 | ! ngrid Number of horizontal gridpoints |
---|
| 27 | ! nlayer Number of layers |
---|
| 28 | ! nq Number of tracers |
---|
| 29 | ! pplev Pressure (Pa) at each layer boundary |
---|
| 30 | ! pq Aerosol mixing ratio |
---|
| 31 | ! reffrad(ngrid,nlayer,naerkind) Aerosol effective radius |
---|
[787] | 32 | ! QREFvis3d(ngrid,nlayermx,naerkind) \ 3d extinction coefficients |
---|
| 33 | ! QREFir3d(ngrid,nlayermx,naerkind) / at reference wavelengths |
---|
[135] | 34 | ! |
---|
| 35 | ! Output |
---|
| 36 | ! ------ |
---|
| 37 | ! aerosol Aerosol optical depth in layer l, grid point ig |
---|
| 38 | ! tau_col Total column optical depth at grid point ig |
---|
| 39 | ! |
---|
| 40 | !======================================================================= |
---|
| 41 | |
---|
| 42 | #include "dimensions.h" |
---|
| 43 | #include "dimphys.h" |
---|
| 44 | #include "callkeys.h" |
---|
| 45 | #include "comcstfi.h" |
---|
[726] | 46 | #include "comvert.h" |
---|
[135] | 47 | |
---|
[858] | 48 | INTEGER,INTENT(IN) :: ngrid ! number of atmospheric columns |
---|
| 49 | INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers |
---|
| 50 | INTEGER,INTENT(IN) :: nq ! number of tracers |
---|
| 51 | REAL,INTENT(IN) :: pplay(ngrid,nlayer) ! mid-layer pressure (Pa) |
---|
| 52 | REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa) |
---|
| 53 | REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) ! tracers (.../kg_of_air) |
---|
| 54 | REAL,INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) ! aerosol optical depth |
---|
| 55 | REAL,INTENT(IN) :: reffrad(ngrid,nlayer,naerkind) ! aerosol effective radius |
---|
| 56 | REAL,INTENT(IN) :: QREFvis3d(ngrid,nlayermx,naerkind) ! extinction coefficient in the visible |
---|
| 57 | REAL,INTENT(IN) :: QREFir3d(ngrid,nlayermx,naerkind) |
---|
| 58 | REAL,INTENT(OUT):: tau_col(ngrid) !column integrated visible optical depth |
---|
| 59 | ! BENJAMIN MODIFS |
---|
| 60 | real,intent(in) :: cloudfrac(ngrid,nlayermx) ! cloud fraction |
---|
| 61 | real,intent(out) :: totcloudfrac(ngrid) ! total cloud fraction |
---|
| 62 | logical,intent(in) :: clearsky |
---|
[135] | 63 | |
---|
[253] | 64 | real aerosol0 |
---|
| 65 | |
---|
[135] | 66 | INTEGER l,ig,iq,iaer |
---|
| 67 | |
---|
[858] | 68 | LOGICAL,SAVE :: firstcall=.true. |
---|
[135] | 69 | REAL CBRT |
---|
| 70 | EXTERNAL CBRT |
---|
| 71 | |
---|
| 72 | INTEGER,SAVE :: i_co2ice=0 ! co2 ice |
---|
| 73 | INTEGER,SAVE :: i_h2oice=0 ! water ice |
---|
| 74 | CHARACTER(LEN=20) :: tracername ! to temporarily store text |
---|
| 75 | |
---|
[253] | 76 | ! for fixed dust profiles |
---|
| 77 | real topdust, expfactor, zp |
---|
[787] | 78 | REAL taudusttmp(ngrid) ! Temporary dust opacity used before scaling |
---|
| 79 | REAL tauh2so4tmp(ngrid) ! Temporary h2so4 opacity used before scaling |
---|
[858] | 80 | |
---|
[728] | 81 | real CLFtot |
---|
[253] | 82 | |
---|
| 83 | ! identify tracers |
---|
[135] | 84 | IF (firstcall) THEN |
---|
| 85 | |
---|
[726] | 86 | write(*,*) "Tracers found in aeropacity:" |
---|
[787] | 87 | do iq=1,nq |
---|
[135] | 88 | tracername=noms(iq) |
---|
| 89 | if (tracername.eq."co2_ice") then |
---|
| 90 | i_co2ice=iq |
---|
[726] | 91 | write(*,*) "i_co2ice=",i_co2ice |
---|
| 92 | |
---|
[135] | 93 | endif |
---|
| 94 | if (tracername.eq."h2o_ice") then |
---|
| 95 | i_h2oice=iq |
---|
[726] | 96 | write(*,*) "i_h2oice=",i_h2oice |
---|
[135] | 97 | endif |
---|
| 98 | enddo |
---|
| 99 | |
---|
[741] | 100 | if (noaero) then |
---|
| 101 | print*, "No active aerosols found in aeropacity" |
---|
| 102 | else |
---|
[726] | 103 | print*, "If you would like to use aerosols, make sure any old" |
---|
| 104 | print*, "start files are updated in newstart using the option" |
---|
| 105 | print*, "q=0" |
---|
[741] | 106 | write(*,*) "Active aerosols found in aeropacity:" |
---|
[726] | 107 | endif |
---|
[253] | 108 | |
---|
[741] | 109 | if ((iaero_co2.ne.0).and.(.not.noaero)) then |
---|
[726] | 110 | print*, 'iaero_co2= ',iaero_co2 |
---|
| 111 | endif |
---|
| 112 | if (iaero_h2o.ne.0) then |
---|
| 113 | print*,'iaero_h2o= ',iaero_h2o |
---|
| 114 | endif |
---|
| 115 | if (iaero_dust.ne.0) then |
---|
| 116 | print*,'iaero_dust= ',iaero_dust |
---|
| 117 | endif |
---|
| 118 | if (iaero_h2so4.ne.0) then |
---|
| 119 | print*,'iaero_h2so4= ',iaero_h2so4 |
---|
| 120 | endif |
---|
[1026] | 121 | if (iaero_back2lay.ne.0) then |
---|
| 122 | print*,'iaero_back2lay= ',iaero_back2lay |
---|
| 123 | endif |
---|
[726] | 124 | |
---|
[135] | 125 | firstcall=.false. |
---|
| 126 | ENDIF ! of IF (firstcall) |
---|
| 127 | |
---|
| 128 | |
---|
| 129 | ! --------------------------------------------------------- |
---|
| 130 | !================================================================== |
---|
[726] | 131 | ! CO2 ice aerosols |
---|
[135] | 132 | !================================================================== |
---|
| 133 | |
---|
[728] | 134 | if (iaero_co2.ne.0) then |
---|
[726] | 135 | iaer=iaero_co2 |
---|
[135] | 136 | ! 1. Initialization |
---|
[787] | 137 | aerosol(1:ngrid,1:nlayermx,iaer)=0.0 |
---|
[135] | 138 | ! 2. Opacity calculation |
---|
[804] | 139 | if (noaero) then ! aerosol set to zero |
---|
| 140 | aerosol(1:ngrid,1:nlayermx,iaer)=0.0 |
---|
[741] | 141 | elseif (aerofixco2.or.(i_co2ice.eq.0)) then ! CO2 ice cloud prescribed |
---|
[787] | 142 | aerosol(1:ngrid,1:nlayermx,iaer)=1.e-9 |
---|
| 143 | !aerosol(1:ngrid,12,iaer)=4.0 ! single cloud layer option |
---|
[253] | 144 | else |
---|
| 145 | DO ig=1, ngrid |
---|
| 146 | DO l=1,nlayer-1 ! to stop the rad tran bug |
---|
[135] | 147 | |
---|
[253] | 148 | aerosol0 = & |
---|
| 149 | ( 0.75 * QREFvis3d(ig,l,iaer) / & |
---|
| 150 | ( rho_co2 * reffrad(ig,l,iaer) ) ) * & |
---|
| 151 | ( pq(ig,l,i_co2ice) + 1.E-9 ) * & |
---|
| 152 | ( pplev(ig,l) - pplev(ig,l+1) ) / g |
---|
| 153 | aerosol0 = max(aerosol0,1.e-9) |
---|
| 154 | aerosol0 = min(aerosol0,L_TAUMAX) |
---|
| 155 | aerosol(ig,l,iaer) = aerosol0 |
---|
| 156 | ! aerosol(ig,l,iaer) = 0.0 |
---|
[726] | 157 | ! print*, aerosol(ig,l,iaer) |
---|
[253] | 158 | ! using cloud fraction |
---|
| 159 | ! aerosol(ig,l,iaer) = -log(1 - CLF + CLF*exp(-aerosol0/CLF)) |
---|
| 160 | ! aerosol(ig,l,iaer) = min(aerosol(ig,l,iaer),L_TAUMAX) |
---|
| 161 | |
---|
| 162 | |
---|
| 163 | ENDDO |
---|
| 164 | ENDDO |
---|
[726] | 165 | end if ! if fixed or varying |
---|
[728] | 166 | end if ! if CO2 aerosols |
---|
[135] | 167 | !================================================================== |
---|
[726] | 168 | ! Water ice / liquid |
---|
[135] | 169 | !================================================================== |
---|
| 170 | |
---|
[728] | 171 | if (iaero_h2o.ne.0) then |
---|
[726] | 172 | iaer=iaero_h2o |
---|
[135] | 173 | ! 1. Initialization |
---|
[787] | 174 | aerosol(1:ngrid,1:nlayermx,iaer)=0.0 |
---|
[135] | 175 | ! 2. Opacity calculation |
---|
[726] | 176 | if (aerofixh2o.or.(i_h2oice.eq.0).or.clearsky) then |
---|
[787] | 177 | aerosol(1:ngrid,1:nlayermx,iaer) =1.e-9 |
---|
[305] | 178 | |
---|
| 179 | ! put cloud at cloudlvl |
---|
| 180 | if(kastprof.and.(cloudlvl.ne.0.0))then |
---|
| 181 | ig=1 |
---|
| 182 | do l=1,nlayer |
---|
| 183 | if(int(cloudlvl).eq.l)then |
---|
| 184 | !if(cloudlvl.gt.(pplay(ig,l)/pplev(ig,1)))then |
---|
| 185 | print*,'Inserting cloud at level ',l |
---|
| 186 | !aerosol(ig,l,iaer)=10.0 |
---|
| 187 | |
---|
| 188 | rho_ice=920.0 |
---|
| 189 | |
---|
| 190 | ! the Kasting approximation |
---|
| 191 | aerosol(ig,l,iaer) = & |
---|
| 192 | ( 0.75 * QREFvis3d(ig,l,iaer) / & |
---|
| 193 | ( rho_ice * reffrad(ig,l,iaer) ) ) * & |
---|
| 194 | !( pq(ig,l,i_h2oice) + 1.E-9 ) * & |
---|
| 195 | ( 4.0e-4 + 1.E-9 ) * & |
---|
| 196 | ( pplev(ig,l) - pplev(ig,l+1) ) / g |
---|
| 197 | |
---|
| 198 | |
---|
| 199 | open(115,file='clouds.out',form='formatted') |
---|
| 200 | write(115,*) l,aerosol(ig,l,iaer) |
---|
| 201 | close(115) |
---|
| 202 | |
---|
| 203 | return |
---|
| 204 | endif |
---|
| 205 | end do |
---|
| 206 | |
---|
| 207 | call abort |
---|
| 208 | endif |
---|
| 209 | |
---|
[253] | 210 | else |
---|
[728] | 211 | |
---|
[253] | 212 | do ig=1, ngrid |
---|
| 213 | do l=1,nlayer-1 ! to stop the rad tran bug |
---|
[135] | 214 | |
---|
[728] | 215 | aerosol(ig,l,iaer) = & !modification by BC |
---|
[253] | 216 | ( 0.75 * QREFvis3d(ig,l,iaer) / & |
---|
| 217 | ( rho_ice * reffrad(ig,l,iaer) ) ) * & |
---|
[716] | 218 | ! pq(ig,l,i_h2oice) * & !JL I dropped the +1e-9 here to have the same |
---|
[728] | 219 | !( pplev(ig,l) - pplev(ig,l+1) ) / g ! opacity in the clearsky=true and the |
---|
| 220 | ! clear=false/pq=0 case |
---|
[716] | 221 | ( pq(ig,l,i_h2oice) + 1.E-9 ) * & ! Doing this makes the code unstable, so I have restored it (RW) |
---|
[728] | 222 | ( pplev(ig,l) - pplev(ig,l+1) ) / g |
---|
[135] | 223 | |
---|
[253] | 224 | enddo |
---|
| 225 | enddo |
---|
[135] | 226 | |
---|
[728] | 227 | if(CLFvarying)then |
---|
[787] | 228 | call totalcloudfrac(ngrid,nq,cloudfrac,totcloudfrac,pplev,pq,aerosol(1:ngrid,1:nlayermx,iaer)) |
---|
[728] | 229 | do ig=1, ngrid |
---|
| 230 | do l=1,nlayer-1 ! to stop the rad tran bug |
---|
| 231 | CLFtot = max(totcloudfrac(ig),0.01) |
---|
| 232 | aerosol(ig,l,iaer)=aerosol(ig,l,iaer)/CLFtot |
---|
| 233 | aerosol(ig,l,iaer) = max(aerosol(ig,l,iaer),1.e-9) |
---|
| 234 | enddo |
---|
| 235 | enddo |
---|
| 236 | else |
---|
| 237 | do ig=1, ngrid |
---|
| 238 | do l=1,nlayer-1 ! to stop the rad tran bug |
---|
| 239 | CLFtot = CLFfixval |
---|
| 240 | aerosol(ig,l,iaer)=aerosol(ig,l,iaer)/CLFtot |
---|
| 241 | aerosol(ig,l,iaer) = max(aerosol(ig,l,iaer),1.e-9) |
---|
| 242 | enddo |
---|
| 243 | enddo |
---|
| 244 | end if!(CLFvarying) |
---|
| 245 | endif !(aerofixed.or.(i_h2oice.eq.0).or.clearsky) |
---|
| 246 | |
---|
| 247 | end if ! End if h2o aerosol |
---|
| 248 | |
---|
[726] | 249 | !================================================================== |
---|
| 250 | ! Dust |
---|
| 251 | !================================================================== |
---|
[728] | 252 | if (iaero_dust.ne.0) then |
---|
[726] | 253 | iaer=iaero_dust |
---|
| 254 | ! 1. Initialization |
---|
[787] | 255 | aerosol(1:ngrid,1:nlayermx,iaer)=0.0 |
---|
[726] | 256 | |
---|
[728] | 257 | topdust=30.0 ! km (used to be 10.0 km) LK |
---|
[726] | 258 | |
---|
| 259 | ! 2. Opacity calculation |
---|
| 260 | |
---|
| 261 | ! expfactor=0. |
---|
| 262 | DO l=1,nlayer-1 |
---|
| 263 | DO ig=1,ngrid |
---|
| 264 | ! Typical mixing ratio profile |
---|
| 265 | |
---|
| 266 | zp=(preff/pplay(ig,l))**(70./topdust) |
---|
| 267 | expfactor=max(exp(0.007*(1.-max(zp,1.))),1.e-3) |
---|
| 268 | |
---|
| 269 | ! Vertical scaling function |
---|
| 270 | aerosol(ig,l,iaer)= (pplev(ig,l)-pplev(ig,l+1)) & |
---|
| 271 | *expfactor |
---|
| 272 | |
---|
| 273 | |
---|
| 274 | ENDDO |
---|
| 275 | ENDDO |
---|
| 276 | |
---|
| 277 | ! Rescaling each layer to reproduce the choosen (or assimilated) |
---|
| 278 | ! dust extinction opacity at visible reference wavelength, which |
---|
| 279 | ! is scaled to the "preff" reference surface pressure available in comvert.h |
---|
| 280 | ! and stored in startfi.nc |
---|
| 281 | |
---|
| 282 | taudusttmp(1:ngrid)=0. |
---|
| 283 | DO l=1,nlayer |
---|
| 284 | DO ig=1,ngrid |
---|
| 285 | taudusttmp(ig) = taudusttmp(ig) & |
---|
| 286 | + aerosol(ig,l,iaer) |
---|
| 287 | ENDDO |
---|
| 288 | ENDDO |
---|
| 289 | DO l=1,nlayer-1 |
---|
| 290 | DO ig=1,ngrid |
---|
| 291 | aerosol(ig,l,iaer) = max(1E-20, & |
---|
| 292 | dusttau & |
---|
| 293 | * pplev(ig,1) / preff & |
---|
| 294 | * aerosol(ig,l,iaer) & |
---|
| 295 | / taudusttmp(ig)) |
---|
| 296 | |
---|
| 297 | ENDDO |
---|
| 298 | ENDDO |
---|
[728] | 299 | end if ! If dust aerosol |
---|
[726] | 300 | |
---|
[135] | 301 | !================================================================== |
---|
[726] | 302 | ! H2SO4 |
---|
[253] | 303 | !================================================================== |
---|
[726] | 304 | ! added by LK |
---|
| 305 | if (iaero_h2so4.ne.0) then |
---|
[728] | 306 | iaer=iaero_h2so4 |
---|
[135] | 307 | |
---|
[253] | 308 | ! 1. Initialization |
---|
[787] | 309 | aerosol(1:ngrid,1:nlayermx,iaer)=0.0 |
---|
[253] | 310 | |
---|
| 311 | |
---|
| 312 | ! 2. Opacity calculation |
---|
| 313 | |
---|
[726] | 314 | ! expfactor=0. |
---|
[728] | 315 | DO l=1,nlayer-1 |
---|
| 316 | DO ig=1,ngrid |
---|
| 317 | ! Typical mixing ratio profile |
---|
[253] | 318 | |
---|
[728] | 319 | zp=(preff/pplay(ig,l))**(70./30) !emulating topdust |
---|
| 320 | expfactor=max(exp(0.007*(1.-max(zp,1.))),1.e-3) |
---|
[135] | 321 | |
---|
[726] | 322 | ! Vertical scaling function |
---|
[728] | 323 | aerosol(ig,l,iaer)= (pplev(ig,l)-pplev(ig,l+1))*expfactor |
---|
[135] | 324 | |
---|
[728] | 325 | ENDDO |
---|
| 326 | ENDDO |
---|
| 327 | tauh2so4tmp(1:ngrid)=0. |
---|
| 328 | DO l=1,nlayer |
---|
| 329 | DO ig=1,ngrid |
---|
| 330 | tauh2so4tmp(ig) = tauh2so4tmp(ig) + aerosol(ig,l,iaer) |
---|
| 331 | ENDDO |
---|
| 332 | ENDDO |
---|
| 333 | DO l=1,nlayer-1 |
---|
| 334 | DO ig=1,ngrid |
---|
| 335 | aerosol(ig,l,iaer) = max(1E-20, & |
---|
[726] | 336 | 1 & |
---|
| 337 | * pplev(ig,1) / preff & |
---|
| 338 | * aerosol(ig,l,iaer) & |
---|
| 339 | / tauh2so4tmp(ig)) |
---|
| 340 | |
---|
| 341 | ENDDO |
---|
[728] | 342 | ENDDO |
---|
[726] | 343 | |
---|
| 344 | ! 1/700. is assuming a "sulfurtau" of 1 |
---|
| 345 | ! Sulfur aerosol routine to be improved. |
---|
| 346 | ! aerosol0 = & |
---|
| 347 | ! ( 0.75 * QREFvis3d(ig,l,iaer) / & |
---|
| 348 | ! ( rho_h2so4 * reffrad(ig,l,iaer) ) ) * & |
---|
| 349 | ! ( pq(ig,l,i_h2so4) + 1.E-9 ) * & |
---|
| 350 | ! ( pplev(ig,l) - pplev(ig,l+1) ) / g |
---|
| 351 | ! aerosol0 = max(aerosol0,1.e-9) |
---|
| 352 | ! aerosol0 = min(aerosol0,L_TAUMAX) |
---|
| 353 | ! aerosol(ig,l,iaer) = aerosol0 |
---|
| 354 | |
---|
| 355 | ! ENDDO |
---|
| 356 | ! ENDDO |
---|
| 357 | end if |
---|
| 358 | |
---|
| 359 | |
---|
[1026] | 360 | ! --------------------------------------------------------- |
---|
| 361 | !================================================================== |
---|
| 362 | ! Two-layer aerosols (unknown composition) |
---|
| 363 | ! S. Guerlet (2013) |
---|
| 364 | !================================================================== |
---|
[726] | 365 | |
---|
[1026] | 366 | if (iaero_back2lay .ne.0) then |
---|
| 367 | iaer=iaero_back2lay |
---|
| 368 | ! 1. Initialization |
---|
| 369 | aerosol(1:ngrid,1:nlayermx,iaer)=0.0 |
---|
| 370 | ! 2. Opacity calculation |
---|
| 371 | DO ig=1,ngrid |
---|
| 372 | DO l=1,nlayer-1 |
---|
| 373 | aerosol(ig,l,iaer) = ( pplev(ig,l) - pplev(ig,l+1) ) |
---|
| 374 | !! 1. below tropospheric layer: no aerosols |
---|
| 375 | IF (pplev(ig,l) .gt. pres_bottom_tropo) THEN |
---|
| 376 | aerosol(ig,l,iaer) = 0.*aerosol(ig,l,iaer) |
---|
| 377 | !! 2. tropo layer |
---|
| 378 | ELSEIF (pplev(ig,l) .le. pres_bottom_tropo .and. pplev(ig,l) .ge. pres_top_tropo) THEN |
---|
| 379 | aerosol(ig,l,iaer) = obs_tau_col_tropo*aerosol(ig,l,iaer) |
---|
| 380 | !! 3. linear transition |
---|
| 381 | ELSEIF (pplev(ig,l) .lt. pres_top_tropo .and. pplev(ig,l) .gt. pres_bottom_strato) THEN |
---|
| 382 | expfactor=log(obs_tau_col_strato/obs_tau_col_tropo)/log(pres_bottom_strato/pres_top_tropo) |
---|
[1132] | 383 | aerosol(ig,l,iaer)= obs_tau_col_tropo*((pplev(ig,l)/pres_top_tropo)**expfactor)*aerosol(ig,l,iaer)/1.5 |
---|
[1026] | 384 | !! 4. strato layer |
---|
| 385 | ELSEIF (pplev(ig,l) .le. pres_bottom_strato .and. pplev(ig,l) .gt. pres_top_strato) THEN |
---|
| 386 | aerosol(ig,l,iaer)= obs_tau_col_strato*aerosol(ig,l,iaer) |
---|
| 387 | !! 5. above strato layer: no aerosols |
---|
| 388 | ELSEIF (pplev(ig,l) .lt. pres_top_strato) THEN |
---|
| 389 | aerosol(ig,l,iaer) = 0.*aerosol(ig,l,iaer) |
---|
| 390 | ENDIF |
---|
| 391 | ENDDO |
---|
| 392 | ENDDO |
---|
| 393 | |
---|
| 394 | ! 3. Re-normalize to observed total column |
---|
| 395 | tau_col(:)=0.0 |
---|
| 396 | DO l=1,nlayer |
---|
| 397 | DO ig=1,ngrid |
---|
| 398 | tau_col(ig) = tau_col(ig) & |
---|
| 399 | + aerosol(ig,l,iaer)/(obs_tau_col_tropo+obs_tau_col_strato) |
---|
| 400 | ENDDO |
---|
| 401 | ENDDO |
---|
| 402 | |
---|
| 403 | DO ig=1,ngrid |
---|
| 404 | DO l=1,nlayer-1 |
---|
| 405 | aerosol(ig,l,iaer)=aerosol(ig,l,iaer)/tau_col(ig) |
---|
| 406 | ENDDO |
---|
| 407 | ENDDO |
---|
| 408 | |
---|
| 409 | |
---|
| 410 | end if ! if Two-layer aerosols |
---|
| 411 | |
---|
| 412 | |
---|
[135] | 413 | ! -------------------------------------------------------------------------- |
---|
| 414 | ! Column integrated visible optical depth in each point (used for diagnostic) |
---|
| 415 | |
---|
[253] | 416 | tau_col(:)=0.0 |
---|
[135] | 417 | do iaer = 1, naerkind |
---|
| 418 | do l=1,nlayer |
---|
| 419 | do ig=1,ngrid |
---|
| 420 | tau_col(ig) = tau_col(ig) + aerosol(ig,l,iaer) |
---|
| 421 | end do |
---|
| 422 | end do |
---|
| 423 | end do |
---|
| 424 | |
---|
[787] | 425 | do ig=1,ngrid |
---|
[253] | 426 | do l=1,nlayer |
---|
| 427 | do iaer = 1, naerkind |
---|
| 428 | if(aerosol(ig,l,iaer).gt.1.e3)then |
---|
| 429 | print*,'WARNING: aerosol=',aerosol(ig,l,iaer) |
---|
| 430 | print*,'at ig=',ig,', l=',l,', iaer=',iaer |
---|
| 431 | print*,'QREFvis3d=',QREFvis3d(ig,l,iaer) |
---|
| 432 | print*,'reffrad=',reffrad(ig,l,iaer) |
---|
| 433 | endif |
---|
| 434 | end do |
---|
| 435 | end do |
---|
| 436 | end do |
---|
| 437 | |
---|
[787] | 438 | do ig=1,ngrid |
---|
[253] | 439 | if(tau_col(ig).gt.1.e3)then |
---|
| 440 | print*,'WARNING: tau_col=',tau_col(ig) |
---|
| 441 | print*,'at ig=',ig |
---|
| 442 | print*,'aerosol=',aerosol(ig,:,:) |
---|
| 443 | print*,'QREFvis3d=',QREFvis3d(ig,:,:) |
---|
| 444 | print*,'reffrad=',reffrad(ig,:,:) |
---|
| 445 | endif |
---|
| 446 | end do |
---|
[135] | 447 | return |
---|
| 448 | end subroutine aeropacity |
---|
| 449 | |
---|