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

Last change on this file since 2166 was 2166, checked in by abierjon, 5 years ago

MARS GCM:
Handle case when injecting without dust scenario (iaervar=1) in compute_dtau_mod
AB

File size: 3.4 KB
Line 
1        MODULE compute_dtau_mod
2
3        IMPLICIT NONE
4
5        include "callkeys.h"
6
7        REAL,SAVE :: ti_injection_sol ! time of beginning injection
8        REAL,SAVE :: tf_injection_sol ! time of end injection
9        REAL,PARAMETER :: t_scenario_sol=14/24.   ! time of tauref_scenario
10
11        CONTAINS
12
13        SUBROUTINE compute_dtau(ngrid,nlayer,                           &
14       &                        zday,pplev,tauref,                      &
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       
23        IMPLICIT NONE
24       
25        include "callkeys.h"
26       
27        INTEGER, INTENT(in) :: ngrid
28        INTEGER, INTENT(in) :: nlayer
29        REAL, INTENT(in) :: zday ! date at lon=0, in fraction of sols
30        REAL, INTENT(in) :: pplev(ngrid,nlayer+1) ! pressure (Pa)
31        REAL, INTENT(in) :: tauref(ngrid) ! Computed dust opacity at 610Pa
32        REAL, INTENT(in) :: ptimestep
33        REAL, INTENT(in) :: local_time(ngrid)
34        REAL, INTENT(out) :: dustliftday(ngrid) ! Dust injection rate (s-1)
35       
36        INTEGER :: ig, l
37        INTEGER, SAVE :: nb_daystep ! nomber of step a day
38        REAL :: tauref_scenario(ngrid) ! from dust scenario
39        REAL :: zday_scenario
40        REAL,ALLOCATABLE,SAVE :: local_time_prev(:)
41        REAL,ALLOCATABLE,SAVE :: dtau(:) ! Dust opacity increment (at 610Pa)
42        REAL,PARAMETER :: odpref=610. !DOD reference pressure (Pa)
43       
44        LOGICAL, SAVE :: firstcall=.TRUE. ! signals first call to physics
45       
46       
47        IF(firstcall)THEN
48                ALLOCATE(local_time_prev(ngrid))
49                ALLOCATE(dtau(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                   dtau(ig)=0.
55                ENDDO
56                nb_daystep=(daysec/dtphys)
57                ! Local time in sol fraction
58                ti_injection_sol=ti_injection/24.
59                tf_injection_sol=tf_injection/24.
60                firstcall=.FALSE.
61        ENDIF
62       
63        ! 1. Obtain tauref_scenario from dust scenario at zday+1
64        if (iaervar.eq.1) then
65          tauref_scenario = tauvis
66        else
67          zday_scenario=zday-modulo(zday,1.) ! integer value of the day: the scenario opacity is measured at 14:00
68          zday_scenario=zday_scenario+1      ! opacity of the dust scenario is read the day after
69          call read_dust_scenario(ngrid,nlayer,zday_scenario,pplev,     &
70       &                          tauref_scenario)
71        endif
72       ! for diagnostics
73        call WRITEDIAGFI(ngrid,"tauref_scenario","tauref_scenario",     &
74       &                  "",2,tauref_scenario)
75
76        ! 2. Compute dtau
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)=tauref_scenario(ig)-tauref(ig)
81         ENDIF
82         IF (dtau(ig).LT.0) THEN
83                  dtau(ig)=0.
84         ENDIF
85         
86        ! 3. Use dtau to compute dustliftday
87         dustliftday(ig)=coeff_injection*                               &
88       &                (dtau(ig)*pplev(ig,1)/odpref)                   &
89       &                /(daysec*(tf_injection_sol-ti_injection_sol))
90        ENDDO
91         
92        ! 4. Save local time
93        local_time_prev(1:ngrid)=local_time(1:ngrid)
94       
95        end subroutine compute_dtau
96
97        end module compute_dtau_mod
Note: See TracBrowser for help on using the repository browser.