source: trunk/LMDZ.COMMON/libf/evolution/surface.F90 @ 4065

Last change on this file since 4065 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: 11.0 KB
Line 
1MODULE surface
2!-----------------------------------------------------------------------
3! NAME
4!     surface
5!
6! DESCRIPTION
7!     Contains global parameters used for the surface.
8!
9! AUTHORS & DATE
10!     JB Clement, 12/2025
11!
12! NOTES
13!
14!-----------------------------------------------------------------------
15
16! DEPENDENCIES
17! ------------
18use numerics, only: dp, di, k4
19
20! DECLARATION
21! -----------
22implicit none
23
24! PARAMETERS
25! ----------
26real(dp), dimension(:),   allocatable, protected :: albedodat_PCM  ! Albedo of bare ground
27real(dp), dimension(:,:), allocatable, protected :: albedo_PCM     ! Surface albedo_PCM
28real(dp), dimension(:,:), allocatable, protected :: emissivity_PCM ! Thermal IR surface emissivity_PCM
29
30contains
31!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32
33!=======================================================================
34SUBROUTINE ini_surface()
35!-----------------------------------------------------------------------
36! NAME
37!     ini_surface
38!
39! DESCRIPTION
40!     Initialize the parameters of module 'surface'.
41!
42! AUTHORS & DATE
43!     JB Clement, 12/2025
44!
45! NOTES
46!
47!-----------------------------------------------------------------------
48
49! DEPENDENCIES
50! ------------
51use geometry, only: ngrid, nslope
52
53! DECLARATION
54! -----------
55implicit none
56
57! CODE
58! ----
59if (.not. allocated(albedo_PCM)) allocate(albedo_PCM(ngrid,nslope))
60if (.not. allocated(albedodat_PCM)) allocate(albedodat_PCM(ngrid))
61if (.not. allocated(emissivity_PCM)) allocate(emissivity_PCM(ngrid,nslope))
62
63END SUBROUTINE ini_surface
64!=======================================================================
65
66!=======================================================================
67SUBROUTINE end_surface()
68!-----------------------------------------------------------------------
69! NAME
70!     end_surface
71!
72! DESCRIPTION
73!     Deallocate surface arrays.
74!
75! AUTHORS & DATE
76!     JB Clement, 12/2025
77!
78! NOTES
79!
80!-----------------------------------------------------------------------
81
82! DECLARATION
83! -----------
84implicit none
85
86! CODE
87! ----
88if (allocated(albedo_PCM)) deallocate(albedo_PCM)
89if (allocated(albedodat_PCM)) deallocate(albedodat_PCM)
90if (allocated(emissivity_PCM)) deallocate(emissivity_PCM )
91
92END SUBROUTINE end_surface
93!=======================================================================
94
95!=======================================================================
96SUBROUTINE set_albedodat_PCM(albedodat_PCM_in)
97!-----------------------------------------------------------------------
98! NAME
99!     set_albedodat_PCM
100!
101! DESCRIPTION
102!     Setter for 'albedodat_PCM'.
103!
104! AUTHORS & DATE
105!     JB Clement, 12/2025
106!
107! NOTES
108!
109!-----------------------------------------------------------------------
110
111! DECLARATION
112! -----------
113implicit none
114
115! ARGUMENTS
116! ---------
117real(dp), dimension(:), intent(in) :: albedodat_PCM_in
118
119! CODE
120! ----
121albedodat_PCM(:) = albedodat_PCM_in(:)
122
123END SUBROUTINE set_albedodat_PCM
124!=======================================================================
125
126!=======================================================================
127SUBROUTINE set_albedo_PCM(albedo_PCM_in)
128!-----------------------------------------------------------------------
129! NAME
130!     set_albedo_PCM
131!
132! DESCRIPTION
133!     Setter for 'albedo_PCM'.
134!
135! AUTHORS & DATE
136!     JB Clement, 12/2025
137!
138! NOTES
139!
140!-----------------------------------------------------------------------
141
142! DECLARATION
143! -----------
144implicit none
145
146! ARGUMENTS
147! ---------
148real(dp), dimension(:,:), intent(in) :: albedo_PCM_in
149
150! CODE
151! ----
152albedo_PCM(:,:) = albedo_PCM_in(:,:)
153
154END SUBROUTINE set_albedo_PCM
155!=======================================================================
156
157!=======================================================================
158SUBROUTINE set_emissivity_PCM(emissivity_PCM_in)
159!-----------------------------------------------------------------------
160! NAME
161!     set_emissivity_PCM
162!
163! DESCRIPTION
164!     Setter for 'emissivity_PCM'.
165!
166! AUTHORS & DATE
167!     JB Clement, 12/2025
168!
169! NOTES
170!
171!-----------------------------------------------------------------------
172
173! DECLARATION
174! -----------
175implicit none
176
177! ARGUMENTS
178! ---------
179real(dp), dimension(:,:), intent(in) :: emissivity_PCM_in
180
181! CODE
182! ----
183emissivity_PCM(:,:) = emissivity_PCM_in(:,:)
184
185END SUBROUTINE set_emissivity_PCM
186!=======================================================================
187
188!=======================================================================
189SUBROUTINE build4PCM_surf_rad_prop(h2o_ice,co2_ice,albedo4PCM,emissivity4PCM)
190!-----------------------------------------------------------------------
191! NAME
192!     build4PCM_surf_rad_prop
193!
194! DESCRIPTION
195!     Reconstructs albedo and emissivity for the PCM.
196!
197! AUTHORS & DATE
198!     JB Clement, 12/2025
199!     C. Metz, 01/2026
200!
201! NOTES
202!
203!-----------------------------------------------------------------------
204
205! DEPENDENCIES
206! ------------
207use geometry, only: ngrid, nslope, latitudes
208use frost,    only: h2o_frost4PCM, co2_frost4PCM
209use display,  only: print_msg
210
211! DECLARATION
212! -----------
213implicit none
214
215! ARGUMENTS
216! ---------
217real(dp), dimension(:,:), intent(in)  :: h2o_ice, co2_ice
218real(dp), dimension(:,:), intent(out) :: albedo4PCM, emissivity4PCM
219
220! LOCAL VARIABLES
221! ---------------
222integer(di)            :: i, islope, icap
223real(dp)               :: albedo_h2oice, albedo_h2ofrost, h2ofrost_albedo_threshold, emissivity_baresoil
224real(dp), dimension(2) :: albedo_co2frost, emissivity_co2ice, albedo_co2ice
225
226! CODE
227! ----
228! Fetch parameters from "callphys.def" and "startfi.nc"
229call read_albedo_emis(albedo_h2oice,albedo_h2ofrost,h2ofrost_albedo_threshold, &
230                      albedo_co2ice,albedo_co2frost,emissivity_co2ice,emissivity_baresoil)
231
232! Reconstruction Loop
233call print_msg('> Building albedo and emmissivity for the PCM')
234do i = 1,ngrid
235    ! Determine hemisphere: 1 = Northern, 2 = Southern
236    if (latitudes(i) < 0._dp) then
237        icap = 2
238    else
239        icap = 1
240    end if
241
242    do islope = 1,nslope
243        ! Bare ground (default)
244        albedo4PCM(i,islope) = albedodat_PCM(i)
245        emissivity4PCM(i,islope) = emissivity_baresoil
246
247        ! H2O ice
248        if (h2o_ice(i,islope) > 0._dp) then
249            albedo4PCM(i,islope) = albedo_h2oice
250            emissivity4PCM(i,islope) = 1._dp
251        end if
252
253        ! CO2 ice (dominant over H2O ice)
254        if (co2_ice(i,islope) > 0._dp) then
255            albedo4PCM(i,islope) = albedo_co2ice(icap)
256            emissivity4PCM(i,islope) = albedo_co2ice(icap)
257        end if
258
259        ! H2O frost (subject to threshold)
260        if (h2o_frost4PCM(i,islope) > h2ofrost_albedo_threshold) then
261            albedo4PCM(i,islope) = albedo_h2ofrost
262            emissivity4PCM(i,islope) = 1._dp
263        end if
264
265        ! CO2 frost (final priority)
266        if (co2_frost4PCM(i,islope) > 0.) then
267            albedo4PCM(i,islope) = albedo_co2frost(icap)
268            emissivity4PCM(i,islope) = emissivity_co2ice(icap)
269        end if
270    end do
271end do
272
273END SUBROUTINE build4PCM_surf_rad_prop
274!=======================================================================
275
276!=======================================================================
277SUBROUTINE read_albedo_emis(albedo_h2oice,albedo_h2ofrost,h2ofrost_albedo_threshold, &
278                            albedo_co2ice,albedo_co2frost,emissivity_co2ice,emissivity_baresoil)
279!-----------------------------------------------------------------------
280! NAME
281!     read_albedo_emis
282!
283! DESCRIPTION
284!     Reads albedo/emissivity parameters from "callphys.def" and
285!     "startfi.nc" ('controle').
286!
287! AUTHORS & DATE
288!     C. Metz, 01/2026
289!
290! NOTES
291!
292!-----------------------------------------------------------------------
293
294! DEPENDENCIES
295! ------------
296#ifdef CPP_IOIPSL
297use IOIPSL,          only: getin
298#else
299use ioipsl_getincom, only: getin
300#endif
301use config,          only: callphys_name
302use io_netcdf,       only: open_nc, close_nc, startfi_name, get_dim_nc, get_var_nc
303use stoppage,        only: stop_clean
304
305! DECLARATION
306! -----------
307implicit none
308
309! ARGUMENTS
310! ---------
311real(dp),               intent(out) :: albedo_h2oice, albedo_h2ofrost, h2ofrost_albedo_threshold, emissivity_baresoil
312real(dp), dimension(2), intent(out) :: albedo_co2frost, emissivity_co2ice, albedo_co2ice
313
314! LOCAL VARIABLES
315! ---------------
316logical(k4)                         :: here, cst_h2oice_albedo
317integer(di)                         :: i, nindex
318real(dp), dimension(:), allocatable :: controle
319
320! CODE
321! ----
322inquire(file = callphys_name,exist = here)
323if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find "'//callphys_name//'"!', 1)
324
325! Read albedos of H2O frost and perennial ice from "callphys.def"
326albedo_h2oice = 0.35_dp ! Default
327call getin('albedo_h2o_cap',albedo_h2oice)
328if (albedo_h2oice < 0._dp .or. albedo_h2oice > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_h2o_cap'' in "'//callphys_name//'" is out of bounds [0,1]!',1)
329
330albedo_h2ofrost = 0.35_dp ! Default
331call getin('albedo_h2ofrost',albedo_h2ofrost)
332if (albedo_h2ofrost < 0._dp .or. albedo_h2ofrost > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_h2ofrost'' in "'//callphys_name//'" is out of bounds [0,1]!',1)
333
334h2ofrost_albedo_threshold = 0.005_dp ! Default [kg/m2]
335call getin('frost_albedo_threshold',h2ofrost_albedo_threshold)
336
337cst_h2oice_albedo = .false. ! Default
338call getin('cst_cap_albedo',cst_h2oice_albedo)
339if (cst_h2oice_albedo) albedo_h2ofrost = albedo_h2oice ! If true, then we don't account for water frost albedo effect
340
341! Read albedos of CO2 perennial ice from "callphys.def"
342albedo_co2ice(1) = 0.6_dp ! Default, north hemisphere
343albedo_co2ice(2) = 0.85_dp ! Default, south hemisphere
344call getin('albedo_co2ice_north',albedo_co2ice(1))
345call getin('albedo_co2ice_south',albedo_co2ice(2))
346if (albedo_co2ice(1) < 0._dp .or. albedo_co2ice(1) > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_co2ice_north'' in "'//callphys_name//'" is out of bounds [0,1]!',1)
347if (albedo_co2ice(2) < 0._dp .or. albedo_co2ice(2) > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_co2ice_south'' in "'//callphys_name//'" is out of bounds [0,1]!',1)
348
349! Read albedos of CO2 frost from "callphys.def"
350call open_nc(startfi_name,'read')
351call get_dim_nc('index',nindex)
352allocate(controle(nindex))
353call get_var_nc('controle',controle)
354call close_nc(startfi_name)
355albedo_co2frost(1)   = controle(22)
356albedo_co2frost(2)   = controle(23)
357emissivity_co2ice(1) = controle(24)
358emissivity_co2ice(2) = controle(25)
359emissivity_baresoil  = controle(26)
360deallocate(controle)
361do i = 1,2
362    if (albedo_co2frost(i) < 0._dp .or. albedo_co2frost(i) > 1._dp) call stop_clean(__FILE__,__LINE__,'''albedo_co2frost'' from ''controle(22:23)'' in "'//startfi_name//'" is out of bounds [0,1]!',1)
363    if (emissivity_co2ice(i) < 0._dp .or. emissivity_co2ice(i) > 1._dp) call stop_clean(__FILE__,__LINE__,'''emissivity_co2ice'' from ''controle(24:25)'' in "'//startfi_name//'" is out of bounds [0,1]!',1)
364end do
365if (emissivity_baresoil < 0._dp .or. emissivity_baresoil > 1._dp) call stop_clean(__FILE__,__LINE__,'''emissivity_baresoil'' from ''controle(26)'' in "'//startfi_name//'" is out of bounds [0,1]!',1)
366
367END SUBROUTINE read_albedo_emis
368!=======================================================================
369
370END MODULE surface
Note: See TracBrowser for help on using the repository browser.