source: trunk/LMDZ.COMMON/libf/evolution/recomp_tend_co2_mod.F90 @ 3584

Last change on this file since 3584 was 3584, checked in by jbclement, 4 days ago

PEM:

  • Albedo is now updated only at the end of the PEM (and not at every iteration) + Correct way to set it taking into account CO2/H2O ice and frost.
  • Cosmetic cleanings.

JBC

File size: 4.3 KB
Line 
1MODULE recomp_tend_co2_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 tendencie 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(*,*) "Update of the CO2 tendency from the current 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(ngrid,nslope,timelen,d_h2oice,PCM_temp,PEM_temp)
64
65implicit none
66
67!=======================================================================
68!
69!  To compute the evolution of the tendencie for h2o ice
70!
71!=======================================================================
72
73! Inputs
74! ------
75integer,            intent(in) :: timelen, ngrid, nslope
76real, dimension(:), intent(in) :: PCM_temp, PEM_temp
77! Outputs
78! -------
79real, dimension(ngrid,nslope), intent(inout) :: d_h2oice ! physical point field: Evolution of perennial ice over one year
80
81! Local:
82! ------
83
84write(*,*) "Update of the H2O tendency due to lag layer"
85
86! Flux correction due to lag layer
87!~ Rz_old = h2oice_depth_old*0.0325/4.e-4              ! resistance from PCM
88!~ Rz_new = h2oice_depth_new*0.0325/4.e-4              ! new resistance based on new depth
89!~ R_dec = (1./Rz_old)/(1./Rz_new)                     ! decrease because of resistance
90!~ soil_psv_old = psv(max(PCM_temp(h2oice_depth_old))) ! the maxmimum annual mean saturation vapor pressure at the temperature of the GCM run temperature at the old ice location
91!~ soil_psv_new = psv(max(PEM_temp(h2oice_depth_new))) ! the maxmimum annual mean saturation vapor pressure at the temperature of the PEM run temperature at the new ice location
92!~ hum_dec = soil_psv_old/soil_psv_new                 ! decrease because of lower water vapor pressure at the new depth
93!~ d_h2oice = d_h2oice*R_dec*hum_dec                   ! decrease of flux
94
95END SUBROUTINE recomp_tend_h2o
96
97END MODULE recomp_tend_co2_mod
Note: See TracBrowser for help on using the repository browser.