source: trunk/LMDZ.MARS/libf/phymars/dust_scaling_mod.F90 @ 2613

Last change on this file since 2613 was 2417, checked in by emillour, 4 years ago

Mars GCM:
Add a new scheme to handle the dust and its radiative impact. Triggered using
a new flag dustscaling_mode=2 (dustscaling_mod=0: no rescaling at all, and
dustscaling_mode=1: full rescaling using tauscaling, GCMv5.3 style). Rescaling
is then only done on the radiative impact (see dust_scaling_mod.F90) of dust.
Moreover the scaling factor "dust_rad_adjust" is evaluated using the target dust
scenario opacity for the next sol and left to evolve linearly until then to not
impose the diurnal evolution of dust.
In practice, main changes or additions in the code are:

  • renamed flag "tauscaling_mode" as "dustscaling_mode"
  • moved parameter "t_scenario_sol" to "dust_param_mod"
  • adapted "compute_dustscaling" routine in "dust_scaling_mod"
  • added module "dust_rad_adjust_mod"
  • 2D fields "dust_rad_adjust_prev" and "dust_rad_adjust_next" required to compute coefficient "dust_rad_adjust" need to be stored in (re)startfi files

EM

File size: 3.3 KB
RevLine 
[2413]1module dust_scaling_mod
2
3implicit none
4
5contains
6
[2417]7  subroutine compute_dustscaling(ngrid,nlayer,naerkind,naerdust, &
8                                 zday,pplev, &
9                                 tau_pref_scenario,tauscaling, &
10                                 dust_rad_adjust,aerosol)
[2413]11   
[2417]12    use dust_param_mod, only: dustscaling_mode, odpref
13    use dust_rad_adjust_mod, only: compute_dust_rad_adjust
[2413]14    use dimradmars_mod, only: iaerdust ! dust aerosol indexes
15   
16    implicit none
17   
18    integer,intent(in) :: ngrid ! number of atmospheric columns
19    integer,intent(in) :: nlayer ! number of atmospheric layers
20    integer,intent(in) :: naerkind ! total number of aerosols
21    integer,intent(in) :: naerdust ! number of dust aerosols
22    real,intent(in) :: zday
23    real,intent(in) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa)
[2415]24    real,intent(in) :: tau_pref_scenario(ngrid) ! prescribed visible dust
25                       ! opacity column at odpref reference pressure
[2413]26    real,intent(out) :: tauscaling(ngrid) ! dust scaling factor
[2417]27    real,intent(out) :: dust_rad_adjust(ngrid) ! Radiative adjustment
28                          ! factor for dust
[2413]29    real,intent(inout) :: aerosol(ngrid,nlayer,naerkind) ! opacities
30   
31    integer :: ig, l , iaer
32    real :: taudust(ngrid)
[2417]33
[2413]34  ! 1. compute/set tauscaling
35
[2417]36    if (dustscaling_mode /= 1) then
37      ! simple "freedust" case, no effective rescaling using tauscaling, ever
[2413]38      tauscaling(:) = 1
39    endif
40   
[2417]41    if (dustscaling_mode == 1) then
42      ! Compute dust column opacity using aerosol() dusts
[2413]43      taudust(:) = 0
44      do iaer=1,naerdust ! loop on all dust aerosols
45        do l=1,nlayer
46          do ig=1,ngrid
47              taudust(ig)=taudust(ig)+aerosol(ig,l,iaerdust(iaer))
48          enddo
49        enddo
50      enddo
[2417]51
52    elseif (dustscaling_mode == 2) then
53      ! Compute dust column opacity using only background dust
54      taudust(:) = 0
55      do l=1,nlayer
56        do ig=1,ngrid
57            taudust(ig)=taudust(ig)+aerosol(ig,l,iaerdust(1))
58        enddo
59      enddo
60
61    endif ! of if (dustscaling_mode == 1) elseif (dustscaling_mode == 2)
[2413]62     
[2417]63  ! 2. compute the scaling factors (tauscaling or dust_rad_adjust)
64    if (dustscaling_mode==1) then
65      ! GCM v5.3 style: tauscaling is computed so that
66      ! aerosol() opacities correspond to the prescribed tau_pref_scenario()
[2415]67      tauscaling(:)=tau_pref_scenario(:)*pplev(:,1)/odpref/taudust(:)
[2417]68    elseif (dustscaling_mode==2) then
69      ! GCM v6 style, compute dust_rad_adjust
70      call compute_dust_rad_adjust(ngrid,nlayer,zday,pplev, &
71                                   taudust,dust_rad_adjust)
72    endif
73
74  ! 3. Apply dust aerosol opacities rescaling
75    if (dustscaling_mode <=1) then
76      do iaer=1,naerdust
77        do l=1,nlayer
78          do ig=1,ngrid
[2413]79            aerosol(ig,l,iaerdust(iaer)) = max(1E-20, &
80                      aerosol(ig,l,iaerdust(iaer))* tauscaling(ig))
[2417]81          enddo
[2413]82        enddo
83      enddo
[2417]84    else ! duscaling_mode==2, use dust_rad_adjust
85      do iaer=1,naerdust
86        do l=1,nlayer
87          do ig=1,ngrid
88            aerosol(ig,l,iaerdust(iaer)) = max(1E-20, &
89                      aerosol(ig,l,iaerdust(iaer))*dust_rad_adjust(ig))
90          enddo
91        enddo
92      enddo
93    endif
94  end subroutine compute_dustscaling
[2413]95
96end module dust_scaling_mod
Note: See TracBrowser for help on using the repository browser.