source: trunk/LMDZ.MARS/libf/phymars/dust_rad_adjust_mod.F90 @ 2437

Last change on this file since 2437 was 2417, checked in by emillour, 5 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: 5.4 KB
Line 
1module dust_rad_adjust_mod
2
3implicit none
4
5real,save,allocatable :: dust_rad_adjust_prev(:) ! adjustment coefficient
6                         ! computed when at current t_scenario
7real,save,allocatable :: dust_rad_adjust_next(:) ! adjustment coefficient
8                         ! computed for t_scenario of the next sol
9
10contains
11
12  subroutine compute_dust_rad_adjust(ngrid,nlayer,zday,pplev, &
13                                     taudust,dust_rad_adjust)
14 
15  use geometry_mod, only: longitude_deg
16  use time_phylmdz_mod, only: dtphys, daysec
17  use dust_param_mod, only: odpref, t_scenario_sol
18 
19  implicit none
20 
21  integer,intent(in) :: ngrid ! number of atmospheric columns
22  integer,intent(in) :: nlayer ! number of atmospheric levels
23  real,intent(in) :: zday ! tim (in sols and fraction thereof)
24  real,intent(in) :: pplev(ngrid,nlayer+1) ! pressure (Pa) at layer boundaries
25  real,intent(in) :: taudust(ngrid) ! visible dust columns opacity in the GCM
26  real,intent(out) :: dust_rad_adjust(ngrid) ! radiative adjustment coefficient
27                      ! for dust
28 
29  real,allocatable,save :: local_time(:) ! LT at current physics time step
30  real,allocatable,save :: local_time_prevdt(:) ! LT at previous physics time step
31  real :: zday_prevdt !value of zday at previous physics time step
32  real,save :: zday_scenario ! to fetch dod values from the scenario
33  real,save :: zday_scenario_next ! to fetch dod values from the scenario the next day
34  logical,save :: firstcall=.true.
35  integer :: ig
36!  real,allocatable,save :: tau_pref_scenario(:)
37  real,allocatable,save :: tau_pref_scenario_next(:)
38  real :: weight ! interpolation weight
39  real,save :: zday_prev_call=-666. ! stored value of zday from previous call
40 
41  ! 0. preliminary stuff
42  ! NB: this routine may be called multiple times per physics
43  ! so we have to save some arrays to store the information and not
44  ! recompute it for each call
45 
46  if (firstcall) then
47    write(*,*) "compute_dust_rad_adjust: dust scenario assumed exact at", &
48               " time(sol)=",t_scenario_sol
49    allocate(local_time(ngrid))
50    allocate(local_time_prevdt(ngrid))
51!    allocate(tau_pref_scenario(ngrid))
52    allocate(tau_pref_scenario_next(ngrid))
53    firstcall=.false.
54  endif ! of if firstcall
55 
56  ! 1. Compute local times (in sol fraction), if not already done
57  if (zday/=zday_prev_call) then
58    local_time(1:ngrid)=modulo(1.+(zday-INT(zday)) + &
59                         (longitude_deg(1:ngrid)/15)/24,1.)
60    zday_prevdt=zday-dtphys/daysec
61    local_time_prevdt(1:ngrid)=modulo(1.+(zday_prevdt-INT(zday_prevdt)) + &
62                         (longitude_deg(1:ngrid)/15)/24,1.)
63   
64    zday_scenario=zday-modulo(zday,1.) ! integer value of the day: the scenario
65    ! opacity is assumed to be measured at 2pm but stored at nidnight
66    zday_scenario_next=zday_scenario+1
67  endif ! of if (zday/=zday_prev_call)
68 
69  ! 2. Load dust opacities for zday_scenario and zday_scenario_next
70  !    if not already done
71  if (zday/=zday_prev_call) then
72!    call read_dust_scenario(ngrid,nlayer,zday_scenario,pplev,     &
73!                                         tau_pref_scenario)
74    call read_dust_scenario(ngrid,nlayer,zday_scenario_next,pplev, &
75                                         tau_pref_scenario_next)
76  endif ! of if (zday/=zday_prev_call)
77
78  ! 3. Update dust_rad_adjust_* for grid points which just reached 2pm
79  ! but only when this routine is called for the first time
80  ! during this time step
81  if (zday/=zday_prev_call) then
82   do ig=1,ngrid
83    if ((local_time(ig).ge.t_scenario_sol).and. &
84                 (local_time_prevdt(ig).lt.(t_scenario_sol))) then
85      ! store previous "next" as "prev" (NB we could also decide to recompute
86      ! it using the current taudust...)
87      dust_rad_adjust_prev(ig)=dust_rad_adjust_next(ig)
88      ! compute new target based on current dust opacity
89      dust_rad_adjust_next(ig)=tau_pref_scenario_next(ig)* &
90                               pplev(ig,1)/odpref/taudust(ig)
91    endif
92   enddo
93  endif ! of if (zday/=zday_prev_call)
94 
95  ! 4. Compute dust_rad_adjust using linear interpolation
96  ! between dust_rad_adjust_prev and dust_rad_adjust_next
97  do ig=1,ngrid
98    ! prev and next are separated by a sol exactly
99    ! we just need the distance (in sol) between current local time
100    ! and 2pm the day before
101    if (local_time(ig).ge.t_scenario_sol) then
102      ! we are between t_scenario_sol and midnight
103      weight=local_time(ig)-t_scenario_sol
104    else
105      ! we are between midnight and t_scenario_sol of the next day
106      weight=(1.-t_scenario_sol)+local_time(ig)
107    endif
108    dust_rad_adjust(ig)=dust_rad_adjust_prev(ig)+ &
109                        weight* &
110                        (dust_rad_adjust_next(ig)-dust_rad_adjust_prev(ig))
111  enddo! of do=ig=1,ngrid
112 
113  ! update zday_prev_call
114  zday_prev_call=zday
115 
116  end subroutine compute_dust_rad_adjust
117
118!=======================================================================
119! Initialization of the module variables
120
121  subroutine ini_dust_rad_adjust_mod(ngrid)
122       
123  implicit none
124       
125  integer, intent(in) :: ngrid
126       
127  allocate(dust_rad_adjust_prev(ngrid))
128  allocate(dust_rad_adjust_next(ngrid))
129       
130  end subroutine ini_dust_rad_adjust_mod
131       
132  subroutine end_dust_rad_adjust_mod
133       
134  implicit none
135       
136  if (allocated(dust_rad_adjust_prev)) deallocate(dust_rad_adjust_prev)
137  if (allocated(dust_rad_adjust_next)) deallocate(dust_rad_adjust_next)
138
139  end subroutine end_dust_rad_adjust_mod
140
141end module dust_rad_adjust_mod
Note: See TracBrowser for help on using the repository browser.