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