source: trunk/LMDZ.COMMON/libf/evolution/frost.F90 @ 4066

Last change on this file since 4066 was 4065, checked in by jbclement, 2 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: 5.0 KB
Line 
1MODULE frost
2!-----------------------------------------------------------------------
3! NAME
4!     frost
5!
6! DESCRIPTION
7!     Module for managing frost variables.
8!
9! AUTHORS & DATE
10!     JB Clement, 12/2025
11!
12! NOTES
13!
14!-----------------------------------------------------------------------
15
16! DEPENDENCIES
17! ------------
18use numerics, only: dp, di
19
20! DECLARATION
21! -----------
22implicit none
23
24! VARIABLES
25! ---------
26! Different types of frost retained by the PEM to give back to the PCM at the end
27real(dp), dimension(:,:), allocatable :: h2o_frost4PCM
28real(dp), dimension(:,:), allocatable :: co2_frost4PCM
29
30! PARAMETERS
31! ----------
32! Different types of frost in the PCM at the beginning [Pa]
33real(dp), dimension(:,:), allocatable, protected :: h2ofrost_PCM
34real(dp), dimension(:,:), allocatable, protected :: co2frost_PCM
35
36contains
37!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
38
39!=======================================================================
40SUBROUTINE ini_frost()
41!-----------------------------------------------------------------------
42! NAME
43!     ini_frost
44!
45! DESCRIPTION
46!     Initialize the parameters of module 'frost'.
47!
48! AUTHORS & DATE
49!     JB Clement, 12/2025
50!
51! NOTES
52!
53!-----------------------------------------------------------------------
54
55! DEPENDENCIES
56! ------------
57use geometry, only: ngrid, nslope
58
59! DECLARATION
60! -----------
61implicit none
62
63! CODE
64! ----
65if (.not. allocated(h2ofrost_PCM)) allocate(h2ofrost_PCM(ngrid,nslope))
66if (.not. allocated(co2frost_PCM)) allocate(co2frost_PCM(ngrid,nslope))
67if (.not. allocated(h2o_frost4PCM)) allocate(h2o_frost4PCM(ngrid,nslope))
68if (.not. allocated(co2_frost4PCM)) allocate(co2_frost4PCM(ngrid,nslope))
69
70END SUBROUTINE ini_frost
71!=======================================================================
72
73!=======================================================================
74SUBROUTINE end_frost()
75!-----------------------------------------------------------------------
76! NAME
77!     end_frost
78!
79! DESCRIPTION
80!     Deallocate frost arrays.
81!
82! AUTHORS & DATE
83!     JB Clement, 12/2025
84!
85! NOTES
86!
87!-----------------------------------------------------------------------
88
89! DECLARATION
90! -----------
91implicit none
92
93! CODE
94! ----
95if (allocated(h2ofrost_PCM)) deallocate(h2ofrost_PCM)
96if (allocated(co2frost_PCM)) deallocate(co2frost_PCM)
97if (allocated(h2o_frost4PCM)) deallocate(h2o_frost4PCM)
98if (allocated(co2_frost4PCM)) deallocate(co2_frost4PCM)
99
100END SUBROUTINE end_frost
101!=======================================================================
102
103!=======================================================================
104SUBROUTINE set_h2ofrost_PCM(h2ofrost_PCM_in)
105!-----------------------------------------------------------------------
106! NAME
107!     set_h2ofrost_PCM
108!
109! DESCRIPTION
110!     Setter for 'h2ofrost_PCM'.
111!
112! AUTHORS & DATE
113!     JB Clement, 12/2025
114!
115! NOTES
116!
117!-----------------------------------------------------------------------
118
119! DECLARATION
120! -----------
121implicit none
122
123! ARGUMENTS
124! ---------
125real(dp), dimension(:,:), intent(in) :: h2ofrost_PCM_in
126
127! CODE
128! ----
129h2ofrost_PCM(:,:) = h2ofrost_PCM_in(:,:)
130
131END SUBROUTINE set_h2ofrost_PCM
132!=======================================================================
133
134!=======================================================================
135SUBROUTINE set_co2frost_PCM(co2frost_PCM_in)
136!-----------------------------------------------------------------------
137! NAME
138!     set_co2frost_PCM
139!
140! DESCRIPTION
141!     Setter for 'co2frost_PCM'.
142!
143! AUTHORS & DATE
144!     JB Clement, 12/2025
145!
146! NOTES
147!
148!-----------------------------------------------------------------------
149
150! DECLARATION
151! -----------
152implicit none
153
154! ARGUMENTS
155! ---------
156real(dp), dimension(:,:), intent(in) :: co2frost_PCM_in
157
158! CODE
159! ----
160co2frost_PCM(:,:) = co2frost_PCM_in(:,:)
161
162END SUBROUTINE set_co2frost_PCM
163!=======================================================================
164
165!=======================================================================
166SUBROUTINE compute_frost4PCM(min_h2ofrost,min_co2frost)
167!-----------------------------------------------------------------------
168! NAME
169!     compute_frost4PCM
170!
171! DESCRIPTION
172!     Compute the frost to give back to the PCM (metamorphism).
173!
174! AUTHORS & DATE
175!     JB Clement, 12/2025
176!
177! NOTES
178!     Frost for the PEM is the extra part of the PCM frost above the
179!     yearly minimum.
180!-----------------------------------------------------------------------
181
182! DEPENDENCIES
183! ------------
184use display, only: print_msg
185
186! DECLARATION
187! -----------
188implicit none
189
190! ARGUMENTS
191! ---------
192real(dp), dimension(:,:), intent(in) :: min_h2ofrost, min_co2frost
193
194! CODE
195! ----
196call print_msg('> Computing frost to give back to the PCM (metamorphism)')
197
198h2o_frost4PCM(:,:) = 0._dp
199co2_frost4PCM(:,:) = 0._dp
200where (h2ofrost_PCM(:,:) > 0._dp) h2o_frost4PCM(:,:) = h2ofrost_PCM(:,:) - min_h2ofrost(:,:)
201where (co2frost_PCM(:,:) > 0._dp) co2_frost4PCM(:,:) = co2frost_PCM(:,:) - min_co2frost(:,:)
202
203END SUBROUTINE compute_frost4PCM
204!=======================================================================
205
206END MODULE frost
Note: See TracBrowser for help on using the repository browser.