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