Changeset 3498 for trunk/LMDZ.COMMON/libf/evolution/evol_ice_mod.F90
- Timestamp:
- Nov 7, 2024, 2:48:08 PM (2 weeks ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.COMMON/libf/evolution/evol_ice_mod.F90
r3368 r3498 7 7 !======================================================================= 8 8 9 SUBROUTINE evol_co2_ice(ngrid,nslope,co2_ice, tend_co2_ice)9 SUBROUTINE evol_co2_ice(ngrid,nslope,co2_ice,d_co2ice) 10 10 11 use time_evol_mod, only: dt _pem11 use time_evol_mod, only: dt 12 12 13 13 implicit none … … 25 25 ! OUTPUT 26 26 real, 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 year27 real, dimension(ngrid,nslope), intent(inout) :: d_co2ice ! Evolution of perennial ice over one year 28 28 29 29 ! local: … … 35 35 36 36 co2_ice_old = co2_ice 37 co2_ice = co2_ice + tend_co2_ice*dt_pem37 co2_ice = co2_ice + d_co2ice*dt 38 38 where (co2_ice < 0.) 39 39 co2_ice = 0. 40 tend_co2_ice = -co2_ice_old/dt_pem40 d_co2ice = -co2_ice_old/dt 41 41 end where 42 42 … … 45 45 !======================================================================= 46 46 47 SUBROUTINE evol_h2o_ice(ngrid,nslope,cell_area,delta_h2o_adsorbded,delta_h2o_icetablesublim,h2o_ice, tend_h2o_ice,stopPEM)47 SUBROUTINE evol_h2o_ice(ngrid,nslope,cell_area,delta_h2o_adsorbded,delta_h2o_icetablesublim,h2o_ice,d_h2oice,stopPEM) 48 48 49 use time_evol_mod, only: dt _pem49 use time_evol_mod, only: dt 50 50 use comslope_mod, only: subslope_dist, def_slope_mean 51 51 … … 73 73 ! OUTPUT 74 74 real, 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)75 real, dimension(ngrid,nslope), intent(inout) :: d_h2oice ! Evolution of perennial ice over one year (kg/m^2/year) 76 76 integer, intent(inout) :: stopPEM ! Stopping criterion code 77 77 … … 80 80 integer :: i, islope ! Loop variables 81 81 real :: pos_tend, neg_tend, real_coefficient, negative_part ! Variable to conserve h2o 82 real, dimension(ngrid,nslope) :: new_tend encies! Tendencies computed in order to conserve h2o ice on the surface, only exchange between surface are done82 real, dimension(ngrid,nslope) :: new_tend ! Tendencies computed in order to conserve h2o ice on the surface, only exchange between surface are done 83 83 !======================================================================= 84 84 if (ngrid /= 1) then ! Not in 1D … … 99 99 do islope = 1,nslope 100 100 if (h2o_ice(i,islope) > 0.) then 101 if ( tend_h2o_ice(i,islope) > 0.) then102 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.) 103 103 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.) 105 105 endif 106 106 endif … … 108 108 enddo 109 109 110 new_tend encies= 0.110 new_tend = 0. 111 111 if (abs(pos_tend) < 1.e-10 .or. abs(neg_tend) < 1.e-10) then 112 112 write(*,*) "Reason of stopping: there is no sublimating or condensing h2o ice!" … … 118 118 ! We adapt the tendencies to conserve h2o and do only exchange between grid points 119 119 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 coefficient121 new_tend encies = tend_h2o_ice*pos_tend/neg_tend120 where (d_h2oice < 0.) ! We lower the sublimating rate by a coefficient 121 new_tend = d_h2oice*pos_tend/neg_tend 122 122 elsewhere ! We don't change the condensing rate 123 new_tend encies = tend_h2o_ice123 new_tend = d_h2oice 124 124 end where 125 125 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 rate127 new_tend encies = tend_h2o_ice126 where (d_h2oice < 0.) ! We don't change the sublimating rate 127 new_tend = d_h2oice 128 128 elsewhere ! We lower the condensing rate by a coefficient 129 new_tend encies = tend_h2o_ice*neg_tend/pos_tend129 new_tend = d_h2oice*neg_tend/pos_tend 130 130 end where 131 131 endif … … 133 133 134 134 ! Evolution of the h2o ice for each physical point 135 h2o_ice = h2o_ice + new_tend encies*dt_pem135 h2o_ice = h2o_ice + new_tend*dt 136 136 137 137 ! We compute the amount of h2o that is sublimated in excess … … 142 142 negative_part = negative_part - h2o_ice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.) 143 143 h2o_ice(i,islope) = 0. 144 tend_h2o_ice(i,islope) = 0.144 d_h2oice(i,islope) = 0. 145 145 endif 146 146 enddo … … 157 157 do islope = 1,nslope 158 158 do i = 1,ngrid 159 if (new_tend encies(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.) 160 160 enddo 161 161 enddo 162 162 else ! ngrid == 1, i.e. in 1D 163 h2o_ice = h2o_ice + tend_h2o_ice*dt_pem163 h2o_ice = h2o_ice + d_h2oice*dt 164 164 endif 165 165
Note: See TracChangeset
for help on using the changeset viewer.