source: trunk/LMDZ.MARS/libf/phymars/dust_rad_adjust_mod.F90 @ 2597

Last change on this file since 2597 was 2584, checked in by romain.vande, 4 years ago

Second stage of implementation of Open_MP in the physic.
Run with callrad=.true.

File size: 5.5 KB
Line 
1module dust_rad_adjust_mod
2
3implicit none
4
5real,save,allocatable :: dust_rad_adjust_prev(:) ! adjustment coefficient
6                         ! computed when at current t_scenario
7real,save,allocatable :: dust_rad_adjust_next(:) ! adjustment coefficient
8                         ! computed for t_scenario of the next sol
9
10!$OMP THREADPRIVATE(dust_rad_adjust_prev,dust_rad_adjust_next)
11
12contains
13
14  subroutine compute_dust_rad_adjust(ngrid,nlayer,zday,pplev, &
15                                     taudust,dust_rad_adjust)
16 
17  use geometry_mod, only: longitude_deg
18  use time_phylmdz_mod, only: dtphys, daysec
19  use dust_param_mod, only: odpref, t_scenario_sol
20 
21  implicit none
22 
23  integer,intent(in) :: ngrid ! number of atmospheric columns
24  integer,intent(in) :: nlayer ! number of atmospheric levels
25  real,intent(in) :: zday ! tim (in sols and fraction thereof)
26  real,intent(in) :: pplev(ngrid,nlayer+1) ! pressure (Pa) at layer boundaries
27  real,intent(in) :: taudust(ngrid) ! visible dust columns opacity in the GCM
28  real,intent(out) :: dust_rad_adjust(ngrid) ! radiative adjustment coefficient
29                      ! for dust
30 
31  real,allocatable,save :: local_time(:) ! LT at current physics time step
32  real,allocatable,save :: local_time_prevdt(:) ! LT at previous physics time step
33  real :: zday_prevdt !value of zday at previous physics time step
34  real,save :: zday_scenario ! to fetch dod values from the scenario
35  real,save :: zday_scenario_next ! to fetch dod values from the scenario the next day
36  logical,save :: firstcall=.true.
37  integer :: ig
38!  real,allocatable,save :: tau_pref_scenario(:)
39  real,allocatable,save :: tau_pref_scenario_next(:)
40  real :: weight ! interpolation weight
41  real,save :: zday_prev_call=-666. ! stored value of zday from previous call
42 
43  ! 0. preliminary stuff
44  ! NB: this routine may be called multiple times per physics
45  ! so we have to save some arrays to store the information and not
46  ! recompute it for each call
47 
48  if (firstcall) then
49    write(*,*) "compute_dust_rad_adjust: dust scenario assumed exact at", &
50               " time(sol)=",t_scenario_sol
51    allocate(local_time(ngrid))
52    allocate(local_time_prevdt(ngrid))
53!    allocate(tau_pref_scenario(ngrid))
54    allocate(tau_pref_scenario_next(ngrid))
55    firstcall=.false.
56  endif ! of if firstcall
57 
58  ! 1. Compute local times (in sol fraction), if not already done
59  if (zday/=zday_prev_call) then
60    local_time(1:ngrid)=modulo(1.+(zday-INT(zday)) + &
61                         (longitude_deg(1:ngrid)/15)/24,1.)
62    zday_prevdt=zday-dtphys/daysec
63    local_time_prevdt(1:ngrid)=modulo(1.+(zday_prevdt-INT(zday_prevdt)) + &
64                         (longitude_deg(1:ngrid)/15)/24,1.)
65   
66    zday_scenario=zday-modulo(zday,1.) ! integer value of the day: the scenario
67    ! opacity is assumed to be measured at 2pm but stored at nidnight
68    zday_scenario_next=zday_scenario+1
69  endif ! of if (zday/=zday_prev_call)
70 
71  ! 2. Load dust opacities for zday_scenario and zday_scenario_next
72  !    if not already done
73  if (zday/=zday_prev_call) then
74!    call read_dust_scenario(ngrid,nlayer,zday_scenario,pplev,     &
75!                                         tau_pref_scenario)
76    call read_dust_scenario(ngrid,nlayer,zday_scenario_next,pplev, &
77                                         tau_pref_scenario_next)
78  endif ! of if (zday/=zday_prev_call)
79
80  ! 3. Update dust_rad_adjust_* for grid points which just reached 2pm
81  ! but only when this routine is called for the first time
82  ! during this time step
83  if (zday/=zday_prev_call) then
84   do ig=1,ngrid
85    if ((local_time(ig).ge.t_scenario_sol).and. &
86                 (local_time_prevdt(ig).lt.(t_scenario_sol))) then
87      ! store previous "next" as "prev" (NB we could also decide to recompute
88      ! it using the current taudust...)
89      dust_rad_adjust_prev(ig)=dust_rad_adjust_next(ig)
90      ! compute new target based on current dust opacity
91      dust_rad_adjust_next(ig)=tau_pref_scenario_next(ig)* &
92                               pplev(ig,1)/odpref/taudust(ig)
93    endif
94   enddo
95  endif ! of if (zday/=zday_prev_call)
96 
97  ! 4. Compute dust_rad_adjust using linear interpolation
98  ! between dust_rad_adjust_prev and dust_rad_adjust_next
99  do ig=1,ngrid
100    ! prev and next are separated by a sol exactly
101    ! we just need the distance (in sol) between current local time
102    ! and 2pm the day before
103    if (local_time(ig).ge.t_scenario_sol) then
104      ! we are between t_scenario_sol and midnight
105      weight=local_time(ig)-t_scenario_sol
106    else
107      ! we are between midnight and t_scenario_sol of the next day
108      weight=(1.-t_scenario_sol)+local_time(ig)
109    endif
110    dust_rad_adjust(ig)=dust_rad_adjust_prev(ig)+ &
111                        weight* &
112                        (dust_rad_adjust_next(ig)-dust_rad_adjust_prev(ig))
113  enddo! of do=ig=1,ngrid
114 
115  ! update zday_prev_call
116  zday_prev_call=zday
117 
118  end subroutine compute_dust_rad_adjust
119
120!=======================================================================
121! Initialization of the module variables
122
123  subroutine ini_dust_rad_adjust_mod(ngrid)
124       
125  implicit none
126       
127  integer, intent(in) :: ngrid
128       
129  allocate(dust_rad_adjust_prev(ngrid))
130  allocate(dust_rad_adjust_next(ngrid))
131       
132  end subroutine ini_dust_rad_adjust_mod
133       
134  subroutine end_dust_rad_adjust_mod
135       
136  implicit none
137       
138  if (allocated(dust_rad_adjust_prev)) deallocate(dust_rad_adjust_prev)
139  if (allocated(dust_rad_adjust_next)) deallocate(dust_rad_adjust_next)
140
141  end subroutine end_dust_rad_adjust_mod
142
143end module dust_rad_adjust_mod
Note: See TracBrowser for help on using the repository browser.