Ignore:
Timestamp:
Nov 7, 2024, 2:48:08 PM (2 weeks ago)
Author:
jbclement
Message:

PEM:

  • Correction of the variable name for the ice table depth in "pemetat0.F90". So it is now got as intended from the "startpem.nc" file;
  • Renaming of the tendencies in the PEM with the prefix 'd_' instead of 'tend_';
  • Modification of the PEM time step type from integer to real. As a consequence, all time variables are now of real type. This change adds the possibility to consider fractions of year as time step.

JBC

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.COMMON/libf/evolution/evol_ice_mod.F90

    r3368 r3498  
    77!=======================================================================
    88
    9 SUBROUTINE evol_co2_ice(ngrid,nslope,co2_ice,tend_co2_ice)
     9SUBROUTINE evol_co2_ice(ngrid,nslope,co2_ice,d_co2ice)
    1010
    11 use time_evol_mod, only: dt_pem
     11use time_evol_mod, only: dt
    1212
    1313implicit none
     
    2525!   OUTPUT
    2626real, dimension(ngrid,nslope), intent(inout) :: co2_ice      ! Previous and actual density of CO2 ice
    27 real, dimension(ngrid,nslope), intent(inout) :: tend_co2_ice ! Evolution of perennial ice over one year
     27real, dimension(ngrid,nslope), intent(inout) :: d_co2ice ! Evolution of perennial ice over one year
    2828
    2929!   local:
     
    3535
    3636co2_ice_old = co2_ice
    37 co2_ice = co2_ice + tend_co2_ice*dt_pem
     37co2_ice = co2_ice + d_co2ice*dt
    3838where (co2_ice < 0.)
    3939    co2_ice = 0.
    40     tend_co2_ice = -co2_ice_old/dt_pem
     40    d_co2ice = -co2_ice_old/dt
    4141end where
    4242
     
    4545!=======================================================================
    4646
    47 SUBROUTINE evol_h2o_ice(ngrid,nslope,cell_area,delta_h2o_adsorbded,delta_h2o_icetablesublim,h2o_ice,tend_h2o_ice,stopPEM)
     47SUBROUTINE evol_h2o_ice(ngrid,nslope,cell_area,delta_h2o_adsorbded,delta_h2o_icetablesublim,h2o_ice,d_h2oice,stopPEM)
    4848
    49 use time_evol_mod,          only: dt_pem
     49use time_evol_mod,          only: dt
    5050use comslope_mod,           only: subslope_dist, def_slope_mean
    5151
     
    7373!   OUTPUT
    7474real, dimension(ngrid,nslope), intent(inout) :: h2o_ice      ! Previous and actual density of h2o ice (kg/m^2)
    75 real, dimension(ngrid,nslope), intent(inout) :: tend_h2o_ice ! Evolution of perennial ice over one year (kg/m^2/year)
     75real, dimension(ngrid,nslope), intent(inout) :: d_h2oice ! Evolution of perennial ice over one year (kg/m^2/year)
    7676integer,                       intent(inout) :: stopPEM      ! Stopping criterion code
    7777
     
    8080integer                       :: i, islope                                           ! Loop variables
    8181real                          :: pos_tend, neg_tend, real_coefficient, negative_part ! Variable to conserve h2o
    82 real, dimension(ngrid,nslope) :: new_tendencies                                      ! Tendencies computed in order to conserve h2o ice on the surface, only exchange between surface are done
     82real, dimension(ngrid,nslope) :: new_tend                                      ! Tendencies computed in order to conserve h2o ice on the surface, only exchange between surface are done
    8383!=======================================================================
    8484if (ngrid /= 1) then ! Not in 1D
     
    9999        do islope = 1,nslope
    100100            if (h2o_ice(i,islope) > 0.) then
    101                 if (tend_h2o_ice(i,islope) > 0.) then
    102                     pos_tend = pos_tend + tend_h2o_ice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.)
     101                if (d_h2oice(i,islope) > 0.) then
     102                    pos_tend = pos_tend + d_h2oice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.)
    103103                else
    104                     neg_tend = neg_tend - tend_h2o_ice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.)
     104                    neg_tend = neg_tend - d_h2oice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.)
    105105                endif
    106106            endif
     
    108108    enddo
    109109
    110     new_tendencies = 0.
     110    new_tend = 0.
    111111    if (abs(pos_tend) < 1.e-10 .or. abs(neg_tend) < 1.e-10) then
    112112        write(*,*) "Reason of stopping: there is no sublimating or condensing h2o ice!"
     
    118118        ! We adapt the tendencies to conserve h2o and do only exchange between grid points
    119119        if (neg_tend > pos_tend .and. pos_tend > 0.) then ! More sublimation on the planet than condensation
    120             where (tend_h2o_ice < 0.) ! We lower the sublimating rate by a coefficient
    121                 new_tendencies = tend_h2o_ice*pos_tend/neg_tend
     120            where (d_h2oice < 0.) ! We lower the sublimating rate by a coefficient
     121                new_tend = d_h2oice*pos_tend/neg_tend
    122122            elsewhere ! We don't change the condensing rate
    123                 new_tendencies = tend_h2o_ice
     123                new_tend = d_h2oice
    124124            end where
    125125        else if (neg_tend < pos_tend .and. neg_tend > 0.) then ! More condensation on the planet than sublimation
    126             where (tend_h2o_ice < 0.) ! We don't change the sublimating rate
    127                 new_tendencies = tend_h2o_ice
     126            where (d_h2oice < 0.) ! We don't change the sublimating rate
     127                new_tend = d_h2oice
    128128            elsewhere ! We lower the condensing rate by a coefficient
    129                 new_tendencies = tend_h2o_ice*neg_tend/pos_tend
     129                new_tend = d_h2oice*neg_tend/pos_tend
    130130            end where
    131131        endif
     
    133133
    134134    ! Evolution of the h2o ice for each physical point
    135     h2o_ice = h2o_ice + new_tendencies*dt_pem
     135    h2o_ice = h2o_ice + new_tend*dt
    136136
    137137    ! We compute the amount of h2o that is sublimated in excess
     
    142142                negative_part = negative_part - h2o_ice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.)
    143143                h2o_ice(i,islope) = 0.
    144                 tend_h2o_ice(i,islope) = 0.
     144                d_h2oice(i,islope) = 0.
    145145            endif
    146146        enddo
     
    157157    do islope = 1,nslope
    158158        do i = 1,ngrid
    159              if (new_tendencies(i,islope) > 0.) h2o_ice(i,islope) = h2o_ice(i,islope) - new_tendencies(i,islope)*real_coefficient*dt_pem*cos(def_slope_mean(islope)*pi/180.)
     159             if (new_tend(i,islope) > 0.) h2o_ice(i,islope) = h2o_ice(i,islope) - new_tend(i,islope)*real_coefficient*dt*cos(def_slope_mean(islope)*pi/180.)
    160160        enddo
    161161     enddo
    162162else ! ngrid == 1, i.e. in 1D
    163     h2o_ice = h2o_ice + tend_h2o_ice*dt_pem
     163    h2o_ice = h2o_ice + d_h2oice*dt
    164164endif
    165165
Note: See TracChangeset for help on using the changeset viewer.