source: trunk/LMDZ.MARS/libf/phymars/compute_dtau_mod.F90 @ 2429

Last change on this file since 2429 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: 4.7 KB
Line 
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       CONTAINS
12
13        SUBROUTINE compute_dtau(ngrid,nlayer,                           &
14                                 zday,pplev,tau_pref_gcm,               &
15                                 ptimestep,dustliftday,local_time)
16
17        USE geometry_mod, only: longitude_deg
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
21        USE dimradmars_mod, only: tauvis
22        USE dust_param_mod, only: odpref, t_scenario_sol
23       
24        IMPLICIT NONE
25       
26        include "callkeys.h"
27       
28        INTEGER, INTENT(in) :: ngrid
29        INTEGER, INTENT(in) :: nlayer
30        REAL, INTENT(in) :: zday ! date at lon=0, in fraction of sols
31        REAL, INTENT(in) :: pplev(ngrid,nlayer+1) ! pressure (Pa)
32        REAL, INTENT(in) :: tau_pref_gcm(ngrid) ! Visible dust opacity column
33                            ! at 610Pa as computed in the GCM
34        REAL, INTENT(in) :: ptimestep
35        REAL, INTENT(in) :: local_time(ngrid)
36        REAL, INTENT(out) :: dustliftday(ngrid) ! Dust injection rate (s-1)
37       
38        INTEGER :: ig, l
39        INTEGER, SAVE :: nb_daystep ! nomber of step a day
40        REAL :: tau_pref_target(ngrid) ! dust opacity column at odpref=610 Pa
41                ! as extracted from dust scenario
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)&
52                                      -INT(zday-ptimestep/daysec)       &
53                                      +(longitude_deg(ig)/15)/24,1.)
54                ENDDO
55                nb_daystep=(daysec/dtphys)
56                ! Local time in sol fraction
57                ti_injection_sol=ti_injection/24.
58                tf_injection_sol=tf_injection/24.
59                firstcall=.FALSE.
60        ENDIF
61       
62        ! 1. Obtain tau_pref_target from dust scenario at zday+1
63        if (iaervar.eq.1) then
64          tau_pref_target = tauvis
65        else
66          zday_scenario=zday-modulo(zday,1.) ! integer value of the day: the scenario opacity is measured at 14:00
67          zday_scenario=zday_scenario+1      ! opacity of the dust scenario is read the day after
68          call read_dust_scenario(ngrid,nlayer,zday_scenario,pplev,     &
69                                         tau_pref_target)
70        endif
71       ! for diagnostics
72        call WRITEDIAGFI(ngrid,"tau_pref_target", &
73                          "target visible dust opacity column at 610Pa", &
74                          "",2,tau_pref_target)
75
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
80                 dtau(ig)=tau_pref_target(ig)-tau_pref_gcm(ig)
81         ENDIF
82
83        ! Use dtau (when positive) to compute dustliftday
84         IF (dtau(ig).LT.0) THEN
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))
90         ENDIF
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)
98         
99        ! 4. Save local time
100        local_time_prev(1:ngrid)=local_time(1:ngrid)
101       
102        end subroutine compute_dtau
103
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
Note: See TracBrowser for help on using the repository browser.