source: trunk/LMDZ.COMMON/libf/evolution/recomp_tend_mod.F90 @ 3620

Last change on this file since 3620 was 3620, checked in by jbclement, 29 hours ago

PEM:
Few small adjustments in the code.
JBC

File size: 4.6 KB
Line 
1MODULE recomp_tend_mod
2
3implicit none
4
5!=======================================================================
6contains
7!=======================================================================
8
9SUBROUTINE recomp_tend_co2(ngrid,nslope,timelen,d_co2ice_phys,d_co2ice_ini,co2ice,emissivity, &
10                           vmr_co2_PCM,vmr_co2_PEM,ps_PCM,ps_avg_global_ini,ps_avg_global)
11
12use constants_marspem_mod, only : alpha_clap_co2, beta_clap_co2, sigmaB, Lco2, sols_per_my, sec_per_sol
13
14implicit none
15
16!=======================================================================
17!
18!  To compute the evolution of the tendency for co2 ice
19!
20!=======================================================================
21
22! Inputs
23! ------
24integer,                        intent(in) :: timelen, ngrid, nslope
25real, dimension(ngrid,timelen), intent(in) :: vmr_co2_PCM       ! physical point field: Volume mixing ratio of co2 in the first layer
26real, dimension(ngrid,timelen), intent(in) :: vmr_co2_PEM       ! physical point field: Volume mixing ratio of co2 in the first layer
27real, dimension(ngrid,timelen), intent(in) :: ps_PCM            ! physical point field: Surface pressure in the PCM
28real,                           intent(in) :: ps_avg_global_ini ! global averaged pressure at previous timestep
29real,                           intent(in) :: ps_avg_global     ! global averaged pressure at current timestep
30real, dimension(ngrid,nslope),  intent(in) :: d_co2ice_ini      ! physical point field: Evolution of perennial ice over one year
31real, dimension(ngrid,nslope),  intent(in) :: co2ice            ! CO2 ice per mesh and sub-grid slope (kg/m^2)
32real, dimension(ngrid,nslope),  intent(in) :: emissivity        ! Emissivity per mesh and sub-grid slope(1)
33! Outputs
34! -------
35real, dimension(ngrid,nslope), intent(inout) :: d_co2ice_phys ! physical point field: Evolution of perennial ice over one year
36
37! Local:
38! ------
39integer :: i, t, islope
40real    :: coef, avg
41
42write(*,*) "> Updating the CO2 ice tendency for the new pressure"
43
44! Evolution of the water ice for each physical point
45do i = 1,ngrid
46    do islope = 1,nslope
47        coef = sols_per_my*sec_per_sol*emissivity(i,islope)*sigmaB/Lco2
48        avg = 0.
49        if (co2ice(i,islope) > 1.e-4 .and. abs(d_co2ice_phys(i,islope)) > 1.e-5) then
50            do t = 1,timelen
51                avg = avg + (beta_clap_co2/(alpha_clap_co2-log(vmr_co2_PCM(i,t)*ps_PCM(i,t)/100.)))**4 &
52                      - (beta_clap_co2/(alpha_clap_co2-log(vmr_co2_PEM(i,t)*ps_PCM(i,t)*(ps_avg_global/ps_avg_global_ini)/100.)))**4
53            enddo
54            if (avg < 1.e-4) avg = 0.
55            d_co2ice_phys(i,islope) = d_co2ice_ini(i,islope) - coef*avg/timelen
56        endif
57    enddo
58enddo
59
60END SUBROUTINE recomp_tend_co2
61!=======================================================================
62
63SUBROUTINE recomp_tend_h2o(icetable_depth_old,icetable_depth_new,tsurf,tsoil_PEM_timeseries_old,tsoil_PEM_timeseries_new,d_h2oice)
64
65use compute_soiltemp_mod, only: itp_tsoil
66use fast_subs_mars,       only: psv
67
68implicit none
69
70!=======================================================================
71!
72!  To compute the evolution of the tendency for h2o ice
73!
74!=======================================================================
75
76! Inputs
77! ------
78real,                 intent(in) :: icetable_depth_old, icetable_depth_new, tsurf
79real, dimension(:,:), intent(in) :: tsoil_PEM_timeseries_old, tsoil_PEM_timeseries_new
80! Outputs
81! -------
82real, intent(inout) :: d_h2oice ! Evolution of perennial ice over one year
83
84! Local:
85! ------
86real            :: Rz_old, Rz_new, R_dec, hum_dec, psv_max_old, psv_max_new
87integer         :: t
88real, parameter :: coef_diff = 4.e-4 ! Diffusion coefficient
89real, parameter :: zcdv = 0.0325     ! Drag coefficient
90
91! Higher resistance due to growing lag layer (higher depth)
92Rz_old = icetable_depth_old*zcdv/coef_diff ! Old resistance from PCM
93Rz_new = icetable_depth_new*zcdv/coef_diff ! New resistance based on new depth
94R_dec = Rz_new/Rz_old ! Decrease because of resistance
95
96! The maxmimum of the daily averages over one year for the saturation vapor pressure at the ice table location
97psv_max_old = 0.
98psv_max_new = 0.
99do t = 1,size(tsoil_PEM_timeseries_old,2)
100    psv_max_old = max(psv_max_old,psv(itp_tsoil(tsoil_PEM_timeseries_old(:,t),tsurf,icetable_depth_old)))
101    psv_max_new = max(psv_max_new,psv(itp_tsoil(tsoil_PEM_timeseries_new(:,t),tsurf,icetable_depth_new)))
102enddo
103
104! Lower humidity due to growing lag layer (higher depth)
105hum_dec = psv_max_old/psv_max_new ! Decrease because of lower water vapor pressure at the new depth
106
107! Flux correction (decrease)
108d_h2oice = d_h2oice*R_dec*hum_dec
109
110END SUBROUTINE recomp_tend_h2o
111
112END MODULE recomp_tend_mod
Note: See TracBrowser for help on using the repository browser.