source: trunk/LMDZ.COMMON/libf/evolution/clim_state_init.F90 @ 4068

Last change on this file since 4068 was 4068, checked in by jbclement, 3 weeks ago

PEM:
Following r4065, some safeguards forgotten in "io_netcdf.F90" + 'tsurf' is not needed anymore in the "start1D.txt" + few forgotten small updates.
JBC

File size: 22.1 KB
Line 
1MODULE clim_state_init
2!-----------------------------------------------------------------------
3! NAME
4!     clim_state_init
5!
6! DESCRIPTION
7!     Read the start files to initialize 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, get_var_nc, get_dim_nc, put_var_nc, start_name, start1D_name, startfi_name, startpem_name
20
21! DECLARATION
22! -----------
23implicit none
24
25contains
26!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27
28!=======================================================================
29SUBROUTINE read_start()
30!-----------------------------------------------------------------------
31! NAME
32!     read_start
33!
34! DESCRIPTION
35!     Read the file "start.nc".
36!
37! AUTHORS & DATE
38!     JB Clement, 12/2025
39!
40! NOTES
41!
42!-----------------------------------------------------------------------
43
44! DEPENDENCIES
45! ------------
46use geometry,   only: nlayer, nlon, nlat, ngrid
47use atmosphere, only: set_ap, set_bp, set_ps_PCM, set_teta_PCM
48use tracers,    only: nq, qnames, set_q_PCM
49use stoppage,   only: stop_clean
50use display,    only: print_msg
51
52! DECLARATION
53! -----------
54implicit none
55
56! LOCAL VARIABLES
57! ---------------
58real(dp), dimension(nlayer + 1)           :: tmp1d
59real(dp), dimension(nlon + 1,nlat)        :: tmp2d
60real(dp), dimension(nlon + 1,nlat,nlayer) :: tmp3d
61logical(k4)                               :: here
62integer(di)                               :: i
63
64! CODE
65! ----
66! In case of 1D
67! ~~~~~~~~~~~~~
68if (ngrid == 1) then
69    call read_start1D()
70    return
71end if
72
73! In case of 3D
74! ~~~~~~~~~~~~~
75! Open
76inquire(file = start_name,exist = here)
77if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//start_name//'"!',1)
78call print_msg('> Reading "'//start_name//'"')
79call open_nc(start_name,'read')
80
81! Get the variables
82call get_var_nc('ap',tmp1d)
83call set_ap(tmp1d)
84
85call get_var_nc('bp',tmp1d)
86call set_bp(tmp1d)
87
88call get_var_nc('ps',tmp2d)
89call set_ps_PCM(tmp2d)
90
91call get_var_nc('teta',tmp3d)
92call set_teta_PCM(tmp3d)
93
94do i = 1,nq
95    call get_var_nc(trim(qnames(i)),tmp3d)
96    call set_q_PCM(tmp3d,i)
97end do
98
99! Close
100call close_nc(start_name)
101
102END SUBROUTINE read_start
103!=======================================================================
104
105!=======================================================================
106SUBROUTINE read_start1D()
107!-----------------------------------------------------------------------
108! NAME
109!     read_start1D
110!
111! DESCRIPTION
112!     Read the file "start1D.txt".
113!
114! AUTHORS & DATE
115!     JB Clement, 12/2025
116!
117! NOTES
118!
119!-----------------------------------------------------------------------
120
121! DEPENDENCIES
122! ------------
123use geometry,   only: nlayer, nslope
124use atmosphere, only: set_ap, set_bp, set_ps_PCM, set_teta_PCM, set_u_PCM, set_v_PCM, compute_alt_coord, set_ap, set_bp
125use tracers,    only: nq, set_q_PCM
126use stoppage,   only: stop_clean
127use config,     only: read_callphys
128use display,    only: print_msg
129
130! DECLARATION
131! -----------
132implicit none
133
134! LOCAL VARIABLES
135! ---------------
136integer(di)                     :: i, j, k, ierr, funit
137character(30)                   :: header
138real(dp), dimension(1,1)        :: tmp
139real(dp), dimension(1,1,nlayer) :: q_tmp, teta_tmp, wind_tmp
140real(dp), dimension(nlayer + 1) :: ap_tmp, bp_tmp
141real(dp)                        :: pa, preff
142logical(k4)                     :: here
143
144! CODE
145! ----
146! Open "start1D.txt"
147inquire(file = start1D_name,exist = here)
148if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//start1D_name//'"!',1)
149call print_msg('> Reading "'//start1D_name//'"')
150open(newunit = funit,file = start1D_name,status = "old",action = "read",iostat = ierr)
151if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening file "'//start1D_name//'"!',ierr)
152
153! Get the variables
154read(funit,*) header, tmp, pa, preff
155call set_ps_PCM(tmp)
156
157do i = 1,nq
158    read(funit,*,iostat = ierr) header, (q_tmp(1,1,k),k = 1,nlayer)
159    if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'not enough atmospheric layers defined in the file "'//start1D_name//'" for the tracer "'//trim(header)//'"!',1)
160    call set_q_PCM(q_tmp,i)
161end do
162read(funit,*,iostat = ierr) header, (wind_tmp(1,1,k),k = 1,nlayer)
163call set_u_PCM(wind_tmp)
164read(funit,*,iostat = ierr) header, (wind_tmp(1,1,k),k = 1,nlayer)
165call set_v_PCM(wind_tmp)
166read(funit,*,iostat = ierr) header, (teta_tmp(1,1,k),k = 1,nlayer)
167call set_teta_PCM(teta_tmp)
168
169! Close
170close(funit)
171
172! Construct altitude coordinates (not stored in "start1D.txt")
173call compute_alt_coord(pa,preff,ap_tmp,bp_tmp)
174call set_ap(ap_tmp)
175call set_bp(bp_tmp)
176
177! Read the "callphys.def"
178call read_callphys() ! To get 'CO2cond_ps'
179
180END SUBROUTINE read_start1D
181!=======================================================================
182
183!=======================================================================
184SUBROUTINE read_startfi()
185!-----------------------------------------------------------------------
186! NAME
187!     read_startfi
188!
189! DESCRIPTION
190!     Read the file "startfi.nc".
191!
192! AUTHORS & DATE
193!     JB Clement, 12/2025
194!
195! NOTES
196!
197!-----------------------------------------------------------------------
198
199! DEPENDENCIES
200! ------------
201use geometry,  only: ngrid, nslope, nsoil_PCM
202use stoppage,  only: stop_clean
203use config,    only: read_controldata
204use slopes,    only: set_def_slope_mean, set_subslope_dist, set_iflat
205use surface,   only: set_albedodat_PCM, set_albedo_PCM, set_emissivity_PCM
206use surf_temp, only: set_tsurf_PCM
207use soil_temp, only: set_tsoil_PCM, set_flux_geo_PCM
208use frost,     only: set_h2ofrost_PCM, set_co2frost_PCM
209use surf_ice,  only: set_is_h2o_perice_PCM, set_co2_perice_PCM
210use soil,      only: set_TI_PCM, set_inertiedat_PCM
211use display,   only: print_msg
212
213! DECLARATION
214! -----------
215implicit none
216
217! LOCAL VARIABLES
218! ---------------
219real(dp), dimension(:),   allocatable       :: tmp1d
220real(dp), dimension(:,:), allocatable       :: tmp2d
221real(dp), dimension(ngrid,nsoil_PCM,nslope) :: tmp3d
222logical(k4)                                 :: here
223
224! CODE
225! ----
226inquire(file = startfi_name,exist = here)
227if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//startfi_name//'"!',1)
228
229! Allocate the array to store the variables
230call print_msg('> Reading "'//startfi_name//'"')
231allocate(tmp1d(nslope + 1),tmp2d(ngrid,nslope))
232
233! Get control data
234call read_controldata()
235
236! Open
237call open_nc(startfi_name,'read')
238
239! Get the variables
240if (nslope > 1) then
241    call get_var_nc('def_slope_mean',tmp1d)
242    call set_def_slope_mean(tmp1d)
243
244    call get_var_nc('subslope_dist',tmp2d)
245    call set_subslope_dist(tmp2d)
246end if
247
248call get_var_nc('flux_geo',tmp2d)
249call set_flux_geo_PCM(tmp2d)
250
251call get_var_nc('h2o_ice',tmp2d)
252where (tmp2d < 0.) tmp2d = 0._dp ! Remove unphysical values of surface tracer
253call set_h2ofrost_PCM(tmp2d)
254
255call get_var_nc('co2',tmp2d)
256where (tmp2d < 0.) tmp2d = 0._dp ! Remove unphysical values of surface tracer
257call set_co2frost_PCM(tmp2d)
258
259call get_var_nc('perennial_co2ice',tmp2d)
260where (tmp2d < 0.) tmp2d = 0._dp ! Remove unphysical values of surface tracer
261call set_co2_perice_PCM(tmp2d)
262
263deallocate(tmp1d); allocate(tmp1d(ngrid))
264call get_var_nc('watercaptag',tmp1d)
265call set_is_h2o_perice_PCM(tmp1d)
266
267call get_var_nc('albedodat',tmp1d)
268call set_albedodat_PCM(tmp1d)
269
270call get_var_nc('albedo',tmp2d)
271call set_albedo_PCM(tmp2d)
272
273call get_var_nc('emis',tmp2d)
274call set_emissivity_PCM(tmp2d)
275
276call get_var_nc('tsurf',tmp2d)
277call set_tsurf_PCM(tmp2d)
278
279call get_var_nc('tsoil',tmp3d)
280call set_tsoil_PCM(tmp3d)
281
282deallocate(tmp2d); allocate(tmp2d(ngrid,nsoil_PCM))
283call get_var_nc('inertiedat',tmp2d)
284call set_inertiedat_PCM(tmp2d)
285
286call get_var_nc('inertiesoil',tmp3d)
287call set_TI_PCM(tmp3d)
288
289! To do?
290!   h2oice_depth
291!   d_coef
292!   lag_co2_ice
293
294! Close/Deallocate
295call close_nc(startfi_name)
296deallocate(tmp1d,tmp2d)
297
298END SUBROUTINE read_startfi
299!=======================================================================
300
301!=======================================================================
302SUBROUTINE read_startpem(tsurf_avg_yr1,tsurf_avg_yr2,ps_avg_glob,ps_ts,q_co2_ts,q_h2o_ts,h2o_surfdensity_avg,d_h2oice,d_co2ice,h2o_ice,co2_ice, &
303                         tsoil_avg,h2o_soildensity_avg,icetable_depth,icetable_thickness,ice_porefilling,layerings_map,                         &
304                         h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
305!-----------------------------------------------------------------------
306! NAME
307!     read_startpem
308!
309! DESCRIPTION
310!     Read the file "startpem.nc" which stores the PEM state.
311!
312! AUTHORS & DATE
313!     L. Lange, 2023
314!     JB Clement, 2023-2026
315!
316! NOTES
317!
318!-----------------------------------------------------------------------
319
320! DEPENDENCIES
321! ------------
322use stoppage,           only: stop_clean
323use geometry,           only: ngrid, nslope, nsoil, nsoil_PCM
324use evolution,          only: dt
325use physics,            only: mugaz, r
326use planet,             only: TI_breccia, TI_bedrock, alpha_clap_h2o, beta_clap_h2o
327use frost,              only: h2ofrost_PCM, h2o_frost4PCM, co2frost_PCM, co2_frost4PCM
328use surf_ice,           only: is_h2o_perice_PCM, h2oice_huge_ini, co2_perice_PCM
329use soil,               only: do_soil, TI, depth_breccia, depth_bedrock, index_breccia, index_bedrock, inertiedat, layer, inertiedat_PCM
330use soil_temp,          only: ini_tsoil_pem, compute_tsoil
331use soil_therm_inertia, only: update_soil_TI
332use ice_table,          only: icetable_equilibrium, icetable_dynamic, computeice_table_equilibrium
333use sorption,           only: do_sorption, evolve_regolith_adsorption
334use tracers,            only: mmol, iPCM_qh2o
335use layered_deposits,   only: layering, do_layering, array2map, ini_layering, add_stratum
336use surf_ice,           only: rho_co2ice, rho_h2oice
337use display,            only: print_msg
338use utility,            only: int2str
339
340! DECLARATION
341! -----------
342implicit none
343
344! ARGUMENTS
345! ---------
346real(dp), dimension(:,:),       intent(in)    :: tsurf_avg_yr1, tsurf_avg_yr2 ! Surface temperature at the last but one PCM run and at the last PCM run
347real(dp),                       intent(in)    :: ps_avg_glob                  ! Mean average pressure on the planet
348real(dp), dimension(:,:),       intent(in)    :: ps_ts                        ! Surface pressure timeseries
349real(dp), dimension(:,:),       intent(in)    :: q_co2_ts, q_h2o_ts           ! MMR of CO2 and H2O
350real(dp), dimension(:,:),       intent(in)    :: h2o_surfdensity_avg          ! Average of surface water density
351real(dp), dimension(:,:),       intent(in)    :: d_h2oice, d_co2ice           ! Tendencies on ice
352real(dp), dimension(:,:),       intent(inout) :: h2o_ice, co2_ice             ! Surface ice
353real(dp), dimension(:,:,:),     intent(inout) :: tsoil_avg                    ! Soil temperature
354real(dp), dimension(:,:,:),     intent(inout) :: h2o_soildensity_avg          ! Average of soil water soil density
355real(dp), dimension(:,:),       intent(inout) :: icetable_depth               ! Ice table depth
356real(dp), dimension(:,:),       intent(inout) :: icetable_thickness           ! Ice table thickness
357real(dp), dimension(:,:,:),     intent(inout) :: ice_porefilling              ! Subsurface ice pore filling
358real(dp), dimension(:,:,:),     intent(inout) :: h2o_ads_reg, co2_ads_reg     ! Mass of H2O and CO2 adsorbed
359type(layering), dimension(:,:), intent(inout) :: layerings_map                ! Layerings
360real(dp), dimension(:),         intent(inout) :: delta_h2o_ads, delta_co2_ads ! Mass of H2O and CO2 exchanged due to adsorption/desorption
361
362! LOCAL VARIABLES
363! ---------------
364logical(k4)                               :: here
365integer(di)                               :: i, islope, k, nb_str_max, nsoil_startpem
366real(dp)                                  :: delta           ! Depth of the interface regolith/breccia, breccia/bedrock [m]
367real(dp), dimension(ngrid,nsoil,nslope)   :: TI_startpem     ! Soil thermal inertia saved in the startpem [SI]
368real(dp), dimension(ngrid,nsoil,nslope)   :: tsoil_startpem  ! Soil temperature saved in the startpem [K]
369real(dp), dimension(:,:,:,:), allocatable :: layerings_array ! Array for layerings
370
371! CODE
372! ----
373! Check if the file exists
374inquire(file = startpem_name,exist = here)
375
376! If the file is not here
377! ~~~~~~~~~~~~~~~~~~~~~~~
378if (.not. here) then
379    ! It is possibly because it is the very first PEM run so everything needs to be initalized
380    call print_msg('> The file "'//startpem_name//'" was not found (possibly because this is the first PEM run)')
381
382    ! H2O ice
383    call print_msg("'h2o_ice' is initialized with default value 'h2oice_huge_ini' where 'watercaptag' is true.")
384    h2o_ice(:,:) = 0._dp
385    do i = 1,ngrid
386        if (is_h2o_perice_PCM(i)) h2o_ice(i,:) = h2oice_huge_ini + h2ofrost_PCM(i,:) - h2o_frost4PCM(i,:)
387    end do
388
389    ! CO2 ice
390    call print_msg("'co2_ice' is initialized with 'perennial_co2ice' found in the PCM.")
391    co2_ice(:,:) = co2_perice_PCM(:,:) + co2frost_PCM(:,:) - co2_frost4PCM(:,:)
392
393    if (do_soil) then
394        ! Thermal inertia
395        do islope = 1,nslope
396            do i = 1,ngrid
397                if (TI(i,index_breccia,islope) < TI_breccia) then
398                    !!! transition
399                    delta = depth_breccia
400                    TI(i,index_breccia + 1,islope) = sqrt((layer(index_breccia + 1) - layer(index_breccia))/                       &
401                                                               (((delta - layer(index_breccia))/(TI(i,index_breccia,islope)**2)) + &
402                                                               ((layer(index_breccia + 1) - delta)/(TI_breccia**2))))
403                    do k = index_breccia + 2,index_bedrock
404                        TI(i,k,islope) = TI_breccia
405                    end do
406                else ! we keep the high ti values
407                    do k = index_breccia + 1,index_bedrock
408                        TI(i,k,islope) = TI(i,index_breccia,islope)
409                    end do
410                end if
411                !! transition
412                delta = depth_bedrock
413                TI(i,index_bedrock + 1,islope) = sqrt((layer(index_bedrock + 1) - layer(index_bedrock))/                   &
414                                                       (((delta - layer(index_bedrock))/(TI(i,index_bedrock,islope)**2)) + &
415                                                       ((layer(index_bedrock + 1) - delta)/(TI_breccia**2))))
416                do k = index_bedrock + 2,nsoil
417                    TI(i,k,islope) = TI_bedrock
418                end do
419            end do
420        end do
421
422        do k = 1,nsoil_PCM
423            inertiedat(:,k) = inertiedat_PCM(:,k)
424        end do
425        !!! zone de transition
426        delta = depth_breccia
427        do i = 1,ngrid
428            if (inertiedat(i,index_breccia) < TI_breccia) then
429                inertiedat(i,index_breccia + 1) = sqrt((layer(index_breccia + 1) - layer(index_breccia))/                    &
430                                                        (((delta - layer(index_breccia))/(inertiedat(i,index_breccia)**2)) + &
431                                                        ((layer(index_breccia + 1) - delta)/(TI_breccia**2))))
432                do k = index_breccia + 2,index_bedrock
433                    inertiedat(i,k) = TI_breccia
434                end do
435            else
436                do k = index_breccia + 1,index_bedrock
437                    inertiedat(i,k) = inertiedat(i,index_breccia)
438                end do
439            end if
440        end do
441        !!! zone de transition
442        delta = depth_bedrock
443        do i = 1,ngrid
444            inertiedat(i,index_bedrock + 1) = sqrt((layer(index_bedrock + 1) - layer(index_bedrock))/                    &
445                                                    (((delta - layer(index_bedrock))/(inertiedat(i,index_bedrock)**2)) + &
446                                                    ((layer(index_bedrock + 1) - delta)/(TI_bedrock**2))))
447        end do
448
449        do k = index_bedrock + 2,nsoil
450            do i = 1,ngrid
451                inertiedat(i,k) = TI_bedrock
452            end do
453        end do
454
455        ! Soil temperature
456        do islope = 1,nslope
457            call ini_tsoil_pem(ngrid,nsoil,TI(:,:,islope),tsurf_avg_yr2(:,islope),tsoil_avg(:,:,islope))
458            call compute_tsoil(ngrid,nsoil,.true.,TI(:,:,islope),dt,tsurf_avg_yr2(:,islope),tsoil_avg(:,:,islope))
459
460            ! First raw initialization
461            h2o_soildensity_avg(:,nsoil_PCM + 1:nsoil,islope) = exp(beta_clap_h2o/tsoil_avg(:,nsoil_PCM + 1:nsoil,islope) + alpha_clap_h2o)/tsoil_avg(:,nsoil_PCM + 1:nsoil,islope)**mmol(iPCM_qh2o)/(mugaz*r)
462        end do
463
464        ! Ice table
465        if (icetable_equilibrium) then
466            call computeice_table_equilibrium(is_h2o_perice_PCM,h2o_surfdensity_avg,h2o_soildensity_avg,TI(:,1,:),icetable_depth,icetable_thickness)
467            call update_soil_TI(d_h2oice,h2o_ice,ps_avg_glob,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
468            do islope = 1,nslope
469                call ini_tsoil_pem(ngrid,nsoil,TI(:,:,islope),tsurf_avg_yr2(:,islope),tsoil_avg(:,:,islope))
470            end do
471        else if (icetable_dynamic) then
472            call update_soil_TI(d_h2oice,h2o_ice,ps_avg_glob,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
473            do islope = 1,nslope
474                call ini_tsoil_pem(ngrid,nsoil,TI(:,:,islope),tsurf_avg_yr2(:,islope),tsoil_avg(:,:,islope))
475            end do
476        end if
477
478        ! Absorption in regolith
479        if (do_sorption) then
480            h2o_ads_reg = 0._dp
481            co2_ads_reg = 0._dp
482            call evolve_regolith_adsorption(d_h2oice,d_co2ice,h2o_ice,co2_ice,tsoil_avg,TI,ps_ts,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
483            delta_co2_ads = 0._dp
484            delta_h2o_ads = 0._dp
485        end if ! do_sorption
486    end if ! do_soil
487
488    ! Layering
489    if (do_layering) then
490        call print_msg('layerings_map is initialized with sub-surface strata.')
491        call print_msg("Ice is added with 'h2oice_huge_ini' where 'watercaptag' is true and otherwise with 'perennial_co2ice' found in the PCM.")
492        do i = 1,ngrid
493            if (is_h2o_perice_PCM(i)) then
494                do islope = 1,nslope
495                    call ini_layering(layerings_map(i,islope))
496                    call add_stratum(layerings_map(i,islope),1.05_dp*h2oice_huge_ini/rho_h2oice,0._dp,h2oice_huge_ini/rho_h2oice,0.05_dp*h2oice_huge_ini/rho_h2oice,0._dp,0._dp)
497                end do
498            else
499                do islope = 1,nslope
500                    call ini_layering(layerings_map(i,islope))
501                    if (co2_perice_PCM(i,islope) > 0.) call add_stratum(layerings_map(i,islope),1.05_dp*co2_perice_PCM(i,islope)/rho_co2ice,co2_perice_PCM(i,islope)/rho_co2ice,0._dp,0.05_dp*co2_perice_PCM(i,islope)/rho_co2ice,0._dp,0._dp)
502                end do
503            end if
504        end do
505    end if ! do_layering
506
507    return
508end if
509
510! If the file is present
511! ~~~~~~~~~~~~~~~~~~~~~~
512call print_msg('> Reading "'//startpem_name//'"')
513call open_nc(startpem_name,'read')
514
515! H2O ice
516h2o_ice(:,:) = 0._dp
517call get_var_nc('h2o_ice',h2o_ice(:,:))
518
519! CO2 ice
520co2_ice(:,:) = 0._dp
521call get_var_nc('co2_ice',co2_ice(:,:))
522
523if (do_soil) then
524    ! Check if the number of soil layers is compatible
525    call get_dim_nc('subsurface_layers',nsoil_startpem)
526    if (nsoil_startpem /= nsoil) then
527        call print_msg('nsoil (PEM) = '//int2str(nsoil)//' | nsoil ("'//startpem_name//'") = '//int2str(nsoil_startpem))
528        call stop_clean(__FILE__,__LINE__,'nsoil defined in the PEM is different from the one read in "'//startpem_name//'"!',1)
529    end if
530
531    do islope = 1,nslope
532        ! Soil temperature
533        call get_var_nc('tsoil',tsoil_startpem(:,:,islope))
534        ! Predictor-corrector: due to changes of surface temperature in the PCM, the tsoil_startpem is adapted firstly
535        !                      for PCM year 1 and then for PCM year 2 in order to build the evolution of the profile at depth
536        call compute_tsoil(ngrid,nsoil,.true.,TI(:,:,islope),dt,tsurf_avg_yr1(:,islope),tsoil_startpem(:,:,islope))
537        call compute_tsoil(ngrid,nsoil,.false.,TI(:,:,islope),dt,tsurf_avg_yr1(:,islope),tsoil_startpem(:,:,islope))
538        call compute_tsoil(ngrid,nsoil,.false.,TI(:,:,islope),dt,tsurf_avg_yr2(:,islope),tsoil_startpem(:,:,islope))
539    end do
540    tsoil_avg(:,nsoil_PCM + 1:nsoil,:) = tsoil_startpem(:,nsoil_PCM + 1:nsoil,:)
541    if (any(isnan(tsoil_avg(:,:,:)))) call stop_clean(__FILE__,__LINE__,"NaN detected in 'tsoil_avg'",1)
542
543    ! Thermal inertia
544    call get_var_nc('TI',TI_startpem(:,:,:))
545    TI(:,nsoil_PCM + 1:nsoil,:) = TI_startpem(:,nsoil_PCM + 1:nsoil,:) ! 1st layers can change because of the presence of ice at the surface, so we don't change it here
546    call get_var_nc('inertiedat',inertiedat(:,:))
547
548    ! Ice table
549    call get_var_nc('icetable_depth',icetable_depth(:,:))
550    if (icetable_dynamic) call get_var_nc('ice_porefilling',ice_porefilling(:,:,:))
551    if (icetable_equilibrium) then
552        call get_var_nc('icetable_thickness',icetable_thickness(:,:))
553    else if (icetable_dynamic) then
554        call get_var_nc('ice_porefilling',ice_porefilling(:,:,:))
555    end if
556
557    ! Absorption in regolith
558    if (do_sorption) then
559        call get_var_nc('co2_ads_reg',co2_ads_reg(:,:,:))
560        call get_var_nc('h2o_ads_reg',h2o_ads_reg(:,:,:))
561        call evolve_regolith_adsorption(d_h2oice,d_co2ice,h2o_ice,co2_ice,tsoil_avg,TI,ps_ts,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
562    end if ! do_sorption
563end if ! do_soil
564
565! Layering
566if (do_layering) then
567    call get_dim_nc('nb_str_max',nb_str_max)
568    allocate(layerings_array(ngrid,nslope,nb_str_max,6))
569    layerings_array = 0.
570    call get_var_nc('stratif_top_elevation',layerings_array(:,:,:,1))
571    call get_var_nc('stratif_h_co2ice',layerings_array(:,:,:,2))
572    call get_var_nc('stratif_h_h2oice',layerings_array(:,:,:,3))
573    call get_var_nc('stratif_h_dust',layerings_array(:,:,:,4))
574    call get_var_nc('stratif_h_pore',layerings_array(:,:,:,5))
575    call get_var_nc('stratif_poreice_volfrac',layerings_array(:,:,:,6))
576    call array2map(layerings_array,layerings_map)
577    deallocate(layerings_array)
578end if
579
580! Close
581call close_nc(startpem_name)
582
583END SUBROUTINE read_startpem
584!=======================================================================
585
586END MODULE clim_state_init
Note: See TracBrowser for help on using the repository browser.