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

Last change on this file was 4090, checked in by jbclement, 30 hours ago

PEM:
Update README and change "startpem.nc" into "startevol.nc".
JBC

File size: 18.9 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 "restartevol.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
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//'"')
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,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
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) :: q4PCM
154
155! LOCAL VARIABLES
156! ---------------
157integer(di) :: funit, ierr, i, l
158
159! CODE
160! ----
161! Write "restart1D.txt"
162call print_msg('> Writing "re'//start1D_name//'"')
163open(newunit = funit,file = 're'//start1D_name,status = "replace",action = "write",iostat = ierr)
164if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening file "re'//start1D_name//'"!',ierr)
165write(funit,*) 'ps', ps4PCM(1), pa4PCM, preff4PCM
166do i = 1,nq
167    write(funit,*) qnames(i), (q4PCM(1,l,i), l = 1,nlayer)
168end do
169write(funit,*) 'u', (u_PCM(1,l), l = 1,nlayer)
170write(funit,*) 'v', (v_PCM(1,l), l = 1,nlayer)
171write(funit,*) 'teta', (teta_PCM(1,l), l = 1,nlayer)
172close(funit)
173
174END SUBROUTINE write_restart1D
175!=======================================================================
176
177!=======================================================================
178SUBROUTINE write_restartfi(is_h2o_perice,h2o_ice4PCM,co2_ice4PCM,tsurf4PCM,tsoil4PCM,inertiesoil4PCM,albedo4PCM,emissivity4PCM,flux_geo4PCM)
179!-----------------------------------------------------------------------
180! NAME
181!    write_restartfi
182!
183! DESCRIPTION
184!    Write the file "restartfi.nc".
185!
186! AUTHORS & DATE
187!    JB Clement, 01/2026
188!
189! NOTES
190!
191!-----------------------------------------------------------------------
192
193! DEPENDENCIES
194! ------------
195use orbit,    only: obliquity, aphelion, perihelion, date_peri
196use frost,    only: h2o_frost4PCM, co2_frost4PCM
197use stoppage, only: stop_clean
198use display,  only: print_msg
199
200! DECLARATION
201! -----------
202implicit none
203
204! ARGUMENTS
205! ---------
206real(dp),    dimension(:,:),   intent(in) :: h2o_ice4PCM, co2_ice4PCM, albedo4PCM, emissivity4PCM, tsurf4PCM, flux_geo4PCM
207real(dp),    dimension(:,:,:), intent(in) :: tsoil4PCM, inertiesoil4PCM
208logical(k4), dimension(:),     intent(in) :: is_h2o_perice
209
210! LOCAL VARIABLES
211! ---------------
212integer(di)                         :: cstat
213integer(di)                         :: nindex ! Size of dimension 'index'
214real(dp), dimension(:), allocatable :: controle
215
216! CODE
217! ----
218! Copy "startfi.nc" into "restartfi.nc"
219call print_msg('> Writing "re'//startfi_name//'"')
220call execute_command_line('cp '//startfi_name//' re'//startfi_name,cmdstat = cstat)
221if (cstat > 0) then
222    call stop_clean(__FILE__,__LINE__,'command execution failed!',1)
223else if (cstat < 0) then
224    call stop_clean(__FILE__,__LINE__,'command execution not supported!',1)
225end if
226
227! Load the variable 'controle' to modify it with new values
228call open_nc('re'//startfi_name,'read')
229call get_dim_nc('index',nindex)
230allocate(controle(nindex))
231call get_var_nc('controle',controle)
232call close_nc('re'//startfi_name)
233
234! Rewrite the variables modified by the PEM
235call open_nc('re'//startfi_name,'write')
236
237! Orbital parameters (controle)
238controle(18) = obliquity ! Obliquity
239controle(15) = perihelion ! Perihelion
240controle(16) = aphelion ! Aphelion
241controle(17) = date_peri ! Date of perihelion
242call put_var_nc('controle',controle)
243deallocate(controle)
244
245! Variables that have been modified
246call put_var_nc('watercaptag',merge(1._dp,0._dp,is_h2o_perice))
247call put_var_nc('watercap',h2o_ice4PCM,1)
248call put_var_nc('h2o_ice',h2o_frost4PCM,1)
249call put_var_nc('co2',co2_frost4PCM,1)
250call put_var_nc('perennial_co2ice',co2_ice4PCM,1)
251call put_var_nc('tsurf',tsurf4PCM,1)
252call put_var_nc('tsoil',tsoil4PCM,1)
253call put_var_nc('inertiesoil',inertiesoil4PCM,1)
254call put_var_nc('albedo',albedo4PCM,1)
255call put_var_nc('emis',emissivity4PCM,1)
256call put_var_nc('flux_geo',flux_geo4PCM,1)
257
258! h2oice_depth
259! lag_co2_ice
260! zdqsdif_ssi_tot
261! d_coef
262
263call close_nc('re'//startfi_name)
264
265END SUBROUTINE write_restartfi
266!=======================================================================
267
268!=======================================================================
269SUBROUTINE create_nc_pem(nb_str_max)
270!-----------------------------------------------------------------------
271! NAME
272!     create_nc_pem
273!
274! DESCRIPTION
275!     Create a NetCDF file for the PEM.
276!
277! AUTHORS & DATE
278!     JB Clement, 01/2026
279!
280! NOTES
281!
282!-----------------------------------------------------------------------
283
284! DEPENDENCIES
285! ------------
286use netcdf,           only: nf90_double, nf90_clobber, nf90_unlimited, nf90_global, &
287                            nf90_create, nf90_enddef, nf90_def_dim, nf90_def_var, nf90_put_att
288use geometry,         only: dim_init, ngrid, nsoil, nslope, longitudes, latitudes, cell_area
289use stoppage,         only: stop_clean
290use soil,             only: do_soil, mlayer
291use sorption,         only: do_sorption
292use layered_deposits, only: do_layering
293use display,          only: print_msg
294use utility,          only: int2str
295
296! DECLARATION
297! -----------
298implicit none
299
300! ARGUMENTS
301! ---------
302integer(di), intent(in) :: nb_str_max
303
304! LOCAL VARIABLES
305! ---------------
306integer(di) :: ncid  ! File ID
307integer(di) :: varid ! Variable ID
308integer(di) :: dim_ngrid, dim_nsoil, dim_nslope, dim_time, dim_nb_str
309
310! CODE
311! ----
312! Check if dimensions are well initialized
313if (.not. dim_init) call stop_clean(__FILE__,__LINE__,'dimensions are not initilized in the PEM!',1)
314
315! Create file
316call print_msg('> Creating "re'//startpem_name//'"')
317call check_nc(nf90_create('re'//startpem_name,nf90_clobber,ncid),'creating re'//startpem_name) ! Enter define mode
318
319! Define dimensions
320call check_nc(nf90_def_dim(ncid,'Time',nf90_unlimited,dim_time),'defining dimension Time')
321call check_nc(nf90_def_dim(ncid,'physical_points',ngrid,dim_ngrid),'defining dimension physical_points')
322call check_nc(nf90_def_dim(ncid,'subsurface_layers',nsoil,dim_nsoil),'defining dimension subsurface_layers')
323call check_nc(nf90_def_dim(ncid,'nslope',nslope,dim_nslope),'defining dimension nslope')
324if (do_layering) then
325    call print_msg('nb_str_max = '//int2str(nb_str_max))
326    call check_nc(nf90_def_dim(ncid,'nb_str_max',nb_str_max,dim_nb_str),'defining dimension nb_str_max')
327end if
328
329! Define variables
330call check_nc(nf90_def_var(ncid,'Time',nf90_double,(/dim_time/),varid),'defining variable Time')
331call check_nc(nf90_put_att(ncid,varid,'title','Year of simulation'),'putting title attribute for Time')
332call check_nc(nf90_put_att(ncid,varid,'units','Planetary year'),'putting units attribute for Time')
333
334call check_nc(nf90_def_var(ncid,'longitude',nf90_double,(/dim_ngrid/),varid),'defining variable longitude')
335call check_nc(nf90_put_att(ncid,varid,'title','Longitudes of the grid'),'putting title attribute for longitude')
336
337call check_nc(nf90_def_var(ncid,'latitude',nf90_double,(/dim_ngrid/),varid),'defining variable latitude')
338call check_nc(nf90_put_att(ncid,varid,'title','Latitudes of the grid'),'putting title attribute for latitude')
339
340call check_nc(nf90_def_var(ncid,'cell_area',nf90_double,(/dim_ngrid/),varid),'defining variable cell_area')
341call check_nc(nf90_put_att(ncid,varid,'title','Cell area'),'putting title attribute for cell_area')
342
343call check_nc(nf90_def_var(ncid,'h2o_ice',nf90_double,(/dim_ngrid,dim_nslope,dim_time/),varid),'defining variable h2o_ice')
344call check_nc(nf90_put_att(ncid,varid,'title','H2O ice'),'putting title attribute for h2o_ice')
345
346call check_nc(nf90_def_var(ncid,'co2_ice',nf90_double,(/dim_ngrid,dim_nslope,dim_time/),varid),'defining variable co2_ice')
347call check_nc(nf90_put_att(ncid,varid,'title','CO2 ice'),'putting title attribute for co2_ice')
348
349if (do_soil) then
350    call check_nc(nf90_def_var(ncid,'soildepth',nf90_double,(/dim_nsoil/),varid),'defining variable soildepth')
351    call check_nc(nf90_put_att(ncid,varid,'title','Depths of soil layers'),'putting title attribute for soildepth')
352
353    call check_nc(nf90_def_var(ncid,'tsoil',nf90_double,(/dim_ngrid,dim_nsoil,dim_nslope,dim_time/),varid),'defining variable tsoil')
354    call check_nc(nf90_put_att(ncid,varid,'title','Soil temperature'),'putting title attribute for tsoil')
355
356    call check_nc(nf90_def_var(ncid,'TI',nf90_double,(/dim_ngrid,dim_nsoil,dim_nslope,dim_time/),varid),'defining variable inertiedat')
357    call check_nc(nf90_put_att(ncid,varid,'title','Thermal inertie of PEM'),'putting title attribute for TI')
358
359    call check_nc(nf90_def_var(ncid,'inertiedat',nf90_double,(/dim_ngrid,dim_nsoil,dim_time/),varid),'defining variable longitude')
360    call check_nc(nf90_put_att(ncid,varid,'title','Soil thermal inertia'),'putting title attribute for inertiedat')
361
362    call check_nc(nf90_def_var(ncid,'icetable_depth',nf90_double,(/dim_ngrid,dim_nslope,dim_time/),varid),'defining variable icetable_depth')
363    call check_nc(nf90_put_att(ncid,varid,'title','Depth of ice table'),'putting title attribute for icetable_depth')
364
365    call check_nc(nf90_def_var(ncid,'icetable_thickness',nf90_double,(/dim_ngrid,dim_nslope,dim_time/),varid),'defining variable icetable_thickness')
366    call check_nc(nf90_put_att(ncid,varid,'title','Thickness of ice table'),'putting title attribute for icetable_thickness')
367
368    call check_nc(nf90_def_var(ncid,'ice_porefilling',nf90_double,(/dim_ngrid,dim_nsoil,dim_nslope,dim_time/),varid),'defining variable ice_porefilling')
369    call check_nc(nf90_put_att(ncid,varid,'title','Subsurface ice pore filling'),'putting title attribute for ice_porefilling')
370
371    if (do_sorption) then
372        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')
373        call check_nc(nf90_put_att(ncid,varid,'title','Mass of H2O adsorbded in the regolith'),'putting title attribute for h2o_ads_reg')
374
375        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')
376        call check_nc(nf90_put_att(ncid,varid,'title','Mass of CO2 adsorbded in the regolith'),'putting title attribute for co2_ads_reg')
377    end if
378end if
379
380if (do_layering) then
381    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')
382    call check_nc(nf90_put_att(ncid,varid,'title','Layering top elevation'),'putting title attribute for stratif_top_elevation')
383
384    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')
385    call check_nc(nf90_put_att(ncid,varid,'title','Layering H2O ice height'),'putting title attribute for stratif_h_h2oice')
386
387    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')
388    call check_nc(nf90_put_att(ncid,varid,'title','Layering CO2 ice height'),'putting title attribute for stratif_h_co2ice')
389
390    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')
391    call check_nc(nf90_put_att(ncid,varid,'title','Layering dust height'),'putting title attribute for stratif_h_dust')
392
393    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')
394    call check_nc(nf90_put_att(ncid,varid,'title','Layering pore height'),'putting title attribute for stratif_h_pore')
395
396    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')
397    call check_nc(nf90_put_att(ncid,varid,'title','Layering ice pore volume fraction'),'putting title attribute for stratif_poreice_volfrac')
398end if
399
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 "restartevol.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.