1 | MODULE compute_dtau_mod |
---|
2 | |
---|
3 | IMPLICIT NONE |
---|
4 | |
---|
5 | REAL,SAVE :: ti_injection_sol ! time of beginning injection |
---|
6 | REAL,SAVE :: tf_injection_sol ! time of end injection |
---|
7 | |
---|
8 | REAL,SAVE,ALLOCATABLE :: dtau(:) ! Dust opacity difference (at 610Pa) |
---|
9 | ! between GCM and dust scenario |
---|
10 | |
---|
11 | !$OMP THREADPRIVATE(ti_injection_sol,tf_injection_sol,dtau) |
---|
12 | |
---|
13 | CONTAINS |
---|
14 | |
---|
15 | SUBROUTINE compute_dtau(ngrid,nlayer, & |
---|
16 | zday,pplev,tau_pref_gcm, & |
---|
17 | ptimestep,local_time,IRtoVIScoef, & |
---|
18 | dustliftday) |
---|
19 | |
---|
20 | USE geometry_mod, only: longitude_deg |
---|
21 | USE time_phylmdz_mod, only: dtphys, daysec |
---|
22 | USE comcstfi_h, only: g |
---|
23 | USE tracer_mod, only: alpha_lift,igcm_dust_mass,igcm_dust_number |
---|
24 | USE dimradmars_mod, only: tauvis |
---|
25 | USE dust_param_mod, only: odpref, t_scenario_sol |
---|
26 | USE read_dust_scenario_mod, only: read_dust_scenario |
---|
27 | |
---|
28 | IMPLICIT NONE |
---|
29 | |
---|
30 | include "callkeys.h" |
---|
31 | |
---|
32 | INTEGER, INTENT(in) :: ngrid |
---|
33 | INTEGER, INTENT(in) :: nlayer |
---|
34 | REAL, INTENT(in) :: zday ! date at lon=0, in fraction of sols |
---|
35 | REAL, INTENT(in) :: pplev(ngrid,nlayer+1) ! pressure (Pa) |
---|
36 | REAL, INTENT(in) :: tau_pref_gcm(ngrid) ! Visible dust opacity column |
---|
37 | ! at 610Pa as computed in the GCM |
---|
38 | REAL, INTENT(in) :: ptimestep |
---|
39 | REAL, INTENT(in) :: local_time(ngrid) |
---|
40 | REAL, INTENT(in) :: IRtoVIScoef(ngrid) ! conversion coefficient to apply on |
---|
41 | ! scenario absorption IR (9.3um) CDOD |
---|
42 | ! = tau_pref_gcm_VIS / tau_pref_gcm_IR |
---|
43 | REAL, INTENT(out) :: dustliftday(ngrid) ! Dust injection rate (s-1) |
---|
44 | |
---|
45 | INTEGER :: ig, l |
---|
46 | INTEGER, SAVE :: nb_daystep ! nomber of step a day |
---|
47 | REAL :: tau_pref_target(ngrid) ! dust opacity column at odpref=610 Pa |
---|
48 | ! as extracted from dust scenario |
---|
49 | REAL :: zday_scenario |
---|
50 | REAL,ALLOCATABLE,SAVE :: local_time_prev(:) |
---|
51 | |
---|
52 | LOGICAL, SAVE :: firstcall=.TRUE. ! signals first call to physics |
---|
53 | |
---|
54 | !$OMP THREADPRIVATE(nb_daystep,local_time_prev,firstcall) |
---|
55 | |
---|
56 | |
---|
57 | IF(firstcall)THEN |
---|
58 | ALLOCATE(local_time_prev(ngrid)) |
---|
59 | DO ig=1,ngrid |
---|
60 | local_time_prev(ig)=modulo(1.+(zday-ptimestep/daysec)& |
---|
61 | -INT(zday-ptimestep/daysec) & |
---|
62 | +(longitude_deg(ig)/15)/24,1.) |
---|
63 | ENDDO |
---|
64 | nb_daystep=(daysec/dtphys) |
---|
65 | ! Local time in sol fraction |
---|
66 | ti_injection_sol=ti_injection/24. |
---|
67 | tf_injection_sol=tf_injection/24. |
---|
68 | firstcall=.FALSE. |
---|
69 | ENDIF |
---|
70 | |
---|
71 | ! 1. Obtain tau_pref_target from dust scenario at zday+1 |
---|
72 | if (iaervar.eq.1) then |
---|
73 | tau_pref_target = tauvis |
---|
74 | else |
---|
75 | zday_scenario=zday-modulo(zday,1.) ! integer value of the day: the scenario opacity is measured at 14:00 |
---|
76 | zday_scenario=zday_scenario+1 ! opacity of the dust scenario is read the day after |
---|
77 | call read_dust_scenario(ngrid,nlayer,zday_scenario,pplev, & |
---|
78 | IRtoVIScoef, & |
---|
79 | tau_pref_target) |
---|
80 | endif |
---|
81 | ! for diagnostics |
---|
82 | call WRITEDIAGFI(ngrid,"tau_pref_target", & |
---|
83 | "target visible dust opacity column at 610Pa", & |
---|
84 | "",2,tau_pref_target) |
---|
85 | |
---|
86 | ! 2. Compute dtau() and dustliftday() |
---|
87 | DO ig=1,ngrid |
---|
88 | IF ((local_time(ig).ge.t_scenario_sol).and. & |
---|
89 | (local_time_prev(ig).lt.(t_scenario_sol)))THEN |
---|
90 | dtau(ig)=tau_pref_target(ig)-tau_pref_gcm(ig) |
---|
91 | ENDIF |
---|
92 | |
---|
93 | ! Use dtau (when positive) to compute dustliftday |
---|
94 | IF (dtau(ig).LT.0) THEN |
---|
95 | dustliftday(ig)=0. |
---|
96 | ELSE |
---|
97 | dustliftday(ig)=coeff_injection* & |
---|
98 | (dtau(ig)*pplev(ig,1)/odpref) & |
---|
99 | /(daysec*(tf_injection_sol-ti_injection_sol)) |
---|
100 | ENDIF |
---|
101 | ENDDO ! of DO ig=1,ngrid |
---|
102 | |
---|
103 | ! for diagnostics |
---|
104 | call WRITEDIAGFI(ngrid,"dtau","opacity difference wrt scenario",& |
---|
105 | "",2,dtau) |
---|
106 | call WRITEDIAGFI(ngrid,"dustliftday","dust injection rate", & |
---|
107 | "s-1",2,dustliftday) |
---|
108 | |
---|
109 | ! 4. Save local time |
---|
110 | local_time_prev(1:ngrid)=local_time(1:ngrid) |
---|
111 | |
---|
112 | end subroutine compute_dtau |
---|
113 | |
---|
114 | !======================================================================= |
---|
115 | ! Initialization of the module variables |
---|
116 | |
---|
117 | subroutine ini_compute_dtau_mod(ngrid) |
---|
118 | |
---|
119 | implicit none |
---|
120 | |
---|
121 | integer, intent(in) :: ngrid |
---|
122 | |
---|
123 | allocate(dtau(ngrid)) |
---|
124 | |
---|
125 | end subroutine ini_compute_dtau_mod |
---|
126 | |
---|
127 | subroutine end_compute_dtau_mod |
---|
128 | |
---|
129 | implicit none |
---|
130 | |
---|
131 | if (allocated(dtau)) deallocate(dtau) |
---|
132 | |
---|
133 | end subroutine end_compute_dtau_mod |
---|
134 | |
---|
135 | END MODULE compute_dtau_mod |
---|