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

Last change on this file since 4068 was 4065, checked in by jbclement, 3 weeks ago

PEM:
Major refactor following the previous ones (r3989 and r3991) completing the large structural reorganization and cleanup of the PEM codebase. This revision introduces newly designed modules, standardizes interfaces with explicit ini/end APIs and adds native NetCDF I/O together with explicit PCM/PEM adapters. In detail:

  • Some PEM models were corrected or improved:
    • Frost/perennial ice semantics are clarified via renaming;
    • Soil temperature remapping clarified, notably by removing the rescaling of temperature deviation;
    • Geothermal flux for the PCM is computed based on the PEM state;
  • New explicit PEM/PCM adapters ("set_*"/"build4PCM_*") to decouple PEM internal representation from PCM file layouts and reconstruct consistent fields returned to the PCM;
  • New explicit build/teardown routines that centralize allocation and initialization ordering, reducing accidental use of uninitialized data and making the model lifecycle explicit;
  • Add native read/write helpers for NetCDF that centralize all low-level NetCDF interactions with major improvements (and more simplicity) compared to legacy PEM/PCM I/O (see the modules "io_netcdf" and "output"). They support reading, creation and writing of "diagevol.nc" (renamed from "diagpem.nc") and start/restart files;
  • Provide well-focused modules ("numerics"/"maths"/"utility"/"display") to host commonly-used primitives:
    • "numerics" defines numerical types and constants for reproducibility, portability across compilers and future transitions (e.g. quadruple precision experiments);
    • "display" provides a single controlled interface for runtime messages, status output and diagnostics, avoiding direct 'print'/'write' to enable silent mode, log redirection, and MPI-safe output in the future.
    • "utility" (new module) hosts generic helpers used throughout the code (e.g. "int2str" or "real2str");
  • Add modules "clim_state_init"/"clim_state_rec" which provide robust read/write logic for "start/startfi/startpem", including 1D fallbacks, mesh conversions and dimension checks. NetCDF file creation is centralized and explicit. Restart files are now self-consistent and future-proof, requiring changes only to affected variables;
  • Add module "atmosphere" which computes pressure fields, reconstructs potential temperature and air mass. It also holds the whole logic to define sigma or hybrid coordinates for altitudes;
  • Add module "geometry" to centrilize dimensions logic and grid conversions routines (including 2 new ones "dyngrd2vect"/"vect2dyngrd");
  • Add module "slopes" to isolate slopes handling;
  • Add module "surface" to isolate surface management. Notably, albedo and emissivity are now fully reconstructed following the PCM settings;
  • Add module "allocation" to check the dimension initialization and centrilize allocation/deallocation;
  • Finalize module decomposition and renaming to consolidate domain-based modules, purpose-based routines and physics/process-based variables;
  • The main program now drives a clearer sequence of conceptual steps (initialization / reading / evolution / update / build / writing) and fails explicitly instead of silently defaulting;
  • Ice table logic is made restart-safe;
  • 'Open'/'read' intrinsic logic is made safe and automatic;
  • Improve discoverability and standardize the data handling (private vs protected vs public);
  • Apply consistent documentation/header style already introduced;
  • Update deftank scripts to reflect new names and launch wrappers;

This revision is a structural milestone aiming to be behavior-preserving where possible. It has been tested via compilation and short integration runs. However, due to extensive renames, moves, and API changes, full validation is still ongoing.
Note: the revision includes one (possibly two) easter egg hidden in the code for future archaeologists of the PEM. No physics were harmed.
JBC

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