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

Last change on this file since 2909 was 2909, checked in by romain.vande, 21 months ago

Mars PEM:
New Boolean options for following orbital parameters of ob_ex_lsp.asc: var_obl, var_ex, var_lsp.
If using evol_orbit_pem=true, you can specify which parameter to follow.
True by default: Do you want to change the parameter XXX during the PEM run as prescribed in ob_ex_lsp.asc.
If false, it is set to constant (to the value of the tab_cntrl in the start)
RV

File size: 4.3 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                                      ! Number of slopes for the statistic  (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
24CONTAINS
25  SUBROUTINE ini_comslope_h(ngrid,nslope_in)
26
27  IMPLICIT NONE
28  INTEGER, INTENT(IN) :: ngrid
29  INTEGER, INTENT(IN) :: nslope_in
30
31    allocate(def_slope(nslope_in+1))
32    allocate(def_slope_mean(nslope_in))
33    allocate(sky_slope(nslope_in))
34    allocate(subslope_dist(ngrid,nslope_in))
35    allocate(major_slope(ngrid))
36  END SUBROUTINE ini_comslope_h
37
38
39  SUBROUTINE end_comslope_h
40
41  IMPLICIT NONE
42
43        IF (allocated(def_slope)) deallocate(def_slope)
44        IF (allocated(def_slope_mean)) deallocate(def_slope_mean)
45        IF (allocated(sky_slope)) deallocate(sky_slope)
46        IF (allocated(subslope_dist)) deallocate(subslope_dist)
47        IF (allocated(major_slope)) deallocate(major_slope)
48
49  END SUBROUTINE end_comslope_h
50
51  SUBROUTINE compute_meshgridavg(ngrid,nq,albedo_slope,emis_slope,tsurf_slope,qsurf_slope,albedo_meshavg,emis_meshavg,tsurf_meshavg, qsurf_meshavg)
52  USE comcstfi_h, only:  pi
53  IMPLICIT NONE
54  INTEGER, INTENT(IN) :: ngrid,nq                  !  # points on the physical grid, tracers  (1)
55  REAL, INTENT(IN) :: albedo_slope(ngrid,2,nslope) !  albedo on each sub-grid slope (1)
56  REAL, INTENT(IN) :: emis_slope(ngrid,nslope)     !  emissivity on each sub-grid slope (1)
57  REAL, INTENT(IN) :: tsurf_slope(ngrid,nslope)    !  Surface Temperature on each sub-grid slope (K)
58  REAL, INTENT(IN) :: qsurf_slope(ngrid,nq,nslope) !  Surface Tracers on each sub-grid slope (kg/m^2 sloped)
59  REAL, INTENT(OUT) :: albedo_meshavg(ngrid,2)     !  grid box average of the albedo (1)
60  REAL, INTENT(OUT) :: emis_meshavg(ngrid)         !  grid box average of the emissivity (1)
61  REAL, INTENT(OUT) :: tsurf_meshavg(ngrid)        !  grid box average of the surface temperature (K)
62  REAL, INTENT(OUT) :: qsurf_meshavg(ngrid,nq)     !  grid box average of the surface tracers (kg/m^2 flat)
63! Local
64  integer :: islope,ig,l,iq
65
66!-------------------
67
68      albedo_meshavg(:,:) = 0.
69      emis_meshavg(:) = 0.
70      tsurf_meshavg(:) = 0.
71      qsurf_meshavg(:,:) = 0.
72
73  if(nslope.eq.1) then
74      albedo_meshavg(:,:) = albedo_slope(:,:,1)
75      emis_meshavg(:) = emis_slope(:,1)
76      tsurf_meshavg(:) = tsurf_slope(:,1)
77      qsurf_meshavg(:,:) = qsurf_slope(:,:,1)
78  else
79
80  DO ig = 1,ngrid
81      DO islope = 1,nslope
82          DO l = 1,2
83            albedo_meshavg(ig,l) = albedo_meshavg(ig,l) + albedo_slope(ig,l,islope)*subslope_dist(ig,islope)
84          ENDDO
85          DO iq = 1,nq
86             qsurf_meshavg(ig,iq) = qsurf_meshavg(ig,iq) + qsurf_slope(ig,iq,islope)*subslope_dist(ig,islope)/cos(def_slope_mean(islope)*pi/180.)
87          ENDDO
88          emis_meshavg(ig) = emis_meshavg(ig) + emis_slope(ig,islope)*subslope_dist(ig,islope)
89          tsurf_meshavg(ig) = tsurf_meshavg(ig) + (emis_slope(ig,islope)*tsurf_slope(ig,islope)**4)*subslope_dist(ig,islope)
90      ENDDO
91      tsurf_meshavg(ig) = (tsurf_meshavg(ig)/emis_meshavg(ig))**(0.25)
92  ENDDO
93
94  endif
95
96  END SUBROUTINE compute_meshgridavg
97
98END MODULE comslope_mod
Note: See TracBrowser for help on using the repository browser.