source: trunk/LMDZ.COMMON/libf/evolution/surf_ice.F90 @ 4110

Last change on this file since 4110 was 4110, checked in by jbclement, 6 days ago

PEM:

  • Introduction of a configurable display/logging system with options 'out2term', 'out2log', 'verbosity_lvl'. All messages now use verbosity levels ('LVL_NFO', 'LVL_WRN', 'LVL_ERR' and 'LVL_DBG').
  • Code encapsulation improvements with systematic privacy/protection of module variables.
  • Addition of workflow safety checks for required executables, dependencies (e.g. 'ncdump'), input files and callphys keys.
  • Renaming of PEM starting and diagnostic files ("startevol.nc" into "startevo.nc", "diagevol.nc" into "diagevo.nc").

JBC

File size: 13.5 KB
RevLine 
[3989]1MODULE surf_ice
[3991]2!-----------------------------------------------------------------------
3! NAME
4!     surf_ice
5!
6! DESCRIPTION
7!     Surface ice management.
8!
9! AUTHORS & DATE
10!     JB Clement, 12/2025
11!
12! NOTES
13!
14!-----------------------------------------------------------------------
[3149]15
[4065]16! DEPENDENCIES
17! ------------
18use numerics, only: dp, qp, di, k4
19
[3991]20! DECLARATION
21! -----------
[3149]22implicit none
23
[4065]24! PARAMETERS
25! ----------
26real(dp),                                 parameter :: rho_co2ice = 1650._dp ! Density of CO2 ice [kg.m-3]
27real(dp),                                 parameter :: rho_h2oice = 920._dp  ! Density of H2O ice [kg.m-3]
28real(dp),                                 protected :: threshold_h2oice_cap  ! Threshold to consider H2O ice as infinite reservoir [kg.m-2]
29real(dp),                                 protected :: h2oice_huge_ini       ! Initial value for huge reservoirs of H2O ice [kg/m^2]
30logical(k4), dimension(:),   allocatable, protected :: is_h2o_perice_PCM     ! Flag for the presence of perennial H2O ice in the PCM at the beginning [kg.m-2]
31real(dp),    dimension(:,:), allocatable, protected :: co2_perice_PCM        ! Perennial CO2 ice in the PCM at the beginning [kg.m-2]
[3989]32
[3149]33contains
[3991]34!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
[3149]35
[3989]36!=======================================================================
[4065]37SUBROUTINE ini_surf_ice()
[3991]38!-----------------------------------------------------------------------
39! NAME
[4065]40!     ini_surf_ice
[3991]41!
42! DESCRIPTION
[4065]43!     Initialize surface ice arrays.
[3991]44!
45! AUTHORS & DATE
[4065]46!     JB Clement 12/2025
47!
48! NOTES
49!
50!-----------------------------------------------------------------------
51
52! DEPENDENCIES
53! ------------
54use geometry, only: ngrid, nslope
55
56! DECLARATION
57! -----------
58implicit none
59
60! CODE
61! ----
62if (.not. allocated(is_h2o_perice_PCM)) allocate(is_h2o_perice_PCM(ngrid))
63if (.not. allocated(co2_perice_PCM)) allocate(co2_perice_PCM(ngrid,nslope))
[4074]64is_h2o_perice_PCM(:) = .false.
65co2_perice_PCM(:,:) = 0._dp
[4065]66
67END SUBROUTINE ini_surf_ice
68!=======================================================================
69
70!=======================================================================
71SUBROUTINE end_surf_ice()
72!-----------------------------------------------------------------------
73! NAME
74!     end_surf_ice
75!
76! DESCRIPTION
77!     Deallocate surface ice arrays.
78!
79! AUTHORS & DATE
[3991]80!     JB Clement, 12/2025
81!
82! NOTES
83!
84!-----------------------------------------------------------------------
[3149]85
[4065]86! DECLARATION
87! -----------
88implicit none
89
90! CODE
91! ----
92if (allocated(is_h2o_perice_PCM)) deallocate(is_h2o_perice_PCM)
93if (allocated(co2_perice_PCM)) deallocate(co2_perice_PCM)
94
95END SUBROUTINE end_surf_ice
96!=======================================================================
97
98!=======================================================================
99SUBROUTINE set_surf_ice_config(threshold_h2oice_cap_in,h2oice_huge_ini_in)
100!-----------------------------------------------------------------------
101! NAME
102!     set_surf_ice_config
103!
104! DESCRIPTION
105!     Setter for 'surf_ice' configuration parameters.
106!
107! AUTHORS & DATE
108!     JB Clement, 02/2026
109!
110! NOTES
111!
112!-----------------------------------------------------------------------
113
[3991]114! DEPENDENCIES
115! ------------
[4110]116use utility,  only: real2str
117use display,  only: print_msg, LVL_NFO
118use stoppage, only: stop_clean
[3149]119
[3991]120! DECLARATION
121! -----------
[3149]122implicit none
123
[3991]124! ARGUMENTS
125! ---------
[4065]126real(dp), intent(in) :: threshold_h2oice_cap_in, h2oice_huge_ini_in
[3989]127
[4065]128! CODE
129! ----
130threshold_h2oice_cap = threshold_h2oice_cap_in
131h2oice_huge_ini = h2oice_huge_ini_in
[4110]132call print_msg('threshold_h2oice_cap = '//real2str(threshold_h2oice_cap),LVL_NFO)
133call print_msg('h2oice_huge_ini      = '//real2str(h2oice_huge_ini),LVL_NFO)
134if (threshold_h2oice_cap < 0._dp) call stop_clean(__FILE__,__LINE__,'''threshold_h2oice_cap'' must be positive or null!',1)
135if (h2oice_huge_ini < 0._dp) call stop_clean(__FILE__,__LINE__,'''h2oice_huge_ini'' must be positive or null!',1)
[3989]136
[4065]137END SUBROUTINE set_surf_ice_config
138!=======================================================================
139
140!=======================================================================
141SUBROUTINE set_co2_perice_PCM(co2_perice_PCM_in)
142!-----------------------------------------------------------------------
143! NAME
144!     set_co2_perice_PCM
145!
146! DESCRIPTION
147!     Setter for 'co2_perice_PCM'.
148!
149! AUTHORS & DATE
150!     JB Clement, 12/2025
151!
152! NOTES
153!
154!-----------------------------------------------------------------------
155
156! DECLARATION
157! -----------
158implicit none
159
160! ARGUMENTS
161! ---------
162real(dp), dimension(:,:), intent(in) :: co2_perice_PCM_in
163
[3991]164! CODE
165! ----
[4065]166co2_perice_PCM(:,:) = co2_perice_PCM_in(:,:)
[3989]167
[4065]168END SUBROUTINE set_co2_perice_PCM
[3149]169!=======================================================================
[3989]170
171!=======================================================================
[4065]172SUBROUTINE set_is_h2o_perice_PCM(is_h2o_perice_PCM_in)
[3991]173!-----------------------------------------------------------------------
174! NAME
[4065]175!     set_is_h2o_perice_PCM
[3991]176!
177! DESCRIPTION
[4065]178!     Setter for 'is_h2o_perice_PCM'.
[3991]179!
180! AUTHORS & DATE
[4065]181!     JB Clement, 12/2025
[3991]182!
183! NOTES
184!
185!-----------------------------------------------------------------------
[3989]186
[3991]187! DECLARATION
188! -----------
[3989]189implicit none
190
[3991]191! ARGUMENTS
192! ---------
[4065]193real(dp), dimension(:), intent(in) :: is_h2o_perice_PCM_in
[3989]194
[3991]195! CODE
196! ----
[4065]197where(is_h2o_perice_PCM_in > 0.5_dp)
198    is_h2o_perice_PCM = .true.
199else where
200    is_h2o_perice_PCM = .false.
201end where
[3989]202
[4065]203END SUBROUTINE set_is_h2o_perice_PCM
[3989]204!=======================================================================
205
206!=======================================================================
[4071]207SUBROUTINE build4PCM_perice(h2o_ice,co2_ice,is_h2o_perice,h2o_ice4PCM,co2_ice4PCM)
[3991]208!-----------------------------------------------------------------------
209! NAME
[4065]210!     build4PCM_perice
[3991]211!
212! DESCRIPTION
[4065]213!     Reconstructs perennial ice and frost for the PCM.
[3991]214!
215! AUTHORS & DATE
216!     JB Clement, 12/2025
217!
218! NOTES
219!
220!-----------------------------------------------------------------------
[3989]221
[4065]222! DEPENDENCIES
223! ------------
224use geometry, only: ngrid
225use frost,    only: h2o_frost4PCM
226use slopes,   only: subslope_dist, def_slope_mean
227use maths,    only: pi
[4110]228use display,  only: print_msg, LVL_NFO
[4065]229
[3991]230! DECLARATION
231! -----------
[3989]232implicit none
233
[4065]234! ARGUMENTS
235! ---------
236real(dp),    dimension(:,:), intent(inout) :: h2o_ice, co2_ice
[4071]237logical(k4), dimension(:),   intent(out)   :: is_h2o_perice            ! H2O perennial ice flag
238real(dp),    dimension(:,:), intent(out)   :: h2o_ice4PCM, co2_ice4PCM ! Ice for PCM
[4065]239
240! LOCAL VARIABLES
241! ---------------
242integer(di) :: i
243
[3991]244! CODE
245! ----
[4110]246call print_msg('> Building surface ice/frost for the PCM',LVL_NFO)
[4071]247co2_ice4PCM(:,:) = co2_ice(:,:)
248h2o_ice4PCM(:,:) = 0._dp ! Because in the Mars PCM, only the variation of perennial H2O ice is monitored, not the absolute quantity
[4065]249do i = 1,ngrid
250    ! Is H2O ice still considered as an infinite reservoir for the PCM?
251    if (sum(h2o_ice(i,:)*subslope_dist(i,:)/cos(pi*def_slope_mean(:)/180._dp)) > threshold_h2oice_cap) then
252        ! There is enough ice to be considered as an infinite reservoir
253        is_h2o_perice(i) = .true.
254    else
255        ! Too little ice to be considered as an infinite reservoir so ice is transferred to the frost
256        is_h2o_perice(i) = .false.
257        h2o_frost4PCM(i,:) = h2o_frost4PCM(i,:) + h2o_ice(i,:)
258        h2o_ice(i,:) = 0._dp
259    end if
260end do
[3989]261
[4065]262END SUBROUTINE build4PCM_perice
[3989]263!=======================================================================
264
265!=======================================================================
[4065]266SUBROUTINE evolve_co2ice(co2_ice,d_co2ice,zshift_surf)
[3991]267!-----------------------------------------------------------------------
268! NAME
[4065]269!     evolve_co2ice
[3991]270!
271! DESCRIPTION
272!     Compute the evolution of CO2 ice over one year.
273!
274! AUTHORS & DATE
275!     R. Vandemeulebrouck
276!     L. Lange
277!     JB Clement, 2023-2025
278!
279! NOTES
280!
281!-----------------------------------------------------------------------
[3989]282
[3991]283! DEPENDENCIES
284! ------------
[4065]285use geometry,  only: ngrid, nslope
[3989]286use evolution, only: dt
[4110]287use display,   only: print_msg, LVL_NFO
[3989]288
[3991]289! DECLARATION
290! -----------
[3989]291implicit none
292
[3991]293! ARGUMENTS
[3989]294! ---------
[4065]295real(dp), dimension(:,:), intent(inout) :: co2_ice
296real(dp), dimension(:,:), intent(inout) :: d_co2ice
297real(dp), dimension(:,:), intent(out)   :: zshift_surf
[3149]298
[3991]299! LOCAL VARIABLES
[3989]300! ---------------
[4065]301real(dp), dimension(ngrid,nslope) :: co2_ice_old ! Old density of CO2 ice
[3989]302
[3991]303! CODE
[3989]304! ----
[3149]305! Evolution of CO2 ice for each physical point
[4110]306call print_msg('> Evolution of CO2 ice',LVL_NFO)
[3367]307
308co2_ice_old = co2_ice
[3498]309co2_ice = co2_ice + d_co2ice*dt
[4065]310where (co2_ice < 0._dp)
311    co2_ice = 0._dp
[3498]312    d_co2ice = -co2_ice_old/dt
[3149]313end where
[3367]314
[3553]315zshift_surf = d_co2ice*dt/rho_co2ice
316
[4065]317END SUBROUTINE evolve_co2ice
[3989]318!=======================================================================
[3149]319
320!=======================================================================
[4065]321SUBROUTINE evolve_h2oice(delta_h2o_ads,delta_icetable,h2o_ice,d_h2oice,zshift_surf,stopcrit)
[3991]322!-----------------------------------------------------------------------
323! NAME
[4065]324!     evolve_h2oice
[3991]325!
326! DESCRIPTION
327!     Compute the evolution of H2O ice over one year.
328!
329! AUTHORS & DATE
330!     R. Vandemeulebrouck
331!     L. Lange
332!     JB Clement, 2023-2025
333!
334! NOTES
335!
336!-----------------------------------------------------------------------
[3149]337
[3991]338! DEPENDENCIES
339! ------------
[4065]340use geometry,      only: ngrid, nslope
[3991]341use evolution,     only: dt
[3989]342use stopping_crit, only: stopping_crit_h2o, stopFlags
[4110]343use display,       only: print_msg, LVL_NFO
[3149]344
[3991]345! DECLARATION
346! -----------
[3149]347implicit none
348
[3991]349! ARGUMENTS
[3989]350! ---------
[4065]351real(dp), dimension(:),   intent(in)    :: delta_h2o_ads, delta_icetable
352real(dp), dimension(:,:), intent(inout) :: h2o_ice, d_h2oice
353type(stopFlags),          intent(inout) :: stopcrit
354real(dp), dimension(:,:), intent(out)   :: zshift_surf
[3149]355
[3991]356! LOCAL VARIABLES
[3989]357! ---------------
[4065]358real(qp)                          :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm ! H2O balance variables
359real(dp), dimension(ngrid,nslope) :: d_h2oice_new                                             ! Tendencies computed to keep balance
[3989]360
[3991]361! CODE
[3989]362! ----
[4110]363call print_msg('> Evolution of H2O ice',LVL_NFO)
[3308]364
[3938]365if (ngrid == 1) then ! In 1D
[4074]366    where (d_h2oice(:,:) < 0._dp .and. abs(d_h2oice(:,:)*dt) > h2o_ice(:,:)) ! Not enough ice to sublimate everything
367        ! We sublimate what we can
368        d_h2oice_new(:,:) = h2o_ice(:,:)/dt
369        ! It means the tendency is zero next time
370        d_h2oice(:,:) = 0._dp
371    else where
372        d_h2oice_new(:,:) = d_h2oice(:,:)
373    end where
[3938]374else ! In 3D
[4065]375    call stopping_crit_h2o(delta_h2o_ads,delta_icetable,h2o_ice,d_h2oice,S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,stopcrit)
376    if (stopcrit%stop_code() > 0) return
377    call balance_h2oice_reservoirs(S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,h2o_ice,d_h2oice,d_h2oice_new)
378end if
[3149]379
[4074]380h2o_ice = h2o_ice + d_h2oice_new*dt
381zshift_surf = d_h2oice_new*dt/rho_h2oice
382
[4065]383END SUBROUTINE evolve_h2oice
[3989]384!=======================================================================
[3938]385
386!=======================================================================
[4065]387SUBROUTINE balance_h2oice_reservoirs(S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,h2o_ice,d_h2oice,d_h2oice_new)
[3991]388!-----------------------------------------------------------------------
389! NAME
390!     balance_h2oice_reservoirs
391!
392! DESCRIPTION
393!     Balance H2O flux from and into atmosphere across reservoirs.
394!
395! AUTHORS & DATE
396!     JB Clement, 2025
397!
398! NOTES
399!
400!-----------------------------------------------------------------------
[3938]401
[3991]402! DEPENDENCIES
403! ------------
[3989]404use evolution, only: dt
[4065]405use geometry,  only: ngrid, nslope
[3938]406
[3991]407! DECLARATION
408! -----------
[3938]409implicit none
410
[3991]411! ARGUMENTS
[3989]412! ---------
[4065]413real(qp),                 intent(in)    :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm
414real(dp), dimension(:,:), intent(in)    :: h2o_ice
415real(dp), dimension(:,:), intent(inout) :: d_h2oice
416real(dp), dimension(:,:), intent(out)   :: d_h2oice_new
[3938]417
[3991]418! LOCAL VARIABLES
[3989]419! ---------------
[4065]420integer(di) :: i, islope
[4074]421real(qp)    :: S_target, S_target_subl_h2oice, S_target_cond_h2oice, S_ghostice ! Balance variables
422real(dp)    :: d_target
[3989]423
[3991]424! CODE
[3989]425! ----
[4065]426S_target = (S_atm_2_h2o + S_h2o_2_atm)/2._dp
[3938]427S_target_cond_h2oice = S_atm_2_h2oice + S_target - S_atm_2_h2o
428S_target_subl_h2oice = S_h2oice_2_atm + S_target - S_h2o_2_atm
429
[4065]430d_h2oice_new = 0._dp
431S_ghostice = 0._qp
[3938]432do i = 1,ngrid
433    do islope = 1,nslope
[4065]434        if (d_h2oice(i,islope) > 0._dp) then ! Condensing
[4074]435            d_h2oice_new(i,islope) = d_h2oice(i,islope)*real(S_target_cond_h2oice/S_atm_2_h2oice,dp)
[4065]436        else if (d_h2oice(i,islope) < 0._dp) then ! Sublimating
437            d_target = d_h2oice(i,islope)*real(S_target_subl_h2oice/S_h2oice_2_atm,dp)
[4074]438            if (abs(d_target*dt) <= h2o_ice(i,islope)) then ! Enough ice to sublimate everything
439                d_h2oice_new(i,islope) = d_target
[3938]440            else ! Not enough ice to sublimate everything
441                ! We sublimate what we can
442                d_h2oice_new(i,islope) = h2o_ice(i,islope)/dt
443                ! It means the tendency is zero next time
[4065]444                d_h2oice(i,islope) = 0._dp
[3938]445                ! We compute the amount of H2O ice that we could not make sublimate
446                S_ghostice = S_ghostice + abs(d_target*dt) - h2o_ice(i,islope)
[4065]447            end if
448        end if
449    end do
450end do
[3149]451
[4074]452! We need to remove this ice unable to sublimate from places where ice condensed in order to keep balance
[4065]453where (d_h2oice_new > 0._dp) d_h2oice_new = d_h2oice_new*real((S_target_cond_h2oice - S_ghostice)/S_target_cond_h2oice,dp)
[3149]454
[3938]455END SUBROUTINE balance_h2oice_reservoirs
[3989]456!=======================================================================
[3553]457
[3989]458END MODULE surf_ice
Note: See TracBrowser for help on using the repository browser.