source: trunk/LMDZ.COMMON/libf/evolution/clim_state_rec.F90

Last change on this file was 4170, checked in by jbclement, 10 days ago

PEM:

  • Deletion of 'flux_ssice' ('zqdsdif_tot') from the "startfi.nc" as it is not needed.
  • Using the yearly average flux for the sublimating subsurface ice instead of the value at last timestep.
  • Keeping a clear separation between the subsurface ice flux and the surface ice tendency.
  • Making sure that subsurface ice depth is well given to the PCM at the end of the PEM (ice table VS layering).

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, startevo_name
20
21! DECLARATION
22! -----------
23implicit none
24
25! VARIABLES
26! ---------
27logical(k4), protected, private :: is_restartevo = .false. ! Flag to know if "restartevo.nc" exists
28
29contains
30!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
31
32!=======================================================================
33SUBROUTINE write_restart(ps4PCM,pa4PCM,preff4PCM,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, LVL_NFO
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) :: 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,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//'"',LVL_NFO)
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(ps4PCM,ps_dyn)
97call put_var_nc('ps',ps_dyn,1)
98
99! Potential temperature
100do l = 1,nlayer
101    call vect2dyngrd(teta4PCM(:,l),var_dyn(:,:,l))
102end do
103call put_var_nc('teta',var_dyn,1)
104
105! Air mass
106do l = 1,nlayer
107    call vect2dyngrd(air_mass4PCM(:,l),var_dyn(:,:,l),extensive = .true.)
108end do
109call put_var_nc('masse',var_dyn,1)
110
111! Tracers
112do i = 1,nq
113    do l = 1,nlayer
114        call vect2dyngrd(q4PCM(:,l,i),var_dyn(:,:,l))
115    end do
116    call put_var_nc(qnames(i),var_dyn,1)
117end do
118
119! Close
120call close_nc('re'//start_name)
121
122END SUBROUTINE write_restart
123!=======================================================================
124
125!=======================================================================
126SUBROUTINE write_restart1D(ps4PCM,pa4PCM,preff4PCM,q4PCM)
127!-----------------------------------------------------------------------
128! NAME
129!    write_restart1D
130!
131! DESCRIPTION
132!    Write the file "restart1D.txt".
133!
134! AUTHORS & DATE
135!    JB Clement, 01/2026
136!
137! NOTES
138!
139!-----------------------------------------------------------------------
140
141! DEPENDENCIES
142! ------------
143use geometry,   only: nlayer
144use atmosphere, only: teta_PCM, u_PCM, v_PCM
145use tracers,    only: nq, qnames
146use stoppage,   only: stop_clean
147use display,    only: print_msg, LVL_NFO
148
149! DECLARATION
150! -----------
151implicit none
152
153! ARGUMENTS
154! ---------
155real(dp),                   intent(in) :: pa4PCM, preff4PCM
156real(dp), dimension(:),     intent(in) :: ps4PCM
157real(dp), dimension(:,:,:), intent(in) :: q4PCM
158
159! LOCAL VARIABLES
160! ---------------
161integer(di) :: funit, ierr, i, l
162
163! CODE
164! ----
165! Write "restart1D.txt"
166call print_msg('> Writing "re'//start1D_name//'"',LVL_NFO)
167open(newunit = funit,file = 're'//start1D_name,status = "replace",action = "write",iostat = ierr)
168if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening file "re'//start1D_name//'"!',ierr)
169write(funit,*) 'ps', ps4PCM(1), pa4PCM, preff4PCM
170do i = 1,nq
171    write(funit,*) qnames(i), (q4PCM(1,l,i), l = 1,nlayer)
172end do
173write(funit,*) 'u', (u_PCM(1,l), l = 1,nlayer)
174write(funit,*) 'v', (v_PCM(1,l), l = 1,nlayer)
175write(funit,*) 'teta', (teta_PCM(1,l), l = 1,nlayer)
176close(funit)
177
178END SUBROUTINE write_restart1D
179!=======================================================================
180
181!=======================================================================
182SUBROUTINE write_restartfi(is_h2o_perice,h2o_ice4PCM,co2_ice4PCM,tsurf4PCM,tsoil4PCM,inertiesoil4PCM,albedo4PCM,emissivity4PCM,flux_geo4PCM,h2oice_depth4PCM)
183!-----------------------------------------------------------------------
184! NAME
185!    write_restartfi
186!
187! DESCRIPTION
188!    Write the file "restartfi.nc".
189!
190! AUTHORS & DATE
191!    JB Clement, 01/2026
192!
193! NOTES
194!
195!-----------------------------------------------------------------------
196
197! DEPENDENCIES
198! ------------
199use orbit,    only: obliquity, aphelion, perihelion, date_peri
200use frost,    only: h2o_frost4PCM, co2_frost4PCM
201use stoppage, only: stop_clean
202use display,  only: print_msg, LVL_NFO
203
204! DECLARATION
205! -----------
206implicit none
207
208! ARGUMENTS
209! ---------
210real(dp),    dimension(:,:),   intent(in) :: h2o_ice4PCM, co2_ice4PCM, albedo4PCM, emissivity4PCM, tsurf4PCM, flux_geo4PCM, h2oice_depth4PCM
211real(dp),    dimension(:,:,:), intent(in) :: tsoil4PCM, inertiesoil4PCM
212logical(k4), dimension(:),     intent(in) :: is_h2o_perice
213
214! LOCAL VARIABLES
215! ---------------
216integer(di)                         :: cstat
217integer(di)                         :: nindex ! Size of dimension 'index'
218real(dp), dimension(:), allocatable :: controle
219
220! CODE
221! ----
222! Copy "startfi.nc" into "restartfi.nc"
223call print_msg('> Writing "re'//startfi_name//'"',LVL_NFO)
224call execute_command_line('cp '//startfi_name//' re'//startfi_name,cmdstat = cstat)
225if (cstat > 0) then
226    call stop_clean(__FILE__,__LINE__,'command execution failed!',1)
227else if (cstat < 0) then
228    call stop_clean(__FILE__,__LINE__,'command execution not supported!',1)
229end if
230
231! Load the variable 'controle' to modify it with new values
232call open_nc('re'//startfi_name,'read')
233call get_dim_nc('index',nindex)
234allocate(controle(nindex))
235call get_var_nc('controle',controle)
236call close_nc('re'//startfi_name)
237
238! Rewrite the variables modified by the PEM
239call open_nc('re'//startfi_name,'write')
240
241! Orbital parameters (controle)
242controle(18) = obliquity  ! Obliquity
243controle(15) = perihelion ! Perihelion
244controle(16) = aphelion   ! Aphelion
245controle(17) = date_peri  ! Date of perihelion
246call put_var_nc('controle',controle)
247deallocate(controle)
248
249! Variables that have been modified
250call put_var_nc('watercaptag',merge(1._dp,0._dp,is_h2o_perice))
251call put_var_nc('watercap',h2o_ice4PCM,1)
252call put_var_nc('h2o_ice',h2o_frost4PCM,1)
253call put_var_nc('co2',co2_frost4PCM,1)
254call put_var_nc('perennial_co2ice',co2_ice4PCM,1)
255call put_var_nc('tsurf',tsurf4PCM,1)
256call put_var_nc('tsoil',tsoil4PCM,1)
257call put_var_nc('inertiesoil',inertiesoil4PCM,1)
258call put_var_nc('albedo',albedo4PCM,1)
259call put_var_nc('emis',emissivity4PCM,1)
260call put_var_nc('flux_geo',flux_geo4PCM,1)
261call put_var_nc('h2oice_depth',h2oice_depth4PCM)
262
263! Close
264call close_nc('re'//startfi_name)
265
266END SUBROUTINE write_restartfi
267!=======================================================================
268
269!=======================================================================
270SUBROUTINE create_startevo(nb_str_max)
271!-----------------------------------------------------------------------
272! NAME
273!     create_startevo
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, LVL_NFO
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'//startevo_name//'"',LVL_NFO)
318call check_nc(nf90_create('re'//startevo_name,nf90_clobber,ncid),'creating re'//startevo_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),LVL_NFO)
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
401! Global attributes
402call check_nc(nf90_put_att(ncid,nf90_global,'title','PEM start file'),'putting global attribute')
403
404! Leave define mode and putting variables defining dimensions
405call check_nc(nf90_enddef(ncid),'leaving define mode')
406call put_var_nc('longitude',longitudes)
407call put_var_nc('latitude',latitudes)
408call put_var_nc('cell_area',cell_area)
409call put_var_nc('soildepth',mlayer)
410call close_nc('re'//startevo_name)
411
412! File creation done
413is_restartevo = .true.
414
415END SUBROUTINE create_startevo
416!=======================================================================
417
418!=======================================================================
419SUBROUTINE write_restartevo(h2o_ice,co2_ice,tsoil,TI,icetable_depth,icetable_thickness,ice_porefilling,h2o_ads_reg,co2_ads_reg,layerings_map)
420!-----------------------------------------------------------------------
421! NAME
422!    write_restartevo
423!
424! DESCRIPTION
425!    Write the file "restartevo.nc".
426!
427! AUTHORS & DATE
428!    JB Clement, 01/2026
429!
430! NOTES
431!
432!-----------------------------------------------------------------------
433
434! DEPENDENCIES
435! ------------
436use geometry,         only: ngrid, nslope
437use evolution,        only: pem_ini_date, n_yr_sim
438use layered_deposits, only: layering, do_layering, map2array, get_nb_str_max
439use soil,             only: do_soil, inertiedat
440use sorption,         only: do_sorption
441use ice_table,        only: icetable_equilibrium, icetable_dynamic
442use stoppage,         only: stop_clean
443use display,          only: print_msg, LVL_NFO
444
445! DECLARATION
446! -----------
447implicit none
448
449! ARGUMENTS
450! ---------
451real(dp),       dimension(:,:),   intent(in) :: h2o_ice, co2_ice, icetable_depth, icetable_thickness
452real(dp),       dimension(:,:,:), intent(in) :: tsoil, TI, ice_porefilling, h2o_ads_reg, co2_ads_reg
453type(layering), dimension(:,:),   intent(in) :: layerings_map
454
455! LOCAL VARIABLES
456! ---------------
457integer(di)                               :: itime ! Time index to record variables
458integer(di)                               :: nb_str_max
459real(dp), dimension(:,:,:,:), allocatable :: layerings_array
460
461! CODE
462! ----
463! Create the file
464nb_str_max = get_nb_str_max(layerings_map)
465call create_startevo(nb_str_max)
466
467call print_msg('> Writing "re'//startevo_name//'"',LVL_NFO)
468if (.not. is_restartevo) call stop_clean(__FILE__,__LINE__,'The file"'//startevo_name//'" has not been created',1)
469
470! Writing time counter
471call open_nc('re'//startevo_name,'write',itime)
472call put_var_nc('Time',pem_ini_date + n_yr_sim,itime)
473
474! Writing other variables
475call put_var_nc('h2o_ice',h2o_ice,itime)
476call put_var_nc('co2_ice',co2_ice,itime)
477
478if (do_soil) then
479    call put_var_nc('tsoil',tsoil,itime)
480    call put_var_nc('TI',TI,itime)
481    call put_var_nc('inertiedat',inertiedat,itime)
482    call put_var_nc('icetable_depth',icetable_depth,itime)
483    if (icetable_equilibrium) then
484        call put_var_nc('icetable_thickness',icetable_thickness,itime)
485    else if (icetable_dynamic) then
486        call put_var_nc('ice_porefilling',ice_porefilling,itime)
487    end if
488    if (do_sorption) then
489        call put_var_nc('h2o_ads_reg',h2o_ads_reg,itime)
490        call put_var_nc('co2_ads_reg',co2_ads_reg,itime)
491    end if
492end if
493
494if (do_layering) then
495    allocate(layerings_array(ngrid,nslope,nb_str_max,6))
496    call map2array(layerings_map,layerings_array)
497    call put_var_nc('stratif_top_elevation',layerings_array(:,:,:,1),itime)
498    call put_var_nc('stratif_h_co2ice',layerings_array(:,:,:,2),itime)
499    call put_var_nc('stratif_h_h2oice',layerings_array(:,:,:,3),itime)
500    call put_var_nc('stratif_h_dust',layerings_array(:,:,:,4),itime)
501    call put_var_nc('stratif_h_pore',layerings_array(:,:,:,5),itime)
502    call put_var_nc('stratif_poreice_volfrac',layerings_array(:,:,:,6),itime)
503    deallocate(layerings_array)
504end if
505
506! Close
507call close_nc('re'//startevo_name)
508
509END SUBROUTINE write_restartevo
510!=======================================================================
511
512END MODULE clim_state_rec
Note: See TracBrowser for help on using the repository browser.