| 1 | module dust_scaling_mod |
|---|
| 2 | |
|---|
| 3 | implicit none |
|---|
| 4 | |
|---|
| 5 | contains |
|---|
| 6 | |
|---|
| 7 | subroutine compute_dustscaling(ngrid,nlayer,naerkind,naerdust, & |
|---|
| 8 | zday,pplev, & |
|---|
| 9 | tau_pref_scenario,IRtoVIScoef, & |
|---|
| 10 | tauscaling, & |
|---|
| 11 | dust_rad_adjust,aerosol) |
|---|
| 12 | |
|---|
| 13 | use dust_param_mod, only: dustscaling_mode, odpref |
|---|
| 14 | use dust_rad_adjust_mod, only: compute_dust_rad_adjust |
|---|
| 15 | use dimradmars_mod, only: iaerdust ! dust aerosol indexes |
|---|
| 16 | |
|---|
| 17 | implicit none |
|---|
| 18 | |
|---|
| 19 | integer,intent(in) :: ngrid ! number of atmospheric columns |
|---|
| 20 | integer,intent(in) :: nlayer ! number of atmospheric layers |
|---|
| 21 | integer,intent(in) :: naerkind ! total number of aerosols |
|---|
| 22 | integer,intent(in) :: naerdust ! number of dust aerosols |
|---|
| 23 | real,intent(in) :: zday |
|---|
| 24 | real,intent(in) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa) |
|---|
| 25 | real,intent(in) :: tau_pref_scenario(ngrid) ! prescribed visible dust |
|---|
| 26 | ! opacity column at odpref reference pressure |
|---|
| 27 | real,intent(in) :: IRtoVIScoef(ngrid) ! conversion coefficient to apply on |
|---|
| 28 | ! scenario absorption IR (9.3um) CDOD |
|---|
| 29 | ! = tau_pref_gcm_VIS / tau_pref_gcm_IR |
|---|
| 30 | real,intent(out) :: tauscaling(ngrid) ! dust scaling factor |
|---|
| 31 | real,intent(inout) :: dust_rad_adjust(ngrid) ! Radiative adjustment |
|---|
| 32 | ! factor for dust |
|---|
| 33 | real,intent(inout) :: aerosol(ngrid,nlayer,naerkind) ! opacities |
|---|
| 34 | |
|---|
| 35 | integer :: ig, l , iaer |
|---|
| 36 | real :: taudust(ngrid) |
|---|
| 37 | real,save :: zday_prev_call=-666. ! stored value of zday from previous call |
|---|
| 38 | !$OMP THREADPRIVATE(zday_prev_call) |
|---|
| 39 | |
|---|
| 40 | ! 1. compute/set tauscaling |
|---|
| 41 | |
|---|
| 42 | if (dustscaling_mode /= 1) then |
|---|
| 43 | ! simple "freedust" case, no effective rescaling using tauscaling, ever |
|---|
| 44 | tauscaling(:) = 1 |
|---|
| 45 | endif |
|---|
| 46 | |
|---|
| 47 | if (dustscaling_mode == 1) then |
|---|
| 48 | ! Compute dust column opacity using aerosol() dusts |
|---|
| 49 | taudust(:) = 0 |
|---|
| 50 | do iaer=1,naerdust ! loop on all dust aerosols |
|---|
| 51 | do l=1,nlayer |
|---|
| 52 | do ig=1,ngrid |
|---|
| 53 | taudust(ig)=taudust(ig)+aerosol(ig,l,iaerdust(iaer)) |
|---|
| 54 | enddo |
|---|
| 55 | enddo |
|---|
| 56 | enddo |
|---|
| 57 | |
|---|
| 58 | elseif (dustscaling_mode == 2) then |
|---|
| 59 | ! Compute dust column opacity using only background dust |
|---|
| 60 | taudust(:) = 0 |
|---|
| 61 | do l=1,nlayer |
|---|
| 62 | do ig=1,ngrid |
|---|
| 63 | taudust(ig)=taudust(ig)+aerosol(ig,l,iaerdust(1)) |
|---|
| 64 | enddo |
|---|
| 65 | enddo |
|---|
| 66 | |
|---|
| 67 | endif ! of if (dustscaling_mode == 1) elseif (dustscaling_mode == 2) |
|---|
| 68 | |
|---|
| 69 | ! 2. compute the scaling factors (tauscaling or dust_rad_adjust) |
|---|
| 70 | if (dustscaling_mode==1) then |
|---|
| 71 | ! GCM v5.3 style: tauscaling is computed so that |
|---|
| 72 | ! aerosol() opacities correspond to the prescribed tau_pref_scenario() |
|---|
| 73 | tauscaling(:)=tau_pref_scenario(:)*pplev(:,1)/odpref/taudust(:) |
|---|
| 74 | elseif (dustscaling_mode==2) then |
|---|
| 75 | ! GCM v6 style, compute dust_rad_adjust |
|---|
| 76 | ! only when this routine is called for the first time |
|---|
| 77 | ! during this time step |
|---|
| 78 | if (zday/=zday_prev_call) then |
|---|
| 79 | call compute_dust_rad_adjust(ngrid,nlayer,zday,pplev, & |
|---|
| 80 | taudust,IRtoVIScoef, & |
|---|
| 81 | dust_rad_adjust) |
|---|
| 82 | endif ! of if (zday/=zday_prev_call) |
|---|
| 83 | |
|---|
| 84 | ! update zday_prev_call |
|---|
| 85 | zday_prev_call=zday |
|---|
| 86 | endif |
|---|
| 87 | |
|---|
| 88 | ! 3. Apply dust aerosol opacities rescaling |
|---|
| 89 | if (dustscaling_mode <=1) then |
|---|
| 90 | do iaer=1,naerdust |
|---|
| 91 | do l=1,nlayer |
|---|
| 92 | do ig=1,ngrid |
|---|
| 93 | aerosol(ig,l,iaerdust(iaer)) = max(1E-20, & |
|---|
| 94 | aerosol(ig,l,iaerdust(iaer))* tauscaling(ig)) |
|---|
| 95 | enddo |
|---|
| 96 | enddo |
|---|
| 97 | enddo |
|---|
| 98 | else ! duscaling_mode==2, use dust_rad_adjust |
|---|
| 99 | do iaer=1,naerdust |
|---|
| 100 | do l=1,nlayer |
|---|
| 101 | do ig=1,ngrid |
|---|
| 102 | aerosol(ig,l,iaerdust(iaer)) = max(1E-20, & |
|---|
| 103 | aerosol(ig,l,iaerdust(iaer))*dust_rad_adjust(ig)) |
|---|
| 104 | enddo |
|---|
| 105 | enddo |
|---|
| 106 | enddo |
|---|
| 107 | endif |
|---|
| 108 | end subroutine compute_dustscaling |
|---|
| 109 | |
|---|
| 110 | end module dust_scaling_mod |
|---|