source: trunk/LMDZ.COMMON/libf/evolution/clim_state_rec.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: 19.0 KB
Line 
1MODULE clim_state_rec
2!-----------------------------------------------------------------------
3! NAME
4!     clim_state_rec
5!
6! DESCRIPTION
7!     Write the restart files to save the climate state.
8!
9! AUTHORS & DATE
10!     JB Clement, 12/2025
11!
12! NOTES
13!
14!-----------------------------------------------------------------------
15
16! DEPENDENCIES
17! ------------
18use numerics,  only: dp, di, k4
19use io_netcdf, only: open_nc, close_nc, check_nc, put_var_nc, get_dim_nc, get_var_nc, start_name, start1D_name, startfi_name, startpem_name
20
21! DECLARATION
22! -----------
23implicit none
24
25! VARIABLES
26! ---------
27logical(k4), private :: is_restartpem = .false. ! Flag to know if "restartpem.nc" exists
28
29contains
30!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
31
32!=======================================================================
33SUBROUTINE write_restart(ps4PCM,pa4PCM,preff4PCM,tsurf4PCM,q4PCM,teta4PCM,air_mass4PCM)
34!-----------------------------------------------------------------------
35! NAME
36!    write_restart
37!
38! DESCRIPTION
39!    Write the file "restart.nc".
40!
41! AUTHORS & DATE
42!    JB Clement, 12/2025
43!
44! NOTES
45!
46!-----------------------------------------------------------------------
47
48! DEPENDENCIES
49! ------------
50use geometry, only: vect2dyngrd, ngrid, nlon, nlat, nlayer
51use tracers,  only: qnames, nq
52use stoppage, only: stop_clean
53use display,  only: print_msg
54
55! DECLARATION
56! -----------
57implicit none
58
59! ARGUMENTS
60! ---------
61real(dp),                   intent(in) :: pa4PCM, preff4PCM
62real(dp), dimension(:),     intent(in) :: ps4PCM
63real(dp), dimension(:,:),   intent(in) :: tsurf4PCM, teta4PCM, air_mass4PCM
64real(dp), dimension(:,:,:), intent(in) :: q4PCM
65
66! LOCAL VARIABLES
67! ---------------
68integer(di)                               :: cstat, l, i
69real(dp), dimension(nlon + 1,nlat)        :: ps_dyn
70real(dp), dimension(nlon + 1,nlat,nlayer) :: var_dyn
71
72! CODE
73! ----
74! In case of 1D
75! ~~~~~~~~~~~~~
76if (ngrid == 1) then
77    call write_restart1D(ps4PCM,pa4PCM,preff4PCM,tsurf4PCM,q4PCM)
78    return
79end if
80
81! In case of 3D
82! ~~~~~~~~~~~~~
83! Copy "start.nc" into "restart.nc"
84call print_msg('> Writing "re'//start_name//'"')
85call execute_command_line('cp '//start_name//' re'//start_name,cmdstat = cstat)
86if (cstat > 0) then
87    call stop_clean(__FILE__,__LINE__,'command execution failed!',1)
88else if (cstat < 0) then
89    call stop_clean(__FILE__,__LINE__,'command execution not supported!',1)
90end if
91
92! Rewrite the variables modified by the PEM
93call open_nc('re'//start_name,'write')
94
95! Surface pressure
96call vect2dyngrd(nlon + 1,nlat,ngrid,ps4PCM,ps_dyn)
97call put_var_nc('ps',ps_dyn)
98
99! Potential temperature
100call vect2dyngrd(nlon + 1,nlat,ngrid,teta4PCM,var_dyn)
101call put_var_nc('teta',var_dyn)
102
103! Air mass
104call vect2dyngrd(nlon + 1,nlat,ngrid,air_mass4PCM,var_dyn)
105call put_var_nc('masse',var_dyn)
106
107! Tracers
108do i = 1,nq
109    do l = 1,nlayer
110        call vect2dyngrd(nlon + 1,nlat,ngrid,q4PCM(:,l,i),var_dyn(:,:,l))
111    end do
112    call put_var_nc(qnames(i),var_dyn)
113end do
114
115! Close
116call close_nc('re'//start_name)
117
118END SUBROUTINE write_restart
119!=======================================================================
120
121!=======================================================================
122SUBROUTINE write_restart1D(ps4PCM,pa4PCM,preff4PCM,tsurf4PCM,q4PCM)
123!-----------------------------------------------------------------------
124! NAME
125!    write_restart1D
126!
127! DESCRIPTION
128!    Write the file "restart1D.txt".
129!
130! AUTHORS & DATE
131!    JB Clement, 01/2026
132!
133! NOTES
134!
135!-----------------------------------------------------------------------
136
137! DEPENDENCIES
138! ------------
139use geometry,   only: nlayer, nslope
140use atmosphere, only: teta_PCM, u_PCM, v_PCM
141use tracers,    only: nq, qnames
142use stoppage,   only: stop_clean
143use display,    only: print_msg
144
145! DECLARATION
146! -----------
147implicit none
148
149! ARGUMENTS
150! ---------
151real(dp),                   intent(in) :: pa4PCM, preff4PCM
152real(dp), dimension(:),     intent(in) :: ps4PCM
153real(dp), dimension(:,:),   intent(in) :: tsurf4PCM
154real(dp), dimension(:,:,:), intent(in) :: q4PCM
155
156! LOCAL VARIABLES
157! ---------------
158integer(di) :: funit, ierr, i, j, l
159
160! CODE
161! ----
162! Write "restart1D.txt"
163call print_msg('> Writing "re'//start1D_name//'"')
164open(newunit = funit,file = 're'//start1D_name,status = "replace",action = "write",iostat = ierr)
165if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening file "re'//start1D_name//'"!',ierr)
166write(funit,*) 'ps', ps4PCM(1), pa4PCM, preff4PCM
167do i = 1,nq
168    write(funit,*) qnames(i), (q4PCM(1,l,i), l = 1,nlayer)
169end do
170write(funit,*) 'u', (u_PCM(1,l), l = 1,nlayer)
171write(funit,*) 'v', (v_PCM(1,l), l = 1,nlayer)
172write(funit,*) 'teta', (tsurf4PCM(1,j), j = 1,nslope), (teta_PCM(1,l), l = 1,nlayer)
173close(funit)
174
175END SUBROUTINE write_restart1D
176!=======================================================================
177
178!=======================================================================
179SUBROUTINE write_restartfi(is_h2o_perice,h2o_ice4PCM,co2_ice4PCM,tsurf4PCM,tsoil4PCM,inertiesoil4PCM,albedo4PCM,emissivity4PCM,flux_geo4PCM)
180!-----------------------------------------------------------------------
181! NAME
182!    write_restartfi
183!
184! DESCRIPTION
185!    Write the file "restartfi.nc".
186!
187! AUTHORS & DATE
188!    JB Clement, 01/2026
189!
190! NOTES
191!
192!-----------------------------------------------------------------------
193
194! DEPENDENCIES
195! ------------
196use orbit,    only: obliquity, aphelion, perihelion, date_peri
197use frost,    only: h2o_frost4PCM, co2_frost4PCM
198use stoppage, only: stop_clean
199use display,  only: print_msg
200
201! DECLARATION
202! -----------
203implicit none
204
205! ARGUMENTS
206! ---------
207real(dp),    dimension(:,:),   intent(in) :: h2o_ice4PCM, co2_ice4PCM, albedo4PCM, emissivity4PCM, tsurf4PCM, flux_geo4PCM
208real(dp),    dimension(:,:,:), intent(in) :: tsoil4PCM, inertiesoil4PCM
209logical(k4), dimension(:),     intent(in) :: is_h2o_perice
210
211! LOCAL VARIABLES
212! ---------------
213integer(di)                         :: cstat
214integer(di)                         :: nindex ! Size of dimension 'index'
215real(dp), dimension(:), allocatable :: controle
216
217! CODE
218! ----
219! Copy "startfi.nc" into "restartfi.nc"
220call print_msg('> Writing "re'//startfi_name//'"')
221call execute_command_line('cp '//startfi_name//' re'//startfi_name,cmdstat = cstat)
222if (cstat > 0) then
223    call stop_clean(__FILE__,__LINE__,'command execution failed!',1)
224else if (cstat < 0) then
225    call stop_clean(__FILE__,__LINE__,'command execution not supported!',1)
226end if
227
228! Load the variable 'controle' to modify it with new values
229call open_nc('re'//startfi_name,'read')
230call get_dim_nc('index',nindex)
231allocate(controle(nindex))
232call get_var_nc('controle',controle)
233call close_nc('re'//startfi_name)
234
235! Rewrite the variables modified by the PEM
236call open_nc('re'//startfi_name,'write')
237
238! Orbital parameters (controle)
239controle(18) = obliquity ! Obliquity
240controle(15) = perihelion ! Perihelion
241controle(16) = aphelion ! Aphelion
242controle(17) = date_peri ! Date of perihelion
243call put_var_nc('controle',controle)
244deallocate(controle)
245
246! Variables that have been modified
247call put_var_nc('watercaptag',merge(1.,0.,is_h2o_perice))
248call put_var_nc('watercap',h2o_ice4PCM,1)
249call put_var_nc('h2o_ice',h2o_frost4PCM,1)
250call put_var_nc('co2',co2_frost4PCM,1)
251call put_var_nc('perennial_co2ice',co2_ice4PCM,1)
252call put_var_nc('tsurf',tsurf4PCM,1)
253call put_var_nc('tsoil',tsoil4PCM,1)
254call put_var_nc('inertiesoil',inertiesoil4PCM,1)
255call put_var_nc('albedo',albedo4PCM,1)
256call put_var_nc('emis',emissivity4PCM,1)
257call put_var_nc('flux_geo',flux_geo4PCM,1)
258
259! h2oice_depth
260! lag_co2_ice
261! zdqsdif_ssi_tot
262! d_coef
263
264call close_nc('re'//startfi_name)
265
266END SUBROUTINE write_restartfi
267!=======================================================================
268
269!=======================================================================
270SUBROUTINE create_nc_pem(nb_str_max)
271!-----------------------------------------------------------------------
272! NAME
273!     create_nc_pem
274!
275! DESCRIPTION
276!     Create a NetCDF file for the PEM.
277!
278! AUTHORS & DATE
279!     JB Clement, 01/2026
280!
281! NOTES
282!
283!-----------------------------------------------------------------------
284
285! DEPENDENCIES
286! ------------
287use netcdf,           only: nf90_double, nf90_clobber, nf90_unlimited, nf90_global, &
288                            nf90_create, nf90_enddef, nf90_def_dim, nf90_def_var, nf90_put_att
289use geometry,         only: dim_init, ngrid, nsoil, nslope, longitudes, latitudes, cell_area
290use stoppage,         only: stop_clean
291use soil,             only: do_soil, mlayer
292use sorption,         only: do_sorption
293use layered_deposits, only: do_layering
294use display,          only: print_msg
295use utility,          only: int2str
296
297! DECLARATION
298! -----------
299implicit none
300
301! ARGUMENTS
302! ---------
303integer(di), intent(in) :: nb_str_max
304
305! LOCAL VARIABLES
306! ---------------
307integer(di) :: ncid  ! File ID
308integer(di) :: varid ! Variable ID
309integer(di) :: dim_ngrid, dim_nsoil, dim_nslope, dim_time, dim_nb_str
310
311! CODE
312! ----
313! Check if dimensions are well initialized
314if (.not. dim_init) call stop_clean(__FILE__,__LINE__,'dimensions are not initilized in the PEM!',1)
315
316! Create file
317call print_msg('> Creating "re'//startpem_name//'"')
318call check_nc(nf90_create('re'//startpem_name,nf90_clobber,ncid),'creating re'//startpem_name) ! Enter define mode
319
320! Define dimensions
321call check_nc(nf90_def_dim(ncid,'Time',nf90_unlimited,dim_time),'defining dimension Time')
322call check_nc(nf90_def_dim(ncid,'physical_points',ngrid,dim_ngrid),'defining dimension physical_points')
323call check_nc(nf90_def_dim(ncid,'subsurface_layers',nsoil,dim_nsoil),'defining dimension subsurface_layers')
324call check_nc(nf90_def_dim(ncid,'nslope',nslope,dim_nslope),'defining dimension nslope')
325if (do_layering) then
326    call print_msg('nb_str_max = '//int2str(nb_str_max))
327    call check_nc(nf90_def_dim(ncid,'nb_str_max',nb_str_max,dim_nb_str),'defining dimension nb_str_max')
328end if
329
330! Define variables
331call check_nc(nf90_def_var(ncid,'Time',nf90_double,(/dim_time/),varid),'defining variable Time')
332call check_nc(nf90_put_att(ncid,varid,'title','Year of simulation'),'putting title attribute for Time')
333call check_nc(nf90_put_att(ncid,varid,'units','Planetary year'),'putting units attribute for Time')
334
335call check_nc(nf90_def_var(ncid,'longitude',nf90_double,(/dim_ngrid/),varid),'defining variable longitude')
336call check_nc(nf90_put_att(ncid,varid,'title','Longitudes of the grid'),'putting title attribute for longitude')
337
338call check_nc(nf90_def_var(ncid,'latitude',nf90_double,(/dim_ngrid/),varid),'defining variable latitude')
339call check_nc(nf90_put_att(ncid,varid,'title','Latitudes of the grid'),'putting title attribute for latitude')
340
341call check_nc(nf90_def_var(ncid,'cell_area',nf90_double,(/dim_ngrid/),varid),'defining variable cell_area')
342call check_nc(nf90_put_att(ncid,varid,'title','Cell area'),'putting title attribute for cell_area')
343
344call check_nc(nf90_def_var(ncid,'h2o_ice',nf90_double,(/dim_ngrid,dim_nslope,dim_time/),varid),'defining variable h2o_ice')
345call check_nc(nf90_put_att(ncid,varid,'title','H2O ice'),'putting title attribute for h2o_ice')
346
347call check_nc(nf90_def_var(ncid,'co2_ice',nf90_double,(/dim_ngrid,dim_nslope,dim_time/),varid),'defining variable co2_ice')
348call check_nc(nf90_put_att(ncid,varid,'title','CO2 ice'),'putting title attribute for co2_ice')
349
350if (do_soil) then
351    call check_nc(nf90_def_var(ncid,'soildepth',nf90_double,(/dim_nsoil/),varid),'defining variable soildepth')
352    call check_nc(nf90_put_att(ncid,varid,'title','Depths of soil layers'),'putting title attribute for soildepth')
353
354    call check_nc(nf90_def_var(ncid,'tsoil',nf90_double,(/dim_ngrid,dim_nsoil,dim_nslope,dim_time/),varid),'defining variable tsoil')
355    call check_nc(nf90_put_att(ncid,varid,'title','Soil temperature'),'putting title attribute for tsoil')
356
357    call check_nc(nf90_def_var(ncid,'TI',nf90_double,(/dim_ngrid,dim_nsoil,dim_nslope,dim_time/),varid),'defining variable inertiedat')
358    call check_nc(nf90_put_att(ncid,varid,'title','Thermal inertie of PEM'),'putting title attribute for TI')
359
360    call check_nc(nf90_def_var(ncid,'inertiedat',nf90_double,(/dim_ngrid,dim_nsoil,dim_time/),varid),'defining variable longitude')
361    call check_nc(nf90_put_att(ncid,varid,'title','Soil thermal inertia'),'putting title attribute for inertiedat')
362
363    call check_nc(nf90_def_var(ncid,'icetable_depth',nf90_double,(/dim_ngrid,dim_nslope,dim_time/),varid),'defining variable icetable_depth')
364    call check_nc(nf90_put_att(ncid,varid,'title','Depth of ice table'),'putting title attribute for icetable_depth')
365
366    call check_nc(nf90_def_var(ncid,'icetable_thickness',nf90_double,(/dim_ngrid,dim_nslope,dim_time/),varid),'defining variable icetable_thickness')
367    call check_nc(nf90_put_att(ncid,varid,'title','Thickness of ice table'),'putting title attribute for icetable_thickness')
368
369    call check_nc(nf90_def_var(ncid,'ice_porefilling',nf90_double,(/dim_ngrid,dim_nsoil,dim_nslope,dim_time/),varid),'defining variable ice_porefilling')
370    call check_nc(nf90_put_att(ncid,varid,'title','Subsurface ice pore filling'),'putting title attribute for ice_porefilling')
371
372    if (do_sorption) then
373        call check_nc(nf90_def_var(ncid,'h2o_ads_reg',nf90_double,(/dim_ngrid,dim_nsoil,dim_nslope,dim_time/),varid),'defining variable h2o_ads_reg')
374        call check_nc(nf90_put_att(ncid,varid,'title','Mass of H2O adsorbded in the regolith'),'putting title attribute for h2o_ads_reg')
375
376        call check_nc(nf90_def_var(ncid,'co2_ads_reg',nf90_double,(/dim_ngrid,dim_nsoil,dim_nslope,dim_time/),varid),'defining variable co2_ads_reg')
377        call check_nc(nf90_put_att(ncid,varid,'title','Mass of CO2 adsorbded in the regolith'),'putting title attribute for co2_ads_reg')
378    end if
379end if
380
381if (do_layering) then
382    call check_nc(nf90_def_var(ncid,'stratif_top_elevation',nf90_double,(/dim_ngrid,dim_nb_str,dim_nslope,dim_time/),varid),'defining variable stratif_top_elevation')
383    call check_nc(nf90_put_att(ncid,varid,'title','Layering top elevation'),'putting title attribute for stratif_top_elevation')
384
385    call check_nc(nf90_def_var(ncid,'stratif_h_h2oice',nf90_double,(/dim_ngrid,dim_nb_str,dim_nslope,dim_time/),varid),'defining variable stratif_h_h2oice')
386    call check_nc(nf90_put_att(ncid,varid,'title','Layering H2O ice height'),'putting title attribute for stratif_h_h2oice')
387
388    call check_nc(nf90_def_var(ncid,'stratif_h_co2ice',nf90_double,(/dim_ngrid,dim_nb_str,dim_nslope,dim_time/),varid),'defining variable longitude')
389    call check_nc(nf90_put_att(ncid,varid,'title','Layering CO2 ice height'),'putting title attribute for stratif_h_co2ice')
390
391    call check_nc(nf90_def_var(ncid,'stratif_h_dust',nf90_double,(/dim_ngrid,dim_nb_str,dim_nslope,dim_time/),varid),'defining variable stratif_h_co2ice')
392    call check_nc(nf90_put_att(ncid,varid,'title','Layering dust height'),'putting title attribute for stratif_h_dust')
393
394    call check_nc(nf90_def_var(ncid,'stratif_h_pore',nf90_double,(/dim_ngrid,dim_nb_str,dim_nslope,dim_time/),varid),'defining variable stratif_h_pore')
395    call check_nc(nf90_put_att(ncid,varid,'title','Layering pore height'),'putting title attribute for stratif_h_pore')
396
397    call check_nc(nf90_def_var(ncid,'stratif_poreice_volfrac',nf90_double,(/dim_ngrid,dim_nb_str,dim_nslope,dim_time/),varid),'defining variable stratif_poreice_volfrac')
398    call check_nc(nf90_put_att(ncid,varid,'title','Layering ice pore volume fraction'),'putting title attribute for stratif_poreice_volfrac')
399end if
400! Global attributes
401call check_nc(nf90_put_att(ncid,nf90_global,'title','PEM start file'),'putting global attribute')
402
403! Leave define mode and putting variables defining dimensions
404call check_nc(nf90_enddef(ncid),'leaving define mode')
405call put_var_nc('longitude',longitudes)
406call put_var_nc('latitude',latitudes)
407call put_var_nc('cell_area',cell_area)
408call put_var_nc('soildepth',mlayer)
409call close_nc('re'//startpem_name)
410
411! File creation done
412is_restartpem = .true.
413
414END SUBROUTINE create_nc_pem
415!=======================================================================
416
417!=======================================================================
418SUBROUTINE write_restartpem(h2o_ice,co2_ice,tsoil,TI,icetable_depth,icetable_thickness,ice_porefilling,h2o_ads_reg,co2_ads_reg,layerings_map)
419!-----------------------------------------------------------------------
420! NAME
421!    write_restartpem
422!
423! DESCRIPTION
424!    Write the file "restartpem.nc".
425!
426! AUTHORS & DATE
427!    JB Clement, 01/2026
428!
429! NOTES
430!
431!-----------------------------------------------------------------------
432
433! DEPENDENCIES
434! ------------
435use geometry,         only: ngrid, nslope
436use evolution,        only: pem_ini_date, n_yr_sim
437use layered_deposits, only: layering, do_layering, map2array, get_nb_str_max
438use soil,             only: do_soil, inertiedat
439use sorption,         only: do_sorption
440use ice_table,        only: icetable_equilibrium, icetable_dynamic
441use stoppage,         only: stop_clean
442use display,          only: print_msg
443
444! DECLARATION
445! -----------
446implicit none
447
448! ARGUMENTS
449! ---------
450real(dp),       dimension(:,:),   intent(in) :: h2o_ice, co2_ice, icetable_depth, icetable_thickness
451real(dp),       dimension(:,:,:), intent(in) :: tsoil, TI, ice_porefilling, h2o_ads_reg, co2_ads_reg
452type(layering), dimension(:,:),   intent(in) :: layerings_map
453
454! LOCAL VARIABLES
455! ---------------
456integer(di)                               :: itime ! Time index to record variables
457integer(di)                               :: nb_str_max
458real(dp), dimension(:,:,:,:), allocatable :: layerings_array
459
460! CODE
461! ----
462! Create the file
463nb_str_max = get_nb_str_max(layerings_map)
464call create_nc_pem(nb_str_max)
465
466call print_msg('> Writing "re'//startpem_name//'"')
467if (.not. is_restartpem) call stop_clean(__FILE__,__LINE__,'The file"'//startpem_name//'" has not been created',1)
468
469! Writing time counter
470call open_nc('re'//startpem_name,'write',itime)
471call put_var_nc('Time',pem_ini_date + n_yr_sim,itime)
472
473! Writing other variables
474call put_var_nc('h2o_ice',h2o_ice,itime)
475call put_var_nc('co2_ice',co2_ice,itime)
476
477if (do_soil) then
478    call put_var_nc('tsoil',tsoil,itime)
479    call put_var_nc('TI',TI,itime)
480    call put_var_nc('inertiedat',inertiedat,itime)
481    call put_var_nc('icetable_depth',icetable_depth,itime)
482    if (icetable_equilibrium) then
483        call put_var_nc('icetable_thickness',icetable_thickness,itime)
484    else if (icetable_dynamic) then
485        call put_var_nc('ice_porefilling',ice_porefilling,itime)
486    end if
487    if (do_sorption) then
488        call put_var_nc('h2o_ads_reg',h2o_ads_reg,itime)
489        call put_var_nc('co2_ads_reg',co2_ads_reg,itime)
490    end if
491end if
492
493if (do_layering) then
494    allocate(layerings_array(ngrid,nslope,nb_str_max,6))
495    call map2array(layerings_map,layerings_array)
496    call put_var_nc('stratif_top_elevation',layerings_array(:,:,:,1),itime)
497    call put_var_nc('stratif_h_co2ice',layerings_array(:,:,:,2),itime)
498    call put_var_nc('stratif_h_h2oice',layerings_array(:,:,:,3),itime)
499    call put_var_nc('stratif_h_dust',layerings_array(:,:,:,4),itime)
500    call put_var_nc('stratif_h_pore',layerings_array(:,:,:,5),itime)
501    call put_var_nc('stratif_poreice_volfrac',layerings_array(:,:,:,6),itime)
502    deallocate(layerings_array)
503end if
504
505! Close
506call close_nc('re'//startpem_name)
507
508END SUBROUTINE write_restartpem
509!=======================================================================
510
511END MODULE clim_state_rec
Note: See TracBrowser for help on using the repository browser.