source: trunk/LMDZ.MARS/libf/phymars/comslope_mod.F90 @ 3325

Last change on this file since 3325 was 3183, checked in by jbclement, 10 months ago

Mars PCM:
Some changes to prepare the introduction of slopes in 1D: in particular, the subroutine "getslopes.F90" and "param_slope.F90" are now inside the module "slope_mod.F90" + Few small cleanings throughout the code.
JBC

File size: 4.5 KB
Line 
1MODULE comslope_mod
2!======================================================================================================================!
3! Subject:
4!---------
5!   Module used for the parametrization of subgrid scale slope
6!----------------------------------------------------------------------------------------------------------------------!
7! Reference:
8!-----------
9!   Lange et al. (2023 in prep.), 'Introduction of Sub-Grid Scale Slope Microclimates in the Mars Planetary Climate Model', JGR
10!
11!======================================================================================================================!
12
13implicit none
14
15integer, save                              :: nslope         ! Number of slopes for the statistic (1)
16integer, save                              :: iflat          ! Flat slope for islope (1)
17real,    save, allocatable, dimension(:)   :: def_slope      ! Bound of the slope statistics / repartitions (°)
18real,    save, allocatable, dimension(:)   :: def_slope_mean ! Mean slope of each bin (°)
19real,    save, allocatable, dimension(:)   :: sky_slope      ! Portion of the sky view by each slopes (1)
20real,    save, allocatable, dimension(:,:) :: subslope_dist  ! Distribution of the slopes (1)
21integer, save, allocatable, dimension(:)   :: major_slope    ! Index of the subslope that occupies most of themesh (1)
22!$OMP THREADPRIVATE(nslope,def_slope,def_slope_mean,sky_slope,subslope_dist,iflat,major_slope)
23
24!=======================================================================
25contains
26!=======================================================================
27
28SUBROUTINE ini_comslope_h(ngrid,nslope_in)
29
30implicit none
31
32integer, intent(in) :: ngrid
33integer, intent(in) :: nslope_in
34
35allocate(def_slope(nslope_in+1))
36allocate(def_slope_mean(nslope_in))
37allocate(sky_slope(nslope_in))
38allocate(subslope_dist(ngrid,nslope_in))
39allocate(major_slope(ngrid))
40
41END SUBROUTINE ini_comslope_h
42
43!=======================================================================
44
45SUBROUTINE end_comslope_h
46
47implicit none
48
49if (allocated(def_slope)) deallocate(def_slope)
50if (allocated(def_slope_mean)) deallocate(def_slope_mean)
51if (allocated(sky_slope)) deallocate(sky_slope)
52if (allocated(subslope_dist)) deallocate(subslope_dist)
53if (allocated(major_slope)) deallocate(major_slope)
54
55END SUBROUTINE end_comslope_h
56
57!=======================================================================
58
59SUBROUTINE compute_meshgridavg(ngrid,nq,albedo_slope,emis_slope,tsurf_slope,qsurf_slope,albedo_meshavg,emis_meshavg,tsurf_meshavg, qsurf_meshavg)
60
61use comcstfi_h, only: pi
62
63implicit none
64
65integer, intent(in) :: ngrid, nq                             ! # points on the physical grid, tracers (1)
66real, dimension(ngrid,2,nslope),  intent(in) :: albedo_slope ! albedo on each sub-grid slope (1)
67real, dimension(ngrid,nslope),    intent(in) :: emis_slope   ! emissivity on each sub-grid slope (1)
68real, dimension(ngrid,nslope),    intent(in) :: tsurf_slope  ! Surface Temperature on each sub-grid slope (K)
69real, dimension(ngrid,nq,nslope), intent(in) :: qsurf_slope  ! Surface Tracers on each sub-grid slope (kg/m^2 sloped)
70real, dimension(ngrid,2),  intent(out) :: albedo_meshavg ! grid box average of the albedo (1)
71real, dimension(ngrid),    intent(out) :: emis_meshavg   ! grid box average of the emissivity (1)
72real, dimension(ngrid),    intent(out) :: tsurf_meshavg  ! grid box average of the surface temperature (K)
73real, dimension(ngrid,nq), intent(out) :: qsurf_meshavg  ! grid box average of the surface tracers (kg/m^2 flat)
74! Local
75integer :: islope, ig, l, iq
76
77!-------------------
78
79albedo_meshavg = 0.
80emis_meshavg = 0.
81tsurf_meshavg = 0.
82qsurf_meshavg = 0.
83
84if (nslope == 1) then
85    albedo_meshavg = albedo_slope(:,:,1)
86    emis_meshavg = emis_slope(:,1)
87    tsurf_meshavg = tsurf_slope(:,1)
88    qsurf_meshavg = qsurf_slope(:,:,1)
89else
90    do ig = 1,ngrid
91        do islope = 1,nslope
92            do l = 1,2
93                albedo_meshavg(ig,l) = albedo_meshavg(ig,l) + albedo_slope(ig,l,islope)*subslope_dist(ig,islope)
94            enddo
95            do iq = 1,nq
96                qsurf_meshavg(ig,iq) = qsurf_meshavg(ig,iq) + qsurf_slope(ig,iq,islope)*subslope_dist(ig,islope)/cos(def_slope_mean(islope)*pi/180.)
97            enddo
98            emis_meshavg(ig) = emis_meshavg(ig) + emis_slope(ig,islope)*subslope_dist(ig,islope)
99            tsurf_meshavg(ig) = tsurf_meshavg(ig) + (emis_slope(ig,islope)*tsurf_slope(ig,islope)**4)*subslope_dist(ig,islope)
100        enddo
101        tsurf_meshavg(ig) = (tsurf_meshavg(ig)/emis_meshavg(ig))**(0.25)
102    enddo
103endif
104
105END SUBROUTINE compute_meshgridavg
106
107END MODULE comslope_mod
Note: See TracBrowser for help on using the repository browser.