source: trunk/LMDZ.COMMON/libf/evolution/slopes.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: 4.5 KB
RevLine 
[4065]1MODULE slopes
2!-----------------------------------------------------------------------
3! NAME
4!     slopes
5!
6! DESCRIPTION
7!     Contains global parameters used for the slopes.
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! PARAMETERS
25! ----------
26real(dp), dimension(:),   allocatable, protected :: def_slope_mean ! Mean slople of each bin [degree]
27real(dp), dimension(:,:), allocatable, protected :: subslope_dist  ! Distribution of the slopes
28integer(di),                           protected :: iflat          ! Index of the flat slope
29
30contains
31!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32
33!=======================================================================
34SUBROUTINE ini_slopes()
35!-----------------------------------------------------------------------
36! NAME
37!     ini_slopes
38!
39! DESCRIPTION
40!     Initialize the parameters of module 'slopes'.
41!
42! AUTHORS & DATE
43!     JB Clement, 12/2025
44!
45! NOTES
46!
47!-----------------------------------------------------------------------
48
49! DEPENDENCIES
50! ------------
51use geometry, only: nslope, ngrid
52
53! DECLARATION
54! -----------
55implicit none
56
57! CODE
58! ----
59if (.not. allocated(def_slope_mean)) allocate(def_slope_mean(nslope))
60if (.not. allocated(subslope_dist)) allocate(subslope_dist(ngrid,nslope))
61
62END SUBROUTINE ini_slopes
63!=======================================================================
64
65!=======================================================================
66SUBROUTINE end_slopes()
67!-----------------------------------------------------------------------
68! NAME
69!     end_slopes
70!
71! DESCRIPTION
72!     Deallocate slopes arrays.
73!
74! AUTHORS & DATE
75!     JB Clement, 12/2025
76!
77! NOTES
78!
79!-----------------------------------------------------------------------
80
81! DECLARATION
82! -----------
83implicit none
84
85! CODE
86! ----
87if (allocated(def_slope_mean)) deallocate(def_slope_mean)
88if (allocated(subslope_dist)) deallocate(subslope_dist)
89
90END SUBROUTINE end_slopes
91!=======================================================================
92
93!=======================================================================
94SUBROUTINE set_def_slope_mean(def_slope_mean_in)
95!-----------------------------------------------------------------------
96! NAME
97!     set_def_slope_mean
98!
99! DESCRIPTION
100!     Setter for 'def_slope_mean'.
101!
102! AUTHORS & DATE
103!     JB Clement, 12/2025
104!
105! NOTES
106!
107!-----------------------------------------------------------------------
108
109! DECLARATION
110! -----------
111implicit none
112
113! ARGUMENTS
114! ---------
115real(dp), dimension(:), intent(in) :: def_slope_mean_in
116
117! CODE
118! ----
119def_slope_mean(:) = def_slope_mean_in(:)
120
121END SUBROUTINE set_def_slope_mean
122!=======================================================================
123
124!=======================================================================
125SUBROUTINE set_subslope_dist(subslope_dist_in)
126!-----------------------------------------------------------------------
127! NAME
128!     set_subslope_dist
129!
130! DESCRIPTION
131!     Setter for 'subslope_dist'.
132!
133! AUTHORS & DATE
134!     JB Clement, 12/2025
135!
136! NOTES
137!
138!-----------------------------------------------------------------------
139
140! DECLARATION
141! -----------
142implicit none
143
144! ARGUMENTS
145! ---------
146real(dp), dimension(:,:), intent(in) :: subslope_dist_in
147
148! CODE
149! ----
150subslope_dist(:,:) = subslope_dist_in(:,:)
151
152END SUBROUTINE set_subslope_dist
153!=======================================================================
154
155!=======================================================================
156SUBROUTINE set_iflat()
157!-----------------------------------------------------------------------
158! NAME
159!     set_iflat
160!
161! DESCRIPTION
162!     Setter for 'iflat'.
163!
164! AUTHORS & DATE
165!     JB Clement, 12/2025
166!
167! NOTES
168!
169!-----------------------------------------------------------------------
170
171! DEPENDENCIES
172! ------------
173use geometry, only: nslope
174use display,  only: print_msg
175use utility,  only: int2str, real2str
176
177! DECLARATION
178! -----------
179implicit none
180
181! LOCAL VARIABLES
182! ---------------
183integer(di) :: islope
184
185! CODE
186! ----
187iflat = 1
188do islope = 2,nslope
189    if (abs(def_slope_mean(islope)) < abs(def_slope_mean(iflat))) iflat = islope
190end do
191call print_msg('Flat slope for islope = '//int2str(iflat))
192call print_msg('Corresponding criterium = '//real2str(def_slope_mean(iflat)))
193
194END SUBROUTINE set_iflat
195!=======================================================================
196
197END MODULE slopes
Note: See TracBrowser for help on using the repository browser.