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