1 | MODULE 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 | |
---|
13 | implicit none |
---|
14 | |
---|
15 | integer, save :: nslope ! Number of slopes for the statistic (1) |
---|
16 | integer, save :: iflat ! Flat slope for islope (1) |
---|
17 | real, save, allocatable, dimension(:) :: def_slope ! Bound of the slope statistics / repartitions (°) |
---|
18 | real, save, allocatable, dimension(:) :: def_slope_mean ! Mean slope of each bin (°) |
---|
19 | real, save, allocatable, dimension(:) :: sky_slope ! Portion of the sky view by each slopes (1) |
---|
20 | real, save, allocatable, dimension(:,:) :: subslope_dist ! Distribution of the slopes (1) |
---|
21 | integer, 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 | !======================================================================= |
---|
25 | contains |
---|
26 | !======================================================================= |
---|
27 | |
---|
28 | SUBROUTINE ini_comslope_h(ngrid,nslope_in) |
---|
29 | |
---|
30 | implicit none |
---|
31 | |
---|
32 | integer, intent(in) :: ngrid |
---|
33 | integer, intent(in) :: nslope_in |
---|
34 | |
---|
35 | allocate(def_slope(nslope_in+1)) |
---|
36 | allocate(def_slope_mean(nslope_in)) |
---|
37 | allocate(sky_slope(nslope_in)) |
---|
38 | allocate(subslope_dist(ngrid,nslope_in)) |
---|
39 | allocate(major_slope(ngrid)) |
---|
40 | |
---|
41 | END SUBROUTINE ini_comslope_h |
---|
42 | |
---|
43 | !======================================================================= |
---|
44 | |
---|
45 | SUBROUTINE end_comslope_h |
---|
46 | |
---|
47 | implicit none |
---|
48 | |
---|
49 | if (allocated(def_slope)) deallocate(def_slope) |
---|
50 | if (allocated(def_slope_mean)) deallocate(def_slope_mean) |
---|
51 | if (allocated(sky_slope)) deallocate(sky_slope) |
---|
52 | if (allocated(subslope_dist)) deallocate(subslope_dist) |
---|
53 | if (allocated(major_slope)) deallocate(major_slope) |
---|
54 | |
---|
55 | END SUBROUTINE end_comslope_h |
---|
56 | |
---|
57 | !======================================================================= |
---|
58 | |
---|
59 | SUBROUTINE compute_meshgridavg(ngrid,nq,albedo_slope,emis_slope,tsurf_slope,qsurf_slope,albedo_meshavg,emis_meshavg,tsurf_meshavg, qsurf_meshavg) |
---|
60 | |
---|
61 | use comcstfi_h, only: pi |
---|
62 | |
---|
63 | implicit none |
---|
64 | |
---|
65 | integer, intent(in) :: ngrid, nq ! # points on the physical grid, tracers (1) |
---|
66 | real, dimension(ngrid,2,nslope), intent(in) :: albedo_slope ! albedo on each sub-grid slope (1) |
---|
67 | real, dimension(ngrid,nslope), intent(in) :: emis_slope ! emissivity on each sub-grid slope (1) |
---|
68 | real, dimension(ngrid,nslope), intent(in) :: tsurf_slope ! Surface Temperature on each sub-grid slope (K) |
---|
69 | real, dimension(ngrid,nq,nslope), intent(in) :: qsurf_slope ! Surface Tracers on each sub-grid slope (kg/m^2 sloped) |
---|
70 | real, dimension(ngrid,2), intent(out) :: albedo_meshavg ! grid box average of the albedo (1) |
---|
71 | real, dimension(ngrid), intent(out) :: emis_meshavg ! grid box average of the emissivity (1) |
---|
72 | real, dimension(ngrid), intent(out) :: tsurf_meshavg ! grid box average of the surface temperature (K) |
---|
73 | real, dimension(ngrid,nq), intent(out) :: qsurf_meshavg ! grid box average of the surface tracers (kg/m^2 flat) |
---|
74 | ! Local |
---|
75 | integer :: islope, ig, l, iq |
---|
76 | |
---|
77 | !------------------- |
---|
78 | |
---|
79 | albedo_meshavg = 0. |
---|
80 | emis_meshavg = 0. |
---|
81 | tsurf_meshavg = 0. |
---|
82 | qsurf_meshavg = 0. |
---|
83 | |
---|
84 | if (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) |
---|
89 | else |
---|
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 |
---|
103 | endif |
---|
104 | |
---|
105 | END SUBROUTINE compute_meshgridavg |
---|
106 | |
---|
107 | END MODULE comslope_mod |
---|