| 1 | module aeropacity_mod |
|---|
| 2 | |
|---|
| 3 | implicit none |
|---|
| 4 | |
|---|
| 5 | contains |
|---|
| 6 | |
|---|
| 7 | Subroutine aeropacity(ngrid,nlayer,nq,pplay,pplev,pt, pq, & |
|---|
| 8 | aerosol,reffrad,nueffrad, QREFvis3d,QREFir3d,tau_col, & |
|---|
| 9 | cloudfrac,totcloudfrac,clearsky) |
|---|
| 10 | |
|---|
| 11 | use radinc_h, only : L_TAUMAX,naerkind |
|---|
| 12 | use aerosol_mod, only: iaero_haze, i_haze, iaero_generic, & |
|---|
| 13 | noaero |
|---|
| 14 | USE tracer_h, only: noms,rho_n2,rho_ice,rho_q,mmol |
|---|
| 15 | use comcstfi_mod, only: g, pi, mugaz, avocado |
|---|
| 16 | use geometry_mod, only: latitude |
|---|
| 17 | use callkeys_mod, only: aerofixn2,kastprof,cloudlvl, & |
|---|
| 18 | pres_bottom_tropo,pres_top_tropo,obs_tau_col_tropo, & |
|---|
| 19 | pres_bottom_strato,pres_top_strato,obs_tau_col_strato, & |
|---|
| 20 | nlayaero, aeronlay_tauref, aeronlay_choice, & |
|---|
| 21 | aeronlay_pbot, aeronlay_ptop, aeronlay_sclhght, & |
|---|
| 22 | aerogeneric |
|---|
| 23 | use generic_tracer_index_mod, only: generic_tracer_index |
|---|
| 24 | implicit none |
|---|
| 25 | |
|---|
| 26 | !================================================================== |
|---|
| 27 | ! |
|---|
| 28 | ! Purpose |
|---|
| 29 | ! ------- |
|---|
| 30 | ! Compute aerosol optical depth in each gridbox. |
|---|
| 31 | ! |
|---|
| 32 | ! Authors |
|---|
| 33 | ! ------- |
|---|
| 34 | ! F. Forget |
|---|
| 35 | ! F. Montmessin (water ice scheme) |
|---|
| 36 | ! update J.-B. Madeleine (2008) |
|---|
| 37 | ! dust removal, simplification by Robin Wordsworth (2009) |
|---|
| 38 | ! Generic n-layer aerosol - J. Vatant d'Ollone (2020) |
|---|
| 39 | ! Radiative Generic Condensable Species - Lucas Teinturier (2022) |
|---|
| 40 | ! |
|---|
| 41 | ! Input |
|---|
| 42 | ! ----- |
|---|
| 43 | ! ngrid Number of horizontal gridpoints |
|---|
| 44 | ! nlayer Number of layers |
|---|
| 45 | ! nq Number of tracers |
|---|
| 46 | ! pplev Pressure (Pa) at each layer boundary |
|---|
| 47 | ! pq Aerosol mixing ratio |
|---|
| 48 | ! reffrad(ngrid,nlayer,naerkind) Aerosol effective radius |
|---|
| 49 | ! QREFvis3d(ngrid,nlayer,naerkind) \ 3d extinction coefficients |
|---|
| 50 | ! QREFir3d(ngrid,nlayer,naerkind) / at reference wavelengths |
|---|
| 51 | ! |
|---|
| 52 | ! Output |
|---|
| 53 | ! ------ |
|---|
| 54 | ! aerosol Aerosol optical depth in layer l, grid point ig |
|---|
| 55 | ! tau_col Total column optical depth at grid point ig |
|---|
| 56 | ! |
|---|
| 57 | !======================================================================= |
|---|
| 58 | |
|---|
| 59 | INTEGER,INTENT(IN) :: ngrid ! number of atmospheric columns |
|---|
| 60 | INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers |
|---|
| 61 | INTEGER,INTENT(IN) :: nq ! number of tracers |
|---|
| 62 | REAL,INTENT(IN) :: pplay(ngrid,nlayer) ! mid-layer pressure (Pa) |
|---|
| 63 | REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa) |
|---|
| 64 | REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) ! tracers (.../kg_of_air) |
|---|
| 65 | REAL,INTENT(IN) :: pt(ngrid,nlayer) ! mid-layer temperature (K) |
|---|
| 66 | REAL,INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) ! aerosol optical depth |
|---|
| 67 | REAL,INTENT(IN) :: reffrad(ngrid,nlayer,naerkind) ! aerosol effective radius |
|---|
| 68 | REAL,INTENT(IN) :: nueffrad(ngrid,nlayer,naerkind) ! aerosol effective variance |
|---|
| 69 | REAL,INTENT(IN) :: QREFvis3d(ngrid,nlayer,naerkind) ! extinction coefficient in the visible |
|---|
| 70 | REAL,INTENT(IN) :: QREFir3d(ngrid,nlayer,naerkind) |
|---|
| 71 | REAL,INTENT(OUT):: tau_col(ngrid) !column integrated visible optical depth |
|---|
| 72 | ! BENJAMIN MODIFS |
|---|
| 73 | real,intent(in) :: cloudfrac(ngrid,nlayer) ! cloud fraction |
|---|
| 74 | real,intent(out) :: totcloudfrac(ngrid) ! total cloud fraction |
|---|
| 75 | logical,intent(in) :: clearsky |
|---|
| 76 | |
|---|
| 77 | real aerosol0, obs_tau_col_aurora, pm |
|---|
| 78 | real pcloud_deck, cloud_slope |
|---|
| 79 | |
|---|
| 80 | real dp_strato(ngrid) |
|---|
| 81 | real dp_tropo(ngrid) |
|---|
| 82 | real dp_layer(ngrid) |
|---|
| 83 | |
|---|
| 84 | INTEGER l,ig,iq,iaer,ia |
|---|
| 85 | |
|---|
| 86 | LOGICAL,SAVE :: firstcall=.true. |
|---|
| 87 | !$OMP THREADPRIVATE(firstcall) |
|---|
| 88 | REAL CBRT |
|---|
| 89 | EXTERNAL CBRT |
|---|
| 90 | |
|---|
| 91 | INTEGER,SAVE :: i_n2ice=0 ! n2 ice |
|---|
| 92 | !$OMP THREADPRIVATE(i_n2ice) |
|---|
| 93 | CHARACTER(LEN=20) :: tracername ! to temporarily store text |
|---|
| 94 | |
|---|
| 95 | real CLFtot |
|---|
| 96 | integer igen_ice,igen_vap ! to store the index of generic tracer |
|---|
| 97 | logical dummy_bool ! dummy boolean just in case we need one |
|---|
| 98 | ! integer i_rgcs_ice(aerogeneric) |
|---|
| 99 | ! for venus clouds |
|---|
| 100 | real :: p_bot,p_top,h_bot,h_top,mode_dens,h_lay |
|---|
| 101 | |
|---|
| 102 | ! identify tracers |
|---|
| 103 | IF (firstcall) THEN |
|---|
| 104 | ia =0 |
|---|
| 105 | write(*,*) "Tracers found in aeropacity:" |
|---|
| 106 | do iq=1,nq |
|---|
| 107 | tracername=noms(iq) |
|---|
| 108 | if (tracername.eq."n2_ice") then |
|---|
| 109 | i_n2ice=iq |
|---|
| 110 | write(*,*) "i_n2ice=",i_n2ice |
|---|
| 111 | |
|---|
| 112 | endif |
|---|
| 113 | enddo |
|---|
| 114 | |
|---|
| 115 | if (noaero) then |
|---|
| 116 | print*, "No active aerosols found in aeropacity" |
|---|
| 117 | else |
|---|
| 118 | print*, "If you would like to use aerosols, make sure any old" |
|---|
| 119 | print*, "start files are updated in newstart using the option" |
|---|
| 120 | print*, "q=0" |
|---|
| 121 | write(*,*) "Active aerosols found in aeropacity:" |
|---|
| 122 | endif |
|---|
| 123 | |
|---|
| 124 | firstcall=.false. |
|---|
| 125 | ENDIF ! of IF (firstcall) |
|---|
| 126 | |
|---|
| 127 | |
|---|
| 128 | ! --------------------------------------------------------- |
|---|
| 129 | !================================================================== |
|---|
| 130 | ! Haze aerosols |
|---|
| 131 | !================================================================== |
|---|
| 132 | |
|---|
| 133 | if (iaero_haze.ne.0) then |
|---|
| 134 | iaer=iaero_haze |
|---|
| 135 | ! 1. Initialization |
|---|
| 136 | aerosol(1:ngrid,1:nlayer,iaer)=0.0 |
|---|
| 137 | ! 2. Opacity calculation |
|---|
| 138 | DO ig=1, ngrid |
|---|
| 139 | DO l=1,nlayer-1 ! to stop the rad tran bug |
|---|
| 140 | ! if fractal, radius doit etre equivalent sphere radius |
|---|
| 141 | aerosol0 = & |
|---|
| 142 | ( 0.75 * QREFvis3d(ig,l,iaer) / & |
|---|
| 143 | ( rho_q(i_haze) * reffrad(ig,l,iaer) ) ) * & |
|---|
| 144 | ( pq(ig,l,i_haze) + 1.E-10 ) * & |
|---|
| 145 | ( pplev(ig,l) - pplev(ig,l+1) ) / g |
|---|
| 146 | aerosol0 = max(aerosol0,1.e-10) |
|---|
| 147 | aerosol0 = min(aerosol0,L_TAUMAX) |
|---|
| 148 | aerosol(ig,l,iaer) = aerosol0 |
|---|
| 149 | ENDDO |
|---|
| 150 | ENDDO |
|---|
| 151 | !QREF est le meme dans toute la colonne par def si size uniforme |
|---|
| 152 | !print*, 'TB17: QREFvis3d=',QREFvis3d(1,:,1) |
|---|
| 153 | !print*, 'TB17: rho_q=',rho_q(i_haze) |
|---|
| 154 | !print*, 'TB17: reffrad=',reffrad(1,:,1) |
|---|
| 155 | !print*, 'TB17: pq=',pq(1,:,i_haze) |
|---|
| 156 | !print*, 'TB17: deltap=',pplev(1,1) - pplev(1,nlayer) |
|---|
| 157 | end if ! if haze aerosols |
|---|
| 158 | |
|---|
| 159 | ! -------------------------------------------------------------------------- |
|---|
| 160 | ! Column integrated visible optical depth in each point (used for diagnostic) |
|---|
| 161 | |
|---|
| 162 | tau_col(:)=0.0 |
|---|
| 163 | do iaer = 1, naerkind |
|---|
| 164 | do l=1,nlayer |
|---|
| 165 | do ig=1,ngrid |
|---|
| 166 | tau_col(ig) = tau_col(ig) + aerosol(ig,l,iaer) |
|---|
| 167 | end do |
|---|
| 168 | end do |
|---|
| 169 | end do |
|---|
| 170 | |
|---|
| 171 | do ig=1,ngrid |
|---|
| 172 | do l=1,nlayer |
|---|
| 173 | do iaer = 1, naerkind |
|---|
| 174 | if(aerosol(ig,l,iaer).gt.1.e3)then |
|---|
| 175 | print*,'WARNING: aerosol=',aerosol(ig,l,iaer) |
|---|
| 176 | print*,'at ig=',ig,', l=',l,', iaer=',iaer |
|---|
| 177 | print*,'QREFvis3d=',QREFvis3d(ig,l,iaer) |
|---|
| 178 | print*,'reffrad=',reffrad(ig,l,iaer) |
|---|
| 179 | endif |
|---|
| 180 | end do |
|---|
| 181 | end do |
|---|
| 182 | end do |
|---|
| 183 | |
|---|
| 184 | do ig=1,ngrid |
|---|
| 185 | if(tau_col(ig).gt.1.e3)then |
|---|
| 186 | print*,'WARNING: tau_col=',tau_col(ig) |
|---|
| 187 | print*,'at ig=',ig |
|---|
| 188 | print*,'aerosol=',aerosol(ig,:,:) |
|---|
| 189 | print*,'QREFvis3d=',QREFvis3d(ig,:,:) |
|---|
| 190 | print*,'reffrad=',reffrad(ig,:,:) |
|---|
| 191 | endif |
|---|
| 192 | end do |
|---|
| 193 | end subroutine aeropacity |
|---|
| 194 | |
|---|
| 195 | end module aeropacity_mod |
|---|